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

Return the LU-Pivoted decomposition of the input square matrix mat(ndim,ndim).
More...

Detailed Description

Return the LU-Pivoted decomposition of the input square matrix mat(ndim,ndim).

The LUP factorization takes the form,

\begin{equation} \text{mat} = P * L * U \end{equation}

where,

  1. \(P\) is the permutation matrix,
  2. \(L\) is lower triangular with unit diagonal elements (that are missing in the output mat),
  3. \(U\) is upper triangular of the decomposition.

The abbreviation LUP stands for the LU factorization with partial Pivoting.

Parameters
[in,out]mat: The input contiguous square 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 matrix whose LUP factorization must be returned.
On output, mat is completely overwritten by its LU Pivoted decomposition.
[out]rperm: The output contiguous vector of size ndim of type integer of default kind IK containing the pivot indices (row permutations).
For 1 <= i <= ndim, the ith row of the matrix is interchanged with row rperm(i).
[out]info: The output scalar integer of default kind IK that is non-zero if a singular matrix is detected, indicating the LUP decomposition failure.


Possible calling interfaces

call setMatLUP(mat(1 : ndim, 1 : ndim), rperm(1 : ndim), info)
Return the LU-Pivoted decomposition of the input square matrix mat(ndim,ndim).
This module contains procedures and generic interfaces relevant to the partially LU Pivoted decomposi...
Warning
The condition size(mat, 1) == size(mat, 2)) must hold for the corresponding input arguments (unless dispatch to LAPACK is enabled).
The condition size(rperm) == shape(mat, 1, IK) must hold for the corresponding input arguments.
These conditions are verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
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
This routine can be used in combination with the generic interfaces of pm_matrixMulTri to iteratively solve systems of linear equations or to invert matrices.
See also
getMatInv
setMatInv
getMatChol
setMatChol
setMatLUP
BLAS/LAPACK equivalent:
The procedures under discussion combine, modernize, and extend the interface and functionalities of Version 3.11 of BLAS/LAPACK routine(s): SGETRF, DGETRF, CGETRF, and ZGETRF.


Example usage

