Generate and return the input scalar string wrapped within the specified width while optionally prefixing each line and respecting the initial indentation in all subsequent lines.
More...
Generate and return the input scalar string wrapped within the specified width while optionally prefixing each line and respecting the initial indentation in all subsequent lines.
- Parameters
-
[in] | str | : The input scalar of type character of kind any supported by the processor (e.g., SK, SKA, SKD , or SKU) that is to be wrapped.
|
[in] | prefix | : The input scalar of the same type and kind as the input str of arbitrary length type parameter, whose contents will be prefixed to each new line (including the first line) before an indentation occurs.
Note that prefix does not count toward filling the specified width or maxwidth of each wrapped line.
(optional, default = "" ) |
[in] | indent | : The input scalar of the same type and kind as the input str of arbitrary length type parameter, whose contents are the pattern (a set of characters) to search for at the beginning of each paragraph in the input str .
The repetitions of the indent at the beginning of each paragraph (including the beginning of str ) collectively construct the indentation to be applied to all wrapped lines in the specific paragraph within str .
For example, the following choice of input arguments, strWrapped = getStrWrapped(str = "~+~+This is a test.", indent = "~+")
will lead to prefixing each wrapped line in the output with ~+~+ .
If indent is as an empty string, then no search for indentation character will be made and no line will be indented.
(optional, default = " " . If indent = "" , then no indentation of the lines will occur.) |
[in] | break | : The input scalar of the same type and kind as the input str of arbitrary length type parameter, whose contents are the set of characters at which breaking the input str is allowed.
Setting break = "" will lead to wrapping the input str at any arbitrary character as soon as the length of a given line reaches the specified width .
(optional, default = " " ) |
[in] | newline | : The input scalar of the same type and kind as the input str of arbitrary length type parameter, whose contents represent the dividing point between the wrapped lines in the input str .
For example, if the string contains the C-style newline character \(\ms{\n}\), then set \(\ms{newline = '\\n'}\).
Note that newline does not count toward the width or maxwidth of the wrapped line.
The presence of two adjacent newline subsrtings signals the beginning of a new paragraph.
In such a case, the indentation of the new paragraph will be computed afresh by couting the number of repetitions of the input indent pattern at the beginning of the new paragraph.
(optional, default = new_line(SKG_"a") where SKG refers to the kind of str .) |
[in] | linefeed | : The input scalar of the same type and kind as the input str of arbitrary length type parameter, whose contents represent the linefeed to use for dividing str into separate wrapped lines.
Each line that is to be wrapped will end with a linefeed .
Additionally, all instances of newline in the input str will be replaced with linefeed in the wrapped str .
Note that linefeed does not count toward the width or maxwidth of the wrapped line.
This argument is extremely useful for fast conversion of C-style newline instances to linefeed characters while wrapping str .
See below for example usages of this input argument.
(optional, default = the input argument newline ) |
[in] | width | : The input scalar of type integer of default kind IK, representing the length of each line in the output wrapping, with the possibility of word overflows.
If a word overflows the specified input line width, it is allowed to continue in the same line until the next break character is encountered.
Note that the lengths of the input arguments prefix , newline , and linefeed do not count toward the width or maxwidth of the wrapped line.
(optional, default = 132_IK ) |
[in] | maxwidth | : The input scalar of type integer of default kind IK, representing the maximum length of each wrapped line, beyond which no character can appear (except newline ).
If a word overflows the specified input maxwidth line limit, it will be broken and continued on the next line.
If break = "" , then the value of maxwidth becomes irrelevant because line wrapping will strictly occur at width .
(optional, default = huge(0_IK) / 2_IK ) |
- Returns
strWrapped
: The output allocatable
scalar character
of the same kind as the input str
, containing the input str
whose contents is divided into lines of maximum width maxwidth
by inserting the input newline
at the appropriate positions specified via break
.
Possible calling interfaces ⛓
character(:,kind(str)), allocatable :: strWrapped
strWrapped
= getStrWrapped(str, prefix
= prefix, indent
= indent, break
= break, newline
= newline, linefeed
= linefeed, width
= width, maxwidth
= maxwidth)
Generate and return the input scalar string wrapped within the specified width while optionally prefi...
This module contains classes and procedures for various string manipulations and inquiries.
- Warning
- The condition
len(indent) < width <= maxwidth
must hold at all times.
Otherwise, appropriate adjustments to the values of width
and maxWidth
will be made to infinite loops.
- See also
- display_type
getCentered
setCentered
isFailedGetShellWidth
isFailedGetShellHeight
Example usage ⛓
12 character(:, SK),
allocatable :: str, strWrapped
14 type(display_type) :: disp
18 str
= "A man asked my friend Jaime Cohen: `What is the human being's funniest characteristic?` &
19‘
&Cohen said: Our contradictoriness. We are in such a hurry to grow up, and then we long for our lost childhood. &
20 &We make ourselves ill earning money, and then spend all our money on getting well again. &
21 &We think so much about the future that we neglect the present, and thus experience neither the present nor the future. &
22 &We live as if we were never going to die, and die as if we had never lived."
25 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
26 call disp%show(
"! Wrap string within the default 132 characters.")
27 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
31 call disp%show( str, deliml
= SK_
"""" )
32 call disp%show(
"strWrapped = getStrWrapped(str)")
35 call disp%show( strWrapped, deliml
= SK_
"""" )
39 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
40 call disp%show(
"! Wrap string with an indent.")
41 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
45 call disp%show(
' '//str, deliml
= SK_
"""" )
46 call disp%show(
"strWrapped = getStrWrapped(' '//str)")
49 call disp%show(
LF//strWrapped, deliml
= SK_
"""" )
53 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
54 call disp%show(
"! Wrap string with an indent and a prefix for a given non-default width.")
55 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
59 call disp%show(
'....'//str, deliml
= SK_
"""" )
60 call disp%show(
"strWrapped = getStrWrapped('....'//str, prefix = SK_'NOTE - ', indent = SK_'.', width = 100_IK)")
61 strWrapped
= getStrWrapped(
'....'//str, prefix
= SK_
'NOTE - ', indent
= SK_
'.', width
= 100_IK)
63 call disp%show(
LF//strWrapped, deliml
= SK_
"""" )
67 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
68 call disp%show(
"! Wrap string with an indent and a prefix for a given non-default width at any character by specifying an empty `break` to ensure a precise width.")
69 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
73 call disp%show(
' '//str, deliml
= SK_
"""" )
74 call disp%show(
"strWrapped = getStrWrapped(' '//str, prefix = SK_'NOTE - ', indent = SK_' ', break = SK_'', width = 50_IK)")
75 strWrapped
= getStrWrapped(
' '//str, prefix
= SK_
'NOTE - ', indent
= SK_
' ', break
= SK_
'', width
= 50_IK)
77 call disp%show(
LF//strWrapped, deliml
= SK_
"""" )
81 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
82 call disp%show(
"! Setting a too small `width` (e.g., less than the initial indentation)can lead to awkward wrapping.")
83 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
87 call disp%show(
'....'//str, deliml
= SK_
"""" )
88 call disp%show(
"strWrapped = getStrWrapped(repeat(SK_'.', 8)//str, prefix = SK_'NOTE - ', indent = repeat(SK_'.', 4), width = 5_IK)")
89 strWrapped
= getStrWrapped(
repeat(SK_
'.',
8)
//str, prefix
= SK_
'NOTE - ', indent
= repeat(SK_
'.',
4), width
= 5_IK)
91 call disp%show(
LF//strWrapped, deliml
= SK_
"""" )
95 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
96 call disp%show(
"! Wrap string while converting all C-style newline instances (\n) to actual linefeed characters.")
97 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
100 str
= "The ParaMonte serial/parallel library contains Monte Carlo routines for \n\n&
101 & + optimization, \n&
104 &of mathematical objective functions of arbitrary-dimensions in particular, the posterior distributions of Bayesian models in \n\n&
105 & + data science, \n&
106 &+ Machine Learning, and \n&
107 &+ scientific inference, \n\n&
108 &with the design goal of unifying the \n&
109 & + automation (of Monte Carlo simulations), \n&
110 & + user-friendliness (of the library), \n&
111 & + accessibility (from multiple programming environments), \n&
112 & + high-performance (at runtime), and \n&
113 & + scalability (across many parallel processors)."
116 call disp%show( str, deliml
= SK_
"""" )
117 call disp%show(
"strWrapped = getStrWrapped(str, prefix = SK_'NOTE - ', newline = '\n', linefeed = LF, width = 80_IK)")
118 strWrapped
= getStrWrapped(str, prefix
= SK_
'NOTE - ', newline
= '\n', linefeed
= LF, width
= 80_IK)
120 call disp%show(
LF//strWrapped, deliml
= SK_
"""" )
122 call disp%show(
"strWrapped = getReplaced(getStrWrapped(str, prefix = SK_'NOTE - ', newline = '\n', width = 80_IK), pattern = '\n', replacement = LF) ! The slower verbose way of achieving the above.")
123 strWrapped
= getReplaced(
getStrWrapped(str, prefix
= SK_
'NOTE - ', newline
= '\n', width
= 80_IK), pattern
= '\n', replacement
= LF)
125 call disp%show(
LF//strWrapped, deliml
= SK_
"""" )
131 call disp%show( str, deliml
= SK_
"""" )
132 call disp%show(
"strWrapped = getStrWrapped(str, prefix = SK_'NOTE - ', newline = '\n', linefeed = LF, width = 80_IK)")
133 strWrapped
= getStrWrapped(str, prefix
= SK_
'NOTE - ', newline
= '\n', linefeed
= LF, width
= 80_IK)
135 call disp%show(
LF//strWrapped, deliml
= SK_
"""" )
137 call disp%show(
"strWrapped = getReplaced(getStrWrapped(str, prefix = SK_'NOTE - ', newline = '\n', width = 80_IK), pattern = '\n', replacement = LF) ! The slower verbose way of achieving the above.")
138 strWrapped
= getReplaced(
getStrWrapped(str, prefix
= SK_
'NOTE - ', newline
= '\n', width
= 80_IK), pattern
= '\n', replacement
= LF)
140 call disp%show(
LF//strWrapped, deliml
= SK_
"""" )
Generate and return an arrayNew of the same type and kind as the input array, in which the requested ...
This is a generic method of the derived type display_type with pass attribute.
This is a generic method of the derived type display_type with pass attribute.
This module contains procedures and generic interfaces for replacing patterns within arrays of variou...
This module contains classes and procedures for input/output (IO) or generic display operations on st...
type(display_type) disp
This is a scalar module variable an object of type display_type for general display.
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
This module contains the uncommon and hardly representable ASCII characters as well as procedures for...
character(1, SK), parameter LF
The scalar character constant of default kind SK of length 1 representing ASCII character: Line Feed ...
Generate and return an object of type display_type.
Example Unix compile command via Intel ifort
compiler ⛓
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
Example Windows Batch compile command via Intel ifort
compiler ⛓
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
Example Unix / MinGW compile command via GNU gfortran
compiler ⛓
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
Example output ⛓
7‘
"A man asked my friend Jaime Cohen: `What is the human being's funniest characteristic?` Cohen said: Our contradictoriness. We are in such a hurry to grow up, and then we long for our lost childhood. We make ourselves ill earning money, and then spend all our money on getting well again. We think so much about the future that we neglect the present, and thus experience neither the present nor the future. We live as if we were never going to die, and die as if we had never lived."
10‘
"A man asked my friend Jaime Cohen: `What is the human being's funniest characteristic?` Cohen said: Our contradictoriness. We are
11in such a hurry to grow up, and then we long for our lost childhood. We make ourselves ill earning money, and then spend all our money
12on getting well again. We think so much about the future that we neglect the present, and thus experience neither the present nor the
13future. We live as if we were never going to die, and die as if we had never lived."
21‘
" A man asked my friend Jaime Cohen: `What is the human being's funniest characteristic?` Cohen said: Our contradictoriness. We are in such a hurry to grow up, and then we long for our lost childhood. We make ourselves ill earning money, and then spend all our money on getting well again. We think so much about the future that we neglect the present, and thus experience neither the present nor the future. We live as if we were never going to die, and die as if we had never lived."
25‘
A man asked my friend Jaime Cohen: `What is the human being's funniest characteristic?` Cohen said: Our contradictoriness. We
26 are in such a hurry to grow up, and then we long for our lost childhood. We make ourselves ill earning money, and then spend all
27 our money on getting well again. We think so much about the future that we neglect the present, and thus experience neither the
28 present nor the future. We live as if we were never going to die, and die as if we had never lived."
36‘
"....A man asked my friend Jaime Cohen: `What is the human being's funniest characteristic?` Cohen said: Our contradictoriness. We are in such a hurry to grow up, and then we long for our lost childhood. We make ourselves ill earning money, and then spend all our money on getting well again. We think so much about the future that we neglect the present, and thus experience neither the present nor the future. We live as if we were never going to die, and die as if we had never lived."
37strWrapped
= getStrWrapped(
'....'//str, prefix
= SK_
'NOTE - ', indent
= SK_
'.', width
= 100_IK)
40NOTE - ....A man asked my friend Jaime Cohen: `What is the human being's funniest characteristic?` Cohen said:
41‘
NOTE - ....Our contradictoriness. We are in such a hurry to grow up, and then we long for our lost childhood.
42NOTE - ....We make ourselves ill earning money, and then spend all our money on getting well again. We think
43NOTE - ....so much about the future that we neglect the present, and thus experience neither the present nor
44NOTE - ....the future. We live as if we were never going to die, and die as if we had never lived."
52‘
" A man asked my friend Jaime Cohen: `What is the human being's funniest characteristic?` Cohen said: Our contradictoriness. We are in such a hurry to grow up, and then we long for our lost childhood. We make ourselves ill earning money, and then spend all our money on getting well again. We think so much about the future that we neglect the present, and thus experience neither the present nor the future. We live as if we were never going to die, and die as if we had never lived."
53strWrapped
= getStrWrapped(
' '//str, prefix
= SK_
'NOTE - ', indent
= SK_
' ', break
= SK_
'', width
= 50_IK)
56NOTE - A man asked my friend Jaime Cohen: `What is th
57NOTE - e human being's funniest characteristic?` Cohe
58‘
NOTE - n said: Our contradictoriness. We are in su
59NOTE - ch a hurry to grow up, and then we long for ou
60NOTE - r lost childhood. We make ourselves ill earnin
61NOTE - g money, and then spend all our money on getti
62NOTE - ng well again. We think so much about the futu
63NOTE - re that we neglect the present, and thus exper
64NOTE - ience neither the present nor the future. We l
65NOTE - ive as if we were never going to die, and die
66NOTE - as if we had never lived."
74‘
"....A man asked my friend Jaime Cohen: `What is the human being's funniest characteristic?` Cohen said: Our contradictoriness. We are in such a hurry to grow up, and then we long for our lost childhood. We make ourselves ill earning money, and then spend all our money on getting well again. We think so much about the future that we neglect the present, and thus experience neither the present nor the future. We live as if we were never going to die, and die as if we had never lived."
75strWrapped
= getStrWrapped(
repeat(SK_
'.',
8)
//str, prefix
= SK_
'NOTE - ', indent
= repeat(SK_
'.',
4), width
= 5_IK)
571"The ParaMonte serial/parallel library contains Monte Carlo routines for \n\n + optimization, \n+ sampling, and \n+ integration \n\nof mathematical objective functions of arbitrary-dimensions in particular, the posterior distributions of Bayesian models in \n\n + data science, \n+ Machine Learning, and \n+ scientific inference, \n\nwith the design goal of unifying the \n + automation (of Monte Carlo simulations), \n + user-friendliness (of the library), \n + accessibility (from multiple programming environments), \n + high-performance (at runtime), and \n + scalability (across many parallel processors)."
572strWrapped
= getStrWrapped(str, prefix
= SK_
'NOTE - ', newline
= '\n', linefeed
= LF, width
= 80_IK)
575NOTE - The ParaMonte serial/parallel library contains Monte Carlo routines for
577NOTE - + optimization,
578NOTE - + sampling, and
581NOTE - of mathematical objective functions of arbitrary-dimensions in particular, the posterior
582NOTE - distributions of Bayesian models in
584NOTE - + data science,
585NOTE - + Machine Learning, and
586NOTE - + scientific inference,
588NOTE - with the design goal of unifying the
589NOTE - + automation (of Monte Carlo simulations),
590NOTE - + user-friendliness (of the library),
591NOTE - + accessibility (from multiple programming environments),
592NOTE - + high-performance (at runtime), and
593NOTE - + scalability (across many parallel processors)."
595strWrapped
= getReplaced(
getStrWrapped(str, prefix
= SK_
'NOTE - ', newline
= '\n', width
= 80_IK), pattern
= '\n', replacement
= LF)
598NOTE - The ParaMonte serial/parallel library contains Monte Carlo routines for
600NOTE - + optimization,
601NOTE - + sampling, and
604NOTE - of mathematical objective functions of arbitrary-dimensions in particular, the posterior
605NOTE - distributions of Bayesian models in
607NOTE - + data science,
608NOTE - + Machine Learning, and
609NOTE - + scientific inference,
611NOTE - with the design goal of unifying the
612NOTE - + automation (of Monte Carlo simulations),
613NOTE - + user-friendliness (of the library),
614NOTE - + accessibility (from multiple programming environments),
615NOTE - + high-performance (at runtime), and
616NOTE - + scalability (across many parallel processors)."
619" The ParaMonte serial/parallel library contains Monte Carlo routines for \n\n + optimization, \n+ sampling, and \n+ integration \n\nof mathematical objective functions of arbitrary-dimensions in particular, the posterior distributions of Bayesian models in \n\n + data science, \n+ Machine Learning, and \n+ scientific inference, \n\nwith the design goal of unifying the \n + automation (of Monte Carlo simulations), \n + user-friendliness (of the library), \n + accessibility (from multiple programming environments), \n + high-performance (at runtime), and \n + scalability (across many parallel processors)."
620strWrapped
= getStrWrapped(str, prefix
= SK_
'NOTE - ', newline
= '\n', linefeed
= LF, width
= 80_IK)
623NOTE - The ParaMonte serial/parallel library contains Monte Carlo routines for
625NOTE - + optimization,
626NOTE - + sampling, and
629NOTE - of mathematical objective functions of arbitrary-dimensions in particular, the posterior
630NOTE - distributions of Bayesian models in
632NOTE - + data science,
633NOTE - + Machine Learning, and
634NOTE - + scientific inference,
636NOTE - with the design goal of unifying the
637NOTE - + automation (of Monte Carlo simulations),
638NOTE - + user-friendliness (of the library),
639NOTE - + accessibility (from multiple programming environments),
640NOTE - + high-performance (at runtime), and
641NOTE - + scalability (across many parallel processors)."
643strWrapped
= getReplaced(
getStrWrapped(str, prefix
= SK_
'NOTE - ', newline
= '\n', width
= 80_IK), pattern
= '\n', replacement
= LF)
646NOTE - The ParaMonte serial/parallel library contains Monte Carlo routines for
648NOTE - + optimization,
649NOTE - + sampling, and
652NOTE - of mathematical objective functions of arbitrary-dimensions in particular, the posterior
653NOTE - distributions of Bayesian models in
655NOTE - + data science,
656NOTE - + Machine Learning, and
657NOTE - + scientific inference,
659NOTE - with the design goal of unifying the
660NOTE - + automation (of Monte Carlo simulations),
661NOTE - + user-friendliness (of the library),
662NOTE - + accessibility (from multiple programming environments),
663NOTE - + high-performance (at runtime), and
664NOTE - + scalability (across many parallel processors)."
- Test:
- test_pm_str
Final Remarks ⛓
If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.
-
If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
-
If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.
This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.
- Copyright
- Computational Data Science Lab
- Author:
- Amir Shahmoradi, Tuesday March 7, 2017, 3:50 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas Austin
Definition at line 1617 of file pm_str.F90.