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 file contains procedure implementations of [pm_str](@ref pm_str).
19 : !>
20 : !> \finmain
21 : !>
22 : !> \author
23 : !> \AmirShahmoradi, September 1, 2017, 12:00 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
24 :
25 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
26 :
27 : !%%%%%%%%%%%%
28 : #if alleq_ENABLED
29 : !%%%%%%%%%%%%
30 :
31 : integer(IK) :: i, lenStr1, lenStr2
32 14 : lenStr1 = len(str1, IK)
33 14 : lenStr2 = len(str2, IK)
34 14 : allEqual = logical(lenStr1 == lenStr2, LK)
35 14 : if (allEqual) then
36 4 : allEqual = logical(str1 == str2, LK)
37 10 : elseif (lenStr1 == 1_IK) then
38 5 : do i = 1, lenStr2
39 4 : allEqual = logical(str1(1:1) == str2(i:i), LK)
40 5 : if (.not. allEqual) return
41 : end do
42 9 : elseif (lenStr2 == 1_IK) then
43 19 : do i = 1, lenStr1
44 15 : allEqual = logical(str1(i:i) == str2(1:1), LK)
45 19 : if (.not. allEqual) return
46 : end do
47 : end if
48 :
49 : !%%%%%%%%%%%%%%%%%
50 : #elif getCharSeq_ENABLED
51 : !%%%%%%%%%%%%%%%%%
52 :
53 : integer(IK) :: i
54 : do concurrent(i = 1_IK : size(charVec, kind = IK))
55 108 : str(i:i) = charVec(i)
56 : end do
57 :
58 : !%%%%%%%%%%%%%%%%%
59 : #elif getCharVec_ENABLED
60 : !%%%%%%%%%%%%%%%%%
61 :
62 : integer(IK) :: i
63 : do concurrent(i = 1_IK : len(str, kind = IK))
64 3424 : charVec(i) = str(i:i)
65 : end do
66 :
67 : !%%%%%%%%%%%%%%%%%%
68 : #elif isEndedWith_ENABLED
69 : !%%%%%%%%%%%%%%%%%%
70 :
71 : #if BSSK_ENABLED || PSSK_ENABLED
72 : #define GET_VAL(X)X%val
73 : #elif SK_ENABLED
74 : #define GET_VAL(X)X
75 : #else
76 : #error "Unrecognized interface."
77 : #endif
78 386606 : endedWith = logical(GET_VAL(str)(max(1_IK, len(GET_VAL(str), IK) - len(GET_VAL(suffix), IK) + 1_IK) : len(GET_VAL(str), IK)) == GET_VAL(suffix), LK)
79 : #undef GET_VAL
80 :
81 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
82 : #elif getMinLoc_ENABLED || getMinVal_ENABLED || getMaxLoc_ENABLED || getMaxVal_ENABLED
83 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
84 :
85 : #if getMinLoc_ENABLED || getMinVal_ENABLED
86 : #define EXTREMUM_LOC minLoc
87 : #define EXTREMUM_VAL minVal
88 : #define COMPARES_WITH <
89 : #elif getMaxLoc_ENABLED || getMaxVal_ENABLED
90 : #define EXTREMUM_LOC maxLoc
91 : #define EXTREMUM_VAL maxVal
92 : #define COMPARES_WITH >
93 : #else
94 : #error "Unrecognized interface."
95 : #endif
96 : integer(IK) :: i
97 : #if getMinLoc_ENABLED || getMaxLoc_ENABLED
98 : #define getLoc_ENABLED 1
99 : character(1,SKC) :: EXTREMUM_VAL
100 : EXTREMUM_LOC = 1_IK
101 : #elif !(getMinVal_ENABLED || getMaxVal_ENABLED)
102 : #error "Unrecognized interface."
103 : #endif
104 41 : CHECK_ASSERTION(__LINE__, 0_IK < len(str, IK), SK_": The condition `0 < len(str)` must hold. len(str) = "//getStr(len(str, IK)))
105 : #if Masked_ENABLED
106 24 : CHECK_ASSERTION(__LINE__, size(mask, 1, IK) == len(str, IK), SK_": The condition `size(mask) == len(str)` must hold. size(mask), len(str) = "//getStr([size(mask, 1, IK), len(str, IK)]))
107 : #endif
108 : #if getLoc_ENABLED
109 27 : if (present(back)) then
110 14 : if (back) then
111 10 : EXTREMUM_VAL = str(len(str, kind = IK) : len(str, kind = IK))
112 186 : do i = len(str, kind = IK) - 1_IK, 1_IK, -1_IK
113 : if ( & ! LCOV_EXCL_LINE
114 : #if Masked_ENABLED
115 : mask(i) .and. & ! LCOV_EXCL_LINE
116 : #endif
117 10 : str(i:i) COMPARES_WITH EXTREMUM_VAL) then
118 : EXTREMUM_VAL = str(i:i)
119 : #if getLoc_ENABLED
120 : EXTREMUM_LOC = i
121 : #endif
122 : end if
123 : end do
124 10 : return
125 : end if
126 : end if
127 : #endif
128 31 : EXTREMUM_VAL = str(1:1)
129 618 : do i = 2_IK, len(str, kind = IK)
130 : if (& ! LCOV_EXCL_LINE
131 : #if Masked_ENABLED
132 : mask(i) .and. & ! LCOV_EXCL_LINE
133 : #endif
134 31 : str(i:i) COMPARES_WITH EXTREMUM_VAL) then
135 30 : EXTREMUM_VAL = str(i:i)
136 : #if getLoc_ENABLED
137 : EXTREMUM_LOC = i
138 : #endif
139 : end if
140 : end do
141 :
142 : !%%%%%%%%%%%%%%%%%%%
143 : #elif getTrimmedTZ_ENABLED
144 : !%%%%%%%%%%%%%%%%%%%
145 :
146 : integer(IK) :: i, j, lenstr, dotpos, tzstart
147 : logical(LK) :: isDigitSeq, leadingDigitIsMissing
148 386 : lenstr = len(str, IK)
149 386 : allocate(character(lenstr,SKC) :: strt)
150 386 : if (lenstr > 0_IK) then
151 : i = 1_IK
152 : j = 1_IK
153 386 : strt(j:j) = str(i:i)
154 : loopOverString: do
155 136212 : if (str(i:i) == SKC_".") then
156 : dotpos = i
157 : loopOverDigitSeq: do
158 2013 : i = i + 1_IK
159 2013 : if (str(i:i) == SKC_"0") then
160 : tzstart = i
161 : loopOverTZ: do
162 16739 : if (i == lenstr) exit loopOverTZ
163 16695 : i = i + 1_IK
164 16695 : if (str(i:i) == SKC_"0") cycle loopOverTZ
165 1878 : isDigitSeq = logical(index(SKC_"123456789", str(i:i)) > 0, LK) ! isCharDigit(str(i:i))
166 1878 : if (isDigitSeq) then ! the zeros are significant, keep them.
167 9 : j = j + 1_IK
168 9 : strt(j:j+i-tzstart) = str(tzstart:i)
169 : j = j + i - tzstart
170 9 : if (i == lenstr) exit loopOverString
171 : cycle loopOverDigitSeq
172 : end if
173 14861 : exit loopOverTZ
174 : end do loopOverTZ
175 1913 : if (tzstart == dotpos + 1_IK) then
176 : ! Ignore the all trailing zeros only if there is digit before the decimal point.
177 : ! Otherwise, keep a zero as the leading digit before the decimal point.
178 1873 : leadingDigitIsMissing = logical(dotpos == 1_IK, LK)
179 1873 : if (.not. leadingDigitIsMissing) leadingDigitIsMissing = logical(index(SKC_"0123456789", str(dotpos-1:dotpos-1)) == 0, LK)
180 1873 : if (leadingDigitIsMissing) then
181 4 : strt(j : j + 1) = SKC_"0."
182 : j = j + 1_IK
183 : !j = j + 1_IK
184 : !strt(j:j+i-tzstart) = str(tzstart:i)
185 : !j = j + i - tzstart
186 : end if
187 : end if
188 1913 : if (i == lenstr) exit loopOverString
189 : else
190 91 : isDigitSeq = logical(index(SKC_"123456789", str(i:i)) > 0, LK) ! isCharDigit(str(i:i))
191 : end if
192 1960 : j = j + 1_IK
193 1960 : strt(j:j) = str(i:i)
194 1960 : if (i == lenstr) exit loopOverString
195 1957 : if (isDigitSeq) cycle loopOverDigitSeq
196 77 : cycle loopOverString
197 : end do loopOverDigitSeq
198 : end if
199 134285 : if (i == lenstr) exit loopOverString
200 133946 : i = i + 1_IK
201 133946 : j = j + 1_IK
202 134332 : strt(j:j) = str(i:i)
203 : end do loopOverString
204 386 : strt = strt(1:j)
205 : else
206 0 : strt = SKC_""
207 : end if
208 :
209 : !%%%%%%%%%%%%%%%%%%%
210 : #elif getLenIndent_ENABLED
211 : !%%%%%%%%%%%%%%%%%%%
212 :
213 : integer(IK) :: lenStr, lenPattern, ibeg, iend
214 4092 : lenPattern = len(pattern, IK)
215 4092 : lenStr = len(str, IK)
216 4092 : if (lenStr == 0_IK .or. lenPattern == 0_IK .or. lenStr < lenPattern) then
217 : lenIndent = 0_IK
218 : return
219 : else
220 : ! Determine how many instances of pattern exist at the beginning of the string, and compute the indentation length.
221 : lenIndent = 0_IK
222 : iend = 0_IK
223 4039 : do
224 : ibeg = iend
225 8118 : iend = ibeg + lenPattern
226 8118 : if (iend > lenStr) return
227 8118 : if (str(ibeg + 1_IK : iend) /= pattern) return
228 4039 : lenIndent = lenIndent + 1_IK
229 : end do
230 : end if
231 :
232 : !%%%%%%%%%%%%%%%%%%%%
233 : #elif getStrWrapped_ENABLED
234 : !%%%%%%%%%%%%%%%%%%%%
235 :
236 : !logical(LK) :: isNewLine ! new line indicator
237 : logical(LK) :: isNewPar ! new paragraph indicator
238 : integer(IK) :: i, j, iwidth
239 : integer(IK) :: def_width, def_maxwidth
240 : integer(IK) :: widthCurrent, maxWidthCurrent
241 : integer(IK) :: lenStr, lenStrWrapped, lenStrWrappedTemp, countIndentPatternOld, countIndentPatternNew
242 : integer(IK) :: lenPrefix, lenIndent, lenIndentPattern, lenBreak, lenNewLine, lenLineFeed, lenPrefixIndent
243 1604 : character(:,SKC), allocatable :: strWrappedTemp, prefix_def, prefixIndent_def, indentPattern_def, break_def, newline_def, linefeed_def
244 :
245 1604 : lenStr = len(str, kind = IK)
246 1604 : if (lenStr == 0_IK) then
247 0 : strWrapped = SKC_""
248 : return
249 : end if
250 :
251 : ! Determine the prefix to each line.
252 :
253 1604 : if (present(prefix)) then
254 1503 : lenPrefix = len(prefix, IK)
255 1503 : prefix_def = prefix
256 : else
257 : lenPrefix = 0_IK
258 101 : prefix_def = SKC_""
259 : end if
260 :
261 : ! Determine the indentation pattern.
262 :
263 1604 : if (present(indent)) then
264 3 : lenIndentPattern = len(indent, IK)
265 3 : indentPattern_def = indent
266 : else
267 : lenIndentPattern = 1_IK
268 1601 : indentPattern_def = SKC_" "
269 : end if
270 :
271 : ! Determine the allowed break characters.
272 :
273 1604 : if (present(break)) then
274 1 : lenBreak = len(break, kind = IK)
275 1 : break_def = break
276 : else
277 : lenBreak = 1_IK
278 1603 : break_def = SKC_" "
279 : end if
280 :
281 : ! Determine the newline.
282 :
283 1604 : if (present(newline)) then
284 36 : lenNewLine = len(newline, IK)
285 36 : newline_def = newline
286 : else
287 : lenNewLine = 1_IK
288 1568 : newline_def = new_line(SKC_"a")
289 : end if
290 :
291 : ! Determine the linefeed.
292 :
293 1604 : if (present(linefeed)) then
294 1498 : lenLineFeed = len(linefeed, IK)
295 1498 : linefeed_def = linefeed
296 : else
297 : lenLineFeed = lenNewLine
298 106 : linefeed_def = newline_def
299 : end if
300 :
301 : ! Set the width of each line.
302 :
303 1604 : if (present(width)) then
304 1534 : def_width = max(width, 1_IK)
305 : else
306 : def_width = 132_IK
307 : end if
308 :
309 : ! Set the max-width of each line.
310 :
311 1604 : if (present(maxwidth)) then
312 0 : def_maxwidth = maxwidth
313 : else
314 : def_maxwidth = (huge(0_IK) - 1_IK) / 2_IK
315 : end if
316 1604 : def_maxwidth = max(def_maxwidth, def_width) ! `def_maxwidth > def_width` must always hold.
317 :
318 : ! Construct the wrapped array.
319 :
320 1604 : lenStrWrapped = lenStr * 3_IK / 2_IK ! best guess
321 1604 : allocate(character(lenStrWrapped,SKC) :: strWrappedTemp)
322 :
323 1604 : if (lenIndentPattern == 0_IK) then
324 0 : lenPrefixIndent = lenPrefix
325 0 : prefixIndent_def = prefix_def
326 : maxWidthCurrent = def_maxwidth
327 : widthCurrent = def_width
328 : lenIndent = 0_IK
329 : end if
330 :
331 : i = 0_IK ! counter for `str`.
332 : j = 0_IK ! counter for `strWrapped`.
333 : isNewPar = .true._LK
334 : !isNewLine = .false._LK
335 : countIndentPatternOld = -1_IK
336 : loopOverLines: do
337 :
338 : ! Set the prefix + indentation for the new paragraph.
339 :
340 12969 : if (isNewPar .and. lenIndentPattern > 0_IK) then
341 4087 : countIndentPatternNew = getLenIndent(str(i + 1 : lenStr), indentPattern_def)
342 4087 : if (countIndentPatternNew /= countIndentPatternOld) then ! This must be true at least the first time reaching here, which requires `countIndentPatternOld < 0_IK` for the first time.
343 3356 : lenIndent = lenIndentPattern * countIndentPatternNew
344 7216 : prefixIndent_def = prefix_def//repeat(indentPattern_def, countIndentPatternNew)
345 3356 : lenPrefixIndent = len(prefixIndent_def, IK)
346 : countIndentPatternOld = countIndentPatternNew
347 3356 : if (lenIndentPattern * countIndentPatternNew > def_width) then
348 1 : widthCurrent = lenIndentPattern * countIndentPatternNew + 1_IK
349 : maxWidthCurrent = widthCurrent
350 : else
351 : maxWidthCurrent = def_maxwidth
352 : widthCurrent = def_width
353 : end if
354 : end if
355 : end if
356 : !write(*,*) "widthCurrent, def_width, maxWidthCurrent, def_maxwidth", widthCurrent, def_width, maxWidthCurrent, def_maxwidth
357 :
358 : ! Indent the line/paragraph.
359 :
360 : iwidth = lenIndent ! width of the current line.
361 :
362 12969 : if (lenPrefixIndent > 0_IK) then
363 12169 : if (j + lenPrefixIndent > lenStrWrapped) call resizeStrWrapped()
364 12169 : strWrappedTemp(j + 1 : j + lenPrefixIndent) = prefixIndent_def
365 12169 : j = j + lenPrefixIndent
366 12169 : if (isNewPar) then
367 3896 : i = i + lenIndent ! The first indent is already taken into account in the above, so skip it.
368 3896 : if (i > lenStr) exit
369 : end if
370 : end if
371 :
372 : ! Wrap the string until reaching a newline.
373 :
374 1604 : loopOverChars: do
375 :
376 : ! Check for the presence of a newline.
377 :
378 1131560 : if (i + lenNewLine <= lenStr) then
379 1129920 : if (str(i + 1 : i + lenNewLine) == newline_def) then
380 3651 : if (j + lenLineFeed > lenStrWrapped) call resizeStrWrapped()
381 3651 : strWrappedTemp(j + 1 : j + lenLineFeed) = linefeed_def
382 : j = j + lenLineFeed
383 : i = i + lenNewLine
384 : !isNewLine = .true._LK
385 3651 : if (i + lenNewLine <= lenStr) then ! Check for new paragraph
386 3638 : if (str(i + 1 : i + lenNewLine) == newline_def) then
387 2470 : if (j + lenPrefix + lenLineFeed > lenStrWrapped) call resizeStrWrapped()
388 2470 : strWrappedTemp(j + 1 : j + lenPrefix) = prefix_def
389 : j = j + lenPrefix
390 2470 : strWrappedTemp(j + 1 : j + lenLineFeed) = linefeed_def
391 : j = j + lenLineFeed
392 : i = i + lenNewLine
393 : isNewPar = .true._LK
394 : cycle loopOverLines
395 : end if
396 : isNewPar = .false._LK
397 : end if
398 : cycle loopOverLines
399 : end if
400 : end if
401 :
402 : ! Copy the next character to `strWrapped`.
403 :
404 1127909 : i = i + 1_IK; if (i > lenStr) exit loopOverLines
405 1126305 : j = j + 1_IK; if (j > lenStrWrapped) call resizeStrWrapped()
406 1126305 : strWrappedTemp(j : j) = str(i : i)
407 1126305 : iwidth = iwidth + 1_IK
408 1126305 : if (iwidth >= widthCurrent) then
409 47798 : if (lenBreak == 0_IK .or. index(break_def, str(i:i), kind = IK) > 0_IK .or. iwidth == maxWidthCurrent) then
410 7714 : if (j + lenLineFeed > lenStrWrapped) call resizeStrWrapped()
411 7714 : strWrappedTemp(j + 1 : j + lenLineFeed) = linefeed_def
412 : j = j + lenLineFeed
413 : !isNewLine = .false._LK
414 : isNewPar = .false._LK
415 : cycle loopOverLines
416 : end if
417 : end if
418 :
419 : end do loopOverChars
420 :
421 : end do loopOverLines
422 :
423 3208 : strWrapped = strWrappedTemp(1 : j)
424 :
425 : contains
426 :
427 71 : subroutine resizeStrWrapped()
428 71 : lenStrWrappedTemp = max(lenStrWrapped * 3_IK / 2_IK, lenStrWrapped + lenPrefixIndent)
429 71 : allocate(character(lenStrWrappedTemp,SKC) :: strWrapped)
430 71 : strWrapped(1_IK:lenStrWrapped) = strWrappedTemp(1_IK:lenStrWrapped)
431 71 : call move_alloc(strWrapped, strWrappedTemp)
432 : #if __GFORTRAN__
433 71 : if (allocated(strWrapped)) deallocate(strWrapped) ! Bypass gfortran <12 bug with Heap allocation.
434 : #endif
435 71 : lenStrWrapped = lenStrWrappedTemp
436 71 : end subroutine
437 : #else
438 : !%%%%%%%%%%%%%%%%%%%%%%%%
439 : #error "Unrecognized interface."
440 : !%%%%%%%%%%%%%%%%%%%%%%%%
441 : #endif
442 : #undef getLoc_ENABLED
443 : #undef COMPARES_WITH
444 : #undef EXTREMUM_LOC
445 : #undef EXTREMUM_VAL
|