linalg 1.6.1
A linear algebra library that provides a user-friendly interface to several BLAS and LAPACK routines.
Loading...
Searching...
No Matches
linalg_immutable.f90
1! linalg_immutable.f90
2
15 use, intrinsic :: iso_fortran_env, only : int32, real64
16 use linalg
17 implicit none
18 private
19 public :: mat_rank1_update
20 public :: mat_mult_diag
21 public :: mat_mult_upper_tri
22 public :: mat_mult_lower_tri
23 public :: mat_det
24 public :: mat_lu
25 public :: mat_qr
26 public :: mat_qr_rank1_update
27 public :: mat_svd
28 public :: mat_cholesky
29 public :: mat_cholesky_rank1_update
30 public :: mat_cholesky_rank1_downdate
31 public :: mat_inverse
32 public :: mat_pinverse
33 public :: mat_solve_upper_tri
34 public :: mat_solve_lower_tri
35 public :: mat_eigen
36 public :: lu_results
37 public :: lu_results_cmplx
38 public :: qr_results
39 public :: qr_results_cmplx
40 public :: svd_results
41 public :: svd_results_cmplx
42 public :: eigen_results
43 public :: identity
44
45! ------------------------------------------------------------------------------
48 interface mat_mult_diag
49 module procedure :: mat_mult_diag_1
50 module procedure :: mat_mult_diag_2
51 module procedure :: mat_mult_diag_3
52 module procedure :: mat_mult_diag_1_cmplx
53 module procedure :: mat_mult_diag_2_cmplx
54 module procedure :: mat_mult_diag_3_cmplx
55 end interface
56
57! ------------------------------------------------------------------------------
61 module procedure :: mat_mult_upper_tri_1
62 module procedure :: mat_mult_upper_tri_2
63 module procedure :: mat_mult_upper_tri_1_cmplx
64 module procedure :: mat_mult_upper_tri_2_cmplx
65 end interface
66
67! ------------------------------------------------------------------------------
71 module procedure :: mat_mult_lower_tri_1
72 module procedure :: mat_mult_lower_tri_2
73 module procedure :: mat_mult_lower_tri_1_cmplx
74 module procedure :: mat_mult_lower_tri_2_cmplx
75 end interface
76
77! ------------------------------------------------------------------------------
81 module procedure :: mat_solve_upper_tri_1
82 module procedure :: mat_solve_upper_tri_2
83 module procedure :: mat_solve_upper_tri_1_cmplx
84 module procedure :: mat_solve_upper_tri_2_cmplx
85 end interface
86
87! ------------------------------------------------------------------------------
91 module procedure :: mat_solve_lower_tri_1
92 module procedure :: mat_solve_lower_tri_2
93 module procedure :: mat_solve_lower_tri_1_cmplx
94 module procedure :: mat_solve_lower_tri_2_cmplx
95 end interface
96
97! ------------------------------------------------------------------------------
100 interface mat_lu
101 module procedure :: mat_lu_dbl
102 module procedure :: mat_lu_cmplx
103 end interface
104
105! ------------------------------------------------------------------------------
108 interface mat_eigen
109 module procedure :: mat_eigen_1
110 module procedure :: mat_eigen_2
111 end interface
112
113! ------------------------------------------------------------------------------
117 real(real64), allocatable, dimension(:,:) :: l
119 real(real64), allocatable, dimension(:,:) :: u
121 real(real64), allocatable, dimension(:,:) :: p
122 end type
123
124! ------------------------------------------------------------------------------
128 complex(real64), allocatable, dimension(:,:) :: l
130 complex(real64), allocatable, dimension(:,:) :: u
132 real(real64), allocatable, dimension(:,:) :: p
133 end type
134
135! ------------------------------------------------------------------------------
139 real(real64), allocatable, dimension(:,:) :: q
141 real(real64), allocatable, dimension(:,:) :: r
144 real(real64), allocatable, dimension(:,:) :: p
145 end type
146
147! ------------------------------------------------------------------------------
151 complex(real64), allocatable, dimension(:,:) :: q
153 complex(real64), allocatable, dimension(:,:) :: r
156 complex(real64), allocatable, dimension(:,:) :: p
157 end type
158
159! ------------------------------------------------------------------------------
164 real(real64), allocatable, dimension(:,:) :: u
166 real(real64), allocatable, dimension(:,:) :: s
168 real(real64), allocatable, dimension(:,:) :: vt
169 end type
170
171! ------------------------------------------------------------------------------
176 complex(real64), allocatable, dimension(:,:) :: u
178 real(real64), allocatable, dimension(:,:) :: s
180 complex(real64), allocatable, dimension(:,:) :: vt
181 end type
182
183! ------------------------------------------------------------------------------
188 complex(real64), allocatable, dimension(:) :: values
191 complex(real64), allocatable, dimension(:,:) :: vectors
192 end type
193
194contains
195! ------------------------------------------------------------------------------
204 function mat_rank1_update(a, x, y) result(b)
205 ! Arguments
206 real(real64), intent(in), dimension(:,:) :: a
207 real(real64), intent(in), dimension(:) :: x, y
208 real(real64), dimension(size(a, 1), size(a, 2)) :: b
209
210 ! Process
211 b = a
212 call rank1_update(1.0d0, x, y, b)
213 end function
214
215! ------------------------------------------------------------------------------
223 function mat_mult_diag_1(a, b) result(c)
224 ! Arguments
225 real(real64), intent(in), dimension(:) :: a
226 real(real64), intent(in), dimension(:,:) :: b
227 real(real64), dimension(size(a), size(b, 2)) :: c
228
229 ! Process
230 if (size(b, 1) > size(a)) then
231 call diag_mtx_mult(.true., .false., 1.0d0, a, b(1:size(a),:), &
232 0.0d0, c)
233 else
234 call diag_mtx_mult(.true., .false., 1.0d0, a, b, 0.0d0, c)
235 end if
236 end function
237
238! ------------------------------------------------------------------------------
246 function mat_mult_diag_2(a, b) result(c)
247 ! Arguments
248 real(real64), intent(in), dimension(:) :: a, b
249 real(real64), dimension(size(a)) :: c
250
251 ! Local Variables
252 real(real64), dimension(size(a), 1) :: bc, cc
253
254 ! Process
255 bc(:,1) = b(1:min(size(a), size(b)))
256 call diag_mtx_mult(.true., .false., 1.0d0, a, bc, 0.0d0, cc)
257 c = cc(:,1)
258 end function
259
260! ------------------------------------------------------------------------------
268 function mat_mult_diag_3(a, b) result(c)
269 ! Arguments
270 real(real64), intent(in), dimension(:,:) :: a
271 real(real64), intent(in), dimension(:) :: b
272 real(real64), dimension(size(a, 1), size(b)) :: c
273
274 ! Process
275 if (size(a, 2) > size(b)) then
276 call diag_mtx_mult(.false., .false., 1.0d0, b, a(:,1:size(b)), &
277 0.0d0, c)
278 else
279 call diag_mtx_mult(.false., .false., 1.0d0, b, a, 0.0d0, c)
280 end if
281 end function
282
283! ------------------------------------------------------------------------------
291 function mat_mult_diag_1_cmplx(a, b) result(c)
292 ! Arguments
293 complex(real64), intent(in), dimension(:) :: a
294 complex(real64), intent(in), dimension(:,:) :: b
295 complex(real64), dimension(size(a), size(b, 2)) :: c
296
297 ! Parameters
298 complex(real64), parameter :: zero = (0.0d0, 0.0d0)
299 complex(real64), parameter :: one = (1.0d0, 0.0d0)
300
301 ! Process
302 if (size(b, 1) > size(a)) then
303 call diag_mtx_mult(.true., la_no_operation, one, a, b(1:size(a),:), &
304 zero, c)
305 else
306 call diag_mtx_mult(.true., la_no_operation, one, a, b, zero, c)
307 end if
308 end function
309
310! ------------------------------------------------------------------------------
318 function mat_mult_diag_2_cmplx(a, b) result(c)
319 ! Arguments
320 complex(real64), intent(in), dimension(:) :: a, b
321 complex(real64), dimension(size(a)) :: c
322
323 ! Parameters
324 complex(real64), parameter :: zero = (0.0d0, 0.0d0)
325 complex(real64), parameter :: one = (1.0d0, 0.0d0)
326
327 ! Local Variables
328 complex(real64), dimension(size(a), 1) :: bc, cc
329
330 ! Process
331 bc(:,1) = b(1:min(size(a), size(b)))
332 call diag_mtx_mult(.true., la_no_operation, one, a, bc, zero, cc)
333 c = cc(:,1)
334 end function
335
336! ------------------------------------------------------------------------------
344 function mat_mult_diag_3_cmplx(a, b) result(c)
345 ! Arguments
346 complex(real64), intent(in), dimension(:,:) :: a
347 complex(real64), intent(in), dimension(:) :: b
348 complex(real64), dimension(size(a, 1), size(b)) :: c
349
350 ! Process
351 if (size(a, 2) > size(b)) then
352 call diag_mtx_mult(.false., la_no_operation, 1.0d0, b, a(:,1:size(b)), &
353 0.0d0, c)
354 else
355 call diag_mtx_mult(.false., la_no_operation, 1.0d0, b, a, 0.0d0, c)
356 end if
357 end function
358
359! ------------------------------------------------------------------------------
366 function mat_mult_upper_tri_1(a, b) result(c)
367 ! Arguments
368 real(real64), intent(in), dimension(:,:) :: a, b
369 real(real64), dimension(size(a, 1), size(b, 2)) :: c
370
371 ! Process
372 c = b
373 call dtrmm('L', 'U', 'N', 'N', size(b, 1), size(b, 2), 1.0d0, &
374 a, size(a, 1), c, size(c, 1))
375 end function
376
377! ------------------------------------------------------------------------------
384 function mat_mult_upper_tri_2(a, b) result(c)
385 ! Arguments
386 real(real64), intent(in), dimension(:,:) :: a
387 real(real64), intent(in), dimension(:) :: b
388 real(real64), dimension(size(a, 1)) :: c
389
390 ! Process
391 c = b
392 call dtrmv('U', 'N', 'N', size(a, 1), a, size(a, 1), c, 1)
393 end function
394
395 ! ------------------------------------------------------------------------------
402 function mat_mult_lower_tri_1(a, b) result(c)
403 ! Arguments
404 real(real64), intent(in), dimension(:,:) :: a, b
405 real(real64), dimension(size(a, 1), size(b, 2)) :: c
406
407 ! Process
408 c = b
409 call dtrmm('L', 'L', 'N', 'N', size(b, 1), size(b, 2), 1.0d0, &
410 a, size(a, 1), c, size(c, 1))
411 end function
412
413 ! ------------------------------------------------------------------------------
420 function mat_mult_lower_tri_2(a, b) result(c)
421 ! Arguments
422 real(real64), intent(in), dimension(:,:) :: a
423 real(real64), intent(in), dimension(:) :: b
424 real(real64), dimension(size(a, 1)) :: c
425
426 ! Process
427 c = b
428 call dtrmv('L', 'N', 'N', size(a, 1), a, size(a, 1), c, 1)
429 end function
430
431! ------------------------------------------------------------------------------
438 function mat_mult_upper_tri_1_cmplx(a, b) result(c)
439 ! Arguments
440 complex(real64), intent(in), dimension(:,:) :: a, b
441 complex(real64), dimension(size(a, 1), size(b, 2)) :: c
442
443 ! Process
444 c = b
445 call ztrmm('L', 'U', 'N', 'N', size(b, 1), size(b, 2), 1.0d0, &
446 a, size(a, 1), c, size(c, 1))
447 end function
448
449! ------------------------------------------------------------------------------
456 function mat_mult_upper_tri_2_cmplx(a, b) result(c)
457 ! Arguments
458 complex(real64), intent(in), dimension(:,:) :: a
459 complex(real64), intent(in), dimension(:) :: b
460 complex(real64), dimension(size(a, 1)) :: c
461
462 ! Process
463 c = b
464 call ztrmv('U', 'N', 'N', size(a, 1), a, size(a, 1), c, 1)
465 end function
466
467 ! ------------------------------------------------------------------------------
474 function mat_mult_lower_tri_1_cmplx(a, b) result(c)
475 ! Arguments
476 complex(real64), intent(in), dimension(:,:) :: a, b
477 complex(real64), dimension(size(a, 1), size(b, 2)) :: c
478
479 ! Process
480 c = b
481 call ztrmm('L', 'L', 'N', 'N', size(b, 1), size(b, 2), 1.0d0, &
482 a, size(a, 1), c, size(c, 1))
483 end function
484
485 ! ------------------------------------------------------------------------------
492 function mat_mult_lower_tri_2_cmplx(a, b) result(c)
493 ! Arguments
494 complex(real64), intent(in), dimension(:,:) :: a
495 complex(real64), intent(in), dimension(:) :: b
496 complex(real64), dimension(size(a, 1)) :: c
497
498 ! Process
499 c = b
500 call ztrmv('L', 'N', 'N', size(a, 1), a, size(a, 1), c, 1)
501 end function
502
503! ------------------------------------------------------------------------------
508 function mat_det(a) result(x)
509 ! Arguments
510 real(real64), intent(in), dimension(:,:) :: a
511 real(real64) :: x
512
513 ! Local Variables
514 real(real64), dimension(size(a, 1), size(a, 2)) :: b
515
516 ! Process
517 b = a
518 x = det(b)
519 end function
520
521! ------------------------------------------------------------------------------
527 function mat_lu_dbl(a) result(x)
528 ! Arguments
529 real(real64), intent(in), dimension(:,:) :: a
530 type(lu_results) :: x
531
532 ! Local Variables
533 integer(int32) :: n
534 integer(int32), allocatable, dimension(:) :: ipvt
535
536 ! Memory Allocation
537 n = size(a, 1)
538 allocate(ipvt(n))
539 allocate(x%l(n,n))
540 allocate(x%u(n,n))
541 allocate(x%p(n,n))
542
543 ! Compute the factorization
544 x%l = a
545 call lu_factor(x%l, ipvt)
546
547 ! Form L, U, and P
548 call form_lu(x%l, ipvt, x%u, x%p)
549 end function
550
551! ------------------------------------------------------------------------------
557 function mat_lu_cmplx(a) result(x)
558 ! Arguments
559 complex(real64), intent(in), dimension(:,:) :: a
560 type(lu_results_cmplx) :: x
561
562 ! Local Variables
563 integer(int32) :: n
564 integer(int32), allocatable, dimension(:) :: ipvt
565
566 ! Memory Allocation
567 n = size(a, 1)
568 allocate(ipvt(n))
569 allocate(x%l(n,n))
570 allocate(x%u(n,n))
571 allocate(x%p(n,n))
572
573 ! Compute the factorization
574 x%l = a
575 call lu_factor(x%l, ipvt)
576
577 ! Form L, U, and P
578 call form_lu(x%l, ipvt, x%u, x%p)
579 end function
580
581! ------------------------------------------------------------------------------
591 function mat_qr(a, pvt) result(x)
592 ! Arguments
593 real(real64), intent(in), dimension(:,:) :: a
594 logical, intent(in), optional :: pvt
595 type(qr_results) :: x
596
597 ! Local Variables
598 logical :: use_pivot
599 integer(int32) :: m, n, mn
600 integer(int32), allocatable, dimension(:) :: jpvt
601 real(real64), allocatable, dimension(:) :: tau
602
603 ! Memory Allocation
604 use_pivot = .false.
605 if (present(pvt)) use_pivot = pvt
606 m = size(a, 1)
607 n = size(a, 2)
608 mn = min(m, n)
609 allocate(tau(mn))
610 allocate(x%q(m,m))
611 allocate(x%r(m,n))
612
613 ! Compute the factorization, and then form Q, R, and P
614 x%r = a
615 if (use_pivot) then
616 allocate(x%p(n,n))
617 allocate(jpvt(n))
618 jpvt = 0 ! Ensure all columns are free columns
619 call qr_factor(x%r, tau, jpvt)
620 call form_qr(x%r, tau, jpvt, x%q, x%p)
621 else
622 call qr_factor(x%r, tau)
623 call form_qr(x%r, tau, x%q)
624 end if
625 end function
626
627! ------------------------------------------------------------------------------
637 function mat_qr_rank1_update(q, r, x, y) result(rst)
638 ! Arguments
639 real(real64), intent(in), dimension(:,:) :: q, r
640 real(real64), intent(in), dimension(:) :: x, y
641 type(qr_results) :: rst
642
643 ! Local Variables
644 integer(int32) :: i, m, n
645 real(real64), allocatable, dimension(:) :: xc, yc
646
647 ! Memory allocation
648 m = size(q, 1)
649 n = size(r, 2)
650 allocate(xc(m))
651 allocate(yc(n))
652 allocate(rst%q(m,m))
653 allocate(rst%r(m,n))
654
655 ! Process
656 do i = 1, m
657 xc(i) = x(i)
658 rst%q(:,i) = q(:,i)
659 end do
660 do i = 1, n
661 yc(i) = y(i)
662 rst%r(:,i) = r(:,i)
663 end do
664 call qr_rank1_update(rst%q, rst%r, xc, yc)
665 end function
666
667! ------------------------------------------------------------------------------
673 function mat_svd(a) result(x)
674 ! Arguments
675 real(real64), intent(in), dimension(:,:) :: a
676 type(svd_results) :: x
677
678 ! Local Variables
679 integer(int32) :: i, m, n, mn
680 real(real64), allocatable, dimension(:) :: s
681 real(real64), allocatable, dimension(:,:) :: ac
682
683 ! Memory Allocation
684 m = size(a, 1)
685 n = size(a, 2)
686 mn = min(m, n)
687 allocate(s(mn))
688 allocate(ac(m,n))
689 allocate(x%u(m,m))
690 allocate(x%s(m,n))
691 allocate(x%vt(n,n))
692
693 ! Process
694 ac = a
695 call svd(ac, s, x%u, x%vt)
696
697 ! Extract the singular values, and populate the results matrix
698 x%s = 0.0d0
699 do i = 1, mn
700 x%s(i,i) = s(i)
701 end do
702 end function
703
704! ------------------------------------------------------------------------------
714 function mat_cholesky(a, upper) result(r)
715 ! Arguments
716 real(real64), intent(in), dimension(:,:) :: a
717 logical, intent(in), optional :: upper
718 real(real64), dimension(size(a, 1), size(a, 2)) :: r
719
720 ! Local Variables
721 logical :: compute_upper
722
723 ! Process
724 compute_upper = .true.
725 if (present(upper)) compute_upper = upper
726 r = a
727 call cholesky_factor(r, compute_upper)
728 end function
729
730! ------------------------------------------------------------------------------
737 function mat_cholesky_rank1_update(a, x) result(r)
738 ! Arguments
739 real(real64), intent(in), dimension(:,:) :: a
740 real(real64), intent(in), dimension(:) :: x
741 real(real64), dimension(size(a, 1), size(a, 2)) :: r
742
743 ! Local Variables
744 real(real64), dimension(size(x)) :: xc
745
746 ! Process
747 r = a
748 xc = x
749 call cholesky_rank1_update(r, xc)
750 end function
751
752! ------------------------------------------------------------------------------
759 function mat_cholesky_rank1_downdate(a, x) result(r)
760 ! Arguments
761 real(real64), intent(in), dimension(:,:) :: a
762 real(real64), intent(in), dimension(:) :: x
763 real(real64), dimension(size(a, 1), size(a, 2)) :: r
764
765 ! Local Variables
766 real(real64), dimension(size(x)) :: xc
767
768 ! Process
769 r = a
770 xc = x
771 call cholesky_rank1_downdate(r, xc)
772 end function
773
774! ------------------------------------------------------------------------------
779 function mat_inverse(a) result(x)
780 ! Arguments
781 real(real64), intent(in), dimension(:,:) :: a
782 real(real64), dimension(size(a, 2), size(a, 1)) :: x
783
784 ! Compute the inverse of A
785 x = a
786 call mtx_inverse(x)
787 end function
788
789! ------------------------------------------------------------------------------
794 function mat_pinverse(a) result(x)
795 ! Arguments
796 real(real64), intent(in), dimension(:,:) :: a
797 real(real64), dimension(size(a, 2), size(a, 1)) :: x
798
799 ! Local Variables
800 real(real64), dimension(size(a, 1), size(a, 2)) :: ac
801
802 ! Compute the inverse of A
803 ac = a
804 call mtx_pinverse(ac, x)
805 end function
806
807! ------------------------------------------------------------------------------
814 function mat_solve_upper_tri_1(a, b) result(x)
815 ! Arguments
816 real(real64), intent(in), dimension(:,:) :: a, b
817 real(real64), dimension(size(b, 1), size(b, 2)) :: x
818
819 ! Process
820 x = b
821 call solve_triangular_system(.true., .true., .false., .true., 1.0d0, &
822 a, x)
823 end function
824
825! ------------------------------------------------------------------------------
832 function mat_solve_upper_tri_2(a, b) result(x)
833 ! Arguments
834 real(real64), intent(in), dimension(:,:) :: a
835 real(real64), intent(in), dimension(:) :: b
836 real(real64), dimension(size(b)) :: x
837
838 ! Process
839 x = b
840 call solve_triangular_system(.true., .false., .true., a, x)
841 end function
842
843! ------------------------------------------------------------------------------
850 function mat_solve_upper_tri_1_cmplx(a, b) result(x)
851 ! Arguments
852 complex(real64), intent(in), dimension(:,:) :: a, b
853 complex(real64), dimension(size(b, 1), size(b, 2)) :: x
854
855 ! Process
856 x = b
857 call solve_triangular_system(.true., .true., .false., .true., &
858 (1.0d0, 0.0d0), a, x)
859 end function
860
861! ------------------------------------------------------------------------------
868 function mat_solve_upper_tri_2_cmplx(a, b) result(x)
869 ! Arguments
870 complex(real64), intent(in), dimension(:,:) :: a
871 complex(real64), intent(in), dimension(:) :: b
872 complex(real64), dimension(size(b)) :: x
873
874 ! Process
875 x = b
876 call solve_triangular_system(.true., .false., .true., a, x)
877 end function
878
879! ------------------------------------------------------------------------------
886 function mat_solve_lower_tri_1(a, b) result(x)
887 ! Arguments
888 real(real64), intent(in), dimension(:,:) :: a, b
889 real(real64), dimension(size(b, 1), size(b, 2)) :: x
890
891 ! Process
892 x = b
893 call solve_triangular_system(.true., .false., .false., .true., 1.0d0, &
894 a, x)
895 end function
896
897! ------------------------------------------------------------------------------
904 function mat_solve_lower_tri_2(a, b) result(x)
905 ! Arguments
906 real(real64), intent(in), dimension(:,:) :: a
907 real(real64), intent(in), dimension(:) :: b
908 real(real64), dimension(size(b)) :: x
909
910 ! Process
911 x = b
912 call solve_triangular_system(.false., .false., .true., a, x)
913 end function
914
915! ------------------------------------------------------------------------------
922 function mat_solve_lower_tri_1_cmplx(a, b) result(x)
923 ! Arguments
924 complex(real64), intent(in), dimension(:,:) :: a, b
925 complex(real64), dimension(size(b, 1), size(b, 2)) :: x
926
927 ! Process
928 x = b
929 call solve_triangular_system(.true., .false., .false., .true., &
930 (1.0d0, 0.0d0), a, x)
931 end function
932
933! ------------------------------------------------------------------------------
940 function mat_solve_lower_tri_2_cmplx(a, b) result(x)
941 ! Arguments
942 complex(real64), intent(in), dimension(:,:) :: a
943 complex(real64), intent(in), dimension(:) :: b
944 complex(real64), dimension(size(b)) :: x
945
946 ! Process
947 x = b
948 call solve_triangular_system(.false., .false., .true., a, x)
949 end function
950
951! ------------------------------------------------------------------------------
958 function mat_eigen_1(a) result(x)
959 ! Arguments
960 real(real64), intent(in), dimension(:,:) :: a
961 type(eigen_results) :: x
962
963 ! Local Variables
964 integer(int32) :: n
965 real(real64), dimension(size(a, 1), size(a, 2)) :: ac
966
967 ! Memory Allocation
968 n = size(a, 1)
969 allocate(x%values(n))
970 allocate(x%vectors(n,n))
971
972 ! Process
973 ac = a
974 call eigen(ac, x%values, x%vectors)
975
976 ! Sort the eigenvalues and eigenvectors.
977 call sort(x%values, x%vectors, .true.)
978 end function
979
980! ------------------------------------------------------------------------------
988 function mat_eigen_2(a, b) result(x)
989 ! Arguments
990 real(real64), intent(in), dimension(:,:) :: a, b
991 type(eigen_results) :: x
992
993 ! Local Variables
994 integer(int32) :: i, j, n
995 real(real64), dimension(size(a, 1), size(a, 2)) :: ac
996 real(real64), dimension(size(b, 1), size(b, 2)) :: bc
997
998 ! Memory Allocation
999 n = size(a, 1)
1000 allocate(x%values(n))
1001 allocate(x%vectors(n,n))
1002
1003 ! Process
1004 do j = 1, n
1005 do i = 1, n
1006 ac(i,j) = a(i,j)
1007 bc(i,j) = b(i,j)
1008 end do
1009 end do
1010 call eigen(ac, bc, x%values, vecs = x%vectors)
1011
1012 ! Sort the eigenvalues and eigenvectors.
1013 call sort(x%values, x%vectors, .true.)
1014 end function
1015
1016! ------------------------------------------------------------------------------
1021 pure function identity(n) result(x)
1022 integer(int32), intent(in) :: n
1023 real(real64), dimension(n, n) :: x
1024 integer(int32) :: i
1025 x = 0.0d0
1026 do i = 1, n
1027 x(i,i) = 1.0d0
1028 end do
1029 end function
1030
1031end module
Computes the Cholesky factorization of a symmetric, positive definite matrix.
Definition: linalg.f90:1433
Computes the rank 1 downdate to a Cholesky factored matrix (upper triangular).
Definition: linalg.f90:1639
Computes the rank 1 update to a Cholesky factored matrix (upper triangular).
Definition: linalg.f90:1532
Computes the determinant of a square matrix.
Definition: linalg.f90:434
Multiplies a diagonal matrix with another matrix or array.
Definition: linalg.f90:329
Computes the eigenvalues, and optionally the eigenvectors, of a matrix.
Definition: linalg.f90:3098
Extracts the L and U matrices from the condensed [L\U] storage format used by the lu_factor.
Definition: linalg.f90:717
Forms the full M-by-M orthogonal matrix Q from the elementary reflectors returned by the base QR fact...
Definition: linalg.f90:1031
Computes the LU factorization of an M-by-N matrix.
Definition: linalg.f90:595
Computes the inverse of a square matrix.
Definition: linalg.f90:2778
Computes the Moore-Penrose pseudo-inverse of a M-by-N matrix using the singular value decomposition o...
Definition: linalg.f90:2884
Computes the QR factorization of an M-by-N matrix.
Definition: linalg.f90:871
Computes the rank 1 update to an M-by-N QR factored matrix A (M >= N) where , and such that .
Definition: linalg.f90:1334
Performs the rank-1 update to matrix A such that: , where is an M-by-N matrix, is a scalar,...
Definition: linalg.f90:194
Solves a triangular system of equations.
Definition: linalg.f90:2061
Sorts an array.
Definition: linalg.f90:3181
Computes the singular value decomposition of a matrix A. The SVD is defined as: , where is an M-by-M...
Definition: linalg.f90:1927
Computes the eigenvalues and eigenvectors (right) of a general N-by-N matrix.
Computes the LU factorization of a square matrix. Notice, partial row pivoting is utilized.
Computes the matrix operation: C = A * B, where A is a diagonal matrix.
Computes the matrix operation C = A * B, where A is a lower triangular matrix.
Computes the matrix operation C = A * B, where A is an upper triangular matrix.
Solves the lower triangular system A X = B, where A is a lower triangular matrix.
Solves the upper triangular system A X = B, where A is an upper triangular matrix.
Provides an immutable interface to many of the core linear algebra routines in this library....
Provides a set of common linear algebra routines.
Definition: linalg.f90:15
Defines a container for the output of an Eigen analysis of a square matrix.
Defines a container for the output of an LU factorization.
Defines a container for the output of an LU factorization.
Defines a container for the output of a QR factorization.
Defines a container for the output of a QR factorization.
Defines a container for the output of a singular value decomposition of a matrix.
Defines a container for the output of a singular value decomposition of a matrix.