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