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_distPiwiPoweto](@ref pm_distPiwiPoweto).
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 : use pm_complexAbs, only: abs, log, operator(<), operator(<=)
29 : #elif !RK_ENABLED
30 : #error "Unrecognized interface."
31 : #endif
32 :
33 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
34 : #if getPiwiPowetoLogPDF_ENABLED || setPiwiPowetoLogPDF_ENABLED
35 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
36 :
37 : integer(IK) :: i, numComp
38 :
39 : real(RKC) , parameter :: TOL = sqrt(epsilon(0._RKC))
40 : real(RKC) , parameter :: LOG_HUGE = log(huge(0._RKC))
41 : real(RKC) , allocatable :: logLimX(:), alpha(:), logPDFNF(:)
42 : real(RKC) :: logx, output, output_ref, diff
43 :
44 8 : assertion = .true._LK
45 :
46 168 : do i = 1_IK, 20_IK
47 :
48 160 : call setUnifRand(numComp, 1_IK, 5_IK)
49 : !if (getUnifRand()) then
50 : ! logLimX = [getUnifRand(sqrt(epsilon(0._RKC)), 10._RKC, numComp), LOG_HUGE]
51 : !else
52 : ! logLimX = getUnifRand(sqrt(epsilon(0._RKC)), 10._RKC, numComp + 1)
53 : !end if
54 : !logLimX = getUnifRand(log(sqrt(epsilon(0._RKC))), sqrt(LOG_HUGE), numComp + 1)
55 : !logLimX = getUnifRand(log(sqrt(epsilon(0._RKC))), 10._RKC, numComp + 1)
56 953 : logLimX = getUnifRand(-5._RKC, 8._RKC, numComp + 1)
57 160 : call setSorted(logLimX)
58 160 : if (logLimX(size(logLimX)) == LOG_HUGE) then
59 0 : alpha = [getUnifRand(-2._RKC, 2._RKC), getUnifRand(-2._RKC, 2._RKC, numComp - 2_IK), getUnifRand(-8._RKC, -2._RKC)]
60 : else
61 793 : alpha = getUnifRand(-2._RKC, 2._RKC, numComp)
62 : end if
63 160 : if (getUnifRand() .and. size(alpha) > 1) alpha(getUnifRand(1, size(alpha) - 1)) = -1._RKC
64 160 : if (getUnifRand() .and. getUnifRand() .and. size(alpha) > 1) alpha(getUnifRand(1, size(alpha) - 1)) = 0._RKC
65 160 : if (getUnifRand()) then
66 76 : alpha(size(alpha)) = -1._RKC
67 84 : elseif (getUnifRand()) then
68 42 : alpha(size(alpha)) = 0._RKC
69 : end if
70 793 : logPDFNF = getPiwiPowetoLogPDFNF(alpha, logLimX)
71 :
72 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
73 :
74 : #if getPiwiPowetoLogPDF_ENABLED
75 :
76 80 : logx = getUnifRand(logLimX(1), logLimX(size(logLimX)))
77 80 : output = getPiwiPowetoLogPDF(logx, alpha, logLimX)
78 40 : call setPiwiPowetoLogPDF(output_ref, logx, alpha, logLimX, logPDFNF)
79 80 : call report()
80 80 : call test%assert(assertion, SK_"@getPiwiPowetoLogPDF(): The PDF of the PiwiPoweto distribution must be computed correctly when `logPDFNF` is missing.", int(__LINE__, IK))
81 :
82 80 : logx = getUnifRand(logLimX(1), logLimX(size(logLimX)))
83 80 : output = getPiwiPowetoLogPDF(logx, alpha, logLimX, logPDFNF)
84 40 : call setPiwiPowetoLogPDF(output_ref, logx, alpha, logLimX, logPDFNF)
85 80 : call report()
86 244 : call test%assert(assertion, SK_"@getPiwiPowetoLogPDF(): The PDF of the PiwiPoweto distribution must be computed correctly when `logPDFNF` is present.", int(__LINE__, IK))
87 :
88 : #elif setPiwiPowetoLogPDF_ENABLED
89 :
90 4 : block
91 :
92 : use pm_distPiwiPoweto, only: getPiwiPowetoCDF
93 : use pm_arraySearch, only: getBin
94 : use pm_val2str, only: getStr
95 : character(255, SK) :: msg
96 : integer(IK) :: j
97 :
98 : !output_ref = 1._RKC
99 80 : msg = SK_" "
100 : !output = 0._RKC
101 248 : do j = 1, numComp - 1
102 168 : output_ref = getPiwiPowetoCDF(logLimX(j + 1), alpha, logLimX) - getPiwiPowetoCDF(logLimX(j), alpha, logLimX)
103 168 : assertion = assertion .and. .not. isFailedQuad(getPiwiPowetoPDF, exp(logLimX(j)), exp(logLimX(j + 1)), output, reltol = TOL, msg = msg)
104 168 : call report()
105 248 : call test%assert(assertion, SK_"@setPiwiPowetoLogPDF(): The integral of the PiwiPoweto distribution must be computed correctly.", int(__LINE__, IK))
106 : !call test%assert(assertion, SK_"@setPiwiPowetoLogPDF(): The integration of the PiwiPoweto distribution must not fail. j, msg: "//getStr(j)//SK_", "//trim(msg), int(__LINE__, IK))
107 : !output = output + logx
108 : end do
109 :
110 : !call report()
111 : !call test%assert(assertion, SK_"@setPiwiPowetoLogPDF(): The integral of the PiwiPoweto distribution must be computed correctly.", int(__LINE__, IK))
112 :
113 80 : logx = getUnifRand(logLimX(1), logLimX(size(logLimX)))
114 40 : call setPiwiPowetoLogPDF(output_ref, logx, alpha, logLimX, logPDFNF)
115 80 : call setPiwiPowetoLogPDF(output, logx, alpha, getBin(logLimX, logx), logPDFNF)
116 80 : call report()
117 80 : call test%assert(assertion, SK_"@setPiwiPowetoLogPDF(): The PDF of the PiwiPoweto distribution must be computed correctly when the value bin is specified.", int(__LINE__, IK))
118 :
119 : end block
120 :
121 : #else
122 : #error "Unrecognized interface."
123 : #endif
124 :
125 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
126 :
127 : end do
128 :
129 : contains
130 :
131 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
132 :
133 : #if setPiwiPowetoLogPDF_ENABLED
134 18522 : function getPiwiPowetoPDF(x) result(pdf)
135 : real(RKC), intent(in) :: x
136 : real(RKC) :: pdf
137 18522 : call setPiwiPowetoLogPDF(pdf, log(x), alpha, logLimX, logPDFNF)
138 18522 : pdf = exp(pdf)
139 : !if (pdf > -huge(pdf)**0.9) then
140 : ! pdf = exp(pdf)
141 : !else
142 : ! pdf = 0._RKC
143 : !end if
144 18522 : end function
145 : #elif !getPiwiPowetoLogPDF_ENABLED
146 : #error "Unrecognized interface."
147 : #endif
148 :
149 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
150 :
151 408 : subroutine report()
152 408 : diff = abs(output - output_ref)
153 408 : assertion = assertion .and. diff <= max(numComp * abs(output_ref) * TOL, TOL)
154 408 : if (test%traceable .and. .not. assertion) then
155 : ! LCOV_EXCL_START
156 : write(test%disp%unit,"(*(g0,:,', '))")
157 : write(test%disp%unit,"(*(g0,:,', '))") "numComp ", numComp
158 : write(test%disp%unit,"(*(g0,:,', '))") "size(alpha) ", size(alpha)
159 : write(test%disp%unit,"(*(g0,:,', '))") "size(logLimX) ", size(logLimX)
160 : write(test%disp%unit,"(*(g0,:,', '))") "logPDFNF ", logPDFNF
161 : write(test%disp%unit,"(*(g0,:,', '))") "logLimX ", logLimX
162 : write(test%disp%unit,"(*(g0,:,', '))") "alpha ", alpha
163 : write(test%disp%unit,"(*(g0,:,', '))") "diff ", diff
164 : write(test%disp%unit,"(*(g0,:,', '))") "TOL ", TOL
165 : write(test%disp%unit,"(*(g0,:,', '))") "output ", output
166 : write(test%disp%unit,"(*(g0,:,', '))") "output_ref ", output_ref
167 : write(test%disp%unit,"(*(g0,:,', '))")
168 : ! LCOV_EXCL_STOP
169 : end if
170 408 : end subroutine
171 :
172 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
173 :
174 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
175 : #elif getPiwiPowetoCDF_ENABLED || setPiwiPowetoCDF_ENABLED
176 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
177 :
178 : integer(IK) :: i, numComp
179 : real(RKC) , parameter :: TOL = sqrt(epsilon(0._RKC))
180 : real(RKC) , allocatable :: logLimX(:), cumSumArea(:), alpha(:), logPDFNF(:)
181 : real(RKC) :: output, output_ref, diff, lb, ub
182 :
183 8 : assertion = .true._LK
184 :
185 808 : do i = 1_IK, 100_IK
186 :
187 800 : call setUnifRand(numComp, 1_IK, 5_IK)
188 : !if (getUnifRand()) then
189 : ! logLimX = [getUnifRand(sqrt(epsilon(0._RKC)), 10._RKC, numComp), log(huge(0._RKC))]
190 : !else
191 : ! logLimX = getUnifRand(sqrt(epsilon(0._RKC)), 10._RKC, numComp + 1)
192 : !end if
193 4765 : logLimX = getUnifRand(-5._RKC, 5._RKC, numComp + 1)
194 800 : call setSorted(logLimX)
195 800 : if (numComp > 1_IK) then
196 4398 : alpha = [getUnifRand(-1.01_RKC, 3._RKC), getUnifRand(-3._RKC, 3._RKC, numComp - 2_IK), getUnifRand(-8._RKC, -2._RKC)]
197 : else
198 498 : alpha = [getUnifRand(-8._RKC, -2._RKC)]
199 : end if
200 800 : if (getUnifRand() .and. size(alpha) > 1) alpha(getUnifRand(1, size(alpha) - 1)) = -1._RKC
201 800 : if (getUnifRand() .and. getUnifRand() .and. size(alpha) > 1) alpha(getUnifRand(1, size(alpha) - 1)) = 0._RKC
202 800 : if (size(logLimX) > size(alpha)) then
203 800 : if (getUnifRand()) then
204 434 : alpha(size(alpha)) = -1._RKC
205 366 : elseif (getUnifRand()) then
206 190 : alpha(size(alpha)) = 0._RKC
207 : end if
208 : end if
209 1600 : if (allocated(cumSumArea)) deallocate(cumSumArea); allocate(cumSumArea, mold = logLimX)
210 3965 : logPDFNF = getPiwiPowetoLogPDFNF(alpha, logLimX, cumSumArea)
211 :
212 : ! ub = merge(logLimX(size(logLimX)), logLimX(size(logLimX)), size(logLimX,1,IK) > numComp)
213 800 : lb = getUnifRand(logLimX(1), logLimX(size(logLimX)))
214 800 : ub = exp(getUnifRand(lb, logLimX(size(logLimX))))
215 800 : lb = exp(lb)
216 :
217 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
218 :
219 : #if getPiwiPowetoCDF_ENABLED
220 :
221 4 : block
222 :
223 : real(RKC) :: logx
224 :
225 400 : logx = merge(getUnifRand(logLimX(1), logLimX(size(logLimX))), getUnifRand(logLimX(1), logLimX(size(logLimX))), getUnifRand())
226 400 : output = getPiwiPowetoCDF(logx, alpha, logLimX)
227 200 : call setPiwiPowetoCDF(output_ref, logx, alpha, logLimX, logPDFNF, cumSumArea)
228 400 : call report()
229 400 : call test%assert(assertion, SK_"@getPiwiPowetoCDF(): The CDF of the PiwiPoweto distribution must be computed correctly when `logPDFNF` is missing.", int(__LINE__, IK))
230 :
231 400 : logx = merge(getUnifRand(logLimX(1), logLimX(size(logLimX))), getUnifRand(logLimX(1), logLimX(size(logLimX))), getUnifRand())
232 400 : output = getPiwiPowetoCDF(logx, alpha, logLimX, logPDFNF, cumSumArea)
233 200 : call setPiwiPowetoCDF(output_ref, logx, alpha, logLimX, logPDFNF, cumSumArea)
234 400 : call report()
235 1200 : call test%assert(assertion, SK_"@getPiwiPowetoCDF(): The CDF of the PiwiPoweto distribution must be computed correctly when `logPDFNF` is present.", int(__LINE__, IK))
236 :
237 : end block
238 :
239 : #elif setPiwiPowetoCDF_ENABLED
240 :
241 4 : block
242 :
243 : use pm_arraySearch, only: getBin
244 : use pm_val2str, only: getStr
245 : real(RKC) :: upperCDF, lowerCDF
246 : character(255, SK) :: msg
247 : integer(IK) :: j
248 :
249 : ! Validate the CDF on the boundary points.
250 :
251 1627 : do j = 1, numComp
252 :
253 1227 : output_ref = 0._RKC
254 1227 : if (j > 1_IK) output_ref = cumSumArea(j)
255 :
256 1227 : call setPiwiPowetoCDF(output, logLimX(j), alpha, logLimX, logPDFNF, cumSumArea)
257 1227 : call report()
258 1227 : call test%assert(assertion, SK_"@setPiwiPowetoCDF(): The CDF at the component boundary values must equal the corresponding elements of `cumSumArea`.", int(__LINE__, IK))
259 :
260 1227 : call setPiwiPowetoCDF(output, logLimX(j), alpha, logLimX, logPDFNF, cumSumArea, bin = j)
261 1227 : call report()
262 4081 : call test%assert(assertion, SK_"@setPiwiPowetoCDF(): The CDF at the component boundary values must equal the corresponding elements of `cumSumArea`.", int(__LINE__, IK))
263 :
264 : end do
265 :
266 : ! The CDF over the entire support must equal unity.
267 :
268 400 : output_ref = 1._RKC
269 :
270 400 : call setPiwiPowetoCDF(output, logLimX(size(logLimX)), alpha, logLimX, logPDFNF, cumSumArea)
271 400 : call report()
272 400 : call test%assert(assertion, SK_"@setPiwiPowetoCDF(): The CDF over the entire support must equal unity.", int(__LINE__, IK))
273 :
274 400 : call setPiwiPowetoCDF(output, logLimX(size(logLimX)), alpha, logLimX, logPDFNF, cumSumArea, bin = numComp)
275 400 : call report()
276 400 : call test%assert(assertion, SK_"@setPiwiPowetoCDF(): The CDF at the component boundary values must equal the corresponding elements of `cumSumArea` with `bin` present.", int(__LINE__, IK))
277 :
278 400 : msg = SK_" "
279 :
280 : ! Segments of CDF must be computed correctly.
281 :
282 400 : j = getBin(logLimX, log(lb))
283 400 : assertion = assertion .and. .not. isFailedQuad(getPiwiPowetoPDF, exp(logLimX(j)), lb, lowerCDF, reltol = TOL, msg = msg)
284 4454 : call test%assert(assertion, SK_"@setPiwiPowetoCDF(): The integration of the PiwiPoweto distribution must not fail. lb, ub, exp(logLimX), msg: "//getStr([lb, ub, exp(logLimX)])//SK_", "//trim(msg), int(__LINE__, IK))
285 400 : output_ref = cumSumArea(j) + lowerCDF
286 :
287 400 : j = getBin(logLimX, log(ub))
288 400 : assertion = assertion .and. .not. isFailedQuad(getPiwiPowetoPDF, exp(logLimX(j)), ub, upperCDF, reltol = TOL, msg = msg)
289 4454 : call test%assert(assertion, SK_"@setPiwiPowetoCDF(): The integration of the PiwiPoweto distribution must not fail. lb, ub, exp(logLimX), msg: "//getStr([lb, ub, exp(logLimX)])//SK_", "//trim(msg), int(__LINE__, IK))
290 400 : output_ref = cumSumArea(j) + upperCDF - output_ref
291 :
292 400 : call setPiwiPowetoCDF(lowerCDF, log(lb), alpha, logLimX, logPDFNF, cumSumArea)
293 400 : call setPiwiPowetoCDF(upperCDF, log(ub), alpha, logLimX, logPDFNF, cumSumArea)
294 400 : output = upperCDF - lowerCDF
295 400 : call report()
296 400 : call test%assert(assertion, SK_"@setPiwiPowetoCDF(): The CDF between two arbitrary lower and upper limits must be computed correctly.", int(__LINE__, IK))
297 :
298 400 : call setPiwiPowetoCDF(lowerCDF, log(lb), alpha, logLimX, logPDFNF, cumSumArea, bin = getBin(logLimX, log(lb)))
299 400 : call setPiwiPowetoCDF(upperCDF, log(ub), alpha, logLimX, logPDFNF, cumSumArea, bin = getBin(logLimX, log(ub)))
300 400 : output = upperCDF - lowerCDF
301 400 : call report()
302 400 : call test%assert(assertion, SK_"@setPiwiPowetoCDF(): The CDF between two arbitrary lower and upper limits must be computed correctly with `bin` present.", int(__LINE__, IK))
303 :
304 : end block
305 :
306 : #else
307 : #error "Unrecognized interface."
308 : #endif
309 :
310 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
311 :
312 : end do
313 :
314 : contains
315 :
316 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
317 :
318 : #if setPiwiPowetoCDF_ENABLED
319 65856 : function getPiwiPowetoPDF(x) result(pdf)
320 : use pm_distPiwiPoweto, only: setPiwiPowetoLogPDF
321 : real(RKC), intent(in) :: x
322 : real(RKC) :: pdf
323 65856 : call setPiwiPowetoLogPDF(pdf, log(x), alpha, logLimX, logPDFNF)
324 65856 : pdf = exp(pdf)
325 65856 : end function
326 : #endif
327 :
328 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
329 :
330 4854 : subroutine report()
331 4854 : diff = abs(output - output_ref)
332 4854 : assertion = assertion .and. diff <= max(numComp * abs(output_ref) * TOL, TOL)
333 4854 : if (test%traceable .and. .not. assertion) then
334 : ! LCOV_EXCL_START
335 : write(test%disp%unit,"(*(g0,:,', '))")
336 : write(test%disp%unit,"(*(g0,:,', '))") "numComp ", numComp
337 : write(test%disp%unit,"(*(g0,:,', '))") "size(alpha) ", size(alpha)
338 : write(test%disp%unit,"(*(g0,:,', '))") "size(logLimX) ", size(logLimX)
339 : write(test%disp%unit,"(*(g0,:,', '))") "cumSumArea ", cumSumArea
340 : write(test%disp%unit,"(*(g0,:,', '))") "logPDFNF ", logPDFNF
341 : write(test%disp%unit,"(*(g0,:,', '))") "logLimX ", logLimX
342 : write(test%disp%unit,"(*(g0,:,', '))") "alpha ", alpha
343 : write(test%disp%unit,"(*(g0,:,', '))") "diff ", diff
344 : write(test%disp%unit,"(*(g0,:,', '))") "TOL ", TOL
345 : write(test%disp%unit,"(*(g0,:,', '))") "output ", output
346 : write(test%disp%unit,"(*(g0,:,', '))") "output_ref ", output_ref
347 : write(test%disp%unit,"(*(g0,:,', '))") "lb, ub ", lb, ub
348 : write(test%disp%unit,"(*(g0,:,', '))")
349 : ! LCOV_EXCL_STOP
350 : end if
351 4854 : end subroutine
352 :
353 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
354 :
355 : #else
356 : !%%%%%%%%%%%%%%%%%%%%%%%%
357 : #error "Unrecognized interface."
358 : !%%%%%%%%%%%%%%%%%%%%%%%%
359 : #endif
|