1program example
2
3 use pm_kind, only: IK, LK, SK, RKS, RKD, RKH
6 use pm_matrixLUP, only: setMatLUP
7 use pm_io, only: display_type
8
9 implicit none
10
11 integer(IK), allocatable :: rperm(:), rperm_ref(:)
12 integer(IK) :: info
13
14 type(display_type) :: disp
15 disp = display_type(file = "main.out.F90")
16
17 call disp%skip
18 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
19 call disp%show("! Compute the LU-Pivoted decomposition of a real square matrix.")
20 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
21 call disp%skip
22
23 block
24 use pm_kind, only: TKG => RKS
25 real(TKG), allocatable :: mat_lup(:,:), lup_ref(:,:)
26 mat_lup = reshape( [ 1.0_TKG, +1.2_TKG, 1.4_TKG, 1.6_TKG, 1.8_TKG, 2.0_TKG, 2.2_TKG, 2.4_TKG, 2.6_TKG &
27 , 1.2_TKG, +1.0_TKG, 1.2_TKG, 1.4_TKG, 1.6_TKG, 1.8_TKG, 2.0_TKG, 2.2_TKG, 2.4_TKG &
28 , 1.4_TKG, +1.2_TKG, 1.0_TKG, 1.2_TKG, 1.4_TKG, 1.6_TKG, 1.8_TKG, 2.0_TKG, 2.2_TKG &
29 , 1.6_TKG, +1.4_TKG, 1.2_TKG, 1.0_TKG, 1.2_TKG, 1.4_TKG, 1.6_TKG, 1.8_TKG, 2.0_TKG &
30 , 1.8_TKG, +1.6_TKG, 1.4_TKG, 1.2_TKG, 1.0_TKG, 1.2_TKG, 1.4_TKG, 1.6_TKG, 1.8_TKG &
31 , 2.0_TKG, +1.8_TKG, 1.6_TKG, 1.4_TKG, 1.2_TKG, 1.0_TKG, 1.2_TKG, 1.4_TKG, 1.6_TKG &
32 , 2.2_TKG, +2.0_TKG, 1.8_TKG, 1.6_TKG, 1.4_TKG, 1.2_TKG, 1.0_TKG, 1.2_TKG, 1.4_TKG &
33 , 2.4_TKG, +2.2_TKG, 2.0_TKG, 1.8_TKG, 1.6_TKG, 1.4_TKG, 1.2_TKG, 1.0_TKG, 1.2_TKG &
34 , 2.6_TKG, +2.4_TKG, 2.2_TKG, 2.0_TKG, 1.8_TKG, 1.6_TKG, 1.4_TKG, 1.2_TKG, 1.0_TKG &
35 ], shape = [9, 9], order = [2, 1])
36 lup_ref = reshape( [ 2.6_TKG, 2.4_TKG, 2.2_TKG, 2.0_TKG, 1.8_TKG, 1.6_TKG, 1.4_TKG, 1.2_TKG, 1.0_TKG &
37 , 0.4_TKG, 0.3_TKG, 0.6_TKG, 0.8_TKG, 1.1_TKG, 1.4_TKG, 1.7_TKG, 1.9_TKG, 2.2_TKG &
38 , 0.5_TKG, -0.4_TKG, 0.4_TKG, 0.8_TKG, 1.2_TKG, 1.6_TKG, 2.0_TKG, 2.4_TKG, 2.8_TKG &
39 , 0.5_TKG, -0.3_TKG, 0.0_TKG, 0.4_TKG, 0.8_TKG, 1.2_TKG, 1.6_TKG, 2.0_TKG, 2.4_TKG &
40 , 0.6_TKG, -0.3_TKG, 0.0_TKG, 0.0_TKG, 0.4_TKG, 0.8_TKG, 1.2_TKG, 1.6_TKG, 2.0_TKG &
41 , 0.7_TKG, -0.2_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.4_TKG, 0.8_TKG, 1.2_TKG, 1.6_TKG &
42 , 0.8_TKG, -0.2_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.4_TKG, 0.8_TKG, 1.2_TKG &
43 , 0.8_TKG, -0.1_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.4_TKG, 0.8_TKG &
44 , 0.9_TKG, -0.1_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.4_TKG &
45 ], shape = [9, 9], order = [2, 1])
46 rperm_ref = [9, 9, 9, 9, 9, 9, 9, 9, 9]
47 call disp%skip
48 call disp%show("mat_lup")
49 call disp%show( mat_lup )
50 call disp%show("call setResized(rperm, size(mat_lup, 1, IK))")
51 call setResized(rperm, size(mat_lup, 1, IK))
52 call disp%show("call setMatLUP(mat_lup, rperm, info); if (info /= 0) error stop")
53 call setMatLUP(mat_lup, rperm, info); if (info /= 0) error stop
54 call disp%show("mat_lup")
55 call disp%show( mat_lup )
56 call disp%show("call setResized(rperm, size(mat_lup, 1, IK)) ! reconstruct the original matrix.")
57 call setResized(rperm, size(mat_lup, 1, IK))
58 call disp%show("lup_ref ! reference matrix rounded to 1 significant digit.")
59 call disp%show( lup_ref )
60 call disp%show("lup_ref - mat_lup, format = SK_'(*(f0.1,:,"", ""))'")
61 call disp%show( lup_ref - mat_lup, format = SK_"(*(f0.1,:,"", ""))")
62 call disp%show("rperm_ref")
63 call disp%show( rperm_ref )
64 call disp%show("rperm")
65 call disp%show( rperm )
66 call disp%skip
67 end block
68
69 block
70 use pm_kind, only: TKG => RKD
71 real(TKG), allocatable :: mat_lup(:,:), lup_ref(:,:)
72 mat_lup = reshape( [ 1.0_TKG, +1.2_TKG, 1.4_TKG, 1.6_TKG, 1.8_TKG, 2.0_TKG, 2.2_TKG, 2.4_TKG, 2.6_TKG &
73 , 1.2_TKG, +1.0_TKG, 1.2_TKG, 1.4_TKG, 1.6_TKG, 1.8_TKG, 2.0_TKG, 2.2_TKG, 2.4_TKG &
74 , 1.4_TKG, +1.2_TKG, 1.0_TKG, 1.2_TKG, 1.4_TKG, 1.6_TKG, 1.8_TKG, 2.0_TKG, 2.2_TKG &
75 , 1.6_TKG, +1.4_TKG, 1.2_TKG, 1.0_TKG, 1.2_TKG, 1.4_TKG, 1.6_TKG, 1.8_TKG, 2.0_TKG &
76 , 1.8_TKG, +1.6_TKG, 1.4_TKG, 1.2_TKG, 1.0_TKG, 1.2_TKG, 1.4_TKG, 1.6_TKG, 1.8_TKG &
77 , 2.0_TKG, +1.8_TKG, 1.6_TKG, 1.4_TKG, 1.2_TKG, 1.0_TKG, 1.2_TKG, 1.4_TKG, 1.6_TKG &
78 , 2.2_TKG, +2.0_TKG, 1.8_TKG, 1.6_TKG, 1.4_TKG, 1.2_TKG, 1.0_TKG, 1.2_TKG, 1.4_TKG &
79 , 2.4_TKG, +2.2_TKG, 2.0_TKG, 1.8_TKG, 1.6_TKG, 1.4_TKG, 1.2_TKG, 1.0_TKG, 1.2_TKG &
80 , 2.6_TKG, +2.4_TKG, 2.2_TKG, 2.0_TKG, 1.8_TKG, 1.6_TKG, 1.4_TKG, 1.2_TKG, 1.0_TKG &
81 ], shape = [9, 9], order = [2, 1])
82 lup_ref = reshape( [ 2.6_TKG, 2.4_TKG, 2.2_TKG, 2.0_TKG, 1.8_TKG, 1.6_TKG, 1.4_TKG, 1.2_TKG, 1.0_TKG &
83 , 0.4_TKG, 0.3_TKG, 0.6_TKG, 0.8_TKG, 1.1_TKG, 1.4_TKG, 1.7_TKG, 1.9_TKG, 2.2_TKG &
84 , 0.5_TKG, -0.4_TKG, 0.4_TKG, 0.8_TKG, 1.2_TKG, 1.6_TKG, 2.0_TKG, 2.4_TKG, 2.8_TKG &
85 , 0.5_TKG, -0.3_TKG, 0.0_TKG, 0.4_TKG, 0.8_TKG, 1.2_TKG, 1.6_TKG, 2.0_TKG, 2.4_TKG &
86 , 0.6_TKG, -0.3_TKG, 0.0_TKG, 0.0_TKG, 0.4_TKG, 0.8_TKG, 1.2_TKG, 1.6_TKG, 2.0_TKG &
87 , 0.7_TKG, -0.2_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.4_TKG, 0.8_TKG, 1.2_TKG, 1.6_TKG &
88 , 0.8_TKG, -0.2_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.4_TKG, 0.8_TKG, 1.2_TKG &
89 , 0.8_TKG, -0.1_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.4_TKG, 0.8_TKG &
90 , 0.9_TKG, -0.1_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.4_TKG &
91 ], shape = [9, 9], order = [2, 1])
92 rperm_ref = [9, 9, 9, 9, 9, 9, 9, 9, 9]
93 call disp%skip
94 call disp%show("mat_lup")
95 call disp%show( mat_lup )
96 call disp%show("call setResized(rperm, size(mat_lup, 1, IK))")
97 call setResized(rperm, size(mat_lup, 1, IK))
98 call disp%show("call setMatLUP(mat_lup, rperm, info); if (info /= 0) error stop")
99 call setMatLUP(mat_lup, rperm, info); if (info /= 0) error stop
100 call disp%show("mat_lup")
101 call disp%show( mat_lup )
102 call disp%show("call setResized(rperm, size(mat_lup, 1, IK)) ! reconstruct the original matrix.")
103 call setResized(rperm, size(mat_lup, 1, IK))
104 call disp%show("lup_ref ! reference matrix rounded to 1 significant digit.")
105 call disp%show( lup_ref )
106 call disp%show("lup_ref - mat_lup, format = SK_'(*(f0.1,:,"", ""))'")
107 call disp%show( lup_ref - mat_lup, format = SK_"(*(f0.1,:,"", ""))")
108 call disp%show("rperm_ref")
109 call disp%show( rperm_ref )
110 call disp%show("rperm")
111 call disp%show( rperm )
112 call disp%skip
113 end block
114
115 block
116 use pm_kind, only: TKG => RKH
117 real(TKG), allocatable :: mat_lup(:,:), lup_ref(:,:)
118 mat_lup = reshape( [ 1.0_TKG, +1.2_TKG, 1.4_TKG, 1.6_TKG, 1.8_TKG, 2.0_TKG, 2.2_TKG, 2.4_TKG, 2.6_TKG &
119 , 1.2_TKG, +1.0_TKG, 1.2_TKG, 1.4_TKG, 1.6_TKG, 1.8_TKG, 2.0_TKG, 2.2_TKG, 2.4_TKG &
120 , 1.4_TKG, +1.2_TKG, 1.0_TKG, 1.2_TKG, 1.4_TKG, 1.6_TKG, 1.8_TKG, 2.0_TKG, 2.2_TKG &
121 , 1.6_TKG, +1.4_TKG, 1.2_TKG, 1.0_TKG, 1.2_TKG, 1.4_TKG, 1.6_TKG, 1.8_TKG, 2.0_TKG &
122 , 1.8_TKG, +1.6_TKG, 1.4_TKG, 1.2_TKG, 1.0_TKG, 1.2_TKG, 1.4_TKG, 1.6_TKG, 1.8_TKG &
123 , 2.0_TKG, +1.8_TKG, 1.6_TKG, 1.4_TKG, 1.2_TKG, 1.0_TKG, 1.2_TKG, 1.4_TKG, 1.6_TKG &
124 , 2.2_TKG, +2.0_TKG, 1.8_TKG, 1.6_TKG, 1.4_TKG, 1.2_TKG, 1.0_TKG, 1.2_TKG, 1.4_TKG &
125 , 2.4_TKG, +2.2_TKG, 2.0_TKG, 1.8_TKG, 1.6_TKG, 1.4_TKG, 1.2_TKG, 1.0_TKG, 1.2_TKG &
126 , 2.6_TKG, +2.4_TKG, 2.2_TKG, 2.0_TKG, 1.8_TKG, 1.6_TKG, 1.4_TKG, 1.2_TKG, 1.0_TKG &
127 ], shape = [9, 9], order = [2, 1])
128 lup_ref = reshape( [ 2.6_TKG, 2.4_TKG, 2.2_TKG, 2.0_TKG, 1.8_TKG, 1.6_TKG, 1.4_TKG, 1.2_TKG, 1.0_TKG &
129 , 0.4_TKG, 0.3_TKG, 0.6_TKG, 0.8_TKG, 1.1_TKG, 1.4_TKG, 1.7_TKG, 1.9_TKG, 2.2_TKG &
130 , 0.5_TKG, -0.4_TKG, 0.4_TKG, 0.8_TKG, 1.2_TKG, 1.6_TKG, 2.0_TKG, 2.4_TKG, 2.8_TKG &
131 , 0.5_TKG, -0.3_TKG, 0.0_TKG, 0.4_TKG, 0.8_TKG, 1.2_TKG, 1.6_TKG, 2.0_TKG, 2.4_TKG &
132 , 0.6_TKG, -0.3_TKG, 0.0_TKG, 0.0_TKG, 0.4_TKG, 0.8_TKG, 1.2_TKG, 1.6_TKG, 2.0_TKG &
133 , 0.7_TKG, -0.2_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.4_TKG, 0.8_TKG, 1.2_TKG, 1.6_TKG &
134 , 0.8_TKG, -0.2_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.4_TKG, 0.8_TKG, 1.2_TKG &
135 , 0.8_TKG, -0.1_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.4_TKG, 0.8_TKG &
136 , 0.9_TKG, -0.1_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.0_TKG, 0.4_TKG &
137 ], shape = [9, 9], order = [2, 1])
138 rperm_ref = [9, 9, 9, 9, 9, 9, 9, 9, 9]
139 call disp%skip
140 call disp%show("mat_lup")
141 call disp%show( mat_lup )
142 call disp%show("call setResized(rperm, size(mat_lup, 1, IK))")
143 call setResized(rperm, size(mat_lup, 1, IK))
144 call disp%show("call setMatLUP(mat_lup, rperm, info); if (info /= 0) error stop")
145 call setMatLUP(mat_lup, rperm, info); if (info /= 0) error stop
146 call disp%show("mat_lup")
147 call disp%show( mat_lup )
148 call disp%show("call setResized(rperm, size(mat_lup, 1, IK)) ! reconstruct the original matrix.")
149 call setResized(rperm, size(mat_lup, 1, IK))
150 call disp%show("lup_ref ! reference matrix rounded to 1 significant digit.")
151 call disp%show( lup_ref )
152 call disp%show("lup_ref - mat_lup, format = SK_'(*(f0.1,:,"", ""))'")
153 call disp%show( lup_ref - mat_lup, format = SK_"(*(f0.1,:,"", ""))")
154 call disp%show("rperm_ref")
155 call disp%show( rperm_ref )
156 call disp%show("rperm")
157 call disp%show( rperm )
158 call disp%skip
159 end block
160
161#if 1 || LAPACK_ENABLED
162 call disp%skip
163 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
164 call disp%show("! Compute the LU-Pivoted decomposition of a complex square matrix.")
165 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
166 call disp%skip
167
168 block
169 use pm_kind, only: TKG => CKS
170 complex(TKG), allocatable :: mat_lup(:,:), lup_ref(:,:)
171 mat_lup = reshape( [ (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) &
172 , (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) &
173 , (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) &
174 , (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) &
175 , (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) &
176 , (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) &
177 , (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) &
178 , (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) &
179 , (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) &
180 ], shape = [9, 9], order = [2, 1])
181 lup_ref = reshape( [ (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) &
182 , (0.4, 0.1), (0.6,-2.0), (1.1,-1.9), (1.7,-1.9), (+2.3,-1.8), (2.8,-1.8), (3.4,-1.7), (3.9,-1.7), (4.5,-1.6) &
183 , (0.5, 0.1), (0.0,-0.1), (0.6,-1.9), (1.2,-1.8), (+1.8,-1.7), (2.5,-1.6), (3.1,-1.5), (3.7,-1.4), (4.3,-1.3) &
184 , (0.6, 0.1), (0.0,-0.1),(-0.1,-0.1), (0.7,-1.9), (+1.3,-1.7), (2.0,-1.6), (2.7,-1.5), (3.4,-1.4), (4.0,-1.2) &
185 , (0.6, 0.1), (0.0,-0.1),(-0.1,-0.1),(-0.1, 0.0), (+0.7,-1.9), (1.5,-1.7), (2.2,-1.6), (2.9,-1.5), (3.7,-1.3) &
186 , (0.7, 0.1), (0.0,-0.1), (0.0, 0.0),(-0.1, 0.0), (-0.1, 0.0), (0.8,-1.9), (1.6,-1.8), (2.4,-1.6), (3.2,-1.5) &
187 , (0.8, 0.0), (0.0, 0.0), (0.0, 0.0), (0.0, 0.0), (+0.0, 0.0), (0.0, 0.0), (0.8,-1.9), (1.7,-1.8), (2.5,-1.8) &
188 , (0.9, 0.0), (0.0, 0.0), (0.0, 0.0), (0.0, 0.0), (+0.0, 0.0), (0.0, 0.0), (0.0, 0.0), (0.8,-2.0), (1.7,-1.9) &
189 , (0.9, 0.0), (0.0, 0.0), (0.0, 0.0), (0.0, 0.0), (+0.0, 0.0), (0.0, 0.0), (0.0, 0.0), (0.0, 0.0), (0.8,-2.0) &
190 ], shape = [9, 9], order = [2, 1])
191 rperm_ref = [9, 9, 9, 9, 9, 9, 9, 9, 9]
192 call disp%skip
193 call disp%show("mat_lup")
194 call disp%show( mat_lup )
195 call disp%show("call setResized(rperm, size(mat_lup, 1, IK))")
196 call setResized(rperm, size(mat_lup, 1, IK))
197 call disp%show("call setMatLUP(mat_lup, rperm, info); if (info /= 0) error stop")
198 call setMatLUP(mat_lup, rperm, info); if (info /= 0) error stop
199 call disp%show("mat_lup")
200 call disp%show( mat_lup )
201 call disp%show("call setResized(rperm, size(mat_lup, 1, IK)) ! reconstruct the original matrix.")
202 call setResized(rperm, size(mat_lup, 1, IK))
203 call disp%show("lup_ref ! reference matrix rounded to 1 significant digit.")
204 call disp%show( lup_ref )
205 call disp%show("lup_ref - mat_lup, format = SK_'(*(""("",f0.1,"","",f0.1,"")"",:,"", ""))'")
206 call disp%show( lup_ref - mat_lup, format = SK_"(*(""("",f0.1,"","",f0.1,"")"",:,"", ""))")
207 call disp%show("rperm_ref")
208 call disp%show( rperm_ref )
209 call disp%show("rperm")
210 call disp%show( rperm )
211 call disp%skip
212 end block
213
214 block
215 use pm_kind, only: TKG => CKD
216 complex(TKG), allocatable :: mat_lup(:,:), lup_ref(:,:)
217 mat_lup = reshape( [ (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) &
218 , (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) &
219 , (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) &
220 , (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) &
221 , (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) &
222 , (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) &
223 , (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) &
224 , (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) &
225 , (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) &
226 ], shape = [9, 9], order = [2, 1])
227 lup_ref = reshape( [ (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) &
228 , (0.4, 0.1), (0.6,-2.0), (1.1,-1.9), (1.7,-1.9), (+2.3,-1.8), (2.8,-1.8), (3.4,-1.7), (3.9,-1.7), (4.5,-1.6) &
229 , (0.5, 0.1), (0.0,-0.1), (0.6,-1.9), (1.2,-1.8), (+1.8,-1.7), (2.5,-1.6), (3.1,-1.5), (3.7,-1.4), (4.3,-1.3) &
230 , (0.6, 0.1), (0.0,-0.1),(-0.1,-0.1), (0.7,-1.9), (+1.3,-1.7), (2.0,-1.6), (2.7,-1.5), (3.4,-1.4), (4.0,-1.2) &
231 , (0.6, 0.1), (0.0,-0.1),(-0.1,-0.1),(-0.1, 0.0), (+0.7,-1.9), (1.5,-1.7), (2.2,-1.6), (2.9,-1.5), (3.7,-1.3) &
232 , (0.7, 0.1), (0.0,-0.1), (0.0, 0.0),(-0.1, 0.0), (-0.1, 0.0), (0.8,-1.9), (1.6,-1.8), (2.4,-1.6), (3.2,-1.5) &
233 , (0.8, 0.0), (0.0, 0.0), (0.0, 0.0), (0.0, 0.0), (+0.0, 0.0), (0.0, 0.0), (0.8,-1.9), (1.7,-1.8), (2.5,-1.8) &
234 , (0.9, 0.0), (0.0, 0.0), (0.0, 0.0), (0.0, 0.0), (+0.0, 0.0), (0.0, 0.0), (0.0, 0.0), (0.8,-2.0), (1.7,-1.9) &
235 , (0.9, 0.0), (0.0, 0.0), (0.0, 0.0), (0.0, 0.0), (+0.0, 0.0), (0.0, 0.0), (0.0, 0.0), (0.0, 0.0), (0.8,-2.0) &
236 ], shape = [9, 9], order = [2, 1])
237 rperm_ref = [9, 9, 9, 9, 9, 9, 9, 9, 9]
238 call disp%skip
239 call disp%show("mat_lup")
240 call disp%show( mat_lup )
241 call disp%show("call setResized(rperm, size(mat_lup, 1, IK))")
242 call setResized(rperm, size(mat_lup, 1, IK))
243 call disp%show("call setMatLUP(mat_lup, rperm, info); if (info /= 0) error stop")
244 call setMatLUP(mat_lup, rperm, info); if (info /= 0) error stop
245 call disp%show("mat_lup")
246 call disp%show( mat_lup )
247 call disp%show("call setResized(rperm, size(mat_lup, 1, IK)) ! reconstruct the original matrix.")
248 call setResized(rperm, size(mat_lup, 1, IK))
249 call disp%show("lup_ref ! reference matrix rounded to 1 significant digit.")
250 call disp%show( lup_ref )
251 call disp%show("lup_ref - mat_lup, format = SK_'(*(""("",f0.1,"","",f0.1,"")"",:,"", ""))'")
252 call disp%show( lup_ref - mat_lup, format = SK_"(*(""("",f0.1,"","",f0.1,"")"",:,"", ""))")
253 call disp%show("rperm_ref")
254 call disp%show( rperm_ref )
255 call disp%show("rperm")
256 call disp%show( rperm )
257 call disp%skip
258 end block
259
260 block
261 use pm_kind, only: TKG => CKH
262 complex(TKG), allocatable :: mat_lup(:,:), lup_ref(:,:)
263 mat_lup = reshape( [ (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) &
264 , (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) &
265 , (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) &
266 , (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) &
267 , (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) &
268 , (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) &
269 , (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) &
270 , (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) &
271 , (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) &
272 ], shape = [9, 9], order = [2, 1])
273 lup_ref = reshape( [ (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) &
274 , (0.4, 0.1), (0.6,-2.0), (1.1,-1.9), (1.7,-1.9), (+2.3,-1.8), (2.8,-1.8), (3.4,-1.7), (3.9,-1.7), (4.5,-1.6) &
275 , (0.5, 0.1), (0.0,-0.1), (0.6,-1.9), (1.2,-1.8), (+1.8,-1.7), (2.5,-1.6), (3.1,-1.5), (3.7,-1.4), (4.3,-1.3) &
276 , (0.6, 0.1), (0.0,-0.1),(-0.1,-0.1), (0.7,-1.9), (+1.3,-1.7), (2.0,-1.6), (2.7,-1.5), (3.4,-1.4), (4.0,-1.2) &
277 , (0.6, 0.1), (0.0,-0.1),(-0.1,-0.1),(-0.1, 0.0), (+0.7,-1.9), (1.5,-1.7), (2.2,-1.6), (2.9,-1.5), (3.7,-1.3) &
278 , (0.7, 0.1), (0.0,-0.1), (0.0, 0.0),(-0.1, 0.0), (-0.1, 0.0), (0.8,-1.9), (1.6,-1.8), (2.4,-1.6), (3.2,-1.5) &
279 , (0.8, 0.0), (0.0, 0.0), (0.0, 0.0), (0.0, 0.0), (+0.0, 0.0), (0.0, 0.0), (0.8,-1.9), (1.7,-1.8), (2.5,-1.8) &
280 , (0.9, 0.0), (0.0, 0.0), (0.0, 0.0), (0.0, 0.0), (+0.0, 0.0), (0.0, 0.0), (0.0, 0.0), (0.8,-2.0), (1.7,-1.9) &
281 , (0.9, 0.0), (0.0, 0.0), (0.0, 0.0), (0.0, 0.0), (+0.0, 0.0), (0.0, 0.0), (0.0, 0.0), (0.0, 0.0), (0.8,-2.0) &
282 ], shape = [9, 9], order = [2, 1])
283 rperm_ref = [9, 9, 9, 9, 9, 9, 9, 9, 9]
284 call disp%skip
285 call disp%show("mat_lup")
286 call disp%show( mat_lup )
287 call disp%show("call setResized(rperm, size(mat_lup, 1, IK))")
288 call setResized(rperm, size(mat_lup, 1, IK))
289 call disp%show("call setMatLUP(mat_lup, rperm, info); if (info /= 0) error stop")
290 call setMatLUP(mat_lup, rperm, info); if (info /= 0) error stop
291 call disp%show("mat_lup")
292 call disp%show( mat_lup )
293 call disp%show("call setResized(rperm, size(mat_lup, 1, IK)) ! reconstruct the original matrix.")
294 call setResized(rperm, size(mat_lup, 1, IK))
295 call disp%show("lup_ref ! reference matrix rounded to 1 significant digit.")
296 call disp%show( lup_ref )
297 call disp%show("lup_ref - mat_lup, format = SK_'(*(""("",f0.1,"","",f0.1,"")"",:,"", ""))'")
298 call disp%show( lup_ref - mat_lup, format = SK_"(*(""("",f0.1,"","",f0.1,"")"",:,"", ""))")
299 call disp%show("rperm_ref")
300 call disp%show( rperm_ref )
301 call disp%show("rperm")
302 call disp%show( rperm )
303 call disp%skip
304 end block
305#endif
306
307end program example
Allocate or resize (shrink or expand) an input allocatable scalar string or array of rank 1....
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
Copy a desired subset of the input source matrix of arbitrary shape (:) or (:,:) to the target subset...
This module contains procedures and generic interfaces for resizing allocatable arrays of various typ...
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 CKH
The scalar integer constant of intrinsic default kind, representing the highest-precision complex kin...
Definition: pm_kind.F90:843
integer, parameter CKS
The single-precision complex kind in Fortran mode. On most platforms, this is a 32-bit real kind.
Definition: pm_kind.F90:570
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 CKD
The double precision complex kind in Fortran mode. On most platforms, this is a 64-bit real kind.
Definition: pm_kind.F90:571
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
This module contains procedures and generic interfaces relevant to copying (diagonal or upper/lower t...
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! Compute the LU-Pivoted decomposition of a real square matrix.
4!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6
7mat_lup
8+1.00000000, +1.20000005, +1.39999998, +1.60000002, +1.79999995, +2.00000000, +2.20000005, +2.40000010, +2.59999990
9+1.20000005, +1.00000000, +1.20000005, +1.39999998, +1.60000002, +1.79999995, +2.00000000, +2.20000005, +2.40000010
10+1.39999998, +1.20000005, +1.00000000, +1.20000005, +1.39999998, +1.60000002, +1.79999995, +2.00000000, +2.20000005
11+1.60000002, +1.39999998, +1.20000005, +1.00000000, +1.20000005, +1.39999998, +1.60000002, +1.79999995, +2.00000000
12+1.79999995, +1.60000002, +1.39999998, +1.20000005, +1.00000000, +1.20000005, +1.39999998, +1.60000002, +1.79999995
13+2.00000000, +1.79999995, +1.60000002, +1.39999998, +1.20000005, +1.00000000, +1.20000005, +1.39999998, +1.60000002
14+2.20000005, +2.00000000, +1.79999995, +1.60000002, +1.39999998, +1.20000005, +1.00000000, +1.20000005, +1.39999998
15+2.40000010, +2.20000005, +2.00000000, +1.79999995, +1.60000002, +1.39999998, +1.20000005, +1.00000000, +1.20000005
16+2.59999990, +2.40000010, +2.20000005, +2.00000000, +1.79999995, +1.60000002, +1.39999998, +1.20000005, +1.00000000
17call setResized(rperm, size(mat_lup, 1, IK))
18call setMatLUP(mat_lup, rperm, info); if (info /= 0) error stop
19mat_lup
20+2.59999990, +2.40000010, +2.20000005, +2.00000000, +1.79999995, +1.60000002, +1.39999998, +1.20000005, +1.00000000
21+0.384615391, +0.276923060, +0.553846121, +0.830769241, +1.10769224, +1.38461542, +1.66153848, +1.93846154, +2.21538448
22+0.461538494, -0.388889551, +0.400000334, +0.800000429, +1.20000076, +1.60000086, +2.00000095, +2.40000129, +2.80000138
23+0.538461566, -0.333333760, +0.298022968E-6, +0.400000125, +0.800000072, +1.20000005, +1.60000002, +2.00000000, +2.40000033
24+0.615384638, -0.277778417, +0.894068876E-6, -0.558792806E-6, +0.400000066, +0.800000012, +1.20000017, +1.60000002, +2.00000024
25+0.692307711, -0.222222656, +0.298022968E-6, +0.335276212E-6, -0.372528973E-6, +0.400000006, +0.799999952, +1.19999993, +1.59999990
26+0.769230783, -0.166667312, +0.894068876E-6, -0.558792806E-6, +0.223516594E-6, -0.223516778E-6, +0.400000095, +0.800000012, +1.20000017
27+0.846153855, -0.111111544, +0.298022968E-6, +0.316649761E-6, -0.335276070E-6, +0.558793545E-7, -0.372527467E-7, +0.399999976, +0.799999952
28+0.923076987, -0.555566326E-1, +0.119209187E-5, -0.586732369E-6, +0.577419087E-6, -0.866129369E-6, +0.875443050E-6, -0.894069899E-6, +0.400000304
29call setResized(rperm, size(mat_lup, 1, IK)) ! reconstruct the original matrix.
30lup_ref ! reference matrix rounded to 1 significant digit.
31+2.59999990, +2.40000010, +2.20000005, +2.00000000, +1.79999995, +1.60000002, +1.39999998, +1.20000005, +1.00000000
32+0.400000006, +0.300000012, +0.600000024, +0.800000012, +1.10000002, +1.39999998, +1.70000005, +1.89999998, +2.20000005
33+0.500000000, -0.400000006, +0.400000006, +0.800000012, +1.20000005, +1.60000002, +2.00000000, +2.40000010, +2.79999995
34+0.500000000, -0.300000012, +0.00000000, +0.400000006, +0.800000012, +1.20000005, +1.60000002, +2.00000000, +2.40000010
35+0.600000024, -0.300000012, +0.00000000, +0.00000000, +0.400000006, +0.800000012, +1.20000005, +1.60000002, +2.00000000
36+0.699999988, -0.200000003, +0.00000000, +0.00000000, +0.00000000, +0.400000006, +0.800000012, +1.20000005, +1.60000002
37+0.800000012, -0.200000003, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.400000006, +0.800000012, +1.20000005
38+0.800000012, -0.100000001, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.400000006, +0.800000012
39+0.899999976, -0.100000001, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.400000006
40lup_ref - mat_lup, format = SK_'(*(f0.1,:,", "))'
41.0, .0, .0, .0, .0, .0, .0, .0, .0
42.0, .0, .0, -.0, -.0, .0, .0, -.0, -.0
43.0, -.0, -.0, -.0, -.0, -.0, -.0, -.0, -.0
44-.0, .0, -.0, -.0, -.0, .0, .0, .0, -.0
45-.0, -.0, -.0, .0, -.0, .0, -.0, .0, -.0
46.0, .0, -.0, -.0, .0, .0, .0, .0, .0
47.0, -.0, -.0, .0, -.0, .0, -.0, .0, -.0
48-.0, .0, -.0, -.0, .0, -.0, .0, .0, .0
49-.0, -.0, -.0, .0, -.0, .0, -.0, .0, -.0
50rperm_ref
51+9, +9, +9, +9, +9, +9, +9, +9, +9
52rperm
53+9, +9, +9, +9, +9, +9, +9, +9, +9
54
55
56mat_lup
57+1.0000000000000000, +1.2000000000000000, +1.3999999999999999, +1.6000000000000001, +1.8000000000000000, +2.0000000000000000, +2.2000000000000002, +2.3999999999999999, +2.6000000000000001
58+1.2000000000000000, +1.0000000000000000, +1.2000000000000000, +1.3999999999999999, +1.6000000000000001, +1.8000000000000000, +2.0000000000000000, +2.2000000000000002, +2.3999999999999999
59+1.3999999999999999, +1.2000000000000000, +1.0000000000000000, +1.2000000000000000, +1.3999999999999999, +1.6000000000000001, +1.8000000000000000, +2.0000000000000000, +2.2000000000000002
60+1.6000000000000001, +1.3999999999999999, +1.2000000000000000, +1.0000000000000000, +1.2000000000000000, +1.3999999999999999, +1.6000000000000001, +1.8000000000000000, +2.0000000000000000
61+1.8000000000000000, +1.6000000000000001, +1.3999999999999999, +1.2000000000000000, +1.0000000000000000, +1.2000000000000000, +1.3999999999999999, +1.6000000000000001, +1.8000000000000000
62+2.0000000000000000, +1.8000000000000000, +1.6000000000000001, +1.3999999999999999, +1.2000000000000000, +1.0000000000000000, +1.2000000000000000, +1.3999999999999999, +1.6000000000000001
63+2.2000000000000002, +2.0000000000000000, +1.8000000000000000, +1.6000000000000001, +1.3999999999999999, +1.2000000000000000, +1.0000000000000000, +1.2000000000000000, +1.3999999999999999
64+2.3999999999999999, +2.2000000000000002, +2.0000000000000000, +1.8000000000000000, +1.6000000000000001, +1.3999999999999999, +1.2000000000000000, +1.0000000000000000, +1.2000000000000000
65+2.6000000000000001, +2.3999999999999999, +2.2000000000000002, +2.0000000000000000, +1.8000000000000000, +1.6000000000000001, +1.3999999999999999, +1.2000000000000000, +1.0000000000000000
66call setResized(rperm, size(mat_lup, 1, IK))
67call setMatLUP(mat_lup, rperm, info); if (info /= 0) error stop
68mat_lup
69+2.6000000000000001, +2.3999999999999999, +2.2000000000000002, +2.0000000000000000, +1.8000000000000000, +1.6000000000000001, +1.3999999999999999, +1.2000000000000000, +1.0000000000000000
70+0.38461538461538458, +0.27692307692307705, +0.55384615384615377, +0.83076923076923093, +1.1076923076923078, +1.3846153846153846, +1.6615384615384619, +1.9384615384615385, +2.2153846153846155
71+0.46153846153846145, -0.38888888888888812, +0.39999999999999958, +0.79999999999999949, +1.1999999999999995, +1.5999999999999992, +1.9999999999999989, +2.3999999999999990, +2.7999999999999985
72+0.53846153846153832, -0.33333333333333198, -0.13877787807814473E-14, +0.40000000000000024, +0.80000000000000038, +1.2000000000000006, +1.6000000000000010, +2.0000000000000009, +2.4000000000000012
73+0.61538461538461542, -0.27777777777777823, +0.27755575615628943E-15, +0.34694469519536098E-15, +0.39999999999999980, +0.79999999999999960, +1.1999999999999997, +1.5999999999999996, +1.9999999999999991
74+0.69230769230769229, -0.22222222222222132, -0.18735013540549536E-14, +0.19428902930940243E-14, -0.76327832942979404E-15, +0.39999999999999991, +0.80000000000000016, +1.2000000000000000, +1.6000000000000003
75+0.76923076923076916, -0.16666666666666519, -0.18041124150158814E-14, +0.69388939039072392E-15, -0.69388939039071088E-16, -0.69388939039075254E-16, +0.40000000000000041, +0.80000000000000027, +1.2000000000000006
76+0.84615384615384615, -0.11111111111111066, -0.65919492087118739E-15, +0.69388939039072294E-15, -0.65919492087118660E-15, -0.73955709864469874E-30, +0.69388939039072234E-15, +0.39999999999999986, +0.79999999999999971
77+0.92307692307692291, -0.55555555555553728E-1, -0.22724877535296197E-14, +0.17347234759768081E-14, -0.60715321659188229E-15, -0.11102230246251597E-14, +0.17347234759768093E-14, -0.12143064331837640E-14, +0.40000000000000013
78call setResized(rperm, size(mat_lup, 1, IK)) ! reconstruct the original matrix.
79lup_ref ! reference matrix rounded to 1 significant digit.
80+2.6000000000000001, +2.3999999999999999, +2.2000000000000002, +2.0000000000000000, +1.8000000000000000, +1.6000000000000001, +1.3999999999999999, +1.2000000000000000, +1.0000000000000000
81+0.40000000000000002, +0.29999999999999999, +0.59999999999999998, +0.80000000000000004, +1.1000000000000001, +1.3999999999999999, +1.7000000000000000, +1.8999999999999999, +2.2000000000000002
82+0.50000000000000000, -0.40000000000000002, +0.40000000000000002, +0.80000000000000004, +1.2000000000000000, +1.6000000000000001, +2.0000000000000000, +2.3999999999999999, +2.7999999999999998
83+0.50000000000000000, -0.29999999999999999, +0.0000000000000000, +0.40000000000000002, +0.80000000000000004, +1.2000000000000000, +1.6000000000000001, +2.0000000000000000, +2.3999999999999999
84+0.59999999999999998, -0.29999999999999999, +0.0000000000000000, +0.0000000000000000, +0.40000000000000002, +0.80000000000000004, +1.2000000000000000, +1.6000000000000001, +2.0000000000000000
85+0.69999999999999996, -0.20000000000000001, +0.0000000000000000, +0.0000000000000000, +0.0000000000000000, +0.40000000000000002, +0.80000000000000004, +1.2000000000000000, +1.6000000000000001
86+0.80000000000000004, -0.20000000000000001, +0.0000000000000000, +0.0000000000000000, +0.0000000000000000, +0.0000000000000000, +0.40000000000000002, +0.80000000000000004, +1.2000000000000000
87+0.80000000000000004, -0.10000000000000001, +0.0000000000000000, +0.0000000000000000, +0.0000000000000000, +0.0000000000000000, +0.0000000000000000, +0.40000000000000002, +0.80000000000000004
88+0.90000000000000002, -0.10000000000000001, +0.0000000000000000, +0.0000000000000000, +0.0000000000000000, +0.0000000000000000, +0.0000000000000000, +0.0000000000000000, +0.40000000000000002
89lup_ref - mat_lup, format = SK_'(*(f0.1,:,", "))'
90.0, .0, .0, .0, .0, .0, .0, .0, .0
91.0, .0, .0, -.0, -.0, .0, .0, -.0, -.0
92.0, -.0, .0, .0, .0, .0, .0, .0, .0
93-.0, .0, .0, -.0, -.0, -.0, -.0, -.0, -.0
94-.0, -.0, -.0, -.0, .0, .0, .0, .0, .0
95.0, .0, .0, -.0, .0, .0, -.0, .0, -.0
96.0, -.0, .0, -.0, .0, .0, -.0, -.0, -.0
97-.0, .0, .0, -.0, .0, .0, -.0, .0, .0
98-.0, -.0, .0, -.0, .0, .0, -.0, .0, -.0
99rperm_ref
100+9, +9, +9, +9, +9, +9, +9, +9, +9
101rperm
102+9, +9, +9, +9, +9, +9, +9, +9, +9
103
104
105mat_lup
106+1.00000000000000000000000000000000000, +1.19999999999999999999999999999999996, +1.39999999999999999999999999999999992, +1.60000000000000000000000000000000008, +1.80000000000000000000000000000000004, +2.00000000000000000000000000000000000, +2.20000000000000000000000000000000015, +2.39999999999999999999999999999999992, +2.60000000000000000000000000000000008
107+1.19999999999999999999999999999999996, +1.00000000000000000000000000000000000, +1.19999999999999999999999999999999996, +1.39999999999999999999999999999999992, +1.60000000000000000000000000000000008, +1.80000000000000000000000000000000004, +2.00000000000000000000000000000000000, +2.20000000000000000000000000000000015, +2.39999999999999999999999999999999992
108+1.39999999999999999999999999999999992, +1.19999999999999999999999999999999996, +1.00000000000000000000000000000000000, +1.19999999999999999999999999999999996, +1.39999999999999999999999999999999992, +1.60000000000000000000000000000000008, +1.80000000000000000000000000000000004, +2.00000000000000000000000000000000000, +2.20000000000000000000000000000000015
109+1.60000000000000000000000000000000008, +1.39999999999999999999999999999999992, +1.19999999999999999999999999999999996, +1.00000000000000000000000000000000000, +1.19999999999999999999999999999999996, +1.39999999999999999999999999999999992, +1.60000000000000000000000000000000008, +1.80000000000000000000000000000000004, +2.00000000000000000000000000000000000
110+1.80000000000000000000000000000000004, +1.60000000000000000000000000000000008, +1.39999999999999999999999999999999992, +1.19999999999999999999999999999999996, +1.00000000000000000000000000000000000, +1.19999999999999999999999999999999996, +1.39999999999999999999999999999999992, +1.60000000000000000000000000000000008, +1.80000000000000000000000000000000004
111+2.00000000000000000000000000000000000, +1.80000000000000000000000000000000004, +1.60000000000000000000000000000000008, +1.39999999999999999999999999999999992, +1.19999999999999999999999999999999996, +1.00000000000000000000000000000000000, +1.19999999999999999999999999999999996, +1.39999999999999999999999999999999992, +1.60000000000000000000000000000000008
112+2.20000000000000000000000000000000015, +2.00000000000000000000000000000000000, +1.80000000000000000000000000000000004, +1.60000000000000000000000000000000008, +1.39999999999999999999999999999999992, +1.19999999999999999999999999999999996, +1.00000000000000000000000000000000000, +1.19999999999999999999999999999999996, +1.39999999999999999999999999999999992
113+2.39999999999999999999999999999999992, +2.20000000000000000000000000000000015, +2.00000000000000000000000000000000000, +1.80000000000000000000000000000000004, +1.60000000000000000000000000000000008, +1.39999999999999999999999999999999992, +1.19999999999999999999999999999999996, +1.00000000000000000000000000000000000, +1.19999999999999999999999999999999996
114+2.60000000000000000000000000000000008, +2.39999999999999999999999999999999992, +2.20000000000000000000000000000000015, +2.00000000000000000000000000000000000, +1.80000000000000000000000000000000004, +1.60000000000000000000000000000000008, +1.39999999999999999999999999999999992, +1.19999999999999999999999999999999996, +1.00000000000000000000000000000000000
115call setResized(rperm, size(mat_lup, 1, IK))
116call setMatLUP(mat_lup, rperm, info); if (info /= 0) error stop
117mat_lup
118+2.60000000000000000000000000000000008, +2.39999999999999999999999999999999992, +2.20000000000000000000000000000000015, +2.00000000000000000000000000000000000, +1.80000000000000000000000000000000004, +1.60000000000000000000000000000000008, +1.39999999999999999999999999999999992, +1.19999999999999999999999999999999996, +1.00000000000000000000000000000000000
119+0.384615384615384615384615384615384586, +0.276923076923076923076923076923077033, +0.553846153846153846153846153846153777, +0.830769230769230769230769230769230906, +1.10769230769230769230769230769230775, +1.38461538461538461538461538461538459, +1.66153846153846153846153846153846181, +1.93846153846153846153846153846153846, +2.21538461538461538461538461538461549
120+0.461538461538461538461538461538461464, -0.388888888888888888888888888888888220, +0.399999999999999999999999999999999634, +0.799999999999999999999999999999999557, +1.19999999999999999999999999999999958, +1.59999999999999999999999999999999931, +1.99999999999999999999999999999999904, +2.39999999999999999999999999999999915, +2.79999999999999999999999999999999869
121+0.538461538461538461538461538461538343, -0.333333333333333333333333333333332162, -0.120370621524202240815998621411558076E-32, +0.400000000000000000000000000000000212, +0.800000000000000000000000000000000327, +1.20000000000000000000000000000000054, +1.60000000000000000000000000000000085, +2.00000000000000000000000000000000077, +2.40000000000000000000000000000000108
122+0.615384615384615384615384615384615414, -0.277777777777777777777777777777778174, +0.240741243048404481631997242823116137E-33, +0.300926553810505602039996553528894560E-33, +0.399999999999999999999999999999999827, +0.799999999999999999999999999999999653, +1.19999999999999999999999999999999977, +1.59999999999999999999999999999999969, +1.99999999999999999999999999999999923
123+0.692307692307692307692307692307692293, -0.222222222222222222222222222222221441, -0.162500339057673025101598138905603391E-32, +0.168518870133883137142398069976181170E-32, -0.662038418383112324487992417763567950E-33, +0.399999999999999999999999999999999923, +0.800000000000000000000000000000000135, +1.19999999999999999999999999999999996, +1.60000000000000000000000000000000027
124+0.769230769230769230769230769230769172, -0.166666666666666666666666666666665383, -0.156481807981462913060798207835025493E-32, +0.601853107621011204079993107057790603E-33, -0.601853107621011204079993107057780792E-34, -0.601853107621011204079993107057812135E-34, +0.400000000000000000000000000000000356, +0.800000000000000000000000000000000231, +1.20000000000000000000000000000000054
125+0.846153846153846153846153846153846146, -0.111111111111111111111111111111110721, -0.571760452239960643875993451704900817E-33, +0.601853107621011204079993107057789861E-33, -0.571760452239960643875993451704900224E-33, -0.556380922603113207859760289232146041E-66, +0.601853107621011204079993107057789416E-33, +0.399999999999999999999999999999999875, +0.799999999999999999999999999999999750
126+0.923076923076923076923076923076922929, -0.555555555555555555555555555555539700E-1, -0.197106892745881169336197742561426333E-32, +0.150463276905252801019998276764447521E-32, -0.526621469168384803569993968675565915E-33, -0.962964972193617926527988971292466033E-33, +0.150463276905252801019998276764447610E-32, -0.105324293833676960713998793735113139E-32, +0.400000000000000000000000000000000116
127call setResized(rperm, size(mat_lup, 1, IK)) ! reconstruct the original matrix.
128lup_ref ! reference matrix rounded to 1 significant digit.
129+2.60000000000000000000000000000000008, +2.39999999999999999999999999999999992, +2.20000000000000000000000000000000015, +2.00000000000000000000000000000000000, +1.80000000000000000000000000000000004, +1.60000000000000000000000000000000008, +1.39999999999999999999999999999999992, +1.19999999999999999999999999999999996, +1.00000000000000000000000000000000000
130+0.400000000000000000000000000000000019, +0.299999999999999999999999999999999990, +0.599999999999999999999999999999999981, +0.800000000000000000000000000000000039, +1.10000000000000000000000000000000008, +1.39999999999999999999999999999999992, +1.69999999999999999999999999999999996, +1.89999999999999999999999999999999992, +2.20000000000000000000000000000000015
131+0.500000000000000000000000000000000000, -0.400000000000000000000000000000000019, +0.400000000000000000000000000000000019, +0.800000000000000000000000000000000039, +1.19999999999999999999999999999999996, +1.60000000000000000000000000000000008, +2.00000000000000000000000000000000000, +2.39999999999999999999999999999999992, +2.79999999999999999999999999999999985
132+0.500000000000000000000000000000000000, -0.299999999999999999999999999999999990, +0.00000000000000000000000000000000000, +0.400000000000000000000000000000000019, +0.800000000000000000000000000000000039, +1.19999999999999999999999999999999996, +1.60000000000000000000000000000000008, +2.00000000000000000000000000000000000, +2.39999999999999999999999999999999992
133+0.599999999999999999999999999999999981, -0.299999999999999999999999999999999990, +0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000, +0.400000000000000000000000000000000019, +0.800000000000000000000000000000000039, +1.19999999999999999999999999999999996, +1.60000000000000000000000000000000008, +2.00000000000000000000000000000000000
134+0.699999999999999999999999999999999961, -0.200000000000000000000000000000000010, +0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000, +0.400000000000000000000000000000000019, +0.800000000000000000000000000000000039, +1.19999999999999999999999999999999996, +1.60000000000000000000000000000000008
135+0.800000000000000000000000000000000039, -0.200000000000000000000000000000000010, +0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000, +0.400000000000000000000000000000000019, +0.800000000000000000000000000000000039, +1.19999999999999999999999999999999996
136+0.800000000000000000000000000000000039, -0.100000000000000000000000000000000005, +0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000, +0.400000000000000000000000000000000019, +0.800000000000000000000000000000000039
137+0.900000000000000000000000000000000019, -0.100000000000000000000000000000000005, +0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000, +0.400000000000000000000000000000000019
138lup_ref - mat_lup, format = SK_'(*(f0.1,:,", "))'
139.0, .0, .0, .0, .0, .0, .0, .0, .0
140.0, .0, .0, -.0, -.0, .0, .0, -.0, -.0
141.0, -.0, .0, .0, .0, .0, .0, .0, .0
142-.0, .0, .0, -.0, -.0, -.0, -.0, -.0, -.0
143-.0, -.0, -.0, -.0, .0, .0, .0, .0, .0
144.0, .0, .0, -.0, .0, .0, -.0, .0, -.0
145.0, -.0, .0, -.0, .0, .0, -.0, -.0, -.0
146-.0, .0, .0, -.0, .0, .0, -.0, .0, .0
147-.0, -.0, .0, -.0, .0, .0, -.0, .0, -.0
148rperm_ref
149+9, +9, +9, +9, +9, +9, +9, +9, +9
150rperm
151+9, +9, +9, +9, +9, +9, +9, +9, +9
152
153
154!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
155! Compute the LU-Pivoted decomposition of a complex square matrix.
156!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
157
158
159mat_lup
160(+2.00000000, +1.00000000), (+2.40000010, -1.00000000), (+2.79999995, -1.00000000), (+3.20000005, -1.00000000), (+3.59999990, -1.00000000), (+4.00000000, -1.00000000), (+4.40000010, -1.00000000), (+4.80000019, -1.00000000), (+5.19999981, -1.00000000)
161(+2.40000010, +1.00000000), (+2.00000000, +1.00000000), (+2.40000010, -1.00000000), (+2.79999995, -1.00000000), (+3.20000005, -1.00000000), (+3.59999990, -1.00000000), (+4.00000000, -1.00000000), (+4.40000010, -1.00000000), (+4.80000019, -1.00000000)
162(+2.79999995, +1.00000000), (+2.40000010, +1.00000000), (+2.00000000, +1.00000000), (+2.40000010, -1.00000000), (+2.79999995, -1.00000000), (+3.20000005, -1.00000000), (+3.59999990, -1.00000000), (+4.00000000, -1.00000000), (+4.40000010, -1.00000000)
163(+3.20000005, +1.00000000), (+2.79999995, +1.00000000), (+2.40000010, +1.00000000), (+2.00000000, +1.00000000), (+2.40000010, -1.00000000), (+2.79999995, -1.00000000), (+3.20000005, -1.00000000), (+3.59999990, -1.00000000), (+4.00000000, -1.00000000)
164(+3.59999990, +1.00000000), (+3.20000005, +1.00000000), (+2.79999995, +1.00000000), (+2.40000010, +1.00000000), (+2.00000000, +1.00000000), (+2.40000010, -1.00000000), (+2.79999995, -1.00000000), (+3.20000005, -1.00000000), (+3.59999990, -1.00000000)
165(+4.00000000, +1.00000000), (+3.59999990, +1.00000000), (+3.20000005, +1.00000000), (+2.79999995, +1.00000000), (+2.40000010, +1.00000000), (+2.00000000, +1.00000000), (+2.40000010, -1.00000000), (+2.79999995, -1.00000000), (+3.20000005, -1.00000000)
166(+4.40000010, +1.00000000), (+4.00000000, +1.00000000), (+3.59999990, +1.00000000), (+3.20000005, +1.00000000), (+2.79999995, +1.00000000), (+2.40000010, +1.00000000), (+2.00000000, +1.00000000), (+2.40000010, -1.00000000), (+2.79999995, -1.00000000)
167(+4.80000019, +1.00000000), (+4.40000010, +1.00000000), (+4.00000000, +1.00000000), (+3.59999990, +1.00000000), (+3.20000005, +1.00000000), (+2.79999995, +1.00000000), (+2.40000010, +1.00000000), (+2.00000000, +1.00000000), (+2.40000010, -1.00000000)
168(+5.19999981, +1.00000000), (+4.80000019, +1.00000000), (+4.40000010, +1.00000000), (+4.00000000, +1.00000000), (+3.59999990, +1.00000000), (+3.20000005, +1.00000000), (+2.79999995, +1.00000000), (+2.40000010, +1.00000000), (+2.00000000, +1.00000000)
169call setResized(rperm, size(mat_lup, 1, IK))
170call setMatLUP(mat_lup, rperm, info); if (info /= 0) error stop
171mat_lup
172(+5.19999981, +1.00000000), (+4.80000019, +1.00000000), (+4.40000010, +1.00000000), (+4.00000000, +1.00000000), (+3.59999990, +1.00000000), (+3.20000005, +1.00000000), (+2.79999995, +1.00000000), (+2.40000010, +1.00000000), (+2.00000000, +1.00000000)
173(+0.406562090, +0.114122696), (+0.562624693, -1.95435107), (+1.12524951, -1.90870190), (+1.68787444, -1.86305285), (+2.25049925, -1.81740379), (+2.81312418, -1.77175474), (+3.37574911, -1.72610569), (+3.93837380, -1.68045664), (+4.50099850, -1.63480747)
174(+0.480741858, +0.998573601E-1), (-0.471276082E-1, -0.927102715E-1), (+0.614579797, -1.90574467), (+1.22915947, -1.81148922), (+1.84373927, -1.71723378), (+2.45831895, -1.62297845), (+3.07289886, -1.52872312), (+3.68747878, -1.43446767), (+4.30205870, -1.34021246)
175(+0.554921627, +0.855920240E-1), (-0.403950885E-1, -0.794659629E-1), (-0.627603084E-1, -0.631566793E-1), (+0.673686862, -1.87447929), (+1.34737337, -1.74895859), (+2.02106023, -1.62343812), (+2.69474673, -1.49791741), (+3.36843348, -1.37239671), (+4.04212046, -1.24687600)
176(+0.629101396, +0.713266879E-1), (-0.336626060E-1, -0.662218928E-1), (-0.523002893E-1, -0.526303090E-1), (-0.672924370E-1, -0.319699794E-1), (+0.736060262, -1.86541510), (+1.47212029, -1.73083019), (+2.20818043, -1.59624541), (+2.94424057, -1.46166050), (+3.68030071, -1.32707548)
177(+0.703281105, +0.570613593E-1), (-0.269300230E-1, -0.529773608E-1), (-0.418402553E-1, -0.421044528E-1), (-0.538339317E-1, -0.255758688E-1), (-0.593044236E-1, -0.402076915E-2), (+0.791958690, -1.88139129), (+1.58391690, -1.76278222), (+2.37587571, -1.64417326), (+3.16783404, -1.52556443)
178(+0.777460873, +0.427960157E-1), (-0.201975685E-1, -0.397332832E-1), (-0.313801952E-1, -0.315780789E-1), (-0.403755009E-1, -0.191820823E-1), (-0.444782786E-1, -0.301529467E-2), (-0.413116254E-1, +0.141840773E-1), (+0.828368783, -1.91737688), (+1.65673697, -1.83475339), (+2.48510551, -1.75213003)
179(+0.851640642, +0.285306722E-1), (-0.134651121E-1, -0.264890790E-1), (-0.209201779E-1, -0.210520662E-1), (-0.269169938E-1, -0.127878133E-1), (-0.296522155E-1, -0.201031473E-2), (-0.275410283E-1, +0.945628528E-2), (-0.206180140E-1, +0.187712964E-1), (+0.837542593, -1.95876396), (+1.67508531, -1.91752791)
180(+0.925820410, +0.142653435E-1), (-0.673253182E-2, -0.132446680E-1), (-0.104601085E-1, -0.105259735E-1), (-0.134585742E-1, -0.639407430E-2), (-0.148261329E-1, -0.100500276E-2), (-0.137705412E-1, +0.472813379E-2), (-0.103089772E-1, +0.938579254E-2), (-0.543472543E-2, +0.119069861E-1), (+0.823814571, -1.98913050)
181call setResized(rperm, size(mat_lup, 1, IK)) ! reconstruct the original matrix.
182lup_ref ! reference matrix rounded to 1 significant digit.
183(+5.19999981, +1.00000000), (+4.80000019, +1.00000000), (+4.40000010, +1.00000000), (+4.00000000, +1.00000000), (+3.59999990, +1.00000000), (+3.20000005, +1.00000000), (+2.79999995, +1.00000000), (+2.40000010, +1.00000000), (+2.00000000, +1.00000000)
184(+0.400000006, +0.100000001), (+0.600000024, -2.00000000), (+1.10000002, -1.89999998), (+1.70000005, -1.89999998), (+2.29999995, -1.79999995), (+2.79999995, -1.79999995), (+3.40000010, -1.70000005), (+3.90000010, -1.70000005), (+4.50000000, -1.60000002)
185(+0.500000000, +0.100000001), (+0.00000000, -0.100000001), (+0.600000024, -1.89999998), (+1.20000005, -1.79999995), (+1.79999995, -1.70000005), (+2.50000000, -1.60000002), (+3.09999990, -1.50000000), (+3.70000005, -1.39999998), (+4.30000019, -1.29999995)
186(+0.600000024, +0.100000001), (+0.00000000, -0.100000001), (-0.100000001, -0.100000001), (+0.699999988, -1.89999998), (+1.29999995, -1.70000005), (+2.00000000, -1.60000002), (+2.70000005, -1.50000000), (+3.40000010, -1.39999998), (+4.00000000, -1.20000005)
187(+0.600000024, +0.100000001), (+0.00000000, -0.100000001), (-0.100000001, -0.100000001), (-0.100000001, +0.00000000), (+0.699999988, -1.89999998), (+1.50000000, -1.70000005), (+2.20000005, -1.60000002), (+2.90000010, -1.50000000), (+3.70000005, -1.29999995)
188(+0.699999988, +0.100000001), (+0.00000000, -0.100000001), (+0.00000000, +0.00000000), (-0.100000001, +0.00000000), (-0.100000001, +0.00000000), (+0.800000012, -1.89999998), (+1.60000002, -1.79999995), (+2.40000010, -1.60000002), (+3.20000005, -1.50000000)
189(+0.800000012, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.800000012, -1.89999998), (+1.70000005, -1.79999995), (+2.50000000, -1.79999995)
190(+0.899999976, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.800000012, -2.00000000), (+1.70000005, -1.89999998)
191(+0.899999976, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.800000012, -2.00000000)
192lup_ref - mat_lup, format = SK_'(*("(",f0.1,",",f0.1,")",:,", "))'
193(.0,.0), (.0,.0), (.0,.0), (.0,.0), (.0,.0), (.0,.0), (.0,.0), (.0,.0), (.0,.0)
194(-.0,-.0), (.0,-.0), (-.0,.0), (.0,-.0), (.0,.0), (-.0,-.0), (.0,.0), (-.0,-.0), (-.0,.0)
195(.0,.0), (.0,-.0), (-.0,.0), (-.0,.0), (-.0,.0), (.0,.0), (.0,.0), (.0,.0), (-.0,.0)
196(.0,.0), (.0,-.0), (-.0,-.0), (.0,-.0), (-.0,.0), (-.0,.0), (.0,-.0), (.0,-.0), (-.0,.0)
197(-.0,.0), (.0,-.0), (-.0,-.0), (-.0,.0), (-.0,-.0), (.0,.0), (-.0,-.0), (-.0,-.0), (.0,.0)
198(-.0,.0), (.0,-.0), (.0,.0), (-.0,.0), (-.0,.0), (.0,-.0), (.0,-.0), (.0,.0), (.0,.0)
199(.0,-.0), (.0,.0), (.0,.0), (.0,.0), (.0,.0), (.0,-.0), (-.0,.0), (.0,.0), (.0,-.0)
200(.0,-.0), (.0,.0), (.0,.0), (.0,.0), (.0,.0), (.0,-.0), (.0,-.0), (-.0,-.0), (.0,.0)
201(-.0,-.0), (.0,.0), (.0,.0), (.0,.0), (.0,.0), (.0,-.0), (.0,-.0), (.0,-.0), (-.0,-.0)
202rperm_ref
203+9, +9, +9, +9, +9, +9, +9, +9, +9
204rperm
205+9, +9, +9, +9, +9, +9, +9, +9, +9
206
207
208mat_lup
209(+2.0000000000000000, +1.0000000000000000), (+2.4000000953674316, -1.0000000000000000), (+2.7999999523162842, -1.0000000000000000), (+3.2000000476837158, -1.0000000000000000), (+3.5999999046325684, -1.0000000000000000), (+4.0000000000000000, -1.0000000000000000), (+4.4000000953674316, -1.0000000000000000), (+4.8000001907348633, -1.0000000000000000), (+5.1999998092651367, -1.0000000000000000)
210(+2.4000000953674316, +1.0000000000000000), (+2.0000000000000000, +1.0000000000000000), (+2.4000000953674316, -1.0000000000000000), (+2.7999999523162842, -1.0000000000000000), (+3.2000000476837158, -1.0000000000000000), (+3.5999999046325684, -1.0000000000000000), (+4.0000000000000000, -1.0000000000000000), (+4.4000000953674316, -1.0000000000000000), (+4.8000001907348633, -1.0000000000000000)
211(+2.7999999523162842, +1.0000000000000000), (+2.4000000953674316, +1.0000000000000000), (+2.0000000000000000, +1.0000000000000000), (+2.4000000953674316, -1.0000000000000000), (+2.7999999523162842, -1.0000000000000000), (+3.2000000476837158, -1.0000000000000000), (+3.5999999046325684, -1.0000000000000000), (+4.0000000000000000, -1.0000000000000000), (+4.4000000953674316, -1.0000000000000000)
212(+3.2000000476837158, +1.0000000000000000), (+2.7999999523162842, +1.0000000000000000), (+2.4000000953674316, +1.0000000000000000), (+2.0000000000000000, +1.0000000000000000), (+2.4000000953674316, -1.0000000000000000), (+2.7999999523162842, -1.0000000000000000), (+3.2000000476837158, -1.0000000000000000), (+3.5999999046325684, -1.0000000000000000), (+4.0000000000000000, -1.0000000000000000)
213(+3.5999999046325684, +1.0000000000000000), (+3.2000000476837158, +1.0000000000000000), (+2.7999999523162842, +1.0000000000000000), (+2.4000000953674316, +1.0000000000000000), (+2.0000000000000000, +1.0000000000000000), (+2.4000000953674316, -1.0000000000000000), (+2.7999999523162842, -1.0000000000000000), (+3.2000000476837158, -1.0000000000000000), (+3.5999999046325684, -1.0000000000000000)
214(+4.0000000000000000, +1.0000000000000000), (+3.5999999046325684, +1.0000000000000000), (+3.2000000476837158, +1.0000000000000000), (+2.7999999523162842, +1.0000000000000000), (+2.4000000953674316, +1.0000000000000000), (+2.0000000000000000, +1.0000000000000000), (+2.4000000953674316, -1.0000000000000000), (+2.7999999523162842, -1.0000000000000000), (+3.2000000476837158, -1.0000000000000000)
215(+4.4000000953674316, +1.0000000000000000), (+4.0000000000000000, +1.0000000000000000), (+3.5999999046325684, +1.0000000000000000), (+3.2000000476837158, +1.0000000000000000), (+2.7999999523162842, +1.0000000000000000), (+2.4000000953674316, +1.0000000000000000), (+2.0000000000000000, +1.0000000000000000), (+2.4000000953674316, -1.0000000000000000), (+2.7999999523162842, -1.0000000000000000)
216(+4.8000001907348633, +1.0000000000000000), (+4.4000000953674316, +1.0000000000000000), (+4.0000000000000000, +1.0000000000000000), (+3.5999999046325684, +1.0000000000000000), (+3.2000000476837158, +1.0000000000000000), (+2.7999999523162842, +1.0000000000000000), (+2.4000000953674316, +1.0000000000000000), (+2.0000000000000000, +1.0000000000000000), (+2.4000000953674316, -1.0000000000000000)
217(+5.1999998092651367, +1.0000000000000000), (+4.8000001907348633, +1.0000000000000000), (+4.4000000953674316, +1.0000000000000000), (+4.0000000000000000, +1.0000000000000000), (+3.5999999046325684, +1.0000000000000000), (+3.2000000476837158, +1.0000000000000000), (+2.7999999523162842, +1.0000000000000000), (+2.4000000953674316, +1.0000000000000000), (+2.0000000000000000, +1.0000000000000000)
218call setResized(rperm, size(mat_lup, 1, IK))
219call setMatLUP(mat_lup, rperm, info); if (info /= 0) error stop
220mat_lup
221(+4.8000001907348633, +1.0000000000000000), (+4.4000000953674316, +1.0000000000000000), (+4.0000000000000000, +1.0000000000000000), (+3.5999999046325684, +1.0000000000000000), (+3.2000000476837158, +1.0000000000000000), (+2.7999999523162842, +1.0000000000000000), (+2.4000000953674316, +1.0000000000000000), (+2.0000000000000000, +1.0000000000000000), (+2.4000000953674316, -1.0000000000000000)
222(+0.44093176264973133, +0.11647254481976956), (+0.57637284247785381, -1.9534109709644047), (+1.1527454465371285, -1.9068219419288095), (+1.7291182890149823, -1.8602329128932142), (+2.3054907879479329, -1.8136439116268377), (+2.8818636304257867, -1.7670548825912424), (+3.4582363677773165, -1.7204658813248657), (+4.0346092102551703, -1.6738768522892704), (+4.0252909920354822, -0.83860235602540301)
223(+0.52079866730101931, +0.99833607011923217E-1), (-0.45439973893843404E-1, -0.84718584222333845E-1), (+0.63056301145646165, -1.9091200724107975), (+1.2611257953280783, -1.8182401246231106), (+1.8916886801920934, -1.7273602310040808), (+2.5222514640637099, -1.6364802832163938), (+3.1528143597614591, -1.5456003693988793), (+3.7833773820516550, -1.4547204216111926), (+3.7042040544544634, -0.41589111428542119)
224(+0.60066552434793374, +0.83194679121654291E-1), (-0.37866614792331624E-1, -0.70598713518433442E-1), (-0.57792589919254116E-1, -0.54871250089531212E-1), (+0.69025768532287168, -1.8844148146499891), (+1.3805149867778459, -1.7688296720294576), (+2.0707726721007176, -1.6532444866794469), (+2.7610299821015896, -1.5376593254650870), (+3.4512876812032882, -1.4220741270327508), (+3.3237325026764735, -0.16735749784694393)
225(+0.68053242899922173, +0.66555741313807953E-1), (-0.30293317433914217E-1, -0.56479058953567654E-1), (-0.46234056670262633E-1, -0.43896872656820095E-1), (-0.58311665170033999E-1, -0.25230067457196938E-1), (+0.74953987047137327, -1.8833767160201977), (+1.4990796878732637, -1.7667533856801301), (+2.2486195652105661, -1.6501300855253955), (+2.9981593941089129, -1.5335067436090046), (+2.8570207711033255, -0.59787015427843371E-1)
226(+0.76039928604613616, +0.49916813423539014E-1), (-0.22719958332402375E-1, -0.42359188249667050E-1), (-0.34675546433463329E-1, -0.32922760731555381E-1), (-0.43733717463510409E-1, -0.18922453442919950E-1), (-0.46995502238813366E-1, -0.13912995938357257E-2), (+0.79721778760510109, -1.9060089586473699), (+1.5944351461833164, -1.8120179412688675), (+2.3916529429349374, -1.7180268906739415), (+2.2771155745492440, -0.43675391451432448E-1)
227(+0.84026619069742414, +0.33277875615692676E-1), (-0.15146660973984916E-1, -0.28239533684801268E-1), (-0.23117013184471846E-1, -0.21948383298844045E-1), (-0.29155820661109833E-1, -0.12615045074678387E-1), (-0.31330312873164096E-1, -0.92745383702202630E-3), (-0.28326727191086122E-1, +0.10874888479857528E-1), (+0.82174974873585849, -1.9433465576141231), (+1.6434994844691331, -1.8866930969081133), (+1.5821054054156254, -0.55117628645948741E-1)
228(+0.92013309534871213, +0.16638937807846338E-1), (-0.75733304869925241E-2, -0.14119766842400957E-1), (-0.11558525187603533E-1, -0.10974247949567641E-1), (-0.14577885737126345E-1, -0.63074015173988163E-2), (-0.15665176722368543E-1, -0.46384351348661285E-3), (-0.14163336206010608E-1, +0.54375606149159785E-2), (-0.10357962156502729E-1, +0.99757332762874310E-2), (+0.81995188925278273, -1.9792840662075553), (+0.80689012042934261, -0.43912511384084786E-1)
229(+1.0798668094425410, -0.16638917972691514E-1), (+0.75733877160311762E-2, +0.14119974565400699E-1), (+0.11558569667465674E-1, +0.10974201101544966E-1), (+0.14577952164102010E-1, +0.63074539850387140E-2), (+0.15665126284363676E-1, +0.46353953016864055E-3), (+0.14163355729338623E-1, -0.54373874031065179E-2), (+0.10357867111659553E-1, -0.99760055630593973E-2), (+0.53689886243285193E-2, -0.12304368340802769E-1), (-0.81068252211711833, +2.0540765456545467)
230call setResized(rperm, size(mat_lup, 1, IK)) ! reconstruct the original matrix.
231lup_ref ! reference matrix rounded to 1 significant digit.
232(+5.1999998092651367, +1.0000000000000000), (+4.8000001907348633, +1.0000000000000000), (+4.4000000953674316, +1.0000000000000000), (+4.0000000000000000, +1.0000000000000000), (+3.5999999046325684, +1.0000000000000000), (+3.2000000476837158, +1.0000000000000000), (+2.7999999523162842, +1.0000000000000000), (+2.4000000953674316, +1.0000000000000000), (+2.0000000000000000, +1.0000000000000000)
233(+0.40000000596046448, +0.10000000149011612), (+0.60000002384185791, -2.0000000000000000), (+1.1000000238418579, -1.8999999761581421), (+1.7000000476837158, -1.8999999761581421), (+2.2999999523162842, -1.7999999523162842), (+2.7999999523162842, -1.7999999523162842), (+3.4000000953674316, -1.7000000476837158), (+3.9000000953674316, -1.7000000476837158), (+4.5000000000000000, -1.6000000238418579)
234(+0.50000000000000000, +0.10000000149011612), (+0.0000000000000000, -0.10000000149011612), (+0.60000002384185791, -1.8999999761581421), (+1.2000000476837158, -1.7999999523162842), (+1.7999999523162842, -1.7000000476837158), (+2.5000000000000000, -1.6000000238418579), (+3.0999999046325684, -1.5000000000000000), (+3.7000000476837158, -1.3999999761581421), (+4.3000001907348633, -1.2999999523162842)
235(+0.60000002384185791, +0.10000000149011612), (+0.0000000000000000, -0.10000000149011612), (-0.10000000149011612, -0.10000000149011612), (+0.69999998807907104, -1.8999999761581421), (+1.2999999523162842, -1.7000000476837158), (+2.0000000000000000, -1.6000000238418579), (+2.7000000476837158, -1.5000000000000000), (+3.4000000953674316, -1.3999999761581421), (+4.0000000000000000, -1.2000000476837158)
236(+0.60000002384185791, +0.10000000149011612), (+0.0000000000000000, -0.10000000149011612), (-0.10000000149011612, -0.10000000149011612), (-0.10000000149011612, +0.0000000000000000), (+0.69999998807907104, -1.8999999761581421), (+1.5000000000000000, -1.7000000476837158), (+2.2000000476837158, -1.6000000238418579), (+2.9000000953674316, -1.5000000000000000), (+3.7000000476837158, -1.2999999523162842)
237(+0.69999998807907104, +0.10000000149011612), (+0.0000000000000000, -0.10000000149011612), (+0.0000000000000000, +0.0000000000000000), (-0.10000000149011612, +0.0000000000000000), (-0.10000000149011612, +0.0000000000000000), (+0.80000001192092896, -1.8999999761581421), (+1.6000000238418579, -1.7999999523162842), (+2.4000000953674316, -1.6000000238418579), (+3.2000000476837158, -1.5000000000000000)
238(+0.80000001192092896, +0.0000000000000000), (+0.0000000000000000, +0.0000000000000000), (+0.0000000000000000, +0.0000000000000000), (+0.0000000000000000, +0.0000000000000000), (+0.0000000000000000, +0.0000000000000000), (+0.0000000000000000, +0.0000000000000000), (+0.80000001192092896, -1.8999999761581421), (+1.7000000476837158, -1.7999999523162842), (+2.5000000000000000, -1.7999999523162842)
239(+0.89999997615814209, +0.0000000000000000), (+0.0000000000000000, +0.0000000000000000), (+0.0000000000000000, +0.0000000000000000), (+0.0000000000000000, +0.0000000000000000), (+0.0000000000000000, +0.0000000000000000), (+0.0000000000000000, +0.0000000000000000), (+0.0000000000000000, +0.0000000000000000), (+0.80000001192092896, -2.0000000000000000), (+1.7000000476837158, -1.8999999761581421)
240(+0.89999997615814209, +0.0000000000000000), (+0.0000000000000000, +0.0000000000000000), (+0.0000000000000000, +0.0000000000000000), (+0.0000000000000000, +0.0000000000000000), (+0.0000000000000000, +0.0000000000000000), (+0.0000000000000000, +0.0000000000000000), (+0.0000000000000000, +0.0000000000000000), (+0.0000000000000000, +0.0000000000000000), (+0.80000001192092896, -2.0000000000000000)
241lup_ref - mat_lup, format = SK_'(*("(",f0.1,",",f0.1,")",:,", "))'
242(.4,.0), (.4,.0), (.4,.0), (.4,.0), (.4,.0), (.4,.0), (.4,.0), (.4,.0), (-.4,2.0)
243(-.0,-.0), (.0,-.0), (-.1,.0), (-.0,-.0), (-.0,.0), (-.1,-.0), (-.1,.0), (-.1,-.0), (.5,-.8)
244(-.0,.0), (.0,-.0), (-.0,.0), (-.1,.0), (-.1,.0), (-.0,.0), (-.1,.0), (-.1,.1), (.6,-.9)
245(-.0,.0), (.0,-.0), (-.0,-.0), (.0,-.0), (-.1,.1), (-.1,.1), (-.1,.0), (-.1,.0), (.7,-1.0)
246(-.1,.0), (.0,-.0), (-.1,-.1), (-.0,.0), (-.0,-.0), (.0,.1), (-.0,.1), (-.1,.0), (.8,-1.2)
247(-.1,.1), (.0,-.1), (.0,.0), (-.1,.0), (-.1,.0), (.0,.0), (.0,.0), (.0,.1), (.9,-1.5)
248(-.0,-.0), (.0,.0), (.0,.0), (.0,.0), (.0,.0), (.0,-.0), (-.0,.0), (.1,.1), (.9,-1.7)
249(-.0,-.0), (.0,.0), (.0,.0), (.0,.0), (.0,.0), (.0,-.0), (.0,-.0), (-.0,-.0), (.9,-1.9)
250(-.2,.0), (-.0,-.0), (-.0,-.0), (-.0,-.0), (-.0,-.0), (-.0,.0), (-.0,.0), (-.0,.0), (1.6,-4.1)
251rperm_ref
252+9, +9, +9, +9, +9, +9, +9, +9, +9
253rperm
254+8, +8, +8, +8, +8, +8, +8, +8, +9
255
256
257mat_lup
258(+2.00000000000000000000000000000000000, +1.00000000000000000000000000000000000), (+2.40000009536743164062500000000000000, -1.00000000000000000000000000000000000), (+2.79999995231628417968750000000000000, -1.00000000000000000000000000000000000), (+3.20000004768371582031250000000000000, -1.00000000000000000000000000000000000), (+3.59999990463256835937500000000000000, -1.00000000000000000000000000000000000), (+4.00000000000000000000000000000000000, -1.00000000000000000000000000000000000), (+4.40000009536743164062500000000000000, -1.00000000000000000000000000000000000), (+4.80000019073486328125000000000000000, -1.00000000000000000000000000000000000), (+5.19999980926513671875000000000000000, -1.00000000000000000000000000000000000)
259(+2.40000009536743164062500000000000000, +1.00000000000000000000000000000000000), (+2.00000000000000000000000000000000000, +1.00000000000000000000000000000000000), (+2.40000009536743164062500000000000000, -1.00000000000000000000000000000000000), (+2.79999995231628417968750000000000000, -1.00000000000000000000000000000000000), (+3.20000004768371582031250000000000000, -1.00000000000000000000000000000000000), (+3.59999990463256835937500000000000000, -1.00000000000000000000000000000000000), (+4.00000000000000000000000000000000000, -1.00000000000000000000000000000000000), (+4.40000009536743164062500000000000000, -1.00000000000000000000000000000000000), (+4.80000019073486328125000000000000000, -1.00000000000000000000000000000000000)
260(+2.79999995231628417968750000000000000, +1.00000000000000000000000000000000000), (+2.40000009536743164062500000000000000, +1.00000000000000000000000000000000000), (+2.00000000000000000000000000000000000, +1.00000000000000000000000000000000000), (+2.40000009536743164062500000000000000, -1.00000000000000000000000000000000000), (+2.79999995231628417968750000000000000, -1.00000000000000000000000000000000000), (+3.20000004768371582031250000000000000, -1.00000000000000000000000000000000000), (+3.59999990463256835937500000000000000, -1.00000000000000000000000000000000000), (+4.00000000000000000000000000000000000, -1.00000000000000000000000000000000000), (+4.40000009536743164062500000000000000, -1.00000000000000000000000000000000000)
261(+3.20000004768371582031250000000000000, +1.00000000000000000000000000000000000), (+2.79999995231628417968750000000000000, +1.00000000000000000000000000000000000), (+2.40000009536743164062500000000000000, +1.00000000000000000000000000000000000), (+2.00000000000000000000000000000000000, +1.00000000000000000000000000000000000), (+2.40000009536743164062500000000000000, -1.00000000000000000000000000000000000), (+2.79999995231628417968750000000000000, -1.00000000000000000000000000000000000), (+3.20000004768371582031250000000000000, -1.00000000000000000000000000000000000), (+3.59999990463256835937500000000000000, -1.00000000000000000000000000000000000), (+4.00000000000000000000000000000000000, -1.00000000000000000000000000000000000)
262(+3.59999990463256835937500000000000000, +1.00000000000000000000000000000000000), (+3.20000004768371582031250000000000000, +1.00000000000000000000000000000000000), (+2.79999995231628417968750000000000000, +1.00000000000000000000000000000000000), (+2.40000009536743164062500000000000000, +1.00000000000000000000000000000000000), (+2.00000000000000000000000000000000000, +1.00000000000000000000000000000000000), (+2.40000009536743164062500000000000000, -1.00000000000000000000000000000000000), (+2.79999995231628417968750000000000000, -1.00000000000000000000000000000000000), (+3.20000004768371582031250000000000000, -1.00000000000000000000000000000000000), (+3.59999990463256835937500000000000000, -1.00000000000000000000000000000000000)
263(+4.00000000000000000000000000000000000, +1.00000000000000000000000000000000000), (+3.59999990463256835937500000000000000, +1.00000000000000000000000000000000000), (+3.20000004768371582031250000000000000, +1.00000000000000000000000000000000000), (+2.79999995231628417968750000000000000, +1.00000000000000000000000000000000000), (+2.40000009536743164062500000000000000, +1.00000000000000000000000000000000000), (+2.00000000000000000000000000000000000, +1.00000000000000000000000000000000000), (+2.40000009536743164062500000000000000, -1.00000000000000000000000000000000000), (+2.79999995231628417968750000000000000, -1.00000000000000000000000000000000000), (+3.20000004768371582031250000000000000, -1.00000000000000000000000000000000000)
264(+4.40000009536743164062500000000000000, +1.00000000000000000000000000000000000), (+4.00000000000000000000000000000000000, +1.00000000000000000000000000000000000), (+3.59999990463256835937500000000000000, +1.00000000000000000000000000000000000), (+3.20000004768371582031250000000000000, +1.00000000000000000000000000000000000), (+2.79999995231628417968750000000000000, +1.00000000000000000000000000000000000), (+2.40000009536743164062500000000000000, +1.00000000000000000000000000000000000), (+2.00000000000000000000000000000000000, +1.00000000000000000000000000000000000), (+2.40000009536743164062500000000000000, -1.00000000000000000000000000000000000), (+2.79999995231628417968750000000000000, -1.00000000000000000000000000000000000)
265(+4.80000019073486328125000000000000000, +1.00000000000000000000000000000000000), (+4.40000009536743164062500000000000000, +1.00000000000000000000000000000000000), (+4.00000000000000000000000000000000000, +1.00000000000000000000000000000000000), (+3.59999990463256835937500000000000000, +1.00000000000000000000000000000000000), (+3.20000004768371582031250000000000000, +1.00000000000000000000000000000000000), (+2.79999995231628417968750000000000000, +1.00000000000000000000000000000000000), (+2.40000009536743164062500000000000000, +1.00000000000000000000000000000000000), (+2.00000000000000000000000000000000000, +1.00000000000000000000000000000000000), (+2.40000009536743164062500000000000000, -1.00000000000000000000000000000000000)
266(+5.19999980926513671875000000000000000, +1.00000000000000000000000000000000000), (+4.80000019073486328125000000000000000, +1.00000000000000000000000000000000000), (+4.40000009536743164062500000000000000, +1.00000000000000000000000000000000000), (+4.00000000000000000000000000000000000, +1.00000000000000000000000000000000000), (+3.59999990463256835937500000000000000, +1.00000000000000000000000000000000000), (+3.20000004768371582031250000000000000, +1.00000000000000000000000000000000000), (+2.79999995231628417968750000000000000, +1.00000000000000000000000000000000000), (+2.40000009536743164062500000000000000, +1.00000000000000000000000000000000000), (+2.00000000000000000000000000000000000, +1.00000000000000000000000000000000000)
267call setResized(rperm, size(mat_lup, 1, IK))
268call setMatLUP(mat_lup, rperm, info); if (info /= 0) error stop
269mat_lup
270(+5.19999980926513671875000000000000000, +1.00000000000000000000000000000000000), (+4.80000019073486328125000000000000000, +1.00000000000000000000000000000000000), (+4.40000009536743164062500000000000000, +1.00000000000000000000000000000000000), (+4.00000000000000000000000000000000000, +1.00000000000000000000000000000000000), (+3.59999990463256835937500000000000000, +1.00000000000000000000000000000000000), (+3.20000004768371582031250000000000000, +1.00000000000000000000000000000000000), (+2.79999995231628417968750000000000000, +1.00000000000000000000000000000000000), (+2.40000009536743164062500000000000000, +1.00000000000000000000000000000000000), (+2.00000000000000000000000000000000000, +1.00000000000000000000000000000000000)
271(+0.406562069365338448307885910489287586, +0.114122683154199215461098838905655892), (+0.562624768022445588468338572918435761, -1.95435097027266905122676021003513422), (+1.12524949149021386472394688493347071, -1.90870188612740218068952073807352266), (+1.68787445337656124254205519694850547, -1.86305280198213531015228126611191111), (+2.25049917684432951879766350896354003, -1.81740371783686843961504179415029975), (+2.81312404179872600194088695143027587, -1.77175466090056952995980216318706098), (+3.37574900368507337975899526344531043, -1.72610557675530265942256269122544962), (+3.93837386863946986290221870591204588, -1.68045651981900374976732306026221085), (+4.50099835368865903759532701792708083, -1.63480743567373687923008358830059949)
272(+0.480741832801958117292577129687623789, +0.998573435085228162373729652872935195E-1), (-0.471276558567860914748836487117786875E-1, -0.927102459830737368039370104475017215E-1), (+0.614579920859427465911694723762965187, -1.90574464223682124948712676866829751), (+1.22915961453638457602325943188654205, -1.81148926236979738353347079866600506), (+1.84373953539581204193495415564950724, -1.71723390460661863302059756733430237), (+2.45831911240935895003054949926379406, -1.62297855881632060856888593362137741), (+3.07289904450489516170461420738737035, -1.52872317894929674261522996361908477), (+3.68747885993702117136270955100165698, -1.43446783315899871816351832990616001), (+4.30205876956033989151203429040401106, -1.34021249749966508309142783724504746)
273(+0.554921552024003836230914432866965713, +0.855920123656494115958241291792419283E-1), (-0.403951121455697039120147223467780034E-1, -0.794658480422347130718860658065729188E-1), (-0.627602857301207443682885519813026064E-1, -0.631566590377966007600423280489925702E-1), (+0.673686866636173558340052587708276332, -1.87447942255446611837934577005883122), (+1.34737350087686453650746793161029105, -1.74895885109422462964276866609622610), (+2.02106022828638003150437123704176893, -1.62343831236404198193722259312836145), (+2.69474687146719263606538071492504245, -1.49791771986078719172324733542136386), (+3.36843361383992627884574091053151772, -1.37239716607288363542438023468767022), (+4.04212047755907473218206212045227362, -1.24687661565565544816520128446945738)
274(+0.629101315460623505215605652065301869, +0.713266727199730123720982555608795434E-1), (-0.336626292856229384379880377705494937E-1, -0.662216675213578649899811814864377975E-1), (-0.523002320715055155823898188827763888E-1, -0.526304138751594806269631989698370126E-1), (-0.672924019125919549643092114853026028E-1, -0.319699768768599811951459288313160240E-1), (+0.736060226553276957549191310734687226, -1.86541520899506366109404033852971182), (+1.47212006122601974191805951002946982, -1.73083044448899078627345177298245764), (+2.20818029515826712985713531376877595, -1.59624563424633089813616582696439779), (+2.94424014282587283163681185313571527, -1.46166085570054558421614426330018272), (+3.68030037804393820712755934953600142, -1.32707607631109633365074541130816035)
275(+0.703281034682669224153942955244643745, +0.570613415770996077305494194528279522E-1), (-0.269300855744065508751191114055488037E-1, -0.529772695805188412579302368455089827E-1), (-0.418402000938279183921447492062127882E-1, -0.421044339682751292558655355090701566E-1), (-0.538338972351025411696959586957909707E-1, -0.255758710441887004966765345983583686E-1), (-0.593043873886219968696752309963618543E-1, -0.402076063076755818085941850815341357E-2), (+0.791958475886349191700948829309255587, -1.88139124956905794717833238870455920), (+1.58391691320452537763642197760745175, -1.76278245823110001708634069034526111), (+2.37587540020097170640325504301979621, -1.64417369568349972739735817982365975), (+3.16783383874811303504275522230175902, -1.52556493085266665857913858026615535)
276(+0.777460798119288893138634174442979996, +0.427960019314232085068235458344655793E-1), (-0.201976027144597854010924268293202458E-1, -0.397330890596419931760253525253738795E-1), (-0.313801464352126896062460161076865405E-1, -0.315781888056380091227864064299146290E-1), (-0.403754439242383030227950331983194926E-1, -0.191819724998771277878334675109860831E-1), (-0.444782589931566153015496414658317552E-1, -0.301543346466736808013204424672655076E-2), (-0.413116146453040020299505299092477523E-1, +0.141840636466419621856831837264744262E-1), (+0.828368526507507596810006089103654252, -1.91737673779223798765799608095070898), (+1.65673661508048797573645073478451095, -1.83475349907119737124140227686285987), (+2.48510514218449015954916482133202790, -1.75213025785630970859552575892376668)
277(+0.851640561555908562123325393641316248, +0.285306622857468092830976722161032064E-1), (-0.134650874224327751257378527300013260E-1, -0.264887958816780474641018892818576325E-1), (-0.209201370441056221305356882034695206E-1, -0.210521864539334850345359668103215481E-1), (-0.269169529667857587460214690805018382E-1, -0.127878377884406372107638608800053681E-1), (-0.296521837717206704407329334349308029E-1, -0.201034335075040916853140300125287678E-2), (-0.275410467123129612507834699188405339E-1, +0.945614108546513173148461045746392651E-2), (-0.206180555192709098685721777582319888E-1, +0.187712773249124954119394791755320377E-1), (+0.837542523430939427095356246190945171, -1.95876389688298615286913953596245227), (+1.67508503041293076305088668543487630, -1.91752780028379045202604934309322035)
278(+0.925820324992528231108016612839652404, +0.142653226400704100593717985977408335E-1), (-0.673257213040576485038327863068239569E-2, -0.132445027037141017521784260383411944E-1), (-0.104600911086233407798025858740813334E-1, -0.105260707821465698001152072463345747E-1), (-0.134585104451901024018939279180789983E-1, -0.639394633652226621726645699311453279E-2), (-0.148260679817866902259470784154405160E-1, -0.100501881481330360696970768247508484E-2), (-0.137705342927948242166573075490478702E-1, +0.472798419457221890381492461113495016E-2), (-0.103089796220610467607581439304703944E-1, +0.938577648888811482770145145899397056E-2), (-0.543474132826119831719717505728610759E-2, +0.119069538582110081459178407213557556E-1), (+0.823814328478359595612029689736234028, -1.98913052044211662862968091128137802)
279call setResized(rperm, size(mat_lup, 1, IK)) ! reconstruct the original matrix.
280lup_ref ! reference matrix rounded to 1 significant digit.
281(+5.19999980926513671875000000000000000, +1.00000000000000000000000000000000000), (+4.80000019073486328125000000000000000, +1.00000000000000000000000000000000000), (+4.40000009536743164062500000000000000, +1.00000000000000000000000000000000000), (+4.00000000000000000000000000000000000, +1.00000000000000000000000000000000000), (+3.59999990463256835937500000000000000, +1.00000000000000000000000000000000000), (+3.20000004768371582031250000000000000, +1.00000000000000000000000000000000000), (+2.79999995231628417968750000000000000, +1.00000000000000000000000000000000000), (+2.40000009536743164062500000000000000, +1.00000000000000000000000000000000000), (+2.00000000000000000000000000000000000, +1.00000000000000000000000000000000000)
282(+0.400000005960464477539062500000000000, +0.100000001490116119384765625000000000), (+0.600000023841857910156250000000000000, -2.00000000000000000000000000000000000), (+1.10000002384185791015625000000000000, -1.89999997615814208984375000000000000), (+1.70000004768371582031250000000000000, -1.89999997615814208984375000000000000), (+2.29999995231628417968750000000000000, -1.79999995231628417968750000000000000), (+2.79999995231628417968750000000000000, -1.79999995231628417968750000000000000), (+3.40000009536743164062500000000000000, -1.70000004768371582031250000000000000), (+3.90000009536743164062500000000000000, -1.70000004768371582031250000000000000), (+4.50000000000000000000000000000000000, -1.60000002384185791015625000000000000)
283(+0.500000000000000000000000000000000000, +0.100000001490116119384765625000000000), (+0.00000000000000000000000000000000000, -0.100000001490116119384765625000000000), (+0.600000023841857910156250000000000000, -1.89999997615814208984375000000000000), (+1.20000004768371582031250000000000000, -1.79999995231628417968750000000000000), (+1.79999995231628417968750000000000000, -1.70000004768371582031250000000000000), (+2.50000000000000000000000000000000000, -1.60000002384185791015625000000000000), (+3.09999990463256835937500000000000000, -1.50000000000000000000000000000000000), (+3.70000004768371582031250000000000000, -1.39999997615814208984375000000000000), (+4.30000019073486328125000000000000000, -1.29999995231628417968750000000000000)
284(+0.600000023841857910156250000000000000, +0.100000001490116119384765625000000000), (+0.00000000000000000000000000000000000, -0.100000001490116119384765625000000000), (-0.100000001490116119384765625000000000, -0.100000001490116119384765625000000000), (+0.699999988079071044921875000000000000, -1.89999997615814208984375000000000000), (+1.29999995231628417968750000000000000, -1.70000004768371582031250000000000000), (+2.00000000000000000000000000000000000, -1.60000002384185791015625000000000000), (+2.70000004768371582031250000000000000, -1.50000000000000000000000000000000000), (+3.40000009536743164062500000000000000, -1.39999997615814208984375000000000000), (+4.00000000000000000000000000000000000, -1.20000004768371582031250000000000000)
285(+0.600000023841857910156250000000000000, +0.100000001490116119384765625000000000), (+0.00000000000000000000000000000000000, -0.100000001490116119384765625000000000), (-0.100000001490116119384765625000000000, -0.100000001490116119384765625000000000), (-0.100000001490116119384765625000000000, +0.00000000000000000000000000000000000), (+0.699999988079071044921875000000000000, -1.89999997615814208984375000000000000), (+1.50000000000000000000000000000000000, -1.70000004768371582031250000000000000), (+2.20000004768371582031250000000000000, -1.60000002384185791015625000000000000), (+2.90000009536743164062500000000000000, -1.50000000000000000000000000000000000), (+3.70000004768371582031250000000000000, -1.29999995231628417968750000000000000)
286(+0.699999988079071044921875000000000000, +0.100000001490116119384765625000000000), (+0.00000000000000000000000000000000000, -0.100000001490116119384765625000000000), (+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000), (-0.100000001490116119384765625000000000, +0.00000000000000000000000000000000000), (-0.100000001490116119384765625000000000, +0.00000000000000000000000000000000000), (+0.800000011920928955078125000000000000, -1.89999997615814208984375000000000000), (+1.60000002384185791015625000000000000, -1.79999995231628417968750000000000000), (+2.40000009536743164062500000000000000, -1.60000002384185791015625000000000000), (+3.20000004768371582031250000000000000, -1.50000000000000000000000000000000000)
287(+0.800000011920928955078125000000000000, +0.00000000000000000000000000000000000), (+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000), (+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000), (+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000), (+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000), (+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000), (+0.800000011920928955078125000000000000, -1.89999997615814208984375000000000000), (+1.70000004768371582031250000000000000, -1.79999995231628417968750000000000000), (+2.50000000000000000000000000000000000, -1.79999995231628417968750000000000000)
288(+0.899999976158142089843750000000000000, +0.00000000000000000000000000000000000), (+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000), (+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000), (+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000), (+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000), (+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000), (+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000), (+0.800000011920928955078125000000000000, -2.00000000000000000000000000000000000), (+1.70000004768371582031250000000000000, -1.89999997615814208984375000000000000)
289(+0.899999976158142089843750000000000000, +0.00000000000000000000000000000000000), (+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000), (+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000), (+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000), (+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000), (+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000), (+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000), (+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000), (+0.800000011920928955078125000000000000, -2.00000000000000000000000000000000000)
290lup_ref - mat_lup, format = SK_'(*("(",f0.1,",",f0.1,")",:,", "))'
291(.0,.0), (.0,.0), (.0,.0), (.0,.0), (.0,.0), (.0,.0), (.0,.0), (.0,.0), (.0,.0)
292(-.0,-.0), (.0,-.0), (-.0,.0), (.0,-.0), (.0,.0), (-.0,-.0), (.0,.0), (-.0,-.0), (-.0,.0)
293(.0,.0), (.0,-.0), (-.0,.0), (-.0,.0), (-.0,.0), (.0,.0), (.0,.0), (.0,.0), (-.0,.0)
294(.0,.0), (.0,-.0), (-.0,-.0), (.0,-.0), (-.0,.0), (-.0,.0), (.0,-.0), (.0,-.0), (-.0,.0)
295(-.0,.0), (.0,-.0), (-.0,-.0), (-.0,.0), (-.0,-.0), (.0,.0), (-.0,-.0), (-.0,-.0), (.0,.0)
296(-.0,.0), (.0,-.0), (.0,.0), (-.0,.0), (-.0,.0), (.0,-.0), (.0,-.0), (.0,.0), (.0,.0)
297(.0,-.0), (.0,.0), (.0,.0), (.0,.0), (.0,.0), (.0,-.0), (-.0,.0), (.0,.0), (.0,-.0)
298(.0,-.0), (.0,.0), (.0,.0), (.0,.0), (.0,.0), (.0,-.0), (.0,-.0), (-.0,-.0), (.0,.0)
299(-.0,-.0), (.0,.0), (.0,.0), (.0,.0), (.0,.0), (.0,-.0), (.0,-.0), (.0,-.0), (-.0,-.0)
300rperm_ref
301+9, +9, +9, +9, +9, +9, +9, +9, +9
302rperm
303+9, +9, +9, +9, +9, +9, +9, +9, +9
304
305
Test:
test_pm_matrixLUP


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:43 PM, Institute for Computational Engineering and Sciences (ICES), The University of Texas Austin

Definition at line 191 of file pm_matrixLUP.F90.


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