chore: append -dev to version number (#1641)
[FMS.git] / test_fms / parser / parser_demo2.F90
blob674ab85fbbc8a4c7c3b8973f7a7bcec12f9d3799
1 !***********************************************************************
2 !*                   GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
14 !* for more details.
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS.  If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
20 program parser_demo
21 !> @brief  This programs demostrates how to use the parser
23 #ifdef use_yaml
24 use FMS_mod, only: fms_init, fms_end
25 use yaml_parser_mod
26 use platform_mod
28 implicit none
30 integer              :: diag_yaml_id   !< Id for the diag_table yaml
31 integer              :: nfiles         !< Number of files in the diag_table yaml
32 integer, allocatable :: file_ids(:)    !< Ids of the files in the diag_table yaml
33 integer              :: nvariables     !< Number of variables in the diag_table yaml
34 integer, allocatable :: var_ids(:)     !< Ids of the variables in the diag_table yaml
35 integer              :: i, j, k        !< For do loops
36 integer              :: nkeys          !< Number of keys
37 integer, allocatable :: key_ids(:)     !< Ids of keys in the diag_table_yaml
38 character(len=255)   :: key_value      !< The value of a key
39 character(len=255)   :: key_name       !< The name of a key
41 call fms_init
43 diag_yaml_id = open_and_parse_file("diag_table.yaml")
44 print *, ""
46 nkeys = get_nkeys(diag_yaml_id, 0)
47 allocate(key_ids(nkeys))
48 call get_key_ids(diag_yaml_id, 0, key_ids)
50 do i = 1, nkeys
51     call get_key_name(diag_yaml_id, key_ids(i), key_name)
52     call get_key_value(diag_yaml_id, key_ids(i), key_value)
53     print *, "Key:", trim(key_name), " Value:", trim(key_value)
54 enddo
56 deallocate(key_ids)
58 nfiles = get_num_blocks(diag_yaml_id, "diag_files")
59 allocate(file_ids(nfiles))
60 call get_block_ids(diag_yaml_id, "diag_files", file_ids)
61 print *, ""
63 do i = 1, nfiles
64    print *, "File number:", i
66    nkeys = get_nkeys(diag_yaml_id, file_ids(i))
67    allocate(key_ids(nkeys))
68    call get_key_ids(diag_yaml_id, file_ids(i), key_ids)
70    do j = 1, nkeys
71       call get_key_name(diag_yaml_id, key_ids(j), key_name)
72       call get_key_value(diag_yaml_id, key_ids(j), key_value)
73       print *, "  Key:", trim(key_name), " Value:", trim(key_value)
74    enddo
76    deallocate(key_ids)
77    print *, ""
78    !< The number of variables that are part of the current file
79    nvariables = get_num_blocks(diag_yaml_id, "varlist", parent_block_id=file_ids(i))
80    allocate(var_ids(nvariables))
81    call get_block_ids(diag_yaml_id, "varlist", var_ids, parent_block_id=file_ids(i))
83    do j = 1, nvariables
84        print *, "  Variable number:", j
86        nkeys = get_nkeys(diag_yaml_id, var_ids(j))
87        allocate(key_ids(nkeys))
88        call get_key_ids(diag_yaml_id, var_ids(j), key_ids)
90        do k = 1, nkeys
91           call get_key_name(diag_yaml_id, key_ids(k), key_name)
92           call get_key_value(diag_yaml_id, key_ids(k), key_value)
93           print *, "     Key:", trim(key_name), " Value:", trim(key_value)
94        enddo
96        deallocate(key_ids)
97        print *, ""
98    end do
100    deallocate(var_ids)
101    print *, ""
102 enddo
103 deallocate(file_ids)
104 call fms_end
106 #endif
108 end program parser_demo