ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_matrixChol::getMatChol Interface Reference

Generate and return the upper or the lower Cholesky factorization of the input symmetric positive-definite matrix represented by the upper-triangle of the input argument \(\ms{mat} = L.L^T\).
More...

Detailed Description

Generate and return the upper or the lower Cholesky factorization of the input symmetric positive-definite matrix represented by the upper-triangle of the input argument \(\ms{mat} = L.L^T\).

This generic interface is merely a wrapper for the lower-level, potentially faster, more comprehensive generic interface setMatChol.

Parameters
[in]mat: The input square positive-definite matrix of shape (ndim, ndim) of
  1. type complex of kind any supported by the processor (e.g., CK, CK32, CK64, or CK128),
  2. type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
containing the specified subset of the positive definite matrix whose Cholesky factorization is to be computed.
[in]subset: The input scalar that can be either,
  1. the constant uppDia implying that the upper-diagonal triangular block of mat should be used for building the Cholesky factor while the lower triangular part of mat is not referenced within the algorithm.
  2. the constant lowDia implying that the lower-diagonal triangular block of mat should be used for building the Cholesky factor while the upper triangular part of mat is not referenced within the algorithm.
This argument is merely a convenience to differentiate the different procedure functionalities within this generic interface.
[in]operation: The input scalar that can be either,
  1. the constant nothing implying that the same subset of the output chol as that of mat specified by the input subset must be populated with the Cholesky factorization.
  2. the constant transHerm implying that the Hermitian transpose of the computed Cholesky factorization must be written to the output chol.
    For example, specifying subset = uppDia, operation = transHerm will lead to using the upper-diagonal subset of the input mat to compute the lower-diagonal triangle of chol.
(optional, default = nothing)
Returns
chol : The output matrix of the same type, kind, and shape as the input mat containing the requested subset of the Cholesky factorization of mat.
  1. If subset = uppDia, operation = nothing, then the upper-diagonal Cholesky factorization of mat is computed.
  2. If subset = lowDia, operation = nothing, then the lower-diagonal Cholesky factorization of mat is computed.
  3. If subset = uppDia, operation = transHerm, then the lower-diagonal Cholesky factorization of mat is computed.
  4. If subset = lowDia, operation = transHerm, then the upper-diagonal Cholesky factorization of mat is computed.
By definition, the opposite-triangular elements of chol that do not contain the Cholesky factorization will be set to zero.
The program will call error stop is the Cholesky factorization fails.


Possible calling interfaces

