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

Read a full record (line) of arbitrary length as a string from the current position of the record-oriented and formatted file connected to the specified input unit. More...

Detailed Description

Read a full record (line) of arbitrary length as a string from the current position of the record-oriented and formatted file connected to the specified input unit.

Parameters
[in]unit: The input scalar integer of default kind IK representing the unit of the connected file to read the line from.
[in,out]record: The input/output allocatable scalar character of default kind SK that will contain the requested line in full.
On input,
  1. The argument record can be preallocated to the best-guess length for the record to be read.
    Preallocating record to the right size can improve the runtime performance of the procedure.
    1. If the optional input argument lb is specified, any old contents record(1 : lb - 1) will remain intact on return.
  2. The argument record can be unallocated, in which case, it will be allocated to the proper size within the algorithm.
On output,
  1. The argument record will potentially be resized to contain the full record from the specified unit.
    1. If the optional output argument ub is specified, then the output record will not be trimmed (i.e., reallocated) on return.
      Specifying ub will eliminate one redundant final reallocation and copy action, which can in return boost the runtime performance.
[out]iostat: The output scalar integer of default kind IK.
  1. If present and no error occurs, it is set to 0 on output.
  2. If present and an end-of-file condition occurs, it is set to iostat_end from Fortran intrinsic module iso_fortran_env.
  3. If present and an error occurs (e.g., if the input argument values are wrong or inconsistent), it is set to a positive non-zero value.
  4. If missing and an error occurs (including end-of-file condition), then the program halts by calling error stop followed by the relevant error message.
(optional. It must be present if and only if the output argument iomsg is also present.)
[in,out]iomsg: The input/output scalar character of default kind SK containing the error message, if any error occurs.
A length type parameter value of LEN_IOMSG is generally sufficient for iomsg to contain the output error messages.
(optional. It must be present if and only if the output argument iostat is also present.)
[in]lb: The input scalar integer of default kind IK, containing the lower starting bound (index) of record from which the writing must begin.
Specifying this argument will keep the segment record(1 : lb - 1) intact on return.
This behavior is extremely useful for consecutively reading a series of lines of a unit.
(optional, default = 1)
[out]ub: The output scalar integer of default kind IK, containing the upper bound (index) of record up to which the record from the specified unit was read.
Specifying this argument leads to faster runtime performance since it prevents an extra final reallocation and string copy.
By definition, if the line in the file is empty, ub is set to 0, unless linefed is set to .true., in which case ub is the length of the linefeed character(s).
(optional. If missing, record will be reallocated to size record(1 : ub), otherwise, record will not be trimmed.)
[in]linefed: The input scalar logical of default kind LK.
If .true., then the output record will end with the new line character(s) as specified by the Fortran intrinsic new_line("a").
If .false., then the output record will not end with the new line character(s).
In either case, only one record will be read from the specified unit, but will or will not end with a new line character if linefed is .true. or .false. respectively.
This behavior is extremely useful for,
  1. reading records from a CSV file potentially containing new line characters in its fields.
  2. consecutively reading a series of lines of a unit into a single string by repeated calls to this generic interface while preserving the linefeed characters as record separators.
(optional, default = .false.)


Possible calling interfaces

