UPP (upp-srw-2.2.0)
Loading...
Searching...
No Matches
SET_LVLSXML.f
1 subroutine set_lvlsxml(param,ifld,irec,kpv,pv,kth,th)
2!
3!$$$ SUBPROGRAM DOCUMENTATION BLOCK
4! . . .
5! SUBPROGRAM: SET_LVLCXML SET field levels from POST xml CONTROL FILE
6! PRGRMMR: J. WANG ORG: NCEP/EMC DATE: 12-01-27
7!
8! ABSTRACT:
9! THIS ROUTINE SET THE LVLS and LVLSXML for contain request field.
10!
11! PROGRAM HISTORY LOG:
12! 01_27_2012 Jun Wang - INITIAL CODE
13! 04_03_2012 Jun Wang - add SPEC_PRES_ABOVE_GRND for different CAPE/CIN
14! 08_06_2013 S Moorthi - fix index out of bound after iloop5
15! 10_03_2013 Jun Wang - add isentropic levels
16! 03_10_2015 Lin Gan - Replace XML file with flat file implementation
17! 07_08_2016 J. Carley - Comment out debug prints
18! 06_01_2017 Y Mao - For MISCLN.f and FDLVL.f, allow FD levels input from control file
19! 01-24-2023 Sam Trahan - IFI flight levels and related unit conversions
20!
21! USAGE: CALL SET_LVLSXML(param,ifld,irec,kpv,pv,kth,th)
22! INPUT ARGUMENT LIST:
23! param: input field
24! ifld : field number in post control file
25! irec : data fields number in output file
26! kpv : total number of potential vorticity levels
27! pv : potential vorticity levels
28! kth : total number of isentropic levels
29! th : isentropic levels
30!
31! OUTPUT ARGUMENT LIST:
32!
33! OUTPUT FILES:
34! NONE
35!
36! SUBPROGRAMS CALLED:
37! UTILITIES:
38!
39! LIBRARY:
40! MODULE: - RQSTFLD_MOD
41! CTLBLK_MOD
42! xml_data_post_t
43! SOIL
44!
45! ATTRIBUTES:
46! LANGUAGE: FORTRAN
47! MACHINE : IBM
48!
49 use xml_perl_data, only: param_t
50 use ctlblk_mod, only: lsm, spl, nsoil, isf_surface_physics, nfd, htfd, &
51 petabnd, nbnd, ifi_nflight, ifi_flight_levels
52 use soil, only: sldpth,sllevel
53 use rqstfld_mod,only : mxlvl,lvls,lvlsxml
54 implicit none
55!
56 type(param_t),intent(inout) :: param
57 integer, intent(in) :: ifld
58 integer, intent(inout) :: irec
59 integer, intent(in) :: kpv
60 real,intent(in) :: pv(1:kpv)
61 integer, intent(in) :: kth
62 real,intent(in) :: th(1:kth)
63!
64 real,parameter :: small=1.e-5
65 real,parameter :: small1=1.e-3
66 real,parameter :: small2=1
67 integer,parameter :: LSIG1=22,lsig2=5
68 integer i,j,l,nlevel,scalef,lvlcape,lvlcin
69 logical READTHK,logrec,found
70 REAL :: SIGO2(LSIG2+1),ASIGO2(LSIG2),DSIGO2(LSIG2)
71 REAL :: SIGO1(LSIG1+1),ASIGO1(LSIG1),DSIGO1(LSIG1)
72!
73 readthk=.false.
74!
75 nlevel=size(param%level)
76!
77
78 if(nlevel<=0) then
79 lvls(1,ifld)=1
80 lvlsxml(1,ifld)=1
81 irec=irec+1
82 return
83 endif
84
85 if(trim(param%fixed_sfc1_type)=='isobaric_sfc') then
86 if(index(param%shortname,"ON_ICAO_STD_SFC")<=0) then
87 do j=1, nlevel
88 iloop: do i=1, lsm
89
90 if(abs(param%level(j)-spl(i))<small1)then
91 lvls(i,ifld)=1
92 lvlsxml(i,ifld)=j
93 irec=irec+1
94 exit iloop
95 endif
96 enddo iloop
97 enddo
98 return
99 endif
100! For standard atmospheric pressures, use levels from control file directly
101 do j=1, nlevel
102 lvls(j,ifld)=1
103 lvlsxml(j,ifld)=j
104 irec=irec+1
105 enddo
106 return
107 endif
108!
109 if(trim(param%fixed_sfc1_type)=='hybrid_lvl') then
110 do j=1, nlevel
111 iloop1: do i=1, mxlvl
112 if(nint(param%level(j))==i)then
113 lvls(i,ifld)=1
114 lvlsxml(i,ifld)=j
115 irec=irec+1
116 exit iloop1
117 endif
118 enddo iloop1
119 enddo
120 return
121 endif
122!
123 if(trim(param%fixed_sfc1_type)=='depth_bel_land_sfc'.and. &
124 trim(param%fixed_sfc2_type)=='depth_bel_land_sfc' ) then
125! if(me==0)print *,'nsoil=',nsoil,'iSF_SURFACE_PHYSICS=',iSF_SURFACE_PHYSICS, &
126! 'level=',param%level,'sldpth=',SLDPTH(1:nsoil),'sum=',sum(SLDPTH(1:nsoil))*100.
127 do j=1, nlevel
128 iloop2: do i=1, nsoil
129 if(isf_surface_physics ==3) then
130 if(nint(param%level(j))==nint(sllevel(i)*100.)) then
131 lvls(i,ifld)=1
132 lvlsxml(i,ifld)=j
133 irec=irec+1
134 exit iloop2
135 endif
136 else
137 if(nint(param%level2(j))==nint(sum(sldpth(1:i))*100.) ) then
138 lvls(i,ifld)=1
139 lvlsxml(i,ifld)=j
140 irec=irec+1
141 exit iloop2
142 endif
143 endif
144 enddo iloop2
145 enddo
146 return
147 endif
148!
149!for pv sfc
150 if(trim(param%fixed_sfc1_type)=='pot_vort_sfc') then
151 do j=1, nlevel
152 scalef=param%scale_fact_fixed_sfc1(j)-6
153 if(param%scale_fact_fixed_sfc1(j)<6) scalef=0
154 iloop3: do i=1, kpv
155 if(pv(i)/=0.and.abs(param%level(j)*10.**(-1*scalef)-pv(i))<=1.e-5) then
156 lvls(i,ifld)=1
157 lvlsxml(i,ifld)=j
158 irec=irec+1
159 exit iloop3
160 endif
161 enddo iloop3
162 enddo
163! print *,'for level type pv,nlevel=',nlevel,'level=', &
164! param%level(1:nlevel)*10.**(-1*scalef), &
165! 'pv=',pv(1:kpv),lvls1(1:kpv),'ifld=',ifld,'var=',trim(param%pname), &
166! 'lvl type=',trim(param%fixed_sfc1_type)
167 return
168 endif
169!
170!for th sfc
171 if(trim(param%fixed_sfc1_type)=='isentropic_lvl') then
172 do j=1, nlevel
173 !print *,'in set_lvl,kth=',kth,'nlevel=',nlevel,'j=',j,param%level(j)
174 iloop3a: do i=1, kth
175 if(th(i)/=0.and.abs(param%level(j)-th(i))<=1.e-5) then
176 lvls(i,ifld)=1
177 lvlsxml(i,ifld)=j
178 irec=irec+1
179 exit iloop3a
180 endif
181 enddo iloop3a
182 enddo
183! print *,'for level type th,nlevel=',nlevel,'level=', &
184! param%level(1:nlevel), &
185! 'th=',th(1:kth),'ifld=',ifld,'var=',trim(param%pname), &
186! 'lvl type=',trim(param%fixed_sfc1_type)
187 return
188 endif
189!
190 if(trim(param%fixed_sfc1_type)=='spec_alt_above_mean_sea_lvl') then
191 if(index(param%shortname,"SPECIFIC_IFI_FLIGHT_LEVEL")>0) then
192 do j=1, nlevel
193 found=.false.
194 iloop411: do i=1, ifi_nflight
195 if(nint(param%level(j)/10)==nint(ifi_flight_levels(i)/10) )then
196 lvls(i,ifld)=1
197 lvlsxml(i,ifld)=j
198 !print *,'SPECIFIC_IFI_FLIGHT_LEVEL ',j,' is ',param%level(j)
199 irec=irec+1
200 found=.true.
201 exit iloop411
202 endif
203 enddo iloop411
204 if(.not.found) then
205 write(*,*) 'ERROR: No such IFI flight level: ',param%level(j)/10
206 lvls(i,ifld)=0
207 endif
208 enddo
209 else if(index(param%shortname,"IFI_FLIGHT_LEVEL")>0) then
210 do j=1, ifi_nflight
211 lvls(j,ifld)=1
212 lvlsxml(j,ifld)=j
213 !print *,'IFI_FLIGHT_LEVEL ',j,' is ',param%level(j)
214 irec=irec+1
215 enddo
216 elseif(index(param%shortname,"GTG_ON_SPEC_ALT_ABOVE_MEAN_SEA_LVL")<=0) then
217 do j=1, nlevel
218 iloop4: do i=1, nfd
219 if(nint(param%level(j))==nint(htfd(i)) )then
220 if(htfd(i)>300.) then
221 lvls(i,ifld)=1
222 else
223 lvls(i,ifld)=2
224 endif
225 lvlsxml(i,ifld)=j
226 irec=irec+1
227 exit iloop4
228 endif
229 enddo iloop4
230 enddo
231 return
232 endif
233! Allow inputs of FD levels from control file. For GTG (EDPARM CATEDR MWTURB)
234! SET LVLS to 1
235 do j=1, nlevel
236 lvls(j,ifld)=1
237 lvlsxml(j,ifld)=j
238 irec=irec+1
239 enddo
240! print *, "GTG levels, n=",nlevel, "irec=",irec
241 return
242 endif
243!
244 if(trim(param%fixed_sfc1_type)=='spec_pres_above_grnd') then
245 logrec=.false.
246 if(trim(param%shortname)=="MIXED_LAYER_CAPE_ON_SPEC_PRES_ABOVE_GRND" .or. &
247 trim(param%shortname)=="MIXED_LAYER_CIN_ON_SPEC_PRES_ABOVE_GRND") then
248 lvlsxml(1,ifld)=1
249 irec=irec+1
250! allocate(param%level(1),param%level2(1))
251 param%level(1)=nint(petabnd(3)+15.)*100
252 param%level2(1)=nint(petabnd(1)-15.)*100
253 else if (trim(param%shortname)=="UNSTABLE_CAPE_ON_SPEC_PRES_ABOVE_GRND" .or. &
254 trim(param%shortname)=="UNSTABLE_CIN_ON_SPEC_PRES_ABOVE_GRND") then
255 lvlsxml(1,ifld)=1
256! allocate(param%level(1),param%level2(1))
257 param%level(1)=25500
258 irec=irec+1
259 param%level2(1)=0
260 else if (trim(param%shortname)=="BEST_CAPE_ON_SPEC_PRES_ABOVE_GRND" .or. &
261 trim(param%shortname)=="BEST_CIN_ON_SPEC_PRES_ABOVE_GRND") then
262 !print *,'in set_vlv,best cape'
263 lvlsxml(1,ifld)=1
264 irec=irec+1
265! allocate(param%level(1),param%level2(1))
266 param%level(1)=nint(petabnd(nbnd)+15.)*100
267 param%level2(1)=nint(petabnd(1)-15.)*100
268 else
269 do j=1, nlevel
270 iloop5: do i=1, nbnd
271 if(nint(param%level(j)/100.)==nint(petabnd(i)+15.))then
272 lvls(i,ifld)=1
273 lvlsxml(i,ifld)=j
274 irec=irec+1
275 logrec=.true.
276 exit iloop5
277 endif
278 enddo iloop5
279 if(nint(param%level(j)/100.) == 255) then
280 lvls(1,ifld) = 1
281 lvlsxml(1,ifld) = j
282 irec = irec+1
283 endif
284 enddo
285 if(.not.logrec.and.nlevel==1) then
286 lvls(1,ifld)=1
287 lvlsxml(1,ifld)=1
288 irec=irec+1
289 endif
290 endif
291 return
292 endif
293!
294 if(trim(param%fixed_sfc1_type)=='spec_hgt_lvl_above_grnd') then
295 if(index(param%shortname,"SPEC_HGT_LVL_ABOVE_GRND_FDHGT")>0) then
296 do j=1, nlevel
297 iloop41: do i=1, nfd
298 if(nint(param%level(j))==nint(htfd(i)) )then
299 lvls(i,ifld)=1
300 lvlsxml(i,ifld)=j
301 irec=irec+1
302 exit iloop41
303 endif
304 enddo iloop41
305 enddo
306 return
307 endif
308 do j=1, nlevel
309 lvls(j,ifld)=1
310 lvlsxml(j,ifld)=j
311 irec=irec+1
312 enddo
313 return
314 endif
315!
316!for hpc tmp at sigma lvl
317 if(trim(param%shortname)=='TMP_ON_SIGMA_LVL_HPC') then
318 IF(readthk)THEN ! EITHER READ DSG THICKNESS
319 READ(41)dsigo2 !DSIGO FROM TOP TO BOTTOM
320!
321 sigo2(1)=0.0
322 DO l=2,lsig2+1
323 sigo2(l)=sigo2(l-1)+dsigo2(lsig2-l+2)
324 END DO
325 sigo2(lsig2+1)=1.0
326 DO l=1,lsig2
327 asigo2(l)=0.5*(sigo2(l)+sigo2(l+1))
328 END DO
329 ELSE ! SPECIFY SIGO
330 asigo2( 1)= 0.7000
331 asigo2( 2)= 0.7500
332 asigo2( 3)= 0.8000
333 asigo2( 4)= 0.8500
334 asigo2( 5)= 0.9000
335 END IF
336!
337 do j=1, nlevel
338 DO i=1,lsig2
339 if(abs(param%level(j)-asigo2(i)*10000)<small1) then
340 lvls(i,ifld)=1
341 lvlsxml(i,ifld)=j
342 irec=irec+1
343 endif
344 enddo
345 enddo
346 return
347!
348 ENDIF
349!
350!for hpc tmp at sigma lvl
351 if(index(trim(param%shortname),'SIGMA_LVLS')>0) then
352 IF(readthk)THEN ! EITHER READ DSG THICKNESS
353 READ(41)dsigo1 !DSIGO FROM TOP TO BOTTOM
354!
355 sigo1(1)=0.0
356 DO l=2,lsig1+1
357 sigo1(l)=sigo1(l-1)+dsigo1(lsig1-l+2)
358 END DO
359 sigo1(lsig1+1)=1.0
360 DO l=1,lsig1
361 asigo1(l)=0.5*(sigo1(l)+sigo1(l+1))
362 END DO
363 ELSE ! SPECIFY SIGO
364 asigo1( 1)= 0.0530
365 asigo1( 2)= 0.1580
366 asigo1( 3)= 0.2605
367 asigo1( 4)= 0.3595
368 asigo1( 5)= 0.4550
369 asigo1( 6)= 0.5470
370 asigo1( 7)= 0.6180
371 asigo1( 8)= 0.6690
372 asigo1( 9)= 0.7185
373 asigo1(10)= 0.7585
374 asigo1(11)= 0.7890
375 asigo1(12)= 0.8190
376 asigo1(13)= 0.8480
377 asigo1(14)= 0.8755
378 asigo1(15)= 0.9015
379 asigo1(16)= 0.9260
380 asigo1(17)= 0.9490
381 asigo1(18)= 0.9650
382 asigo1(19)= 0.9745
383 asigo1(20)= 0.9835
384 asigo1(21)= 0.9915
385 asigo1(22)= 0.9975
386!
387 sigo1( 1)= 0.0
388 sigo1( 2)= 0.1060
389 sigo1( 3)= 0.2100
390 sigo1( 4)= 0.3110
391 sigo1( 5)= 0.4080
392 sigo1( 6)= 0.5020
393 sigo1( 7)= 0.5920
394 sigo1( 8)= 0.6440
395 sigo1( 9)= 0.6940
396 sigo1(10)= 0.7430
397 sigo1(11)= 0.7740
398 sigo1(12)= 0.8040
399 sigo1(13)= 0.8340
400 sigo1(14)= 0.8620
401 sigo1(15)= 0.8890
402 sigo1(16)= 0.9140
403 sigo1(17)= 0.9380
404 sigo1(18)= 0.9600
405 sigo1(19)= 0.9700
406 sigo1(20)= 0.9790
407 sigo1(21)= 0.9880
408 sigo1(22)= 0.9950
409 sigo1(23)= 1.0
410 END IF
411!
412!
413 do j=1, nlevel
414 DO i=1,lsig1
415 if(abs(param%level(j)-asigo1(i)*10000)<small1) then
416 lvls(i,ifld)=1
417 lvlsxml(i,ifld)=j
418 irec=irec+1
419 endif
420 enddo
421 enddo
422 return
423!
424 ENDIF
425!
426!other
427 if(nlevel==1) then
428 lvls(1,ifld)=1
429 lvlsxml(1,ifld)=1
430 irec=irec+1
431 return
432 endif
433!
434
435 end