23 integer(IK):: i, info, ndim, ntry
= 10
24 character(:, SK),
allocatable :: cform, rform
26 type(display_type) :: disp
29 cform
= getFormat(
mold = [(
0._TKG,
0._TKG)], ed
= SK_
"f", ndigit
= 2_IK, signed
= .true._LK)
30 rform
= getFormat(mold
= [
0._TKG], ed
= SK_
"f", ndigit
= 2_IK, signed
= .true._LK)
33 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
34 call disp%show(
"! Compute inverse of a general real matrix.")
35 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
40 real(TKG),
allocatable :: mat(:,:), inv(:,:), mul(:,:)
41 real(TKG) :: invDetSqrt
44 call disp%show(
"mat = reshape([1, 0, 2, -1, 5, 0, 0, 3, -9], shape = [3,3], order = [2, 1])")
45 mat
= reshape([
1,
0,
2,
-1,
5,
0,
0,
3,
-9], shape
= [
3,
3], order
= [
2,
1])
47 call disp%show( mat ,
format = rform )
48 call disp%show(
"inv = getMatInv(mat)")
51 call disp%show( inv ,
format = rform )
52 call disp%show(
"inv = getMatInv(mat, info) ! gracefully catch inversion errors.")
54 call disp%show(
"if (info /= 0) error stop 'inversion failed.'")
55 if (info
/= 0)
error stop 'inversion failed.'
57 call disp%show( inv ,
format = rform )
58 call disp%show(
"inv = getMatInv(mat, invDetSqrt) ! compute the sqrt of the determinant of the inverse of the general matrix along with the inverse.")
63 call disp%show( inv ,
format = rform )
64 call disp%show(
"mul = matmul(mat, inv)")
65 mul
= matmul(mat, inv)
67 call disp%show( mul,
format = rform )
73 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
74 call disp%show(
"! Compute inverse of a general complex matrix.")
75 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
80 complex(TKG) :: invDetSqrt
81 complex(TKG),
allocatable :: inv(:,:), mul(:,:)
82 complex(TKG),
parameter :: mat(
*,
*)
= reshape(
&
83 [ (
2.0,
1.0), (
2.4,
-1.0), (
2.8,
-1.0), (
3.2,
-1.0), (
3.6,
-1.0), (
4.0,
-1.0), (
4.4,
-1.0), (
4.8,
-1.0), (
5.2,
-1.0)
&
84 , (
2.4,
1.0), (
2.0,
1.0), (
2.4,
-1.0), (
2.8,
-1.0), (
3.2,
-1.0), (
3.6,
-1.0), (
4.0,
-1.0), (
4.4,
-1.0), (
4.8,
-1.0)
&
85 , (
2.8,
1.0), (
2.4,
1.0), (
2.0,
1.0), (
2.4,
-1.0), (
2.8,
-1.0), (
3.2,
-1.0), (
3.6,
-1.0), (
4.0,
-1.0), (
4.4,
-1.0)
&
86 , (
3.2,
1.0), (
2.8,
1.0), (
2.4,
1.0), (
2.0,
1.0), (
2.4,
-1.0), (
2.8,
-1.0), (
3.2,
-1.0), (
3.6,
-1.0), (
4.0,
-1.0)
&
87 , (
3.6,
1.0), (
3.2,
1.0), (
2.8,
1.0), (
2.4,
1.0), (
2.0,
1.0), (
2.4,
-1.0), (
2.8,
-1.0), (
3.2,
-1.0), (
3.6,
-1.0)
&
88 , (
4.0,
1.0), (
3.6,
1.0), (
3.2,
1.0), (
2.8,
1.0), (
2.4,
1.0), (
2.0,
1.0), (
2.4,
-1.0), (
2.8,
-1.0), (
3.2,
-1.0)
&
89 , (
4.4,
1.0), (
4.0,
1.0), (
3.6,
1.0), (
3.2,
1.0), (
2.8,
1.0), (
2.4,
1.0), (
2.0,
1.0), (
2.4,
-1.0), (
2.8,
-1.0)
&
90 , (
4.8,
1.0), (
4.4,
1.0), (
4.0,
1.0), (
3.6,
1.0), (
3.2,
1.0), (
2.8,
1.0), (
2.4,
1.0), (
2.0,
1.0), (
2.4,
-1.0)
&
91 , (
5.2,
1.0), (
4.8,
1.0), (
4.4,
1.0), (
4.0,
1.0), (
3.6,
1.0), (
3.2,
1.0), (
2.8,
1.0), (
2.4,
1.0), (
2.0,
1.0)
&
92 ], shape
= [
9,
9], order
= [
2,
1])
95 call disp%show( mat ,
format = cform )
96 call disp%show(
"inv = getMatInv(mat)")
99 call disp%show( inv ,
format = cform )
100 call disp%show(
"inv = getMatInv(mat, info) ! gracefully catch inversion errors.")
102 call disp%show(
"if (info /= 0) error stop 'inversion failed.'")
103 if (info
/= 0)
error stop 'inversion failed.'
105 call disp%show( inv ,
format = cform )
106 call disp%show(
"inv = getMatInv(mat, invDetSqrt) ! compute the sqrt of the determinant of the inverse of the general matrix along with the inverse.")
111 call disp%show( inv ,
format = cform )
112 call disp%show(
"mul = matmul(mat, inv)")
113 mul
= matmul(mat, inv)
115 call disp%show( mul,
format = cform )
120 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
121 call disp%show(
"! Compute inverse of an upperDiag triangular real matrix.")
122 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
127 real(TKG),
allocatable :: mat(:,:), inv(:,:), mul(:,:)
130 call disp%show(
"ndim = getUnifRand(1_IK, 6_IK)")
132 call disp%show(
"ndim ! matrix rank")
134 call disp%show(
"mat = getUnifRand(1_IK, 9_IK, ndim, ndim)")
136 call disp%show(
"call setMatInit(mat(2 : ndim, 1 : ndim - 1), lowDia, 0._TKG, 0._TKG)")
137 call setMatInit(mat(
2 : ndim,
1 : ndim
- 1), lowDia,
0._TKG,
0._TKG)
139 call disp%show( mat ,
format = rform )
140 call disp%show(
"inv = getMatInv(mat, upperDiag)")
143 call disp%show( inv ,
format = rform )
144 call disp%show(
"mul = matmul(mat, inv)")
145 mul
= matmul(mat, inv)
147 call disp%show( mul ,
format = rform )
153 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
154 call disp%show(
"! Compute inverse of a upperDiag triangular complex matrix.")
155 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
160 complex(TKG),
allocatable :: mat(:,:), inv(:,:), mul(:,:)
163 call disp%show(
"ndim = getUnifRand(1_IK, 6_IK)")
165 call disp%show(
"ndim ! matrix rank")
167 call disp%show(
"mat = getUnifRand((1., 1.), (2., 2.), ndim, ndim)")
169 call disp%show(
"call setMatInit(mat(2 : ndim, 1 : ndim - 1), lowDia, (0._TKG, 0._TKG), (0._TKG, 0._TKG))")
170 call setMatInit(mat(
2 : ndim,
1 : ndim
- 1), lowDia, (
0._TKG,
0._TKG), (
0._TKG,
0._TKG))
172 call disp%show( mat ,
format = cform )
173 call disp%show(
"inv = getMatInv(mat, upperDiag)")
176 call disp%show( inv ,
format = cform )
177 call disp%show(
"mul = matmul(mat, inv)")
178 mul
= matmul(mat, inv)
180 call disp%show( mul ,
format = cform )
186 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
187 call disp%show(
"! Compute inverse of a lowerDiag triangular real matrix.")
188 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
193 real(TKG),
allocatable :: mat(:,:), inv(:,:), mul(:,:)
196 call disp%show(
"ndim = getUnifRand(1_IK, 6_IK)")
198 call disp%show(
"ndim ! matrix rank")
200 call disp%show(
"mat = getUnifRand(1_IK, 9_IK, ndim, ndim)")
202 call disp%show(
"call setMatInit(mat(1 : ndim - 1, 2 : ndim), uppDia, 0._TKG, 0._TKG)")
203 call setMatInit(mat(
1 : ndim
- 1,
2 : ndim), uppDia,
0._TKG,
0._TKG)
205 call disp%show( mat ,
format = rform )
206 call disp%show(
"inv = getMatInv(mat, lowerDiag)")
209 call disp%show( inv ,
format = rform )
210 call disp%show(
"mul = matmul(mat, inv)")
211 mul
= matmul(mat, inv)
213 call disp%show( mul ,
format = rform )
219 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
220 call disp%show(
"! Compute inverse of a lowerDiag triangular complex matrix.")
221 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
226 complex(TKG),
allocatable :: mat(:,:), inv(:,:), mul(:,:)
229 call disp%show(
"ndim = getUnifRand(1_IK, 6_IK)")
231 call disp%show(
"ndim ! matrix rank")
233 call disp%show(
"mat = getUnifRand((1., 1.), (2., 2.), ndim, ndim)")
235 call disp%show(
"call setMatInit(mat(1 : ndim - 1, 2 : ndim), uppDia, (0._TKG, 0._TKG), (0._TKG, 0._TKG))")
236 call setMatInit(mat(
1 : ndim
- 1,
2 : ndim), uppDia, (
0._TKG,
0._TKG), (
0._TKG,
0._TKG))
238 call disp%show( mat ,
format = cform )
239 call disp%show(
"inv = getMatInv(mat, lowerDiag)")
242 call disp%show( inv ,
format = cform )
243 call disp%show(
"mul = matmul(mat, inv)")
244 mul
= matmul(mat, inv)
246 call disp%show( mul ,
format = cform )
252 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
253 call disp%show(
"! Compute inverse of a upperUnit triangular real matrix.")
254 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
259 real(TKG),
allocatable :: mat(:,:), inv(:,:), mul(:,:)
262 call disp%show(
"ndim = getUnifRand(1_IK, 6_IK)")
264 call disp%show(
"ndim ! matrix rank")
266 call disp%show(
"mat = getUnifRand(1_IK, 9_IK, ndim, ndim)")
268 call disp%show(
"call setMatInit(mat, lowDia, 0._TKG, 1._TKG)")
271 call disp%show( mat ,
format = rform )
272 call disp%show(
"inv = getMatInv(mat, upperUnit)")
275 call disp%show( inv ,
format = rform )
276 call disp%show(
"mul = matmul(mat, inv)")
277 mul
= matmul(mat, inv)
279 call disp%show( mul ,
format = rform )
285 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
286 call disp%show(
"! Compute inverse of a upperUnit triangular complex matrix.")
287 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
292 complex(TKG),
allocatable :: mat(:,:), inv(:,:), mul(:,:)
295 call disp%show(
"ndim = getUnifRand(1_IK, 6_IK)")
297 call disp%show(
"ndim ! matrix rank")
299 call disp%show(
"mat = getUnifRand((-1., -1.), (+1., +1.), ndim, ndim)")
300 mat
= getUnifRand((
-1.,
-1.), (
+1.,
+1.), ndim, ndim)
301 call disp%show(
"call setMatInit(mat, lowDia, (0._TKG, 0._TKG), (1._TKG, 0._TKG))")
302 call setMatInit(mat, lowDia, (
0._TKG,
0._TKG), (
1._TKG,
0._TKG))
304 call disp%show( mat ,
format = cform )
305 call disp%show(
"inv = getMatInv(mat, upperUnit)")
308 call disp%show( inv ,
format = cform )
309 call disp%show(
"mul = matmul(mat, inv)")
310 mul
= matmul(mat, inv)
312 call disp%show( mul ,
format = cform )
318 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
319 call disp%show(
"! Compute inverse of a lowerUnit triangular real matrix.")
320 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
325 real(TKG),
allocatable :: mat(:,:), inv(:,:), mul(:,:)
328 call disp%show(
"ndim = getUnifRand(1_IK, 6_IK)")
330 call disp%show(
"ndim ! matrix rank")
332 call disp%show(
"mat = getUnifRand(1_IK, 9_IK, ndim, ndim)")
334 call disp%show(
"call setMatInit(mat, uppDia, 0._TKG, 1._TKG)")
337 call disp%show( mat ,
format = rform )
338 call disp%show(
"inv = getMatInv(mat, lowerUnit)")
341 call disp%show( inv ,
format = rform )
342 call disp%show(
"mul = matmul(mat, inv)")
343 mul
= matmul(mat, inv)
345 call disp%show( mul ,
format = rform )
351 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
352 call disp%show(
"! Compute inverse of a lowerUnit triangular complex matrix.")
353 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
358 complex(TKG),
allocatable :: mat(:,:), inv(:,:), mul(:,:)
361 call disp%show(
"ndim = getUnifRand(1_IK, 6_IK)")
363 call disp%show(
"ndim ! matrix rank")
365 call disp%show(
"mat = getUnifRand((-1., -1.), (+1., +1.), ndim, ndim)")
366 mat
= getUnifRand((
-1.,
-1.), (
+1.,
+1.), ndim, ndim)
367 call disp%show(
"call setMatInit(mat, uppDia, (0._TKG, 0._TKG), (1._TKG, 0._TKG))")
368 call setMatInit(mat, uppDia, (
0._TKG,
0._TKG), (
1._TKG,
0._TKG))
370 call disp%show( mat ,
format = cform )
371 call disp%show(
"inv = getMatInv(mat, lowerUnit)")
374 call disp%show( inv ,
format = cform )
375 call disp%show(
"mul = matmul(mat, inv)")
376 mul
= matmul(mat, inv)
378 call disp%show( mul ,
format = cform )
384 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
385 call disp%show(
"! Compute inverse of a general real matrix by passing its LUP factorization.")
386 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
391 integer(IK),
allocatable :: rperm(:)
392 real(TKG),
allocatable :: mat(:,:), lup(:,:), inv(:,:), mul(:,:)
395 call disp%show(
"mat = reshape([1, 0, 2, -1, 5, 0, 0, 3, -9], shape = [3,3], order = [2, 1])")
396 mat
= reshape([
1,
0,
2,
-1,
5,
0,
0,
3,
-9], shape
= [
3,
3], order
= [
2,
1])
398 call disp%show( mat ,
format = rform )
399 call disp%show(
"call setResized(rperm, size(mat, 1, IK))")
404 call disp%show( lup ,
format = rform )
405 call disp%show(
"call setMatLUP(lup, rperm, info) ! compute the LUP factorization of the matrix.")
407 call disp%show(
"if (info /= 0) error stop 'LUP factorization failed.'")
408 if (info
/= 0)
error stop 'LUP factorization failed.'
410 call disp%show( lup ,
format = rform )
411 call disp%show(
"inv = getMatInv(lup, rperm) ! compute the inverse of the matrix by passing its LUP factorization.")
414 call disp%show( inv ,
format = rform )
415 call disp%show(
"mul = matmul(mat, inv)")
416 mul
= matmul(mat, inv)
418 call disp%show( mul ,
format = rform )
424 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
425 call disp%show(
"! Compute inverse of a general complex matrix by passing its LUP factorization.")
426 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
431 integer(IK),
allocatable :: rperm(:)
432 complex(TKG),
parameter :: mat(
*,
*)
= reshape(
&
433 [ (
2.0,
1.0), (
2.4,
-1.0), (
2.8,
-1.0), (
3.2,
-1.0), (
3.6,
-1.0), (
4.0,
-1.0), (
4.4,
-1.0), (
4.8,
-1.0), (
5.2,
-1.0)
&
434 , (
2.4,
1.0), (
2.0,
1.0), (
2.4,
-1.0), (
2.8,
-1.0), (
3.2,
-1.0), (
3.6,
-1.0), (
4.0,
-1.0), (
4.4,
-1.0), (
4.8,
-1.0)
&
435 , (
2.8,
1.0), (
2.4,
1.0), (
2.0,
1.0), (
2.4,
-1.0), (
2.8,
-1.0), (
3.2,
-1.0), (
3.6,
-1.0), (
4.0,
-1.0), (
4.4,
-1.0)
&
436 , (
3.2,
1.0), (
2.8,
1.0), (
2.4,
1.0), (
2.0,
1.0), (
2.4,
-1.0), (
2.8,
-1.0), (
3.2,
-1.0), (
3.6,
-1.0), (
4.0,
-1.0)
&
437 , (
3.6,
1.0), (
3.2,
1.0), (
2.8,
1.0), (
2.4,
1.0), (
2.0,
1.0), (
2.4,
-1.0), (
2.8,
-1.0), (
3.2,
-1.0), (
3.6,
-1.0)
&
438 , (
4.0,
1.0), (
3.6,
1.0), (
3.2,
1.0), (
2.8,
1.0), (
2.4,
1.0), (
2.0,
1.0), (
2.4,
-1.0), (
2.8,
-1.0), (
3.2,
-1.0)
&
439 , (
4.4,
1.0), (
4.0,
1.0), (
3.6,
1.0), (
3.2,
1.0), (
2.8,
1.0), (
2.4,
1.0), (
2.0,
1.0), (
2.4,
-1.0), (
2.8,
-1.0)
&
440 , (
4.8,
1.0), (
4.4,
1.0), (
4.0,
1.0), (
3.6,
1.0), (
3.2,
1.0), (
2.8,
1.0), (
2.4,
1.0), (
2.0,
1.0), (
2.4,
-1.0)
&
441 , (
5.2,
1.0), (
4.8,
1.0), (
4.4,
1.0), (
4.0,
1.0), (
3.6,
1.0), (
3.2,
1.0), (
2.8,
1.0), (
2.4,
1.0), (
2.0,
1.0)
&
442 ], shape
= [
9,
9], order
= [
2,
1])
443 complex(TKG), dimension(
size(mat,
1),
size(mat,
2)) :: inv, lup, mul
445 call disp%show(
"call setResized(rperm, size(mat, 1, IK))")
450 call disp%show( lup ,
format = cform )
451 call disp%show(
"call setMatLUP(lup, rperm, info) ! compute the LUP factorization of the matrix.")
453 call disp%show(
"if (info /= 0) error stop 'LUP factorization failed.'")
454 if (info
/= 0)
error stop 'LUP factorization failed.'
456 call disp%show( lup ,
format = cform )
457 call disp%show(
"inv = getMatInv(lup, rperm) ! compute the inverse of the matrix by passing its LUP factorization.")
460 call disp%show( inv ,
format = cform )
461 call disp%show(
"mul = matmul(mat, inv)")
462 mul
= matmul(mat, inv)
464 call disp%show( mul ,
format = cform )
469 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
470 call disp%show(
"! Compute the inverse of a positive-definite real matrix by passing its Cholesky factorization.")
471 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
476 real(TKG),
allocatable :: mat(:,:), chol(:,:), inv(:,:), mul(:,:)
479 call disp%show(
"ndim = getUnifRand(1_IK, 6_IK)")
481 call disp%show(
"ndim ! matrix rank")
483 call disp%show(
"mat = getCovRand(mold = 1._TKG, ndim = ndim)")
486 call disp%show( mat ,
format = rform )
487 call disp%show(
"chol = getMatChol(mat, subset = lowDia)")
490 call disp%show( chol ,
format = rform )
491 call disp%show(
"inv = getMatInv(chol, auxil = choLow)")
494 call disp%show( inv ,
format = rform )
495 call disp%show(
"mul = matmul(mat, inv)")
496 mul
= matmul(mat, inv)
498 call disp%show( mul ,
format = rform )
505 real(TKG),
allocatable :: mat(:,:), chol(:,:), inv(:,:), mul(:,:)
508 call disp%show(
"ndim = getUnifRand(1_IK, 6_IK)")
510 call disp%show(
"ndim ! matrix rank")
512 call disp%show(
"mat = getCovRand(mold = 1._TKG, ndim = ndim)")
515 call disp%show( mat ,
format = rform )
516 call disp%show(
"chol = getMatChol(mat, subset = uppDia)")
519 call disp%show( chol ,
format = rform )
520 call disp%show(
"inv = getMatInv(chol, auxil = choUpp)")
523 call disp%show( inv ,
format = rform )
524 call disp%show(
"mul = matmul(mat, inv)")
525 mul
= matmul(mat, inv)
527 call disp%show( mul ,
format = rform )
533 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
534 call disp%show(
"! Compute the inverse of a positive-definite complex matrix by passing its Cholesky factorization.")
535 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
540 integer(IK),
parameter :: ndim
= 3_IK
541 complex(TKG),
allocatable :: chol(:,:), inv(:,:), mul(:,:)
542 complex(TKG),
parameter :: mat(
*,
*)
= reshape( [ (
9.0,
0.0), (
3.0,
3.0), (
3.0,
-3.0)
&
543 , (
3.0,
-3.0),(
18.0,
0.0), (
8.0,
-6.0)
&
544 , (
3.0,
3.0), (
8.0,
6.0),(
43.0,
0.0)
&
545 ], shape
= [ndim, ndim], order
= [
2,
1])
548 call disp%show( mat ,
format = cform )
549 call disp%show(
"chol = getMatChol(mat, subset = lowDia)")
552 call disp%show( chol ,
format = cform )
553 call disp%show(
"inv = getMatInv(chol, auxil = choLow)")
556 call disp%show( inv ,
format = cform )
557 call disp%show(
"mul = matmul(mat, inv)")
558 mul
= matmul(mat, inv)
560 call disp%show( mul ,
format = cform )
566 integer(IK),
parameter :: ndim
= 3_IK
567 complex(TKG),
allocatable :: chol(:,:), inv(:,:), mul(:,:)
568 complex(TKG),
parameter :: mat(
*,
*)
= reshape( [ (
9.0,
0.0), (
3.0,
3.0), (
3.0,
-3.0)
&
569 , (
3.0,
-3.0),(
18.0,
0.0), (
8.0,
-6.0)
&
570 , (
3.0,
3.0), (
8.0,
6.0),(
43.0,
0.0)
&
571 ], shape
= [ndim, ndim], order
= [
2,
1])
574 call disp%show( mat ,
format = cform )
575 call disp%show(
"chol = getMatChol(mat, subset = uppDia)")
578 call disp%show( chol ,
format = cform )
579 call disp%show(
"inv = getMatInv(chol, auxil = choUpp)")
582 call disp%show( inv ,
format = cform )
583 call disp%show(
"mul = matmul(mat, inv)")
584 mul
= matmul(mat, inv)
586 call disp%show( mul ,
format = cform )
Generate and return an array of the specified rank and shape of arbitrary intrinsic type and kind wit...
Allocate or resize (shrink or expand) an input allocatable scalar string or array of rank 1....
Generate and return a random positive-definite (correlation or covariance) matrix using the Gram meth...
Generate and return a scalar or a contiguous array of rank 1 of length s1 of randomly uniformly distr...
Verify the input assertion holds and if it does not, print the (optional) input message on stdout and...
Generate and return an object of type stop_type with the user-specified input attributes.
This is a generic method of the derived type display_type with pass attribute.
Generate and return the upper or the lower Cholesky factorization of the input symmetric positive-def...
Generate and return a matrix of shape (shape(1), shape(2)) with the upper/lower triangle and the diag...
Set the upper/lower triangle and the diagonal elements of the input matrix of arbitrary shape (:,...
Return the LU-Pivoted decomposition of the input square matrix mat(ndim,ndim).
Generate and return the conversion of the input value to an output Fortran string,...
This module contains procedures and generic interfaces for convenient allocation and filling of array...
This module contains procedures and generic interfaces for resizing allocatable arrays of various typ...
This module contains classes and procedures for generating random matrices distributed on the space o...
This module contains classes and procedures for computing various statistical quantities related to t...
This module contains classes and procedures for reporting and handling errors.
This module contains classes and procedures for input/output (IO) or generic display operations on st...
type(display_type) disp
This is a scalar module variable an object of type display_type for general display.
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
integer, parameter RKS
The single-precision real kind in Fortran mode. On most platforms, this is an 32-bit real kind.
This module contains procedures and generic interfaces for computing the Cholesky factorization of po...
This module contains procedures and generic interfaces relevant to generating and initializing matric...
This module contains procedures and generic interfaces relevant to the partially LU Pivoted decomposi...
This module contains the generic procedures for converting values of different types and kinds to For...
Generate and return an object of type display_type.