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

[LEGACY code]
Return the lower-triangle of the Cholesky factorization \(L\) of the symmetric positive-definite real-valued matrix represented by the upper-triangle of the input argument \(\ms{mat} = L.L^T\).
More...

Detailed Description

[LEGACY code]
Return the lower-triangle of the Cholesky factorization \(L\) of the symmetric positive-definite real-valued matrix represented by the upper-triangle of the input argument \(\ms{mat} = L.L^T\).

On input, the upper triangle and diagonal of mat must be specified, which remains intact on output.

Parameters
[in,out]mat: The input array of shape (ndim, ndim) of type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128).
  • On input, the upper triangle and diagonals of mat contain the square symmetric positive-definite (covariance) matrix whose Cholesky factorization is to be computed.
  • On output, the lower triangle of mat contains the lower triangle of the Cholesky factorization of the matrix, while the upper triangle and diagonals remain intact.
[out]dia: The output vector of the same type and kind as mat containing the diagonal elements of the lower triangle of the Cholesky factorization.
If the Cholesky factorization fails, dia(1) = -idim will be set, where idim is the column index causing the singularity.
[in]ndim: The input integer of default kind IK representing the rank of the input square matrix (ndim,ndim).


Possible calling interfaces

call setChoLow(mat, dia, ndim) ! Explicit-shape dummy argument unsafe interface(for benchmarking purposes)
[LEGACY code] Return the lower-triangle of the Cholesky factorization of the symmetric positive-def...
This module contains procedures and generic interfaces for computing the Cholesky factorization of po...
Warning
This generic interface with explicit-shape dummy arguments is not recommended for general usage.
It exists merely for benchmarking purposes against the alternative modern interface setMatChol.
However, if used, and the Cholesky factorization fails, dia(1) will be set to -idim where idim is the index of the column causing the singularity.
The pure procedure(s) documented herein become impure when the ParaMonte library is compiled with preprocessor macro CHECK_ENABLED=1.
By default, these procedures are pure in release build and impure in debug and testing builds.
Remarks
If ndim = 1, mat will not be touched, only sqrt(mat) will be written to the output dia.
See also
getMatChol
setMatChol


Example usage

1program example
2
3 use pm_kind, only: SK
4 use pm_kind, only: IK, LK, RKS, RKD, RKH
5 use pm_io, only: display_type
6 use pm_matrixChol, only: setChoLow
7
8 implicit none
9
10 real(RKH), allocatable :: matrix_RKH(:,:), chodia_RKH(:)
11 real(RKD), allocatable :: matrix_RKD(:,:), chodia_RKD(:)
12 real(RKS), allocatable :: matrix_RKS(:,:), chodia_RKS(:)
13
14 logical(LK) :: failed
15 type(display_type) :: disp
16
17 disp = display_type(file = "main.out.F90")
18
19 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
20 ! Construct the Cholesky lower-triangle: RKH
21 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
22
23 matrix_RKH = reshape( [ 1._RKH, 0._RKH, 0._RKH &
24 , 0._RKH, 4._RKH, 0._RKH &
25 , 2._RKH, 0._RKH, 8._RKH ], shape = [3,3])
26
27 call disp%show("matrix Upper Triangle")
28 call disp%show( matrix_RKH )
29
30 allocate(chodia_RKH(size(matrix_RKH, dim = 1)))
31 call setChoLow(matrix_RKH, chodia_RKH, size(chodia_RKH, 1, IK))
32 if (chodia_RKH(1) <= 0) error stop
33
34 call disp%show("matrix Upper Triangle / Cholesky Lower Triangle")
35 call disp%show( matrix_RKH )
36
37 call disp%show("Cholesky Diagonal")
38 call disp%show( chodia_RKH )
39
40 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
41 ! Construct the Cholesky lower-triangle: RKD
42 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
43
44 matrix_RKD = real(matrix_RKH, kind = RKD)
45
46 call disp%skip()
47 call disp%show("matrix Upper Triangle")
48 call disp%show( matrix_RKD )
49
50 allocate(chodia_RKD(size(matrix_RKD, dim = 1)))
51 call setChoLow(matrix_RKD, chodia_RKD, size(chodia_RKD, 1, IK))
52 if (chodia_RKD(1) <= 0) error stop
53
54 call disp%show("matrix Upper Triangle / Cholesky Lower Triangle")
55 call disp%show( matrix_RKD )
56
57 call disp%show("Cholesky Diagonal")
58 call disp%show( chodia_RKD )
59
60 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
61 ! Construct the Cholesky lower-triangle: RKS
62 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
63
64 matrix_RKS = real(matrix_RKH, kind = RKS)
65
66 call disp%skip()
67 call disp%show("matrix Upper Triangle")
68 call disp%show( matrix_RKS )
69
70 allocate(chodia_RKS(size(matrix_RKS, dim = 1)))
71 call setChoLow(matrix_RKS, chodia_RKS, size(chodia_RKS, 1, IK))
72 if (chodia_RKS(1) <= 0) error stop
73
74 call disp%show("matrix Upper Triangle / Cholesky Lower Triangle")
75 call disp%show( matrix_RKS )
76
77 call disp%show("Cholesky Diagonal")
78 call disp%show( chodia_RKS )
79
80end program example
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 RKD
The double precision real kind in Fortran mode. On most platforms, this is an 64-bit real kind.
Definition: pm_kind.F90:568
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 RKH
The scalar integer constant of intrinsic default kind, representing the highest-precision real kind t...
Definition: pm_kind.F90:858
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
1matrix Upper Triangle
2+1.00000000000000000000000000000000000, +0.00000000000000000000000000000000000, +2.00000000000000000000000000000000000
3+0.00000000000000000000000000000000000, +4.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
4+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000, +8.00000000000000000000000000000000000
5matrix Upper Triangle / Cholesky Lower Triangle
6+1.00000000000000000000000000000000000, +0.00000000000000000000000000000000000, +2.00000000000000000000000000000000000
7+0.00000000000000000000000000000000000, +4.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
8+2.00000000000000000000000000000000000, +0.00000000000000000000000000000000000, +8.00000000000000000000000000000000000
9Cholesky Diagonal
10+1.00000000000000000000000000000000000, +2.00000000000000000000000000000000000, +2.00000000000000000000000000000000000
11
12matrix Upper Triangle
13+1.0000000000000000, +0.0000000000000000, +2.0000000000000000
14+0.0000000000000000, +4.0000000000000000, +0.0000000000000000
15+2.0000000000000000, +0.0000000000000000, +8.0000000000000000
16matrix Upper Triangle / Cholesky Lower Triangle
17+1.0000000000000000, +0.0000000000000000, +2.0000000000000000
18+0.0000000000000000, +4.0000000000000000, +0.0000000000000000
19+2.0000000000000000, +0.0000000000000000, +8.0000000000000000
20Cholesky Diagonal
21+1.0000000000000000, +2.0000000000000000, +2.0000000000000000
22
23matrix Upper Triangle
24+1.00000000, +0.00000000, +2.00000000
25+0.00000000, +4.00000000, +0.00000000
26+2.00000000, +0.00000000, +8.00000000
27matrix Upper Triangle / Cholesky Lower Triangle
28+1.00000000, +0.00000000, +2.00000000
29+0.00000000, +4.00000000, +0.00000000
30+2.00000000, +0.00000000, +8.00000000
31Cholesky Diagonal
32+1.00000000, +2.00000000, +2.00000000
33
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 Austin

Definition at line 209 of file pm_matrixChol.F90.


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