ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_io::wrap Interface Reference

This is a generic method of the derived type display_type with pass attribute.
More...

Detailed Description

This is a generic method of the derived type display_type with pass attribute.

This method centers the input string str within the specified width, adds margins as requested on all sides, then displays the wrapped string in the output display.

Note
  • The procedures under this generic interface do not perform any wrapping of the input text as done by getStrWrapped.
  • If needed, the input str must have been already wrapped via getStrWrapped with the desired input newline.
  • The procedures under this generic interface only sandwich each line in the input str with the specified left and right margins, then display the wrapped (sandwiched) lines.
  • If there are multiple lines within the input str, then they are all lef/right sandwiched in order, then displayed.
  • If top/bottom margins are also specified, then the entire set of output sandwiched lines is also collectively sandwiched by a top and a bottom margin, then displayed.
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 (sandwiched).
[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 separate lines in the input str to display.
For example, if the string contains the C-style newline character \(\ms{\n}\), then set \(\ms{newline = '\\n'}\).
See also the documentation of the corresponding argument in getStrWrapped.
(optional, default = new_line(SKG_"a") where SKG refers to the kind of str.)
[in]fill: The input scalar of the same type and kind as the input str of length type parameter len = 1, containing the value to fill the area between the content and left/right wrap margin in each line of the wrapped str in the display.
For more details, see the documentation of the corresponding argument in getCentered.
(optional, default = SKG_" " where SKG is the kind of the input argument str.)
[in]lwfill: The input scalar of the same type and kind as the input str of length type parameter len = 1, containing the value to fill the left wrap margin of each line of the wrapped str in the display.
For more details, see the documentation of the corresponding argument in getCentered.
(optional, default = MFILL)
[in]rwfill: The input scalar of the same type and kind as the input str of length type parameter len = 1, containing the value to fill the right wrap margin of each line of the wrapped str in the display.
For more details, see the documentation of the corresponding argument in getCentered.
(optional, default = MFILL)
[in]twfill: The input scalar of the same type and kind as the input str of length type parameter len = 1, containing the value to fill the top wrap margin of each line of the wrapped str in the display.
(optional, default = MFILL)
[in]bwfill: The input scalar of the same type and kind as the input str of length type parameter len = 1, containing the value to fill the bottom wrap margin of each line of the wrapped str in the display.
(optional, default = MFILL)
[in]lwsize: The input scalar integer of default kind IK representing the size of the left wrap margin of each line of the wrapped str in the display.
For more details, see the documentation of the corresponding argument in getCentered.
(optional, default = 2_IK)
[in]rwsize: The input scalar integer of default kind IK representing the size of the right wrap margin of each line of the wrapped str in the display.
For more details, see the documentation of the corresponding argument in getCentered.
(optional, default = 2_IK)
[in]twsize: The input scalar integer of default kind IK representing the size of the top wrap margin of the wrapped str in the display.
For more details, see the documentation of the corresponding argument in getCentered.
(optional, default = 1_IK)
[in]bwsize: The input scalar integer of default kind IK representing the size of the bottom wrap margin of the wrapped str in the display.
For more details, see the documentation of the corresponding argument in getCentered.
(optional, default = 1_IK)
[in]tmsize: The input scalar integer of default kind IK representing the number of rows in the display to skip before printing the wrapped str (i.e., the top margin).
(optional, default = 0_IK)
[in]bmsize: The input scalar integer of default kind IK representing the number of rows in the display to skip after printing the wrapped str (i.e., the bottom margin).
(optional, default = 0_IK)
[in]width: The input scalar integer of default kind IK representing the width of the contents of each wrapped line in the display.
For more details, see the documentation of the size argument in getCentered.
Note that the full length of each line in the display is lwsize + width + rwsize. (optional, default = 92_IK)
[in]unit: The input scalar integer of default kind IK representing the output display unit number where the text will be displayed.
(optional, default = display_type%unit)
[in]sticky: The input scalar of type logical of default kind LK.
If .true., the global properties (i.e., components) of the object of class wrap_type can be overridden by the corresponding arguments passed to the dynamic methods of the object.
(optional, default = .false.)


Possible calling interfaces

use pm_io, only: wrap_type
type(wrap_type) :: text
call text%wrap( str &
, newline = newline &
, fill = fill &
, lwfill = lwfill &
, rwfill = rwfill &
, twfill = twfill &
, bwfill = bwfill &
, lwsize = lwsize &
, rwsize = rwsize &
, twsize = twsize &
, bwsize = bwsize &
, tmsize = tmsize &
, bmsize = bmsize &
, width = width &
, unit = unit &
, sticky = sticky &
)
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:9995
This module contains classes and procedures for input/output (IO) or generic display operations on st...
Definition: pm_io.F90:252
This is the derived type for constructing objects that contain the specifications of the getStrWrappe...
Definition: pm_io.F90:9830
Warning
Any overflow of a single line beyond the specified width is handled as described in the documentation of getCentered.
Remarks
The procedures under discussion are impure.
Note
If the input str is potentially longer than the specified width, the overflowed text will be symmetrically cut from the beginning and the end of the str.
To avoid this behavior, str can be first passed to getStrWrapped to wrap it within the desired width.
The resulting strWrapped can be then passed to the procedures under this generic interface for centering, sandwiching, and displaying it.
This avoids the usually unnecessary and costly wrapping of the input str with the specified width by this interface.
See also
getStrWrapped
getCentered
setCentered
isFailedGetShellWidth
isFailedGetShellHeight


Example usage

1program example
2
3 use pm_io, only: wrap_type
4 use pm_io, only: display_type
5 use pm_str, only: getStrWrapped
6 use pm_kind, only: SK, IK, LK ! All kinds are supported.
7 use pm_strASCII, only: LF ! linefeed for better display.
8
9 implicit none
10
11 character(:, SK), allocatable :: str, strWrapped
12
13 type(wrap_type) :: text
14 type(display_type) :: disp
15
16 disp = display_type(file = SK_"main.out.F90")
17
18 str = "ParaMonte is a serial/parallel library of Monte Carlo routines for sampling mathematical density functions of arbitrary-dimensions&
19 &and Machine Learning algorithms for scientific inference, with the design goal of unifying **automation** (of simulations and tasks),&
20 &**user-friendliness** (of algorithms), **accessibility** (from any platform or programming environment),&
21 &**high-performance** (at runtime), and **scalability** (across many parallel processors)."
22
23 call disp%skip()
24 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
25 call disp%show("! Decorate the input 72-characters wrapped string within the default width.")
26 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
27 call disp%skip()
28
29 call disp%skip()
30 call disp%show("str")
31 call disp%show( str, deliml = SK_"""" )
32 call disp%show("strWrapped = getStrWrapped(str, width = 72_IK)")
33 strWrapped = getStrWrapped(str, width = 72_IK)
34 call disp%show("text = wrap_type(unit = disp%unit)")
35 text = wrap_type(unit = disp%unit)
36 call disp%show("call text%wrap(LF//strWrapped//LF)")
37 call text%wrap(LF//strWrapped//LF)
38 call disp%skip()
39
40 call disp%skip()
41 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
42 call disp%show("! Decorate the input 72-characters wrapped string within non-default settings.")
43 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
44 call disp%skip()
45
46 call disp%skip()
47 call disp%show("strWrapped")
48 call disp%show( strWrapped, deliml = SK_"""" )
49 call disp%show("call text%wrap(LF//strWrapped//LF, lwsize = 1_IK, rwsize = 1_IK, lwfill = SK_'|', rwfill = SK_'|', twfill = SK_'-', bwfill = SK_'-')")
50 call text%wrap(LF//strWrapped//LF, lwsize = 1_IK, rwsize = 1_IK, lwfill = SK_'|', rwfill = SK_'|', twfill = SK_'-', bwfill = SK_'-')
51 call disp%skip()
52
53 call disp%skip()
54 call disp%show("strWrapped")
55 call disp%show( strWrapped, deliml = SK_"""" )
56 call disp%show("call text%wrap(LF//strWrapped//LF, fill = SK_'.', lwsize = 1_IK, rwsize = 1_IK, lwfill = SK_'\', rwfill = SK_'/', twfill = SK_'_', bwfill = SK_'_')")
57 call text%wrap(LF//strWrapped//LF, fill = SK_'.', lwsize = 1_IK, rwsize = 1_IK, lwfill = SK_'\', rwfill = SK_'/', twfill = SK_'_', bwfill = SK_'_')
58 call disp%skip()
59
60 str = "ParaMonte is a serial/parallel library of Monte Carlo routines for sampling mathematical objective functions of arbitrary-dimensions."
61
62 call disp%skip()
63 call disp%show("str")
64 call disp%show( str, deliml = SK_"""" )
65 call disp%show("strWrapped = getStrWrapped(str)")
66 strWrapped = getStrWrapped(str)
67 call disp%show("call text%wrap(strWrapped) ! Note what happens to the overflowed text.")
68 call text%wrap(strWrapped) ! Note what happens to the overflowed text.
69 call disp%skip()
70
71 call disp%skip()
72 call disp%show("call text%wrap(LF//SK_'ParaMonte'//LF)")
73 call text%wrap(LF//SK_'ParaMonte'//LF)
74 call disp%skip()
75
76 call disp%skip()
77 call disp%show("call text%wrap('ParaMonte')")
78 call text%wrap('ParaMonte')
79 call disp%skip()
80
81 call disp%skip()
82 call disp%show("call text%wrap(LF//SK_''//LF)")
83 call text%wrap(LF//SK_''//LF)
84 call disp%skip()
85
86 call disp%skip()
87 call disp%show("call text%wrap(LF)")
88 call text%wrap(LF)
89 call disp%skip()
90
91 call disp%skip()
92 call disp%show("call text%wrap('')")
93 call text%wrap('')
94 call disp%skip()
95
96end program example
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11726
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11508
Generate and return the input scalar string wrapped within the specified width while optionally prefi...
Definition: pm_str.F90:1617
type(display_type) disp
This is a scalar module variable an object of type display_type for general display.
Definition: pm_io.F90:11393
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
Definition: pm_kind.F90:268
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
Definition: pm_kind.F90:541
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
Definition: pm_kind.F90:540
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
Definition: pm_kind.F90:539
This module contains the uncommon and hardly representable ASCII characters as well as procedures for...
Definition: pm_strASCII.F90:61
character(1, SK), parameter LF
The scalar character constant of default kind SK of length 1 representing ASCII character: Line Feed ...
Definition: pm_strASCII.F90:93
This module contains classes and procedures for various string manipulations and inquiries.
Definition: pm_str.F90:49
Generate and return an object of type display_type.
Definition: pm_io.F90:10282

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1
2!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
3! Decorate the input 72-characters wrapped string within the default width.
4!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6
7str
8"ParaMonte is a serial/parallel library of Monte Carlo routines for sampling mathematical density functions of arbitrary-dimensionsand Machine Learning algorithms for scientific inference, with the design goal of unifying **automation** (of simulations and tasks),**user-friendliness** (of algorithms), **accessibility** (from any platform or programming environment),**high-performance** (at runtime), and **scalability** (across many parallel processors)."
9strWrapped = getStrWrapped(str, width = 72_IK)
10text = wrap_type(unit = disp%unit)
11call text%wrap(LF//strWrapped//LF)
12%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
13%% %%
14%% ParaMonte is a serial/parallel library of Monte Carlo routines for sampling %%
15%% mathematical density functions of arbitrary-dimensionsand Machine Learning %%
16%% algorithms for scientific inference, with the design goal of unifying **automation** %%
17%% (of simulations and tasks),**user-friendliness** (of algorithms), **accessibility** %%
18%% (from any platform or programming environment),**high-performance** (at %%
19%% runtime), and **scalability** (across many parallel processors). %%
20%% %%
21%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
22
23
24!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
25! Decorate the input 72-characters wrapped string within non-default settings.
26!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
27
28
29strWrapped
30"ParaMonte is a serial/parallel library of Monte Carlo routines for sampling
31mathematical density functions of arbitrary-dimensionsand Machine Learning
32algorithms for scientific inference, with the design goal of unifying **automation**
33(of simulations and tasks),**user-friendliness** (of algorithms), **accessibility**
34(from any platform or programming environment),**high-performance** (at
35runtime), and **scalability** (across many parallel processors)."
36call text%wrap(LF//strWrapped//LF, lwsize = 1_IK, rwsize = 1_IK, lwfill = SK_'|', rwfill = SK_'|', twfill = SK_'-', bwfill = SK_'-')
37--------------------------------------------------------------------------------------------------
38| |
39| ParaMonte is a serial/parallel library of Monte Carlo routines for sampling |
40| mathematical density functions of arbitrary-dimensionsand Machine Learning |
41| algorithms for scientific inference, with the design goal of unifying **automation** |
42| (of simulations and tasks),**user-friendliness** (of algorithms), **accessibility** |
43| (from any platform or programming environment),**high-performance** (at |
44| runtime), and **scalability** (across many parallel processors). |
45| |
46--------------------------------------------------------------------------------------------------
47
48
49strWrapped
50"ParaMonte is a serial/parallel library of Monte Carlo routines for sampling
51mathematical density functions of arbitrary-dimensionsand Machine Learning
52algorithms for scientific inference, with the design goal of unifying **automation**
53(of simulations and tasks),**user-friendliness** (of algorithms), **accessibility**
54(from any platform or programming environment),**high-performance** (at
55runtime), and **scalability** (across many parallel processors)."
56call text%wrap(LF//strWrapped//LF, fill = SK_'.', lwsize = 1_IK, rwsize = 1_IK, lwfill = SK_'\', rwfill = SK_'/', twfill = SK_'_', bwfill = SK_'_')
57__________________________________________________________________________________________________
58\................................................................................................/
59\..........ParaMonte is a serial/parallel library of Monte Carlo routines for sampling ........../
60\..........mathematical density functions of arbitrary-dimensionsand Machine Learning .........../
61\.....algorithms for scientific inference, with the design goal of unifying **automation** ....../
62\......(of simulations and tasks),**user-friendliness** (of algorithms), **accessibility** ....../
63\............(from any platform or programming environment),**high-performance** (at ............/
64\................runtime), and **scalability** (across many parallel processors)................./
65\................................................................................................/
66__________________________________________________________________________________________________
67
68
69str
70"ParaMonte is a serial/parallel library of Monte Carlo routines for sampling mathematical objective functions of arbitrary-dimensions."
71strWrapped = getStrWrapped(str)
72call text%wrap(strWrapped) ! Note what happens to the overflowed text.
73%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
74%%ial/parallel library of Monte Carlo routines for sampling mathematical objective functions of ar%%
75%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
76
77
78call text%wrap(LF//SK_'ParaMonte'//LF)
79%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
80%% %%
81%% ParaMonte %%
82%% %%
83%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
84
85
86call text%wrap('ParaMonte')
87%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
88%% ParaMonte %%
89%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
90
91
92call text%wrap(LF//SK_''//LF)
93%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
94%% %%
95%% %%
96%% %%
97%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
98
99
100call text%wrap(LF)
101%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
102%% %%
103%% %%
104%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
105
106
107call text%wrap('')
108%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
109%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
110
111
Test:
test_pm_io


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.

  1. 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.
  2. 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.

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 9995 of file pm_io.F90.


The documentation for this interface was generated from the following file: