source: trunk/LMDZ.PLUTO/libf/muphypluto/swift_asciiread.F90 @ 3590

Last change on this file since 3590 was 3560, checked in by debatzbr, 5 weeks ago

Addition of the microphysics model in moments.

File size: 14.2 KB
RevLine 
[3560]1! Copyright (c) Jeremie Burgalat (2013-2022)
2! Contributor: Jeremie Burgalat (jeremie.burgalat@univ-reims.fr).
3!
4! This file is part of SWIFT
5!
6! Permission is hereby granted, free of charge, to any person obtaining a copy of
7! this software and associated documentation files (the "Software"), to deal in
8! the Software without restriction, including without limitation the rights to
9! use, copy, modify, merge, publish, distribute, sublicense, and/or sell copies of
10! the Software, and to permit persons to whom the Software is furnished to do so,
11! subject to the following conditions:
12!
13! The above copyright notice and this permission notice shall be included in all
14! copies or substantial portions of the Software.
15!
16! THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
17! IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY, FITNESS
18! FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS OR
19! COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER
20! IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN
21! CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.
22
23!! File:    swift_asciiread.F90
24!! Summary: ASCII data file reader source file
25!! Author:  J. Burgalat
26!! Date:    2013-2022
27
28MODULE SWIFT_ASCIIREAD
29    !! ASCII data file reader module
30    !!
31    !! This module provides a single generic method that can be used to read 2D/3D
32    !! data array from ASCII file.
33    !!
34    !! ```
35    !! FUNCTION read_data(path,data) RESULT(err)
36    !! ```
37    !!
38    !! Where:
39    !!
40    !! - __path__ is a string with the path of the data file.
41    !! - __data__ is an output __allocatable__ 2D/3D array of real(kind=8) values.
42    !!
43    !! ## Expected Format of the data file
44    !!
45    !! The input file:
46    !!   - must use blank space(s) as value delimiter.
47    !!   - must have a regular number of columns, that is each data line must
48    !!     have the same number of columns.
49    !!   - can contains any number of empty lines and/or comment line (i.e. line
50    !!     where first non-blank character is "#"). All other lines are assumed
51    !!     to be data.
52    !! Moreover, in the case of 3D data, it must use a SINGLE empty line for "depth" block separator.
53    !!
54    !! Error occured when:
55    !!
56    !! - path does not refer to a existing file (1)
57    !! - No free logical unit available (1)
58    !! - the file does not have regular data-columns number (5)
59    !! - at least a value cannot be cast in double precision (5)
60    !!
61    !! On error, __data__ array is not allocated.
62    !!
63    !! ## 3D data structure
64    !!
65    !! In the case of 3D data the method assumes that the input files consists in _D_ blocks
66    !! of _R_ lines with _C_ columns. Each block must be separated by a single empty line
67    !! and each columns must separated by one or more blank spaces (no tabulation ALLOWED).
68    !!
69    !! On success, the shape of the 3D output array will be _data(R,C,D)_.
70    USE, INTRINSIC :: ISO_FORTRAN_ENV, ONLY: IOSTAT_END, IOSTAT_EOR
71    !USE STRING_OP, ONLY: tokenize,from_string, st_slen
72    USE SWIFT_STRING_OP
73    USE SWIFT_ERRORS
74    IMPLICIT NONE
75 
76    PRIVATE
77    PUBLIC :: noerror,error, error_to_string,aborting
78    PUBLIC :: readline,read_data, OPERATOR(/=), OPERATOR(==)
79 
80    !! Global interface to reading methods
81    INTERFACE read_data
82      MODULE PROCEDURE read_data_2d, read_data_3d
83    END INTERFACE
84 
85    CONTAINS
86 
87 
88    FUNCTION read_data_3d(path,data3d,delimiter) RESULT(err)
89      !! Read an input ASCII data file and stores its content in a 3D array
90      !!
91      !! The function reads an ASCII file and saves its values in a real(kind=8) 3D-array.
92      !!
93      !! The input file:
94      !!
95      !! - must have a regular number of columns, that is each data line must have the same number
96      !!   of columns (according to the delimiter used).
97      !! - must use a SINGLE empty line for "depth" block separator.
98      !! - can contains any number of comment lines (i.e. line where first non-blank character is "#").
99      !!   All other lines (except empty lines) are assumed to be data.
100      !!
101      !! Error occured when:
102      !!    - Path does not refer to a existing file (-11)
103      !!    - No free logical unit available (-1)
104      !!    - The file does not have regular data-columns number (-5)
105      !!    - At least a value cannot be cast in double precision (-10)
106      !!
107      !! The method assumes the input files consists in _D_ block of _R_ lines
108      !! with _C_ columns. Each block must be separated by a single empty line and
109      !! each columns must separated by one or more blank spaces (no tabulation ALLOWED).
110      !!
111      !! On success, the shape of the 3D output array will be _output(R,C,D)_.
112      !! On error, the 3D output array is __not allocated__.
113      CHARACTER(len=*), INTENT(in)                             :: path      !! Path of the input data file
114      REAL(kind=8), INTENT(out), DIMENSION(:,:,:), ALLOCATABLE :: data3d    !! 3D-array with the output values (double precision)
115      CHARACTER(len=*), INTENT(in), OPTIONAL                   :: delimiter !! Optional column delimiter(s)
116      TYPE(error) :: err                                                    !! Error status of the function
117      LOGICAL                                           :: ok
118      INTEGER                                           :: i,lc,tlc
119      INTEGER                                           :: ndr,ndc,ndd
120      INTEGER                                           :: ir,jc,kd,lu
121      REAL(kind=8), DIMENSION(:), ALLOCATABLE           :: tmp
122      CHARACTER(len=5)                                  :: slc
123      CHARACTER(len=15)                                 :: i2s
124      CHARACTER(len=:), ALLOCATABLE                     :: line,lm1,zdelim
125      CHARACTER(len=st_slen), DIMENSION(:), ALLOCATABLE :: wrds
126 
127      zdelim = CHAR(9)//CHAR(32)
128      IF (PRESENT(delimiter)) zdelim = delimiter
129      ! Check file status
130      INQUIRE(FILE=TRIM(path),EXIST=ok)
131      IF (.NOT.ok) THEN
132        err = error(trim(path)//": no such file",-1) ; RETURN
133      ENDIF
134      lu = free_lun()
135      IF (lu == -1) THEN ; err = error("Cannot find available logical unit...",-1) ; RETURN ; ENDIF
136      ! Open file
137      OPEN(lu,FILE=TRIM(path),STATUS='OLD',ACTION='READ')
138 
139      ! First pass :
140      ! ------------
141      !   - get size (rows, columns, depth)
142      !   - check size consistendcy
143      !   - check value type
144      lc = 0 ; tlc = 0
145      ndr = -1 ; ndc = -1 ; ndd = 1
146      DO WHILE(readline(lu,line))
147        lm1 = line
148        ! Read the line
149        lc = lc + 1 ; WRITE(slc,'(I5)') lc ; slc = ADJUSTL(slc)
150        ! skip comment line
151        IF (INDEX(TRIM(ADJUSTL(line)),"#") == 1) CYCLE
152        ! An empty line: new 2D block
153        IF (LEN_TRIM(line) == 0) THEN
154          ndd = ndd + 1
155          IF (ndr < 0) THEN
156            ndr = tlc
157          ELSEIF (tlc /= ndr) THEN
158            WRITE(i2s,'(I15)') ndd ; i2s = ADJUSTL(i2s)
159            err = error(trim(path)//":Invalid number of lines in block #"//i2s//"(line "//TRIM(slc)//")",-5)
160            RETURN
161          ENDIF
162          tlc = 0
163          CYCLE
164        ENDIF
165        tlc = tlc + 1
166        ! Splits line in words
167        IF (.NOT.tokenize(line,wrds,zdelim,.true.)) THEN
168          ! cannot tokenize
169          err = error(trim(path)//": Cannot parse line "//TRIM(slc),-5)
170          RETURN
171        ELSEIF (.NOT.from_string(wrds,tmp)) THEN
172          ! cannot cast values
173          err = error(trim(path)//": Cannot cast values at line "//TRIM(slc),-5)
174          RETURN
175        ELSEIF (ndc > 0 .AND. ndc /= SIZE(wrds)) THEN
176          ! current number of columns not equal to last one
177          err = error(trim(path)//": Invalid number of columns (line "//TRIM(slc)//")",-5)
178          RETURN
179        ENDIF
180        IF (ndc == -1) ndc = SIZE(wrds)
181        IF (ALLOCATED(wrds)) DEALLOCATE(wrds)
182        IF (ALLOCATED(tmp)) DEALLOCATE(tmp)
183      ENDDO
184 
185      ! NOTE:
186      ! there is a possible bug if data file ends by an empty line:
187      ! we will have an extra empty bloc !
188      !   current patch: we save the last line of the file and check it:
189      !     - if we have empty line, we reduce ndd by one.
190      IF (LEN_TRIM(lm1) == 0) ndd = ndd-1
191 
192      ! Rewind input data file
193      REWIND(lu)
194      ! Allocate memory
195      ALLOCATE(data3d(ndr,ndc,ndd))
196      ir = 0 ; kd = 1 ;
197      DO WHILE(readline(lu,line))
198        IF (INDEX(TRIM(ADJUSTL(line)),"#") == 1) CYCLE
199        ir = ir + 1
200        ! empty line update block subscripts
201        IF (LEN_TRIM(line) == 0) THEN
202          kd = kd + 1 ; ir = 0 ; CYCLE
203        ENDIF
204        ok = tokenize(line,wrds,zdelim,.true.)
205        DO jc = 1,ndc ; ok = from_string(wrds(jc),data3d(ir,jc,kd)) ; ENDDO
206        IF (ALLOCATED(wrds)) DEALLOCATE(wrds)
207      ENDDO
208      CLOSE(lu)
209    END FUNCTION read_data_3d
210 
211    FUNCTION read_data_2d(path,data2d,delimiter) RESULT(err)
212      !! Read an input ASCII data file and stores its content in a 2D array
213      !!
214      !! The function reads an ASCII file and saves its values in a real(kind=8) 2D-array.
215      !!
216      !! The input file:
217      !!
218      !! - can contains any number of empty lines and/or comment line (i.e. line where first
219      !!   non-blank character is "#"). All other lines are assumed to be data.
220      !! - must have a regular number of columns, that is each data line must have the same
221      !!   number of columns.
222      !! - must use blank space(s) as value delimiter.
223      !!
224      !! Error occured when:
225      !!
226      !! - Path does not refer to a existing file (-1)
227      !! - No free logical unit available (-1)
228      !! - The file does not have regular data-columns number (-5)
229      !! - At least a value cannot be cast in double precision (-5)
230      !!
231      !! On error, the 2D output array is __not allocated__.
232      USE SWIFT_FSYSTEM
233      CHARACTER(len=*), INTENT(in)                           :: path      !! Path of the input data file
234      REAL(kind=8), INTENT(out), DIMENSION(:,:), ALLOCATABLE :: data2d    !! 2D-array with the output values (double precision)
235      CHARACTER(len=*), INTENT(in), OPTIONAL                 :: delimiter !! Optional column delimiter(s)
236      TYPE(error) :: err                                                  !! Error status of function.
237      LOGICAL                                           :: ok
238      INTEGER                                           :: i,e,vc,lc
239      INTEGER                                           :: nl,nc,lu
240      REAL(kind=8), DIMENSION(:), ALLOCATABLE           :: tmp
241      CHARACTER(len=5)                                  :: slc
242      CHARACTER(len=:), ALLOCATABLE                     :: line,zdelim
243      CHARACTER(len=st_slen), DIMENSION(:), ALLOCATABLE :: wrds
244      zdelim = CHAR(9)//CHAR(32)
245      IF (PRESENT(delimiter)) zdelim = delimiter
246      INQUIRE(FILE=TRIM(path),EXIST=ok)
247      IF (.NOT.ok) THEN
248        err = error(trim(path)//": no such file",-1) ; RETURN
249      ENDIF
250      lu = free_lun()
251      IF (lu == -1) THEN ; err = error("Cannot find available logical unit...",-1) ; RETURN ; ENDIF
252      OPEN(lu,FILE=TRIM(path),STATUS='OLD',ACTION='READ')
253      vc = 0 ; lc=0 ; ok = .true.
254      ! Read the file twice :)
255      lc=0 ; vc = 0 ; nc=-1
256      ! First pass : get number of row values and checks everything !
257      DO
258        ! Read the line
259        IF (.NOT.readline(lu,line)) EXIT
260        lc = lc + 1 ; WRITE(slc,'(I5)') lc ; slc = ADJUSTL(slc)
261        ! skip empty/comment line
262        IF (INDEX(TRIM(ADJUSTL(line)),"#") == 1.OR.LEN_TRIM(line) == 0) CYCLE
263        ! update row counter
264        vc = vc + 1
265        IF (.NOT.tokenize(line,wrds,zdelim,.true.)) THEN
266          ! cannot tokenize
267          err = error(trim(path)//": Cannot parse line "//TRIM(slc),-5)
268          RETURN
269        ELSEIF (.NOT.from_string(wrds,tmp)) THEN
270          ! cannot cast values
271          err = error(trim(path)//": Cannot cast values at line "//TRIM(slc),-5)
272          RETURN
273        ELSEIF (nc > 0 .AND. nc /= SIZE(tmp)) THEN
274          ! current number of columns not equal to last one
275          err = error(trim(path)//": Invalid number of columns (line "//TRIM(slc)//")",-5)
276          RETURN
277        ENDIF
278        IF (nc == -1) nc = SIZE(wrds)
279        IF (ALLOCATED(wrds)) DEALLOCATE(wrds)
280        IF (ALLOCATED(tmp)) DEALLOCATE(tmp)
281      ENDDO
282      ! Rewind input data file
283      REWIND(lu)
284      nl = vc
285      ! allocate memory
286      ALLOCATE(data2d(nl,nc))
287      ! Second pass : saves values :)
288      vc = 0
289      DO WHILE(vc <= nl)
290        ! Reads the line
291        IF (.NOT.readline(lu,line)) EXIT
292        ! Check if we have comment or null string
293        IF (INDEX(TRIM(ADJUSTL(line)),"#") == 1.OR.LEN_TRIM(line) == 0) CYCLE
294        vc = vc + 1
295        ok = tokenize(line,wrds,zdelim,.true.)
296        DO i = 1,nc
297          ok = from_string(wrds(i),data2d(vc,i))
298        ENDDO
299        IF (ALLOCATED(wrds)) DEALLOCATE(wrds)
300      ENDDO
301      CLOSE(lu)
302      RETURN
303    END FUNCTION read_data_2d
304 
305    FUNCTION readline(lun,line) RESULT(not_eof)
306      !! Read a complete line
307      !!
308      !! Each time, it is called, the function reads a complete of the file opened in __lun__
309      !! logical unit and returns .false. if EOF has been reached, .true. otherwise.
310      !!
311      !! The function is intended to read a file line by line:
312      !!
313      !! ```
314      !! lu = 1
315      !! open(lu,file="path/to/the/file/to/read")
316      !! l=0   ! A line counter
317      !! DO WHILE(readline(lu,line))
318      !!   l = l + 1
319      !!   WRITE(*,'("L",I2.2,": ",(a))') il,line
320      !! ENDDO
321      !! CLOSE(1)
322      !! ```
323      INTEGER, INTENT(in)                        :: lun  !! Logical unit with the opened file to read.
324      CHARACTER(len=:), ALLOCATABLE, INTENT(out) :: line !! Processed line
325      LOGICAL  :: not_eof                                !! .true. if EOF has NOT been reached yet, .false. otherwise
326      CHARACTER(len=50) :: buf
327      INTEGER           :: e,sz
328      not_eof = .true. ; line = ''
329      DO
330        READ(lun,'(a)',ADVANCE="no",SIZE=sz,IOSTAT=e) buf
331        IF (e == IOSTAT_END) THEN
332          not_eof = .false.
333          IF (sz > 0) line=line//buf(1:sz)
334          EXIT
335        ELSE IF (e == IOSTAT_EOR) THEN
336          line = line//buf(1:sz)
337          EXIT
338        ELSE
339          line = line//buf
340        ENDIF
341      ENDDO
342    END FUNCTION readline
343 
344  END MODULE SWIFT_ASCIIREAD 
Note: See TracBrowser for help on using the repository browser.