UPP (upp-srw-2.2.0)
Loading...
Searching...
No Matches
xml_perl_data.f
1 module xml_perl_data
2!------------------------------------------------------------------------
3!
4! This module read in Perl XML processed flat file and
5! handle parameter marshalling for existing POST program
6!
7! program log:
8! March, 2015 Lin Gan Initial Code
9! July, 2016 J. Carley Clean up prints
10!
11!------------------------------------------------------------------------
12 implicit none
14 integer :: post_avblfldidx=-9999
15 character(len=80) :: shortname=''
16 character(len=300) :: longname=''
17 integer :: mass_windpoint=1
18 character(len=30) :: pdstmpl='tmpl4_0'
19 character(len=30) :: pname=''
20 character(len=10) :: table_info=''
21 character(len=80) :: stats_proc=''
22 character(len=80) :: fixed_sfc1_type=''
23 integer, dimension(:), pointer :: scale_fact_fixed_sfc1 => null()
24 real, dimension(:), pointer :: level => null()
25 character(len=80) :: fixed_sfc2_type=''
26 integer, dimension(:), pointer :: scale_fact_fixed_sfc2 => null()
27 real, dimension(:), pointer :: level2 => null()
28 character(len=80) :: aerosol_type=''
29 character(len=80) :: typ_intvl_size=''
30 integer :: scale_fact_1st_size=0
31 real :: scale_val_1st_size=0.0
32 integer :: scale_fact_2nd_size=0
33 real :: scale_val_2nd_size=0.0
34 character(len=80) :: typ_intvl_wvlen=''
35 integer :: scale_fact_1st_wvlen=0
36 real :: scale_val_1st_wvlen=0.0
37 integer :: scale_fact_2nd_wvlen=0
38 real :: scale_val_2nd_wvlen=0.0
39 real, dimension(:), pointer :: scale => null()
40 integer :: stat_miss_val=0
41 integer :: leng_time_range_prev=0
42 integer :: time_inc_betwn_succ_fld=0
43 character(len=80) :: type_of_time_inc=''
44 character(len=20) :: stat_unit_time_key_succ=''
45 character(len=20) :: bit_map_flag=''
46 end type param_t
47
49 character(len=6) :: datset=''
50 integer :: grid_num=255
51 character(len=20) :: sub_center=''
52 character(len=20) :: version_no=''
53 character(len=20) :: local_table_vers_no=''
54 character(len=20) :: sigreftime=''
55 character(len=20) :: prod_status=''
56 character(len=20) :: data_type=''
57 character(len=20) :: gen_proc_type=''
58 character(len=30) :: time_range_unit=''
59 character(len=50) :: orig_center=''
60 character(len=30) :: gen_proc=''
61 character(len=50) :: packing_method=''
62 character(len=30) :: order_of_sptdiff='1st_ord_sptdiff'
63 character(len=20) :: field_datatype=''
64 character(len=30) :: comprs_type=''
65 character(len=50) :: type_ens_fcst=''
66 character(len=50) :: type_derived_fcst=''
67 type(param_t), dimension(:), pointer :: param => null()
68 end type paramset_t
69
71 type(param_t), dimension(:), pointer :: param => null()
72 end type post_avblfld_t
73
74 type (paramset_t), dimension(:), pointer :: paramset
75 type (post_avblfld_t),save :: post_avblflds
76
77 contains
78 subroutine read_postxconfig()
79
80 use rqstfld_mod,only: num_post_afld,mxlvl,lvlsxml
81 use ctlblk_mod, only:tprec,tclod,trdlw,trdsw,tsrfc &
82 ,tmaxmin,td3d,me,filenameflat
83 implicit none
84
85! Read in the flat file postxconfig-NT.txt
86! for current working parameters and param
87 integer paramset_count, param_count
88
89! temp array count
90 integer cc
91 integer level_array_count
92 integer cv
93 integer level2_array_count
94 integer scale_array_count
95 integer i,j
96
97! evil for empty default char "?"
98 character(len=80) testcharname
99 character dummy_char
100 integer testintname
101
102! open the Post flat file
103! open(UNIT=22,file="postxconfig-NT.txt", &
104 open(unit=22,file=trim(filenameflat), &
105 form="formatted", access="sequential", &
106 status="old", position="rewind")
107
108! Take the first line as paramset_count
109 read(22,*)paramset_count
110
111! Allocate paramset array size
112 allocate(paramset(paramset_count))
113
114! Take the second line as param_count (on n..1 down loop)
115! stored as FILO
116
117! Initialize num_post_afld here
118 num_post_afld = 0
119
120 do i = paramset_count, 1, -1
121 read(22,*)param_count
122
123 allocate(paramset(i)%param(param_count))
124
125! LinGan lvlsxml is now a sum of flat file read out
126! Also allocate lvlsxml for rqstfld_mod
127 num_post_afld = num_post_afld + param_count
128
129 end do
130
131 if(allocated(lvlsxml)) deallocate(lvlsxml)
132 allocate(lvlsxml(mxlvl,num_post_afld))
133
134! For each paramset_count to read in all 16 control contain
135 do i = 1, paramset_count
136! allocate array size from param for current paramset
137! filter_char_inp is to check if "?" is found
138! then replace to empty string because it means no input.
139
140 read(22,*)paramset(i)%datset
141 call filter_char_inp(paramset(i)%datset)
142
143 param_count = size (paramset(i)%param)
144
145 read(22,*)paramset(i)%grid_num
146 read(22,*)paramset(i)%sub_center
147 call filter_char_inp(paramset(i)%sub_center)
148 read(22,*)paramset(i)%version_no
149 call filter_char_inp(paramset(i)%version_no)
150 read(22,*)paramset(i)%local_table_vers_no
151 call filter_char_inp(paramset(i)%local_table_vers_no)
152 read(22,*)paramset(i)%sigreftime
153 call filter_char_inp(paramset(i)%sigreftime)
154 read(22,*)paramset(i)%prod_status
155 call filter_char_inp(paramset(i)%prod_status)
156 read(22,*)paramset(i)%data_type
157 call filter_char_inp(paramset(i)%data_type)
158 read(22,*)paramset(i)%gen_proc_type
159 call filter_char_inp(paramset(i)%gen_proc_type)
160 read(22,*)paramset(i)%time_range_unit
161 call filter_char_inp(paramset(i)%time_range_unit)
162 read(22,*)paramset(i)%orig_center
163 call filter_char_inp(paramset(i)%orig_center)
164 read(22,*)paramset(i)%gen_proc
165 call filter_char_inp(paramset(i)%gen_proc)
166 read(22,*)paramset(i)%packing_method
167 call filter_char_inp(paramset(i)%packing_method)
168 read(22,*)paramset(i)%order_of_sptdiff
169 read(22,*)paramset(i)%field_datatype
170 call filter_char_inp(paramset(i)%field_datatype)
171 read(22,*)paramset(i)%comprs_type
172 call filter_char_inp(paramset(i)%comprs_type)
173 if(paramset(i)%gen_proc_type=='ens_fcst')then
174 read(22,*)paramset(i)%type_ens_fcst
175 call filter_char_inp(paramset(i)%type_ens_fcst)
176 tprec = 6 ! always 6 hr bucket for gefs
177 tclod = tprec
178 trdlw = tprec
179 trdsw = tprec
180 tsrfc = tprec
181 tmaxmin = tprec
182 td3d = tprec
183 end if
184! Loop param_count (param datas 161) for gfsprs
185 do j = 1, param_count
186 read(22,*)paramset(i)%param(j)%post_avblfldidx
187 read(22,*)paramset(i)%param(j)%shortname
188 read(22,'(A300)')paramset(i)%param(j)%longname
189 call filter_char_inp(paramset(i)%param(j)%longname)
190
191 read(22,*)paramset(i)%param(j)%mass_windpoint
192 read(22,*)paramset(i)%param(j)%pdstmpl
193 read(22,*)paramset(i)%param(j)%pname
194 call filter_char_inp(paramset(i)%param(j)%pname)
195
196 read(22,*)paramset(i)%param(j)%table_info
197 call filter_char_inp(paramset(i)%param(j)%table_info)
198 read(22,*)paramset(i)%param(j)%stats_proc
199 call filter_char_inp(paramset(i)%param(j)%stats_proc)
200 read(22,*)paramset(i)%param(j)%fixed_sfc1_type
201 call filter_char_inp(paramset(i)%param(j)%fixed_sfc1_type)
202! Read array count for scale_fact_fixed_sfc1
203 read(22,*)cc
204!
205 allocate( paramset(i)%param(j)%scale_fact_fixed_sfc1(1))
206
207 if (cc > 0) then
208!
209 deallocate( paramset(i)%param(j)%scale_fact_fixed_sfc1)
210
211 allocate( paramset(i)%param(j)%scale_fact_fixed_sfc1(cc))
212 read(22,*)paramset(i)%param(j)%scale_fact_fixed_sfc1
213 else
214! If array count is zero dummy out the line
215!
216 paramset(i)%param(j)%scale_fact_fixed_sfc1(1)=0
217
218 read(22,*)dummy_char
219 endif
220
221 read(22,*)level_array_count
222 allocate( paramset(i)%param(j)%level(1))
223 if (level_array_count > 0) then
224 deallocate( paramset(i)%param(j)%level)
225 allocate( paramset(i)%param(j)%level(level_array_count))
226 read(22,*)paramset(i)%param(j)%level
227 else
228 paramset(i)%param(j)%level(1)=0
229 read(22,*)dummy_char
230 endif
231
232 read(22,*)paramset(i)%param(j)%fixed_sfc2_type
233 call filter_char_inp(paramset(i)%param(j)%fixed_sfc2_type)
234 read(22,*)cv
235 allocate( paramset(i)%param(j)%scale_fact_fixed_sfc2(1))
236 if (cv > 0) then
237 deallocate(paramset(i)%param(j)%scale_fact_fixed_sfc2)
238 allocate(paramset(i)%param(j)%scale_fact_fixed_sfc2(cv))
239 read(22,*)paramset(i)%param(j)%scale_fact_fixed_sfc2
240 else
241 paramset(i)%param(j)%scale_fact_fixed_sfc2(1)=0
242 read(22,*)dummy_char
243 endif
244
245 read(22,*)level2_array_count
246 if (level2_array_count > 0) then
247 allocate(paramset(i)%param(j)%level2(level2_array_count))
248 read(22,*)paramset(i)%param(j)%level2
249 else
250 read(22,*)dummy_char
251 endif
252
253 read(22,*)paramset(i)%param(j)%aerosol_type
254 call filter_char_inp(paramset(i)%param(j)%aerosol_type)
255 read(22,*)paramset(i)%param(j)%typ_intvl_size
256 call filter_char_inp(paramset(i)%param(j)%typ_intvl_size)
257
258 read(22,*)paramset(i)%param(j)%scale_fact_1st_size
259 read(22,*)paramset(i)%param(j)%scale_val_1st_size
260 read(22,*)paramset(i)%param(j)%scale_fact_2nd_size
261 read(22,*)paramset(i)%param(j)%scale_val_2nd_size
262 read(22,*)paramset(i)%param(j)%typ_intvl_wvlen
263 call filter_char_inp(paramset(i)%param(j)%typ_intvl_wvlen)
264
265 read(22,*)paramset(i)%param(j)%scale_fact_1st_wvlen
266 read(22,*)paramset(i)%param(j)%scale_val_1st_wvlen
267 read(22,*)paramset(i)%param(j)%scale_fact_2nd_wvlen
268 read(22,*)paramset(i)%param(j)%scale_val_2nd_wvlen
269 read(22,*)scale_array_count
270 allocate(paramset(i)%param(j)%scale(1))
271 if (scale_array_count > 0) then
272 deallocate(paramset(i)%param(j)%scale)
273 allocate(paramset(i)%param(j)%scale(scale_array_count))
274 read(22,*)paramset(i)%param(j)%scale
275 else
276 paramset(i)%param(j)%scale(1)=0
277 read(22,*)dummy_char
278 endif
279 read(22,*)paramset(i)%param(j)%stat_miss_val
280 read(22,*)paramset(i)%param(j)%leng_time_range_prev
281 read(22,*)paramset(i)%param(j)%time_inc_betwn_succ_fld
282 read(22,*)paramset(i)%param(j)%type_of_time_inc
283
284 call filter_char_inp(paramset(i)%param(j)%type_of_time_inc)
285 read(22,*)paramset(i)%param(j)%stat_unit_time_key_succ
286 call filter_char_inp(paramset(i)%param(j)%stat_unit_time_key_succ)
287 read(22,*)paramset(i)%param(j)%bit_map_flag
288 call filter_char_inp(paramset(i)%param(j)%bit_map_flag)
289
290! End of reading param
291 end do
292
293 post_avblflds%param => paramset(i)%param
294
295! End of reading paramset
296 end do
297 close (unit=22)
298
299 end subroutine read_postxconfig
300
301
302 subroutine filter_char_inp (inpchar)
303 implicit none
304 character, intent(inout) :: inpchar
305 if (inpchar == "?") then
306 inpchar = ""
307 endif
308 end subroutine filter_char_inp
309
310 end module
311