chol(1:ndim, 1:ndim) = getMatChol(mat(1:ndim, 1:ndim), subset, operation = operation)
Generate and return the upper or the lower Cholesky factorization of the input symmetric positive-def...
This module contains procedures and generic interfaces for computing the Cholesky factorization of po...
Note
This generic interface is a convenience wrapper around the more efficient generic subroutine interface setMatChol.
Remarks
The procedures under discussion are pure. The procedures of this generic interface are impure when the output argument failed is present.
See also
pm_matrixDet
pm_matrixInv
pm_matrixLUP
transHerm
iteration
recursion
lowDia
uppDia


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK
4 use pm_kind, only: TKG => RKS ! all processor type kinds are supported.
5 use pm_matrixChol, only: getMatChol, lowDia, uppDia
6 use pm_io, only: display_type
7 use pm_io, only: getFormat
8
9 implicit none
10
11 type(display_type) :: disp
12
13 character(:, SK), allocatable :: cform, gform
14 real(TKG) , parameter :: DUM = -huge(0._TKG)
15 complex(TKG) , parameter :: CMPLX_DUMM = cmplx(-huge(0._TKG), -huge(0._TKG), TKG)
16 cform = getFormat([CMPLX_DUMM], ed = SK_'f', signed = .true.)
17 gform = getFormat([DUM], ed = SK_'f', signed = .true.)
18
19 disp = display_type(file = "main.out.F90")
20
21 call disp%skip
22 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
23 call disp%show("!Cholesky factorization of real positive-definite square matrix update.")
24 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
25 call disp%skip
26
27 block
28 real(TKG), allocatable :: cholref(:,:), cholmat(:,:)
29
30 cholmat = reshape([real(TKG) :: 1.0, DUM, DUM, DUM, DUM, DUM, DUM, DUM, DUM &
31 , 1.0, 2.0, DUM, DUM, DUM, DUM, DUM, DUM, DUM &
32 , 1.0, 2.0, 3.0, DUM, DUM, DUM, DUM, DUM, DUM &
33 , 1.0, 2.0, 3.0, 4.0, DUM, DUM, DUM, DUM, DUM &
34 , 1.0, 2.0, 3.0, 4.0, 5.0, DUM, DUM, DUM, DUM &
35 , 1.0, 2.0, 3.0, 4.0, 5.0, 6.0, DUM, DUM, DUM &
36 , 1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, DUM, DUM &
37 , 1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0, DUM &
38 , 1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0, 9.0 &
39 ], shape = [9, 9], order = [2, 1])
40 cholref = reshape([real(TKG) :: 1.0, DUM, DUM, DUM, DUM, DUM, DUM, DUM, DUM &
41 , 1.0, 1.0, DUM, DUM, DUM, DUM, DUM, DUM, DUM &
42 , 1.0, 1.0, 1.0, DUM, DUM, DUM, DUM, DUM, DUM &
43 , 1.0, 1.0, 1.0, 1.0, DUM, DUM, DUM, DUM, DUM &
44 , 1.0, 1.0, 1.0, 1.0, 1.0, DUM, DUM, DUM, DUM &
45 , 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, DUM, DUM, DUM &
46 , 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, DUM, DUM &
47 , 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, DUM &
48 , 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0 &
49 ], shape = [9, 9], order = [2, 1])
50 call disp%skip()
51 call disp%show("cholmat")
52 call disp%show( cholmat , format = gform )
53 call disp%show("cholmat = getMatChol(cholmat, lowDia)")
54 cholmat = getMatChol(cholmat, lowDia)
55 call disp%show("cholmat")
56 call disp%show( cholmat , format = gform )
57 call disp%show("cholref")
58 call disp%show( cholref , format = gform )
59 call disp%skip()
60
61
62 cholmat = reshape([real(TKG) :: DUM, DUM, DUM, DUM, DUM, DUM, DUM, DUM, DUM, DUM &
63 , DUM, DUM, DUM, DUM, DUM, DUM, DUM, DUM, DUM, DUM &
64 , DUM, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0 &
65 , DUM, DUM, 2.0, 2.0, 2.0, 2.0, 2.0, 2.0, 2.0, 2.0 &
66 , DUM, DUM, DUM, 3.0, 3.0, 3.0, 3.0, 3.0, 3.0, 3.0 &
67 , DUM, DUM, DUM, DUM, 4.0, 4.0, 4.0, 4.0, 4.0, 4.0 &
68 , DUM, DUM, DUM, DUM, DUM, 5.0, 5.0, 5.0, 5.0, 5.0 &
69 , DUM, DUM, DUM, DUM, DUM, DUM, 6.0, 6.0, 6.0, 6.0 &
70 , DUM, DUM, DUM, DUM, DUM, DUM, DUM, 7.0, 7.0, 7.0 &
71 , DUM, DUM, DUM, DUM, DUM, DUM, DUM, DUM, 8.0, 8.0 &
72 , DUM, DUM, DUM, DUM, DUM, DUM, DUM, DUM, DUM, 9.0 &
73 ], shape = [11, 10], order = [2, 1])
74 cholref = reshape([real(TKG) :: DUM, DUM, DUM, DUM, DUM, DUM, DUM, DUM, DUM, DUM &
75 , DUM, DUM, DUM, DUM, DUM, DUM, DUM, DUM, DUM, DUM &
76 , DUM, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0 &
77 , DUM, DUM, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0 &
78 , DUM, DUM, DUM, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0 &
79 , DUM, DUM, DUM, DUM, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0 &
80 , DUM, DUM, DUM, DUM, DUM, 1.0, 1.0, 1.0, 1.0, 1.0 &
81 , DUM, DUM, DUM, DUM, DUM, DUM, 1.0, 1.0, 1.0, 1.0 &
82 , DUM, DUM, DUM, DUM, DUM, DUM, DUM, 1.0, 1.0, 1.0 &
83 , DUM, DUM, DUM, DUM, DUM, DUM, DUM, DUM, 1.0, 1.0 &
84 , DUM, DUM, DUM, DUM, DUM, DUM, DUM, DUM, DUM, 1.0 &
85 ], shape = [11, 10], order = [2, 1])
86 call disp%skip()
87 call disp%show("cholmat")
88 call disp%show( cholmat , format = gform )
89 call disp%show("cholmat = getMatChol(cholmat(3:,2:), uppDia)")
90 cholmat = getMatChol(cholmat(3:,2:), uppDia)
91 call disp%show("cholmat")
92 call disp%show( cholmat , format = gform )
93 call disp%show("cholref")
94 call disp%show( cholref , format = gform )
95 call disp%skip()
96
97 end block
98
99 call disp%skip
100 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
101 call disp%show("!Cholesky factorization of complex positive-definite square matrix update.")
102 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
103 call disp%skip
104
105 block
106
107 complex(TKG), allocatable :: cholref(:,:), cholmat(:,:)
108
109 cholmat = reshape( [complex(TKG) :: CMPLX_DUMM, CMPLX_DUMM, CMPLX_DUMM, CMPLX_DUMM, CMPLX_DUMM &
110 , CMPLX_DUMM, CMPLX_DUMM, (25.0, DUM), CMPLX_DUMM, CMPLX_DUMM &
111 , CMPLX_DUMM, CMPLX_DUMM, (-5.0, 5.0), (51.0, DUM), CMPLX_DUMM &
112 , CMPLX_DUMM, CMPLX_DUMM, (10.0, -5.0), ( 4.0, 6.0), (71.0, DUM) &
113 ], shape = [4, 5], order = [2, 1])
114 cholref = reshape( [complex(TKG) :: CMPLX_DUMM, CMPLX_DUMM, CMPLX_DUMM, CMPLX_DUMM, CMPLX_DUMM &
115 , CMPLX_DUMM, CMPLX_DUMM, ( 5.0, 0.0), CMPLX_DUMM, CMPLX_DUMM &
116 , CMPLX_DUMM, CMPLX_DUMM, (-1.0, 1.0), (7.0, 0.0), CMPLX_DUMM &
117 , CMPLX_DUMM, CMPLX_DUMM, ( 2.0, -1.0), (1.0, 1.0), (8.0, 0.0) &
118 ], shape = [4, 5], order = [2, 1])
119 call disp%skip()
120 call disp%show("cholmat")
121 call disp%show( cholmat , format = cform )
122 call disp%show("cholmat = getMatChol(cholmat(2:,3:), lowDia)")
123 cholmat = getMatChol(cholmat(2:,3:), lowDia)
124 call disp%show("cholmat")
125 call disp%show( cholmat , format = cform )
126 call disp%show("cholref")
127 call disp%show( cholref , format = cform )
128 call disp%skip()
129
130 cholmat = reshape( [complex(TKG) :: CMPLX_DUMM, CMPLX_DUMM, CMPLX_DUMM, CMPLX_DUMM &
131 , (9.0, DUM), ( 3.0,3.0), ( 3.0,-3.0), CMPLX_DUMM &
132 , CMPLX_DUMM, (18.0,DUM), ( 8.0,-6.0), CMPLX_DUMM &
133 , CMPLX_DUMM, CMPLX_DUMM, (43.0, DUM), CMPLX_DUMM &
134 ], shape = [4, 4], order = [2, 1])
135 cholref = reshape( [complex(TKG) :: CMPLX_DUMM, CMPLX_DUMM, CMPLX_DUMM, CMPLX_DUMM &
136 , (3.0, 0.0), (1.0, 1.0), (1.0, -1.0), CMPLX_DUMM &
137 , CMPLX_DUMM, (4.0, 0.0), (2.0, -1.0), CMPLX_DUMM &
138 , CMPLX_DUMM, CMPLX_DUMM, (6.0, 0.0), CMPLX_DUMM &
139 ], shape = [4, 4], order = [2, 1])
140 call disp%skip()
141 call disp%show("cholmat")
142 call disp%show( cholmat , format = cform )
143 call disp%show("cholmat = getMatChol(cholmat(2:,:3), uppDia)")
144 cholmat = getMatChol(cholmat(2:,:3), uppDia)
145 call disp%show("cholmat")
146 call disp%show( cholmat , format = cform )
147 call disp%show("cholref")
148 call disp%show( cholref , format = cform )
149 call disp%skip()
150
151 end block
152
153end program example
Generate and return a generic or type/kind-specific IO format with the requested specifications that ...
Definition: pm_io.F90:18485
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11726
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11508
This module contains classes and procedures for input/output (IO) or generic display operations on st...
Definition: pm_io.F90:252
type(display_type) disp
This is a scalar module variable an object of type display_type for general display.
Definition: pm_io.F90:11393
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
Definition: pm_kind.F90:268
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
Definition: pm_kind.F90:541
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
Definition: pm_kind.F90:540
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
Definition: pm_kind.F90:539
integer, parameter RKS
The single-precision real kind in Fortran mode. On most platforms, this is an 32-bit real kind.
Definition: pm_kind.F90:567
Generate and return an object of type display_type.
Definition: pm_io.F90:10282

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1
2!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
3!Cholesky factorization of real positive-definite square matrix update.
4!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6
7cholmat
8 +1.000000, ***************, ***************, ***************, ***************, ***************, ***************, ***************, ***************
9 +1.000000, +2.000000, ***************, ***************, ***************, ***************, ***************, ***************, ***************
10 +1.000000, +2.000000, +3.000000, ***************, ***************, ***************, ***************, ***************, ***************
11 +1.000000, +2.000000, +3.000000, +4.000000, ***************, ***************, ***************, ***************, ***************
12 +1.000000, +2.000000, +3.000000, +4.000000, +5.000000, ***************, ***************, ***************, ***************
13 +1.000000, +2.000000, +3.000000, +4.000000, +5.000000, +6.000000, ***************, ***************, ***************
14 +1.000000, +2.000000, +3.000000, +4.000000, +5.000000, +6.000000, +7.000000, ***************, ***************
15 +1.000000, +2.000000, +3.000000, +4.000000, +5.000000, +6.000000, +7.000000, +8.000000, ***************
16 +1.000000, +2.000000, +3.000000, +4.000000, +5.000000, +6.000000, +7.000000, +8.000000, +9.000000
17cholmat = getMatChol(cholmat, lowDia)
18cholmat
19 +1.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
20 +1.000000, +1.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
21 +1.000000, +1.000000, +1.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
22 +1.000000, +1.000000, +1.000000, +1.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
23 +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +0.000000, +0.000000, +0.000000, +0.000000
24 +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +0.000000, +0.000000, +0.000000
25 +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +0.000000, +0.000000
26 +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +0.000000
27 +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000
28cholref
29 +1.000000, ***************, ***************, ***************, ***************, ***************, ***************, ***************, ***************
30 +1.000000, +1.000000, ***************, ***************, ***************, ***************, ***************, ***************, ***************
31 +1.000000, +1.000000, +1.000000, ***************, ***************, ***************, ***************, ***************, ***************
32 +1.000000, +1.000000, +1.000000, +1.000000, ***************, ***************, ***************, ***************, ***************
33 +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, ***************, ***************, ***************, ***************
34 +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, ***************, ***************, ***************
35 +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, ***************, ***************
36 +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, ***************
37 +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000
38
39
40cholmat
41***************, ***************, ***************, ***************, ***************, ***************, ***************, ***************, ***************, ***************
42***************, ***************, ***************, ***************, ***************, ***************, ***************, ***************, ***************, ***************
43***************, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000
44***************, ***************, +2.000000, +2.000000, +2.000000, +2.000000, +2.000000, +2.000000, +2.000000, +2.000000
45***************, ***************, ***************, +3.000000, +3.000000, +3.000000, +3.000000, +3.000000, +3.000000, +3.000000
46***************, ***************, ***************, ***************, +4.000000, +4.000000, +4.000000, +4.000000, +4.000000, +4.000000
47***************, ***************, ***************, ***************, ***************, +5.000000, +5.000000, +5.000000, +5.000000, +5.000000
48***************, ***************, ***************, ***************, ***************, ***************, +6.000000, +6.000000, +6.000000, +6.000000
49***************, ***************, ***************, ***************, ***************, ***************, ***************, +7.000000, +7.000000, +7.000000
50***************, ***************, ***************, ***************, ***************, ***************, ***************, ***************, +8.000000, +8.000000
51***************, ***************, ***************, ***************, ***************, ***************, ***************, ***************, ***************, +9.000000
52cholmat = getMatChol(cholmat(3:,2:), uppDia)
53cholmat
54 +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000
55 +0.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000
56 +0.000000, +0.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000
57 +0.000000, +0.000000, +0.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000
58 +0.000000, +0.000000, +0.000000, +0.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000
59 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +1.000000, +1.000000, +1.000000, +1.000000
60 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +1.000000, +1.000000, +1.000000
61 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +1.000000, +1.000000
62 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +1.000000
63cholref
64***************, ***************, ***************, ***************, ***************, ***************, ***************, ***************, ***************, ***************
65***************, ***************, ***************, ***************, ***************, ***************, ***************, ***************, ***************, ***************
66***************, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000
67***************, ***************, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000
68***************, ***************, ***************, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000
69***************, ***************, ***************, ***************, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000
70***************, ***************, ***************, ***************, ***************, +1.000000, +1.000000, +1.000000, +1.000000, +1.000000
71***************, ***************, ***************, ***************, ***************, ***************, +1.000000, +1.000000, +1.000000, +1.000000
72***************, ***************, ***************, ***************, ***************, ***************, ***************, +1.000000, +1.000000, +1.000000
73***************, ***************, ***************, ***************, ***************, ***************, ***************, ***************, +1.000000, +1.000000
74***************, ***************, ***************, ***************, ***************, ***************, ***************, ***************, ***************, +1.000000
75
76
77!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
78!Cholesky factorization of complex positive-definite square matrix update.
79!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
80
81
82cholmat
83(***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************)
84(***************, ***************), (***************, ***************), ( +25.000000, ***************), (***************, ***************), (***************, ***************)
85(***************, ***************), (***************, ***************), ( -5.000000, +5.000000), ( +51.000000, ***************), (***************, ***************)
86(***************, ***************), (***************, ***************), ( +10.000000, -5.000000), ( +4.000000, +6.000000), ( +71.000000, ***************)
87cholmat = getMatChol(cholmat(2:,3:), lowDia)
88cholmat
89( +5.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000)
90( -1.000000, +1.000000), ( +7.000000, +0.000000), ( +0.000000, +0.000000)
91( +2.000000, -1.000000), ( +1.000000, +1.000000), ( +8.000000, +0.000000)
92cholref
93(***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************)
94(***************, ***************), (***************, ***************), ( +5.000000, +0.000000), (***************, ***************), (***************, ***************)
95(***************, ***************), (***************, ***************), ( -1.000000, +1.000000), ( +7.000000, +0.000000), (***************, ***************)
96(***************, ***************), (***************, ***************), ( +2.000000, -1.000000), ( +1.000000, +1.000000), ( +8.000000, +0.000000)
97
98
99cholmat
100(***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************)
101( +9.000000, ***************), ( +3.000000, +3.000000), ( +3.000000, -3.000000), (***************, ***************)
102(***************, ***************), ( +18.000000, ***************), ( +8.000000, -6.000000), (***************, ***************)
103(***************, ***************), (***************, ***************), ( +43.000000, ***************), (***************, ***************)
104cholmat = getMatChol(cholmat(2:,:3), uppDia)
105cholmat
106( +3.000000, +0.000000), ( +1.000000, +1.000000), ( +1.000000, -1.000000)
107( +0.000000, +0.000000), ( +4.000000, +0.000000), ( +2.000000, -1.000000)
108( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +6.000000, +0.000000)
109cholref
110(***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************)
111( +3.000000, +0.000000), ( +1.000000, +1.000000), ( +1.000000, -1.000000), (***************, ***************)
112(***************, ***************), ( +4.000000, +0.000000), ( +2.000000, -1.000000), (***************, ***************)
113(***************, ***************), (***************, ***************), ( +6.000000, +0.000000), (***************, ***************)
114
115
Test:
test_pm_matrixChol


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, Apr 21, 2017, 1:54 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin

Definition at line 368 of file pm_matrixChol.F90.


The documentation for this interface was generated from the following file: