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

Generate and return a matrix of shape (shape(1), shape(2)) with the upper/lower triangle and the diagonal elements of the matrix set to the corresponding requested input values.
More...

Detailed Description

Generate and return a matrix of shape (shape(1), shape(2)) with the upper/lower triangle and the diagonal elements of the matrix set to the corresponding requested input values.

See the documentation of pm_matrixInit for illustrations of possible initialization formats.

Parameters
[in]shape: The input contiguous vector of size 2 of type integer of default kind IK, containing the shape of the output matrix (as returned by the Fortran intrinsic function shape()).
By definition, both elements of shape must be non-negative.
[in]subset: The input argument that can be either,
  1. the constant upp, signifying the initialization of the upper-triangle elements of the requested subset of the matrix.
  2. the constant low, signifying the initialization of the lower-triangle elements of the requested subset of the matrix.
  3. the constant dia, signifying the initialization of the diagonal elements of the requested subset of the matrix.
  4. the constant lowDia, signifying the initialization of the lower-triangle and diagonal elements of the requested subset of the matrix.
  5. the constant uppDia, signifying the initialization of the upper-triangle and diagonal elements of the requested subset of the matrix.
  6. the constant uppLow, signifying the initialization of the upper/lower-triangle elements of the requested subset of the matrix.
  7. the constant uppLowDia, signifying the initialization of the upper/lower-triangle and diagonal elements of the requested subset of the matrix.
This input argument is merely serves to resolve the different procedures of this generic interface from each other at compile-time.
Note that setting only the lower/upper triangles of a matrix can be readily converted to the problem of setting the upper/lower triangle and diagonal elements of a matrix.
[in]vupp: The input scalar of the same type and kind as mat,
containing the value for all upper-triangular elements of the output matrix.
(optional. It can be present only if the requested input subset format is compatible.)
[in]vlow: The input scalar of the same type and kind as mat,
containing the value for all lower-triangular elements of the output matrix.
(optional. It can be present only if the requested input subset format is compatible.)
[in]vdia: The input scalar or contiguous array of rank 1 of the same type and kind as mat containing the diagonal elements of the output matrix.
(optional. It can be present only if the requested input subset format is compatible.)
[in]nrow: The input non-negative integer of default kind IK representing the number of rows of the initialized block of the matrix.
Setting nrow = size(mat,1), roff = 0 is equivalent to initializing all rows of the matrix.
(optional, default = size(mat,1))
[in]ncol: The input non-negative integer of default kind IK representing the number of columns of the initialized block of the matrix.
Setting ncol = size(mat,2), roff = 0 is equivalent to initializing all columns of the matrix.
(optional, default = size(mat,2))
[in]ndia: The input non-negative integer of default kind IK representing the number of diagonal elements to write to the initialized block of the matrix.
It can be present only if the input argument subset is set to dia and the input argument vdia is a scalar.
If vdia is a vector of rank 1, then ndia is automatically set to the size of the input vector vdia.
(optional. default = min(shape(1) - roff, shape(2) - coff).)
[in]roff: The input non-negative integer of default kind IK representing the offset of the top-left corner of the initialization block from the first row of the matrix.
The initialization of the matrix will start at row 1 + roff.
(optional, default = 0)
[in]coff: The input non-negative integer of default kind IK representing the offset of the top-left corner of the initialization block from the first column of the matrix.
The initialization of the matrix will start at column 1 + coff.
(optional, default = 0)
[in]doff: The input integer of default kind IK representing the offset of the diagonal of the initialization block from the top-left corner of the initialization block.
  • Setting doff > 0 implies a diagonal start with column offset doff from the top-left corner of the initialization block.
  • Setting doff < 0 implies a diagonal start with row offset -doff from the top-left corner of the initialization block.
  • Setting doff = 0 implies the same diagonal start as the top-left corner of the initialization block.
(optional, default = 0. It can be present only if the input arguments vupp or vlow are present.)
Returns
mat : The output contiguous matrix of arbitrary shape (:, :) of,
  1. type character of kind any supported by the processor (e.g., SK, SKA, SKD , or SKU) with arbitrary len type parameter, or
  2. type integer of kind any supported by the processor (e.g., IK, IK8, IK16, IK32, or IK64), or
  3. type logical of kind any supported by the processor (e.g., LK), or
  4. type complex of kind any supported by the processor (e.g., CK, CK32, CK64, or CK128), or
  5. type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128).
On output, the upper/lower triangle and diagonal elements of mat (as specified by the input argument subset) are set to the corresponding input values.
All other matrix elements remain uninitialized on output.


Possible calling interfaces

use pm_matrixInit, only: getMatInit, dia, uppDia, lowDia, uppLow, uppLowDia
mat(@shape) = getMatInit(shape, subset, vupp, roff = roff, coff = coff, doff = doff) ! subset = upp: Initialize upper-triangular subset of matrix.
mat(@shape) = getMatInit(shape, subset, vlow, roff = roff, coff = coff, doff = doff) ! subset = low: Initialize lower-triangular subset of matrix.
mat(@shape) = getMatInit(shape, subset, vdia , ndia = ndia, roff = roff, coff = coff) ! Initialize diagonal subset of matrix.
mat(@shape) = getMatInit(shape, subset, vdia(:) , roff = roff, coff = coff) ! Initialize diagonal subset of matrix.
mat(@shape) = getMatInit(shape, subset, vlow, vdia , nrow = nrow, ncol = ncol, roff = roff, coff = coff, doff = doff) ! Initialize lower-diagonal subset of matrix.
mat(@shape) = getMatInit(shape, subset, vlow, vdia(:), nrow = nrow, ncol = ncol, roff = roff, coff = coff, doff = doff) ! Initialize lower-diagonal subset of matrix.
mat(@shape) = getMatInit(shape, subset, vupp, vdia , nrow = nrow, ncol = ncol, roff = roff, coff = coff, doff = doff) ! Initialize upper-diagonal subset of matrix.
mat(@shape) = getMatInit(shape, subset, vupp, vdia(:), nrow = nrow, ncol = ncol, roff = roff, coff = coff, doff = doff) ! Initialize upper-diagonal subset of matrix.
mat(@shape) = getMatInit(shape, subset, vupp, vlow, nrow = nrow, ncol = ncol, roff = roff, coff = coff, doff = doff) ! Initialize upper-lower subset of matrix.
mat(@shape) = getMatInit(shape, subset, vupp, vlow, nrow = nrow, ncol = ncol, roff = roff, coff = coff, doff = doff) ! Initialize upper-lower subset of matrix.
mat(@shape) = getMatInit(shape, subset, vupp, vlow, vdia , nrow = nrow, ncol = ncol, roff = roff, coff = coff, doff = doff) ! Initialize upper-lower-diagonal subset of matrix.
mat(@shape) = getMatInit(shape, subset, vupp, vlow, vdia(:), nrow = nrow, ncol = ncol, roff = roff, coff = coff, doff = doff) ! Initialize upper-lower-diagonal subset of matrix.
!
Generate and return a matrix of shape (shape(1), shape(2)) with the upper/lower triangle and the diag...
This module contains procedures and generic interfaces relevant to generating and initializing matric...
Warning
All warnings associated with setMatInit also apply to this generic interface.
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.
BLAS/LAPACK equivalent:
The procedures under discussion combine, modernize, and extend the interface and functionalities of Version 3.11 of BLAS/LAPACK routine(s): SLASET, DLASET, CLASET, ZLASET.
See also
setMatInit
pm_arrayInit
pm_matrixCopy
pm_arrayCopy
pm_arrayCopy


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK, CK, RK
4 use pm_matrixInit, only: getMatInit, dia, uppDia, lowDia, uppLow, uppLowDia
5 use pm_io, only: display_type
6 use pm_arrayRange, only: getRange
7
8 implicit none
9
10 type(display_type) :: disp
11 disp = display_type(file = "main.out.F90")
12
13 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
14
15 call disp%skip()
16 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
17 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
18 call disp%show("! Initialize diagonal matrices.")
19 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
20 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
21 call disp%skip()
22
23 call disp%skip()
24 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
25 call disp%show("! Construct diagonal string matrix.")
26 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
27 call disp%skip()
28
29 block
30
31 character(2) :: mat(10,10)
32
33 call disp%skip()
34 call disp%show("mat = getMatInit(shape(mat, IK), dia, vdia = 'OO')")
35 mat = getMatInit(shape(mat, IK), dia, vdia = 'OO')
36 call disp%show("where(mat /= 'OO'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.")
37 where(mat /= 'OO'); mat = ' '; endwhere
38 call disp%show("mat")
39 call disp%show( mat , deliml = SK_"""" )
40 call disp%skip()
41
42 call disp%skip()
43 call disp%show("mat = getMatInit(shape(mat, IK), dia, vdia = 'AA', ndia = 5_IK, roff = 3_IK)")
44 mat = getMatInit(shape(mat, IK), dia, vdia = 'AA', ndia = 5_IK, roff = 3_IK)
45 call disp%show("where(mat /= 'AA'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.")
46 where(mat /= 'AA'); mat = ' '; endwhere
47 call disp%show("mat")
48 call disp%show( mat , deliml = SK_"""" )
49 call disp%skip()
50
51 call disp%skip()
52 call disp%show("mat = getMatInit(shape(mat, IK), dia, vdia = ['QQ', 'YY', 'WW', 'XX'], roff = 3_IK, coff = 2_IK)")
53 mat = getMatInit(shape(mat, IK), dia, vdia = ['QQ', 'YY', 'WW', 'XX'], roff = 3_IK, coff = 2_IK)
54 call disp%show(.and..and..and."where(mat /= 'QQ' mat /= 'YY' mat /= 'WW' mat /= 'XX'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.")
55 where(mat /= 'QQ' .and. mat /= 'YY' .and. mat /= 'WW' .and. mat /= 'XX'); mat = ' '; endwhere
56 call disp%show("mat")
57 call disp%show( mat , deliml = SK_"""" )
58 call disp%skip()
59
60 end block
61
62 call disp%skip()
63 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
64 call disp%show("! Construct diagonal integer matrix.")
65 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
66 call disp%skip()
67
68 block
69
70 integer :: mat(10,10)
71
72 call disp%skip()
73 call disp%show("mat = getMatInit(shape(mat, IK), dia, vdia = -2)")
74 mat = getMatInit(shape(mat, IK), dia, vdia = -2)
75 call disp%show("where(mat /= -2); mat = 0; endwhere ! get rid of uninitialized values for better illustration.")
76 where(mat /= -2); mat = 0; endwhere
77 call disp%show("mat")
78 call disp%show( mat )
79 call disp%skip()
80
81 call disp%skip()
82 call disp%show("mat = getMatInit(shape(mat, IK), dia, vdia = getRange(1, 7), roff = 3_IK, coff = 3_IK)")
83 mat = getMatInit(shape(mat, IK), dia, vdia = getRange(1, 7), roff = 3_IK, coff = 3_IK)
84 call disp%show("mat")
85 call disp%show( mat )
86 call disp%skip()
87
88 end block
89
90 call disp%skip()
91 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
92 call disp%show("! Construct diagonal logical matrix.")
93 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
94 call disp%skip()
95
96 block
97
98 logical :: mat(10,10)
99
100 call disp%skip()
101 call disp%show("mat = getMatInit(shape(mat, IK), dia, vdia = .true.)")
102 mat = getMatInit(shape(mat, IK), dia, vdia = .true.)
103 call disp%show("mat")
104 call disp%show( mat )
105 call disp%skip()
106
107 call disp%skip()
108 call disp%show("mat = getMatInit(shape(mat, IK), dia, vdia = [.true., .false., .true., .false.], roff = 3_IK, coff = 3_IK)")
109 mat = getMatInit(shape(mat, IK), dia, vdia = [.true., .false., .true., .false.], roff = 3_IK, coff = 3_IK)
110 call disp%show("mat")
111 call disp%show( mat )
112 call disp%skip()
113
114 end block
115
116 call disp%skip()
117 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
118 call disp%show("! Construct diagonal real matrix.")
119 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
120 call disp%skip()
121
122 block
123
124 real :: mat(10,10)
125
126 call disp%skip()
127 call disp%show("mat = getMatInit(shape(mat, IK), dia, vdia = -2., ndia = size(mat,1,IK))")
128 mat = getMatInit(shape(mat, IK), dia, vdia = -2., ndia = size(mat,1,IK))
129 call disp%show("where(mat /= -2.); mat = 0.; endwhere ! get rid of uninitialized values for better illustration.")
130 where(mat /= -2.); mat = 0.; endwhere
131 call disp%show("mat")
132 call disp%show( mat )
133 call disp%skip()
134
135 call disp%skip()
136 call disp%show("mat = getMatInit(shape(mat, IK), dia, vdia = -4., ndia = 7_IK, coff = 3_IK)")
137 mat = getMatInit(shape(mat, IK), dia, vdia = -4., ndia = 7_IK, coff = 3_IK)
138 call disp%show("where(mat /= -4.); mat = 0.; endwhere ! get rid of uninitialized values for better illustration.")
139 where(mat /= -4.); mat = 0.; endwhere
140 call disp%show("mat")
141 call disp%show( mat )
142 call disp%skip()
143
144 end block
145
146 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
147
148 call disp%skip()
149 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
150 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
151 call disp%show("! Initialize upper-diagonal matrices.")
152 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
153 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
154 call disp%skip()
155
156 call disp%skip()
157 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
158 call disp%show("! Construct upper-diagonal string matrix.")
159 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
160 call disp%skip()
161
162 block
163
164 character(2) :: mat(10,10)
165
166 call disp%skip()
167 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = 'AA', vdia = 'ZZ', nrow = 4_IK, ncol = 7_IK, roff = 3_IK)")
168 mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = 'AA', vdia = 'ZZ', nrow = 4_IK, ncol = 7_IK, roff = 3_IK)
169 call disp%show(.and."where(mat /= 'AA' mat /= 'ZZ'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.")
170 where(mat /= 'AA' .and. mat /= 'ZZ'); mat = ' '; endwhere
171 call disp%show("mat")
172 call disp%show( mat , deliml = SK_"""" )
173 call disp%skip()
174 mat = ""
175
176 call disp%skip()
177 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = 'AA', vdia = ['ZZ', 'YY', 'WW', 'XX'], nrow = 4_IK, ncol = 7_IK, roff = 3_IK)")
178 mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = 'AA', vdia = ['ZZ', 'YY', 'WW', 'XX'], nrow = 4_IK, ncol = 7_IK, roff = 3_IK)
179 call disp%show(.and..and..and..and."where(mat /= 'AA' mat /= 'ZZ' mat /= 'YY' mat /= 'WW' mat /= 'XX'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.")
180 where(mat /= 'AA' .and. mat /= 'ZZ' .and. mat /= 'YY' .and. mat /= 'WW' .and. mat /= 'XX'); mat = ' '; endwhere
181 call disp%show("mat")
182 call disp%show( mat , deliml = SK_"""" )
183 call disp%skip()
184 mat = ""
185
186 call disp%skip()
187 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = 'BB', vdia = 'YY', nrow = 4_IK, ncol = 7_IK, roff = 2_IK, coff = 3_IK, doff = -1_IK)")
188 mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = 'BB', vdia = 'YY', nrow = 4_IK, ncol = 7_IK, roff = 2_IK, coff = 3_IK, doff = -1_IK)
189 call disp%show(.and."where(mat /= 'BB' mat /= 'YY'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.")
190 where(mat /= 'BB' .and. mat /= 'YY'); mat = ' '; endwhere
191 call disp%show("mat")
192 call disp%show( mat , deliml = SK_"""" )
193 call disp%skip()
194 mat = ""
195
196 call disp%skip()
197 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = 'BB', vdia = 'YY', nrow = 5_IK, ncol = 7_IK, roff = 2_IK, coff = 3_IK, doff = -4_IK)")
198 mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = 'BB', vdia = 'YY', nrow = 5_IK, ncol = 7_IK, roff = 2_IK, coff = 3_IK, doff = -4_IK)
199 call disp%show(.and."where(mat /= 'BB' mat /= 'YY'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.")
200 where(mat /= 'BB' .and. mat /= 'YY'); mat = ' '; endwhere
201 call disp%show("mat")
202 call disp%show( mat , deliml = SK_"""" )
203 call disp%skip()
204 mat = ""
205
206 end block
207
208 call disp%skip()
209 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
210 call disp%show("! Construct upper-diagonal integer matrix.")
211 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
212 call disp%skip()
213
214 block
215
216 integer :: mat(10,10)
217
218 call disp%skip()
219 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = -1, vdia = -2)")
220 mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = -1, vdia = -2)
221 call disp%show(.and."where(mat /= -1 mat /= -2); mat = 0; endwhere ! get rid of uninitialized values for better illustration.")
222 where(mat /= -1 .and. mat /= -2); mat = 0; endwhere
223 call disp%show("mat")
224 call disp%show( mat )
225 call disp%skip()
226 mat = 0
227
228 call disp%skip()
229 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = -3, vdia = -4, nrow = 4_IK, ncol = 7_IK, roff = 3_IK, coff = 3_IK)")
230 mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = -3, vdia = -4, nrow = 4_IK, ncol = 7_IK, roff = 3_IK, coff = 3_IK)
231 call disp%show(.and."where(mat /= -3 mat /= -4); mat = 0; endwhere ! get rid of uninitialized values for better illustration.")
232 where(mat /= -3 .and. mat /= -4); mat = 0; endwhere
233 call disp%show("mat")
234 call disp%show( mat )
235 call disp%skip()
236 mat = 0
237
238 end block
239
240 call disp%skip()
241 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
242 call disp%show("! Construct upper-diagonal logical matrix.")
243 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
244 call disp%skip()
245
246 block
247
248 logical :: mat(10,10)
249
250 mat = .false.
251 call disp%skip()
252 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = .true., vdia = .true.)")
253 mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = .true., vdia = .true.)
254 call disp%show("mat")
255 call disp%show( mat )
256 call disp%skip()
257 mat = .false.
258
259 mat = .false.
260 call disp%skip()
261 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = .true., vdia = .true., nrow = 4_IK, ncol = 7_IK, roff = 3_IK, coff = 3_IK)")
262 mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = .true., vdia = .true., nrow = 4_IK, ncol = 7_IK, roff = 3_IK, coff = 3_IK)
263 call disp%show("mat")
264 call disp%show( mat )
265 call disp%skip()
266 mat = .false.
267
268 end block
269
270 call disp%skip()
271 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
272 call disp%show("! Construct upper-diagonal real matrix.")
273 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
274 call disp%skip()
275
276 block
277
278 real :: mat(10,10)
279
280 call disp%skip()
281 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = -1., vdia = -2.)")
282 mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = -1., vdia = -2.)
283 call disp%show(.and."where(mat /= -1 mat /= -2); mat = 0; endwhere ! get rid of uninitialized values for better illustration.")
284 where(mat /= -1 .and. mat /= -2); mat = 0; endwhere
285 call disp%show("mat")
286 call disp%show( mat )
287 call disp%skip()
288 mat = 0
289
290 call disp%skip()
291 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = -3., vdia = -4., nrow = 3_IK, ncol = 7_IK, coff = 3_IK)")
292 mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = -3., vdia = -4., nrow = 3_IK, ncol = 7_IK, coff = 3_IK)
293 call disp%show(.and."where(mat /= -3. mat /= -4.); mat = 0; endwhere ! get rid of uninitialized values for better illustration.")
294 where(mat /= -3. .and. mat /= -4.); mat = 0; endwhere
295 call disp%show("mat")
296 call disp%show( mat )
297 call disp%skip()
298 mat = 0
299
300 end block
301
302 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
303
304 call disp%skip()
305 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
306 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
307 call disp%show("! Initialize lower-diagonal matrices.")
308 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
309 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
310 call disp%skip()
311
312 call disp%skip()
313 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
314 call disp%show("! Construct lower-diagonal string matrix.")
315 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
316 call disp%skip()
317
318 block
319
320 character(2) :: mat(10,10)
321
322 call disp%skip()
323 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = 'AA', vdia = 'ZZ', nrow = 7_IK, ncol = 4_IK, roff = 3_IK)")
324 mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = 'AA', vdia = 'ZZ', nrow = 7_IK, ncol = 4_IK, roff = 3_IK)
325 call disp%show(.and."where(mat /= 'AA' mat /= 'ZZ'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.")
326 where(mat /= 'AA' .and. mat /= 'ZZ'); mat = ' '; endwhere
327 call disp%show("mat")
328 call disp%show( mat , deliml = SK_"""" )
329 call disp%skip()
330 mat = ""
331
332 call disp%skip()
333 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = 'AA', vdia = ['ZZ', 'YY', 'WW', 'XX'], nrow = 7_IK, ncol = 4_IK, roff = 3_IK)")
334 mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = 'AA', vdia = ['ZZ', 'YY', 'WW', 'XX'], nrow = 7_IK, ncol = 4_IK, roff = 3_IK)
335 call disp%show(.and..and..and..and."where(mat /= 'AA' mat /= 'ZZ' mat /= 'YY' mat /= 'WW' mat /= 'XX'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.")
336 where(mat /= 'AA' .and. mat /= 'ZZ' .and. mat /= 'YY' .and. mat /= 'WW' .and. mat /= 'XX'); mat = ' '; endwhere
337 call disp%show("mat")
338 call disp%show( mat , deliml = SK_"""" )
339 call disp%skip()
340 mat = ""
341
342 call disp%skip()
343 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = 'BB', vdia = 'YY', nrow = 4_IK, ncol = 6_IK, roff = 2_IK, coff = 3_IK, doff = +2_IK)")
344 mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = 'BB', vdia = 'YY', nrow = 4_IK, ncol = 6_IK, roff = 2_IK, coff = 3_IK, doff = +2_IK)
345 call disp%show(.and."where(mat /= 'BB' mat /= 'YY'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.")
346 where(mat /= 'BB' .and. mat /= 'YY'); mat = ' '; endwhere
347 call disp%show("mat")
348 call disp%show( mat , deliml = SK_"""" )
349 call disp%skip()
350 mat = ""
351
352 call disp%skip()
353 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = 'BB', vdia = 'YY', nrow = 4_IK, ncol = 7_IK, roff = 2_IK, coff = 3_IK, doff = +4_IK)")
354 mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = 'BB', vdia = 'YY', nrow = 4_IK, ncol = 7_IK, roff = 2_IK, coff = 3_IK, doff = +4_IK)
355 call disp%show(.and."where(mat /= 'BB' mat /= 'YY'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.")
356 where(mat /= 'BB' .and. mat /= 'YY'); mat = ' '; endwhere
357 call disp%show("mat")
358 call disp%show( mat , deliml = SK_"""" )
359 call disp%skip()
360 mat = ""
361
362 end block
363
364 call disp%skip()
365 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
366 call disp%show("! Construct lower-diagonal integer matrix.")
367 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
368 call disp%skip()
369
370 block
371
372 integer :: mat(10,10)
373
374 call disp%skip()
375 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = -1, vdia = -2)")
376 mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = -1, vdia = -2)
377 call disp%show(.and."where(mat /= -1 mat /= -2); mat = 0; endwhere ! get rid of uninitialized values for better illustration.")
378 where(mat /= -1 .and. mat /= -2); mat = 0; endwhere
379 call disp%show("mat")
380 call disp%show( mat )
381 call disp%skip()
382 mat = 0
383
384 call disp%skip()
385 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = -3, vdia = -4, nrow = 7_IK, ncol = 4_IK, roff = 3_IK, coff = 3_IK)")
386 mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = -3, vdia = -4, nrow = 7_IK, ncol = 4_IK, roff = 3_IK, coff = 3_IK)
387 call disp%show(.and."where(mat /= -3 mat /= -4); mat = 0; endwhere ! get rid of uninitialized values for better illustration.")
388 where(mat /= -3 .and. mat /= -4); mat = 0; endwhere
389 call disp%show("mat")
390 call disp%show( mat )
391 call disp%skip()
392 mat = 0
393
394 end block
395
396 call disp%skip()
397 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
398 call disp%show("! Construct lower-diagonal logical matrix.")
399 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
400 call disp%skip()
401
402 block
403
404 logical :: mat(10,10)
405
406 mat = .false.
407 call disp%skip()
408 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = .true., vdia = .true.)")
409 mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = .true., vdia = .true.)
410 call disp%show("mat")
411 call disp%show( mat )
412 call disp%skip()
413 mat = .false.
414
415 mat = .false.
416 call disp%skip()
417 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = .true., vdia = .true., nrow = 7_IK, ncol = 4_IK, roff = 3_IK, coff = 3_IK)")
418 mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = .true., vdia = .true., nrow = 7_IK, ncol = 4_IK, roff = 3_IK, coff = 3_IK)
419 call disp%show("mat")
420 call disp%show( mat )
421 call disp%skip()
422 mat = .false.
423
424 end block
425
426 call disp%skip()
427 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
428 call disp%show("! Construct lower-diagonal real matrix.")
429 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
430 call disp%skip()
431
432 block
433
434 real :: mat(10,10)
435
436 call disp%skip()
437 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = -1., vdia = -2.)")
438 mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = -1., vdia = -2.)
439 call disp%show(.and."where(mat /= -1 mat /= -2); mat = 0; endwhere ! get rid of uninitialized values for better illustration.")
440 where(mat /= -1 .and. mat /= -2); mat = 0; endwhere
441 call disp%show("mat")
442 call disp%show( mat )
443 call disp%skip()
444 mat = 0
445
446 call disp%skip()
447 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = -3., vdia = -4., nrow = 7_IK, ncol = 3_IK, coff = 3_IK)")
448 mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = -3., vdia = -4., nrow = 7_IK, ncol = 3_IK, coff = 3_IK)
449 call disp%show(.and."where(mat /= -3. mat /= -4.); mat = 0; endwhere ! get rid of uninitialized values for better illustration.")
450 where(mat /= -3. .and. mat /= -4.); mat = 0; endwhere
451 call disp%show("mat")
452 call disp%show( mat )
453 call disp%skip()
454 mat = 0
455
456 end block
457
458 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
459
460 call disp%skip()
461 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
462 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
463 call disp%show("! Initialize upper-lower matrices.")
464 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
465 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
466 call disp%skip()
467
468 call disp%skip()
469 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
470 call disp%show("! Construct upper-lower string matrix.")
471 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
472 call disp%skip()
473
474 block
475
476 character(2) :: mat(10,10)
477
478 mat = ""
479 call disp%skip()
480 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = 'AA', vlow = 'HH', nrow = 7_IK, ncol = 4_IK, roff = 3_IK)")
481 mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = 'AA', vlow = 'HH', nrow = 7_IK, ncol = 4_IK, roff = 3_IK)
482 call disp%show(.and..and."where(mat /= 'AA' mat /= 'HH' mat /= 'OO'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.")
483 where(mat /= 'AA' .and. mat /= 'HH' .and. mat /= 'OO'); mat = ' '; endwhere
484 call disp%show("mat")
485 call disp%show( mat , deliml = SK_"""" )
486 call disp%skip()
487
488 mat = ""
489 call disp%skip()
490 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = 'AA', vlow = 'HH', nrow = 7_IK, ncol = 4_IK, roff = 3_IK)")
491 mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = 'AA', vlow = 'HH', nrow = 7_IK, ncol = 4_IK, roff = 3_IK)
492 call disp%show(.and."where(mat /= 'AA' mat /= 'HH'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.")
493 where(mat /= 'AA' .and. mat /= 'HH'); mat = ' '; endwhere
494 call disp%show("mat")
495 call disp%show( mat , deliml = SK_"""" )
496 call disp%skip()
497
498 mat = ""
499 call disp%skip()
500 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = 'AA', vlow = 'HH', nrow = 4_IK, ncol = 7_IK, roff = 3_IK)")
501 mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = 'AA', vlow = 'HH', nrow = 4_IK, ncol = 7_IK, roff = 3_IK)
502 call disp%show(.and."where(mat /= 'AA' mat /= 'HH'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.")
503 where(mat /= 'AA' .and. mat /= 'HH'); mat = ' '; endwhere
504 call disp%show("mat")
505 call disp%show( mat , deliml = SK_"""" )
506 call disp%skip()
507
508 mat = ""
509 call disp%skip()
510 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = 'AA', vlow = 'BB', nrow = 4_IK, ncol = 6_IK, roff = 2_IK, coff = 3_IK, doff = +2_IK)")
511 mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = 'AA', vlow = 'BB', nrow = 4_IK, ncol = 6_IK, roff = 2_IK, coff = 3_IK, doff = +2_IK)
512 call disp%show(.and..and."where(mat /= 'AA' mat /= 'BB' mat /= 'YY'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.")
513 where(mat /= 'AA' .and. mat /= 'BB' .and. mat /= 'YY'); mat = ' '; endwhere
514 call disp%show("mat")
515 call disp%show( mat , deliml = SK_"""" )
516 call disp%skip()
517
518 mat = ""
519 call disp%skip()
520 call disp%show("mat = getMatInit(vupp = 'AA', vlow = 'BB', Shape = [10_IK, 10_IK], subset = uppLow, nrow = 4_IK, ncol = 7_IK, roff = 2_IK, coff = 3_IK, doff = +4_IK)")
521 mat = getMatInit(vupp = 'AA', vlow = 'BB', Shape = [10_IK, 10_IK], subset = uppLow, nrow = 4_IK, ncol = 7_IK, roff = 2_IK, coff = 3_IK, doff = +4_IK)
522 call disp%show(.and..and."where(mat /= 'AA' mat /= 'BB' mat /= 'YY'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.")
523 where(mat /= 'AA' .and. mat /= 'BB' .and. mat /= 'YY'); mat = ' '; endwhere
524 call disp%show("mat")
525 call disp%show( mat , deliml = SK_"""" )
526 call disp%skip()
527
528 end block
529
530 call disp%skip()
531 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
532 call disp%show("! Construct upper-lower integer matrix.")
533 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
534 call disp%skip()
535
536 block
537
538 integer :: mat(10,10)
539
540 call disp%skip()
541 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = +1, vlow = -1)")
542 mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = +1, vlow = -1)
543 call disp%show(.and..and."where(mat /= +1 mat /= -1 mat /= -2); mat = 0; endwhere ! get rid of uninitialized values for better illustration.")
544 where(mat /= +1 .and. mat /= -1 .and. mat /= -2); mat = 0; endwhere
545 call disp%show("mat")
546 call disp%show( mat )
547 call disp%skip()
548 mat = 0
549
550 call disp%skip()
551 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = +1, vlow = -3, nrow = 7_IK, ncol = 4_IK, roff = 3_IK, coff = 3_IK)")
552 mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = +1, vlow = -3, nrow = 7_IK, ncol = 4_IK, roff = 3_IK, coff = 3_IK)
553 call disp%show(.and..and."where(mat /= +1 mat /= -3 mat /= -4); mat = 0; endwhere ! get rid of uninitialized values for better illustration.")
554 where(mat /= +1 .and. mat /= -3 .and. mat /= -4); mat = 0; endwhere
555 call disp%show("mat")
556 call disp%show( mat )
557 call disp%skip()
558 mat = 0
559
560 end block
561
562 call disp%skip()
563 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
564 call disp%show("! Construct upper-lower logical matrix.")
565 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
566 call disp%skip()
567
568 block
569
570 logical :: mat(10,10)
571
572 call disp%skip()
573 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = .true., vlow = .true.)")
574 mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = .true., vlow = .true.)
575 call disp%show("mat")
576 call disp%show( mat )
577 call disp%skip()
578
579 call disp%skip()
580 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = .true., vlow = .true., nrow = 7_IK, ncol = 4_IK, roff = 3_IK, coff = 3_IK, doff = +1_IK)")
581 mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = .true., vlow = .true., nrow = 7_IK, ncol = 4_IK, roff = 3_IK, coff = 3_IK, doff = +1_IK)
582 call disp%show("mat")
583 call disp%show( mat )
584 call disp%skip()
585 mat = .false.
586
587 end block
588
589 call disp%skip()
590 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
591 call disp%show("! Construct upper-lower real matrix.")
592 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
593 call disp%skip()
594
595 block
596
597 real :: mat(10,10)
598
599 call disp%skip()
600 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = +1., vlow = -1.)")
601 mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = +1., vlow = -1.)
602 call disp%show(.and..and."where(mat /= +1 mat /= -1 mat /= -2); mat = 0; endwhere ! get rid of uninitialized values for better illustration.")
603 where(mat /= +1 .and. mat /= -1 .and. mat /= -2); mat = 0; endwhere
604 call disp%show("mat")
605 call disp%show( mat )
606 call disp%skip()
607 mat = 0
608
609 call disp%skip()
610 call disp%show("mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = +3., vlow = -3., nrow = 7_IK, ncol = 3_IK, coff = 3_IK)")
611 mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = +3., vlow = -3., nrow = 7_IK, ncol = 3_IK, coff = 3_IK)
612 call disp%show(.and..and."where(mat /= +3. mat /= -3. mat /= -4.); mat = 0; endwhere ! get rid of uninitialized values for better illustration.")
613 where(mat /= +3. .and. mat /= -3. .and. mat /= -4.); mat = 0; endwhere
614 call disp%show("mat")
615 call disp%show( mat )
616 call disp%skip()
617 mat = 0
618
619 end block
620
621 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
622
623 call disp%skip()
624 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
625 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
626 call disp%show("! Initialize upper-lower-diagonal matrices.")
627 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
628 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
629 call disp%skip()
630
631 call disp%skip()
632 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
633 call disp%show("! Construct upper-lower-diagonal string matrix.")
634 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
635 call disp%skip()
636
637 block
638
639 character(2) :: mat(10,10)
640
641 mat = ""
642 call disp%skip()
643 call disp%show("mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = 'AA', vlow = 'HH', vdia = 'OO', nrow = 7_IK, ncol = 4_IK, roff = 3_IK)")
644 mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = 'AA', vlow = 'HH', vdia = 'OO', nrow = 7_IK, ncol = 4_IK, roff = 3_IK)
645 call disp%show(.and..and."where(mat /= 'AA' mat /= 'HH' mat /= 'OO'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.")
646 where(mat /= 'AA' .and. mat /= 'HH' .and. mat /= 'OO'); mat = ' '; endwhere
647 call disp%show("mat")
648 call disp%show( mat , deliml = SK_"""" )
649 call disp%skip()
650
651 mat = ""
652 call disp%skip()
653 call disp%show("mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = 'AA', vlow = 'HH', vdia = ['OO', 'YY', 'WW', 'XX'], nrow = 7_IK, ncol = 4_IK, roff = 3_IK)")
654 mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = 'AA', vlow = 'HH', vdia = ['OO', 'YY', 'WW', 'XX'], nrow = 7_IK, ncol = 4_IK, roff = 3_IK)
655 call disp%show(.and..and..and..and..and."where(mat /= 'AA' mat /= 'HH' mat /= 'OO' mat /= 'YY' mat /= 'WW' mat /= 'XX'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.")
656 where(mat /= 'AA' .and. mat /= 'HH' .and. mat /= 'OO' .and. mat /= 'YY' .and. mat /= 'WW' .and. mat /= 'XX'); mat = ' '; endwhere
657 call disp%show("mat")
658 call disp%show( mat , deliml = SK_"""" )
659 call disp%skip()
660
661 mat = ""
662 call disp%skip()
663 call disp%show("mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = 'AA', vlow = 'HH', vdia = ['OO', 'YY', 'WW', 'XX'], nrow = 4_IK, ncol = 7_IK, roff = 3_IK)")
664 mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = 'AA', vlow = 'HH', vdia = ['OO', 'YY', 'WW', 'XX'], nrow = 4_IK, ncol = 7_IK, roff = 3_IK)
665 call disp%show(.and..and..and..and..and."where(mat /= 'AA' mat /= 'HH' mat /= 'OO' mat /= 'YY' mat /= 'WW' mat /= 'XX'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.")
666 where(mat /= 'AA' .and. mat /= 'HH' .and. mat /= 'OO' .and. mat /= 'YY' .and. mat /= 'WW' .and. mat /= 'XX'); mat = ' '; endwhere
667 call disp%show("mat")
668 call disp%show( mat , deliml = SK_"""" )
669 call disp%skip()
670
671 mat = ""
672 call disp%skip()
673 call disp%show("mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = 'AA', vlow = 'BB', vdia = 'YY', nrow = 4_IK, ncol = 6_IK, roff = 2_IK, coff = 3_IK, doff = +2_IK)")
674 mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = 'AA', vlow = 'BB', vdia = 'YY', nrow = 4_IK, ncol = 6_IK, roff = 2_IK, coff = 3_IK, doff = +2_IK)
675 call disp%show(.and..and."where(mat /= 'AA' mat /= 'BB' mat /= 'YY'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.")
676 where(mat /= 'AA' .and. mat /= 'BB' .and. mat /= 'YY'); mat = ' '; endwhere
677 call disp%show("mat")
678 call disp%show( mat , deliml = SK_"""" )
679 call disp%skip()
680
681 mat = ""
682 call disp%skip()
683 call disp%show("mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = 'AA', vlow = 'BB', vdia = 'YY', nrow = 4_IK, ncol = 7_IK, roff = 2_IK, coff = 3_IK, doff = +4_IK)")
684 mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = 'AA', vlow = 'BB', vdia = 'YY', nrow = 4_IK, ncol = 7_IK, roff = 2_IK, coff = 3_IK, doff = +4_IK)
685 call disp%show(.and..and."where(mat /= 'AA' mat /= 'BB' mat /= 'YY'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.")
686 where(mat /= 'AA' .and. mat /= 'BB' .and. mat /= 'YY'); mat = ' '; endwhere
687 call disp%show("mat")
688 call disp%show( mat , deliml = SK_"""" )
689 call disp%skip()
690
691 end block
692
693 call disp%skip()
694 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
695 call disp%show("! Construct upper-lower-diagonal integer matrix.")
696 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
697 call disp%skip()
698
699 block
700
701 integer :: mat(10,10)
702
703 call disp%skip()
704 call disp%show("mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = +1, vlow = -1, vdia = -2)")
705 mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = +1, vlow = -1, vdia = -2)
706 call disp%show(.and..and."where(mat /= +1 mat /= -1 mat /= -2); mat = 0; endwhere ! get rid of uninitialized values for better illustration.")
707 where(mat /= +1 .and. mat /= -1 .and. mat /= -2); mat = 0; endwhere
708 call disp%show("mat")
709 call disp%show( mat )
710 call disp%skip()
711 mat = 0
712
713 call disp%skip()
714 call disp%show("mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = +1, vlow = -3, vdia = -4, nrow = 7_IK, ncol = 4_IK, roff = 3_IK, coff = 3_IK)")
715 mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = +1, vlow = -3, vdia = -4, nrow = 7_IK, ncol = 4_IK, roff = 3_IK, coff = 3_IK)
716 call disp%show(.and..and."where(mat /= +1 mat /= -3 mat /= -4); mat = 0; endwhere ! get rid of uninitialized values for better illustration.")
717 where(mat /= +1 .and. mat /= -3 .and. mat /= -4); mat = 0; endwhere
718 call disp%show("mat")
719 call disp%show( mat )
720 call disp%skip()
721 mat = 0
722
723 end block
724
725 call disp%skip()
726 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
727 call disp%show("! Construct upper-lower-diagonal logical matrix.")
728 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
729 call disp%skip()
730
731 block
732
733 logical :: mat(10,10)
734
735 call disp%skip()
736 call disp%show("mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = .true., vlow = .true., vdia = .false.)")
737 mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = .true., vlow = .true., vdia = .false.)
738 call disp%show("mat")
739 call disp%show( mat )
740 call disp%skip()
741
742 call disp%skip()
743 call disp%show("mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = .true., vlow = .true., vdia = .false., nrow = 7_IK, ncol = 4_IK, roff = 3_IK, coff = 3_IK, doff = +1_IK)")
744 mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = .true., vlow = .true., vdia = .false., nrow = 7_IK, ncol = 4_IK, roff = 3_IK, coff = 3_IK, doff = +1_IK)
745 call disp%show("mat")
746 call disp%show( mat )
747 call disp%skip()
748 mat = .false.
749
750 end block
751
752 call disp%skip()
753 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
754 call disp%show("! Construct upper-lower-diagonal real matrix.")
755 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
756 call disp%skip()
757
758 block
759
760 real :: mat(10,10)
761
762 call disp%skip()
763 call disp%show("mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = +1., vlow = -1., vdia = -2.)")
764 mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = +1., vlow = -1., vdia = -2.)
765 call disp%show(.and..and."where(mat /= +1 mat /= -1 mat /= -2); mat = 0; endwhere ! get rid of uninitialized values for better illustration.")
766 where(mat /= +1 .and. mat /= -1 .and. mat /= -2); mat = 0; endwhere
767 call disp%show("mat")
768 call disp%show( mat )
769 call disp%skip()
770 mat = 0
771
772 call disp%skip()
773 call disp%show("mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = +3., vlow = -3., vdia = -4., nrow = 7_IK, ncol = 3_IK, coff = 3_IK)")
774 mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = +3., vlow = -3., vdia = -4., nrow = 7_IK, ncol = 3_IK, coff = 3_IK)
775 call disp%show(.and..and."where(mat /= +3. mat /= -3. mat /= -4.); mat = 0; endwhere ! get rid of uninitialized values for better illustration.")
776 where(mat /= +3. .and. mat /= -3. .and. mat /= -4.); mat = 0; endwhere
777 call disp%show("mat")
778 call disp%show( mat )
779 call disp%skip()
780 mat = 0
781
782 end block
783
784 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
785
786end program example
Generate minimally-spaced character, integer, real sequences or sequences at fixed intervals of size ...
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 procedures and generic interfaces for generating ranges of discrete character,...
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 RK
The default real kind in the ParaMonte library: real64 in Fortran, c_double in C-Fortran Interoperati...
Definition: pm_kind.F90:543
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
Definition: pm_kind.F90:541
integer, parameter CK
The default complex kind in the ParaMonte library: real64 in Fortran, c_double_complex in C-Fortran I...
Definition: pm_kind.F90:542
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
Definition: pm_kind.F90:540
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
Definition: pm_kind.F90:539
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!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
4! Initialize diagonal matrices.
5!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
6!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
7
8
9!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
10! Construct diagonal string matrix.
11!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
12
13
14mat = getMatInit(shape(mat, IK), dia, vdia = 'OO')
15where(mat /= 'OO'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.
16mat
17"OO", " ", " ", " ", " ", " ", " ", " ", " ", " "
18" ", "OO", " ", " ", " ", " ", " ", " ", " ", " "
19" ", " ", "OO", " ", " ", " ", " ", " ", " ", " "
20" ", " ", " ", "OO", " ", " ", " ", " ", " ", " "
21" ", " ", " ", " ", "OO", " ", " ", " ", " ", " "
22" ", " ", " ", " ", " ", "OO", " ", " ", " ", " "
23" ", " ", " ", " ", " ", " ", "OO", " ", " ", " "
24" ", " ", " ", " ", " ", " ", " ", "OO", " ", " "
25" ", " ", " ", " ", " ", " ", " ", " ", "OO", " "
26" ", " ", " ", " ", " ", " ", " ", " ", " ", "OO"
27
28
29mat = getMatInit(shape(mat, IK), dia, vdia = 'AA', ndia = 5_IK, roff = 3_IK)
30where(mat /= 'AA'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.
31mat
32" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
33" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
34" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
35"AA", " ", " ", " ", " ", " ", " ", " ", " ", " "
36" ", "AA", " ", " ", " ", " ", " ", " ", " ", " "
37" ", " ", "AA", " ", " ", " ", " ", " ", " ", " "
38" ", " ", " ", "AA", " ", " ", " ", " ", " ", " "
39" ", " ", " ", " ", "AA", " ", " ", " ", " ", " "
40" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
41" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
42
43
44mat = getMatInit(shape(mat, IK), dia, vdia = ['QQ', 'YY', 'WW', 'XX'], roff = 3_IK, coff = 2_IK)
45where(mat /= 'QQ' .and. mat /= 'YY' .and. mat /= 'WW' .and. mat /= 'XX'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.
46mat
47" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
48" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
49" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
50" ", " ", "QQ", " ", " ", " ", " ", " ", " ", " "
51" ", " ", " ", "YY", " ", " ", " ", " ", " ", " "
52" ", " ", " ", " ", "WW", " ", " ", " ", " ", " "
53" ", " ", " ", " ", " ", "XX", " ", " ", " ", " "
54" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
55" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
56" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
57
58
59!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
60! Construct diagonal integer matrix.
61!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
62
63
64mat = getMatInit(shape(mat, IK), dia, vdia = -2)
65where(mat /= -2); mat = 0; endwhere ! get rid of uninitialized values for better illustration.
66mat
67-2, +0, +0, +0, +0, +0, +0, +0, +0, +0
68+0, -2, +0, +0, +0, +0, +0, +0, +0, +0
69+0, +0, -2, +0, +0, +0, +0, +0, +0, +0
70+0, +0, +0, -2, +0, +0, +0, +0, +0, +0
71+0, +0, +0, +0, -2, +0, +0, +0, +0, +0
72+0, +0, +0, +0, +0, -2, +0, +0, +0, +0
73+0, +0, +0, +0, +0, +0, -2, +0, +0, +0
74+0, +0, +0, +0, +0, +0, +0, -2, +0, +0
75+0, +0, +0, +0, +0, +0, +0, +0, -2, +0
76+0, +0, +0, +0, +0, +0, +0, +0, +0, -2
77
78
79mat = getMatInit(shape(mat, IK), dia, vdia = getRange(1, 7), roff = 3_IK, coff = 3_IK)
80mat
81-2, +0, +0, +0, +0, +0, +0, +0, +0, +0
82+0, -2, +0, +0, +0, +0, +0, +0, +0, +0
83+0, +0, -2, +0, +0, +0, +0, +0, +0, +0
84+0, +0, +0, +1, +0, +0, +0, +0, +0, +0
85+0, +0, +0, +0, +2, +0, +0, +0, +0, +0
86+0, +0, +0, +0, +0, +3, +0, +0, +0, +0
87+0, +0, +0, +0, +0, +0, +4, +0, +0, +0
88+0, +0, +0, +0, +0, +0, +0, +5, +0, +0
89+0, +0, +0, +0, +0, +0, +0, +0, +6, +0
90+0, +0, +0, +0, +0, +0, +0, +0, +0, +7
91
92
93!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
94! Construct diagonal logical matrix.
95!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
96
97
98mat = getMatInit(shape(mat, IK), dia, vdia = .true.)
99mat
100T, F, F, F, F, F, F, F, F, F
101F, T, F, F, F, F, F, F, F, F
102F, F, T, F, F, F, F, F, F, F
103F, F, F, T, F, F, F, F, F, F
104F, F, F, F, T, F, F, F, F, F
105F, F, F, F, F, T, F, F, F, F
106F, F, F, F, F, F, T, F, F, F
107F, F, F, F, F, F, F, T, F, F
108F, F, F, F, F, F, F, F, T, F
109F, F, F, F, F, F, F, F, F, T
110
111
112mat = getMatInit(shape(mat, IK), dia, vdia = [.true., .false., .true., .false.], roff = 3_IK, coff = 3_IK)
113mat
114T, F, F, F, F, F, F, F, F, F
115F, T, F, F, F, F, F, F, F, F
116F, F, T, F, F, F, F, F, F, F
117F, F, F, T, F, F, F, F, F, F
118F, F, F, F, F, F, F, F, F, F
119F, F, F, F, F, T, F, F, F, F
120F, F, F, F, F, F, F, F, F, F
121F, F, F, F, F, F, F, T, F, F
122F, F, F, F, F, F, F, F, T, F
123F, F, F, F, F, F, F, F, F, T
124
125
126!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
127! Construct diagonal real matrix.
128!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
129
130
131mat = getMatInit(shape(mat, IK), dia, vdia = -2., ndia = size(mat,1,IK))
132where(mat /= -2.); mat = 0.; endwhere ! get rid of uninitialized values for better illustration.
133mat
134-2.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
135+0.00000000, -2.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
136+0.00000000, +0.00000000, -2.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
137+0.00000000, +0.00000000, +0.00000000, -2.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
138+0.00000000, +0.00000000, +0.00000000, +0.00000000, -2.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
139+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, -2.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
140+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, -2.00000000, +0.00000000, +0.00000000, +0.00000000
141+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, -2.00000000, +0.00000000, +0.00000000
142+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, -2.00000000, +0.00000000
143+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, -2.00000000
144
145
146mat = getMatInit(shape(mat, IK), dia, vdia = -4., ndia = 7_IK, coff = 3_IK)
147where(mat /= -4.); mat = 0.; endwhere ! get rid of uninitialized values for better illustration.
148mat
149+0.00000000, +0.00000000, +0.00000000, -4.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
150+0.00000000, +0.00000000, +0.00000000, +0.00000000, -4.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
151+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, -4.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
152+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, -4.00000000, +0.00000000, +0.00000000, +0.00000000
153+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, -4.00000000, +0.00000000, +0.00000000
154+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, -4.00000000, +0.00000000
155+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, -4.00000000
156+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
157+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
158+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
159
160
161!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
162!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
163! Initialize upper-diagonal matrices.
164!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
165!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
166
167
168!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
169! Construct upper-diagonal string matrix.
170!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
171
172
173mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = 'AA', vdia = 'ZZ', nrow = 4_IK, ncol = 7_IK, roff = 3_IK)
174where(mat /= 'AA' .and. mat /= 'ZZ'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.
175mat
176" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
177" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
178" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
179"ZZ", "AA", "AA", "AA", "AA", "AA", "AA", " ", " ", " "
180" ", "ZZ", "AA", "AA", "AA", "AA", "AA", " ", " ", " "
181" ", " ", "ZZ", "AA", "AA", "AA", "AA", " ", " ", " "
182" ", " ", " ", "ZZ", "AA", "AA", "AA", " ", " ", " "
183" ", " ", " ", " ", "AA", " ", " ", " ", " ", " "
184" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
185" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
186
187
188mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = 'AA', vdia = ['ZZ', 'YY', 'WW', 'XX'], nrow = 4_IK, ncol = 7_IK, roff = 3_IK)
189where(mat /= 'AA' .and. mat /= 'ZZ' .and. mat /= 'YY' .and. mat /= 'WW' .and. mat /= 'XX'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.
190mat
191" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
192" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
193" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
194"ZZ", "AA", "AA", "AA", "AA", "AA", "AA", " ", " ", " "
195" ", "YY", "AA", "AA", "AA", "AA", "AA", " ", " ", " "
196" ", " ", "WW", "AA", "AA", "AA", "AA", " ", " ", " "
197" ", " ", " ", "XX", "AA", "AA", "AA", " ", " ", " "
198" ", " ", " ", " ", "AA", " ", " ", " ", " ", " "
199" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
200" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
201
202
203mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = 'BB', vdia = 'YY', nrow = 4_IK, ncol = 7_IK, roff = 2_IK, coff = 3_IK, doff = -1_IK)
204where(mat /= 'BB' .and. mat /= 'YY'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.
205mat
206" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
207" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
208" ", " ", " ", "BB", "BB", "BB", "BB", "BB", "BB", "BB"
209" ", " ", " ", "YY", "BB", "BB", "BB", "BB", "BB", "BB"
210" ", " ", " ", " ", "YY", "BB", "BB", "BB", "BB", "BB"
211" ", " ", " ", " ", " ", "YY", "BB", "BB", "BB", "BB"
212" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
213" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
214" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
215" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
216
217
218mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = 'BB', vdia = 'YY', nrow = 5_IK, ncol = 7_IK, roff = 2_IK, coff = 3_IK, doff = -4_IK)
219where(mat /= 'BB' .and. mat /= 'YY'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.
220mat
221" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
222" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
223" ", " ", " ", "BB", "BB", "BB", "BB", "BB", "BB", "BB"
224" ", " ", " ", "BB", "BB", "BB", "BB", "BB", "BB", "BB"
225" ", " ", " ", "BB", "BB", "BB", "BB", "BB", "BB", "BB"
226" ", " ", " ", "BB", "BB", "BB", "BB", "BB", "BB", "BB"
227" ", " ", " ", "YY", "BB", "BB", "BB", "BB", "BB", "BB"
228" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
229" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
230" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
231
232
233!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
234! Construct upper-diagonal integer matrix.
235!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
236
237
238mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = -1, vdia = -2)
239where(mat /= -1 .and. mat /= -2); mat = 0; endwhere ! get rid of uninitialized values for better illustration.
240mat
241-2, -1, -1, -1, -1, -1, -1, -1, -1, -1
242+0, -2, -1, -1, -1, -1, -1, -1, -1, -1
243+0, +0, -2, -1, -1, -1, -1, -1, -1, -1
244+0, +0, +0, -2, -1, -1, -1, -1, -1, -1
245+0, +0, +0, +0, -2, -1, -1, -1, -1, -1
246+0, +0, +0, +0, +0, -2, -1, -1, -1, -1
247+0, +0, +0, +0, +0, +0, -2, -1, -1, -1
248+0, +0, +0, +0, +0, +0, +0, -2, -1, -1
249+0, +0, +0, +0, +0, +0, +0, +0, -2, -1
250+0, +0, +0, +0, +0, +0, +0, +0, +0, -2
251
252
253mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = -3, vdia = -4, nrow = 4_IK, ncol = 7_IK, roff = 3_IK, coff = 3_IK)
254where(mat /= -3 .and. mat /= -4); mat = 0; endwhere ! get rid of uninitialized values for better illustration.
255mat
256+0, +0, +0, +0, +0, +0, +0, +0, +0, +0
257+0, +0, +0, +0, +0, +0, +0, +0, +0, +0
258+0, +0, +0, +0, +0, +0, +0, +0, +0, +0
259+0, +0, +0, -4, -3, -3, -3, -3, -3, -3
260+0, +0, +0, +0, -4, -3, -3, -3, -3, -3
261+0, +0, +0, +0, +0, -4, -3, -3, -3, -3
262+0, +0, +0, +0, +0, +0, -4, -3, -3, -3
263+0, +0, +0, +0, +0, +0, +0, +0, +0, +0
264+0, +0, +0, +0, +0, +0, +0, +0, +0, +0
265+0, +0, +0, +0, +0, +0, +0, +0, +0, +0
266
267
268!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
269! Construct upper-diagonal logical matrix.
270!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
271
272
273mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = .true., vdia = .true.)
274mat
275T, T, T, T, T, T, T, T, T, T
276F, T, T, T, T, T, T, T, T, T
277F, F, T, T, T, T, T, T, T, T
278F, F, F, T, T, T, T, T, T, T
279F, F, F, F, T, T, T, T, T, T
280F, F, F, F, F, T, T, T, T, T
281F, F, F, F, F, F, T, T, T, T
282F, F, F, F, F, F, F, T, T, T
283F, F, F, F, F, F, F, F, T, T
284F, F, F, F, F, F, F, F, F, T
285
286
287mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = .true., vdia = .true., nrow = 4_IK, ncol = 7_IK, roff = 3_IK, coff = 3_IK)
288mat
289F, F, F, F, F, F, F, F, F, F
290F, F, F, F, F, F, F, F, F, F
291F, F, F, F, F, F, F, F, F, F
292F, F, F, T, T, T, T, T, T, T
293F, F, F, F, T, T, T, T, T, T
294F, F, F, F, F, T, T, T, T, T
295F, F, F, F, F, F, T, T, T, T
296F, F, F, F, F, F, F, F, F, F
297F, F, F, F, F, F, F, F, F, F
298F, F, F, F, F, F, F, F, F, F
299
300
301!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
302! Construct upper-diagonal real matrix.
303!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
304
305
306mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = -1., vdia = -2.)
307where(mat /= -1 .and. mat /= -2); mat = 0; endwhere ! get rid of uninitialized values for better illustration.
308mat
309-2.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000
310+0.00000000, -2.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000
311+0.00000000, +0.00000000, -2.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000
312+0.00000000, +0.00000000, +0.00000000, -2.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000
313+0.00000000, +0.00000000, +0.00000000, +0.00000000, -2.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000
314+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, -2.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000
315+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, -2.00000000, -1.00000000, -1.00000000, -1.00000000
316+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, -2.00000000, -1.00000000, -1.00000000
317+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, -2.00000000, -1.00000000
318+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, -2.00000000
319
320
321mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppDia, vupp = -3., vdia = -4., nrow = 3_IK, ncol = 7_IK, coff = 3_IK)
322where(mat /= -3. .and. mat /= -4.); mat = 0; endwhere ! get rid of uninitialized values for better illustration.
323mat
324+0.00000000, +0.00000000, +0.00000000, -4.00000000, -3.00000000, -3.00000000, -3.00000000, -3.00000000, -3.00000000, -3.00000000
325+0.00000000, +0.00000000, +0.00000000, +0.00000000, -4.00000000, -3.00000000, -3.00000000, -3.00000000, -3.00000000, -3.00000000
326+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, -4.00000000, -3.00000000, -3.00000000, -3.00000000, -3.00000000
327+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
328+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
329+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
330+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
331+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
332+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
333+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
334
335
336!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
337!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
338! Initialize lower-diagonal matrices.
339!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
340!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
341
342
343!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
344! Construct lower-diagonal string matrix.
345!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
346
347
348mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = 'AA', vdia = 'ZZ', nrow = 7_IK, ncol = 4_IK, roff = 3_IK)
349where(mat /= 'AA' .and. mat /= 'ZZ'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.
350mat
351" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
352" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
353" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
354"ZZ", " ", " ", " ", " ", " ", " ", " ", " ", " "
355"AA", "ZZ", " ", " ", " ", " ", " ", " ", " ", " "
356"AA", "AA", "ZZ", " ", " ", " ", " ", " ", " ", " "
357"AA", "AA", "AA", "ZZ", " ", " ", " ", " ", " ", " "
358"AA", "AA", "AA", "AA", "AA", " ", " ", " ", " ", " "
359"AA", "AA", "AA", "AA", " ", " ", " ", " ", " ", " "
360"AA", "AA", "AA", "AA", " ", " ", " ", " ", " ", " "
361
362
363mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = 'AA', vdia = ['ZZ', 'YY', 'WW', 'XX'], nrow = 7_IK, ncol = 4_IK, roff = 3_IK)
364where(mat /= 'AA' .and. mat /= 'ZZ' .and. mat /= 'YY' .and. mat /= 'WW' .and. mat /= 'XX'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.
365mat
366" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
367" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
368" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
369"ZZ", " ", " ", " ", " ", " ", " ", " ", " ", " "
370"AA", "YY", " ", " ", " ", " ", " ", " ", " ", " "
371"AA", "AA", "WW", " ", " ", " ", " ", " ", " ", " "
372"AA", "AA", "AA", "XX", " ", " ", " ", " ", " ", " "
373"AA", "AA", "AA", "AA", "AA", " ", " ", " ", " ", " "
374"AA", "AA", "AA", "AA", " ", " ", " ", " ", " ", " "
375"AA", "AA", "AA", "AA", " ", " ", " ", " ", " ", " "
376
377
378mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = 'BB', vdia = 'YY', nrow = 4_IK, ncol = 6_IK, roff = 2_IK, coff = 3_IK, doff = +2_IK)
379where(mat /= 'BB' .and. mat /= 'YY'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.
380mat
381" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
382" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
383" ", " ", " ", "BB", "BB", "YY", "BB", "BB", "BB", "BB"
384" ", " ", " ", "BB", "BB", "BB", "YY", "BB", "BB", "BB"
385" ", " ", " ", "BB", "BB", "BB", "BB", "YY", "BB", "BB"
386" ", " ", " ", "BB", "BB", "BB", "BB", "BB", "YY", "BB"
387" ", " ", " ", " ", "BB", "BB", "BB", "BB", "BB", "BB"
388" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
389" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
390" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
391
392
393mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = 'BB', vdia = 'YY', nrow = 4_IK, ncol = 7_IK, roff = 2_IK, coff = 3_IK, doff = +4_IK)
394where(mat /= 'BB' .and. mat /= 'YY'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.
395mat
396" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
397" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
398" ", " ", " ", "BB", "BB", "BB", "BB", "YY", "BB", "BB"
399" ", " ", " ", "BB", "BB", "BB", "BB", "BB", "YY", "BB"
400" ", " ", " ", "BB", "BB", "BB", "BB", "BB", "BB", "YY"
401" ", " ", " ", "BB", "BB", "BB", "BB", "BB", "BB", "BB"
402" ", " ", " ", " ", "BB", "BB", "BB", "BB", "BB", "BB"
403" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
404" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
405" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
406
407
408!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
409! Construct lower-diagonal integer matrix.
410!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
411
412
413mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = -1, vdia = -2)
414where(mat /= -1 .and. mat /= -2); mat = 0; endwhere ! get rid of uninitialized values for better illustration.
415mat
416-2, +0, +0, +0, +0, +0, +0, +0, +0, +0
417-1, -2, +0, +0, +0, +0, +0, +0, +0, +0
418-1, -1, -2, +0, +0, +0, +0, +0, +0, +0
419-1, -1, -1, -2, +0, +0, +0, +0, +0, +0
420-1, -1, -1, -1, -2, +0, +0, +0, +0, +0
421-1, -1, -1, -1, -1, -2, +0, +0, +0, +0
422-1, -1, -1, -1, -1, -1, -2, +0, +0, +0
423-1, -1, -1, -1, -1, -1, -1, -2, +0, +0
424-1, -1, -1, -1, -1, -1, -1, -1, -2, +0
425-1, -1, -1, -1, -1, -1, -1, -1, -1, -2
426
427
428mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = -3, vdia = -4, nrow = 7_IK, ncol = 4_IK, roff = 3_IK, coff = 3_IK)
429where(mat /= -3 .and. mat /= -4); mat = 0; endwhere ! get rid of uninitialized values for better illustration.
430mat
431+0, +0, +0, +0, +0, +0, +0, +0, +0, +0
432+0, +0, +0, +0, +0, +0, +0, +0, +0, +0
433+0, +0, +0, +0, +0, +0, +0, +0, +0, +0
434+0, +0, +0, -4, +0, +0, +0, +0, +0, +0
435+0, +0, +0, -3, -4, +0, +0, +0, +0, +0
436+0, +0, +0, -3, -3, -4, +0, +0, +0, +0
437+0, +0, +0, -3, -3, -3, -4, +0, +0, +0
438+0, +0, +0, -3, -3, -3, -3, +0, +0, +0
439+0, +0, +0, -3, -3, -3, -3, +0, +0, +0
440+0, +0, +0, -3, -3, -3, -3, +0, +0, +0
441
442
443!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
444! Construct lower-diagonal logical matrix.
445!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
446
447
448mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = .true., vdia = .true.)
449mat
450T, F, F, F, F, F, F, F, F, F
451T, T, F, F, F, F, F, F, F, F
452T, T, T, F, F, F, F, F, F, F
453T, T, T, T, F, F, F, F, F, F
454T, T, T, T, T, F, F, F, F, F
455T, T, T, T, T, T, F, F, F, F
456T, T, T, T, T, T, T, F, F, F
457T, T, T, T, T, T, T, T, F, F
458T, T, T, T, T, T, T, T, T, F
459T, T, T, T, T, T, T, T, T, T
460
461
462mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = .true., vdia = .true., nrow = 7_IK, ncol = 4_IK, roff = 3_IK, coff = 3_IK)
463mat
464F, F, F, F, F, F, F, F, F, F
465F, F, F, F, F, F, F, F, F, F
466F, F, F, F, F, F, F, F, F, F
467F, F, F, T, F, F, F, F, F, F
468F, F, F, T, T, F, F, F, F, F
469F, F, F, T, T, T, F, F, F, F
470F, F, F, T, T, T, T, F, F, F
471F, F, F, T, T, T, T, F, F, F
472F, F, F, T, T, T, T, F, F, F
473F, F, F, T, T, T, T, F, F, F
474
475
476!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
477! Construct lower-diagonal real matrix.
478!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
479
480
481mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = -1., vdia = -2.)
482where(mat /= -1 .and. mat /= -2); mat = 0; endwhere ! get rid of uninitialized values for better illustration.
483mat
484-2.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
485-1.00000000, -2.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
486-1.00000000, -1.00000000, -2.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
487-1.00000000, -1.00000000, -1.00000000, -2.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
488-1.00000000, -1.00000000, -1.00000000, -1.00000000, -2.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
489-1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -2.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
490-1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -2.00000000, +0.00000000, +0.00000000, +0.00000000
491-1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -2.00000000, +0.00000000, +0.00000000
492-1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -2.00000000, +0.00000000
493-1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -2.00000000
494
495
496mat = getMatInit(Shape = [10_IK, 10_IK], subset = lowDia, vlow = -3., vdia = -4., nrow = 7_IK, ncol = 3_IK, coff = 3_IK)
497where(mat /= -3. .and. mat /= -4.); mat = 0; endwhere ! get rid of uninitialized values for better illustration.
498mat
499+0.00000000, +0.00000000, +0.00000000, -4.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
500+0.00000000, +0.00000000, +0.00000000, -3.00000000, -4.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
501+0.00000000, +0.00000000, +0.00000000, -3.00000000, -3.00000000, -4.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
502+0.00000000, +0.00000000, +0.00000000, -3.00000000, -3.00000000, -3.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
503+0.00000000, +0.00000000, +0.00000000, -3.00000000, -3.00000000, -3.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
504+0.00000000, +0.00000000, +0.00000000, -3.00000000, -3.00000000, -3.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
505+0.00000000, +0.00000000, +0.00000000, -3.00000000, -3.00000000, -3.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
506+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
507+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
508+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
509
510
511!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
512!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
513! Initialize upper-lower matrices.
514!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
515!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
516
517
518!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
519! Construct upper-lower string matrix.
520!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
521
522
523mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = 'AA', vlow = 'HH', nrow = 7_IK, ncol = 4_IK, roff = 3_IK)
524where(mat /= 'AA' .and. mat /= 'HH' .and. mat /= 'OO'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.
525mat
526" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
527" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
528" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
529" ", "AA", "AA", "AA", " ", " ", " ", " ", " ", " "
530"HH", " ", "AA", "AA", " ", " ", " ", " ", " ", " "
531"HH", "HH", " ", "AA", " ", " ", " ", " ", " ", " "
532"HH", "HH", "HH", " ", " ", " ", " ", " ", " ", " "
533"HH", "HH", "HH", "HH", "AA", " ", " ", "OO", " ", " "
534"HH", "HH", "HH", "HH", " ", " ", " ", " ", "OO", " "
535"HH", "HH", "HH", "HH", " ", " ", " ", " ", " ", "OO"
536
537
538mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = 'AA', vlow = 'HH', nrow = 7_IK, ncol = 4_IK, roff = 3_IK)
539where(mat /= 'AA' .and. mat /= 'HH'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.
540mat
541" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
542" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
543" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
544" ", "AA", "AA", "AA", " ", " ", " ", " ", " ", " "
545"HH", " ", "AA", "AA", " ", " ", " ", " ", " ", " "
546"HH", "HH", " ", "AA", " ", " ", " ", " ", " ", " "
547"HH", "HH", "HH", " ", " ", " ", " ", " ", " ", " "
548"HH", "HH", "HH", "HH", "AA", " ", " ", " ", " ", " "
549"HH", "HH", "HH", "HH", " ", " ", " ", " ", " ", " "
550"HH", "HH", "HH", "HH", " ", " ", " ", " ", " ", " "
551
552
553mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = 'AA', vlow = 'HH', nrow = 4_IK, ncol = 7_IK, roff = 3_IK)
554where(mat /= 'AA' .and. mat /= 'HH'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.
555mat
556" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
557" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
558" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
559" ", "AA", "AA", "AA", "AA", "AA", "AA", " ", " ", " "
560"HH", " ", "AA", "AA", "AA", "AA", "AA", " ", " ", " "
561"HH", "HH", " ", "AA", "AA", "AA", "AA", " ", " ", " "
562"HH", "HH", "HH", " ", "AA", "AA", "AA", " ", " ", " "
563" ", " ", " ", " ", "AA", " ", " ", " ", " ", " "
564" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
565" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
566
567
568mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = 'AA', vlow = 'BB', nrow = 4_IK, ncol = 6_IK, roff = 2_IK, coff = 3_IK, doff = +2_IK)
569where(mat /= 'AA' .and. mat /= 'BB' .and. mat /= 'YY'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.
570mat
571" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
572" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
573" ", " ", " ", "BB", "BB", "BB", "AA", "AA", "AA", "BB"
574" ", " ", " ", "BB", "BB", "BB", "AA", "AA", "AA", "BB"
575" ", " ", " ", "BB", "BB", "BB", "BB", "BB", "AA", "YY"
576" ", " ", " ", "BB", "BB", "BB", "BB", "BB", "BB", "BB"
577" ", " ", " ", " ", "AA", "AA", "AA", "BB", "BB", "BB"
578" ", " ", " ", " ", "AA", " ", " ", " ", " ", " "
579" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
580" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
581
582
583mat = getMatInit(vupp = 'AA', vlow = 'BB', Shape = [10_IK, 10_IK], subset = uppLow, nrow = 4_IK, ncol = 7_IK, roff = 2_IK, coff = 3_IK, doff = +4_IK)
584where(mat /= 'AA' .and. mat /= 'BB' .and. mat /= 'YY'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.
585mat
586" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
587" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
588" ", " ", " ", "BB", "BB", "BB", "BB", "AA", "AA", "AA"
589" ", " ", " ", "BB", "BB", "BB", "BB", "BB", "AA", "AA"
590" ", " ", " ", "BB", "BB", "BB", "BB", "BB", "BB", "YY"
591" ", " ", " ", "BB", "BB", "BB", "BB", "BB", "BB", "BB"
592" ", " ", " ", " ", "AA", "AA", "AA", "BB", "BB", "BB"
593" ", " ", " ", " ", "AA", " ", " ", " ", " ", " "
594" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
595" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
596
597
598!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
599! Construct upper-lower integer matrix.
600!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
601
602
603mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = +1, vlow = -1)
604where(mat /= +1 .and. mat /= -1 .and. mat /= -2); mat = 0; endwhere ! get rid of uninitialized values for better illustration.
605mat
606+0, +1, +1, +1, +1, +1, +1, +1, +1, +1
607-1, +0, +1, +1, +1, +1, +1, +1, +1, +1
608-1, -1, +0, +1, +1, +1, +1, +1, +1, +1
609-1, -1, -1, +0, +1, +1, +1, +1, +1, +1
610-1, -1, -1, -1, +0, +1, +1, +1, +1, +1
611-1, -1, -1, -1, -1, +0, +1, +1, +1, +1
612-1, -1, -1, -1, -1, -1, +0, +1, +1, +1
613-1, -1, -1, -1, -1, -1, -1, +0, +1, +1
614-1, -1, -1, -1, -1, -1, -1, -1, +0, +1
615-1, -1, -1, -1, -1, -1, -1, -1, -1, +0
616
617
618mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = +1, vlow = -3, nrow = 7_IK, ncol = 4_IK, roff = 3_IK, coff = 3_IK)
619where(mat /= +1 .and. mat /= -3 .and. mat /= -4); mat = 0; endwhere ! get rid of uninitialized values for better illustration.
620mat
621+0, +0, +0, +0, +0, +0, +0, +0, +0, +0
622+0, +0, +0, +0, +0, +0, +0, +0, +0, +0
623+0, +0, +0, +0, +0, +0, +0, +0, +0, +0
624+0, +0, +0, +0, +1, +1, +1, +0, +0, +0
625+0, +0, +0, -3, +0, +1, +1, +0, +0, +0
626+0, +0, +0, -3, -3, +0, +1, +0, +0, +0
627+0, +0, +0, -3, -3, -3, +0, +0, +0, +0
628+0, +0, +0, -3, -3, -3, -3, +0, +0, +0
629+0, +0, +0, -3, -3, -3, -3, +0, +0, +0
630+0, +0, +0, -3, -3, -3, -3, +0, +0, +0
631
632
633!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
634! Construct upper-lower logical matrix.
635!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
636
637
638mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = .true., vlow = .true.)
639mat
640T, T, T, T, T, T, T, T, T, T
641T, T, T, T, T, T, T, T, T, T
642T, T, F, T, T, T, T, T, T, T
643T, T, T, F, T, T, T, T, T, T
644T, T, T, T, T, T, T, T, T, T
645T, T, T, T, T, F, T, T, T, T
646T, T, T, T, T, T, T, T, T, T
647T, T, T, T, T, T, T, T, T, T
648T, T, T, T, T, T, T, T, T, T
649T, T, T, T, T, T, T, T, T, T
650
651
652mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = .true., vlow = .true., nrow = 7_IK, ncol = 4_IK, roff = 3_IK, coff = 3_IK, doff = +1_IK)
653mat
654T, T, T, T, T, T, T, T, T, T
655T, T, T, T, T, T, T, T, T, T
656T, T, F, T, T, T, T, T, T, T
657T, T, T, T, T, T, T, T, T, T
658T, T, T, T, T, T, T, T, T, T
659T, T, T, T, T, T, T, T, T, T
660T, T, T, T, T, T, T, T, T, T
661T, T, T, T, T, T, T, T, T, T
662T, T, T, T, T, T, T, T, T, T
663T, T, T, T, T, T, T, T, T, T
664
665
666!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
667! Construct upper-lower real matrix.
668!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
669
670
671mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = +1., vlow = -1.)
672where(mat /= +1 .and. mat /= -1 .and. mat /= -2); mat = 0; endwhere ! get rid of uninitialized values for better illustration.
673mat
674+0.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000
675-1.00000000, +0.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000
676-1.00000000, -1.00000000, +0.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000
677-1.00000000, -1.00000000, -1.00000000, +0.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000
678-1.00000000, -1.00000000, -1.00000000, -1.00000000, +0.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000
679-1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, +0.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000
680-1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, +0.00000000, +1.00000000, +1.00000000, +1.00000000
681-1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, +0.00000000, +1.00000000, +1.00000000
682-1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, +0.00000000, +1.00000000
683-1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, +0.00000000
684
685
686mat = getMatInit(Shape = [10_IK, 10_IK], subset = uppLow, vupp = +3., vlow = -3., nrow = 7_IK, ncol = 3_IK, coff = 3_IK)
687where(mat /= +3. .and. mat /= -3. .and. mat /= -4.); mat = 0; endwhere ! get rid of uninitialized values for better illustration.
688mat
689+0.00000000, +0.00000000, +0.00000000, +0.00000000, +3.00000000, +3.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
690+0.00000000, +0.00000000, +0.00000000, -3.00000000, +0.00000000, +3.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
691+0.00000000, +0.00000000, +0.00000000, -3.00000000, -3.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
692+0.00000000, +0.00000000, +0.00000000, -3.00000000, -3.00000000, -3.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
693+0.00000000, +0.00000000, +0.00000000, -3.00000000, -3.00000000, -3.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
694+0.00000000, +0.00000000, +0.00000000, -3.00000000, -3.00000000, -3.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
695+0.00000000, +0.00000000, +0.00000000, -3.00000000, -3.00000000, -3.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
696+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
697+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
698+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
699
700
701!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
702!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
703! Initialize upper-lower-diagonal matrices.
704!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
705!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
706
707
708!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
709! Construct upper-lower-diagonal string matrix.
710!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
711
712
713mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = 'AA', vlow = 'HH', vdia = 'OO', nrow = 7_IK, ncol = 4_IK, roff = 3_IK)
714where(mat /= 'AA' .and. mat /= 'HH' .and. mat /= 'OO'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.
715mat
716" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
717" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
718" ", " ", " ", " ", " ", " ", " ", "AA", "AA", "AA"
719"OO", "AA", "AA", "AA", " ", " ", " ", " ", "AA", "AA"
720"HH", "OO", "AA", "AA", " ", " ", " ", " ", " ", " "
721"HH", "HH", "OO", "AA", " ", " ", " ", " ", " ", " "
722"HH", "HH", "HH", "OO", "AA", "AA", "AA", " ", " ", " "
723"HH", "HH", "HH", "HH", "AA", " ", " ", "OO", " ", " "
724"HH", "HH", "HH", "HH", " ", " ", " ", " ", "OO", " "
725"HH", "HH", "HH", "HH", " ", " ", " ", " ", " ", "OO"
726
727
728mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = 'AA', vlow = 'HH', vdia = ['OO', 'YY', 'WW', 'XX'], nrow = 7_IK, ncol = 4_IK, roff = 3_IK)
729where(mat /= 'AA' .and. mat /= 'HH' .and. mat /= 'OO' .and. mat /= 'YY' .and. mat /= 'WW' .and. mat /= 'XX'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.
730mat
731" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
732" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
733" ", " ", " ", " ", " ", " ", " ", "AA", "AA", "AA"
734"OO", "AA", "AA", "AA", " ", " ", " ", " ", "AA", "AA"
735"HH", "YY", "AA", "AA", " ", " ", " ", " ", " ", "YY"
736"HH", "HH", "WW", "AA", " ", " ", " ", " ", " ", " "
737"HH", "HH", "HH", "XX", "AA", "AA", "AA", " ", " ", " "
738"HH", "HH", "HH", "HH", "AA", " ", " ", "OO", " ", " "
739"HH", "HH", "HH", "HH", " ", " ", " ", " ", "OO", " "
740"HH", "HH", "HH", "HH", " ", " ", " ", " ", " ", "OO"
741
742
743mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = 'AA', vlow = 'HH', vdia = ['OO', 'YY', 'WW', 'XX'], nrow = 4_IK, ncol = 7_IK, roff = 3_IK)
744where(mat /= 'AA' .and. mat /= 'HH' .and. mat /= 'OO' .and. mat /= 'YY' .and. mat /= 'WW' .and. mat /= 'XX'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.
745mat
746" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
747" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
748" ", " ", " ", " ", " ", " ", " ", "AA", "AA", "AA"
749"OO", "AA", "AA", "AA", "AA", "AA", "AA", " ", "AA", "AA"
750"HH", "YY", "AA", "AA", "AA", "AA", "AA", " ", " ", "YY"
751"HH", "HH", "WW", "AA", "AA", "AA", "AA", " ", " ", " "
752"HH", "HH", "HH", "XX", "AA", "AA", "AA", " ", " ", " "
753" ", " ", " ", " ", "AA", " ", " ", "OO", " ", " "
754" ", " ", " ", " ", " ", " ", " ", " ", "OO", " "
755" ", " ", " ", " ", " ", " ", " ", " ", " ", "OO"
756
757
758mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = 'AA', vlow = 'BB', vdia = 'YY', nrow = 4_IK, ncol = 6_IK, roff = 2_IK, coff = 3_IK, doff = +2_IK)
759where(mat /= 'AA' .and. mat /= 'BB' .and. mat /= 'YY'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.
760mat
761" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
762" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
763" ", " ", " ", "BB", "BB", "YY", "AA", "AA", "AA", "AA"
764" ", " ", " ", "BB", "BB", "BB", "YY", "AA", "AA", "AA"
765" ", " ", " ", "BB", "BB", "BB", "BB", "YY", "AA", "YY"
766" ", " ", " ", "BB", "BB", "BB", "BB", "BB", "YY", "BB"
767" ", " ", " ", " ", "AA", "AA", "AA", "BB", "BB", "BB"
768" ", " ", " ", " ", "AA", " ", " ", " ", " ", " "
769" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
770" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
771
772
773mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = 'AA', vlow = 'BB', vdia = 'YY', nrow = 4_IK, ncol = 7_IK, roff = 2_IK, coff = 3_IK, doff = +4_IK)
774where(mat /= 'AA' .and. mat /= 'BB' .and. mat /= 'YY'); mat = ' '; endwhere ! get rid of uninitialized values for better illustration.
775mat
776" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
777" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
778" ", " ", " ", "BB", "BB", "BB", "BB", "YY", "AA", "AA"
779" ", " ", " ", "BB", "BB", "BB", "BB", "BB", "YY", "AA"
780" ", " ", " ", "BB", "BB", "BB", "BB", "BB", "BB", "YY"
781" ", " ", " ", "BB", "BB", "BB", "BB", "BB", "BB", "BB"
782" ", " ", " ", " ", "AA", "AA", "AA", "BB", "BB", "BB"
783" ", " ", " ", " ", "AA", " ", " ", " ", " ", " "
784" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
785" ", " ", " ", " ", " ", " ", " ", " ", " ", " "
786
787
788!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
789! Construct upper-lower-diagonal integer matrix.
790!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
791
792
793mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = +1, vlow = -1, vdia = -2)
794where(mat /= +1 .and. mat /= -1 .and. mat /= -2); mat = 0; endwhere ! get rid of uninitialized values for better illustration.
795mat
796-2, +1, +1, +1, +1, +1, +1, +1, +1, +1
797-1, -2, +1, +1, +1, +1, +1, +1, +1, +1
798-1, -1, -2, +1, +1, +1, +1, +1, +1, +1
799-1, -1, -1, -2, +1, +1, +1, +1, +1, +1
800-1, -1, -1, -1, -2, +1, +1, +1, +1, +1
801-1, -1, -1, -1, -1, -2, +1, +1, +1, +1
802-1, -1, -1, -1, -1, -1, -2, +1, +1, +1
803-1, -1, -1, -1, -1, -1, -1, -2, +1, +1
804-1, -1, -1, -1, -1, -1, -1, -1, -2, +1
805-1, -1, -1, -1, -1, -1, -1, -1, -1, -2
806
807
808mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = +1, vlow = -3, vdia = -4, nrow = 7_IK, ncol = 4_IK, roff = 3_IK, coff = 3_IK)
809where(mat /= +1 .and. mat /= -3 .and. mat /= -4); mat = 0; endwhere ! get rid of uninitialized values for better illustration.
810mat
811+0, +0, +0, +0, +0, +0, +0, +0, +0, +0
812+0, +0, +0, +0, +0, +0, +0, +0, +0, +0
813+0, +0, +0, +0, +0, +0, +0, +0, +0, +0
814+0, +0, +0, -4, +1, +1, +1, +0, +0, +0
815+0, +0, +0, -3, -4, +1, +1, +0, +0, +0
816+0, +0, +0, -3, -3, -4, +1, +0, +0, +0
817+0, +0, +0, -3, -3, -3, -4, +0, +0, +0
818+0, +0, +0, -3, -3, -3, -3, +0, +0, +0
819+0, +0, +0, -3, -3, -3, -3, +0, +0, +0
820+0, +0, +0, -3, -3, -3, -3, +0, +0, +0
821
822
823!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
824! Construct upper-lower-diagonal logical matrix.
825!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
826
827
828mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = .true., vlow = .true., vdia = .false.)
829mat
830F, T, T, T, T, T, T, T, T, T
831T, F, T, T, T, T, T, T, T, T
832T, T, F, T, T, T, T, T, T, T
833T, T, T, F, T, T, T, T, T, T
834T, T, T, T, F, T, T, T, T, T
835T, T, T, T, T, F, T, T, T, T
836T, T, T, T, T, T, F, T, T, T
837T, T, T, T, T, T, T, F, T, T
838T, T, T, T, T, T, T, T, F, T
839T, T, T, T, T, T, T, T, T, F
840
841
842mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = .true., vlow = .true., vdia = .false., nrow = 7_IK, ncol = 4_IK, roff = 3_IK, coff = 3_IK, doff = +1_IK)
843mat
844F, T, T, T, T, T, T, T, T, T
845T, F, T, T, T, T, T, T, T, T
846T, T, F, T, T, T, T, T, T, T
847T, T, T, T, F, T, T, T, T, T
848T, T, T, T, T, F, T, T, T, T
849T, T, T, T, T, T, F, T, T, T
850T, T, T, T, T, T, T, T, T, T
851T, T, T, T, T, T, T, F, T, T
852T, T, T, T, T, T, T, T, F, T
853T, T, T, T, T, T, T, T, T, F
854
855
856!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
857! Construct upper-lower-diagonal real matrix.
858!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
859
860
861mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = +1., vlow = -1., vdia = -2.)
862where(mat /= +1 .and. mat /= -1 .and. mat /= -2); mat = 0; endwhere ! get rid of uninitialized values for better illustration.
863mat
864-2.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000
865-1.00000000, -2.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000
866-1.00000000, -1.00000000, -2.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000
867-1.00000000, -1.00000000, -1.00000000, -2.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000
868-1.00000000, -1.00000000, -1.00000000, -1.00000000, -2.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000
869-1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -2.00000000, +1.00000000, +1.00000000, +1.00000000, +1.00000000
870-1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -2.00000000, +1.00000000, +1.00000000, +1.00000000
871-1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -2.00000000, +1.00000000, +1.00000000
872-1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -2.00000000, +1.00000000
873-1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -1.00000000, -2.00000000
874
875
876mat = getMatInit([10_IK, 10_IK], uppLowDia, vupp = +3., vlow = -3., vdia = -4., nrow = 7_IK, ncol = 3_IK, coff = 3_IK)
877where(mat /= +3. .and. mat /= -3. .and. mat /= -4.); mat = 0; endwhere ! get rid of uninitialized values for better illustration.
878mat
879+0.00000000, +0.00000000, +0.00000000, -4.00000000, +3.00000000, +3.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
880+0.00000000, +0.00000000, +0.00000000, -3.00000000, -4.00000000, +3.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
881+0.00000000, +0.00000000, +0.00000000, -3.00000000, -3.00000000, -4.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
882+0.00000000, +0.00000000, +0.00000000, -3.00000000, -3.00000000, -3.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
883+0.00000000, +0.00000000, +0.00000000, -3.00000000, -3.00000000, -3.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
884+0.00000000, +0.00000000, +0.00000000, -3.00000000, -3.00000000, -3.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
885+0.00000000, +0.00000000, +0.00000000, -3.00000000, -3.00000000, -3.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
886+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
887+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
888+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
889
890
Test:
test_pm_matrixInit
Internal naming convention:
The following illustrates the internal naming convention used for the procedures within this generic interface.
getMatInitXXD_D2XX0_SK5()
||| ||||| |||
||| ||||| Output type/kind
||| ||||The rank of `vdia` argument. `X` implies missing argument.
||| |||The rank of `vlow` argument. `X` implies missing argument.
||| ||The rank of `vupp` argument. `X` implies missing argument.
||| |The rank of the output `mat`.
||| `D` stands for dimension.
||The `vdia` argument presence: `D` implies presence. `X` implies missing argument.
|The `vlow` argument presence: `L` implies presence. `X` implies missing argument.
The `vupp` argument presence: `U` implies presence. `X` implies missing argument.
Bug:

Status: Unresolved
Source: Intel Classic Fortran Compiler ifort version 2021.8.0 20221119
Description: Intel Classic Fortran Compiler ifort thows a strange error for compiling this module when dia constant is used within this module.
Apparently, Intel Classic Fortran Compiler ifort confuses this imported constant with the dummy procedure argument of the same type in procedures named getMatInitXXD_D2XX0* or getMatInitXXD_D2XX1.

Remedy (as of ParaMonte Library version 2.0.0): There does not appear to exist any way to resolve this intel error currently, other than not importing the constant dia in this module.
Instead, one has to use dia_type() to call the procedures when resolving the generic interface.
As such, the argument names upp, low, dia were renamed to vupp, vlow, vdia to emphasize the value attribute of these arguments and to resolve the Intel name-conflict bug.
Todo:
Critical Priority: This generic interface should be extended to matrices of different packing formats besides the default.


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 289 of file pm_matrixInit.F90.


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