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