Line data Source code
1 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3 : !!!! !!!!
4 : !!!! ParaMonte: Parallel Monte Carlo and Machine Learning Library. !!!!
5 : !!!! !!!!
6 : !!!! Copyright (C) 2012-present, The Computational Data Science Lab !!!!
7 : !!!! !!!!
8 : !!!! This file is part of the ParaMonte library. !!!!
9 : !!!! !!!!
10 : !!!! LICENSE !!!!
11 : !!!! !!!!
12 : !!!! https://github.com/cdslaborg/paramonte/blob/main/LICENSE.md !!!!
13 : !!!! !!!!
14 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
15 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
16 :
17 : !> \brief
18 : !> This include file contains procedure implementations of the tests of [pm_matrixMulTri](@ref pm_matrixMulTri).
19 : !>
20 : !> \fintest
21 : !>
22 : !> \author
23 : !> \FatemehBagheri, 12:27 AM Tuesday, February 22, 2022, Dallas, TX
24 :
25 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
26 :
27 : #if CK_ENABLED
28 : #define TYPE_KIND complex(TKC)
29 : complex(TKC), parameter :: lb = (1._TKC, -2._TKC), ub = (2._TKC, -1._TKC), ONE = (1._TKC, 0._TKC), ZERO = (0._TKC, 0._TKC)
30 : #elif RK_ENABLED
31 : #define TYPE_KIND real(TKC)
32 : real(TKC), parameter :: lb = 1._TKC, ub = 2._TKC, ONE = 1._TKC, ZERO = 0._TKC
33 : #else
34 : #error "Unrecognized interface."
35 : #endif
36 : logical(LK) :: simplified
37 : integer(IK) :: itry, begG, incG, endG, irow, nrow, ncol, ndim, roffT, coffT, roffG, coffG
38 : integer(IK), parameter :: increments(*) = [(irow, irow = -3, -1), (irow, irow = 1, 3)]
39 : TYPE_KIND, allocatable :: solmat(:,:), trimat(:,:), genmat(:,:), genvec(:), dumvec(:), choices(:)
40 : real(TKC), parameter :: TOL = epsilon(0._TKC)**.5
41 : TYPE_KIND :: alpha, alinv
42 :
43 8 : assertion = .true.
44 :
45 808 : do itry = 1, 100
46 :
47 : ! Build the matrices.
48 :
49 800 : incG = getChoice(increments)
50 800 : nrow = getUnifRand(0_IK, 7_IK)
51 800 : ncol = getUnifRand(0_IK, 7_IK)
52 800 : roffT = getUnifRand(0_IK, 3_IK)
53 800 : coffT = getUnifRand(0_IK, 3_IK)
54 800 : roffG = getUnifRand(0_IK, 3_IK)
55 800 : coffG = getUnifRand(0_IK, 3_IK)
56 4000 : choices = [ZERO, ONE, getUnifRand(lb, ub)]
57 800 : alpha = getChoice(choices)
58 800 : if (alpha == ZERO) then
59 263 : alinv = alpha
60 : else
61 537 : alinv = 1 / alpha
62 : end if
63 :
64 800 : ndim = nrow
65 11234 : genvec = getChoice([-1, 1]) * getUnifRand(lb, ub, max(0, 1 + (nrow - 1) * abs(incG)))
66 79638 : genmat = getChoice([-1, 1]) * getUnifRand(lb, ub, 2 * roffG + nrow, 2 * coffG + ncol)
67 83444 : trimat = getChoice([-1, 1]) * getUnifRand(lb, ub, 2 * roffT + ndim, 2 * coffT + ndim)
68 800 : if (incG < 0) then
69 393 : begG = size(genvec)
70 393 : endG = 1
71 : else
72 407 : begG = 1
73 407 : endG = size(genvec)
74 : end if
75 :
76 : ! BLAS - LEVEL 2: ?TRMV / ?TRSV - simplified interface
77 :
78 800 : simplified = .true.
79 :
80 : ! upperDiag, nothing/inversion
81 7036 : dumvec = genvec
82 22970 : call setMatMulTri(trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim), upperDiag, nothing, dumvec(begG:endG:incG))
83 22970 : call setMatMulTri(trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim), upperDiag, inversion, dumvec(begG:endG:incG))
84 800 : call reportTRXV(simplified, "nothing/inversion")
85 :
86 : ! upperDiag, transSymm/transOrth
87 6417 : dumvec = genvec
88 22970 : call setMatMulTri(trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim), upperDiag, transSymm, dumvec(begG:endG:incG))
89 22970 : call setMatMulTri(trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim), upperDiag, transOrth, dumvec(begG:endG:incG))
90 800 : call reportTRXV(simplified, "transSymm/transOrth")
91 :
92 : ! upperDiag, transHerm/transUnit
93 6417 : dumvec = genvec
94 22970 : call setMatMulTri(trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim), upperDiag, transHerm, dumvec(begG:endG:incG))
95 22970 : call setMatMulTri(trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim), upperDiag, transUnit, dumvec(begG:endG:incG))
96 800 : call reportTRXV(simplified, "transHerm/transUnit")
97 :
98 : ! lowerDiag, nothing/inversion
99 6417 : dumvec = genvec
100 22970 : call setMatMulTri(trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim), lowerDiag, nothing, dumvec(begG:endG:incG))
101 22970 : call setMatMulTri(trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim), lowerDiag, inversion, dumvec(begG:endG:incG))
102 800 : call reportTRXV(simplified, "nothing/inversion")
103 :
104 : ! lowerDiag, transSymm/transOrth
105 6417 : dumvec = genvec
106 22970 : call setMatMulTri(trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim), lowerDiag, transSymm, dumvec(begG:endG:incG))
107 22970 : call setMatMulTri(trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim), lowerDiag, transOrth, dumvec(begG:endG:incG))
108 800 : call reportTRXV(simplified, "transSymm/transOrth")
109 :
110 : ! lowerDiag, transHerm/transUnit
111 6417 : dumvec = genvec
112 22970 : call setMatMulTri(trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim), lowerDiag, transHerm, dumvec(begG:endG:incG))
113 22970 : call setMatMulTri(trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim), lowerDiag, transUnit, dumvec(begG:endG:incG))
114 800 : call reportTRXV(simplified, "transHerm/transUnit")
115 :
116 : ! upperUnit, nothing/inversion
117 6417 : dumvec = genvec
118 22970 : call setMatMulTri(trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim), upperUnit, nothing, dumvec(begG:endG:incG))
119 22970 : call setMatMulTri(trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim), upperUnit, inversion, dumvec(begG:endG:incG))
120 800 : call reportTRXV(simplified, "nothing/inversion")
121 :
122 : ! upperUnit, transSymm/transOrth
123 6417 : dumvec = genvec
124 22970 : call setMatMulTri(trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim), upperUnit, transSymm, dumvec(begG:endG:incG))
125 22970 : call setMatMulTri(trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim), upperUnit, transOrth, dumvec(begG:endG:incG))
126 800 : call reportTRXV(simplified, "transSymm/transOrth")
127 :
128 : ! upperUnit, transHerm/transUnit
129 6417 : dumvec = genvec
130 22970 : call setMatMulTri(trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim), upperUnit, transHerm, dumvec(begG:endG:incG))
131 22970 : call setMatMulTri(trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim), upperUnit, transUnit, dumvec(begG:endG:incG))
132 800 : call reportTRXV(simplified, "transHerm/transUnit")
133 :
134 : ! lowerUnit, nothing/inversion
135 6417 : dumvec = genvec
136 22970 : call setMatMulTri(trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim), lowerUnit, nothing, dumvec(begG:endG:incG))
137 22970 : call setMatMulTri(trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim), lowerUnit, inversion, dumvec(begG:endG:incG))
138 800 : call reportTRXV(simplified, "nothing/inversion")
139 :
140 : ! lowerUnit, transSymm/transOrth
141 6417 : dumvec = genvec
142 22970 : call setMatMulTri(trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim), lowerUnit, transSymm, dumvec(begG:endG:incG))
143 22970 : call setMatMulTri(trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim), lowerUnit, transOrth, dumvec(begG:endG:incG))
144 800 : call reportTRXV(simplified, "transSymm/transOrth")
145 :
146 : ! lowerUnit, transHerm/transUnit
147 6417 : dumvec = genvec
148 22970 : call setMatMulTri(trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim), lowerUnit, transHerm, dumvec(begG:endG:incG))
149 22970 : call setMatMulTri(trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim), lowerUnit, transUnit, dumvec(begG:endG:incG))
150 800 : call reportTRXV(simplified, "transHerm/transUnit")
151 :
152 : ! BLAS - LEVEL 2: ?TRMV / ?TRSV - contiguous interface
153 :
154 800 : simplified = .false.
155 :
156 : ! upperDiag, nothing/inversion
157 6417 : dumvec = genvec
158 800 : call setMatMulTri(trimat, upperDiag, nothing, dumvec, ndim, roffT, coffT, incG)
159 800 : call setMatMulTri(trimat, upperDiag, inversion, dumvec, ndim, roffT, coffT, incG)
160 800 : call reportTRXV(simplified, "nothing/inversion")
161 :
162 : ! upperDiag, transSymm/transOrth
163 6417 : dumvec = genvec
164 800 : call setMatMulTri(trimat, upperDiag, transSymm, dumvec, ndim, roffT, coffT, incG)
165 800 : call setMatMulTri(trimat, upperDiag, transOrth, dumvec, ndim, roffT, coffT, incG)
166 800 : call reportTRXV(simplified, "transSymm/transOrth")
167 :
168 : ! upperDiag, transHerm/transUnit
169 6417 : dumvec = genvec
170 800 : call setMatMulTri(trimat, upperDiag, transHerm, dumvec, ndim, roffT, coffT, incG)
171 800 : call setMatMulTri(trimat, upperDiag, transUnit, dumvec, ndim, roffT, coffT, incG)
172 800 : call reportTRXV(simplified, "transHerm/transUnit")
173 :
174 : ! lowerDiag, nothing/inversion
175 6417 : dumvec = genvec
176 800 : call setMatMulTri(trimat, lowerDiag, nothing, dumvec, ndim, roffT, coffT, incG)
177 800 : call setMatMulTri(trimat, lowerDiag, inversion, dumvec, ndim, roffT, coffT, incG)
178 800 : call reportTRXV(simplified, "nothing/inversion")
179 :
180 : ! lowerDiag, transSymm/transOrth
181 6417 : dumvec = genvec
182 800 : call setMatMulTri(trimat, lowerDiag, transSymm, dumvec, ndim, roffT, coffT, incG)
183 800 : call setMatMulTri(trimat, lowerDiag, transOrth, dumvec, ndim, roffT, coffT, incG)
184 800 : call reportTRXV(simplified, "transSymm/transOrth")
185 :
186 : ! lowerDiag, transHerm/transUnit
187 6417 : dumvec = genvec
188 800 : call setMatMulTri(trimat, lowerDiag, transHerm, dumvec, ndim, roffT, coffT, incG)
189 800 : call setMatMulTri(trimat, lowerDiag, transUnit, dumvec, ndim, roffT, coffT, incG)
190 800 : call reportTRXV(simplified, "transHerm/transUnit")
191 :
192 : ! upperUnit, nothing/inversion
193 6417 : dumvec = genvec
194 800 : call setMatMulTri(trimat, upperUnit, nothing, dumvec, ndim, roffT, coffT, incG)
195 800 : call setMatMulTri(trimat, upperUnit, inversion, dumvec, ndim, roffT, coffT, incG)
196 800 : call reportTRXV(simplified, "nothing/inversion")
197 :
198 : ! upperUnit, transSymm/transOrth
199 6417 : dumvec = genvec
200 800 : call setMatMulTri(trimat, upperUnit, transSymm, dumvec, ndim, roffT, coffT, incG)
201 800 : call setMatMulTri(trimat, upperUnit, transOrth, dumvec, ndim, roffT, coffT, incG)
202 800 : call reportTRXV(simplified, "transSymm/transOrth")
203 :
204 : ! upperUnit, transHerm/transUnit
205 6417 : dumvec = genvec
206 800 : call setMatMulTri(trimat, upperUnit, transHerm, dumvec, ndim, roffT, coffT, incG)
207 800 : call setMatMulTri(trimat, upperUnit, transUnit, dumvec, ndim, roffT, coffT, incG)
208 800 : call reportTRXV(simplified, "transHerm/transUnit")
209 :
210 : ! lowerUnit, nothing/inversion
211 6417 : dumvec = genvec
212 800 : call setMatMulTri(trimat, lowerUnit, nothing, dumvec, ndim, roffT, coffT, incG)
213 800 : call setMatMulTri(trimat, lowerUnit, inversion, dumvec, ndim, roffT, coffT, incG)
214 800 : call reportTRXV(simplified, "nothing/inversion")
215 :
216 : ! lowerUnit, transSymm/transOrth
217 6417 : dumvec = genvec
218 800 : call setMatMulTri(trimat, lowerUnit, transSymm, dumvec, ndim, roffT, coffT, incG)
219 800 : call setMatMulTri(trimat, lowerUnit, transOrth, dumvec, ndim, roffT, coffT, incG)
220 800 : call reportTRXV(simplified, "transSymm/transOrth")
221 :
222 : ! lowerUnit, transHerm/transUnit
223 6417 : dumvec = genvec
224 800 : call setMatMulTri(trimat, lowerUnit, transHerm, dumvec, ndim, roffT, coffT, incG)
225 800 : call setMatMulTri(trimat, lowerUnit, transUnit, dumvec, ndim, roffT, coffT, incG)
226 800 : call reportTRXV(simplified, "transHerm/transUnit")
227 :
228 : ! BLAS - LEVEL 3: ?TRMM / ?TRSM
229 :
230 800 : call testTRXM()
231 808 : call testTRXM(alpha, alinv)
232 :
233 : end do
234 :
235 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
236 :
237 : contains
238 :
239 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
240 :
241 19200 : subroutine reportTRXV(simplified, operationA)
242 : character(*), intent(in) :: operationA
243 : logical(LK) , intent(in) :: simplified
244 : logical(LK), allocatable :: tolerable(:)
245 19200 : type(display_type) :: disp
246 269616 : tolerable = isClose(genvec, dumvec, abstol = TOL)
247 134808 : assertion = assertion .and. all(tolerable)
248 19200 : if (test%traceable .and. .not. assertion) then
249 : ! LCOV_EXCL_START
250 : call disp%show("TRMV/TRSV")
251 : call disp%show("simplified")
252 : call disp%show( simplified )
253 : call disp%show("operationA")
254 : call disp%show( operationA )
255 : call disp%show("shape(trimat)")
256 : call disp%show( shape(trimat) )
257 : call disp%show("shape(genvec)")
258 : call disp%show( shape(genvec) )
259 : call disp%show("shape(dumvec)")
260 : call disp%show( shape(dumvec) )
261 : call disp%show("[nrow, ndim]")
262 : call disp%show( [nrow, ndim] )
263 : call disp%show("[roffT, coffT, begG, incG, endG]")
264 : call disp%show( [roffT, coffT, begG, incG, endG] )
265 : call disp%show("trimat(roffT + 1 : roffT + nrow, coffT + 1 : coffT + nrow)")
266 : call disp%show( trimat(roffT + 1 : roffT + nrow, coffT + 1 : coffT + nrow) )
267 : call disp%show("genvec(begG:endG:incG)")
268 : call disp%show( genvec(begG:endG:incG) )
269 : call disp%show("dumvec(begG:endG:incG)")
270 : call disp%show( dumvec(begG:endG:incG) )
271 : call disp%show("tolerable")
272 : call disp%show( tolerable )
273 : call disp%show("pack(dumvec, .not. tolerable)")
274 : call disp%show( pack(dumvec, .not. tolerable) )
275 : call disp%show("pack(genvec, .not. tolerable)")
276 : call disp%show( pack(genvec, .not. tolerable) )
277 : call disp%show("pack(abs(genvec - dumvec), .not. tolerable)")
278 : call disp%show( pack(abs(genvec - dumvec), .not. tolerable) )
279 : call disp%show("TOL")
280 : call disp%show( TOL )
281 : ! LCOV_EXCL_STOP
282 : end if
283 19200 : call test%assert(assertion, SK_"TRMV/TRSV test with the above specifications must successfully pass.")
284 19200 : end subroutine
285 :
286 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
287 :
288 1600 : subroutine testTRXM(alpha, alinv)
289 : TYPE_KIND, intent(in), optional :: alpha, alinv
290 : character(:, SK), allocatable :: subset
291 : if (present(alpha) .neqv. present(alinv)) error stop "`alpha` and `alinv` must be both present or both missing." ! LCOV_EXCL_LINE
292 : #define sliceT trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim)
293 : #define sliceS solmat(roffG + 1 : roffG + nrow, coffG + 1 : coffG + ncol)
294 :
295 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
296 :
297 1600 : ndim = nrow
298 166888 : trimat = getChoice([-1, 1]) * getUnifRand(lb, ub, 2 * roffT + ndim, 2 * coffT + ndim)
299 159276 : genmat = getChoice([-1, 1]) * getUnifRand(lb, ub, 2 * roffG + nrow, 2 * coffG + ncol)
300 :
301 1600 : subset = "upperDiag"
302 :
303 : ! upperDiag, nothing/inversion
304 82774 : solmat = genmat
305 84336 : call setMatMulTri(sliceT, upperDiag, nothing, sliceS, alpha)
306 84336 : call setMatMulTri(sliceT, upperDiag, inversion, sliceS, alinv)
307 1600 : call reportTRXM(simplified, "nothing/inversion", subset, alpha)
308 :
309 : ! upperDiag, transSymm/transOrth
310 81238 : solmat = genmat
311 84336 : call setMatMulTri(sliceT, upperDiag, transSymm, sliceS, alpha)
312 84336 : call setMatMulTri(sliceT, upperDiag, transOrth, sliceS, alinv)
313 1600 : call reportTRXM(simplified, "transSymm/transOrth", subset, alpha)
314 :
315 : ! upperDiag, transHerm/transUnit
316 81238 : solmat = genmat
317 84336 : call setMatMulTri(sliceT, upperDiag, transHerm, sliceS, alpha)
318 84336 : call setMatMulTri(sliceT, upperDiag, transUnit, sliceS, alinv)
319 1600 : call reportTRXM(simplified, "transHerm/transUnit", subset, alpha)
320 :
321 1600 : subset = "upperDiag"
322 :
323 : ! lowerDiag, nothing/inversion
324 81238 : solmat = genmat
325 84336 : call setMatMulTri(sliceT, lowerDiag, nothing, sliceS, alpha)
326 84336 : call setMatMulTri(sliceT, lowerDiag, inversion, sliceS, alinv)
327 1600 : call reportTRXM(simplified, "nothing/inversion", subset, alpha)
328 :
329 : ! lowerDiag, transSymm/transOrth
330 81238 : solmat = genmat
331 84336 : call setMatMulTri(sliceT, lowerDiag, transSymm, sliceS, alpha)
332 84336 : call setMatMulTri(sliceT, lowerDiag, transOrth, sliceS, alinv)
333 1600 : call reportTRXM(simplified, "transSymm/transOrth", subset, alpha)
334 :
335 : ! lowerDiag, transHerm/transUnit
336 81238 : solmat = genmat
337 84336 : call setMatMulTri(sliceT, lowerDiag, transHerm, sliceS, alpha)
338 84336 : call setMatMulTri(sliceT, lowerDiag, transUnit, sliceS, alinv)
339 1600 : call reportTRXM(simplified, "transHerm/transUnit", subset, alpha)
340 :
341 1600 : subset = "upperUnit"
342 :
343 : ! upperUnit, nothing/inversion
344 81238 : solmat = genmat
345 84336 : call setMatMulTri(sliceT, upperUnit, nothing, sliceS, alpha)
346 84336 : call setMatMulTri(sliceT, upperUnit, inversion, sliceS, alinv)
347 1600 : call reportTRXM(simplified, "nothing/inversion", subset, alpha)
348 :
349 : ! upperUnit, transSymm/transOrth
350 81238 : solmat = genmat
351 84336 : call setMatMulTri(sliceT, upperUnit, transSymm, sliceS, alpha)
352 84336 : call setMatMulTri(sliceT, upperUnit, transOrth, sliceS, alinv)
353 1600 : call reportTRXM(simplified, "transSymm/transOrth", subset, alpha)
354 :
355 : ! upperUnit, transHerm/transUnit
356 81238 : solmat = genmat
357 84336 : call setMatMulTri(sliceT, upperUnit, transHerm, sliceS, alpha)
358 84336 : call setMatMulTri(sliceT, upperUnit, transUnit, sliceS, alinv)
359 1600 : call reportTRXM(simplified, "transHerm/transUnit", subset, alpha)
360 :
361 1600 : subset = "upperUnit"
362 :
363 : ! lowerUnit, nothing/inversion
364 81238 : solmat = genmat
365 84336 : call setMatMulTri(sliceT, lowerUnit, nothing, sliceS, alpha)
366 84336 : call setMatMulTri(sliceT, lowerUnit, inversion, sliceS, alinv)
367 1600 : call reportTRXM(simplified, "nothing/inversion", subset, alpha)
368 :
369 : ! lowerUnit, transSymm/transOrth
370 81238 : solmat = genmat
371 84336 : call setMatMulTri(sliceT, lowerUnit, transSymm, sliceS, alpha)
372 84336 : call setMatMulTri(sliceT, lowerUnit, transOrth, sliceS, alinv)
373 1600 : call reportTRXM(simplified, "transSymm/transOrth", subset, alpha)
374 :
375 : ! lowerUnit, transHerm/transUnit
376 81238 : solmat = genmat
377 84336 : call setMatMulTri(sliceT, lowerUnit, transHerm, sliceS, alpha)
378 84336 : call setMatMulTri(sliceT, lowerUnit, transUnit, sliceS, alinv)
379 1600 : call reportTRXM(simplified, "transHerm/transUnit", subset, alpha)
380 :
381 1600 : if (present(alpha)) then
382 :
383 800 : subset = "upperDiag"
384 :
385 : ! upperDiag, nothing/inversion
386 40619 : solmat = genmat
387 800 : call setMatMulTri(trimat, upperDiag, nothing, solmat, alpha, nrow, ncol, roffT, coffT, roffG, coffG)
388 800 : call setMatMulTri(trimat, upperDiag, inversion, solmat, alinv, nrow, ncol, roffT, coffT, roffG, coffG)
389 800 : call reportTRXM(simplified, "nothing/inversion", subset, alpha)
390 :
391 : ! upperDiag, transSymm/transOrth
392 40619 : solmat = genmat
393 800 : call setMatMulTri(trimat, upperDiag, transSymm, solmat, alpha, nrow, ncol, roffT, coffT, roffG, coffG)
394 800 : call setMatMulTri(trimat, upperDiag, transOrth, solmat, alinv, nrow, ncol, roffT, coffT, roffG, coffG)
395 800 : call reportTRXM(simplified, "transSymm/transOrth", subset, alpha)
396 :
397 : ! upperDiag, transHerm/transUnit
398 40619 : solmat = genmat
399 800 : call setMatMulTri(trimat, upperDiag, transHerm, solmat, alpha, nrow, ncol, roffT, coffT, roffG, coffG)
400 800 : call setMatMulTri(trimat, upperDiag, transUnit, solmat, alinv, nrow, ncol, roffT, coffT, roffG, coffG)
401 800 : call reportTRXM(simplified, "transHerm/transUnit", subset, alpha)
402 :
403 800 : subset = "upperDiag"
404 :
405 : ! lowerDiag, nothing/inversion
406 40619 : solmat = genmat
407 800 : call setMatMulTri(trimat, lowerDiag, nothing, solmat, alpha, nrow, ncol, roffT, coffT, roffG, coffG)
408 800 : call setMatMulTri(trimat, lowerDiag, inversion, solmat, alinv, nrow, ncol, roffT, coffT, roffG, coffG)
409 800 : call reportTRXM(simplified, "nothing/inversion", subset, alpha)
410 :
411 : ! lowerDiag, transSymm/transOrth
412 40619 : solmat = genmat
413 800 : call setMatMulTri(trimat, lowerDiag, transSymm, solmat, alpha, nrow, ncol, roffT, coffT, roffG, coffG)
414 800 : call setMatMulTri(trimat, lowerDiag, transOrth, solmat, alinv, nrow, ncol, roffT, coffT, roffG, coffG)
415 800 : call reportTRXM(simplified, "transSymm/transOrth", subset, alpha)
416 :
417 : ! lowerDiag, transHerm/transUnit
418 40619 : solmat = genmat
419 800 : call setMatMulTri(trimat, lowerDiag, transHerm, solmat, alpha, nrow, ncol, roffT, coffT, roffG, coffG)
420 800 : call setMatMulTri(trimat, lowerDiag, transUnit, solmat, alinv, nrow, ncol, roffT, coffT, roffG, coffG)
421 800 : call reportTRXM(simplified, "transHerm/transUnit", subset, alpha)
422 :
423 800 : subset = "upperUnit"
424 :
425 : ! upperUnit, nothing/inversion
426 40619 : solmat = genmat
427 800 : call setMatMulTri(trimat, upperUnit, nothing, solmat, alpha, nrow, ncol, roffT, coffT, roffG, coffG)
428 800 : call setMatMulTri(trimat, upperUnit, inversion, solmat, alinv, nrow, ncol, roffT, coffT, roffG, coffG)
429 800 : call reportTRXM(simplified, "nothing/inversion", subset, alpha)
430 :
431 : ! upperUnit, transSymm/transOrth
432 40619 : solmat = genmat
433 800 : call setMatMulTri(trimat, upperUnit, transSymm, solmat, alpha, nrow, ncol, roffT, coffT, roffG, coffG)
434 800 : call setMatMulTri(trimat, upperUnit, transOrth, solmat, alinv, nrow, ncol, roffT, coffT, roffG, coffG)
435 800 : call reportTRXM(simplified, "transSymm/transOrth", subset, alpha)
436 :
437 : ! upperUnit, transHerm/transUnit
438 40619 : solmat = genmat
439 800 : call setMatMulTri(trimat, upperUnit, transHerm, solmat, alpha, nrow, ncol, roffT, coffT, roffG, coffG)
440 800 : call setMatMulTri(trimat, upperUnit, transUnit, solmat, alinv, nrow, ncol, roffT, coffT, roffG, coffG)
441 800 : call reportTRXM(simplified, "transHerm/transUnit", subset, alpha)
442 :
443 800 : subset = "upperUnit"
444 :
445 : ! lowerUnit, nothing/inversion
446 40619 : solmat = genmat
447 800 : call setMatMulTri(trimat, lowerUnit, nothing, solmat, alpha, nrow, ncol, roffT, coffT, roffG, coffG)
448 800 : call setMatMulTri(trimat, lowerUnit, inversion, solmat, alinv, nrow, ncol, roffT, coffT, roffG, coffG)
449 800 : call reportTRXM(simplified, "nothing/inversion", subset, alpha)
450 :
451 : ! lowerUnit, transSymm/transOrth
452 40619 : solmat = genmat
453 800 : call setMatMulTri(trimat, lowerUnit, transSymm, solmat, alpha, nrow, ncol, roffT, coffT, roffG, coffG)
454 800 : call setMatMulTri(trimat, lowerUnit, transOrth, solmat, alinv, nrow, ncol, roffT, coffT, roffG, coffG)
455 800 : call reportTRXM(simplified, "transSymm/transOrth", subset, alpha)
456 :
457 : ! lowerUnit, transHerm/transUnit
458 40619 : solmat = genmat
459 800 : call setMatMulTri(trimat, lowerUnit, transHerm, solmat, alpha, nrow, ncol, roffT, coffT, roffG, coffG)
460 800 : call setMatMulTri(trimat, lowerUnit, transUnit, solmat, alinv, nrow, ncol, roffT, coffT, roffG, coffG)
461 800 : call reportTRXM(simplified, "transHerm/transUnit", subset, alpha)
462 :
463 : end if
464 :
465 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
466 :
467 1600 : ndim = ncol
468 159276 : genmat = getChoice([-1, 1]) * getUnifRand(lb, ub, 2 * roffG + nrow, 2 * coffG + ncol)
469 168272 : trimat = getChoice([-1, 1]) * getUnifRand(lb, ub, 2 * roffT + ndim, 2 * coffT + ndim)
470 81238 : solmat = genmat
471 :
472 1600 : subset = "upperDiag"
473 :
474 : ! upperDiag, nothing/inversion
475 81238 : solmat = genmat
476 84764 : call setMatMulTri(sliceS, sliceT, upperDiag, nothing, alpha)
477 84764 : call setMatMulTri(sliceS, sliceT, upperDiag, inversion, alinv)
478 1600 : call reportTRXM(simplified, "nothing/inversion", subset, alpha)
479 :
480 : ! upperDiag, transSymm/transOrth
481 81238 : solmat = genmat
482 84764 : call setMatMulTri(sliceS, sliceT, upperDiag, transSymm, alpha)
483 84764 : call setMatMulTri(sliceS, sliceT, upperDiag, transOrth, alinv)
484 1600 : call reportTRXM(simplified, "transSymm/transOrth", subset, alpha)
485 :
486 : ! upperDiag, transHerm/transUnit
487 81238 : solmat = genmat
488 84764 : call setMatMulTri(sliceS, sliceT, upperDiag, transHerm, alpha)
489 84764 : call setMatMulTri(sliceS, sliceT, upperDiag, transUnit, alinv)
490 1600 : call reportTRXM(simplified, "transHerm/transUnit", subset, alpha)
491 :
492 1600 : subset = "lowerDiag"
493 :
494 : ! lowerDiag, nothing/inversion
495 81238 : solmat = genmat
496 84764 : call setMatMulTri(sliceS, sliceT, lowerDiag, nothing, alpha)
497 84764 : call setMatMulTri(sliceS, sliceT, lowerDiag, inversion, alinv)
498 1600 : call reportTRXM(simplified, "nothing/inversion", subset, alpha)
499 :
500 : ! lowerDiag, transSymm/transOrth
501 81238 : solmat = genmat
502 84764 : call setMatMulTri(sliceS, sliceT, lowerDiag, transSymm, alpha)
503 84764 : call setMatMulTri(sliceS, sliceT, lowerDiag, transOrth, alinv)
504 1600 : call reportTRXM(simplified, "transSymm/transOrth", subset, alpha)
505 :
506 : ! lowerDiag, transHerm/transUnit
507 81238 : solmat = genmat
508 84764 : call setMatMulTri(sliceS, sliceT, lowerDiag, transHerm, alpha)
509 84764 : call setMatMulTri(sliceS, sliceT, lowerDiag, transUnit, alinv)
510 1600 : call reportTRXM(simplified, "transHerm/transUnit", subset, alpha)
511 :
512 1600 : subset = "upperUnit"
513 :
514 : ! upperUnit, nothing/inversion
515 81238 : solmat = genmat
516 84764 : call setMatMulTri(sliceS, sliceT, upperUnit, nothing, alpha)
517 84764 : call setMatMulTri(sliceS, sliceT, upperUnit, inversion, alinv)
518 1600 : call reportTRXM(simplified, "nothing/inversion", subset, alpha)
519 :
520 : ! upperUnit, transSymm/transOrth
521 81238 : solmat = genmat
522 84764 : call setMatMulTri(sliceS, sliceT, upperUnit, transSymm, alpha)
523 84764 : call setMatMulTri(sliceS, sliceT, upperUnit, transOrth, alinv)
524 1600 : call reportTRXM(simplified, "transSymm/transOrth", subset, alpha)
525 :
526 : ! upperUnit, transHerm/transUnit
527 81238 : solmat = genmat
528 84764 : call setMatMulTri(sliceS, sliceT, upperUnit, transHerm, alpha)
529 84764 : call setMatMulTri(sliceS, sliceT, upperUnit, transUnit, alinv)
530 1600 : call reportTRXM(simplified, "transHerm/transUnit", subset, alpha)
531 :
532 1600 : subset = "lowerUnit"
533 :
534 : ! lowerUnit, nothing/inversion
535 81238 : solmat = genmat
536 84764 : call setMatMulTri(sliceS, sliceT, lowerUnit, nothing, alpha)
537 84764 : call setMatMulTri(sliceS, sliceT, lowerUnit, inversion, alinv)
538 1600 : call reportTRXM(simplified, "nothing/inversion", subset, alpha)
539 :
540 : ! lowerUnit, transSymm/transOrth
541 81238 : solmat = genmat
542 84764 : call setMatMulTri(sliceS, sliceT, lowerUnit, transSymm, alpha)
543 84764 : call setMatMulTri(sliceS, sliceT, lowerUnit, transOrth, alinv)
544 1600 : call reportTRXM(simplified, "transSymm/transOrth", subset, alpha)
545 :
546 : ! lowerUnit, transHerm/transUnit
547 81238 : solmat = genmat
548 84764 : call setMatMulTri(sliceS, sliceT, lowerUnit, transHerm, alpha)
549 84764 : call setMatMulTri(sliceS, sliceT, lowerUnit, transUnit, alinv)
550 1600 : call reportTRXM(simplified, "transHerm/transUnit", subset, alpha)
551 :
552 1600 : if (present(alpha)) then
553 :
554 800 : subset = "upperDiag"
555 :
556 : ! upperDiag, nothing/inversion
557 40619 : solmat = genmat
558 800 : call setMatMulTri(solmat, trimat, upperDiag, nothing, alpha, nrow, ncol, roffG, coffG, roffT, coffT)
559 800 : call setMatMulTri(solmat, trimat, upperDiag, inversion, alinv, nrow, ncol, roffG, coffG, roffT, coffT)
560 800 : call reportTRXM(simplified, "nothing/inversion", subset, alpha)
561 :
562 : ! upperDiag, transSymm/transOrth
563 40619 : solmat = genmat
564 800 : call setMatMulTri(solmat, trimat, upperDiag, transSymm, alpha, nrow, ncol, roffG, coffG, roffT, coffT)
565 800 : call setMatMulTri(solmat, trimat, upperDiag, transOrth, alinv, nrow, ncol, roffG, coffG, roffT, coffT)
566 800 : call reportTRXM(simplified, "transSymm/transOrth", subset, alpha)
567 :
568 : ! upperDiag, transHerm/transUnit
569 40619 : solmat = genmat
570 800 : call setMatMulTri(solmat, trimat, upperDiag, transHerm, alpha, nrow, ncol, roffG, coffG, roffT, coffT)
571 800 : call setMatMulTri(solmat, trimat, upperDiag, transUnit, alinv, nrow, ncol, roffG, coffG, roffT, coffT)
572 800 : call reportTRXM(simplified, "transHerm/transUnit", subset, alpha)
573 :
574 800 : subset = "lowerDiag"
575 :
576 : ! lowerDiag, nothing/inversion
577 40619 : solmat = genmat
578 800 : call setMatMulTri(solmat, trimat, lowerDiag, nothing, alpha, nrow, ncol, roffG, coffG, roffT, coffT)
579 800 : call setMatMulTri(solmat, trimat, lowerDiag, inversion, alinv, nrow, ncol, roffG, coffG, roffT, coffT)
580 800 : call reportTRXM(simplified, "nothing/inversion", subset, alpha)
581 :
582 : ! lowerDiag, transSymm/transOrth
583 40619 : solmat = genmat
584 800 : call setMatMulTri(solmat, trimat, lowerDiag, transSymm, alpha, nrow, ncol, roffG, coffG, roffT, coffT)
585 800 : call setMatMulTri(solmat, trimat, lowerDiag, transOrth, alinv, nrow, ncol, roffG, coffG, roffT, coffT)
586 800 : call reportTRXM(simplified, "transSymm/transOrth", subset, alpha)
587 :
588 : ! lowerDiag, transHerm/transUnit
589 40619 : solmat = genmat
590 800 : call setMatMulTri(solmat, trimat, lowerDiag, transHerm, alpha, nrow, ncol, roffG, coffG, roffT, coffT)
591 800 : call setMatMulTri(solmat, trimat, lowerDiag, transUnit, alinv, nrow, ncol, roffG, coffG, roffT, coffT)
592 800 : call reportTRXM(simplified, "transHerm/transUnit", subset, alpha)
593 :
594 800 : subset = "upperUnit"
595 :
596 : ! upperUnit, nothing/inversion
597 40619 : solmat = genmat
598 800 : call setMatMulTri(solmat, trimat, upperUnit, nothing, alpha, nrow, ncol, roffG, coffG, roffT, coffT)
599 800 : call setMatMulTri(solmat, trimat, upperUnit, inversion, alinv, nrow, ncol, roffG, coffG, roffT, coffT)
600 800 : call reportTRXM(simplified, "nothing/inversion", subset, alpha)
601 :
602 : ! upperUnit, transSymm/transOrth
603 40619 : solmat = genmat
604 800 : call setMatMulTri(solmat, trimat, upperUnit, transSymm, alpha, nrow, ncol, roffG, coffG, roffT, coffT)
605 800 : call setMatMulTri(solmat, trimat, upperUnit, transOrth, alinv, nrow, ncol, roffG, coffG, roffT, coffT)
606 800 : call reportTRXM(simplified, "transSymm/transOrth", subset, alpha)
607 :
608 : ! upperUnit, transHerm/transUnit
609 40619 : solmat = genmat
610 800 : call setMatMulTri(solmat, trimat, upperUnit, transHerm, alpha, nrow, ncol, roffG, coffG, roffT, coffT)
611 800 : call setMatMulTri(solmat, trimat, upperUnit, transUnit, alinv, nrow, ncol, roffG, coffG, roffT, coffT)
612 800 : call reportTRXM(simplified, "transHerm/transUnit", subset, alpha)
613 :
614 800 : subset = "lowerUnit"
615 :
616 : ! lowerUnit, nothing/inversion
617 40619 : solmat = genmat
618 800 : call setMatMulTri(solmat, trimat, lowerUnit, nothing, alpha, nrow, ncol, roffG, coffG, roffT, coffT)
619 800 : call setMatMulTri(solmat, trimat, lowerUnit, inversion, alinv, nrow, ncol, roffG, coffG, roffT, coffT)
620 800 : call reportTRXM(simplified, "nothing/inversion", subset, alpha)
621 :
622 : ! lowerUnit, transSymm/transOrth
623 40619 : solmat = genmat
624 800 : call setMatMulTri(solmat, trimat, lowerUnit, transSymm, alpha, nrow, ncol, roffG, coffG, roffT, coffT)
625 800 : call setMatMulTri(solmat, trimat, lowerUnit, transOrth, alinv, nrow, ncol, roffG, coffG, roffT, coffT)
626 800 : call reportTRXM(simplified, "transSymm/transOrth", subset, alpha)
627 :
628 : ! lowerUnit, transHerm/transUnit
629 40619 : solmat = genmat
630 800 : call setMatMulTri(solmat, trimat, lowerUnit, transHerm, alpha, nrow, ncol, roffG, coffG, roffT, coffT)
631 800 : call setMatMulTri(solmat, trimat, lowerUnit, transUnit, alinv, nrow, ncol, roffG, coffG, roffT, coffT)
632 800 : call reportTRXM(simplified, "transHerm/transUnit", subset, alpha)
633 :
634 : end if
635 :
636 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
637 :
638 : #undef sliceT
639 : #undef sliceS
640 1600 : end subroutine
641 :
642 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
643 :
644 57600 : subroutine reportTRXM(simplified, operation, subset, alpha)
645 : logical(LK) , intent(in) :: simplified
646 : character(*), intent(in) :: operation, subset
647 : TYPE_KIND, intent(in), optional :: alpha
648 : logical(LK), allocatable :: tolerable(:,:)
649 57600 : type(display_type) :: disp
650 5733936 : tolerable = isClose(genmat, solmat, abstol = TOL)
651 2684808 : assertion = assertion .and. (all(tolerable) .or. (all(solmat(roffG + 1 : roffG + nrow, coffG + 1 : coffG + ncol) == ZERO) .and. getOption(ONE, alpha) == ZERO))
652 57600 : if (test%traceable .and. .not. assertion) then
653 : ! LCOV_EXCL_START
654 : call disp%show("TRMM/TRSM")
655 : call disp%show("simplified")
656 : call disp%show( simplified )
657 : call disp%show("subset")
658 : call disp%show( subset )
659 : call disp%show("operation")
660 : call disp%show( operation )
661 : call disp%show("shape(trimat)")
662 : call disp%show( shape(trimat) )
663 : call disp%show("shape(genmat)")
664 : call disp%show( shape(genmat) )
665 : call disp%show("shape(solmat)")
666 : call disp%show( shape(solmat) )
667 : call disp%show("present(alpha)")
668 : call disp%show( present(alpha) )
669 : call disp%show("getOption(ONE, alpha)")
670 : call disp%show( getOption(ONE, alpha) )
671 : call disp%show("[nrow, ncol, ndim]")
672 : call disp%show( [nrow, ncol, ndim] )
673 : call disp%show("[roffT, coffT, roffG, coffG]")
674 : call disp%show( [roffT, coffT, roffG, coffG] )
675 : call disp%show("trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim)")
676 : call disp%show( trimat(roffT + 1 : roffT + ndim, coffT + 1 : coffT + ndim) )
677 : call disp%show("genmat(roffG + 1 : roffG + nrow, coffG + 1 : coffG + ncol)")
678 : call disp%show( genmat(roffG + 1 : roffG + nrow, coffG + 1 : coffG + ncol) )
679 : call disp%show("solmat(roffG + 1 : roffG + nrow, coffG + 1 : coffG + ncol)")
680 : call disp%show( solmat(roffG + 1 : roffG + nrow, coffG + 1 : coffG + ncol) )
681 : call disp%show("tolerable")
682 : call disp%show( tolerable )
683 : call disp%show("pack(solmat, .not. tolerable)")
684 : call disp%show( pack(solmat, .not. tolerable) )
685 : call disp%show("pack(genmat, .not. tolerable)")
686 : call disp%show( pack(genmat, .not. tolerable) )
687 : call disp%show("pack(abs(genmat - solmat), .not. tolerable)")
688 : call disp%show( pack(abs(genmat - solmat), .not. tolerable) )
689 : call disp%show("TOL")
690 : call disp%show( TOL )
691 : ! LCOV_EXCL_STOP
692 : end if
693 57600 : call test%assert(assertion, SK_"TRMM/TRSM test with the above specifications must successfully pass.")
694 57600 : end subroutine
695 :
696 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
697 :
698 : #undef TYPE_KIND
|