use pm_io, only: setRecordFrom
call setRecordFrom(unit, record, lb = lb, ub = ub, linefed = linefed)
call setRecordFrom(unit, record, iostat, iomsg, lb = lb, ub = ub, linefed = linefed)
Read a full record (line) of arbitrary length as a string from the current position of the record-ori...
Definition: pm_io.F90:2114
This module contains classes and procedures for input/output (IO) or generic display operations on st...
Definition: pm_io.F90:252
Warning
The condition 0 < lb must hold for the corresponding input arguments. This condition is verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
Remarks
Allowing the input record to be preallocated with intent(inout) particularly helps with minimizing the number of reallocations of record within the procedure.
Specifying the output argument ub particularly helps with the runtime performance by removing a redundant reallocation and data copy.
Note
This generic interface can be readily used to read a certain number or all of the lines of a specified unit.
For example, the following code snippet,
use iso_fortran_env, only: iostat_end
open(unit, file, status = "old")
lb = 1_IK
do
call setRecordFrom(unit, record, iostat, iomsg, lb, ub, linefed = .true._LK)
if (iostat /= iostat_end) exit
if (iostat /= 0_IK) error stop trim(iomsg)
lb = ub + 1_IK
end do
close(unit)
record = record(1:ub)
will store the entire contents of the specified file associated with unit in the allocatable string record, within which lines are separated by new_line("a") instances.
Note that setting linefed = .true. leads to each record being padded with a newline character, including the last line in the file.
As such, reading an entire file using the above code snippet will lead to a final output that an has extra final newline character compared to the output of getContentsFrom() or setContentsFrom().
See also
getRecordFrom
getContentsFrom
setContentsFrom
isPreconnected
getFileUnit


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK
4 use pm_io, only: display_type
5 use pm_io, only: setRecordFrom
6
7 implicit none
8
9 integer(IK) :: ub
10 integer(IK) :: unit
11 integer(IK) :: iostat
12 character(255, SK) :: iomsg
13 character( :, SK), allocatable :: record
14
15 type(display_type) :: disp
16
17 disp = display_type(file = "main.out.F90")
18
19 call disp%skip
20 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
21 call disp%show("!Read the record-oriented (sequential access) `main.F90` file, line by line to the end.")
22 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
23 call disp%skip
24
25 iomsg = repeat(" ", len(iomsg))
26
27 open(newunit = unit, file = "main.F90", status = "old")
28
29 do
30 call disp%skip
31 call disp%show("call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)")
32 call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
33 if (iostat == 0_IK) then
34 call disp%show("record")
35 call disp%show( record , deliml = SK_"""" )
36 else
37 call disp%show("iostat")
38 call disp%show( iostat )
39 call disp%show("trim(adjustl(iomsg))")
40 call disp%show( trim(adjustl(iomsg)) , deliml = SK_"""" )
41 exit
42 end if
43 call disp%skip
44 end do
45
46 close(unit)
47
48 call disp%skip
49 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
50 call disp%show("!Read a record-oriented (sequential access) `main.F90` file, line by line to the end, faster without redundant reallocations and copies.")
51 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
52 call disp%skip
53
54 iomsg = repeat(" ", len(iomsg))
55
56 call disp%show('open(newunit = unit, file = "main.F90", status = "old")')
57 open(newunit = unit, file = "main.F90", status = "old")
58 do
59 call disp%skip
60 call disp%show("call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)")
61 call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
62 if (iostat == 0_IK) then
63 call disp%show("record(1:ub)")
64 call disp%show( record(1:ub) , deliml = SK_"""" )
65 else
66 call disp%show("iostat")
67 call disp%show( iostat )
68 call disp%show("trim(adjustl(iomsg))")
69 call disp%show( trim(adjustl(iomsg)) , deliml = SK_"""" )
70 exit
71 end if
72 call disp%skip
73 end do
74
75 close(unit)
76
77 call disp%skip
78 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
79 call disp%show("!Read a record-oriented (sequential access) `main.F90` file line by line to the end in a single string, each line separated by linefeed (see the end of contents for the source code).")
80 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
81 call disp%skip
82
83 block
84 use iso_fortran_env, only: iostat_end
86 integer(IK) :: lb, ub
87 open(newunit = unit, file = "main.F90", status = "old")
88 iomsg = repeat(" ", len(iomsg))
89 lb = 1_IK
90 do
91 call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, lb = lb, ub = ub, linefed = .true._LK)
92 if (iostat == 0_IK) then
93 lb = ub + 1_IK
94 cycle
95 elseif (iostat == iostat_end) then
96 exit
97 else
98 error stop SK_"Unknown IO error occurred: "//trim(iomsg)
99 end if
100 end do
101 call disp%show( getCentered(SK_" records ", size = 132_IK, fill = SK_"_") )
102 call disp%show( record(1:ub) )
103 call disp%show( getCentered(SK_" end of records ", size = 132_IK, fill = SK_"_") )
104 close(unit)
105 end block
106
107end program example
Generate and return a resized and centered copy of the input array within the newly allocated arrayCe...
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
This module contains procedures and generic interfaces for resizing an input array and centering the ...
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
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!Read the record-oriented (sequential access) `main.F90` file, line by line to the end.
4!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6
7call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
8record
9"program example"
10
11
12call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
13record
14""
15
16
17call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
18record
19" use pm_kind, only: SK, IK, LK"
20
21
22call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
23record
24" use pm_io, only: display_type"
25
26
27call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
28record
29" use pm_io, only: setRecordFrom"
30
31
32call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
33record
34""
35
36
37call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
38record
39" implicit none"
40
41
42call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
43record
44""
45
46
47call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
48record
49" integer(IK) :: ub"
50
51
52call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
53record
54" integer(IK) :: unit"
55
56
57call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
58record
59" integer(IK) :: iostat"
60
61
62call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
63record
64" character(255, SK) :: iomsg"
65
66
67call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
68record
69" character( :, SK), allocatable :: record"
70
71
72call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
73record
74""
75
76
77call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
78record
79" type(display_type) :: disp"
80
81
82call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
83record
84""
85
86
87call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
88record
89" disp = display_type(file = "main.out.F90")"
90
91
92call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
93record
94""
95
96
97call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
98record
99" call disp%skip"
100
101
102call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
103record
104" call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")"
105
106
107call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
108record
109" call disp%show("!Read the record-oriented (sequential access) `main.F90` file, line by line to the end.")"
110
111
112call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
113record
114" call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")"
115
116
117call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
118record
119" call disp%skip"
120
121
122call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
123record
124""
125
126
127call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
128record
129" iomsg = repeat(" ", len(iomsg))"
130
131
132call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
133record
134""
135
136
137call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
138record
139" open(newunit = unit, file = "main.F90", status = "old")"
140
141
142call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
143record
144""
145
146
147call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
148record
149" do"
150
151
152call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
153record
154" call disp%skip"
155
156
157call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
158record
159" call disp%show("call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)")"
160
161
162call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
163record
164" call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)"
165
166
167call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
168record
169" if (iostat == 0_IK) then"
170
171
172call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
173record
174" call disp%show("record")"
175
176
177call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
178record
179" call disp%show( record , deliml = SK_"""" )"
180
181
182call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
183record
184" else"
185
186
187call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
188record
189" call disp%show("iostat")"
190
191
192call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
193record
194" call disp%show( iostat )"
195
196
197call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
198record
199" call disp%show("trim(adjustl(iomsg))")"
200
201
202call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
203record
204" call disp%show( trim(adjustl(iomsg)) , deliml = SK_"""" )"
205
206
207call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
208record
209" exit"
210
211
212call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
213record
214" end if"
215
216
217call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
218record
219" call disp%skip"
220
221
222call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
223record
224" end do"
225
226
227call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
228record
229""
230
231
232call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
233record
234" close(unit)"
235
236
237call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
238record
239""
240
241
242call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
243record
244" call disp%skip"
245
246
247call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
248record
249" call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")"
250
251
252call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
253record
254" call disp%show("!Read a record-oriented (sequential access) `main.F90` file, line by line to the end, faster without redundant reallocations and copies.")"
255
256
257call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
258record
259" call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")"
260
261
262call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
263record
264" call disp%skip"
265
266
267call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
268record
269""
270
271
272call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
273record
274" iomsg = repeat(" ", len(iomsg))"
275
276
277call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
278record
279""
280
281
282call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
283record
284" call disp%show('open(newunit = unit, file = "main.F90", status = "old")')"
285
286
287call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
288record
289" open(newunit = unit, file = "main.F90", status = "old")"
290
291
292call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
293record
294" do"
295
296
297call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
298record
299" call disp%skip"
300
301
302call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
303record
304" call disp%show("call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)")"
305
306
307call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
308record
309" call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)"
310
311
312call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
313record
314" if (iostat == 0_IK) then"
315
316
317call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
318record
319" call disp%show("record(1:ub)")"
320
321
322call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
323record
324" call disp%show( record(1:ub) , deliml = SK_"""" )"
325
326
327call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
328record
329" else"
330
331
332call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
333record
334" call disp%show("iostat")"
335
336
337call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
338record
339" call disp%show( iostat )"
340
341
342call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
343record
344" call disp%show("trim(adjustl(iomsg))")"
345
346
347call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
348record
349" call disp%show( trim(adjustl(iomsg)) , deliml = SK_"""" )"
350
351
352call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
353record
354" exit"
355
356
357call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
358record
359" end if"
360
361
362call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
363record
364" call disp%skip"
365
366
367call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
368record
369" end do"
370
371
372call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
373record
374""
375
376
377call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
378record
379" close(unit)"
380
381
382call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
383record
384""
385
386
387call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
388record
389" call disp%skip"
390
391
392call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
393record
394" call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")"
395
396
397call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
398record
399" call disp%show("!Read a record-oriented (sequential access) `main.F90` file line by line to the end in a single string, each line separated by linefeed (see the end of contents for the source code).")"
400
401
402call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
403record
404" call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")"
405
406
407call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
408record
409" call disp%skip"
410
411
412call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
413record
414""
415
416
417call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
418record
419" block"
420
421
422call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
423record
424" use iso_fortran_env, only: iostat_end"
425
426
427call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
428record
429" use pm_arrayCenter, only: getCentered"
430
431
432call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
433record
434" integer(IK) :: lb, ub"
435
436
437call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
438record
439" open(newunit = unit, file = "main.F90", status = "old")"
440
441
442call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
443record
444" iomsg = repeat(" ", len(iomsg))"
445
446
447call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
448record
449" lb = 1_IK"
450
451
452call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
453record
454" do"
455
456
457call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
458record
459" call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, lb = lb, ub = ub, linefed = .true._LK)"
460
461
462call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
463record
464" if (iostat == 0_IK) then"
465
466
467call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
468record
469" lb = ub + 1_IK"
470
471
472call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
473record
474" cycle"
475
476
477call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
478record
479" elseif (iostat == iostat_end) then"
480
481
482call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
483record
484" exit"
485
486
487call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
488record
489" else"
490
491
492call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
493record
494" error stop SK_"Unknown IO error occurred: "//trim(iomsg)"
495
496
497call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
498record
499" end if"
500
501
502call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
503record
504" end do"
505
506
507call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
508record
509" call disp%show( getCentered(SK_" records ", size = 132_IK, fill = SK_"_") )"
510
511
512call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
513record
514" call disp%show( record(1:ub) )"
515
516
517call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
518record
519" call disp%show( getCentered(SK_" end of records ", size = 132_IK, fill = SK_"_") )"
520
521
522call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
523record
524" close(unit)"
525
526
527call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
528record
529" end block"
530
531
532call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
533record
534""
535
536
537call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
538record
539"end program example"
540
541
542call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
543iostat
544-1
545trim(adjustl(iomsg))
546"End of file"
547
548!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
549!Read a record-oriented (sequential access) `main.F90` file, line by line to the end, faster without redundant reallocations and copies.
550!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
551
552open(newunit = unit, file = "main.F90", status = "old")
553
554call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
555record(1:ub)
556"program example"
557
558
559call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
560record(1:ub)
561""
562
563
564call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
565record(1:ub)
566" use pm_kind, only: SK, IK, LK"
567
568
569call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
570record(1:ub)
571" use pm_io, only: display_type"
572
573
574call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
575record(1:ub)
576" use pm_io, only: setRecordFrom"
577
578
579call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
580record(1:ub)
581""
582
583
584call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
585record(1:ub)
586" implicit none"
587
588
589call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
590record(1:ub)
591""
592
593
594call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
595record(1:ub)
596" integer(IK) :: ub"
597
598
599call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
600record(1:ub)
601" integer(IK) :: unit"
602
603
604call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
605record(1:ub)
606" integer(IK) :: iostat"
607
608
609call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
610record(1:ub)
611" character(255, SK) :: iomsg"
612
613
614call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
615record(1:ub)
616" character( :, SK), allocatable :: record"
617
618
619call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
620record(1:ub)
621""
622
623
624call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
625record(1:ub)
626" type(display_type) :: disp"
627
628
629call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
630record(1:ub)
631""
632
633
634call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
635record(1:ub)
636" disp = display_type(file = "main.out.F90")"
637
638
639call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
640record(1:ub)
641""
642
643
644call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
645record(1:ub)
646" call disp%skip"
647
648
649call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
650record(1:ub)
651" call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")"
652
653
654call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
655record(1:ub)
656" call disp%show("!Read the record-oriented (sequential access) `main.F90` file, line by line to the end.")"
657
658
659call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
660record(1:ub)
661" call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")"
662
663
664call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
665record(1:ub)
666" call disp%skip"
667
668
669call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
670record(1:ub)
671""
672
673
674call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
675record(1:ub)
676" iomsg = repeat(" ", len(iomsg))"
677
678
679call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
680record(1:ub)
681""
682
683
684call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
685record(1:ub)
686" open(newunit = unit, file = "main.F90", status = "old")"
687
688
689call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
690record(1:ub)
691""
692
693
694call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
695record(1:ub)
696" do"
697
698
699call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
700record(1:ub)
701" call disp%skip"
702
703
704call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
705record(1:ub)
706" call disp%show("call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)")"
707
708
709call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
710record(1:ub)
711" call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)"
712
713
714call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
715record(1:ub)
716" if (iostat == 0_IK) then"
717
718
719call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
720record(1:ub)
721" call disp%show("record")"
722
723
724call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
725record(1:ub)
726" call disp%show( record , deliml = SK_"""" )"
727
728
729call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
730record(1:ub)
731" else"
732
733
734call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
735record(1:ub)
736" call disp%show("iostat")"
737
738
739call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
740record(1:ub)
741" call disp%show( iostat )"
742
743
744call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
745record(1:ub)
746" call disp%show("trim(adjustl(iomsg))")"
747
748
749call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
750record(1:ub)
751" call disp%show( trim(adjustl(iomsg)) , deliml = SK_"""" )"
752
753
754call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
755record(1:ub)
756" exit"
757
758
759call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
760record(1:ub)
761" end if"
762
763
764call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
765record(1:ub)
766" call disp%skip"
767
768
769call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
770record(1:ub)
771" end do"
772
773
774call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
775record(1:ub)
776""
777
778
779call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
780record(1:ub)
781" close(unit)"
782
783
784call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
785record(1:ub)
786""
787
788
789call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
790record(1:ub)
791" call disp%skip"
792
793
794call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
795record(1:ub)
796" call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")"
797
798
799call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
800record(1:ub)
801" call disp%show("!Read a record-oriented (sequential access) `main.F90` file, line by line to the end, faster without redundant reallocations and copies.")"
802
803
804call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
805record(1:ub)
806" call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")"
807
808
809call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
810record(1:ub)
811" call disp%skip"
812
813
814call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
815record(1:ub)
816""
817
818
819call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
820record(1:ub)
821" iomsg = repeat(" ", len(iomsg))"
822
823
824call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
825record(1:ub)
826""
827
828
829call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
830record(1:ub)
831" call disp%show('open(newunit = unit, file = "main.F90", status = "old")')"
832
833
834call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
835record(1:ub)
836" open(newunit = unit, file = "main.F90", status = "old")"
837
838
839call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
840record(1:ub)
841" do"
842
843
844call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
845record(1:ub)
846" call disp%skip"
847
848
849call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
850record(1:ub)
851" call disp%show("call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)")"
852
853
854call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
855record(1:ub)
856" call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)"
857
858
859call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
860record(1:ub)
861" if (iostat == 0_IK) then"
862
863
864call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
865record(1:ub)
866" call disp%show("record(1:ub)")"
867
868
869call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
870record(1:ub)
871" call disp%show( record(1:ub) , deliml = SK_"""" )"
872
873
874call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
875record(1:ub)
876" else"
877
878
879call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
880record(1:ub)
881" call disp%show("iostat")"
882
883
884call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
885record(1:ub)
886" call disp%show( iostat )"
887
888
889call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
890record(1:ub)
891" call disp%show("trim(adjustl(iomsg))")"
892
893
894call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
895record(1:ub)
896" call disp%show( trim(adjustl(iomsg)) , deliml = SK_"""" )"
897
898
899call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
900record(1:ub)
901" exit"
902
903
904call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
905record(1:ub)
906" end if"
907
908
909call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
910record(1:ub)
911" call disp%skip"
912
913
914call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
915record(1:ub)
916" end do"
917
918
919call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
920record(1:ub)
921""
922
923
924call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
925record(1:ub)
926" close(unit)"
927
928
929call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
930record(1:ub)
931""
932
933
934call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
935record(1:ub)
936" call disp%skip"
937
938
939call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
940record(1:ub)
941" call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")"
942
943
944call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
945record(1:ub)
946" call disp%show("!Read a record-oriented (sequential access) `main.F90` file line by line to the end in a single string, each line separated by linefeed (see the end of contents for the source code).")"
947
948
949call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
950record(1:ub)
951" call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")"
952
953
954call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
955record(1:ub)
956" call disp%skip"
957
958
959call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
960record(1:ub)
961""
962
963
964call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
965record(1:ub)
966" block"
967
968
969call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
970record(1:ub)
971" use iso_fortran_env, only: iostat_end"
972
973
974call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
975record(1:ub)
976" use pm_arrayCenter, only: getCentered"
977
978
979call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
980record(1:ub)
981" integer(IK) :: lb, ub"
982
983
984call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
985record(1:ub)
986" open(newunit = unit, file = "main.F90", status = "old")"
987
988
989call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
990record(1:ub)
991" iomsg = repeat(" ", len(iomsg))"
992
993
994call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
995record(1:ub)
996" lb = 1_IK"
997
998
999call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
1000record(1:ub)
1001" do"
1002
1003
1004call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
1005record(1:ub)
1006" call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, lb = lb, ub = ub, linefed = .true._LK)"
1007
1008
1009call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
1010record(1:ub)
1011" if (iostat == 0_IK) then"
1012
1013
1014call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
1015record(1:ub)
1016" lb = ub + 1_IK"
1017
1018
1019call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
1020record(1:ub)
1021" cycle"
1022
1023
1024call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
1025record(1:ub)
1026" elseif (iostat == iostat_end) then"
1027
1028
1029call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
1030record(1:ub)
1031" exit"
1032
1033
1034call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
1035record(1:ub)
1036" else"
1037
1038
1039call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
1040record(1:ub)
1041" error stop SK_"Unknown IO error occurred: "//trim(iomsg)"
1042
1043
1044call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
1045record(1:ub)
1046" end if"
1047
1048
1049call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
1050record(1:ub)
1051" end do"
1052
1053
1054call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
1055record(1:ub)
1056" call disp%show( getCentered(SK_" records ", size = 132_IK, fill = SK_"_") )"
1057
1058
1059call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
1060record(1:ub)
1061" call disp%show( record(1:ub) )"
1062
1063
1064call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
1065record(1:ub)
1066" call disp%show( getCentered(SK_" end of records ", size = 132_IK, fill = SK_"_") )"
1067
1068
1069call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
1070record(1:ub)
1071" close(unit)"
1072
1073
1074call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
1075record(1:ub)
1076" end block"
1077
1078
1079call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
1080record(1:ub)
1081""
1082
1083
1084call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
1085record(1:ub)
1086"end program example"
1087
1088
1089call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
1090iostat
1091-1
1092trim(adjustl(iomsg))
1093"End of file"
1094
1095!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1096!Read a record-oriented (sequential access) `main.F90` file line by line to the end in a single string, each line separated by linefeed (see the end of contents for the source code).
1097!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1098
1099_____________________________________________________________ records ______________________________________________________________
1100program example
1101
1102 use pm_kind, only: SK, IK, LK
1103 use pm_io, only: display_type
1104 use pm_io, only: setRecordFrom
1105
1106 implicit none
1107
1108 integer(IK) :: ub
1109 integer(IK) :: unit
1110 integer(IK) :: iostat
1111 character(255, SK) :: iomsg
1112 character( :, SK), allocatable :: record
1113
1114 type(display_type) :: disp
1115
1116 disp = display_type(file = "main.out.F90")
1117
1118 call disp%skip
1119 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
1120 call disp%show("!Read the record-oriented (sequential access) `main.F90` file, line by line to the end.")
1121 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
1122 call disp%skip
1123
1124 iomsg = repeat(" ", len(iomsg))
1125
1126 open(newunit = unit, file = "main.F90", status = "old")
1127
1128 do
1129 call disp%skip
1130 call disp%show("call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)")
1131 call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg)
1132 if (iostat == 0_IK) then
1133 call disp%show("record")
1134 call disp%show( record , deliml = SK_"""" )
1135 else
1136 call disp%show("iostat")
1137 call disp%show( iostat )
1138 call disp%show("trim(adjustl(iomsg))")
1139 call disp%show( trim(adjustl(iomsg)) , deliml = SK_"""" )
1140 exit
1141 end if
1142 call disp%skip
1143 end do
1144
1145 close(unit)
1146
1147 call disp%skip
1148 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
1149 call disp%show("!Read a record-oriented (sequential access) `main.F90` file, line by line to the end, faster without redundant reallocations and copies.")
1150 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
1151 call disp%skip
1152
1153 iomsg = repeat(" ", len(iomsg))
1154
1155 call disp%show('open(newunit = unit, file = "main.F90", status = "old")')
1156 open(newunit = unit, file = "main.F90", status = "old")
1157 do
1158 call disp%skip
1159 call disp%show("call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)")
1160 call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, ub = ub)
1161 if (iostat == 0_IK) then
1162 call disp%show("record(1:ub)")
1163 call disp%show( record(1:ub) , deliml = SK_"""" )
1164 else
1165 call disp%show("iostat")
1166 call disp%show( iostat )
1167 call disp%show("trim(adjustl(iomsg))")
1168 call disp%show( trim(adjustl(iomsg)) , deliml = SK_"""" )
1169 exit
1170 end if
1171 call disp%skip
1172 end do
1173
1174 close(unit)
1175
1176 call disp%skip
1177 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
1178 call disp%show("!Read a record-oriented (sequential access) `main.F90` file line by line to the end in a single string, each line separated by linefeed (see the end of contents for the source code).")
1179 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
1180 call disp%skip
1181
1182 block
1183 use iso_fortran_env, only: iostat_end
1184 use pm_arrayCenter, only: getCentered
1185 integer(IK) :: lb, ub
1186 open(newunit = unit, file = "main.F90", status = "old")
1187 iomsg = repeat(" ", len(iomsg))
1188 lb = 1_IK
1189 do
1190 call setRecordFrom(unit, record, iostat = iostat, iomsg = iomsg, lb = lb, ub = ub, linefed = .true._LK)
1191 if (iostat == 0_IK) then
1192 lb = ub + 1_IK
1193 cycle
1194 elseif (iostat == iostat_end) then
1195 exit
1196 else
1197 error stop SK_"Unknown IO error occurred: "//trim(iomsg)
1198 end if
1199 end do
1200 call disp%show( getCentered(SK_" records ", size = 132_IK, fill = SK_"_") )
1201 call disp%show( record(1:ub) )
1202 call disp%show( getCentered(SK_" end of records ", size = 132_IK, fill = SK_"_") )
1203 close(unit)
1204 end block
1205
1206end program example
1207
1208__________________________________________________________ end of records __________________________________________________________
1209
program main
This is main entry to the tests of the ParaMonte kernel library.
Definition: main.F90:27
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 at Austin

Definition at line 2114 of file pm_io.F90.


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