ROMS
Loading...
Searching...
No Matches
yaml_parser_test.F
Go to the documentation of this file.
2!
3!git $Id$
4!================================================== Hernan G. Arango ===
5! Copyright (c) 2002-2025 The ROMS Group !
6! Licensed under a MIT/X style license !
7! See License_ROMS.md !
8!=======================================================================
9! !
10! This module contains several routines to process input YAML files. !
11! !
12! Notice that several Fortran parsers exist for complex and simple !
13! YAML files coded with Object-Oriented Programming (OOP) principles. !
14! For example, check: !
15! !
16! * FCKit (https://github.com/ecmwf/fckit) !
17! * Fortran-YAML (https://github.com/BoldingBruggeman/fortran-yaml) !
18! * yaFyaml (https://github.com/Goddard-Fortran-Ecosystem/yaFyaml) !
19! !
20! However, this YAML parser is more uncomplicated with substantial !
21! capabilities. It is a hybrid between standard and OOP principles !
22! but without the need for recurrency, inheritance, polymorphism, !
23! and containers. !
24! !
25! The only constraint is that the YAML file is read twice for !
26! simplicity. The first read determines the number indentation of !
27! blanks policy and the length of the collection list(:) pairs object !
28! (CLASS yaml_pair). It supports: !
29! !
30! * Single or multiple line comments start with a hash '#'. Also, !
31! comment after a key/value pair is allowed. All comments are !
32! skipped during processing. !
33! * Unlimited nested structure (lists, mappings, hierarchies). !
34! Indentation of whitespace is used to denote structure. !
35! * Unrestricted schema indentation. However, some schema validators !
36! recommend or impose two whitespace indentations. !
37! * A key is followed by a colon to denote a mapping value (like !
38! ocean_model: ROMS). !
39! * Aliases and Anchors. !
40! * Blocking lists: members are denoted by a leading hyphen and !
41! space, which is considered as part of the indentation. !
42! * Flow sequence: a vector list with values enclosed in square !
43! brackets and separated by a comma-and-space: [val1, ..., valN]. !
44! * Keyword values are processed and stored as strings but converted !
45! to a logical, integer, or floating-point type when appropriate !
46! during extraction. If derived-type values are needed, the caller !
47! can process such structure outside this module, as shown below. !
48! * Remove unwanted control characters like tabs and separators !
49! (ASCII character code 0-31) !
50! * English uppercase and lowercase alphabet, but it can be expanded !
51! to other characters (see yaml_ValueType routine) !
52! * Module is self contained, but it has very minimal association to !
53! four ROMS modules. !
54! * Multiple or continuation lines are supported, for example we can !
55! have: !
56! !
57! state variables: [sea_surface_height_above_geopotential_datum, !
58! barotropic_sea_water_x_velocity, !
59! barotropic_sea_water_y_velocity, !
60! sea_water_x_velocity, !
61! sea_water_y_velocity, !
62! sea_water_potential_temperature, !
63! sea_water_practical_salinity] !
64! !
65! Usage: !
66! !
67! USE yaml_parser_mod, ONLY : yaml_initialize, yaml_get, yaml_tree !
68! !
69! TYPE (yaml_tree) :: YML !
70! !
71! CALL yaml_initialize (YML, 'ocn_coupling.yaml', report) !
72! status = yaml_get(YML, 'Nimport', Nimport) !
73! status = yaml_get(YML, 'import.standard_name', Sstandard) !
74! status = yaml_get(YML, 'import.standard_name.short_name', S%short)!
75! status = yaml_get(YML, 'import.standard_name.unit', S%units) !
76! ... !
77! !
78! and so on for logical, integer, floating-point, and string !
79! key/value pairs. !
80! !
81! Here, 'Sstandard(1:Nimport) will contain all 'standard_name' !
82! values for the YAML block 'import:'. Notice that nested keywords !
83! are concatenated with a period: 'key1.key2.key3' for a three- !
84! level nested block, similar to how Matlab build structures. The !
85! key can have more than one word separated by one space. For !
86! example, we can have 'bulk_flux import.standard_name'. Similarly, !
87! any otherkey/value pair can be extrated from the YML object. !
88! !
89!=======================================================================
90!
91! Define local parameters.
92!
93 logical :: yaml_master = .true. ! Master process
94 logical :: lreport = .true. ! dumps YAML dictionaty
95!
96 integer, parameter :: kind_real = selected_real_kind(12,300)
97 integer, parameter :: dp = selected_real_kind(12,300)
98
99 integer, parameter :: yaml_stdout = 6 ! standard output
100 integer, parameter :: stdout = 6 ! standard output
101
102 integer, parameter :: noerror = 0 ! No error flag
103
104 integer :: yaml_errflag = 0 ! error flag
105!
106!-----------------------------------------------------------------------
107! Structures/Objects to hold YAML dictionary lists with theirs keys and
108! values.
109!-----------------------------------------------------------------------
110!
111! YAML file key/value pair.
112!
113 TYPE, PUBLIC :: yaml_pair
114
115 logical :: has_alias ! alias '*' token
116 logical :: has_anchor ! anchor '&' token
117 logical :: is_block ! block '-' list
118 logical :: is_sequence ! sequence '[]' tokens
119!
120 logical :: is_logical ! logical value
121 logical :: is_integer ! integer value
122 logical :: is_real ! numerical value
123 logical :: is_string ! string value
124!
125 integer :: id ! key/value ID
126 integer :: parent_id ! parent ID
127 integer :: left_padding ! indent level: 0,1,..
128!
129 character (len=:), allocatable :: line ! YAML line
130 character (len=:), allocatable :: key ! YAML keyword:
131 character (len=:), allocatable :: value ! YAML value(s)
132 character (len=:), allocatable :: anchor ! anchor keyword
133
134 END TYPE yaml_pair
135!
136! YAML file dictionary tree.
137!
138 TYPE, PUBLIC :: yaml_tree
139
140 logical :: iscreated = .false. ! Object creation switch
141
142 integer :: nbranches ! total number of branches
143 integer :: npairs ! total number of pairs
144 integer :: indent ! blank indentation policy
145!
146 character (len=:), allocatable :: filename ! YAML file name
147!
148! YAML file collection pairs, [1:Npairs].
149!
150 TYPE (yaml_pair), pointer :: list(:)
151!
152 CONTAINS ! CLASS objects
153!
154 PROCEDURE :: create => yaml_tree_create
155 PROCEDURE :: destroy => yaml_tree_destroy
156 PROCEDURE :: dump => yaml_tree_dump
157 PROCEDURE :: extract => yaml_tree_extract
158 PROCEDURE :: fill => yaml_tree_fill
159 PROCEDURE :: fill_aliases => yaml_tree_fill_aliases
160 PROCEDURE :: has => yaml_tree_has
161 PROCEDURE :: read_line => yaml_tree_read_line
162
163 END TYPE yaml_tree
164!
165! Public structures that can be used in applications to extract block
166! list YAML constructs. The key may represents a sequence flow [...]
167! with a vector of values. The values can be integers, logicals, reals,
168! or strings. For example,
169!
170! import:
171! - standard_name: surface_downward_heat_flux_in_sea_water
172! long_name: surface net heat flux
173! short_name: shflux
174! data_variables: [shf, shf_time]
175! - standard_name: surface_wind_x_stress
176! long_name: surface zonal wind stress component
177! short_name: sustr
178! data_variables: [taux, atm_time]
179! - standard_name: surface_wind_y_stress
180! long_name: surface meridional wind stress component
181! short_name: svstr
182! data_variables: [tauy, atm_time]
183!
184! The extraction is loaded into V, which is a TYPE "yaml_Svec"
185! structure:
186!
187! status = yaml_get(YML, 'import.data_variables', V)
188!
189! or altenatively
190!
191! IF (yaml_Error(yaml_get(YML, 'import.data_variables', V), &
192! & NoErr, __LINE__, MyFile)) RETURN
193!
194! yielding the following block-list, string structure in a single
195! invocation of the overloaded function "yaml_get":
196!
197! V(1)%vector(1)%value = 'shf'
198! V(1)%vector(2)%value = 'shf_time'
199! V(2)%vector(1)%value = 'taux'
200! V(2)%vector(2)%value = 'atm_time'
201! V(3)%vector(1)%value = 'tauy'
202! V(3)%vector(2)%value = 'atm_time'
203!
204! It is a compact way to extract similar data blocks.
205!
206 TYPE, PUBLIC :: yaml_ivec ! integer structure
207 integer, allocatable :: vector(:) ! vector values
208 END TYPE yaml_ivec
209!
210 TYPE, PUBLIC :: yaml_lvec ! logical structure
211 logical, allocatable :: vector(:) ! vector values
212 END TYPE yaml_lvec
213!
214 TYPE, PUBLIC :: yaml_rvec ! real structure
215 real (kind_real), allocatable :: vector(:) ! vector values
216 END TYPE yaml_rvec
217!
218 TYPE, PUBLIC :: yaml_svec ! string structure
219 character (len=:), allocatable :: value ! scalar value
220 TYPE (yaml_svec), pointer :: vector(:) ! recursive vector
221 END TYPE yaml_svec ! values
222!
223! Derived-type structure, extended/inherited from parent "yaml_Svec",
224! is used to extract hierarchies of keys and associated values from
225! YAML dictionary object. The calling program specifies a key-string
226! that may be generated by aggregating nested keys with a period.
227! Also, it can extract flow sequence string element values that are
228! separated by commas.
229!
230 TYPE, PRIVATE, EXTENDS(yaml_svec) :: yaml_extract
231 logical :: has_vector ! true if loaded vector values
232 END TYPE yaml_extract
233!
234!-----------------------------------------------------------------------
235! Define public overloading API to extract key/value pairs from YAML
236! tree dictionary object accounting to variable type.
237!-----------------------------------------------------------------------
238!
239 INTERFACE yaml_get
240 MODULE PROCEDURE yaml_get_i_struc ! Gets integer structure
241 MODULE PROCEDURE yaml_get_l_struc ! Gets logical structure
242 MODULE PROCEDURE yaml_get_r_struc ! Gets real structure
243 MODULE PROCEDURE yaml_get_s_struc ! Gets string structure
244!
245 MODULE PROCEDURE yaml_get_ivar_0d ! Gets integer value
246 MODULE PROCEDURE yaml_get_ivar_1d ! Gest integer values
247 MODULE PROCEDURE yaml_get_lvar_0d ! Gets logical value
248 MODULE PROCEDURE yaml_get_lvar_1d ! Gets logical values
249 MODULE PROCEDURE yaml_get_rvar_0d ! Gets real value
250 MODULE PROCEDURE yaml_get_rvar_1d ! Gets real values
251 MODULE PROCEDURE yaml_get_svar_0d ! Gets string value
252 MODULE PROCEDURE yaml_get_svar_1d ! Gets string values
253 END INTERFACE yaml_get
254!
255!
256!-----------------------------------------------------------------------
257! Define generic coupling field to process import and export states.
258!-----------------------------------------------------------------------
259!
261 logical :: connected
262 logical :: debug_write
263!
264 real(dp) :: add_offset
265 real(dp) :: scale
266!
267 character (len=:), allocatable :: connected_to
268 character (len=:), allocatable :: data_netcdf_vname
269 character (len=:), allocatable :: data_netcdf_tname
270 character (len=:), allocatable :: destination_grid
271 character (len=:), allocatable :: destination_units
272 character (len=:), allocatable :: extrapolate_method
273 character (len=:), allocatable :: long_name
274 character (len=:), allocatable :: map_norm
275 character (len=:), allocatable :: map_type
276 character (len=:), allocatable :: regrid_method
277 character (len=:), allocatable :: source_units
278 character (len=:), allocatable :: source_grid
279 character (len=:), allocatable :: short_name
280 character (len=:), allocatable :: standard_name
281!
282 END TYPE couplingfield
283!
284! <><><><><><><><><><><><><><><><><><><><><><><><><><><><><><><><><><><>
285! YAML Tree Post-proccesing Variables and Structures:
286!
287! Define generic YAML dictionary, containers, and counters used
288! during processing.
289!-----------------------------------------------------------------------
290!
291 logical :: ldebugmetadata = .true. ! Debugging switch
292!
293! Counters.
294!
295 integer :: ientry ! entry counter
296 integer :: nentries ! number of entries
297!
298! logical scalar dummy values.
299!
300 logical, allocatable :: ylogical1(:)
301!
302! Real scalar dummy values.
303!
304 real(kind_real), allocatable :: yreal1(:)
305 real(kind_real), allocatable :: yreal2(:)
306!
307! Derived-type dummy structures for processing string value or set
308! of values from a sequence flow, [val1, ..., valN].
309!
310 TYPE (yaml_svec), allocatable :: ystring1 (:)
311 TYPE (yaml_svec), allocatable :: ystring2 (:)
312 TYPE (yaml_svec), allocatable :: ystring3 (:)
313 TYPE (yaml_svec), allocatable :: ystring4 (:)
314 TYPE (yaml_svec), allocatable :: ystring5 (:)
315 TYPE (yaml_svec), allocatable :: ystring6 (:)
316 TYPE (yaml_svec), allocatable :: ystring7 (:)
317 TYPE (yaml_svec), allocatable :: ystring8 (:)
318 TYPE (yaml_svec), allocatable :: ystring9 (:)
319 TYPE (yaml_svec), allocatable :: ystring10(:)
320 TYPE (yaml_svec), allocatable :: ystring11(:)
321 TYPE (yaml_svec), allocatable :: ystring12(:)
322!
323 PUBLIC :: cmeps_metadata
324 PUBLIC :: coupling_metadata
325 PUBLIC :: process_yaml
326!
327! <><><><><><><><><><><><><><><><><><><><><><><><><><><><><><><><><><><>
328!
329 PUBLIC :: yaml_assignstring
330 PUBLIC :: yaml_error
331 PUBLIC :: yaml_get
332 PUBLIC :: yaml_initialize
333!
334 PRIVATE :: yaml_countkeys
335 PRIVATE :: yaml_lowercase
336 PRIVATE :: yaml_uppercase
337 PRIVATE :: yaml_valuetype
338!
339 PRIVATE
340!
341!-----------------------------------------------------------------------
342! Local module parameters.
343!-----------------------------------------------------------------------
344!
345 logical, parameter :: ldebugyaml = .false. ! debugging switch
346!
347 integer, parameter :: ldim = 8 ! logical Lswitch dimension
348 integer, parameter :: lkey = 254 ! Maximum characters per keyword
349 integer, parameter :: lmax = 2048 ! Maximum characters per line
350 integer, parameter :: noerr = 0 ! no error flag
351 integer, parameter :: iunit = 222 ! Fortan unit for reading
352!
353 character (len=55), dimension(7) :: yaml_errmsg = &
354 & (/ ' YAML_PARSER - Blows up ................ yaml_ErrFlag: ', &
355 & ' YAML_PARSER - Input error ............. yaml_ErrFlag: ', &
356 & ' YAML_PARSER - Output error ............ yaml_ErrFlag: ', &
357 & ' YAML_PARSER - I/O error ............... yaml_ErrFlag: ', &
358 & ' YAML_PARSER - Configuration error ..... yaml_ErrFlag: ', &
359 & ' YAML_PARSER - Partition error ......... yaml_ErrFlag: ', &
360 & ' YAML_PARSER - Illegal input parameter . yaml_ErrFlag: ' /)
361!
362!-----------------------------------------------------------------------
363 CONTAINS
364!-----------------------------------------------------------------------
365!
366 FUNCTION yaml_initialize (self, filename, report) RESULT (status)
367!
368!***********************************************************************
369! !
370! It creates a YAML tree dictionary object. First, it reads the YAML !
371! file to determine the indentation policy and length of list oject !
372! (TYPE yaml_tree). !
373! !
374! After the object is allocated, the Fortran unit is rewinded and the !
375! YAML file is read again to populate the keyword values. !
376! !
377! On Input: !
378! !
379! self YAML tree dictionary object (CLASS yaml_tree) !
380! filename YAML filename (string) !
381! report Switch to dump to standard output (logical, OPTIONAL) !
382! !
383! On Ouptut: !
384! !
385! self Allocated and populated YAML object. !
386! status Error flag (integer) !
387! !
388!***********************************************************************
389!
390! Imported variable declarations.
391!
392 class(yaml_tree), intent(inout) :: self
393 character (len=*), intent(in ) :: filename
394 logical, optional, intent(in ) :: report
395!
396! Local variable declarations.
397!
398 logical :: ldump
399!
400 integer :: lenstr, status
401!
402 character (len=*), parameter :: myfile = &
403 & __FILE__//", yaml_initialize"
404!
405!-----------------------------------------------------------------------
406! Initialize YAML object.
407!-----------------------------------------------------------------------
408!
409 status=noerr
410!
411! Set switch to print the processed YAML key/value pairs to standard
412! ouput.
413!
414 IF (PRESENT(report)) THEN
415 ldump = report
416 ELSE
417 ldump = .false.
418 END IF
419!
420! Set YAML file path and name.
421!
422 IF (yaml_error(yaml_assignstring(self%filename, &
423 & filename, lenstr), &
424 & noerr, __line__, myfile)) THEN
425 status=yaml_errflag
426 RETURN
427 END IF
428!
429! Create and populate YAML object.
430!
431 IF (self%IsCreated) CALL self%destroy ()
432!
433 CALL self%create ()
434 IF (yaml_error(yaml_errflag, noerr, __line__, myfile)) THEN
435 status=yaml_errflag
436 RETURN
437 END IF
438!
439! Report YAML tree dictionary, for debugging.
440!
441 IF (ldump) CALL self%dump ()
442!
443 RETURN
444 END FUNCTION yaml_initialize
445!
446 SUBROUTINE yaml_tree_create (self)
447!
448!***********************************************************************
449! !
450! It creates a YAML tree dictionary object. First, it reads the YAML !
451! file to determine the dimensions of parent and children structures. !
452! After the structures are allocate, the Fortran unit is rewinded and !
453! YAML file is read again to populate the keyword values. !
454! !
455! On Input: !
456! !
457! self YAML tree dictionary object (CLASS yaml_tree) !
458! !
459! On Ouptut: !
460! !
461! self Allocated and populated YAML object. !
462! !
463!***********************************************************************
464!
465! Imported variable declarations.
466!
467 class(yaml_tree), intent(inout) :: self
468!
469! Local variable declarations.
470!
471 logical :: FirstPass, Lswitch(Ldim)
472!
473 integer :: Nblanks, Nbranches, Npairs
474 integer :: i, io_err, status
475!
476 character (len=Lkey) :: anchor, io_errmsg, key
477 character (len=Lmax) :: line, value
478!
479 character (len=*), parameter :: MyFile = &
480 & __FILE__//", yaml_tree_create"
481!
482!-----------------------------------------------------------------------
483! Open YAML file for reading.
484!-----------------------------------------------------------------------
485!
486 OPEN (unit=iunit, file=self%filename, form='formatted', &
487 & status='old', iostat=io_err, iomsg=io_errmsg)
488 IF (io_err.ne.0) THEN
490 IF (yaml_error(yaml_errflag, noerr, __line__, myfile)) THEN
491 IF (yaml_master) WRITE (yaml_stdout,10) self%filename, &
492 & trim(io_errmsg)
493 RETURN
494 END IF
495 END IF
496!
497! Determine the total number of YAML key/value pairs.
498!
499 nblanks = 0
500 nbranches = 0
501 npairs = 0
502!
503 firstpass = .true.
504!
505 yaml_line : DO WHILE (.true.)
506!
507 status=self%read_line(nblanks, line, key, value, anchor, &
508 & lswitch)
509!
510 SELECT CASE (status)
511 CASE (-1) ! error condition during reading
513 IF (yaml_error(yaml_errflag, noerr, __line__, myfile)) THEN
514 IF (yaml_master) WRITE (yaml_stdout,20) self%filename, &
515 & trim(line)
516 RETURN
517 END IF
518 CASE (0) ! end-of-file
519 EXIT
520 CASE (1) ! blank or comment line, move to the next line
521 cycle
522 CASE (2) ! processed viable line
523 npairs=npairs+1
524 END SELECT
525!
526! If no leading blanks, advance YAML tree branch counter.
527!
528 IF (nblanks.eq.0) THEN
529 nbranches = nbranches+1 ! hierarchy start
530 END IF
531!
532! On first pass, establish indentation policy: number of blanks.
533!
534 IF (firstpass.and.(nblanks.gt.0)) THEN
535 firstpass=.false.
536 self%indent=nblanks
537 END IF
538!
539! Check for consitent identation policy. Some YAML validators impose
540! a two-blank spacing policy.
541!
542 IF (nblanks.gt.0) THEN
543 IF (mod(nblanks, self%indent).ne.0) THEN
545 IF (yaml_error(yaml_errflag, noerr, __line__, myfile)) THEN
546 IF (yaml_master) WRITE (yaml_stdout,30) nblanks, &
547 & self%indent
548 RETURN
549 END IF
550 END IF
551 END IF
552
553 END DO yaml_line
554!
555! Rewind unit since we need to reprocess the file again to load the
556! data into the allocated "self%list(1:Npairs)" container.
557!
558 rewind(iunit)
559!
560!-----------------------------------------------------------------------
561! Allocate YAML dictionary container.
562!-----------------------------------------------------------------------
563!
564! Set number of main branches and number of key/values in YAML 'list'.
565!
566 self%Nbranches=nbranches
567 self%Npairs=npairs
568!
569! Allocate YAML key/value pair list (TYPE 'yaml_pair') object.
570!
571 IF (.not.self%IsCreated) THEN
572 ALLOCATE ( self%list(npairs) )
573 self%IsCreated=.true.
574 END IF
575!
576!-----------------------------------------------------------------------
577! Re-read YAML file again and populate dictionary.
578!-----------------------------------------------------------------------
579!
580 CALL self%fill ()
581!
582 10 FORMAT (/,' YAML_TREE_CREATE - Unable to open input YAML file: ', &
583 & a,/,20x,'ERROR: ',a)
584 20 FORMAT (/,' YAML_TREE_CREATE - Error while reading YAML file: ', &
585 & a,/,20x,'LINE: ',a)
586 30 FORMAT (/,' YAML_TREE_CREATE - Inconsistent indentation, ', &
587 & 'self%indent = ',i0,', nblanks = ',i0)
588 40 FORMAT (/,' YAML_TREE_CREATE - Inconsistent indentation, ', &
589 & 'nblanks = ',i0,', indent blank policy = ',i0,/,20x, &
590 & 'Number of nested collections = ',i0, &
591 & ' is greater than 3',/,20x,'Line: ',a)
592 50 FORMAT (/,' YAML_TREE_CREATE - Cannot find branches in YAML ', &
593 & 'file: ',a)
594!
595 RETURN
596 END SUBROUTINE yaml_tree_create
597!
598 SUBROUTINE yaml_tree_destroy (self)
599!
600!***********************************************************************
601! !
602! It deallocates/destroys the YAML dictionary object. !
603! !
604! Arguments: !
605! !
606! self YAML tree dictionary object (CLASS yaml_tree) !
607! !
608!***********************************************************************
609!
610! Imported variable declarations.
611!
612 class(yaml_tree), intent(inout) :: self
613!
614! Local variable declarations.
615!
616 logical :: Lopened
617!
618!-----------------------------------------------------------------------
619! Deallocate YAML dictionary object.
620!-----------------------------------------------------------------------
621!
622! If applicable, close YAML Fortran unit.
623!
624 INQUIRE (unit=iunit, opened=lopened)
625 IF (lopened) THEN
626 CLOSE (iunit)
627 END IF
628!
629! Recall that Fortran 2003 standard allows deallocating just the
630! parent object to deallocate all array variables within its scope
631! automatically.
632!
633 IF (ASSOCIATED(self%list)) THEN
634 self%IsCreated=.false.
635 DEALLOCATE (self%list)
636 END IF
637!
638 RETURN
639 END SUBROUTINE yaml_tree_destroy
640!
641 SUBROUTINE yaml_tree_dump (self)
642!
643!***********************************************************************
644! !
645! It prints the YAML dictionary to standard output for debugging !
646! purposes. !
647! !
648! Arguments: !
649! !
650! self YAML tree dictionary object (CLASS yaml_tree) !
651! !
652!***********************************************************************
653!
654! Imported variable declarations.
655!
656 class(yaml_tree), intent(in) :: self
657!
658! Local variable declarations.
659!
660 integer :: Lstr, Nblanks, i, indent, padding
661!
662 character (len=Lkey) :: key
663 character (len=Lmax) :: string, value
664!
665!-----------------------------------------------------------------------
666! Report the contents of the YAML tree directory.
667!-----------------------------------------------------------------------
668!
669 IF (yaml_master) THEN
670 WRITE (yaml_stdout,10) self%filename
671 indent=self%indent
672!
673 DO i=1,self%Npairs
674 padding=self%list(i)%left_padding
675 nblanks=indent*padding
676 key=self%list(i)%key
677 IF (ALLOCATED(self%list(i)%value)) THEN
678 value=self%list(i)%value
679 lstr=len_trim(value)
680 ELSE
681 lstr=0
682 END IF
683 IF (self%list(i)%is_block) THEN
684 IF (lstr.gt.0) THEN
685 IF (self%list(i)%is_sequence) THEN
686 WRITE (string,20) i, repeat(char(32),nblanks), '- ', &
687 & trim(key), ': [', trim(value), ']'
688 ELSE
689 WRITE (string,30) i, repeat(char(32),nblanks), '- ', &
690 & trim(key), ': ', trim(value)
691 END IF
692 ELSE
693 WRITE (string,40) i, repeat(char(32),nblanks), '- ', &
694 & trim(key), ':'
695 END IF
696 ELSE
697 IF (lstr.gt.0) THEN
698 IF (self%list(i)%is_sequence) THEN
699 WRITE (string,30) i, repeat(char(32),nblanks), &
700 & trim(key), ': [', trim(value), ']'
701 ELSE
702 WRITE (string,40) i, repeat(char(32),nblanks), &
703 & trim(key), ': ', trim(value)
704 END IF
705 ELSE
706 WRITE (string,50) i, repeat(char(32),nblanks), &
707 & trim(key), ':'
708 END IF
709 END IF
710 WRITE (yaml_stdout,60) trim(string)
711 END DO
712 END IF
713!
714 10 FORMAT (/,"YAML Tree Dictinary, file: '",a,"'",/, &
715 & '==========================',/)
716 20 FORMAT ('L=',i4.4,1x,'% ',6a)
717 30 FORMAT ('L=',i4.4,1x,'% ',5a)
718 40 FORMAT ('L=',i4.4,1x,'% ',4a)
719 50 FORMAT ('L=',i4.4,1x,'% ',3a)
720 60 FORMAT (a)
721!
722 RETURN
723 END SUBROUTINE yaml_tree_dump
724!
725 FUNCTION yaml_tree_extract (self, keystring, S) RESULT (status)
726!
727!***********************************************************************
728! !
729! It extracts YAML value(s) from processing key-string. The key !
730! string may be a set of keys aggregated with a period, CHAR(46). !
731! !
732! On Input: !
733! !
734! self YAML tree dictionary object (CLASS yaml_tree) !
735! keystring Aggregated YAML keys (string) !
736! !
737! On Output: !
738! !
739! S Separated YAML key/pair value (TYPE yaml_extract) !
740! nkeys Number of extracted keys (integer) !
741! !
742!***********************************************************************
743!
744! Imported variable declarations.
745!
746 class(yaml_tree), intent(in) :: self
747 character (len=*), intent(in) :: keystring
748 TYPE (yaml_extract), allocatable, intent(inout) :: s(:)
749!
750! Local variable declarations.
751!
752 TYPE (yaml_extract), allocatable :: k(:)
753!
754 logical :: blockflow
755!
756 integer :: i, ib, ic, ie, ik, ipair, is, j, li, pid
757 integer :: lstr, lenstr, nkeys, npairs, nvalues
758 integer :: icomma, idot
759 integer :: status
760!
761 integer, allocatable :: p(:) ! pair index
762!
763 character (len=:), allocatable :: kstring ! key string
764 character (len=:), allocatable :: vstring ! value string
765!
766 character (len=*), parameter :: myfile = &
767 & __FILE__//", yaml_tree_extract"
768!
769!-----------------------------------------------------------------------
770! Extract YAML key(s) from key-string.
771!-----------------------------------------------------------------------
772!
773 status=noerr
774!
775! Count the number keys in key-string, separated by period.
776!
777 lstr=len_trim(keystring)
778 IF (yaml_error(yaml_assignstring(kstring, &
779 & keystring, lenstr), &
780 & noerr, __line__, myfile)) RETURN
781!
782 nkeys=yaml_countkeys(kstring, char(46))
783!
784! Allocate key structure.
785!
786 ALLOCATE ( k(nkeys) )
787!
788! Extract keys.
789!
790 is=1
791 DO i=1,nkeys
792 idot=index(kstring,char(46),back=.false.)
793 ie=idot
794 IF (idot.eq.0) THEN
795 ie=len_trim(kstring)
796 ELSE
797 ie=ie-1
798 END IF
799 IF (yaml_error(yaml_assignstring(k(i)%value, &
800 & kstring(is:ie), lenstr), &
801 & noerr, __line__, myfile)) RETURN
802 IF (idot.gt.0) kstring(is:ie+1) = repeat(char(32), ie-is+2)
803 kstring=trim(adjustl(kstring))
804 END DO
805!
806!-----------------------------------------------------------------------
807! Determine the number of YAML tree pairs to process and assoicated
808! list array element.
809!-----------------------------------------------------------------------
810!
811! Check if key-string is a blocking list where any of the members has
812! a leading hyphen as indentation. If blocking, all the pairs with the
813! same key-string are extracted.
814!
815 blockflow=.false. ! Switch to process all block membert
816 ib=0 ! block list(s) counter
817 ic=0 ! key/value pair counter
818 ik=1 ! key counter to avoid double DO-loops
819!
820 DO i=1,self%Npairs
821 lstr=len_trim(self%list(i)%key)
822 IF ((self%list(i)%key).eq.(k(ik)%value)) THEN
823
824 IF (yaml_master.and.ldebugyaml) THEN
825 print '(2(a,i0,2a))', 'key ',ik,' = ', trim(k(ik)%value), &
826 & ', YAML list ',i,' = ', &
827 & trim(self%list(i)%key)
828 END IF
829
830 pid=self%list(i)%parent_id
831 IF (self%list(i)%is_block.or.self%list(pid)%is_block) THEN
832 ib=ib+1 ! advance block counter
833 END IF
834 IF (ik.eq.nkeys) THEN
835 ic=ic+1 ! advance pair counter
836 IF (ib.eq.0) THEN
837 li=i ! list index to extract
838 EXIT ! no blocking list found
839 ELSE
840 blockflow=.true.
841 END IF
842 ELSE
843 ik=ik+1 ! advance key counter
844 END IF
845 END IF
846 IF (blockflow.and.(self%list(i)%left_padding.eq.0)) THEN
847 EXIT ! processed all blocks
848 END IF
849 END DO
850 npairs=ic ! pairs to process
851!
852! Allocate pair index array, P.
853!
854 IF (npairs.ne.0) THEN
855 IF (.not.ALLOCATED(p)) ALLOCATE ( p(npairs) )
856 ELSE
858 status=yaml_errflag
859 IF (yaml_error(yaml_errflag, noerr, __line__, myfile)) THEN
860 IF (yaml_master) WRITE (yaml_stdout,10) keystring, &
861 & self%filename
862 RETURN
863 END IF
864 END IF
865!
866! Set pair element(s) to extract.
867!
868 IF (blockflow) THEN
869 ic=0
870 ik=1
871 DO i=1,self%Npairs
872 IF ((self%list(i)%key).eq.(k(ik)%value)) THEN
873 IF (ik.eq.nkeys) THEN
874 ic=ic+1
875 p(ic)=i ! multiple pair index
876 ELSE
877 ik=ik+1
878 END IF
879 END IF
880 IF ((ic.gt.0).and.(self%list(i)%left_padding.eq.0)) THEN
881 EXIT ! processed all blocks
882 END IF
883 END DO
884 ELSE
885 p(1)=li ! single pair index
886 END IF
887!
888!-----------------------------------------------------------------------
889! Extract pair(s) value(s).
890!-----------------------------------------------------------------------
891!
892 DO i=1,npairs
893 ipair=p(i)
894!
895! Get key-string value.
896!
897 IF (yaml_error(yaml_assignstring(vstring, &
898 & self%list(ipair)%value, &
899 & lenstr), &
900 & noerr, __line__, myfile)) RETURN
901!
902!-----------------------------------------------------------------------
903! Extract/load keys-tring value(s). In a sequence, values are separated
904! by commas.
905!-----------------------------------------------------------------------
906!
907! Extract pair from a blocking list. Currently, scalar value are
908! allowed. Nested blocking is not supported since it is complicated
909! to process. For example, the following nested blocking lists are
910! not supported.
911!
912! branch:
913! - blocklist1: value
914! blocklist1_key1: value
915! - blocklist1A: value ! not supported
916!
917 IF (blockflow) THEN
918!
919! Process a block list pair with vector values (flow sequence):
920! S(i)%vector(1)%value to S(i)%vector(nvalues)%value
921!
922 IF (self%list(ipair)%is_sequence) THEN
923 lstr=len_trim(vstring)
924 nvalues=count((/(vstring(j:j), j=1,lstr)/) == char(44)) + 1
925!
926 IF (.not.ALLOCATED(s)) THEN
927 ALLOCATE ( s(npairs) ) ! main structure
928 END IF
929 ALLOCATE ( s(i)%vector(nvalues) ) ! sub-structure
930 s(i)%has_vector=.true.
931!
932 is=1
933 DO j=1,nvalues
934 icomma=index(vstring,char(44),back=.false.)
935 ie=icomma
936 IF (icomma.eq.0) THEN
937 ie=len_trim(vstring)
938 ELSE
939 ie=ie-1
940 END IF
941 IF (yaml_error(yaml_assignstring(s(i)%vector(j)%value, &
942 & vstring(is:ie), &
943 & lenstr), &
944 & noerr, __line__, myfile)) RETURN
945
946 IF (yaml_master.and.ldebugyaml) THEN
947 print '(3a,2(i0,a),a)', 'keystring = ',trim(keystring), &
948 & ', S(', i, ')%vector(', j, ') = ', &
949 & trim(s(i)%vector(j)%value)
950 END IF
951
952 IF (icomma.gt.0) vstring(is:ie+1)=repeat(char(32),ie-is+2)
953 vstring=trim(adjustl(vstring))
954 END DO
955!
956! Process a block list pair with a scalar value, S(i)%value
957!
958 ELSE
959!
960 IF (.not.ALLOCATED(s)) THEN
961 ALLOCATE ( s(npairs) )
962 END IF
963 s(i)%has_vector=.false.
964!
965 IF (yaml_error(yaml_assignstring(s(i)%value, &
966 & vstring, lenstr), &
967 & noerr, __line__, myfile)) RETURN
968
969 IF (yaml_master.and.ldebugyaml) THEN
970 print '(a,i0,4a)', 'keystring ',i,' = ', trim(keystring), &
971 & ', value = ', trim(s(i)%value)
972 END IF
973 END IF
974!
975! Otherwise, process a non-block list.
976!
977 ELSE
978!
979! Process a flow sequence: S(1)%value to S(nvalues)%value.
980!
981 IF (self%list(ipair)%is_sequence) THEN
982 lstr=len_trim(vstring)
983 nvalues=count((/(vstring(j:j), j=1,lstr)/) == char(44)) + 1
984!
985 ALLOCATE ( s(nvalues) )
986 s(i)%has_vector=.false.
987!
988 is=1
989 DO j=1,nvalues
990 icomma=index(vstring,char(44),back=.false.)
991 ie=icomma
992 IF (icomma.eq.0) THEN
993 ie=len_trim(vstring)
994 ELSE
995 ie=ie-1
996 END IF
997 IF (yaml_error(yaml_assignstring(s(j)%value, &
998 & vstring(is:ie), &
999 & lenstr), &
1000 & noerr, __line__, myfile)) RETURN
1001 IF (icomma.gt.0) vstring(is:ie+1)=repeat(char(32),ie-is+2)
1002 vstring=trim(adjustl(vstring))
1003
1004 IF (yaml_master.and.ldebugyaml) THEN
1005 print '(a,i0,4a)', 'keystring ',j,' = ', &
1006 & trim(keystring), &
1007 & ', value = ', trim(s(j)%value)
1008 END IF
1009 END DO
1010!
1011! Process a single scalar value, S(1)%value.
1012!
1013 ELSE
1014!
1015 ALLOCATE ( s(1) )
1016 s(1)%has_vector=.false.
1017!
1018 IF (yaml_error(yaml_assignstring(s(1)%value, &
1019 & vstring, lenstr), &
1020 & noerr, __line__, myfile)) RETURN
1021
1022 IF (yaml_master.and.ldebugyaml) THEN
1023 print '(4a)', 'keystring = ', trim(keystring), &
1024 & ', value = ', trim(s(1)%value)
1025 END IF
1026 END IF
1027 END IF
1028 END DO
1029!
1030!-----------------------------------------------------------------------
1031! Deallocate local variables.
1032!-----------------------------------------------------------------------
1033!
1034 IF (ALLOCATED(k)) DEALLOCATE (k)
1035 IF (ALLOCATED(p)) DEALLOCATE (p)
1036 IF (ALLOCATED(kstring)) DEALLOCATE (kstring)
1037 IF (ALLOCATED(vstring)) DEALLOCATE (vstring)
1038!
1039 10 FORMAT (/," YAML_TREE_EXTRACT - Cannot find key-string: '",a, &
1040 & "'",/,21x,'File: ',a)
1041 20 FORMAT (/," YAML_TREE_EXTRACT - Not supported key-string: '",a, &
1042 & "'",/,21x,'nested sub-blocking in a leading blocking ', &
1043 & 'list',/,21x,'File: ',a)
1044!
1045 RETURN
1046 END FUNCTION yaml_tree_extract
1047!
1048 SUBROUTINE yaml_tree_fill (self)
1049!
1050!***********************************************************************
1051! !
1052! It reads YAML file and fills structure with the key/value pairs. !
1053! Both the key and value pairs are strings. The numercal convertions !
1054! are done elsewhere when needed. !
1055! !
1056! Arguments: !
1057! !
1058! self YAML tree dictionary object (CLASS yaml_tree) !
1059! !
1060!***********************************************************************
1061!
1062! Imported variable declarations.
1063!
1064 class(yaml_tree), intent(inout) :: self
1065!
1066! Local variable declarations.
1067!
1068 logical :: Lswitch(Ldim)
1069!
1070 integer :: Nblanks, LenStr, left_padding, new_parent, old_parent
1071 integer :: i, ibranch, icount, ic_alias, ic_anchor, status
1072!
1073 character (len=Lkey) :: anchor, key
1074 character (len=Lmax) :: line, value
1075!
1076 character (len=*), parameter :: MyFile = &
1077 & __FILE__//", yaml_tree_fill"
1078!
1079!-----------------------------------------------------------------------
1080! Read and populate YAML structure.
1081!-----------------------------------------------------------------------
1082!
1083 ibranch=0
1084 icount=0
1085 ic_alias=0
1086 ic_anchor=0
1087 new_parent=0
1088 old_parent=0
1089!
1090 yaml_line : DO WHILE (.true.)
1091!
1092 status=self%read_line(nblanks, line, key, value, anchor, &
1093 & lswitch)
1094!
1095 SELECT CASE (status)
1096 CASE (-1) ! error condition during reading
1097 yaml_errflag=4
1098 IF (yaml_error(yaml_errflag, noerr, __line__, myfile)) THEN
1099 IF (yaml_master) WRITE (yaml_stdout,10) self%filename, &
1100 & trim(line)
1101 RETURN
1102 END IF
1103 CASE (0) ! end-of-file
1104 EXIT
1105 CASE (1) ! blank or comment line, move to the next line
1106 cycle
1107 CASE (2) ! processed viable line
1108 icount=icount+1
1109 END SELECT
1110!
1111! Determine structure indices according to the nested levels counters.
1112!
1113 IF (nblanks.eq.0) THEN
1114 ibranch=ibranch+1
1115 new_parent=icount
1116 old_parent=icount
1117 END IF
1118!
1119! Load YAML pair switch identifiers.
1120!
1121 self%list(icount)%has_alias = lswitch(1)
1122 self%list(icount)%has_anchor = lswitch(2)
1123 self%list(icount)%is_block = lswitch(3)
1124 self%list(icount)%is_sequence = lswitch(4)
1125 self%list(icount)%is_logical = lswitch(5)
1126 self%list(icount)%is_integer = lswitch(6)
1127 self%list(icount)%is_real = lswitch(7)
1128 self%list(icount)%is_string = lswitch(8)
1129!
1130 IF (lswitch(1)) ic_alias=ic_alias+1
1131 IF (lswitch(2)) ic_anchor=ic_anchor+1
1132!
1133! Set left-padding indentation level: 0, 1, 2, ...
1134!
1135 left_padding=nblanks/self%indent
1136 self%list(icount)%left_padding=left_padding
1137!
1138! Load key/value ID and parent ID.
1139!
1140 IF (nblanks.gt.0) THEN
1141 IF (left_padding.gt.self%list(icount-1)%left_padding) THEN
1142 new_parent=old_parent
1143 old_parent=icount
1144 END IF
1145 END IF
1146 self%list(icount)%id=icount
1147 self%list(icount)%parent_id=new_parent
1148!
1149! Allocate and load line, key, and value strings. If applicable, loal
1150! anchor value. Notice that it is possible to have keyword without
1151! value in the main branches.
1152!
1153 IF (yaml_error(yaml_assignstring(self%list(icount)%line, &
1154 & line, lenstr), &
1155 & noerr, __line__, myfile)) RETURN
1156!
1157 IF (yaml_error(yaml_assignstring(self%list(icount)%key, &
1158 & key, lenstr), &
1159 & noerr, __line__, myfile)) RETURN
1160!
1161 IF (len_trim(value).gt.0) THEN
1162 IF (yaml_error(yaml_assignstring(self%list(icount)%value, &
1163 & value, lenstr), &
1164 & noerr, __line__, myfile)) RETURN
1165 END IF
1166!
1167 IF (lswitch(2).and.len_trim(anchor).gt.0) THEN
1168 IF (yaml_error(yaml_assignstring(self%list(icount)%anchor, &
1169 & anchor, lenstr), &
1170 & noerr, __line__, myfile)) RETURN
1171 END IF
1172!
1173 END DO yaml_line
1174!
1175! Substitute 'Alias' with 'Anchor' values. Notice that Anchors and
1176! Aliases let you identify an item with an 'anchor' in a YAML pair,
1177! and then refer to that item with an alias later in the same tile.
1178! It is done to avoid repetitions and errors.
1179!
1180 IF ((ic_alias.gt.0).and.(ic_anchor.gt.0)) THEN
1181 CALL self%fill_aliases (ic_alias, ic_anchor)
1182 END IF
1183!
1184! Close YAML file.
1185!
1186 CLOSE (iunit)
1187!
1188 10 FORMAT (/,' YAML_TREE_FILL - Error while reading YAML file: ', &
1189 & a,/,18x,'LINE: ',a)
1190!
1191 RETURN
1192 END SUBROUTINE yaml_tree_fill
1193!
1194 SUBROUTINE yaml_tree_fill_aliases (self, Nalias, Nanchor)
1195!
1196!***********************************************************************
1197! !
1198! It replaces the Aliases items with its respective Anchors values in !
1199! YAML dictionary. !
1200! !
1201! Arguments: !
1202! !
1203! self YAML tree dictionary object (CLASS yaml_tree) !
1204! Nalias Number of Aliases items in YAML file (integer) !
1205! Nanchors Number of Anchors in YAML file (integer) !
1206! !
1207!***********************************************************************
1208!
1209! Imported variable declarations.
1210!
1211 class(yaml_tree), intent(inout) :: self
1212 integer, intent(in) :: Nalias
1213 integer, intent(in) :: Nanchor
1214!
1215! Local variable declarations.
1216!
1217 logical :: Lswitch(Ldim)
1218!
1219 integer :: Ialias, LenStr, i, j, ic
1220!
1221 character (len=Lkey), dimension(Nanchor) :: AnchorKey, AnchorVal
1222 character (len=Lmax) :: AliasVal
1223!
1224 character (len=*), parameter :: MyFile = &
1225 & __FILE__//", yaml_tree_fill_aliases"
1226!
1227!-----------------------------------------------------------------------
1228! Extract Anchors keyword and value.
1229!-----------------------------------------------------------------------
1230!
1231 ic=0
1232 DO i=1,self%Npairs
1233 IF (self%list(i)%has_anchor) THEN
1234 ic=ic+1
1235 anchorkey(ic)=self%list(i)%anchor
1236 anchorval(ic)=self%list(i)%value
1237 END IF
1238 END DO
1239!
1240!-----------------------------------------------------------------------
1241! Replace Aliases with associate Anchors values.
1242!-----------------------------------------------------------------------
1243!
1244 DO j=1,self%Npairs
1245!
1246! Get Aliases keyword.
1247!
1248 IF (self%list(j)%has_alias) THEN
1249 aliasval=self%list(j)%value
1250 ialias=index(aliasval,char(42),back=.false.) ! alias '*'
1251 IF (ialias.gt.0) THEN
1252 aliasval(ialias:ialias)=char(32) ! blank
1253 aliasval=trim(adjustl(aliasval))
1254 END IF
1255!
1256! Replace Aliases with anchor values and update value type.
1257!
1258 DO i=1,ic
1259 IF (trim(aliasval).eq.trim(anchorkey(i))) THEN
1260 DEALLOCATE (self%list(j)%value)
1261 IF (yaml_error(yaml_assignstring(self%list(j)%value, &
1262 & trim(anchorval(i)), &
1263 & lenstr), &
1264 & noerr, __line__, myfile)) RETURN
1265!
1266 lswitch=.false.
1267 CALL yaml_valuetype (self%list(j)%value, lswitch)
1268 self%list(j)%is_logical = lswitch(5)
1269 self%list(j)%is_integer = lswitch(6)
1270 self%list(j)%is_real = lswitch(7)
1271 self%list(j)%is_string = lswitch(8)
1272 END IF
1273 END DO
1274 END IF
1275 END DO
1276!
1277 RETURN
1278 END SUBROUTINE yaml_tree_fill_aliases
1279!
1280 FUNCTION yaml_tree_has (self, keystring) RESULT (FoundIt)
1281!
1282!***********************************************************************
1283! !
1284! It checks if YAML dictionary has requested key string. !
1285! !
1286! On Input: !
1287! !
1288! self YAML tree dictionary object (CLASS yaml_tree) !
1289! !
1290! keystring Requested YAML key-string (string) !
1291! !
1292! On Output: !
1293! !
1294! FoundIt Switch indicating if the key string was found or not !
1295! in the YAML dictionary. !
1296! !
1297!***********************************************************************
1298!
1299! Imported variable declarations.
1300!
1301 class(yaml_tree), intent(inout) :: self
1302 character (len=*), intent(in ) :: keystring
1303!
1304! Local variable declarations.
1305!
1306 TYPE (yaml_extract), allocatable :: k(:)
1307!
1308 logical :: foundit
1309!
1310 integer :: lstr, lenstr
1311 integer :: i, idot, ie, is, j, nkeys
1312!
1313 character (len=:), allocatable :: kstring
1314
1315 character (len=*), parameter :: myfile = &
1316 & __FILE__//", yaml_tree_has"
1317!
1318!-----------------------------------------------------------------------
1319! Check if requested key-string is available in YAML dictionary
1320!-----------------------------------------------------------------------
1321!
1322 foundit=.false.
1323!
1324! Count number of keys in key-string, separated by period, CHAR(46).
1325!
1326 lstr=len_trim(keystring)
1327 IF (yaml_error(yaml_assignstring(kstring, &
1328 & keystring, lenstr), &
1329 & noerr, __line__, myfile)) RETURN
1330!
1331 nkeys=yaml_countkeys(kstring, char(46))
1332!
1333! Check single key.
1334!
1335 IF (nkeys.eq.1) THEN
1336 DO i=1,SIZE(self%list)
1337 IF (self%list(i)%key.eq.kstring) THEN
1338 foundit=.true.
1339 EXIT
1340 END IF
1341 END DO
1342 DEALLOCATE (kstring)
1343 RETURN
1344!
1345! Otherwise, check for compound key separated by period.
1346!
1347 ELSE
1348 ALLOCATE ( k(nkeys) )
1349!
1350! Extract keys.
1351!
1352 is=1
1353 DO j=1,nkeys
1354 idot=index(kstring,char(46),back=.false.)
1355 ie=idot
1356 IF (idot.eq.0) THEN
1357 ie=len_trim(kstring)
1358 ELSE
1359 ie=ie-1
1360 END IF
1361 IF (yaml_error(yaml_assignstring(k(j)%value, &
1362 & kstring(is:ie), lenstr), &
1363 & noerr, __line__, myfile)) RETURN
1364 IF (idot.gt.0) kstring(is:ie+1) = repeat(char(32), ie-is+2)
1365 kstring=trim(adjustl(kstring))
1366 END DO
1367!
1368! Check for compound key: val1.val2 ...
1369! (It fails if any of the compound keys is not found)
1370!
1371 is=1
1372 DO j=1,nkeys
1373 foundit=.false.
1374 DO i=is,SIZE(self%list)
1375 IF (self%list(i)%key.eq.k(j)%value) THEN
1376 foundit=.true.
1377 is=i+1
1378 EXIT
1379 END IF
1380 END DO
1381 IF (.not.foundit) EXIT
1382 END DO
1383 DEALLOCATE (k, kstring)
1384 RETURN
1385 END IF
1386!
1387 END FUNCTION yaml_tree_has
1388!
1389 FUNCTION yaml_tree_read_line (self, Nblanks, line, key, value, &
1390 & anchor, Lswitch) &
1391 & result(status)
1392!
1393!***********************************************************************
1394! !
1395! It reads a single text line from current YAML file. It uses the !
1396! ASCII character set extensively. !
1397! !
1398! On Input: !
1399! !
1400! self YAML tree dictionary object (CLASS yaml_tree) !
1401! !
1402! On Output: !
1403! !
1404! Nblanks YAML leading blanks identation (integer) !
1405! !
1406! line Current YAML file line (string) !
1407! !
1408! key YAML line keyword (string) !
1409! !
1410! value YAML line keyword value (string) !
1411! !
1412! anchor YAML value ancher keyword (string) !
1413! !
1414! Lswitch YAML key/value switches (logical, vector) !
1415! !
1416! Lswitch(1) = T/F, value has an alias '*' token !
1417! Lswitch(2) = T/F, value has an anchor '&' token !
1418! Lswitch(3) = T/F, key/value starts a block '-' !
1419! Lswitch(4) = T/F, value has sequence '[...]' tokens !
1420! Lswitch(5) = T/F, logical value(s) !
1421! Lswitch(6) = T/F, integer value(s) !
1422! Lswitch(7) = T/F, floating-point value(s) !
1423! Lswitch(8) = T/F, string value(s) !
1424! !
1425!***********************************************************************
1426!
1427! Imported variable declarations.
1428!
1429 class(yaml_tree), intent(inout) :: self
1430!
1431 logical, intent(out) :: lswitch(:)
1432!
1433 integer, intent(out) :: nblanks
1434!
1435 character (len=*), intent(out) :: line
1436 character (len=*), intent(out) :: key
1437 character (len=*), intent(out) :: value
1438 character (len=*), intent(out) :: anchor
1439!
1440! Local variable declarations.
1441!
1442 logical :: lbracket, rbracket
1443!
1444 integer :: ialias, ianchor, iblank, icolon, idash, ihash, ispace
1445 integer :: ibracketl, ibracketr
1446 integer :: lstr, lstrnext, lstrval, i, j, k
1447 integer :: status
1448!
1449 character (len=Lmax) :: linecopy, next_line
1450
1451 character (len=*), parameter :: myfile = &
1452 & __FILE__//", yaml_tree_read_line"
1453!
1454!-----------------------------------------------------------------------
1455! Read a single YAML file line
1456!-----------------------------------------------------------------------
1457!
1458! Initialize.
1459!
1460 status=-1 ! error condition
1461!
1462 DO i=1,len(key)
1463 key(i:i)=char(32)
1464 END DO
1465
1466 DO i=1,len(value)
1467 value(i:i)=char(32)
1468 END DO
1469
1470 DO i=1,len(anchor)
1471 anchor(i:i)=char(32)
1472 END DO
1473!
1474 lswitch=.false.
1475 lbracket=.false.
1476 rbracket=.false.
1477!
1478! Read in next YAML file line.
1479!
1480 READ (iunit,'(a)',err=10,END=20) line
1481!
1482! Replace illegal control characters with a blank space CHAR(32)
1483!
1484 lstr=len_trim(line)
1485
1486 DO i=1,len_trim(line)
1487 j=ichar(line(i:i))
1488 IF (j.lt.32) THEN
1489 line(i:i)=char(32)
1490 END IF
1491 END DO
1492 linecopy=trim(adjustl(line))
1493!
1494! Get length of "line". Remove comment after the KEYWORD, if any.
1495! In YAML, a comment starts with a hash (#), CHAR(35).
1496!
1497 IF ((lstr.gt.0).and.(linecopy(1:1).ne.char(35))) THEN
1498 ihash=index(line,char(35),back=.false.)
1499 IF (ihash.gt.0) THEN
1500 lstr=ihash-1
1501 line=trim(line(1:lstr)) ! remove trailing comment
1502 lstr=len_trim(line)
1503 END IF
1504 status=2 ! viable line
1505 ELSE
1506 status=1 ! blank or commented line,
1507 RETURN ! move to the next line
1508 END IF
1509!
1510! Find colon sign CHAR(58) and determine the number of leading blank
1511! spaces. YAML uses indentations with one or more spaces to describe
1512! nested collections (lists, mappings). YAML also uses dash plus space,
1513! '- ', as enumerator of block lists.
1514!
1515! It checks if the KEYWORD is a multi-word separated by space.
1516!
1517 icolon=index(line,char(58),back=.false.)
1518 IF (icolon.eq.0) THEN
1519 status=-1
1520 yaml_errflag=2
1521 IF (yaml_master) THEN
1522 WRITE (yaml_stdout,30) trim(line)
1523 END IF
1524 IF (yaml_error(yaml_errflag, noerr, __line__, myfile)) RETURN
1525 END IF
1526!
1527 idash =index(line,char(45),back=.false.)
1528 IF ((idash.gt.0).and.(idash.lt.icolon)) THEN
1529 iblank=index(line(1:idash),char(32),back=.true.)
1530 ELSE
1531 iblank=index(line(1:icolon),char(32),back=.true.)
1532 IF (iblank.gt.0) THEN
1533 k=iblank-1
1534 IF ((65.le.ichar(line(1:1))).and. &
1535 & (ichar(line(1:1)).le.122)) THEN ! multi-word branch key
1536 iblank=0
1537 ELSE IF ((65.le.ichar(line(k:k))).and. &
1538 & (ichar(line(k:k)).le.122)) THEN ! multi-word pair key
1539 iblank=index(line(1:k),char(32),back=.true.)
1540 END IF
1541 END IF
1542 END IF
1543!
1544! Set number of YAML line leading blacks.
1545!
1546 nblanks=iblank
1547!
1548! Extract key and value pair and return.
1549!
1550 IF ((idash.gt.0).and.(idash.lt.icolon)) THEN
1551 key=trim(adjustl(line(idash+1:icolon-1)))
1552 ELSE
1553 key=trim(adjustl(line(1:icolon-1)))
1554 END IF
1555 value=trim(adjustl(line(icolon+1:lstr)))
1556!
1557! Check for special YAML tokens in value string and replace with blank.
1558!
1559 ialias=index(value,char(42),back=.false.) ! alias '*'
1560 IF (ialias.gt.0) THEN
1561 lswitch(1)=.true.
1562 END IF
1563!
1564 ianchor=index(value,char(38),back=.false.) ! anchor '&'
1565 IF (ianchor.gt.0) THEN
1566 ispace=index(value(ianchor+1:),char(32),back=.false.)
1567 anchor=value(ianchor+1:ispace) ! anchor value
1568 lswitch(2)=.true.
1569 value(ianchor:ispace)=repeat(char(32),ispace-ianchor+1)
1570 value=trim(adjustl(value))
1571 END IF
1572!
1573 IF ((idash.gt.0).and.(idash.lt.icolon)) THEN
1574 lswitch(3)=.true. ! block pair '- '
1575 END IF
1576!
1577 ibracketl=index(value,char(91),back=.false.) ! left bracket '['
1578 IF (ibracketl.gt.0) THEN
1579 lbracket=.true.
1580 value(ibracketl:ibracketl)=char(32)
1581 value=trim(adjustl(value))
1582 END IF
1583!
1584 ibracketr=index(value,char(93),back=.false.) ! right bracket ']'
1585 IF (ibracketr.gt.0) THEN
1586 rbracket=.true.
1587 value(ibracketr:ibracketr)=char(32)
1588 value=trim(adjustl(value))
1589 END IF
1590!
1591!-----------------------------------------------------------------------
1592! If right square bracket is not found, the key values are broken into
1593! multiple lines. Process the necessary lines.
1594!-----------------------------------------------------------------------
1595!
1596 IF (.not.rbracket.and.lbracket) THEN
1597 DO WHILE (.not.rbracket)
1598 READ (iunit,'(a)',err=10,END=20) next_line
1599!
1600! Replace illegal control characters with a blank space CHAR(32)
1601!
1602 DO i=1,len_trim(next_line)
1603 j=ichar(next_line(i:i))
1604 IF (j.lt.32) THEN
1605 next_line(i:i)=char(32)
1606 END IF
1607 END DO
1608 next_line=trim(adjustl(next_line))
1609!
1610! If applicable, remove trailing comments starting with a hash (#),
1611! CHAR(35).
1612!
1613 ihash=index(next_line,char(35),back=.false.)
1614 lstrnext=len_trim(next_line)
1615 IF ((lstrnext.gt.0).and.(ihash.gt.0)) THEN
1616 lstrnext=ihash-1
1617 next_line=trim(next_line(1:lstrnext))
1618 lstrnext=len_trim(next_line)
1619 END IF
1620!
1621! Aggregate new 'next_line' to previous 'line' and 'value'.
1622!
1623 lstr=len_trim(line)
1624 lstrval=len_trim(value)
1625 line=line(1:lstr)//char(32)//next_line(1:lstrnext)
1626 value=value(1:lstrval)//char(32)//next_line(1:lstrnext)
1627!
1628 ibracketr=index(value,char(93),back=.false.)
1629 IF (ibracketr.gt.0) THEN
1630 rbracket=.true.
1631 value(ibracketr:ibracketr)=char(32)
1632 value=trim(adjustl(value))
1633 END IF
1634 END DO
1635 END IF
1636 IF (lbracket.and.rbracket) lswitch(4)=.true.
1637!
1638!-----------------------------------------------------------------------
1639! Determine value representation: logical (boolean), numerical
1640! (integer/reals), or string(s). Others can be added as needed.
1641!-----------------------------------------------------------------------
1642!
1643 CALL yaml_valuetype (value, lswitch)
1644!
1645 RETURN
1646!
1647! Read flow control determines the status flag.
1648!
1649 10 status=-1 ! error during reading
1650 20 status=0 ! end-of-file encoutered
1651!
1652 30 FORMAT (/,' YAML_TREE_READ_LINE - Unable to find colon token ', &
1653 & 'after key word',/,23x,'LINE: ',a)
1654
1655 END FUNCTION yaml_tree_read_line
1656!
1657 FUNCTION yaml_assignstring (OutString, InpString, Lstr) &
1658 & result(errflag)
1659!
1660!=======================================================================
1661! !
1662! It assigns allocatable strings. It allocates/reallocates output !
1663! string variable. !
1664! !
1665! On Input: !
1666! !
1667! InpString String to be assigned (character) !
1668! !
1669! On Output: !
1670! !
1671! OutString Assigned allocatable string (character) !
1672! Lstr Length allocated string (integer) !
1673! ErrFlag Error flag (integer) !
1674! !
1675!=======================================================================
1676!
1677! Imported variable declarations.
1678!
1679 character (len=:), allocatable, intent(inout) :: outstring
1680 character (len=*), intent(in) :: inpstring
1681 integer, intent(out) :: lstr
1682!
1683! Local variable declarations.
1684!
1685 integer :: errflag
1686!
1687!-----------------------------------------------------------------------
1688! Allocate output string to the size of input string.
1689!-----------------------------------------------------------------------
1690!
1691 errflag = -1
1692!
1693 lstr=len_trim(inpstring)
1694 IF (.not.allocated(outstring)) THEN
1695 allocate ( character(LEN=Lstr) :: outstring, stat=errflag)
1696 ELSE
1697 deallocate (outstring)
1698 allocate ( character(LEN=Lstr) :: outstring, stat=errflag)
1699 END IF
1700!
1701! Assign requested value.
1702!
1703 outstring = inpstring
1704!
1705 RETURN
1706 END FUNCTION yaml_assignstring
1707!
1708 FUNCTION yaml_countkeys (string, token) RESULT (nkeys)
1709!
1710!=======================================================================
1711! !
1712! It counts the number of concatenated key separated by the specified !
1713! token during processing extraction from YAML object. The same task !
1714! can be done elegantly as: !
1715! !
1716! nkeys=COUNT((/ (string(i:i), i=1,Lstr) /) == token) + 1 !
1717! !
1718! But compilier like 'gfortran' cannot handle such abstraction. !
1719! !
1720! On Input: !
1721! !
1722! string Aggregated YAML keys (string) !
1723! token Key separator (string) !
1724! !
1725! On Output: !
1726! !
1727! nkeys Number of aggregated keys (integer) !
1728! !
1729!=======================================================================
1730!
1731! Imported variable declarations.
1732!
1733 character (len=*), intent(in) :: string
1734 character (len=*), intent(in) :: token
1735!
1736! Local variable declarations.
1737!
1738 integer :: nkeys
1739 integer :: lstr, i
1740!
1741!-----------------------------------------------------------------------
1742! Count number of concatenated keys in input string.
1743!-----------------------------------------------------------------------
1744!
1745 nkeys=1
1746 lstr=len_trim(string)
1747 DO i=1,lstr
1748 IF (string(i:i).eq.token) THEN
1749 nkeys=nkeys+1
1750 END IF
1751 END DO
1752!
1753 RETURN
1754 END FUNCTION yaml_countkeys
1755!
1756 FUNCTION yaml_error (flag, NoErr, line, routine) RESULT (foundit)
1757!
1758!=======================================================================
1759! !
1760! It checks YAML execution flag against no-error code and issue a !
1761! message if they are not equal. !
1762! !
1763! On Input: !
1764! !
1765! flag YAML execution flag (integer) !
1766! NoErr No Error code (integer) !
1767! line Calling model routine line (integer) !
1768! routine Calling model routine (string) !
1769! !
1770! On Output: !
1771! !
1772! foundit The value of the result is TRUE/FALSE if the !
1773! execution flag is in error. !
1774! !
1775!=======================================================================
1776!
1777! Imported variable declarations.
1778!
1779 integer, intent(in) :: flag, noerr, line
1780
1781 character (len=*), intent(in) :: routine
1782!
1783! Local variable declarations.
1784!
1785 logical :: foundit
1786!
1787!-----------------------------------------------------------------------
1788! Scan array for requested string.
1789!-----------------------------------------------------------------------
1790!
1791 foundit=.false.
1792 IF (flag.ne.noerr) THEN
1793 foundit=.true.
1794 IF (yaml_master) THEN
1795 WRITE (yaml_stdout,10) flag, line, trim(routine)
1796 10 FORMAT (' Found Error: ', i2.2, t20, 'Line: ', i0, &
1797 & t35, 'Source: ', a)
1798 END IF
1799 FLUSH (yaml_stdout)
1800 END IF
1801
1802 RETURN
1803 END FUNCTION yaml_error
1804!
1805 FUNCTION yaml_get_i_struc (self, keystring, V) RESULT (status)
1806!
1807!***********************************************************************
1808! !
1809! It loads a vector set of integers in YAML block-list structure, !
1810! V(1:Nitems)%vector(1:Nvalues). If the dummy argument V structure !
1811! is allocated, it deallocates/allocates the required Nitems and !
1812! Nvalues dimensions. !
1813! !
1814! On Input: !
1815! !
1816! self YAML tree dictionary object (CLASS yaml_tree) !
1817! keystring YAML tree aggregated keys (string) !
1818! !
1819! On Output: !
1820! !
1821! V Vector of integers in block list (TYPE yaml_Ivec) !
1822! status Error code (integer) !
1823! !
1824!***********************************************************************
1825!
1826! Imported variable declarations.
1827!
1828 class(yaml_tree), intent(in) :: self
1829 character (len=*), intent(in) :: keystring
1830 TYPE (yaml_ivec), allocatable, intent(out) :: v(:)
1831!
1832! Local variable declarations.
1833!
1834 TYPE (yaml_extract), allocatable :: s(:)
1835!
1836 integer :: lenstr, nitems, nvalues, i, n
1837 integer :: status
1838!
1839 character (len=Lmax) :: msg
1840
1841 character (len=*), parameter :: myfile = &
1842 & __FILE__//", yaml_Get_i_struc"
1843!
1844!-----------------------------------------------------------------------
1845! Extract requested key-string values.
1846!-----------------------------------------------------------------------
1847!
1848 status=noerr
1849!
1850 IF (yaml_error(self%extract(keystring, s), &
1851 & noerr, __line__, myfile)) RETURN
1852!
1853! Allocate output structure.
1854!
1855 nitems=SIZE(s, dim=1)
1856 IF (ALLOCATED(v)) DEALLOCATE (v)
1857 ALLOCATE ( v(nitems) )
1858!
1859! Convert string vector values to integers.
1860!
1861 DO n=1,nitems
1862 IF (s(n)%has_vector) THEN
1863 nvalues=SIZE(s(n)%vector)
1864 ALLOCATE ( v(n)%vector(nvalues) )
1865 DO i=1,nvalues
1866 READ (s(n)%vector(i)%value, * ,iostat=status, iomsg=msg) &
1867 & v(n)%vector(i)
1868 IF (yaml_error(status, noerr, __line__, myfile)) THEN
1869 yaml_errflag=5
1870 IF (yaml_master) WRITE (yaml_stdout,10) trim(keystring), &
1871 & s(n)%vector(i)%value, &
1872 & trim(msg)
1873 RETURN
1874 END IF
1875 END DO
1876 END IF
1877 END DO
1878!
1879! Deallocate local extraction structure.
1880!
1881 IF (ALLOCATED(s)) DEALLOCATE (s)
1882!
1883 10 FORMAT (/,' YAML_GET_I_STRUC - Error while converting string to', &
1884 & ' integer, key: ',a,', value = ',a,/,20x,'ErrMsg: ',a)
1885!
1886 RETURN
1887 END FUNCTION yaml_get_i_struc
1888!
1889 FUNCTION yaml_get_l_struc (self, keystring, V) RESULT (status)
1890!
1891!***********************************************************************
1892! !
1893! It loads a vector set of logicals in YAML block-list structure, !
1894! V(1:Nitems)%vector(1:Nvalues). If the dummy argument V structure !
1895! is allocated, it deallocates/allocates the required Nitems and !
1896! Nvalues dimensions. !
1897! !
1898! On Input: !
1899! !
1900! self YAML tree dictionary object (CLASS yaml_tree) !
1901! keystring YAML tree aggregated keys (string) !
1902! !
1903! On Output: !
1904! !
1905! V Vector of logicals in block list (TYPE yaml_Lvec) !
1906! status Error code (integer) !
1907! !
1908!***********************************************************************
1909!
1910! Imported variable declarations.
1911!
1912 class(yaml_tree), intent(in) :: self
1913 character (len=*), intent(in) :: keystring
1914 TYPE (yaml_lvec), allocatable, intent(out) :: v(:)
1915!
1916! Local variable declarations.
1917!
1918 TYPE (yaml_extract), allocatable :: s(:)
1919!
1920 integer :: nitems, nvalues, n, i
1921 integer :: status
1922!
1923 character (len=*), parameter :: myfile = &
1924 & __FILE__//", yaml_Get_l_struc"
1925!
1926!-----------------------------------------------------------------------
1927! Extract requested key-string values.
1928!-----------------------------------------------------------------------
1929!
1930 status=noerr
1931!
1932 IF (yaml_error(self%extract(keystring, s), &
1933 & noerr, __line__, myfile)) RETURN
1934!
1935! Allocate output structure.
1936!
1937 nitems=SIZE(s, dim=1)
1938 IF (ALLOCATED(v)) DEALLOCATE (v)
1939 ALLOCATE ( v(nitems) )
1940!
1941! Convert string vector values to logicals.
1942!
1943 DO n=1,nitems
1944 IF (s(n)%has_vector) THEN
1945 nvalues=SIZE(s(n)%vector)
1946 ALLOCATE ( v(n)%vector(nvalues) )
1947 DO i=1,nvalues
1948 IF (yaml_uppercase(s(n)%vector(i)%value(1:1)).eq.'T') THEN
1949 v(n)%vector(i)=.true.
1950 ELSE
1951 v(n)%vector(i)=.false.
1952 END IF
1953 END DO
1954 END IF
1955 END DO
1956!
1957! Deallocate local extraction structure.
1958!
1959 IF (ALLOCATED(s)) DEALLOCATE (s)
1960!
1961 RETURN
1962 END FUNCTION yaml_get_l_struc
1963!
1964 FUNCTION yaml_get_r_struc (self, keystring, V) RESULT (status)
1965!
1966!***********************************************************************
1967! !
1968! It loads a vector set of real values in YAML block-list structure, !
1969! V(1:Nitems)%vector(1:Nvalues). If the dummy argument V structure !
1970! is allocated, it deallocates/allocates the required Nitems and !
1971! Nvalues dimensions. !
1972! !
1973! On Input: !
1974! !
1975! self YAML tree dictionary object (CLASS yaml_tree) !
1976! keystring YAML tree aggregated keys (string) !
1977! !
1978! On Output: !
1979! !
1980! V Vector of reals in block list (TYPE yaml_Rvec) !
1981! status Error code (integer) !
1982! !
1983!***********************************************************************
1984!
1985! Imported variable declarations.
1986!
1987 class(yaml_tree), intent(in) :: self
1988 character (len=*), intent(in) :: keystring
1989 TYPE (yaml_rvec), allocatable, intent(out) :: v(:)
1990!
1991! Local variable declarations.
1992!
1993 TYPE (yaml_extract), allocatable :: s(:)
1994!
1995 integer :: nitems, nvalues, i, n
1996 integer :: status
1997!
1998 character (len=Lmax) :: msg
1999
2000 character (len=*), parameter :: myfile = &
2001 & __FILE__//", yaml_Get_r_struc"
2002!
2003!-----------------------------------------------------------------------
2004! Extract requested key-string values.
2005!-----------------------------------------------------------------------
2006!
2007 status=noerr
2008!
2009 IF (yaml_error(self%extract(keystring, s), &
2010 & noerr, __line__, myfile)) RETURN
2011!
2012! Allocate output structure.
2013!
2014 nitems=SIZE(s, dim=1)
2015 IF (ALLOCATED(v)) DEALLOCATE (v)
2016 ALLOCATE ( v(nitems) )
2017!
2018! Convert string vector values to floating-point.
2019!
2020 DO n=1,nitems
2021 IF (s(n)%has_vector) THEN
2022 nvalues=SIZE(s(n)%vector)
2023 ALLOCATE ( v(n)%vector(nvalues) )
2024 DO i=1,nvalues
2025 READ (s(n)%vector(i)%value, * ,iostat=status, iomsg=msg) &
2026 & v(n)%vector(i)
2027 IF (yaml_error(status, noerr, __line__, myfile)) THEN
2028 yaml_errflag=5
2029 IF (yaml_master) WRITE (yaml_stdout,10) trim(keystring), &
2030 & s(n)%vector(i)%value, &
2031 & trim(msg)
2032 RETURN
2033 END IF
2034 END DO
2035 ELSE
2036 READ (s(n)%value, * ,iostat=status, iomsg=msg) &
2037 & v(n)%vector(i)
2038 IF (yaml_error(status, noerr, __line__, myfile)) THEN
2039 yaml_errflag=5
2040 IF (yaml_master) WRITE (yaml_stdout,10) trim(keystring), &
2041 & s(n)%value, &
2042 & trim(msg)
2043 RETURN
2044 END IF
2045 END IF
2046 END DO
2047!
2048! Deallocate local extraction structure.
2049!
2050 IF (ALLOCATED(s)) DEALLOCATE (s)
2051!
2052 10 FORMAT (/,' YAML_GET_R_STRUC - Error while converting string to', &
2053 & ' real, key: ',a,', value = ',a,/,20x,'ErrMsg: ',a)
2054!
2055 RETURN
2056 END FUNCTION yaml_get_r_struc
2057!
2058 FUNCTION yaml_get_s_struc (self, keystring, V) RESULT (status)
2059!
2060!***********************************************************************
2061! !
2062! It loads a vector set of strings in YAML block-list structure, !
2063! V(1:Nitems)%vector(1:Nvalues)%value. If the dummy argument V !
2064! structure is allocated, it deallocates/allocates the required !
2065! Nitems and Nvalues dimensions. !
2066! !
2067! On Input: !
2068! !
2069! self YAML tree dictionary object (CLASS yaml_tree) !
2070! keystring YAML tree aggregated keys (string) !
2071! !
2072! On Output: !
2073! !
2074! V Vector of strings in block list (TYPE yaml_Svec) !
2075! status Error code (integer) !
2076! !
2077!***********************************************************************
2078!
2079! Imported variable declarations.
2080!
2081 class(yaml_tree), intent(in) :: self
2082 character (len=*), intent(in) :: keystring
2083 TYPE (yaml_svec), allocatable, intent(out) :: v(:)
2084!
2085! Local variable declarations.
2086!
2087 TYPE (yaml_extract), allocatable :: s(:)
2088!
2089 integer :: nitems, nvalues, n, i
2090 integer :: status
2091!
2092 character (len=*), parameter :: myfile = &
2093 & __FILE__//", yaml_Get_s_struc"
2094!
2095!-----------------------------------------------------------------------
2096! Extract requested key-string values.
2097!-----------------------------------------------------------------------
2098!
2099 status=noerr
2100!
2101 IF (yaml_error(self%extract(keystring, s), &
2102 & noerr, __line__, myfile)) RETURN
2103!
2104! Allocate output structure.
2105!
2106 nitems=SIZE(s, dim=1)
2107 IF (ALLOCATED(v)) DEALLOCATE (v)
2108 ALLOCATE ( v(nitems) )
2109!
2110! Load string vector values to output structure.
2111!
2112 DO n=1,nitems
2113 IF (s(n)%has_vector) THEN
2114 nvalues=SIZE(s(n)%vector)
2115 ALLOCATE ( v(n)%vector(nvalues) )
2116 DO i=1,nvalues
2117 v(n)%vector(i)%value=s(n)%vector(i)%value
2118 END DO
2119 ELSE
2120 v(n)%value=s(n)%value
2121 END IF
2122 END DO
2123!
2124! Deallocate local extraction structure.
2125!
2126 IF (ALLOCATED(s)) DEALLOCATE (s)
2127!
2128 RETURN
2129 END FUNCTION yaml_get_s_struc
2130!
2131 FUNCTION yaml_get_ivar_0d (self, keystring, value) RESULT (status)
2132!
2133!***********************************************************************
2134! !
2135! It gets scalar integer data from YAML dictionary object. !
2136! !
2137! On Input: !
2138! !
2139! self YAML tree dictionary object (CLASS yaml_tree) !
2140! keystring YAML tree aggregated keys (string) !
2141! !
2142! On Output: !
2143! !
2144! value YAML value (integer; scalar) !
2145! status Error code (integer) !
2146! !
2147!***********************************************************************
2148!
2149! Imported variable declarations.
2150!
2151 class(yaml_tree), intent(in) :: self
2152 character (len=*), intent(in) :: keystring
2153 integer, intent(out) :: value
2154!
2155! Local variable declarations.
2156!
2157 TYPE (yaml_extract), allocatable :: s(:)
2158!
2159 integer :: nvalues, status
2160!
2161 character (len=Lmax) :: msg
2162
2163 character (len=*), parameter :: myfile = &
2164 & __FILE__//", yaml_Get_ivar_0d"
2165!
2166!-----------------------------------------------------------------------
2167! Extract requested key-string value.
2168!-----------------------------------------------------------------------
2169!
2170 status=noerr
2171!
2172 IF (yaml_error(self%extract(keystring, s), &
2173 & noerr, __line__, myfile)) RETURN
2174!
2175! Make sure that extracted value is a scalar.
2176!
2177 nvalues=SIZE(s)
2178!
2179 IF (nvalues.gt.1) THEN
2180 status=7
2181 yaml_errflag=status
2182 IF (yaml_error(status, noerr, __line__, myfile)) THEN
2183 IF (yaml_master) WRITE (yaml_stdout,10) keystring, nvalues, &
2184 & self%filename
2185 RETURN
2186 END IF
2187!
2188! Convert string value to integer.
2189!
2190 ELSE
2191 READ (s(1)%value, *, iostat=status, iomsg=msg) value
2192 IF (yaml_error(status, noerr, __line__, myfile)) THEN
2193 yaml_errflag=5
2194 IF (yaml_master) WRITE (yaml_stdout,20) trim(keystring), &
2195 & s(1)%value, trim(msg)
2196 RETURN
2197 END IF
2198 END IF
2199!
2200! Deallocate.
2201!
2202 IF (ALLOCATED(s)) DEALLOCATE (s)
2203!
2204 10 FORMAT (/,' YAML_GET_IVAR_0D - Extracted value is not a scalar,', &
2205 & ' key-string: ',a,/,20x,'Nvalues = ',i0,/,20x,'File: ',a)
2206 20 FORMAT (/,' YAML_GET_IVAR_0D - Error while converting string to', &
2207 & ' integer, key: ',a,', value = ',a,/,20x,'ErrMsg: ',a)
2208!
2209 RETURN
2210 END FUNCTION yaml_get_ivar_0d
2211!
2212 FUNCTION yaml_get_ivar_1d (self, keystring, value) RESULT (status)
2213!
2214!***********************************************************************
2215! !
2216! It gets 1D integer data from YAML dictionary object. !
2217! !
2218! On Input: !
2219! !
2220! self YAML tree dictionary object (CLASS yaml_tree) !
2221! keystring YAML tree aggregated keys (string) !
2222! !
2223! On Output: !
2224! !
2225! value YAML value (integer; 1D-array) !
2226! status Error code (integer) !
2227! !
2228!***********************************************************************
2229!
2230! Imported variable declarations.
2231!
2232 class(yaml_tree), intent(in) :: self
2233 character (len=*), intent(in) :: keystring
2234 integer, intent(out) :: value(:)
2235!
2236! Local variable declarations.
2237!
2238 TYPE (yaml_extract), allocatable :: s(:)
2239!
2240 integer :: nvalues, i, status
2241!
2242 character (len=Lmax) :: msg
2243
2244 character (len=*), parameter :: myfile = &
2245 & __FILE__//", yaml_Get_ivar_1d"
2246!
2247!-----------------------------------------------------------------------
2248! Extract requested key-string values.
2249!-----------------------------------------------------------------------
2250!
2251 status=noerr
2252!
2253 IF (yaml_error(self%extract(keystring, s), &
2254 & noerr, __line__, myfile)) RETURN
2255!
2256! Make sure that extracted value is a scalar.
2257!
2258 nvalues=SIZE(s, dim=1)
2259!
2260 IF (SIZE(value, dim=1).lt.nvalues) THEN
2261 status=7
2262 yaml_errflag=status
2263 IF (yaml_error(status, noerr, __line__, myfile)) THEN
2264 IF (yaml_master) WRITE (yaml_stdout,10) keystring, nvalues, &
2265 & SIZE(value, dim=1), &
2266 & self%filename
2267 RETURN
2268 END IF
2269 END IF
2270!
2271! Convert string values to integers.
2272!
2273 DO i=1,nvalues
2274 READ (s(i)%value, *, iostat=status, iomsg=msg) value(i)
2275 IF (yaml_error(status, noerr, __line__, myfile)) THEN
2276 yaml_errflag=5
2277 IF (yaml_master) WRITE (yaml_stdout,20) trim(keystring), &
2278 & s(i)%value, trim(msg)
2279 RETURN
2280 END IF
2281 END DO
2282!
2283! Deallocate.
2284!
2285 IF (ALLOCATED(s)) DEALLOCATE (s)
2286!
2287 10 FORMAT (/,' YAML_GET_IVAR_1D - Inconsistent number of values,', &
2288 & ' key-string: ',a,/,20x,'YAML size = ',i0, &
2289 & ', Variable size = ',i0,/,20x,'File: ',a)
2290 20 FORMAT (/,' YAML_GET_IVAR_1D - Error while converting string to', &
2291 & ' integer, key: ',a,', value = ',a,/,20x,'ErrMsg: ',a)
2292!
2293 RETURN
2294 END FUNCTION yaml_get_ivar_1d
2295!
2296 FUNCTION yaml_get_lvar_0d (self, keystring, value) RESULT (status)
2297!
2298!***********************************************************************
2299! !
2300! It gets scalar logical data from YAML disctionary object. !
2301! !
2302! On Input: !
2303! !
2304! self YAML tree dictionary object (CLASS yaml_tree) !
2305! keystring YAML tree aggregated keys (string) !
2306! !
2307! On Output: !
2308! !
2309! value YAML value (logical; scalar) !
2310! !
2311!***********************************************************************
2312!
2313! Imported variable declarations.
2314!
2315 class(yaml_tree), intent(in) :: self
2316 character (len=*), intent(in) :: keystring
2317 logical, intent(out) :: value
2318!
2319! Local variable declarations.
2320!
2321 TYPE (yaml_extract), allocatable :: s(:)
2322!
2323 integer :: nvalues, status
2324!
2325 character (len=*), parameter :: myfile = &
2326 & __FILE__//", yaml_Get_lvar_0d"
2327!
2328!-----------------------------------------------------------------------
2329! Extract requested key-string value.
2330!-----------------------------------------------------------------------
2331!
2332 status=noerr
2333!
2334 IF (yaml_error(self%extract(keystring, s), &
2335 & noerr, __line__, myfile)) RETURN
2336!
2337! Make sure that extracted value is a scalar.
2338!
2339 nvalues=SIZE(s)
2340!
2341 IF (nvalues.gt.1) THEN
2342 status=7
2343 yaml_errflag=status
2344 IF (yaml_error(status, noerr, __line__, myfile)) THEN
2345 IF (yaml_master) WRITE (yaml_stdout,10) keystring, nvalues, &
2346 & self%filename
2347 RETURN
2348 END IF
2349!
2350! Convert string value to logical.
2351!
2352 ELSE
2353 IF (yaml_uppercase(s(1)%value(1:1)).eq.'T') THEN
2354 value=.true.
2355 ELSE
2356 value=.false.
2357 END IF
2358 END IF
2359!
2360! Deallocate.
2361!
2362 IF (ALLOCATED(s)) DEALLOCATE (s)
2363!
2364 10 FORMAT (/,' YAML_GET_LVAR_0D - Extracted value is not a scalar,', &
2365 & ' key-string: ',a,/,20x,'Nvalues = ',i0,/,20x,'File: ',a)
2366!
2367 RETURN
2368 END FUNCTION yaml_get_lvar_0d
2369!
2370 FUNCTION yaml_get_lvar_1d (self, keystring, value) RESULT (status)
2371!
2372!***********************************************************************
2373! !
2374! It gets 1D logical data from YAML dictionary object. !
2375! !
2376! On Input: !
2377! !
2378! self YAML tree dictionary object (CLASS yaml_tree) !
2379! keystring YAML tree aggregated keys (string) !
2380! !
2381! On Output: !
2382! !
2383! value YAML value (logical; 1D-array) !
2384! status Error code (integer) !
2385! !
2386!***********************************************************************
2387!
2388! Imported variable declarations.
2389!
2390 class(yaml_tree), intent(in) :: self
2391 character (len=*), intent(in) :: keystring
2392 logical, intent(out) :: value(:)
2393!
2394! Local variable declarations.
2395!
2396 TYPE (yaml_extract), allocatable :: s(:)
2397!
2398 integer :: nvalues, i, status
2399!
2400 character (len=*), parameter :: myfile = &
2401 & __FILE__//", yaml_Get_lvar_1d"
2402!
2403!-----------------------------------------------------------------------
2404! Extract requested key-string values.
2405!-----------------------------------------------------------------------
2406!
2407 status=noerr
2408!
2409 IF (yaml_error(self%extract(keystring, s), &
2410 & noerr, __line__, myfile)) RETURN
2411!
2412! Make sure that extracted value is a scalar.
2413!
2414 nvalues=SIZE(s, dim=1)
2415!
2416 IF (SIZE(value, dim=1).lt.nvalues) THEN
2417 yaml_errflag=7
2418 status=yaml_errflag
2419 IF (yaml_error(yaml_errflag, noerr, __line__, myfile)) THEN
2420 IF (yaml_master) WRITE (yaml_stdout,10) keystring, nvalues, &
2421 & SIZE(value, dim=1), &
2422 & self%filename
2423 RETURN
2424 END IF
2425 END IF
2426!
2427! Convert string values to logicals.
2428!
2429 DO i=1,nvalues
2430 IF (yaml_uppercase(s(i)%value(1:1)).eq.'T') THEN
2431 value(i)=.true.
2432 ELSE
2433 value(i)=.false.
2434 END IF
2435 END DO
2436!
2437! Deallocate.
2438!
2439 IF (ALLOCATED(s)) DEALLOCATE (s)
2440!
2441 10 FORMAT (/,' YAML_GET_LVAR_1D - Inconsistent number of values,', &
2442 & ' key-string: ',a,/,20x,'YAML size = ',i0, &
2443 & ', Variable size = ',i0,/,20x,'File: ',a)
2444!
2445 RETURN
2446 END FUNCTION yaml_get_lvar_1d
2447!
2448 FUNCTION yaml_get_rvar_0d (self, keystring, value) RESULT (status)
2449!
2450!***********************************************************************
2451! !
2452! It gets scalar floating-point data from YAML dictionary object. !
2453! !
2454! On Input: !
2455! !
2456! self YAML tree dictionary object (CLASS yaml_tree) !
2457! keystring YAML tree aggregated keys (string) !
2458! !
2459! On Output: !
2460! !
2461! value YAML value (real; scalar) !
2462! status Error code (integer) !
2463! !
2464!***********************************************************************
2465!
2466! Imported variable declarations.
2467!
2468 class(yaml_tree), intent(in) :: self
2469 character (len=*), intent(in) :: keystring
2470 real (kind_real), intent(out) :: value
2471!
2472! Local variable declarations.
2473!
2474 TYPE (yaml_extract), allocatable :: s(:)
2475!
2476 integer :: nvalues, ie, status
2477!
2478 character (len=Lmax) :: msg
2479
2480 character (len=*), parameter :: myfile = &
2481 & __FILE__//", yaml_Get_rvar_0d"
2482!
2483!-----------------------------------------------------------------------
2484! Extract requested key-string value.
2485!-----------------------------------------------------------------------
2486!
2487 status=noerr
2488!
2489 IF (yaml_error(self%extract(keystring, s), &
2490 & noerr, __line__, myfile)) RETURN
2491!
2492! Make sure that extracted value is a scalar.
2493!
2494 nvalues=SIZE(s)
2495!
2496 IF (nvalues.gt.1) THEN
2497 yaml_errflag=7
2498 status=yaml_errflag
2499 IF (yaml_error(yaml_errflag, noerr, __line__, myfile)) THEN
2500 IF (yaml_master) WRITE (yaml_stdout,10) keystring, nvalues, &
2501 & self%filename
2502 RETURN
2503 END IF
2504!
2505! Convert string value to real.
2506!
2507 ELSE
2508 s(1)%value=adjustl(s(1)%value)
2509 ie=len_trim(s(1)%value)
2510 READ (s(1)%value(1:ie), *, iostat=status, iomsg=msg) value
2511 IF (yaml_error(status, noerr, __line__, myfile)) THEN
2512 yaml_errflag=5
2513 IF (yaml_master) WRITE (yaml_stdout,20) trim(keystring), &
2514 & s(1)%value, trim(msg)
2515 RETURN
2516 END IF
2517 END IF
2518!
2519! Deallocate.
2520!
2521 IF (ALLOCATED(s)) DEALLOCATE (s)
2522!
2523 10 FORMAT (/,' YAML_GET_RVAR_0D - Extracted value is not a scalar,', &
2524 & ' key-string: ',a,/,20x,'Nvalues = ',i0,/,20x,'File: ',a)
2525 20 FORMAT (/,' YAML_GET_RVAR_0D - Error while converting string to', &
2526 & ' integer, key: ',a,', value = ',a,/,20x,'ErrMsg: ',a)
2527!
2528 RETURN
2529 END FUNCTION yaml_get_rvar_0d
2530!
2531 FUNCTION yaml_get_rvar_1d (self, keystring, value) RESULT (status)
2532!
2533!***********************************************************************
2534! !
2535! It gets 1D floating-point data from YAML dictionary object. !
2536! !
2537! On Input: !
2538! !
2539! self YAML tree dictionary object (CLASS yaml_tree) !
2540! keystring YAML tree aggregated keys (string) !
2541! !
2542! On Output: !
2543! !
2544! value YAML value (real; 1D-array) !
2545! status Error code (integer) !
2546! !
2547!***********************************************************************
2548!
2549! Imported variable declarations.
2550!
2551 class(yaml_tree), intent(in) :: self
2552 character (len=*), intent(in) :: keystring
2553 real (kind_real), intent(out) :: value(:)
2554!
2555! Local variable declarations.
2556!
2557 TYPE (yaml_extract), allocatable :: s(:)
2558!
2559 integer :: nvalues, i, ie, status
2560!
2561 character (len=Lmax) :: msg
2562
2563 character (len=*), parameter :: myfile = &
2564 & __FILE__//", yaml_Get_rvar_1d"
2565!
2566!-----------------------------------------------------------------------
2567! Extract requested key-string values.
2568!-----------------------------------------------------------------------
2569!
2570 status=noerr
2571!
2572 IF (yaml_error(self%extract(keystring, s), &
2573 & noerr, __line__, myfile)) RETURN
2574!
2575! Make sure that extracted value is a scalar.
2576!
2577 nvalues=SIZE(s, dim=1)
2578!
2579 IF (SIZE(value, dim=1).lt.nvalues) THEN
2580 yaml_errflag=7
2581 status=yaml_errflag
2582 IF (yaml_error(yaml_errflag, noerr, __line__, myfile)) THEN
2583 IF (yaml_master) WRITE (yaml_stdout,10) keystring, nvalues, &
2584 & SIZE(value, dim=1), &
2585 & self%filename
2586 RETURN
2587 END IF
2588 END IF
2589!
2590! Convert string values to reals.
2591!
2592 DO i=1,nvalues
2593 s(i)%value=adjustl(s(i)%value)
2594 ie=len_trim(s(i)%value)
2595 READ (s(i)%value(1:ie), *, iostat=status, iomsg=msg) value(i)
2596 IF (yaml_error(status, noerr, __line__, myfile)) THEN
2597 yaml_errflag=5
2598 IF (yaml_master) WRITE (yaml_stdout,20) trim(keystring), &
2599 & s(i)%value, trim(msg)
2600 RETURN
2601 END IF
2602 END DO
2603!
2604! Deallocate.
2605!
2606 IF (ALLOCATED(s)) DEALLOCATE (s)
2607!
2608 10 FORMAT (/,' YAML_GET_RVAR_1D - Inconsistent number of values,', &
2609 & ' key-string: ',a,/,20x,'YAML size = ',i0, &
2610 & ', Variable size = ',i0,/,20x,'File: ',a)
2611 20 FORMAT (/,' YAML_GET_RVAR_1D - Error while converting string to', &
2612 & ' real, key: ',a,', value = ',a,/,20x,'ErrMsg: ',a)
2613!
2614 RETURN
2615 END FUNCTION yaml_get_rvar_1d
2616!
2617 FUNCTION yaml_get_svar_0d (self, keystring, value) RESULT (status)
2618!
2619!***********************************************************************
2620! !
2621! It gets scalar string data from YAML dictionary object. !
2622! !
2623! On Input: !
2624! !
2625! self YAML tree dictionary object (CLASS yaml_tree) !
2626! keystring YAML tree aggregated keys (string) !
2627! !
2628! On Output: !
2629! !
2630! value YAML value (string; scalar) !
2631! status Error code (integer) !
2632! !
2633!***********************************************************************
2634!
2635! Imported variable declarations.
2636!
2637 class(yaml_tree), intent(in ) :: self
2638 character (len=*), intent(in ) :: keystring
2639 character (len=*), intent(out) :: value
2640!
2641! Local variable declarations.
2642!
2643 TYPE (yaml_extract), allocatable :: s(:)
2644!
2645 integer :: nvalues, status
2646!
2647 character (len=*), parameter :: myfile = &
2648 & __FILE__//", yaml_Get_lvar_0d"
2649!
2650!-----------------------------------------------------------------------
2651! Extract requested key-string value.
2652!-----------------------------------------------------------------------
2653!
2654 status=noerr
2655!
2656 IF (yaml_error(self%extract(keystring, s), &
2657 & noerr, __line__, myfile)) RETURN
2658!
2659! Make sure that extracted value is a scalar.
2660!
2661 nvalues=SIZE(s)
2662!
2663 IF (nvalues.gt.1) THEN
2664 yaml_errflag=7
2665 status=yaml_errflag
2666 IF (yaml_error(yaml_errflag, noerr, __line__, myfile)) THEN
2667 IF (yaml_master) WRITE (yaml_stdout,10) keystring, nvalues, &
2668 & self%filename
2669 RETURN
2670 END IF
2671!
2672! Load string value.
2673!
2674 ELSE
2675 value=s(1)%value
2676 END IF
2677!
2678! Deallocate.
2679!
2680 IF (ALLOCATED(s)) DEALLOCATE (s)
2681!
2682 10 FORMAT (/,' YAML_GET_SVAR_0D - Extracted value is not a scalar,', &
2683 & ' key-string: ',a,/,20x,'Nvalues = ',i0,/,20x,'File: ',a)
2684!
2685 RETURN
2686 END FUNCTION yaml_get_svar_0d
2687!
2688 FUNCTION yaml_get_svar_1d (self, keystring, value) RESULT (status)
2689!
2690!***********************************************************************
2691! !
2692! It gets 1D string data from YAML dictionary object. !
2693! !
2694! On Input: !
2695! !
2696! self YAML tree dictionary object (CLASS yaml_tree) !
2697! keystring YAML tree aggregated keys (string) !
2698! !
2699! On Output: !
2700! !
2701! value YAML value (string; 1D-array) !
2702! status Error code (integer) !
2703! !
2704!***********************************************************************
2705!
2706! Imported variable declarations.
2707!
2708 class(yaml_tree), intent(in ) :: self
2709 character (len=*), intent(in ) :: keystring
2710 character (len=*), intent(out) :: value(:)
2711!
2712! Local variable declarations.
2713!
2714 TYPE (yaml_extract), allocatable :: s(:)
2715!
2716 integer :: nvalues, i, status
2717!
2718 character (len=*), parameter :: myfile = &
2719 & __FILE__//", yaml_Get_lvar_1d"
2720!
2721!-----------------------------------------------------------------------
2722! Extract requested key-string values.
2723!-----------------------------------------------------------------------
2724!
2725 status=noerr
2726!
2727 IF (yaml_error(self%extract(keystring, s), &
2728 & noerr, __line__, myfile)) RETURN
2729!
2730! Make sure that extracted value is a scalar.
2731!
2732 nvalues=SIZE(s, dim=1)
2733!
2734 IF (SIZE(value, dim=1).lt.nvalues) THEN
2735 yaml_errflag=7
2736 status=yaml_errflag
2737 IF (yaml_error(yaml_errflag, noerr, __line__, myfile)) THEN
2738 IF (yaml_master) WRITE (yaml_stdout,10) keystring, nvalues, &
2739 & SIZE(value, dim=1), &
2740 & self%filename
2741 RETURN
2742 END IF
2743 END IF
2744!
2745! Load string values.
2746!
2747 DO i=1,nvalues
2748 value(i)=s(i)%value
2749 END DO
2750!
2751! Deallocate.
2752!
2753 IF (ALLOCATED(s)) DEALLOCATE (s)
2754!
2755 10 FORMAT (/,' YAML_GET_SVAR_1D - Inconsistent number of values,', &
2756 & ' key-string: ',a,/,20x,'YAML size = ',i0, &
2757 & ', Variable size = ',i0,/,20x,'File: ',a)
2758!
2759 RETURN
2760 END FUNCTION yaml_get_svar_1d
2761!
2762 FUNCTION yaml_lowercase (Sinp) RESULT (Sout)
2763!
2764!=======================================================================
2765! !
2766! It converts all string elements to lowercase. !
2767! !
2768! Reference: !
2769! !
2770! Cooper Redwine, 1995: "Upgrading to Fortran 90", Springer- !
2771! Verlag, New York, pp 416. !
2772! !
2773!=======================================================================
2774!
2775! Imported variable declarations.
2776!
2777 character(len=*), intent(in) :: sinp
2778!
2779! Local variable definitions.
2780!
2781 integer :: lstr, i, j
2782
2783 character (len=LEN(Sinp)) :: sout
2784
2785 character (len=26), parameter :: l = 'abcdefghijklmnopqrstuvwxyz'
2786 character (len=26), parameter :: u = 'ABCDEFGHIJKLMNOPQRSTUVWXYZ'
2787!
2788!-----------------------------------------------------------------------
2789! Convert input string to lowercase.
2790!-----------------------------------------------------------------------
2791!
2792 lstr=len(sinp)
2793 sout=sinp
2794 DO i=1,lstr
2795 j=index(u, sout(i:i))
2796 IF (j.ne.0) THEN
2797 sout(i:i)=l(j:j)
2798 END IF
2799 END DO
2800!
2801 RETURN
2802 END FUNCTION yaml_lowercase
2803!
2804 FUNCTION yaml_uppercase (Sinp) RESULT (Sout)
2805!
2806!=======================================================================
2807! !
2808! It converts all string elements to uppercase. !
2809! !
2810! Reference: !
2811! !
2812! Cooper Redwine, 1995: "Upgrading to Fortran 90", Springer- !
2813! Verlag, New York, pp 416. !
2814! !
2815!=======================================================================
2816!
2817! Imported variable declarations.
2818!
2819 character (len=*), intent(in) :: sinp
2820!
2821! Local variable definitions.
2822!
2823 integer :: lstr, i, j
2824!
2825 character (len=LEN(Sinp)) :: sout
2826!
2827 character (len=26), parameter :: l = 'abcdefghijklmnopqrstuvwxyz'
2828 character (len=26), parameter :: u = 'ABCDEFGHIJKLMNOPQRSTUVWXYZ'
2829!
2830!-----------------------------------------------------------------------
2831! Convert input string to uppercase.
2832!-----------------------------------------------------------------------
2833!
2834 lstr=len(sinp)
2835 sout=sinp
2836 DO i=1,lstr
2837 j=index(l, sout(i:i))
2838 IF (j.ne.0) THEN
2839 sout(i:i)=u(j:j)
2840 END IF
2841 END DO
2842!
2843 RETURN
2844 END FUNCTION yaml_uppercase
2845!
2846 SUBROUTINE yaml_valuetype (value, Lswitch)
2847!
2848!***********************************************************************
2849! !
2850! It determines the YAML value type as logical, integer, real, or !
2851! string. !
2852! !
2853! On Input: !
2854! !
2855! self YAML pair value (string) !
2856! !
2857! On Output: !
2858! !
2859! Lswitch YAML key/value switches (logical, vector) !
2860! !
2861! Lswitch(1) = T/F, value has an alias '*' token !
2862! Lswitch(2) = T/F, value has an anchor '&' token !
2863! Lswitch(3) = T/F, key/value starts a block '-' !
2864! Lswitch(4) = T/F, value has sequence '[...]' tokens !
2865! Lswitch(5) = T/F, logical value(s) !
2866! Lswitch(6) = T/F, integer value(s) !
2867! Lswitch(7) = T/F, floating-point value(s) !
2868! Lswitch(8) = T/F, string value(s) !
2869! !
2870!***********************************************************************
2871!
2872! Imported variable declarations.
2873!
2874 logical, intent(inout) :: lswitch(:)
2875!
2876 character (len=*), intent(inout) :: value
2877!
2878! Local variable declarations.
2879!
2880 integer :: lstr, schar
2881 integer :: colon, dot, exponent, letter, numeric, precision
2882 integer :: i, multiple
2883!
2884!-----------------------------------------------------------------------
2885! Set keyword value type.
2886!-----------------------------------------------------------------------
2887!
2888! Initialize.
2889!
2890 colon=0
2891 dot=0
2892 exponent=0
2893 letter=0
2894 multiple=0
2895 numeric=0
2896 precision=0
2897!
2898! Check input value string.
2899!
2900 lstr=len_trim(value)
2901 IF (lstr.gt.0) THEN
2902 DO i=1,lstr
2903 schar=ichar(value(i:i))
2904!
2905! Check for numbers, plus, and minus signs. For example, value=-1.0E+37
2906! is a floating-point numerical value.
2907!
2908 IF (((48.le.schar).and.(schar.le.57)).or. &
2909 & (schar.eq.43).or.(schar.eq.45)) numeric=numeric+1
2910!
2911! Check for dot/period, CHAR(46). If period is not found in a numerical
2912! expression, identify value as an integer.
2913!
2914 IF (schar.eq.46) dot=dot+1
2915!
2916! Check for precision character: D=ICHAR(68) and d=ICHAR(100). For
2917! example, value=0.0d0 and others in the future.
2918!
2919 IF ((schar.eq.69).or.(schar.eq.101)) precision=precision+1
2920!
2921! Check for exponent character: E=CHAR(69) and e=CHAR(101).
2922!
2923 IF ((schar.eq.69).or.(schar.eq.101)) exponent=exponent+1
2924!
2925! Check for multiple values separate by comma, CHAR(44), in a sequence
2926! of values: [val1, val2, ..., valN].
2927!
2928 IF (lswitch(4)) multiple=multiple+1
2929!
2930! Check for colon, CHAR(58). We can have value=2020-01-01T00:00:00Z,
2931! which has numbers, hyphen, and letters. It is the colon that makes
2932! it a string variable (https://www.myroms.org).
2933!
2934 IF (schar.eq.58) colon=colon+1
2935!
2936! English uppercase and lowercase alphabet.
2937!
2938 IF (((65.le.schar).and.(schar.le.90)).or. &
2939 & (schar.eq.97).or.(schar.eq.122)) letter=letter+1
2940 END DO
2941!
2942! Set integer or floating-point type.
2943!
2944 IF ((numeric.gt.0).and.(colon.eq.0)) THEN
2945 IF ((dot.gt.0).or.(exponent.gt.0).or.(precision.gt.0)) THEN
2946 lswitch(7)=.true. ! floating-point
2947 ELSE
2948 lswitch(6)=.true. ! integer
2949 END IF
2950 ELSE
2951 IF ((value.eq.'true').or.(value.eq.'false').or. &
2952 & (value.eq.'TRUE').or.(value.eq.'FALSE')) THEN
2953 lswitch(5)=.true. ! logical
2954 ELSE IF (letter.gt.0) THEN
2955 lswitch(8)=.true. ! string
2956 END IF
2957 END IF
2958 END IF
2959!
2960 RETURN
2961 END SUBROUTINE yaml_valuetype
2962!
2963! <><><><><><><><><><><><><><><><><><><><><><><><><><><><><><><><><><><>
2964!
2965 SUBROUTINE process_yaml (yaml_file)
2966!
2967! Imported variable declarations.
2968!
2969 character (len=*), intent(in) :: yaml_file
2970!
2971! Local variable declarations.
2972!
2973 character (len=*), parameter :: myfile = &
2974 & __FILE__//', process_yaml'
2975
2976 TYPE (yaml_tree) :: yml ! YAML tree structure
2977
2978 TYPE (couplingfield), allocatable :: e(:), i(:), s(:)
2979!
2980! Create YAML object and report its content.
2981!
2982 IF (yaml_error(yaml_initialize(yml, trim(yaml_file), lreport), &
2983 & noerror, __line__, myfile)) THEN
2984 WRITE (6,10) trim(yaml_file)
2985 10 FORMAT (/,' PROCESS_YAML - Unable to create YAML object for:', &
2986 & /,13x,'FileName = ',a)
2987 RETURN
2988 END IF
2989!
2990! If CMEPS file, extract export and import states metadata.
2991!
2992 IF (yml%has('export').and.yml%has('import')) THEN
2993 CALL cmeps_metadata (yml, trim(yaml_file), 'export', e)
2994 CALL cmeps_metadata (yml, trim(yaml_file), 'import', i)
2995 CALL yml%destroy ()
2996!
2997! If ROMS coupling file, extract export nas import state metadata.
2998!
2999 ELSE IF (.not.yml%has('svn_repository').and. &
3000 & yml%has('metadata')) THEN
3001 CALL coupling_metadata (yml, trim(yaml_file), s)
3002!
3003! IF ROMS 'varingo.yaml', extract input and output variables metadata.
3004!
3005 ELSE IF (yml%has('svn_repository').and. &
3006 & yml%has('metadata')) THEN
3007 CALL io_metadata (yml, trim(yaml_file))
3008 END IF
3009!
3010 RETURN
3011!
3012 END SUBROUTINE process_yaml
3013!
3014 SUBROUTINE cmeps_metadata (self, filename, key, S)
3015!
3016!=======================================================================
3017! !
3018! It process either import or export fields which are stored as block !
3019! lists (leading key/value is hyphenated) in the YAML file. The YAML !
3020! file is used to configure ROMS ESMF/NUOPC 'cap' module to be run by !
3021! the Community Mediator for Earth Prediction Systems (CMEPS). !
3022! !
3023! On Input: !
3024! !
3025! self YAML tree dictionary (TYPE yaml_tree) !
3026! !
3027! filename ROMS YAML configuration filename for CMEPS (string) !
3028! !
3029! key Leading blocking key to process (string), for example: !
3030! 'export', 'import', or 'bulk_flux import' !
3031! !
3032! On Output: !
3033! !
3034! S Import or Export coupling fields (TYPE CouplingField) !
3035! !
3036!=======================================================================
3037!
3038! Imported variable declarations.
3039!
3040 TYPE (yaml_tree), intent(inout) :: self
3041 character (len=*), intent(in ) :: filename
3042 character (len=*), intent(in ) :: key
3043 TYPE (couplingfield), allocatable, intent(out) :: s(:)
3044!
3045! Local variable declarations.
3046!
3047 integer :: i, lenstr
3048!
3049 character (len=*), parameter :: myfile = &
3050 & __FILE__//", cmeps_metadata"
3051!
3052!-----------------------------------------------------------------------
3053! Process coupling import or export metadata for CMEPS.
3054!-----------------------------------------------------------------------
3055!
3056! If applicable, create YAML tree dictionary.
3057!
3058 IF (.not.ASSOCIATED(self%list)) THEN
3059 IF (yaml_error(yaml_initialize(self, trim(filename), &
3060 & .false.), &
3061 & noerror, __line__, myfile)) THEN
3062 IF (yaml_master) WRITE (stdout,10) trim(filename)
3063 RETURN
3064 END IF
3065 END IF
3066!
3067! Extract requested blocking list.
3068!
3069 IF (yaml_error(yaml_get(self, trim(key)//'.standard_name', &
3070 & ystring1), &
3071 & noerror, __line__, myfile)) RETURN
3072 nentries=SIZE(ystring1)
3073!
3074 IF (yaml_error(yaml_get(self, trim(key)//'.long_name', &
3075 & ystring2), &
3076 & noerror, __line__, myfile)) RETURN
3077!
3078 IF (yaml_error(yaml_get(self, trim(key)//'.short_name', &
3079 & ystring3), &
3080 & noerror, __line__, myfile)) RETURN
3081!
3082 IF (yaml_error(yaml_get(self, trim(key)//'.data_variables', &
3083 & ystring4), &
3084 & noerror, __line__, myfile)) RETURN
3085!
3086 IF (yaml_error(yaml_get(self, trim(key)//'.source_units', &
3087 & ystring5), &
3088 & noerror, __line__, myfile)) RETURN
3089!
3090 IF (yaml_error(yaml_get(self, trim(key)//'.destination_units', &
3091 & ystring6), &
3092 & noerror, __line__, myfile)) RETURN
3093!
3094 IF (yaml_error(yaml_get(self, trim(key)//'.source_grid', &
3095 & ystring7), &
3096 & noerror, __line__, myfile)) RETURN
3097!
3098 IF (yaml_error(yaml_get(self, trim(key)//'.destination_grid', &
3099 & ystring8), &
3100 & noerror, __line__, myfile)) RETURN
3101!
3102 IF (yaml_error(yaml_get(self, trim(key)//'.connected_to', &
3103 & ystring9), &
3104 & noerror, __line__, myfile)) RETURN
3105!
3106 IF (yaml_error(yaml_get(self, trim(key)//'.map_type', &
3107 & ystring10), &
3108 & noerror, __line__, myfile)) RETURN
3109!
3110 IF (yaml_error(yaml_get(self, trim(key)//'.map_norm', &
3111 & ystring11), &
3112 & noerror, __line__, myfile)) RETURN
3113!
3114 IF (.not.allocated(yreal1)) THEN
3115 allocate ( yreal1(nentries) )
3116 END IF
3117 IF (yaml_error(yaml_get(self, trim(key)//'.add_offset', &
3118 & yreal1), &
3119 & noerror, __line__, myfile)) RETURN
3120!
3121 IF (.not.allocated(yreal2)) THEN
3122 allocate ( yreal2(nentries) )
3123 END IF
3124 IF (yaml_error(yaml_get(self, trim(key)//'.scale', &
3125 & yreal2), &
3126 & noerror, __line__, myfile)) RETURN
3127 IF (.not.allocated(yreal1)) THEN
3128 allocate ( yreal1(nentries) )
3129 END IF
3130!
3131 IF (.not.allocated(ylogical1)) THEN
3132 allocate ( ylogical1(nentries) )
3133 END IF
3134 IF (yaml_error(yaml_get(self, trim(key)//'.debug_write', &
3135 & ylogical1), &
3136 & noerror, __line__, myfile)) RETURN
3137!
3138! Load metadata into output structure.
3139!
3140 IF (.not.allocated(s)) allocate ( s(nentries) )
3141!
3142 DO i=1,nentries
3143 s(i)%debug_write = ylogical1(i)
3144 s(i)%add_offset = yreal1(i)
3145 s(i)%scale = yreal2(i)
3146!
3147 IF (yaml_error(yaml_assignstring(s(i)%standard_name, &
3148 & ystring1(i)%value, lenstr), &
3149 & noerror, __line__, myfile)) RETURN
3150!
3151 IF (yaml_error(yaml_assignstring(s(i)%long_name, &
3152 & ystring2(i)%value, lenstr), &
3153 & noerror, __line__, myfile)) RETURN
3154!
3155 IF (yaml_error(yaml_assignstring(s(i)%short_name, &
3156 & ystring3(i)%value, lenstr), &
3157 & noerror, __line__, myfile)) RETURN
3158!
3159 IF (yaml_error(yaml_assignstring(s(i)%data_netcdf_vname, &
3160 & ystring4(i)%vector(1)%value, &
3161 & lenstr), &
3162 & noerror, __line__, myfile)) RETURN
3163!
3164 IF (yaml_error(yaml_assignstring(s(i)%data_netcdf_tname, &
3165 & ystring4(i)%vector(2)%value, &
3166 & lenstr), &
3167 & noerror, __line__, myfile)) RETURN
3168!
3169 IF (yaml_error(yaml_assignstring(s(i)%source_units, &
3170 & ystring5(i)%value, lenstr), &
3171 & noerror, __line__, myfile)) RETURN
3172!
3173 IF (yaml_error(yaml_assignstring(s(i)%destination_units, &
3174 & ystring6(i)%value, lenstr), &
3175 & noerror, __line__, myfile)) RETURN
3176!
3177 IF (yaml_error(yaml_assignstring(s(i)%source_grid, &
3178 & ystring7(i)%value, lenstr), &
3179 & noerror, __line__, myfile)) RETURN
3180!
3181 IF (yaml_error(yaml_assignstring(s(i)%destination_grid, &
3182 & ystring8(i)%value, lenstr), &
3183 & noerror, __line__, myfile)) RETURN
3184!
3185 IF (yaml_error(yaml_assignstring(s(i)%connected_to, &
3186 & ystring9(i)%value, lenstr), &
3187 & noerror, __line__, myfile)) RETURN
3188 IF (yaml_lowercase(s(i)%connected_to).eq.'false') THEN
3189 s(i)%connected=.false.
3190 ELSE
3191 s(i)%connected=.true.
3192 END IF
3193!
3194 IF (yaml_error(yaml_assignstring(s(i)%map_type, &
3195 & ystring10(i)%value, lenstr), &
3196 & noerror, __line__, myfile)) RETURN
3197!
3198 IF (yaml_error(yaml_assignstring(s(i)%map_norm, &
3199 & ystring11(i)%value, lenstr), &
3200 & noerror, __line__, myfile)) RETURN
3201 END DO
3202!
3203! Deallocate generic structures.
3204!
3205 IF (allocated(ystring1 )) deallocate (ystring1 )
3206 IF (allocated(ystring2 )) deallocate (ystring2 )
3207 IF (allocated(ystring3 )) deallocate (ystring3 )
3208 IF (allocated(ystring4 )) deallocate (ystring4 )
3209 IF (allocated(ystring5 )) deallocate (ystring5 )
3210 IF (allocated(ystring6 )) deallocate (ystring6 )
3211 IF (allocated(ystring7 )) deallocate (ystring7 )
3212 IF (allocated(ystring8 )) deallocate (ystring8 )
3213 IF (allocated(ystring9 )) deallocate (ystring9 )
3214 IF (allocated(ystring10)) deallocate (ystring10)
3215 IF (allocated(ystring11)) deallocate (ystring11)
3216 IF (allocated(ylogical1)) deallocate (ylogical1)
3217 IF (allocated(yreal1)) deallocate (yreal1)
3218 IF (allocated(yreal2)) deallocate (yreal2)
3219!
3220! Report.
3221!
3222 IF (yaml_master.and.ldebugmetadata) THEN
3223 WRITE (stdout,'(/,3a,/,3a)') &
3224 & "Coupling Metadata Dictionary, key: '", trim(key), "',", &
3225 & repeat('=',28), ' File: ', trim(filename)
3226 DO i=1,SIZE(s)
3227 WRITE (stdout,'(/,a,a)') ' - standard_name: ', &
3228 & trim(s(i)%standard_name)
3229 WRITE (stdout,'(a,a)') ' long_name: ', &
3230 & trim(s(i)%long_name)
3231 WRITE (stdout,'(a,a)') ' short_name: ', &
3232 & trim(s(i)%short_name)
3233 WRITE (stdout,'(a,a)') ' data_netcdf_variable: ', &
3234 & trim(s(i)%data_netcdf_vname)
3235 WRITE (stdout,'(a,a)') ' data_netcdf_time: ', &
3236 & trim(s(i)%data_netcdf_tname)
3237 WRITE (stdout,'(a,a)') ' source_units: ', &
3238 & trim(s(i)%source_units)
3239 WRITE (stdout,'(a,a)') ' destination_units: ', &
3240 & trim(s(i)%destination_units)
3241 WRITE (stdout,'(a,a)') ' source_grid: ', &
3242 & trim(s(i)%source_grid)
3243 WRITE (stdout,'(a,a)') ' destination_grid: ', &
3244 & trim(s(i)%destination_grid)
3245 WRITE (stdout,'(a,1p,e15.8)') ' add_offset: ', &
3246 & s(i)%add_offset
3247 WRITE (stdout,'(a,1p,e15.8)') ' scale: ', &
3248 & s(i)%scale
3249 WRITE (stdout,'(a,l1)') ' debug_write: ', &
3250 & s(i)%debug_write
3251 WRITE (stdout,'(a,l1)') ' connected: ', &
3252 & s(i)%connected
3253 WRITE (stdout,'(a,a)') ' connected_to: ', &
3254 & trim(s(i)%connected_to)
3255 WRITE (stdout,'(a,a)') ' map_type: ', &
3256 & trim(s(i)%map_type)
3257 WRITE (stdout,'(a,a)') ' map_norm: ', &
3258 & trim(s(i)%map_norm)
3259 END DO
3260 FLUSH (stdout)
3261 END IF
3262!
3263 10 FORMAT (/,' CMEPS_METADATA - Unable to create YAML object for', &
3264 & ' ROMS/CMEPS configuration metadata file: ',/,21x,a,/, &
3265 & 21x,'Default file is located in source directory.')
3266!
3267 RETURN
3268 END SUBROUTINE cmeps_metadata
3269!
3270 SUBROUTINE coupling_metadata (YML, filename, S)
3271!
3272!=======================================================================
3273! !
3274! It processes import and export field dictionary for ROMS coupling !
3275! system with the ESMF/NUOPC library. If processes field metadata !
3276! entry-by-entry from 'coupling_*.yaml'. !
3277! !
3278! On Input: !
3279! !
3280! YML YAML tree dictionary (TYPE yaml_tree) !
3281! !
3282! filename Coupling metadata filename (string) !
3283! !
3284! On Output: !
3285! !
3286! S Import/Export coupling fields (TYPE CouplingField) !
3287! !
3288!=======================================================================
3289!
3290! Imported variable declarations.
3291!
3292 character (len=*), intent(in) :: filename
3293!
3294 TYPE (yaml_tree), intent(inout) :: yml
3295 TYPE (couplingfield), allocatable, intent(out) :: s(:)
3296!
3297! Local variable declarations.
3298!
3299 integer :: i
3300!
3301 character (len=*), parameter :: myfile = &
3302 & __FILE__//", coupling_metadata"
3303!
3304!-----------------------------------------------------------------------
3305! Process coupling import/export metadata.
3306!-----------------------------------------------------------------------
3307!
3308! Create YAML dictionary.
3309!
3310 ientry=0
3311!
3312! Extract key/value pair (blocking list).
3313!
3314 IF (yaml_error(yaml_get(yml, 'metadata.standard_name', &
3315 & ystring1), &
3316 & noerror, __line__, myfile)) RETURN
3317 nentries=SIZE(ystring1)
3318!
3319 IF (yaml_error(yaml_get(yml, 'metadata.long_name', &
3320 & ystring2), &
3321 & noerror, __line__, myfile)) RETURN
3322!
3323 IF (yaml_error(yaml_get(yml, 'metadata.short_name', &
3324 & ystring3), &
3325 & noerror, __line__, myfile)) RETURN
3326!
3327 IF (yaml_error(yaml_get(yml, 'metadata.data_variables', &
3328 & ystring4), &
3329 & noerror, __line__, myfile)) RETURN
3330!
3331 IF (yaml_error(yaml_get(yml, 'metadata.source_units', &
3332 & ystring5), &
3333 & noerror, __line__, myfile)) RETURN
3334!
3335 IF (yaml_error(yaml_get(yml, 'metadata.destination_units', &
3336 & ystring6), &
3337 & noerror, __line__, myfile)) RETURN
3338!
3339 IF (yaml_error(yaml_get(yml, 'metadata.source_grid', &
3340 & ystring7), &
3341 & noerror, __line__, myfile)) RETURN
3342!
3343 IF (yaml_error(yaml_get(yml, 'metadata.destination_grid', &
3344 & ystring8), &
3345 & noerror, __line__, myfile)) RETURN
3346!
3347 IF (yaml_error(yaml_get(yml, 'metadata.connected_to', &
3348 & ystring9), &
3349 & noerror, __line__, myfile)) RETURN
3350!
3351 IF (yaml_error(yaml_get(yml, 'metadata.regrid_method', &
3352 & ystring10), &
3353 & noerror, __line__, myfile)) RETURN
3354!
3355 IF (yaml_error(yaml_get(yml, 'metadata.extrapolate_method', &
3356 & ystring11), &
3357 & noerror, __line__, myfile)) RETURN
3358!
3359 IF (allocated(yreal1)) deallocate (yreal1)
3360 allocate ( yreal1(nentries) )
3361 IF (yaml_error(yaml_get(yml, 'metadata.add_offset', &
3362 & yreal1), &
3363 & noerror, __line__, myfile)) RETURN
3364!
3365 IF (allocated(yreal2)) deallocate (yreal2)
3366 allocate ( yreal2(nentries) )
3367 IF (yaml_error(yaml_get(yml, 'metadata.scale', &
3368 & yreal2), &
3369 & noerror, __line__, myfile)) RETURN
3370!
3371 IF (allocated(ylogical1)) deallocate (ylogical1)
3372 allocate ( ylogical1(nentries) )
3373 IF (yaml_error(yaml_get(yml, 'metadata.debug_write', &
3374 & ylogical1), &
3375 & noerror, __line__, myfile)) RETURN
3376!
3377!-----------------------------------------------------------------------
3378! Load metadata information from YAML structures.
3379!-----------------------------------------------------------------------
3380!
3381 IF (.not.allocated(s)) allocate ( s(nentries) )
3382!
3383 DO i=1,nentries
3384 s(i)%debug_write = ylogical1(i)
3385 s(i)%add_offset = yreal1(i)
3386 s(i)%scale = yreal2(i)
3387!
3388 IF (yaml_error(yaml_assignstring(s(i)%standard_name, &
3389 & ystring1(i)%value, lenstr), &
3390 & noerror, __line__, myfile)) RETURN
3391!
3392 IF (yaml_error(yaml_assignstring(s(i)%long_name, &
3393 & ystring2(i)%value, lenstr), &
3394 & noerror, __line__, myfile)) RETURN
3395!
3396 IF (yaml_error(yaml_assignstring(s(i)%short_name, &
3397 & ystring3(i)%value, lenstr), &
3398 & noerror, __line__, myfile)) RETURN
3399!
3400 IF (yaml_error(yaml_assignstring(s(i)%data_netcdf_vname, &
3401 & ystring4(i)%vector(1)%value, &
3402 & lenstr), &
3403 & noerror, __line__, myfile)) RETURN
3404!
3405 IF (yaml_error(yaml_assignstring(s(i)%data_netcdf_tname, &
3406 & ystring4(i)%vector(2)%value, &
3407 & lenstr),
3408 & noerror, __line__, myfile)) RETURN
3409!
3410 IF (yaml_error(yaml_assignstring(s(i)%source_units, &
3411 & ystring5(i)%value, lenstr), &
3412 & noerror, __line__, myfile)) RETURN
3413!
3414 IF (yaml_error(yaml_assignstring(s(i)%destination_units, &
3415 & ystring6(i)%value, lenstr), &
3416 & noerror, __line__, myfile)) RETURN
3417!
3418 IF (yaml_error(yaml_assignstring(s(i)%source_grid, &
3419 & ystring7(i)%value, lenstr), &
3420 & noerror, __line__, myfile)) RETURN
3421!
3422 IF (yaml_error(yaml_assignstring(s(i)%destination_grid, &
3423 & ystring8(i)%value, lenstr), &
3424 & noerror, __line__, myfile)) RETURN
3425!
3426 IF (yaml_error(yaml_assignstring(s(i)%connected_to, &
3427 & ystring9(i)%value, lenstr), &
3428 & noerror, __line__, myfile)) RETURN
3429 IF (yaml_lowercase(s(i)%connected_to).eq.'false') THEN
3430 s(i)%connected=.false.
3431 ELSE
3432 s(i)%connected=.true.
3433 END IF
3434!
3435 IF (yaml_error(yaml_assignstring(s(i)%regrid_method, &
3436 & ystring10(i)%value, lenstr), &
3437 & noerror, __line__, myfile)) RETURN
3438!
3439 IF (yaml_error(yaml_assignstring(s(i)%extrapolate_method, &
3440 & ystring11(i)%value, lenstr), &
3441 & noerror, __line__, myfile)) RETURN
3442 END DO
3443!
3444! Deallocate generic structures.
3445!
3446 CALL yml%destroy ()
3447 IF (allocated(ystring1 )) deallocate (ystring1 )
3448 IF (allocated(ystring2 )) deallocate (ystring2 )
3449 IF (allocated(ystring3 )) deallocate (ystring3 )
3450 IF (allocated(ystring4 )) deallocate (ystring4 )
3451 IF (allocated(ystring5 )) deallocate (ystring5 )
3452 IF (allocated(ystring6 )) deallocate (ystring6 )
3453 IF (allocated(ystring7 )) deallocate (ystring7 )
3454 IF (allocated(ystring8 )) deallocate (ystring8 )
3455 IF (allocated(ystring9 )) deallocate (ystring9 )
3456 IF (allocated(ystring10)) deallocate (ystring10)
3457 IF (allocated(ystring11)) deallocate (ystring11)
3458 IF (allocated(ylogical1)) deallocate (ylogical1)
3459 IF (allocated(yreal1)) deallocate (yreal1)
3460 IF (allocated(yreal2)) deallocate (yreal2)
3461!
3462! Report.
3463!
3464 IF (yaml_master.and.ldebugmetadata) THEN
3465 WRITE (stdout,'(/,2a,/,a)') &
3466 & 'Coupling Metadata Dictionary, File: ', &
3467 & trim(filename), repeat('=',28)
3468 DO i=1,SIZE(s)
3469 WRITE (stdout,'(/,a,a)') ' - standard_name: ', &
3470 & trim(s(i)%standard_name)
3471 WRITE (stdout,'(a,a)') ' long_name: ', &
3472 & trim(s(i)%long_name)
3473 WRITE (stdout,'(a,a)') ' short_name: ', &
3474 & trim(s(i)%short_name)
3475 WRITE (stdout,'(a,a)') ' data_netcdf_variable: ', &
3476 & trim(s(i)%data_netcdf_vname)
3477 WRITE (stdout,'(a,a)') ' data_netcdf_time: ', &
3478 & trim(s(i)%data_netcdf_tname)
3479 WRITE (stdout,'(a,a)') ' source_units: ', &
3480 & trim(s(i)%source_units)
3481 WRITE (stdout,'(a,a)') ' destination_units: ', &
3482 & trim(s(i)%destination_units)
3483 WRITE (stdout,'(a,a)') ' source_grid: ', &
3484 & trim(s(i)%source_grid)
3485 WRITE (stdout,'(a,a)') ' destination_grid: ', &
3486 & trim(s(i)%destination_grid)
3487 WRITE (stdout,'(a,1p,e15.8)') ' add_offset: ', &
3488 & s(i)%add_offset
3489 WRITE (stdout,'(a,1p,e15.8)') ' scale: ', &
3490 & s(i)%scale
3491 WRITE (stdout,'(a,l1)') ' debug_write: ', &
3492 & s(i)%debug_write
3493 WRITE (stdout,'(a,l1)') ' connected: ', &
3494 & s(i)%connected
3495 WRITE (stdout,'(a,a)') ' connected_to: ', &
3496 & trim(s(i)%connected_to)
3497 WRITE (stdout,'(a,a)') ' regrid_method: ', &
3498 & trim(s(i)%regrid_method)
3499 WRITE (stdout,'(a,a)') ' extrapolate_method: ', &
3500 & trim(s(i)%extrapolate_method)
3501 END DO
3502 FLUSH (stdout)
3503 END IF
3504!
3505 RETURN
3506!
3507 10 FORMAT (/,' COUPLING_METADATA - Unable to create YAML object', &
3508 & ' for ROMS I/O metadata file: ',/,21x,a,/, &
3509 & 21x,'Default file is located in source directory.')
3510!
3511 END SUBROUTINE coupling_metadata
3512!
3513 SUBROUTINE io_metadata (YML, filename)
3514!
3515!=======================================================================
3516! !
3517! It processes ROMS input/output fields metadata entry-by-entry from !
3518! YAML tree dictionary. !
3519! !
3520! On Input: !
3521! !
3522! YML YAML tree dictionary (TYPE yaml_tree) !
3523! !
3524! filename ROMS I/O metadata filename (string) !
3525! !
3526! It reports internal variables metadata: !
3527! !
3528! Vinfo I/O Variable information (string array) !
3529! Vinfo(1): Field variable short-name !
3530! Vinfo(2): Long-name attribute !
3531! Vinfo(3): Units attribute !
3532! Vinfo(4): Field attribute !
3533! Vinfo(5): Associated time variable name !
3534! Vinfo(6): Standard-name attribute !
3535! Vinfo(7): Staggered C-grid variable type: !
3536! 'nulvar' => non-grided variable !
3537! 'p2dvar' => 2D PHI-variable !
3538! 'r2dvar' => 2D RHO-variable !
3539! 'u2dvar' => 2D U-variable !
3540! 'v2dvar' => 2D V-variable !
3541! 'p3dvar' => 3D PSI-variable !
3542! 'r3dvar' => 3D RHO-variable !
3543! 'u3dvar' => 3D U-variable !
3544! 'v3dvar' => 3D V-variable !
3545! 'w3dvar' => 3D W-variable !
3546! 'b3dvar' => 3D BED-sediment !
3547! 'l3dvar' => 3D spectral light variable !
3548! 'l4dvar' => 4D spectral light variable !
3549! Vinfo(8): Index code for information arrays !
3550! !
3551! scale Scale to convert input data to model units (real) !
3552! !
3553! offeset Value to add to input data (real) !
3554! !
3555! Ldone True if end-of-file or dictionary found !
3556! !
3557!=======================================================================
3558!
3559! Imported variable declarations.
3560!
3561 character (len=*), intent(in) :: filename
3562!
3563 TYPE (yaml_tree), intent(inout) :: YML
3564!
3565! Local variable declarations.
3566!
3567 integer :: i, j, n
3568!
3569 real (dp) :: offset, scale
3570!
3571 character (len=160) :: Vinfo(8)
3572
3573 character (len=*), parameter :: MyFile = &
3574 & __FILE__//", io_metadata"
3575!
3576!-----------------------------------------------------------------------
3577! On first pass, initialize metadata processing.
3578!-----------------------------------------------------------------------
3579!
3580! Extract values from YML tree.
3581!
3582 IF (yaml_error(yaml_get(yml, 'metadata.variable', &
3583 & ystring1), &
3584 & noerror, __line__, myfile)) RETURN
3585 nentries=SIZE(ystring1, dim=1)
3586!
3587 IF (yaml_error(yaml_get(yml, 'metadata.long_name', &
3588 & ystring2), &
3589 & noerror, __line__, myfile)) RETURN
3590!
3591 IF (yaml_error(yaml_get(yml, 'metadata.units', &
3592 & ystring3), &
3593 & noerror, __line__, myfile)) RETURN
3594!
3595 IF (yaml_error(yaml_get(yml, 'metadata.field', &
3596 & ystring4), &
3597 & noerror, __line__, myfile)) RETURN
3598!
3599 IF (yaml_error(yaml_get(yml, 'metadata.time', &
3600 & ystring5), &
3601 & noerror, __line__, myfile)) RETURN
3602!
3603 IF (yaml_error(yaml_get(yml, 'metadata.standard_name', &
3604 & ystring6), &
3605 & noerror, __line__, myfile)) RETURN
3606!
3607 IF (yaml_error(yaml_get(yml, 'metadata.type', &
3608 & ystring7), &
3609 & noerror, __line__, myfile)) RETURN
3610!
3611 IF (yaml_error(yaml_get(yml, 'metadata.index_code', &
3612 & ystring8), &
3613 & noerror, __line__, myfile)) RETURN
3614!
3615 IF (allocated(yreal1)) deallocate (yreal1)
3616 allocate ( yreal1(nentries) )
3617 IF (yaml_error(yaml_get(yml, 'metadata.add_offset', &
3618 & yreal1), &
3619 & noerror, __line__, myfile)) RETURN
3620!
3621 IF (allocated(yreal2)) deallocate (yreal2)
3622 allocate ( yreal2(nentries) )
3623 IF (yaml_error(yaml_get(yml, 'metadata.scale', &
3624 & yreal2), &
3625 & noerror, __line__, myfile)) RETURN
3626!
3627!-----------------------------------------------------------------------
3628! Process metadata entries.
3629!-----------------------------------------------------------------------
3630!
3631! Extract metadata information from YAML structures.
3632!
3633 IF (yaml_master.and.ldebugmetadata) THEN
3634 WRITE (stdout,'(/,2a,/,a)') &
3635 & 'ROMS I/O Metadata Dictionary, File: ', &
3636 & trim(filename), repeat('=',28)
3637!
3638 DO n=1,nentries
3639 DO j=1,SIZE(vinfo)
3640 DO i=1,len(vinfo(1))
3641 vinfo(j)(i:i)=char(32)
3642 END DO
3643 END DO
3644!
3645 vinfo(1)=ystring1(n)%value ! 'variable' key
3646 vinfo(2)=ystring2(n)%value ! 'long_name' key
3647 vinfo(3)=ystring3(n)%value ! 'units' key
3648 vinfo(4)=ystring4(n)%value ! 'field' key
3649 vinfo(5)=ystring5(n)%value ! 'time' key
3650 vinfo(6)=ystring6(n)%value ! 'standard_name' key
3651 vinfo(7)=ystring7(n)%value ! 'type' key
3652 vinfo(8)=ystring8(n)%value ! 'index_code' key
3653 offset =yreal1(n) ! 'add_offset' key
3654 scale =yreal2(n) ! 'scale' key
3655!
3656 WRITE (stdout,'(/,a,a)') ' - variable: ', &
3657 & trim(vinfo(1))
3658 WRITE (stdout,'(a,a)') ' standard_name: ', &
3659 & trim(vinfo(6))
3660 WRITE (stdout,'(a,a)') ' long_name: ', &
3661 & trim(vinfo(2))
3662 WRITE (stdout,'(a,a)') ' units: ', &
3663 & trim(vinfo(3))
3664 WRITE (stdout,'(a,a)') ' field: ', &
3665 & trim(vinfo(4))
3666 WRITE (stdout,'(a,a)') ' time: ', &
3667 & trim(vinfo(5))
3668 WRITE (stdout,'(a,a)') ' index_code: ', &
3669 & trim(vinfo(8))
3670 WRITE (stdout,'(a,a)') ' type: ', &
3671 & trim(vinfo(7))
3672 WRITE (stdout,'(a,1p,e15.8)') ' add_offset: ', &
3673 & offset
3674 WRITE (stdout,'(a,1p,e15.8)') ' scale: ', &
3675 & scale
3676 END DO
3677 FLUSH (stdout)
3678 END IF
3679!
3680! Destroy and deallocate.
3681!
3682 CALL yml%destroy ()
3683 IF (allocated(ystring1)) deallocate (ystring1)
3684 IF (allocated(ystring2)) deallocate (ystring2)
3685 IF (allocated(ystring3)) deallocate (ystring3)
3686 IF (allocated(ystring4)) deallocate (ystring4)
3687 IF (allocated(ystring5)) deallocate (ystring5)
3688 IF (allocated(ystring6)) deallocate (ystring6)
3689 IF (allocated(ystring7)) deallocate (ystring7)
3690 IF (allocated(ystring8)) deallocate (ystring8)
3691 IF (allocated(yreal1)) deallocate (yreal1)
3692 IF (allocated(yreal2)) deallocate (yreal2)
3693!
3694 RETURN
3695 END SUBROUTINE io_metadata
3696!
3697!-----------------------------------------------------------------------
3698 END MODULE yaml_parser
3699!-----------------------------------------------------------------------
3700!
3702!
3703!=======================================================================
3704! !
3705! This program can be used to test the YAML parser available in ROMS. !
3706! !
3707!=======================================================================
3708!
3709 USE yaml_parser, ONLY : process_yaml
3710!
3711! Local variable definitions.
3712!
3713 character (len=120) :: yaml_file ! YAML filename
3714
3715 character (len=*), parameter :: myfile = &
3716 & __FILE__//', read_yaml'
3717!
3718! Get YAML file name.
3719!
3720 DO WHILE (.true.)
3721 WRITE (6,10)
3722 10 FORMAT (/,' Enter YAML filename: ',$)
3723 READ(5,'(a)',err=20) yaml_file
3724!
3725! Create and process YAML object.
3726!
3727 CALL process_yaml (trim(yaml_file))
3728 END DO
3729
3730 20 stop
3731
3732 END PROGRAM yaml_parser_test
3733
integer, parameter iunit
type(yaml_svec), dimension(:), allocatable ystring8
logical function yaml_tree_has(self, keystring)
integer function, private yaml_countkeys(string, token)
integer function yaml_get_s_struc(self, keystring, v)
integer function, public yaml_initialize(self, filename, report)
integer function yaml_tree_read_line(self, nblanks, line, key, value, anchor, lswitch)
integer, parameter lmax
type(yaml_svec), dimension(:), allocatable ystring4
integer, parameter ldim
subroutine yaml_tree_fill_aliases(self, nalias, nanchor)
subroutine io_metadata(yml, filename)
logical, dimension(:), allocatable ylogical1
logical yaml_master
integer function yaml_get_ivar_0d(self, keystring, value)
integer function yaml_tree_extract(self, keystring, s)
integer, parameter lkey
character(len=len(sinp)) function, private yaml_uppercase(sinp)
type(yaml_svec), dimension(:), allocatable ystring6
type(yaml_svec), dimension(:), allocatable ystring12
logical, parameter ldebugyaml
integer function, public yaml_assignstring(outstring, inpstring, lstr)
integer function yaml_get_lvar_1d(self, keystring, value)
subroutine yaml_tree_destroy(self)
integer, parameter dp
subroutine yaml_tree_fill(self)
real(kind_real), dimension(:), allocatable yreal1
type(yaml_svec), dimension(:), allocatable ystring9
character(len=55), dimension(7) yaml_errmsg
type(yaml_svec), dimension(:), allocatable ystring7
type(yaml_svec), dimension(:), allocatable ystring3
integer function yaml_get_rvar_0d(self, keystring, value)
integer function yaml_get_i_struc(self, keystring, v)
subroutine, public process_yaml(yaml_file)
subroutine yaml_tree_dump(self)
integer function yaml_get_svar_0d(self, keystring, value)
integer function yaml_get_rvar_1d(self, keystring, value)
integer, parameter noerror
integer function yaml_get_l_struc(self, keystring, v)
integer, parameter kind_real
subroutine, public coupling_metadata(yml, filename, s)
integer, parameter stdout
logical function, public yaml_error(flag, noerr, line, routine)
logical ldebugmetadata
type(yaml_svec), dimension(:), allocatable ystring2
integer function yaml_get_svar_1d(self, keystring, value)
integer, parameter yaml_stdout
subroutine, private yaml_valuetype(value, lswitch)
subroutine, public cmeps_metadata(self, filename, key, s)
subroutine yaml_tree_create(self)
type(yaml_svec), dimension(:), allocatable ystring1
type(yaml_svec), dimension(:), allocatable ystring11
type(yaml_svec), dimension(:), allocatable ystring5
integer function yaml_get_r_struc(self, keystring, v)
type(yaml_svec), dimension(:), allocatable ystring10
real(kind_real), dimension(:), allocatable yreal2
character(len=len(sinp)) function, private yaml_lowercase(sinp)
integer function yaml_get_lvar_0d(self, keystring, value)
integer function yaml_get_ivar_1d(self, keystring, value)
integer, parameter noerr
program yaml_parser_test