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 module contains implementations of the tests of the procedures under the generic interfaces [pm_arrayCenter](@ref pm_arrayCenter).
19 : !>
20 : !> \fintest
21 : !>
22 : !> \author
23 : !> \AmirShahmoradi, September 1, 2017, 11:35 PM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
24 :
25 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
26 :
27 : #if LK_ENABLED
28 : #define IS_EQUAL .eqv.
29 : #else
30 : #define IS_EQUAL ==
31 : #endif
32 :
33 : #if SK_ENABLED && D0_ENABLED
34 : #define GET_LBOUND(Array) 1_IK
35 : #define GEN_UBOUND(Array) len(Array, kind = IK)
36 : #define GET_SIZE(Array) len(Array, kind = IK)
37 : #define ALLOC_ARRAY(Array,lb,ub) allocate(character(ub-lb+1,SKC) :: Array)
38 : #define GEN_UBOLD(ub) ub - lbold + 1_IK
39 : #define GEN_UBNEW(ub) ub - lbnew + 1_IK + getOption(0_IK,lmsize) + getOption(0_IK,rmsize)
40 : #define GEN_LBOLD(lb) 1_IK
41 : #define GEN_LBNEW(lb) 1_IK
42 : #elif D1_ENABLED
43 : #define ALLOC_ARRAY(Array,lb,ub) allocate(Array(lb:ub))
44 : #define GET_LBOUND(Array) lbound(Array, dim = 1, kind = IK)
45 : #define GEN_UBOUND(Array) ubound(Array, dim = 1, kind = IK)
46 : #define GET_SIZE(Array) size(Array, kind = IK)
47 : #define GEN_UBOLD(ub) ub
48 : #define GEN_UBNEW(ub) ub + getOption(0_IK,lmsize) + getOption(0_IK,rmsize)
49 : #define GEN_LBOLD(lb) lb
50 : #define GEN_LBNEW(lb) lb
51 : #else
52 : #error "Unrecognized interface."
53 : #endif
54 :
55 : #if SK_ENABLED && D0_ENABLED
56 : #define ALL
57 2 : character(:,SKC), allocatable :: Array, ArrayCentered, ArrayCentered_ref
58 : character(1,SKC), parameter :: fill = SKC_"*"
59 : character(1,SKC), parameter :: lmfill = SKC_"-"
60 : character(1,SKC), parameter :: rmfill = SKC_"+"
61 : character(1,SKC), parameter :: fill_def = fill
62 : character(1,SKC), parameter :: lmfill_def = lmfill
63 : character(1,SKC), parameter :: rmfill_def = rmfill
64 : #elif SK_ENABLED && D1_ENABLED
65 : character(2,SKC), dimension(:), allocatable :: Array, ArrayCentered, ArrayCentered_ref
66 : character(2,SKC), parameter :: fill = SKC_"**"
67 : character(2,SKC), parameter :: lmfill = SKC_"--"
68 : character(2,SKC), parameter :: rmfill = SKC_"++"
69 : character(2,SKC), parameter :: fill_def = fill
70 : character(2,SKC), parameter :: lmfill_def = lmfill
71 : character(2,SKC), parameter :: rmfill_def = rmfill
72 : #elif LK_ENABLED && D1_ENABLED
73 : logical(LKC) , dimension(:), allocatable :: Array, ArrayCentered, ArrayCentered_ref
74 : logical(LKC) , parameter :: fill = .false._LKC
75 : logical(LKC) , parameter :: lmfill = .false._LKC
76 : logical(LKC) , parameter :: rmfill = .false._LKC
77 : logical(LKC) , parameter :: fill_def = fill
78 : logical(LKC) , parameter :: lmfill_def = lmfill
79 : logical(LKC) , parameter :: rmfill_def = rmfill
80 : #elif IK_ENABLED && D1_ENABLED
81 : integer(IKC) , dimension(:), allocatable :: Array, ArrayCentered, ArrayCentered_ref
82 : integer(IKC) , parameter :: fill = huge(1_IKC)
83 : integer(IKC) , parameter :: lmfill = huge(1_IKC)
84 : integer(IKC) , parameter :: rmfill = huge(1_IKC)
85 : integer(IKC) , parameter :: fill_def = fill
86 : integer(IKC) , parameter :: lmfill_def = lmfill
87 : integer(IKC) , parameter :: rmfill_def = rmfill
88 : #elif CK_ENABLED && D1_ENABLED
89 : complex(CKC) , dimension(:), allocatable :: Array, ArrayCentered, ArrayCentered_ref
90 : complex(CKC) , parameter :: fill = cmplx(huge(0._CKC), huge(0._CKC), kind = CKC)
91 : complex(CKC) , parameter :: lmfill = cmplx(huge(0._CKC), huge(0._CKC), kind = CKC)
92 : complex(CKC) , parameter :: rmfill = cmplx(huge(0._CKC), huge(0._CKC), kind = CKC)
93 : complex(CKC) , parameter :: fill_def = fill
94 : complex(CKC) , parameter :: lmfill_def = lmfill
95 : complex(CKC) , parameter :: rmfill_def = rmfill
96 : #elif RK_ENABLED && D1_ENABLED
97 : real(RKC) , dimension(:), allocatable :: Array, ArrayCentered, ArrayCentered_ref
98 : real(RKC) , parameter :: fill = huge(0._RKC)
99 : real(RKC) , parameter :: lmfill = huge(0._RKC)
100 : real(RKC) , parameter :: rmfill = huge(0._RKC)
101 : real(RKC) , parameter :: fill_def = fill
102 : real(RKC) , parameter :: lmfill_def = lmfill
103 : real(RKC) , parameter :: rmfill_def = rmfill
104 : #else
105 : #error "Unrecognized interface."
106 : #endif
107 : integer(IK) :: cbeg
108 : integer(IK) :: cend
109 : integer(IK) :: sizeold
110 : integer(IK) :: sizenew
111 : integer(IK) :: sizecentered
112 : integer(IK) :: lbold, ubold
113 : integer(IK) :: lbnew, ubnew
114 : logical(LK) :: menabled
115 :
116 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
117 :
118 40 : assertion = .true._LK
119 40 : call runTestsWith()
120 40 : call runTestsWith(fill = fill)
121 40 : call runTestsWith(lmsize = 2_IK, rmsize = 3_IK)
122 40 : call runTestsWith(lmsize = 2_IK, rmsize = 3_IK, fill = fill)
123 40 : call runTestsWith(lmsize = 2_IK, rmsize = 3_IK, lmfill = lmfill)
124 40 : call runTestsWith(lmsize = 2_IK, rmsize = 3_IK, lmfill = lmfill, fill = fill)
125 40 : call runTestsWith(lmsize = 2_IK, rmsize = 3_IK, rmfill = rmfill)
126 40 : call runTestsWith(lmsize = 2_IK, rmsize = 3_IK, rmfill = rmfill, fill = fill)
127 40 : call runTestsWith(lmsize = 2_IK, rmsize = 3_IK, lmfill = lmfill, rmfill = rmfill)
128 40 : call runTestsWith(lmsize = 2_IK, rmsize = 3_IK, lmfill = lmfill, rmfill = rmfill, fill = fill)
129 :
130 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
131 :
132 : contains
133 :
134 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
135 :
136 400 : subroutine runTestsWith (lmsize, rmsize, fill, lmfill, rmfill)
137 :
138 : integer(IK) , intent(in), optional :: lmsize, rmsize
139 : #if SK_ENABLED && D0_ENABLED
140 : character(1,SKC), intent(in), optional :: fill, lmfill, rmfill
141 : #elif SK_ENABLED && D1_ENABLED
142 : character(2,SKC), intent(in), optional :: fill, lmfill, rmfill
143 : #elif LK_ENABLED && D1_ENABLED
144 : logical(LKC) , intent(in), optional :: fill, lmfill, rmfill
145 : #elif IK_ENABLED && D1_ENABLED
146 : integer(IKC) , intent(in), optional :: fill, lmfill, rmfill
147 : #elif CK_ENABLED && D1_ENABLED
148 : complex(CKC) , intent(in), optional :: fill, lmfill, rmfill
149 : #elif RK_ENABLED && D1_ENABLED
150 : real(RKC) , intent(in), optional :: fill, lmfill, rmfill
151 : #else
152 : #error "Unrecognized interface."
153 : #endif
154 :
155 400 : if (present(lmsize) .neqv. present(rmsize)) error stop "Internal ParaMonte Testing error occurred: present(lmsize) .neqv. present(rmsize) must be .true."
156 400 : menabled = present(lmsize) .and. present(rmsize)
157 :
158 : !> \bug
159 : !> GNU Fortran 10.3 cannot concatenate empty character array of length 2 with a non-empty character array of the same length.
160 : !> Fortran runtime error: Different CHARACTER lengths (0/2) in array constructor
161 400 : if (present(lmsize) .and. present(rmsize)) then
162 320 : if (lmsize == 0_IK .and. rmsize == 0_IK) error stop "Internal ParaMonte Testing error occurred: GNU bug exception."
163 : end if
164 :
165 400 : assertion = .true._LK
166 :
167 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
168 : ! Enlarge and center array
169 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
170 :
171 400 : call reset()
172 :
173 400 : lbold = GEN_LBOLD(1_IK)
174 400 : ubold = GEN_UBOLD(3_IK)
175 400 : lbnew = GEN_LBNEW(1_IK)
176 400 : ubnew = GEN_UBNEW(4_IK)
177 400 : sizeold = ubold - lbold + 1_IK
178 400 : sizenew = ubnew - lbnew + 1_IK
179 400 : sizecentered = sizenew - getOption(0_IK, lmsize) - getOption(0_IK, rmsize)
180 400 : ALLOC_ARRAY(Array,lbold,ubold)
181 440 : ALLOC_ARRAY(ArrayCentered,lbnew,ubnew)
182 400 : ALLOC_ARRAY(ArrayCentered_ref,lbnew,ubnew)
183 : #if SK_ENABLED && D0_ENABLED
184 140 : call setUnifRand(Array, repeat(SKC_"A", len(Array)), repeat(SKC_"Z",len(Array)))
185 68 : ArrayCentered_ref = getRepeat(lmfill, lmsize)//Array//fill_def//getRepeat(rmfill, rmsize)
186 : #elif SK_ENABLED && D1_ENABLED
187 80 : Array(:) = [SKC_"AA", SKC_"BB", SKC_"CC"]
188 20 : if (getOption(0_IK,lmsize) > 0_IK) then
189 320 : ArrayCentered_ref(:) = [getRepeat(lmfill,lmsize), Array, fill_def, getRepeat(rmfill,rmsize)]
190 : else
191 32 : ArrayCentered_ref(:) = [Array, fill_def]
192 : end if
193 : #elif LK_ENABLED && D1_ENABLED
194 400 : Array(:) = [.true._LKC, .true._LKC, .true._LKC]
195 1600 : ArrayCentered_ref(:) = [getRepeat(lmfill,lmsize), Array, fill_def, getRepeat(rmfill,rmsize)]
196 : #elif IK_ENABLED && D1_ENABLED
197 400 : Array(:) = [1_IKC, 2_IKC, 3_IKC]
198 1600 : ArrayCentered_ref(:) = [getRepeat(lmfill,lmsize), Array, fill_def, getRepeat(rmfill,rmsize)]
199 : #elif CK_ENABLED && D1_ENABLED
200 320 : Array(:) = [(1._CKC,-1._CKC), (2._CKC,-2._CKC), (3._CKC,-3._CKC)]
201 1280 : ArrayCentered_ref(:) = [getRepeat(lmfill,lmsize), Array, fill_def, getRepeat(rmfill,rmsize)]
202 : #elif RK_ENABLED && D1_ENABLED
203 320 : Array(:) = [1._RKC, 2._RKC, 3._RKC]
204 1280 : ArrayCentered_ref(:) = [getRepeat(lmfill,lmsize), Array, fill_def, getRepeat(rmfill,rmsize)]
205 : #endif
206 400 : cbeg = lbold + getOption(0_IK,lmsize) + (sizecentered - sizeold) / 2_IK
207 400 : cend = cbeg + sizeold - 1_IK
208 :
209 : #if setCentered_ENABLED
210 200 : if (menabled) then
211 192 : call setCentered(ArrayCentered, Array, lmsize, rmsize, fill, lmfill, rmfill)
212 : else
213 44 : call setCentered(ArrayCentered, Array, fill)
214 : end if
215 :
216 390 : assertion = assertion .and. GET_LBOUND(ArrayCentered) == lbnew
217 200 : call report()
218 200 : call test%assert(assertion, SK_"Call to setCentered()/getCentered() must properly set the lower bound of the output array, with present(fill) = "//getStr(present(fill)), int(__LINE__, IK))
219 :
220 390 : assertion = assertion .and. GEN_UBOUND(ArrayCentered) == ubnew
221 200 : call report()
222 200 : call test%assert(assertion, SK_"Call to setCentered()/getCentered() must properly set the upper bound of the output array, with present(fill) = "//getStr(present(fill)), int(__LINE__, IK))
223 : #elif getCentered_ENABLED
224 200 : if (menabled) then
225 1712 : ArrayCentered = getCentered(Array, sizecentered, lmsize, rmsize, fill, lmfill, rmfill)
226 : else
227 234 : ArrayCentered = getCentered(Array, sizecentered, fill)
228 : end if
229 : #else
230 : #error "Unrecognized interface."
231 : #endif
232 :
233 400 : assertion = assertion .and. GET_SIZE(ArrayCentered) == sizenew
234 400 : call report()
235 400 : call test%assert(assertion, SK_"Call to setCentered()/getCentered() must yield an array of proper size, with present(fill) = "//getStr(present(fill)), int(__LINE__, IK))
236 :
237 : #if setCentered_ENABLED
238 770 : assertion = assertion .and. ALL(ArrayCentered(cbeg:cend) IS_EQUAL ArrayCentered_ref(cbeg:cend))
239 : #elif getCentered_ENABLED
240 770 : assertion = assertion .and. ALL(ArrayCentered(cbeg-lbold+1_IK:cend-lbold+1_IK) IS_EQUAL ArrayCentered_ref(cbeg:cend))
241 : #endif
242 400 : call report()
243 400 : call test%assert(assertion, SK_"Call to setCentered()/getCentered() must properly set the upper bound of the output array, with present(fill) = "//getStr(present(fill)), int(__LINE__, IK))
244 :
245 400 : if (present(fill)) then
246 1340 : assertion = assertion .and. ALL(ArrayCentered(GET_LBOUND(ArrayCentered)+getOption(0_IK,lmsize):GEN_UBOUND(ArrayCentered)-getOption(0_IK,rmsize)) IS_EQUAL ArrayCentered_ref(lbnew+getOption(0_IK,lmsize):ubnew-getOption(0_IK,rmsize)))
247 200 : call report()
248 200 : call test%assert(assertion, SK_"Call to setCentered()/getCentered() must properly fill the new elements with `fill` = "//getStr(present(fill)), int(__LINE__, IK))
249 : end if
250 :
251 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
252 : ! Enlarge and center array perfectly
253 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
254 :
255 400 : call reset()
256 :
257 400 : lbold = GEN_LBOLD(1_IK)
258 400 : ubold = GEN_UBOLD(3_IK)
259 400 : lbnew = GEN_LBNEW(1_IK)
260 400 : ubnew = GEN_UBNEW(5_IK)
261 400 : sizeold = ubold - lbold + 1_IK
262 400 : sizenew = ubnew - lbnew + 1_IK
263 400 : sizecentered = sizenew - getOption(0_IK,lmsize) - getOption(0_IK,rmsize)
264 400 : ALLOC_ARRAY(Array,lbold,ubold)
265 440 : ALLOC_ARRAY(ArrayCentered,lbnew,ubnew)
266 400 : ALLOC_ARRAY(ArrayCentered_ref,lbnew,ubnew)
267 : #if SK_ENABLED && D0_ENABLED
268 140 : call setUnifRand(Array, repeat(SKC_"A",len(Array)), repeat(SKC_"Z",len(Array)))
269 20 : ArrayCentered_ref = getRepeat(lmfill,lmsize)//fill_def//Array//fill_def//getRepeat(rmfill,rmsize)
270 : #elif SK_ENABLED && D1_ENABLED
271 80 : Array(:) = [SKC_"AA", SKC_"BB", SKC_"CC"]
272 20 : if (getOption(0_IK,lmsize) > 0_IK) then
273 336 : ArrayCentered_ref(:) = [getRepeat(lmfill,lmsize), fill_def, Array, fill_def, getRepeat(rmfill,rmsize)]
274 : else
275 36 : ArrayCentered_ref(:) = [fill_def, Array, fill_def]
276 : end if
277 : #elif LK_ENABLED && D1_ENABLED
278 400 : Array(:) = [.true._LKC, .true._LKC, .true._LKC]
279 1700 : ArrayCentered_ref(:) = [getRepeat(lmfill,lmsize), fill_def, Array, fill_def, getRepeat(rmfill,rmsize)]
280 : #elif IK_ENABLED && D1_ENABLED
281 400 : Array(:) = [1_IKC, 2_IKC, 3_IKC]
282 1700 : ArrayCentered_ref(:) = [getRepeat(lmfill,lmsize), fill_def, Array, fill_def, getRepeat(rmfill,rmsize)]
283 : #elif CK_ENABLED && D1_ENABLED
284 320 : Array(:) = [(1._CKC,-1._CKC), (2._CKC,-2._CKC), (3._CKC,-3._CKC)]
285 1360 : ArrayCentered_ref(:) = [getRepeat(lmfill,lmsize), fill_def, Array, fill_def, getRepeat(rmfill,rmsize)]
286 : #elif RK_ENABLED && D1_ENABLED
287 320 : Array(:) = [1._RKC, 2._RKC, 3._RKC]
288 1360 : ArrayCentered_ref(:) = [getRepeat(lmfill,lmsize), fill_def, Array, fill_def, getRepeat(rmfill,rmsize)]
289 : #endif
290 400 : cbeg = lbold + getOption(0_IK,lmsize) + (sizecentered - sizeold) / 2_IK
291 400 : cend = cbeg + sizeold - 1_IK
292 :
293 : #if setCentered_ENABLED
294 200 : if (menabled) then
295 192 : call setCentered(ArrayCentered, Array, lmsize, rmsize, fill, lmfill, rmfill)
296 : else
297 44 : call setCentered(ArrayCentered, Array, fill)
298 : end if
299 :
300 390 : assertion = assertion .and. GET_LBOUND(ArrayCentered) == lbnew
301 200 : call report()
302 200 : call test%assert(assertion, SK_"Call to setCentered()/getCentered() must properly set the lower bound of the output array, with present(fill) = "//getStr(present(fill)), int(__LINE__, IK))
303 :
304 390 : assertion = assertion .and. GEN_UBOUND(ArrayCentered) == ubnew
305 200 : call report()
306 200 : call test%assert(assertion, SK_"Call to setCentered()/getCentered() must properly set the upper bound of the output array, with present(fill) = "//getStr(present(fill)), int(__LINE__, IK))
307 : #elif getCentered_ENABLED
308 200 : if (menabled) then
309 1864 : ArrayCentered = getCentered(Array, sizecentered, lmsize, rmsize, fill, lmfill, rmfill)
310 : else
311 272 : ArrayCentered = getCentered(Array, sizecentered, fill)
312 : end if
313 : #else
314 : #error "Unrecognized interface."
315 : #endif
316 :
317 400 : assertion = assertion .and. GET_SIZE(ArrayCentered) == sizenew
318 400 : call report()
319 400 : call test%assert(assertion, SK_"Call to setCentered()/getCentered() must yield an array of proper size, with present(fill) = "//getStr(present(fill)), int(__LINE__, IK))
320 :
321 : #if setCentered_ENABLED
322 770 : assertion = assertion .and. ALL(ArrayCentered(cbeg:cend) IS_EQUAL ArrayCentered_ref(cbeg:cend))
323 : #elif getCentered_ENABLED
324 770 : assertion = assertion .and. ALL(ArrayCentered(cbeg-lbold+1_IK:cend-lbold+1_IK) IS_EQUAL ArrayCentered_ref(cbeg:cend))
325 : #endif
326 400 : call report()
327 400 : call test%assert(assertion, SK_"Call to setCentered()/getCentered() must properly set the upper bound of the output array, with present(fill) = "//getStr(present(fill)), int(__LINE__, IK))
328 :
329 400 : if (present(fill)) then
330 1530 : assertion = assertion .and. ALL(ArrayCentered(GET_LBOUND(ArrayCentered)+getOption(0_IK,lmsize):GEN_UBOUND(ArrayCentered)-getOption(0_IK,rmsize)) IS_EQUAL ArrayCentered_ref(lbnew+getOption(0_IK,lmsize):ubnew-getOption(0_IK,rmsize)))
331 200 : call report()
332 200 : call test%assert(assertion, SK_"Call to setCentered()/getCentered() must properly fill the new elements with `fill` = "//getStr(present(fill)), int(__LINE__, IK))
333 : end if
334 :
335 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
336 : ! Shrink and center array
337 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
338 :
339 400 : call reset()
340 :
341 400 : lbold = GEN_LBOLD(1_IK)
342 400 : ubold = GEN_UBOLD(4_IK)
343 400 : lbnew = GEN_LBNEW(1_IK)
344 400 : ubnew = GEN_UBNEW(3_IK)
345 400 : sizeold = ubold - lbold + 1_IK
346 400 : sizenew = ubnew - lbnew + 1_IK
347 400 : sizecentered = sizenew - getOption(0_IK,lmsize) - getOption(0_IK,rmsize)
348 400 : ALLOC_ARRAY(Array,lbold,ubold)
349 440 : ALLOC_ARRAY(ArrayCentered,lbnew,ubnew)
350 400 : ALLOC_ARRAY(ArrayCentered_ref,lbnew,ubnew)
351 : #if SK_ENABLED && D0_ENABLED
352 180 : call setUnifRand(Array, repeat(SKC_"A",len(Array)), repeat(SKC_"Z",len(Array)))
353 20 : ArrayCentered_ref(lbnew:ubnew) = getRepeat(lmfill,lmsize)//Array(lbold:ubold-1)//getRepeat(rmfill,rmsize)
354 : #elif SK_ENABLED && D1_ENABLED
355 100 : Array(:) = [SKC_"AA", SKC_"BB", SKC_"CC", SKC_"DD"]
356 20 : if (getOption(0_IK,lmsize) > 0_IK) then
357 304 : ArrayCentered_ref(lbnew:ubnew) = [getRepeat(lmfill,lmsize), Array(lbold:ubold-1), getRepeat(rmfill,rmsize)]
358 : else
359 28 : ArrayCentered_ref(lbnew:ubnew) = [Array(lbold:ubold-1)]
360 : end if
361 : #elif LK_ENABLED && D1_ENABLED
362 500 : Array(:) = [.true._LKC, .true._LKC, .true._LKC, .true._LKC]
363 1500 : ArrayCentered_ref(lbnew:ubnew) = [getRepeat(lmfill,lmsize), Array(lbold:ubold-1), getRepeat(rmfill,rmsize)]
364 : #elif IK_ENABLED && D1_ENABLED
365 500 : Array(:) = [1_IKC, 2_IKC, 3_IKC, 4_IKC]
366 1500 : ArrayCentered_ref(lbnew:ubnew) = [getRepeat(lmfill,lmsize), Array(lbold:ubold-1), getRepeat(rmfill,rmsize)]
367 : #elif CK_ENABLED && D1_ENABLED
368 400 : Array(:) = [(1._CKC,-1._CKC), (2._CKC,-2._CKC), (3._CKC,-3._CKC), (4._CKC,-4._CKC)]
369 1200 : ArrayCentered_ref(lbnew:ubnew) = [getRepeat(lmfill,lmsize), Array(lbold:ubold-1), getRepeat(rmfill,rmsize)]
370 : #elif RK_ENABLED && D1_ENABLED
371 400 : Array(:) = [1._RKC, 2._RKC, 3._RKC, 4._RKC]
372 1200 : ArrayCentered_ref(lbnew:ubnew) = [getRepeat(lmfill,lmsize), Array(lbold:ubold-1), getRepeat(rmfill,rmsize)]
373 : #endif
374 400 : cbeg = lbold + getOption(0_IK,lmsize) + abs(sizecentered - sizeold) / 2_IK
375 400 : cend = cbeg + min(sizeold,sizecentered) - 1_IK
376 :
377 : #if setCentered_ENABLED
378 200 : if (menabled) then
379 192 : call setCentered(ArrayCentered, Array, lmsize, rmsize, fill, lmfill, rmfill)
380 : else
381 44 : call setCentered(ArrayCentered, Array, fill)
382 : end if
383 :
384 390 : assertion = assertion .and. GET_LBOUND(ArrayCentered) == lbnew
385 200 : call report()
386 200 : call test%assert(assertion, SK_"Call to setCentered()/getCentered() must properly set the lower bound of the output array, with present(fill) = "//getStr(present(fill)), int(__LINE__, IK))
387 :
388 390 : assertion = assertion .and. GEN_UBOUND(ArrayCentered) == ubnew
389 200 : call report()
390 200 : call test%assert(assertion, SK_"Call to setCentered()/getCentered() must properly set the upper bound of the output array, with present(fill) = "//getStr(present(fill)), int(__LINE__, IK))
391 : #elif getCentered_ENABLED
392 200 : if (menabled) then
393 1560 : ArrayCentered = getCentered(Array, sizecentered, lmsize, rmsize, fill, lmfill, rmfill)
394 : else
395 196 : ArrayCentered = getCentered(Array, sizecentered, fill)
396 : end if
397 : #else
398 : #error "Unrecognized interface."
399 : #endif
400 :
401 400 : assertion = assertion .and. GET_SIZE(ArrayCentered) == sizenew
402 400 : call report()
403 400 : call test%assert(assertion, SK_"Call to setCentered()/getCentered() must yield an array of proper size, with present(fill) = "//getStr(present(fill)), int(__LINE__, IK))
404 :
405 2300 : assertion = assertion .and. ALL(ArrayCentered(GET_LBOUND(ArrayCentered)+getOption(0_IK,lmsize):GEN_UBOUND(ArrayCentered)-getOption(0_IK,rmsize)) IS_EQUAL ArrayCentered_ref(lbnew+getOption(0_IK,lmsize):ubnew-getOption(0_IK,rmsize)))
406 400 : call report()
407 400 : call test%assert(assertion, SK_"Call to setCentered()/getCentered() must properly set the upper bound of the output array, with present(fill) = "//getStr(present(fill)), int(__LINE__, IK))
408 :
409 400 : if (present(fill)) then
410 1150 : assertion = assertion .and. ALL(ArrayCentered(GET_LBOUND(ArrayCentered)+getOption(0_IK,lmsize):GEN_UBOUND(ArrayCentered)-getOption(0_IK,rmsize)) IS_EQUAL ArrayCentered_ref(lbnew+getOption(0_IK,lmsize):ubnew-getOption(0_IK,rmsize)))
411 200 : call report()
412 200 : call test%assert(assertion, SK_"Call to setCentered()/getCentered() must properly fill the new elements with `fill` = "//getStr(present(fill)), int(__LINE__, IK))
413 : end if
414 :
415 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
416 : ! Shrink and center array perfectly
417 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
418 :
419 400 : call reset()
420 :
421 400 : lbold = GEN_LBOLD(1_IK)
422 400 : ubold = GEN_UBOLD(4_IK)
423 400 : lbnew = GEN_LBNEW(1_IK)
424 400 : ubnew = GEN_UBNEW(2_IK)
425 400 : sizeold = ubold - lbold + 1_IK
426 400 : sizenew = ubnew - lbnew + 1_IK
427 400 : sizecentered = sizenew - getOption(0_IK,lmsize) - getOption(0_IK,rmsize)
428 400 : ALLOC_ARRAY(Array,lbold,ubold)
429 440 : ALLOC_ARRAY(ArrayCentered,lbnew,ubnew)
430 400 : ALLOC_ARRAY(ArrayCentered_ref,lbnew,ubnew)
431 : #if SK_ENABLED && D0_ENABLED
432 180 : call setUnifRand(Array, repeat(SKC_"A",len(Array)), repeat(SKC_"Z",len(Array)))
433 20 : ArrayCentered_ref(lbnew:ubnew) = getRepeat(lmfill,lmsize)//Array(lbold+1:ubold-1)//getRepeat(rmfill,rmsize)
434 : #elif SK_ENABLED && D1_ENABLED
435 100 : Array(:) = [SKC_"AA", SKC_"BB", SKC_"CC", SKC_"DD"]
436 20 : if (getOption(0_IK,lmsize) > 0_IK) then
437 272 : ArrayCentered_ref(lbnew:ubnew) = [getRepeat(lmfill,lmsize), Array(lbold+1:ubold-1), getRepeat(rmfill,rmsize)]
438 : else
439 20 : ArrayCentered_ref(lbnew:ubnew) = [Array(lbold+1:ubold-1)]
440 : end if
441 : #elif LK_ENABLED && D1_ENABLED
442 400 : Array(:) = [.true._LKC, .true._LKC, .true._LKC]
443 1300 : ArrayCentered_ref(lbnew:ubnew) = [getRepeat(lmfill,lmsize), Array(lbold+1:ubold-1), getRepeat(rmfill,rmsize)]
444 : #elif IK_ENABLED && D1_ENABLED
445 500 : Array(:) = [1_IKC, 2_IKC, 3_IKC, 4_IKC]
446 1300 : ArrayCentered_ref(lbnew:ubnew) = [getRepeat(lmfill,lmsize), Array(lbold+1:ubold-1), getRepeat(rmfill,rmsize)]
447 : #elif CK_ENABLED && D1_ENABLED
448 400 : Array(:) = [(1._CKC,-1._CKC), (2._CKC,-2._CKC), (3._CKC,-3._CKC), (4._CKC,-4._CKC)]
449 1040 : ArrayCentered_ref(lbnew:ubnew) = [getRepeat(lmfill,lmsize), Array(lbold+1:ubold-1), getRepeat(rmfill,rmsize)]
450 : #elif RK_ENABLED && D1_ENABLED
451 400 : Array(:) = [1._RKC, 2._RKC, 3._RKC, 4._RKC]
452 1040 : ArrayCentered_ref(lbnew:ubnew) = [getRepeat(lmfill,lmsize), Array(lbold+1:ubold-1), getRepeat(rmfill,rmsize)]
453 : #endif
454 400 : cbeg = lbold + abs(sizecentered - sizeold) / 2_IK
455 400 : cend = cbeg + min(sizeold,sizecentered) - 1_IK
456 :
457 : #if setCentered_ENABLED
458 200 : if (menabled) then
459 192 : call setCentered(ArrayCentered, Array, lmsize, rmsize, fill, lmfill, rmfill)
460 : else
461 44 : call setCentered(ArrayCentered, Array, fill)
462 : end if
463 :
464 :
465 390 : assertion = assertion .and. GET_LBOUND(ArrayCentered) == lbnew
466 200 : call report()
467 200 : call test%assert(assertion, SK_"Call to setCentered()/getCentered() must properly set the lower bound of the output array, with present(fill) = "//getStr(present(fill)), int(__LINE__, IK))
468 :
469 390 : assertion = assertion .and. GEN_UBOUND(ArrayCentered) == ubnew
470 200 : call report()
471 200 : call test%assert(assertion, SK_"Call to setCentered()/getCentered() must properly set the upper bound of the output array, with present(fill) = "//getStr(present(fill)), int(__LINE__, IK))
472 : #elif getCentered_ENABLED
473 200 : if (menabled) then
474 1408 : ArrayCentered = getCentered(Array, sizecentered, lmsize, rmsize, fill, lmfill, rmfill)
475 : else
476 158 : ArrayCentered = getCentered(Array, sizecentered, fill)
477 : end if
478 : #else
479 : #error "Unrecognized interface."
480 : #endif
481 :
482 400 : assertion = assertion .and. GET_SIZE(ArrayCentered) == sizenew
483 400 : call report()
484 400 : call test%assert(assertion, SK_"Call to setCentered()/getCentered() must yield an array of proper size, with present(fill) = "//getStr(present(fill)), int(__LINE__, IK))
485 :
486 1920 : assertion = assertion .and. ALL(ArrayCentered(GET_LBOUND(ArrayCentered)+getOption(0_IK,lmsize):GEN_UBOUND(ArrayCentered)-getOption(0_IK,rmsize)) IS_EQUAL ArrayCentered_ref(lbnew+getOption(0_IK,lmsize):ubnew-getOption(0_IK,rmsize)))
487 400 : call report()
488 400 : call test%assert(assertion, SK_"Call to setCentered()/getCentered() must properly set the upper bound of the output array, with present(fill) = "//getStr(present(fill)), int(__LINE__, IK))
489 :
490 400 : if (present(fill)) then
491 960 : assertion = assertion .and. ALL(ArrayCentered(GET_LBOUND(ArrayCentered)+getOption(0_IK,lmsize):GEN_UBOUND(ArrayCentered)-getOption(0_IK,rmsize)) IS_EQUAL ArrayCentered_ref(lbnew+getOption(0_IK,lmsize):ubnew-getOption(0_IK,rmsize)))
492 200 : call report()
493 200 : call test%assert(assertion, SK_"Call to setCentered()/getCentered() must properly fill the new elements with `fill` = "//getStr(present(fill)), int(__LINE__, IK))
494 : end if
495 :
496 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
497 :
498 400 : end subroutine
499 :
500 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
501 :
502 1600 : subroutine reset()
503 1600 : if (allocated(Array)) deallocate(Array)
504 1600 : if (allocated(ArrayCentered)) deallocate(ArrayCentered)
505 1600 : if (allocated(ArrayCentered_ref)) deallocate(ArrayCentered_ref)
506 1600 : end subroutine
507 :
508 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
509 :
510 : #if SK_ENABLED && D0_ENABLED
511 160 : PURE function getRepeat(fill, count) result(Array)
512 : integer(IK) , intent(in), optional :: count
513 : character(1,SKC), intent(in), optional :: fill
514 : character(:,SKC), allocatable :: Array
515 672 : Array = repeat(getOption(SKC_" ",fill), getOption(0_IK,count))
516 160 : end function
517 : #else
518 3008 : PURE function getRepeat(fill, count) result(Array)
519 : integer(IK) , intent(in), optional :: count
520 : #if SK_ENABLED && D1_ENABLED
521 : character(2,SKC), intent(in), optional :: fill
522 : character(2,SKC), allocatable :: Array(:)
523 : #elif LK_ENABLED && D1_ENABLED
524 : logical(LKC) , intent(in), optional :: fill
525 : logical(LKC) , allocatable :: Array(:)
526 : #elif IK_ENABLED && D1_ENABLED
527 : integer(IKC) , intent(in), optional :: fill
528 : integer(IKC) , allocatable :: Array(:)
529 : #elif CK_ENABLED && D1_ENABLED
530 : complex(CKC) , intent(in), optional :: fill
531 : complex(CKC) , allocatable :: Array(:)
532 : #elif RK_ENABLED && D1_ENABLED
533 : real(RKC) , intent(in), optional :: fill
534 : real(RKC) , allocatable :: Array(:)
535 : #else
536 : #error "Unrecognized interface."
537 : #endif
538 : integer(IK) :: i
539 24384 : Array = [( getOption(fill_def,fill), i = 1_IK, getOption(0_IK,count) )]
540 3008 : end function
541 : #endif
542 :
543 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
544 :
545 5600 : subroutine report()
546 5600 : if (test%traceable .and. .not. assertion) then
547 : ! LCOV_EXCL_START
548 : write(test%disp%unit,"(*(g0,:,', '))")
549 : write(test%disp%unit,"(*(g0,:,', '))") "Array ", Array
550 : write(test%disp%unit,"(*(g0,:,', '))") "ArrayCentered ", ArrayCentered
551 : write(test%disp%unit,"(*(g0,:,', '))")
552 : ! LCOV_EXCL_STOP
553 : end if
554 5600 : end subroutine
555 :
556 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
557 :
558 : #undef ALLOC_ARRAY_CENTERED
559 : #undef ALLOC_ARRAY
560 : #undef GEN_UBOUND
561 : #undef GET_LBOUND
562 : #undef GEN_UBOLD
563 : #undef GEN_UBNEW
564 : #undef GEN_LBOLD
565 : #undef GEN_LBNEW
566 : #undef IS_EQUAL
567 : #undef GET_SIZE
568 : #undef ALL
|