UPP  001
 All Data Structures Files Functions Pages
MISCLN.f
Go to the documentation of this file.
1 
2 !
85  SUBROUTINE miscln
86 
87 !
88 !
89  use vrbls3d, only: pmid, uh, vh, t, zmid, zint, pint, alpint, q, omga
90  use vrbls3d, only: catedr,mwt,gtg
91  use vrbls2d, only: pblh, cprate, fis, t500, t700, z500, z700,&
92  teql,ieql
93  use masks, only: lmh
94  use params_mod, only: d00, d50, h99999, h100, h1, h1m12, pq0, a2, a3, a4, &
95  rhmin, rgamog, tfrz, small, g
96  use ctlblk_mod, only: grib, cfld, fld_info, datapd, im, jsta, jend, jm, jsta_m, jend_m, &
97  nbnd, nbin_du, lm, htfd, spval, pthresh, nfd, petabnd, me,&
98  jsta_2l, jend_2u, modelname, submodelname
99  use rqstfld_mod, only: iget, lvls, id, iavblfld, lvlsxml
100  use grib2_module, only: pset
101  use upp_physics, only: fpvsnew,calrh_pw,calcape,calcape2,tvirtual
102  use gridspec_mod, only: gridtype
103 !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104  implicit none
105 !
106 ! SET LOCAL PARAMETERS. MAKE SURE NFD AND NBND AGREE
107 ! WITH THE VALUES SET IN SUBROUTINES FDLVL AND BNDLYR,
108 ! RESPECTIVELY.
109  real,PARAMETER :: c2k=273.15
110  real,parameter :: con_rd =2.8705e+2 ! gas constant air (J/kg/K)
111  real,parameter :: con_rv =4.6150e+2 ! gas constant H2O
112  real,parameter :: con_eps =con_rd/con_rv
113  real,parameter :: con_epsm1 =con_rd/con_rv-1
114  real,parameter :: cpthresh =0.000004
115  real,PARAMETER :: d1000=1000
116  real,PARAMETER :: d1500=1500
117  real,PARAMETER :: d2000=2000
118  real,PARAMETER :: hconst=42000000.
119  real,PARAMETER :: k2c=273.16
120  REAL,PARAMETER :: dm9999=-9999.0
121 
122 !
123 ! DECLARE VARIABLES.
124 !
125  LOGICAL north, field1,field2
126  LOGICAL, dimension(IM,JSTA:JEND) :: done, done1
127 
128  INTEGER, allocatable :: lvlbnd(:,:,:),lb2(:,:)
129 ! INTEGER LVLBND(IM,JM,NBND),LB2(IM,JM),LPBL(IM,JM)
130 
131  real,dimension(im,jm) :: grid1, grid2
132  real,dimension(im,jsta:jend) :: p1d, t1d, q1d, u1d, v1d, shr1d, z1d, &
133  rh1d, egrid1, egrid2, egrid3, egrid4, &
134  egrid5, egrid6, egrid7, egrid8, &
135  mlcape,mlcin,mllcl,mucape,mucin,mumixr, &
136  freezelvl,muq1d,slcl,the,maxthe
137  integer,dimension(im,jsta:jend) :: maxthepos
138  real, dimension(:,:,:),allocatable :: omgbnd, pwtbnd, qcnvbnd, &
139  pbnd, tbnd, qbnd, &
140  ubnd, vbnd, rhbnd, &
141  wbnd, t7d, q7d, &
142  u7d, v6d, p7d, &
143  icingfd,gtgfd,catfd,mwtfd
144  real, dimension(:,:,:,:),allocatable :: aerfd
145 
146  real, dimension(:,:),allocatable :: qm8510, rh4710, rh8498, &
147  rh4796, rh1847, ust, vst, &
148  rh3310, rh6610, rh3366, &
149  pw3310, rh4410, rh7294, &
150  rh4472, &
151  t78483, t89671, p78483, p89671
152 
153  REAL, dimension(:,:,:),allocatable :: heli
154  real, dimension(:,:), allocatable :: ushr1, vshr1, ushr6, vshr6, &
155  maxwp, maxwz, maxwu, maxwv, &
156  maxwt
157  INTEGER,dimension(:,:),allocatable :: llow, lupp
158  REAL, dimension(:,:),allocatable :: cangle,eshr,uvect,vvect,&
159  effust,effvst,fshr,htsfc,&
160  esrh
161 !
162  integer i,j,jj,l,itype,isvalue,lbnd,ilvl,ifd,itypefdlvl(nfd), &
163  iget1, iget2, iget3, llmh,imax,jmax,lmax
164  real dpbnd,pkl1,pku1,fac1,fac2,pl,tl,ql,qsat,rhl,tvrl,tvrblo, &
165  es1,es2,qs1,qs2,rh1,rh2,zsf,depth(2),work1,work2,work3, &
166  scintmp,mucapetmp,mucintmp,mllcltmp,eshrtmp,mlcapetmp,stp,&
167  fshrtmp,mlcintmp,slcltmp,lapse,ship
168 
169 ! Variables introduced to allow FD levels from control file - Y Mao
170  integer :: n,nfdctl
171  REAL, allocatable :: htfdctl(:)
172  integer, allocatable :: itypefdlvlctl(:)
173  integer ie,iw,jn,js,ive(jm),ivw(jm),jvn,jvs
174  integer istart,istop,jstart,jstop,midcal
175  real dummy(im,jsta:jend)
176  integer idummy(im,jsta:jend)
177 ! NEW VARIABLES USED FOR EFFECTIVE LAYER
178  INTEGER,dimension(:,:),allocatable :: el_base, el_tops
179  LOGICAL,dimension(:,:),allocatable :: found_base, found_tops
180  INTEGER,dimension(:,:),allocatable :: l_thetae_max
181  INTEGER,dimension(:,:),allocatable :: cape9, cins9
182  CHARACTER(LEN=5) :: im_ch, jsta_ch, jend_ch, me_ch
183  CHARACTER(LEN=60) :: effl_fname
184  CHARACTER(LEN=60) :: effl_fname2
185  INTEGER :: irec, iunit
186  INTEGER :: irec2, iunit2
187  LOGICAL :: debugprint
188  INTEGER :: lll
189  INTEGER :: llcl_par, leql_par
190  REAL :: lmask, psfc, cape_par, cins_par, lpar0
191  REAL, DIMENSION(4) :: parcel0
192  REAL, DIMENSION(:), ALLOCATABLE :: tpar_b, tpar_t
193  REAL, DIMENSION(:), ALLOCATABLE :: tpar_tmp
194  REAL, DIMENSION(:), ALLOCATABLE :: p_amb, t_amb, q_amb, zint_amb
195  REAL, DIMENSION(:,:,:), ALLOCATABLE :: tpar_base, tpar_tops
196 
197 !
198 !****************************************************************************
199 ! START MISCLN HERE.
200 !
201  debugprint = .false.
202 
203 
204  allocate(ushr1(im,jsta_2l:jend_2u),vshr1(im,jsta_2l:jend_2u), &
205  ushr6(im,jsta_2l:jend_2u),vshr6(im,jsta_2l:jend_2u))
206  allocate(ust(im,jsta_2l:jend_2u),vst(im,jsta_2l:jend_2u), &
207  heli(im,jsta_2l:jend_2u,2),fshr(im,jsta_2l:jend_2u))
208 !
209 ! HELICITY AND STORM MOTION.
210  iget1 = iget(162)
211  iget2 = -1
212  iget3 = -1
213  if (iget1 > 0) then
214  iget2 = lvls(1,iget1)
215  iget3 = lvls(2,iget1)
216  endif
217  IF (iget1 > 0 .OR. iget(163) > 0 .OR. iget(164) > 0) THEN
218  depth(1) = 3000.0
219  depth(2) = 1000.0
220  CALL calhel(depth,ust,vst,heli,ushr1,vshr1,ushr6,vshr6)
221  IF (iget2 > 0) then
222 !$omp parallel do private(i,j)
223  DO j=jsta,jend
224  DO i=1,im
225  grid1(i,j) = heli(i,j,1)
226  ENDDO
227  ENDDO
228  if(grib=='grib2') then
229  cfld=cfld+1
230  fld_info(cfld)%ifld=iavblfld(iget1)
231  fld_info(cfld)%lvl=lvlsxml(1,iget1)
232 !$omp parallel do private(i,j,jj)
233  do j=1,jend-jsta+1
234  jj = jsta+j-1
235  do i=1,im
236  datapd(i,j,cfld) = grid1(i,jj)
237  enddo
238  enddo
239  endif
240  ENDIF
241 
242  IF (iget3 > 0) then
243 !$omp parallel do private(i,j)
244  DO j=jsta,jend
245  DO i=1,im
246  grid1(i,j) = heli(i,j,2)
247  ENDDO
248  ENDDO
249  if(grib=='grib2') then
250  cfld=cfld+1
251  fld_info(cfld)%ifld=iavblfld(iget1)
252  fld_info(cfld)%lvl=lvlsxml(2,iget1)
253 !$omp parallel do private(i,j,jj)
254  do j=1,jend-jsta+1
255  jj = jsta+j-1
256  do i=1,im
257  datapd(i,j,cfld) = grid1(i,jj)
258  enddo
259  enddo
260  endif
261  ENDIF
262 
263  IF (iget(163) > 0) THEN
264 !$omp parallel do private(i,j)
265  DO j=jsta,jend
266  DO i=1,im
267  grid1(i,j) = ust(i,j)
268  ENDDO
269  ENDDO
270  if(grib=='grib2') then
271  cfld=cfld+1
272  fld_info(cfld)%ifld=iavblfld(iget(163))
273 !$omp parallel do private(i,j,jj)
274  do j=1,jend-jsta+1
275  jj = jsta+j-1
276  do i=1,im
277  datapd(i,j,cfld) = grid1(i,jj)
278  enddo
279  enddo
280  endif
281  ENDIF
282  IF (iget(164) > 0) THEN
283 !$omp parallel do private(i,j)
284  DO j=jsta,jend
285  DO i=1,im
286  grid1(i,j) = vst(i,j)
287  ENDDO
288  ENDDO
289  if(grib=='grib2') then
290  cfld=cfld+1
291  fld_info(cfld)%ifld=iavblfld(iget(164))
292 !$omp parallel do private(i,j,jj)
293  do j=1,jend-jsta+1
294  jj = jsta+j-1
295  do i=1,im
296  datapd(i,j,cfld) = grid1(i,jj)
297  enddo
298  enddo
299  endif
300  ENDIF
301  ENDIF
302 
303 ! UPDRAFT HELICITY
304 
305  if (iget(427) > 0) THEN
306  CALL calupdhel(grid1(1,jsta_2l))
307  if(grib=='grib2') then
308  cfld=cfld+1
309  fld_info(cfld)%ifld=iavblfld(iget(427))
310 !$omp parallel do private(i,j,jj)
311  do j=1,jend-jsta+1
312  jj = jsta+j-1
313  do i=1,im
314  datapd(i,j,cfld) = grid1(i,jj)
315  enddo
316  enddo
317  endif
318 
319  ENDIF
320 
321 ! CRA 0-1 KM AND 0-6 KM SHEAR
322 
323  IF(iget(430) > 0 .OR. iget(431) > 0 .OR. iget(432) > 0 &
324  .OR. iget(433) > 0) THEN
325 
326  depth = 6000.0
327  CALL calhel(depth,ust,vst,heli,ushr1,vshr1,ushr6,vshr6)
328 ! 0-6 km shear magnitude
329 !$omp parallel do private(i,j)
330  DO j=jsta,jend
331  DO i=1,im
332  fshr(i,j) = sqrt(ushr6(i,j)**2+vshr6(i,j)**2)
333  ENDDO
334  ENDDO
335  IF(iget(430) > 0) THEN
336 !$omp parallel do private(i,j)
337  DO j=jsta,jend
338  DO i=1,im
339  grid1(i,j) = ushr1(i,j)
340  ENDDO
341  ENDDO
342  if(grib=='grib2') then
343  cfld=cfld+1
344  fld_info(cfld)%ifld=iavblfld(iget(430))
345 !$omp parallel do private(i,j,jj)
346  do j=1,jend-jsta+1
347  jj = jsta+j-1
348  do i=1,im
349  datapd(i,j,cfld) = grid1(i,jj)
350  enddo
351  enddo
352  endif
353  ENDIF
354  IF(iget(431) > 0) THEN
355 !$omp parallel do private(i,j)
356  DO j=jsta,jend
357  DO i=1,im
358  grid1(i,j) = vshr1(i,j)
359  ENDDO
360  ENDDO
361  if(grib=='grib2') then
362  cfld=cfld+1
363  fld_info(cfld)%ifld=iavblfld(iget(431))
364 !$omp parallel do private(i,j,jj)
365  do j=1,jend-jsta+1
366  jj = jsta+j-1
367  do i=1,im
368  datapd(i,j,cfld) = grid1(i,jj)
369  enddo
370  enddo
371  endif
372  ENDIF
373  IF(iget(432) > 0) THEN
374 !$omp parallel do private(i,j)
375  DO j=jsta,jend
376  DO i=1,im
377  grid1(i,j) = ushr6(i,j)
378  ENDDO
379  ENDDO
380  if(grib=='grib2') then
381  cfld=cfld+1
382  fld_info(cfld)%ifld=iavblfld(iget(432))
383 !$omp parallel do private(i,j,jj)
384  do j=1,jend-jsta+1
385  jj = jsta+j-1
386  do i=1,im
387  datapd(i,j,cfld) = grid1(i,jj)
388  enddo
389  enddo
390  endif
391  ENDIF
392  IF(iget(433) > 0) THEN
393 !$omp parallel do private(i,j)
394  DO j=jsta,jend
395  DO i=1,im
396  grid1(i,j) = vshr6(i,j)
397  ENDDO
398  ENDDO
399  if(grib=='grib2') then
400  cfld=cfld+1
401  fld_info(cfld)%ifld=iavblfld(iget(433))
402 !$omp parallel do private(i,j,jj)
403  do j=1,jend-jsta+1
404  jj = jsta+j-1
405  do i=1,im
406  datapd(i,j,cfld) = grid1(i,jj)
407  enddo
408  enddo
409  endif
410  ENDIF
411  ENDIF
412 !
413  if (allocated(ushr1)) deallocate(ushr1)
414  if (allocated(vshr1)) deallocate(vshr1)
415  if (allocated(ushr6)) deallocate(ushr6)
416  if (allocated(vshr6)) deallocate(vshr6)
417  if (allocated(ust)) deallocate(ust)
418  if (allocated(vst)) deallocate(vst)
419  if (allocated(heli)) deallocate(heli)
420 ! CRA
421 !
422 !
423 ! ***BLOCK 1: TROPOPAUSE P, Z, T, U, V, AND WIND SHEAR.
424 !
425  IF ((iget(054)>0).OR.(iget(055)>0).OR. &
426  (iget(056)>0).OR.(iget(057)>0).OR. &
427  (iget(177)>0).OR. &
428  (iget(058)>0).OR.(iget(108)>0) ) THEN
429 ! Chuang: Use GFS algorithm per Iredell's and DiMego's decision on unification
430 !$omp parallel do private(i,j)
431  DO j=jsta,jend
432  DO i=1,im
433 
434  if(pmid(i,j,1)<spval) then
435 ! INPUT
436  CALL tpause(lm,pmid(i,j,1:lm),uh(i,j,1:lm) &
437 ! INPUT
438  ,vh(i,j,1:lm),t(i,j,1:lm),zmid(i,j,1:lm) &
439 ! OUTPUT
440  ,p1d(i,j),u1d(i,j),v1d(i,j),t1d(i,j) &
441 ! OUTPUT
442  ,z1d(i,j),shr1d(i,j)) ! OUTPUT
443  else
444  p1d(i,j) = spval
445  u1d(i,j) = spval
446  v1d(i,j) = spval
447  t1d(i,j) = spval
448  z1d(i,j) = spval
449  shr1d(i,j) = spval
450  endif
451 
452  END DO
453  END DO
454 !
455 ! TROPOPAUSE PRESSURE.
456  IF (iget(054) > 0) THEN
457 !$omp parallel do private(i,j)
458  DO j=jsta,jend
459  DO i=1,im
460  grid1(i,j) = p1d(i,j)
461  ENDDO
462  ENDDO
463  if(grib=='grib2') then
464  cfld=cfld+1
465  fld_info(cfld)%ifld=iavblfld(iget(054))
466 !$omp parallel do private(i,j,jj)
467  do j=1,jend-jsta+1
468  jj = jsta+j-1
469  do i=1,im
470  datapd(i,j,cfld) = grid1(i,jj)
471  enddo
472  enddo
473  endif
474  ENDIF
475 
476 ! ICAO HEIGHT OF TROPOPAUSE
477  IF (iget(399)>0) THEN
478  CALL icaoheight(p1d, grid1(1,jsta))
479 ! print*,'sample TROPOPAUSE ICAO HEIGHTS',GRID1(im/2,(jsta+jend)/2)
480  if(grib=='grib2') then
481  cfld=cfld+1
482  fld_info(cfld)%ifld=iavblfld(iget(399))
483 !$omp parallel do private(i,j,jj)
484  do j=1,jend-jsta+1
485  jj = jsta+j-1
486  do i=1,im
487  datapd(i,j,cfld) = grid1(i,jj)
488  enddo
489  enddo
490  endif
491  ENDIF
492 
493 ! TROPOPAUSE HEIGHT.
494  IF (iget(177) > 0) THEN
495 !$omp parallel do private(i,j)
496  DO j=jsta,jend
497  DO i=1,im
498  grid1(i,j) = z1d(i,j)
499  ENDDO
500  ENDDO
501  if(grib=='grib2') then
502  cfld=cfld+1
503  fld_info(cfld)%ifld=iavblfld(iget(177))
504 !$omp parallel do private(i,j,jj)
505  do j=1,jend-jsta+1
506  jj = jsta+j-1
507  do i=1,im
508  datapd(i,j,cfld) = grid1(i,jj)
509  enddo
510  enddo
511  endif
512  ENDIF
513 !
514 ! TROPOPAUSE TEMPERATURE.
515  IF (iget(055) > 0) THEN
516 !$omp parallel do private(i,j)
517  DO j=jsta,jend
518  DO i=1,im
519  grid1(i,j) = t1d(i,j)
520  ENDDO
521  ENDDO
522  if(grib=='grib2') then
523  cfld=cfld+1
524  fld_info(cfld)%ifld=iavblfld(iget(055))
525 !$omp parallel do private(i,j,jj)
526  do j=1,jend-jsta+1
527  jj = jsta+j-1
528  do i=1,im
529  datapd(i,j,cfld) = grid1(i,jj)
530  enddo
531  enddo
532  endif
533  ENDIF
534 !
535 ! TROPOPAUSE POTENTIAL TEMPERATURE.
536  IF (iget(108) > 0) THEN
537  CALL calpot(p1d,t1d,grid1(1,jsta))
538  if(grib=='grib2') then
539  cfld=cfld+1
540  fld_info(cfld)%ifld=iavblfld(iget(108))
541 !$omp parallel do private(i,j,jj)
542  do j=1,jend-jsta+1
543  jj = jsta+j-1
544  do i=1,im
545  datapd(i,j,cfld) = grid1(i,jj)
546  enddo
547  enddo
548  endif
549  ENDIF
550 !
551 ! TROPOPAUSE U WIND AND/OR V WIND.
552  IF ((iget(056) > 0).OR.(iget(057) > 0)) THEN
553 !$omp parallel do private(i,j)
554  DO j=jsta,jend
555  DO i=1,im
556  grid1(i,j)=u1d(i,j)
557  grid2(i,j)=v1d(i,j)
558  ENDDO
559  ENDDO
560  if(grib=='grib2') then
561  if(iget(056)>0) then
562  cfld=cfld+1
563  fld_info(cfld)%ifld=iavblfld(iget(056))
564 !$omp parallel do private(i,j,jj)
565  do j=1,jend-jsta+1
566  jj = jsta+j-1
567  do i=1,im
568  datapd(i,j,cfld) = grid1(i,jj)
569  enddo
570  enddo
571  endif
572  if(iget(057)>0) then
573  cfld=cfld+1
574  fld_info(cfld)%ifld=iavblfld(iget(057))
575 !$omp parallel do private(i,j,jj)
576  do j=1,jend-jsta+1
577  jj = jsta+j-1
578  do i=1,im
579  datapd(i,j,cfld) = grid2(i,jj)
580  enddo
581  enddo
582  endif
583  endif
584  ENDIF
585 !
586 ! TROPOPAUSE WIND SHEAR.
587  IF (iget(058) > 0) THEN
588 !$omp parallel do private(i,j)
589  DO j=jsta,jend
590  DO i=1,im
591  grid1(i,j) = shr1d(i,j)
592  ENDDO
593  ENDDO
594  if(grib=='grib2') then
595  cfld=cfld+1
596  fld_info(cfld)%ifld=iavblfld(iget(058))
597 !$omp parallel do private(i,j,jj)
598  do j=1,jend-jsta+1
599  jj = jsta+j-1
600  do i=1,im
601  datapd(i,j,cfld) = grid1(i,jj)
602  enddo
603  enddo
604  endif
605  ENDIF
606  ENDIF
607 !
608 !
609 ! ***BLOCK 2: MAX WIND LEVEL P, Z, U, AND V
610 !
611 ! MAX WIND LEVEL CALCULATIONS
612  IF ((iget(173)>0) .OR. (iget(174)>0) .OR. &
613  (iget(175)>0) .OR. (iget(176)>0)) THEN
614 
615  allocate(maxwp(im,jsta:jend), maxwz(im,jsta:jend), &
616  maxwu(im,jsta:jend), maxwv(im,jsta:jend),maxwt(im,jsta:jend))
617 !$omp parallel do private(i,j)
618  DO j=jsta,jend
619  DO i=1,im
620  maxwp(i,j)=spval
621  maxwz(i,j)=spval
622  maxwu(i,j)=spval
623  maxwv(i,j)=spval
624  ENDDO
625  ENDDO
626 
627 ! CALL CALMXW(MAXWP,MAXWZ,MAXWU,MAXWV,MAXWT)
628 ! Chuang: Use GFS algorithm per Iredell's and DiMego's decision on unification
629 !$omp parallel do private(i,j)
630  DO j=jsta,jend
631  loopi:DO i=1,im
632  DO l=1,lm
633  IF (abs(pmid(i,j,l)-spval)<=small .OR. &
634  abs(uh(i,j,l)-spval)<=small .OR. &
635  abs(uh(i,j,l)-spval)<=small .OR. &
636  abs(vh(i,j,l)-spval)<=small .OR. &
637  abs(t(i,j,l)-spval)<=small .OR. &
638  abs(zmid(i,j,l)-spval)<=small) cycle loopi
639  ENDDO
640 ! INPUT
641  CALL mxwind(lm,pmid(i,j,1:lm),uh(i,j,1:lm) &
642 ! INPUT
643  ,vh(i,j,1:lm),t(i,j,1:lm),zmid(i,j,1:lm) &
644 ! OUTPUT
645  ,maxwp(i,j),maxwu(i,j),maxwv(i,j) &
646 ! OUTPUT
647  ,maxwt(i,j),maxwz(i,j))
648  ENDDO loopi
649  END DO
650 ! PRESSURE OF MAX WIND LEVEL
651  IF (iget(173) > 0) THEN
652 !$omp parallel do private(i,j)
653  DO j=jsta,jend
654  DO i=1,im
655  grid1(i,j) = maxwp(i,j)
656  ENDDO
657  ENDDO
658  if(grib=='grib2') then
659  cfld=cfld+1
660  fld_info(cfld)%ifld=iavblfld(iget(173))
661 !$omp parallel do private(i,j,jj)
662  do j=1,jend-jsta+1
663  jj = jsta+j-1
664  do i=1,im
665  datapd(i,j,cfld) = grid1(i,jj)
666  enddo
667  enddo
668  endif
669  ENDIF
670 ! ICAO HEIGHT OF MAX WIND LEVEL
671  IF (iget(398)>0) THEN
672  CALL icaoheight(maxwp, grid1(1,jsta))
673 ! print*,'sample MAX WIND ICAO HEIGHTS',GRID1(im/2,(jsta+jend)/2)
674  if(grib=='grib2') then
675  cfld=cfld+1
676  fld_info(cfld)%ifld=iavblfld(iget(398))
677 !$omp parallel do private(i,j,jj)
678  do j=1,jend-jsta+1
679  jj = jsta+j-1
680  do i=1,im
681  datapd(i,j,cfld) = grid1(i,jj)
682  enddo
683  enddo
684  endif
685  ENDIF
686 ! HEIGHT OF MAX WIND LEVEL
687  IF (iget(174) > 0) THEN
688 !$omp parallel do private(i,j)
689  DO j=jsta,jend
690  DO i=1,im
691  grid1(i,j) = maxwz(i,j)
692  ENDDO
693  ENDDO
694  if(grib=='grib2') then
695  cfld=cfld+1
696  fld_info(cfld)%ifld=iavblfld(iget(174))
697 !$omp parallel do private(i,j,jj)
698  do j=1,jend-jsta+1
699  jj = jsta+j-1
700  do i=1,im
701  datapd(i,j,cfld) = grid1(i,jj)
702  enddo
703  enddo
704  endif
705  ENDIF
706 
707 ! MAX WIND LEVEL U WIND AND/OR V WIND.
708  IF ((iget(175) > 0).OR.(iget(176) > 0)) THEN
709 !$omp parallel do private(i,j)
710  DO j=jsta,jend
711  DO i=1,im
712  grid1(i,j) = maxwu(i,j)
713  grid2(i,j) = maxwv(i,j)
714  ENDDO
715  ENDDO
716  if(grib=='grib2') then
717  cfld=cfld+1
718  fld_info(cfld)%ifld=iavblfld(iget(175))
719 !$omp parallel do private(i,j,jj)
720  do j=1,jend-jsta+1
721  jj = jsta+j-1
722  do i=1,im
723  datapd(i,j,cfld) = grid1(i,jj)
724  enddo
725  enddo
726  cfld=cfld+1
727  fld_info(cfld)%ifld=iavblfld(iget(176))
728 !$omp parallel do private(i,j,jj)
729  do j=1,jend-jsta+1
730  jj = jsta+j-1
731  do i=1,im
732  datapd(i,j,cfld) = grid2(i,jj)
733  enddo
734  enddo
735  endif
736  ENDIF
737 ! TEMPERATURE OF MAX WIND LEVEL
738  IF (iget(314) > 0) THEN
739 !$omp parallel do private(i,j)
740  DO j=jsta,jend
741  DO i=1,im
742  grid1(i,j)=maxwt(i,j)
743  ENDDO
744  ENDDO
745  if(grib=='grib2') then
746  cfld=cfld+1
747  fld_info(cfld)%ifld=iavblfld(iget(314))
748 !$omp parallel do private(i,j,jj)
749  do j=1,jend-jsta+1
750  jj = jsta+j-1
751  do i=1,im
752  datapd(i,j,cfld) = grid1(i,jj)
753  enddo
754  enddo
755  endif
756  ENDIF
757  deallocate(maxwp,maxwz,maxwu,maxwv,maxwt)
758  ENDIF
759 
760 !
761 ! ***BLOCK 3-1: FD LEVEL (selected) T, Q, U, AND V.
762 !
763  IF ( (iget(059)>0.or.iget(586)>0).OR.iget(911)>0.OR. &
764  (iget(060)>0.or.iget(576)>0).OR. &
765  (iget(061)>0.or.iget(577)>0).OR. &
766  (iget(601)>0.or.iget(602)>0.or.iget(603)>0).OR. &
767  (iget(604)>0.or.iget(605)>0).OR. &
768  (iget(451)>0.or.iget(578)>0).OR.iget(580)>0 ) THEN
769 
770  ALLOCATE(t7d(im,jsta:jend,nfd), q7d(im,jsta:jend,nfd), &
771  u7d(im,jsta:jend,nfd), v6d(im,jsta:jend,nfd), &
772  p7d(im,jsta:jend,nfd), icingfd(im,jsta:jend,nfd) &
773  ,aerfd(im,jsta:jend,nfd,nbin_du))
774 
775 !
776 ! DETERMINE WHETHER TO DO MSL OR AGL FD LEVELS
777 !
778  itypefdlvl=1
779  DO ifd = 1,nfd
780  IF (iget(059)>0) THEN
781  IF (lvls(ifd,iget(059))>1) itypefdlvl(ifd)=2
782  ENDIF
783  IF (iget(911)>0) THEN
784  IF (lvls(ifd,iget(911))>1) itypefdlvl(ifd)=2
785  ENDIF
786 !for grib2, spec hgt only
787  IF (iget(586)>0) THEN
788  IF(lvls(ifd,iget(586))>0) itypefdlvl(ifd)=2
789  ENDIF
790  IF (iget(060)>0) THEN
791  IF (lvls(ifd,iget(060))>1) itypefdlvl(ifd)=2
792  ENDIF
793  IF (iget(576)>0) THEN
794  IF(lvls(ifd,iget(576))>0) itypefdlvl(ifd)=2
795  ENDIF
796  IF (iget(061)>0) THEN
797  IF (lvls(ifd,iget(061))>1) itypefdlvl(ifd)=2
798  ENDIF
799  IF (iget(577)>0) then
800  if(lvls(ifd,iget(577))>0) itypefdlvl(ifd)=2
801  ENDIF
802  IF (iget(451)>0) THEN
803  IF (lvls(ifd,iget(451))>1) itypefdlvl(ifd)=2
804  ENDIF
805  IF (iget(578)>0) then
806  if(lvls(ifd,iget(578))>0) itypefdlvl(ifd)=2
807  ENDIF
808 
809  IF (iget(580)>0) then
810  if(lvls(ifd,iget(580))>1) itypefdlvl(ifd)=2
811  ENDIF
812  IF (iget(587)>0) then
813  if(lvls(ifd,iget(587))>0) itypefdlvl(ifd)=2
814  ENDIF
815 
816  IF (iget(601)>0) THEN
817  IF (lvls(ifd,iget(601))>1) itypefdlvl(ifd)=2
818  ENDIF
819  IF (iget(602)>0) THEN
820  IF (lvls(ifd,iget(602))>1) itypefdlvl(ifd)=2
821  ENDIF
822  IF (iget(603)>0) THEN
823  IF (lvls(ifd,iget(603))>1) itypefdlvl(ifd)=2
824  ENDIF
825  IF (iget(604)>0) THEN
826  IF (lvls(ifd,iget(604))>1) itypefdlvl(ifd)=2
827  ENDIF
828  IF (iget(605)>0) THEN
829  IF (lvls(ifd,iget(605))>1) itypefdlvl(ifd)=2
830  ENDIF
831 
832  ENDDO
833 ! print *,'call FDLVL with ITYPEFDLVL: ', ITYPEFDLVL,'for tmp,lvls=',LVLS(1:15,iget(59)), &
834 ! 'grib2tmp lvs=',LVLS(1:15,iget(586))
835 
836  CALL fdlvl(itypefdlvl,t7d,q7d,u7d,v6d,p7d,icingfd,aerfd)
837 !
838  DO 10 ifd = 1,nfd
839 !
840 ! FD LEVEL TEMPERATURE.
841  iget1 = iget(059)
842  iget2 = iget(586)
843  if (iget1 > 0) then
844  work1 = lvls(ifd,iget1)
845  else
846  work1 = 0.0
847  endif
848  if (iget2 > 0) then
849  work2 = lvls(ifd,iget2)
850  else
851  work2 = 0.0
852  endif
853  IF (iget1 > 0 .or. iget2 > 0) THEN
854  IF (work1 > 0 .or. work2 > 0) THEN
855 
856 !$omp parallel do private(i,j)
857  DO j=jsta,jend
858  DO i=1,im
859  grid1(i,j) = t7d(i,j,ifd)
860  ENDDO
861  ENDDO
862  IF(work1 > 0) then
863  if(grib == 'grib2') then
864  cfld = cfld + 1
865  fld_info(cfld)%ifld = iavblfld(iget1)
866  fld_info(cfld)%lvl = lvlsxml(ifd,iget1)
867 !$omp parallel do private(i,j,jj)
868  do j=1,jend-jsta+1
869  jj = jsta+j-1
870  do i=1,im
871  datapd(i,j,cfld) = grid1(i,jj)
872  enddo
873  enddo
874  endif
875  ENDIF
876  IF (work2 > 0) THEN
877  if(grib == 'grib2') then
878  cfld = cfld + 1
879  fld_info(cfld)%ifld = iavblfld(iget2)
880  fld_info(cfld)%lvl = lvlsxml(ifd,iget2)
881 !$omp parallel do private(i,j,jj)
882  do j=1,jend-jsta+1
883  jj = jsta+j-1
884  do i=1,im
885  datapd(i,j,cfld) = grid1(i,jj)
886  enddo
887  enddo
888  endif
889  ENDIF
890  ENDIF
891  ENDIF
892 
893 ! FD LEVEL VIRTUAL TEMPERATURE.
894  IF (iget(911)>0) THEN
895  IF (lvls(ifd,iget(911))>0) THEN
896  DO j=jsta,jend
897  DO i=1,im
898  if ( t7d(i,j,ifd) > 600 ) then
899  grid1(i,j)=spval
900  else
901  grid1(i,j)=t7d(i,j,ifd)*(1.+0.608*q7d(i,j,ifd))
902  endif
903  !print *, "grid value ",T7D(I,J,IFD),Q7D(I,J,IFD),T7D(I,J,IFD)*(1.+0.608*Q7D(I,J,IFD)),GRID1(I,J)
904  ENDDO
905  ENDDO
906  IF(lvls(ifd,iget(911))>0) then
907  if(grib=='grib2') then
908  cfld=cfld+1
909  fld_info(cfld)%ifld=iavblfld(iget(911))
910  fld_info(cfld)%lvl=lvlsxml(ifd,iget(911))
911  datapd(1:im,1:jend-jsta+1,cfld)=grid1(1:im,jsta:jend)
912  endif
913  ENDIF
914  ENDIF
915  ENDIF
916 
917 !
918 ! FD LEVEL SPEC HUMIDITY.
919  iget1 = iget(451)
920  iget2 = iget(578)
921  if (iget1 > 0) then
922  work1 = lvls(ifd,iget1)
923  else
924  work1 = 0.0
925  endif
926  if (iget2 > 0) then
927  work2 = lvls(ifd,iget2)
928  else
929  work2 = 0.0
930  endif
931  IF (iget1 > 0 .or. iget2 > 0) THEN
932  IF (work1 > 0 .or. work2 > 0)THEN
933 !$omp parallel do private(i,j)
934  DO j=jsta,jend
935  DO i=1,im
936  grid1(i,j) = q7d(i,j,ifd)
937  ENDDO
938  ENDDO
939  if(work1 > 0) then
940  if(grib == 'grib2') then
941  cfld = cfld + 1
942  fld_info(cfld)%ifld = iavblfld(iget1)
943  fld_info(cfld)%lvl = lvlsxml(ifd,iget1)
944 !$omp parallel do private(i,j,jj)
945  do j=1,jend-jsta+1
946  jj = jsta+j-1
947  do i=1,im
948  datapd(i,j,cfld) = grid1(i,jj)
949  enddo
950  enddo
951  endif
952  endif
953  if(work2 > 0) then
954  if(grib == 'grib2') then
955  cfld = cfld + 1
956  fld_info(cfld)%ifld = iavblfld(iget2)
957  fld_info(cfld)%lvl = lvlsxml(ifd,iget2)
958 !$omp parallel do private(i,j,jj)
959  do j=1,jend-jsta+1
960  jj = jsta+j-1
961  do i=1,im
962  datapd(i,j,cfld) = grid1(i,jj)
963  enddo
964  enddo
965  endif
966  endif
967  ENDIF
968  ENDIF
969 !
970 ! FD LEVEL PRESSURE
971  iget1 = iget(482)
972  iget2 = iget(579)
973  if (iget1 > 0) then
974  work1 = lvls(ifd,iget1)
975  else
976  work1 = 0.0
977  endif
978  if (iget2 > 0) then
979  work2 = lvls(ifd,iget2)
980  else
981  work2 = 0.0
982  endif
983  IF (iget1 > 0 .or. iget2 > 0) THEN
984  IF (work1 > 0 .or. work2 > 0) THEN
985 !$omp parallel do private(i,j)
986  DO j=jsta,jend
987  DO i=1,im
988  grid1(i,j) = p7d(i,j,ifd)
989  ENDDO
990  ENDDO
991  if(work1 > 0) then
992  if(grib == 'grib2') then
993  cfld = cfld + 1
994  fld_info(cfld)%ifld = iavblfld(iget1)
995  fld_info(cfld)%lvl = lvlsxml(ifd,iget1)
996 !$omp parallel do private(i,j,jj)
997  do j=1,jend-jsta+1
998  jj = jsta+j-1
999  do i=1,im
1000  datapd(i,j,cfld) = grid1(i,jj)
1001  enddo
1002  enddo
1003  endif
1004  endif
1005  if(work2 > 0) then
1006  if(grib == 'grib2') then
1007  cfld = cfld + 1
1008  fld_info(cfld)%ifld = iavblfld(iget2)
1009  fld_info(cfld)%lvl = lvlsxml(ifd,iget2)
1010 !$omp parallel do private(i,j,jj)
1011  do j=1,jend-jsta+1
1012  jj = jsta+j-1
1013  do i=1,im
1014  datapd(i,j,cfld) = grid1(i,jj)
1015  enddo
1016  enddo
1017  endif
1018  endif
1019  ENDIF
1020  ENDIF
1021 !
1022 ! FD LEVEL ICING
1023  iget1 = iget(580)
1024  iget2 = iget(587)
1025  if (iget1 > 0) then
1026  work1 = lvls(ifd,iget1)
1027  else
1028  work1 = 0.0
1029  endif
1030  if (iget2 > 0) then
1031  work2 = lvls(ifd,iget2)
1032  else
1033  work2 = 0.0
1034  endif
1035  IF (iget1 > 0 .or. iget2 > 0) THEN
1036  IF (work1 > 0 .or. work2 > 0) THEN
1037 !$omp parallel do private(i,j)
1038  DO j=jsta,jend
1039  DO i=1,im
1040  grid1(i,j) = icingfd(i,j,ifd)
1041  ENDDO
1042  ENDDO
1043  if(iget1 > 0) then
1044  if(grib == 'grib2') then
1045  cfld = cfld + 1
1046  fld_info(cfld)%ifld = iavblfld(iget1)
1047  fld_info(cfld)%lvl = lvlsxml(ifd,iget1)
1048 !$omp parallel do private(i,j,jj)
1049  do j=1,jend-jsta+1
1050  jj = jsta+j-1
1051  do i=1,im
1052  datapd(i,j,cfld) = grid1(i,jj)
1053  enddo
1054  enddo
1055  endif
1056  endif
1057  if(work2 > 0) then
1058  if(grib == 'grib2') then
1059  cfld = cfld + 1
1060  fld_info(cfld)%ifld = iavblfld(iget2)
1061  fld_info(cfld)%lvl = lvlsxml(ifd,iget2)
1062 !$omp parallel do private(i,j,jj)
1063  do j=1,jend-jsta+1
1064  jj = jsta+j-1
1065  do i=1,im
1066  datapd(i,j,cfld) = grid1(i,jj)
1067  enddo
1068  enddo
1069  endif
1070  endif
1071 
1072  ENDIF
1073  ENDIF
1074 !
1075 ! ADD FD LEVEL DUST/ASH (GOCART)
1076  IF (iget(601)>0) THEN ! DUST 1
1077  IF (lvls(ifd,iget(601))>0) THEN
1078 !$omp parallel do private(i,j)
1079  DO j=jsta,jend
1080  DO i=1,im
1081  grid1(i,j)=aerfd(i,j,ifd,1)
1082  ENDDO
1083  ENDDO
1084  if(iget(601)>0) then
1085  if(grib=='grib2') then
1086  cfld=cfld+1
1087  fld_info(cfld)%ifld=iavblfld(iget(601))
1088  fld_info(cfld)%lvl=lvlsxml(ifd,iget(601))
1089 !$omp parallel do private(i,j,jj)
1090  do j=1,jend-jsta+1
1091  jj = jsta+j-1
1092  do i=1,im
1093  datapd(i,j,cfld) = grid1(i,jj)
1094  enddo
1095  enddo
1096  endif
1097  endif
1098  ENDIF
1099  ENDIF
1100 
1101  IF (iget(602)>0) THEN ! DUST 2
1102  IF (lvls(ifd,iget(602))>0) THEN
1103 !$omp parallel do private(i,j)
1104  DO j=jsta,jend
1105  DO i=1,im
1106  grid1(i,j)=aerfd(i,j,ifd,2)
1107  ENDDO
1108  ENDDO
1109  if(iget(602)>0) then
1110  if(grib=='grib2') then
1111  cfld=cfld+1
1112  fld_info(cfld)%ifld=iavblfld(iget(602))
1113  fld_info(cfld)%lvl=lvlsxml(ifd,iget(602))
1114 !$omp parallel do private(i,j,jj)
1115  do j=1,jend-jsta+1
1116  jj = jsta+j-1
1117  do i=1,im
1118  datapd(i,j,cfld) = grid1(i,jj)
1119  enddo
1120  enddo
1121  endif
1122  endif
1123  ENDIF
1124  ENDIF
1125 
1126  IF (iget(603)>0) THEN ! DUST 3
1127  IF (lvls(ifd,iget(603))>0) THEN
1128 !$omp parallel do private(i,j)
1129  DO j=jsta,jend
1130  DO i=1,im
1131  grid1(i,j)=aerfd(i,j,ifd,3)
1132  ENDDO
1133  ENDDO
1134  if(iget(603)>0) then
1135  if(grib=='grib2') then
1136  cfld=cfld+1
1137  fld_info(cfld)%ifld=iavblfld(iget(603))
1138  fld_info(cfld)%lvl=lvlsxml(ifd,iget(603))
1139 !$omp parallel do private(i,j,jj)
1140  do j=1,jend-jsta+1
1141  jj = jsta+j-1
1142  do i=1,im
1143  datapd(i,j,cfld) = grid1(i,jj)
1144  enddo
1145  enddo
1146  endif
1147  endif
1148  ENDIF
1149  ENDIF
1150 
1151  IF (iget(604)>0) THEN ! DUST 4
1152  IF (lvls(ifd,iget(604))>0) THEN
1153 !$omp parallel do private(i,j)
1154  DO j=jsta,jend
1155  DO i=1,im
1156  grid1(i,j)=aerfd(i,j,ifd,4)
1157  ENDDO
1158  ENDDO
1159  if(iget(604)>0) then
1160  if(grib=='grib2') then
1161  cfld=cfld+1
1162  fld_info(cfld)%ifld=iavblfld(iget(604))
1163  fld_info(cfld)%lvl=lvlsxml(ifd,iget(604))
1164 !$omp parallel do private(i,j,jj)
1165  do j=1,jend-jsta+1
1166  jj = jsta+j-1
1167  do i=1,im
1168  datapd(i,j,cfld) = grid1(i,jj)
1169  enddo
1170  enddo
1171  endif
1172  endif
1173  ENDIF
1174  ENDIF
1175 
1176  IF (iget(605)>0) THEN ! DUST 5
1177  IF (lvls(ifd,iget(605))>0) THEN
1178 !$omp parallel do private(i,j)
1179  DO j=jsta,jend
1180  DO i=1,im
1181  grid1(i,j)=aerfd(i,j,ifd,5)
1182  ENDDO
1183  ENDDO
1184  if(iget(605)>0) then
1185  if(grib=='grib2') then
1186  cfld=cfld+1
1187  fld_info(cfld)%ifld=iavblfld(iget(605))
1188  fld_info(cfld)%lvl=lvlsxml(ifd,iget(605))
1189 !$omp parallel do private(i,j,jj)
1190  do j=1,jend-jsta+1
1191  jj = jsta+j-1
1192  do i=1,im
1193  datapd(i,j,cfld) = grid1(i,jj)
1194  enddo
1195  enddo
1196  endif
1197  endif
1198  ENDIF
1199  ENDIF
1200 
1201 !
1202 !
1203 ! FD LEVEL U WIND AND/OR V WIND.
1204  IF ((iget(060)>0).OR.(iget(061)>0)) THEN
1205 !$omp parallel do private(i,j)
1206  DO j=jsta,jend
1207  DO i=1,im
1208  grid1(i,j)=u7d(i,j,ifd)
1209  grid2(i,j)=v6d(i,j,ifd)
1210  ENDDO
1211  ENDDO
1212  IF (iget(060)>0) THEN
1213  IF (lvls(ifd,iget(060))>0) then
1214  if(grib=='grib2') then
1215  cfld=cfld+1
1216  fld_info(cfld)%ifld=iavblfld(iget(060))
1217  fld_info(cfld)%lvl=lvlsxml(ifd,iget(060))
1218 !$omp parallel do private(i,j,jj)
1219  do j=1,jend-jsta+1
1220  jj = jsta+j-1
1221  do i=1,im
1222  datapd(i,j,cfld) = grid1(i,jj)
1223  enddo
1224  enddo
1225  endif
1226  ENDIF
1227  ENDIF
1228  IF (iget(061)>0) THEN
1229  IF (lvls(ifd,iget(061))>0) THEN
1230  if(grib=='grib2') then
1231  cfld=cfld+1
1232  fld_info(cfld)%ifld=iavblfld(iget(061))
1233  fld_info(cfld)%lvl=lvlsxml(ifd,iget(061))
1234 !$omp parallel do private(i,j,jj)
1235  do j=1,jend-jsta+1
1236  jj = jsta+j-1
1237  do i=1,im
1238  datapd(i,j,cfld) = grid2(i,jj)
1239  enddo
1240  enddo
1241  endif
1242  ENDIF
1243  ENDIF
1244  ENDIF
1245 !
1246 ! FD LEVEL U WIND AND/OR V WIND.
1247  IF ((iget(576)>0).OR.(iget(577)>0)) THEN
1248 !$omp parallel do private(i,j)
1249  DO j=jsta,jend
1250  DO i=1,im
1251  grid1(i,j) = u7d(i,j,ifd)
1252  grid2(i,j) = v6d(i,j,ifd)
1253  ENDDO
1254  ENDDO
1255  IF (iget(576)>0) THEN
1256  IF (lvls(ifd,iget(576))>0) then
1257  if(grib=='grib2') then
1258  cfld=cfld+1
1259  fld_info(cfld)%ifld=iavblfld(iget(576))
1260  fld_info(cfld)%lvl=lvlsxml(ifd,iget(576))
1261 !$omp parallel do private(i,j,jj)
1262  do j=1,jend-jsta+1
1263  jj = jsta+j-1
1264  do i=1,im
1265  datapd(i,j,cfld) = grid1(i,jj)
1266  enddo
1267  enddo
1268  endif
1269  ENDIF
1270  ENDIF
1271  IF (iget(577)>0) THEN
1272  IF (lvls(ifd,iget(577))>0) THEN
1273  if(grib=='grib2') then
1274  cfld=cfld+1
1275  fld_info(cfld)%ifld=iavblfld(iget(577))
1276  fld_info(cfld)%lvl=lvlsxml(ifd,iget(577))
1277 !$omp parallel do private(i,j,jj)
1278  do j=1,jend-jsta+1
1279  jj = jsta+j-1
1280  do i=1,im
1281  datapd(i,j,cfld) = grid2(i,jj)
1282  enddo
1283  enddo
1284  endif
1285  ENDIF
1286  ENDIF
1287  ENDIF
1288 
1289  10 CONTINUE
1290  DEALLOCATE(t7d,q7d,u7d,v6d,p7d,icingfd,aerfd)
1291  ENDIF
1292 
1293 !
1294 ! ***BLOCK 3-2: FD LEVEL (from control file) GTG
1295 !
1296  IF(iget(467)>0.or.iget(468)>0.or.iget(469)>0) THEN
1297  if(iget(467)>0) THEN ! GTG
1298  n=iavblfld(iget(467))
1299  nfdctl=size(pset%param(n)%level)
1300  if(allocated(itypefdlvlctl)) deallocate(itypefdlvlctl)
1301  allocate(itypefdlvlctl(nfdctl))
1302  DO ifd = 1,nfdctl
1303  itypefdlvlctl(ifd)=lvls(ifd,iget(467))
1304  enddo
1305  if(allocated(htfdctl)) deallocate(htfdctl)
1306  allocate(htfdctl(nfdctl))
1307  htfdctl=pset%param(n)%level
1308 ! print *, "GTG 467 levels=",pset%param(N)%level
1309  allocate(gtgfd(im,jsta:jend,nfdctl))
1310  call fdlvl_mass(itypefdlvlctl,nfdctl,htfdctl,gtg,gtgfd)
1311 ! print *, "GTG 467 Done GTGFD=",me,GTGFD(IM/2,jend,1:NFDCTL)
1312  DO ifd = 1,nfdctl
1313  IF (lvls(ifd,iget(467))>0) THEN
1314 !$omp parallel do private(i,j)
1315  DO j=jsta,jend
1316  DO i=1,im
1317  grid1(i,j)=gtgfd(i,j,ifd)
1318  ENDDO
1319  ENDDO
1320  if(grib=='grib2') then
1321  cfld=cfld+1
1322  fld_info(cfld)%ifld=iavblfld(iget(467))
1323  fld_info(cfld)%lvl=lvlsxml(ifd,iget(467))
1324 !$omp parallel do private(i,j,jj)
1325  do j=1,jend-jsta+1
1326  jj = jsta+j-1
1327  do i=1,im
1328  datapd(i,j,cfld) = grid1(i,jj)
1329  enddo
1330  enddo
1331  endif
1332  ENDIF
1333  ENDDO
1334  endif
1335 
1336  if(iget(468)>0) THEN ! CAT
1337  n=iavblfld(iget(468))
1338  nfdctl=size(pset%param(n)%level)
1339  if(allocated(itypefdlvlctl)) deallocate(itypefdlvlctl)
1340  allocate(itypefdlvlctl(nfdctl))
1341  DO ifd = 1,nfdctl
1342  itypefdlvlctl(ifd)=lvls(ifd,iget(468))
1343  enddo
1344  if(allocated(htfdctl)) deallocate(htfdctl)
1345  allocate(htfdctl(nfdctl))
1346  htfdctl=pset%param(n)%level
1347  allocate(catfd(im,jsta:jend,nfdctl))
1348  call fdlvl_mass(itypefdlvlctl,nfdctl,htfdctl,catedr,catfd)
1349  DO ifd = 1,nfdctl
1350  IF (lvls(ifd,iget(468))>0) THEN
1351 !$omp parallel do private(i,j)
1352  DO j=jsta,jend
1353  DO i=1,im
1354  grid1(i,j)=catfd(i,j,ifd)
1355  ENDDO
1356  ENDDO
1357  if(grib=='grib2') then
1358  cfld=cfld+1
1359  fld_info(cfld)%ifld=iavblfld(iget(468))
1360  fld_info(cfld)%lvl=lvlsxml(ifd,iget(468))
1361 !$omp parallel do private(i,j,jj)
1362  do j=1,jend-jsta+1
1363  jj = jsta+j-1
1364  do i=1,im
1365  datapd(i,j,cfld) = grid1(i,jj)
1366  enddo
1367  enddo
1368  endif
1369  ENDIF
1370  ENDDO
1371  endif
1372 
1373  if(iget(469)>0) THEN ! MWT
1374  n=iavblfld(iget(469))
1375  nfdctl=size(pset%param(n)%level)
1376  if(allocated(itypefdlvlctl)) deallocate(itypefdlvlctl)
1377  allocate(itypefdlvlctl(nfdctl))
1378  DO ifd = 1,nfdctl
1379  itypefdlvlctl(ifd)=lvls(ifd,iget(469))
1380  enddo
1381  if(allocated(htfdctl)) deallocate(htfdctl)
1382  allocate(htfdctl(nfdctl))
1383  htfdctl=pset%param(n)%level
1384  allocate(mwtfd(im,jsta:jend,nfdctl))
1385  call fdlvl_mass(itypefdlvlctl,nfdctl,htfdctl,mwt,mwtfd)
1386  DO ifd = 1,nfdctl
1387  IF (lvls(ifd,iget(469))>0) THEN
1388 !$omp parallel do private(i,j)
1389  DO j=jsta,jend
1390  DO i=1,im
1391  grid1(i,j)=mwtfd(i,j,ifd)
1392  ENDDO
1393  ENDDO
1394  if(grib=='grib2') then
1395  cfld=cfld+1
1396  fld_info(cfld)%ifld=iavblfld(iget(469))
1397  fld_info(cfld)%lvl=lvlsxml(ifd,iget(469))
1398 !$omp parallel do private(i,j,jj)
1399  do j=1,jend-jsta+1
1400  jj = jsta+j-1
1401  do i=1,im
1402  datapd(i,j,cfld) = grid1(i,jj)
1403  enddo
1404  enddo
1405  endif
1406  ENDIF
1407  ENDDO
1408  endif
1409 
1410  if(allocated(gtgfd)) deallocate(gtgfd)
1411  if(allocated(catfd)) deallocate(catfd)
1412  if(allocated(mwtfd)) deallocate(mwtfd)
1413 
1414  if(allocated(itypefdlvlctl)) deallocate(itypefdlvlctl)
1415  if(allocated(htfdctl)) deallocate(htfdctl)
1416 
1417  ENDIF
1418 !
1419 !
1420 ! ***BLOCK 4: FREEZING LEVEL Z, RH AND P.
1421 !
1422  IF ( (iget(062)>0).OR.(iget(063)>0) ) THEN
1423  CALL frzlvl(z1d,rh1d,p1d)
1424 !
1425 ! FREEZING LEVEL HEIGHT.
1426  IF (iget(062)>0) THEN
1427 !$omp parallel do private(i,j)
1428  DO j=jsta,jend
1429  DO i=1,im
1430  grid1(i,j)=z1d(i,j)
1431  IF (submodelname == 'RTMA') THEN
1432  freezelvl(i,j)=grid1(i,j)
1433  ENDIF
1434  ENDDO
1435  ENDDO
1436  CALL bound(grid1,d00,h99999)
1437  if(grib=='grib2') then
1438  cfld=cfld+1
1439  fld_info(cfld)%ifld=iavblfld(iget(062))
1440 !$omp parallel do private(i,j,jj)
1441  do j=1,jend-jsta+1
1442  jj = jsta+j-1
1443  do i=1,im
1444  datapd(i,j,cfld) = grid1(i,jj)
1445  enddo
1446  enddo
1447  endif
1448  ENDIF
1449 !
1450 ! FREEZING LEVEL RELATIVE HUMIDITY.
1451  IF (iget(063)>0) THEN
1452 !$omp parallel do private(i,j)
1453  DO j=jsta,jend
1454  DO i=1,im
1455  grid1(i,j) = rh1d(i,j)
1456  ENDDO
1457  ENDDO
1458  CALL sclfld(grid1,h100,im,jm)
1459  CALL bound(grid1,h1,h100)
1460  if(grib=='grib2') then
1461  cfld=cfld+1
1462  fld_info(cfld)%ifld=iavblfld(iget(063))
1463 !$omp parallel do private(i,j,jj)
1464  do j=1,jend-jsta+1
1465  jj = jsta+j-1
1466  do i=1,im
1467  datapd(i,j,cfld) = grid1(i,jj)
1468  enddo
1469  enddo
1470  endif
1471  ENDIF
1472 !
1473 ! FREEZING LEVEL PRESSURE
1474  IF (iget(753)>0) THEN
1475 !$omp parallel do private(i,j)
1476  DO j=jsta,jend
1477  DO i=1,im
1478  grid1(i,j) = p1d(i,j)
1479  ENDDO
1480  ENDDO
1481  if(grib=='grib2') then
1482  cfld=cfld+1
1483  fld_info(cfld)%ifld=iavblfld(iget(753))
1484 !$omp parallel do private(i,j,jj)
1485  do j=1,jend-jsta+1
1486  jj = jsta+j-1
1487  do i=1,im
1488  datapd(i,j,cfld) = grid1(i,jj)
1489  enddo
1490  enddo
1491  endif
1492 
1493  ENDIF
1494  ENDIF
1495 
1496  IF (iget(165)>0 .OR. iget(350)>0.OR. iget(756)>0) THEN
1497  CALL frzlvl2(tfrz,z1d,rh1d,p1d)
1498 !
1499 ! HIGHEST FREEZING LEVEL HEIGHT.
1500  IF (iget(165)>0)THEN
1501 !$omp parallel do private(i,j)
1502  DO j=jsta,jend
1503  DO i=1,im
1504  grid1(i,j)=z1d(i,j)
1505  ENDDO
1506  ENDDO
1507  CALL bound(grid1,d00,h99999)
1508  if(grib=='grib2') then
1509  cfld=cfld+1
1510  fld_info(cfld)%ifld=iavblfld(iget(165))
1511 !$omp parallel do private(i,j,jj)
1512  do j=1,jend-jsta+1
1513  jj = jsta+j-1
1514  do i=1,im
1515  datapd(i,j,cfld) = grid1(i,jj)
1516  enddo
1517  enddo
1518  endif
1519  END IF
1520 
1521 ! HIGHEST FREEZING LEVEL RELATIVE HUMIDITY
1522  IF (iget(350)>0)THEN
1523  grid1=spval
1524 !$omp parallel do private(i,j)
1525  DO j=jsta,jend
1526  DO i=1,im
1527  IF(rh1d(i,j) < spval) grid1(i,j)=rh1d(i,j)*100.
1528  ENDDO
1529  ENDDO
1530  CALL bound(grid1,h1,h100)
1531  if(grib=='grib2') then
1532  cfld=cfld+1
1533  fld_info(cfld)%ifld=iavblfld(iget(350))
1534 !$omp parallel do private(i,j,jj)
1535  do j=1,jend-jsta+1
1536  jj = jsta+j-1
1537  do i=1,im
1538  datapd(i,j,cfld) = grid1(i,jj)
1539  enddo
1540  enddo
1541  endif
1542  END IF
1543 
1544 ! HIGHEST FREEZING LEVEL PRESSURE
1545  IF (iget(756)>0) THEN
1546 !$omp parallel do private(i,j)
1547  DO j=jsta,jend
1548  DO i=1,im
1549  grid1(i,j) = p1d(i,j)
1550  ENDDO
1551  ENDDO
1552  if(grib=='grib2') then
1553  cfld=cfld+1
1554  fld_info(cfld)%ifld=iavblfld(iget(756))
1555 !$omp parallel do private(i,j,jj)
1556  do j=1,jend-jsta+1
1557  jj = jsta+j-1
1558  do i=1,im
1559  datapd(i,j,cfld) = grid1(i,jj)
1560  enddo
1561  enddo
1562  endif
1563  ENDIF
1564 
1565  ENDIF
1566 !
1567 ! HIGHEST -10 C ISOTHERM VALUES
1568 !
1569  IF (iget(776)>0 .OR. iget(777)>0.OR. iget(778)>0) THEN
1570  CALL frzlvl2(263.15,z1d,rh1d,p1d)
1571 !
1572 ! HIGHEST -10C ISOTHERM LEVEL HEIGHT.
1573  IF (iget(776)>0)THEN
1574 !$omp parallel do private(i,j)
1575  DO j=jsta,jend
1576  DO i=1,im
1577  grid1(i,j)=z1d(i,j)
1578  ENDDO
1579  ENDDO
1580  CALL bound(grid1,d00,h99999)
1581  if(grib=='grib2') then
1582  cfld=cfld+1
1583  fld_info(cfld)%ifld=iavblfld(iget(776))
1584 !$omp parallel do private(i,j,jj)
1585  do j=1,jend-jsta+1
1586  jj = jsta+j-1
1587  do i=1,im
1588  datapd(i,j,cfld) = grid1(i,jj)
1589  enddo
1590  enddo
1591  endif
1592  END IF
1593 
1594 ! HIGHEST -10C ISOTHERM RELATIVE HUMIDITY
1595  IF (iget(777)>0)THEN
1596  grid1=spval
1597 !$omp parallel do private(i,j)
1598  DO j=jsta,jend
1599  DO i=1,im
1600  IF(rh1d(i,j) < spval) grid1(i,j)=rh1d(i,j)*100.
1601  ENDDO
1602  ENDDO
1603  CALL bound(grid1,h1,h100)
1604  if(grib=='grib2') then
1605  cfld=cfld+1
1606  fld_info(cfld)%ifld=iavblfld(iget(777))
1607 !$omp parallel do private(i,j,jj)
1608  do j=1,jend-jsta+1
1609  jj = jsta+j-1
1610  do i=1,im
1611  datapd(i,j,cfld) = grid1(i,jj)
1612  enddo
1613  enddo
1614  endif
1615  END IF
1616 
1617 ! HIGHEST -10C ISOTHERM LEVEL PRESSURE
1618  IF (iget(778)>0) THEN
1619 !$omp parallel do private(i,j)
1620  DO j=jsta,jend
1621  DO i=1,im
1622  grid1(i,j)=p1d(i,j)
1623  ENDDO
1624  ENDDO
1625  if(grib=='grib2') then
1626  cfld=cfld+1
1627  fld_info(cfld)%ifld=iavblfld(iget(778))
1628 !$omp parallel do private(i,j,jj)
1629  do j=1,jend-jsta+1
1630  jj = jsta+j-1
1631  do i=1,im
1632  datapd(i,j,cfld) = grid1(i,jj)
1633  enddo
1634  enddo
1635  endif
1636  ENDIF
1637 
1638  ENDIF
1639 !
1640 ! HIGHEST -20 C ISOTHERM VALUES
1641 !
1642  IF (iget(779)>0 .OR. iget(780)>0.OR. iget(781)>0) THEN
1643  CALL frzlvl2(253.15,z1d,rh1d,p1d)
1644 !
1645 ! HIGHEST -20C ISOTHERM LEVEL HEIGHT.
1646  IF (iget(779)>0)THEN
1647 !$omp parallel do private(i,j)
1648  DO j=jsta,jend
1649  DO i=1,im
1650  grid1(i,j)=z1d(i,j)
1651  ENDDO
1652  ENDDO
1653  CALL bound(grid1,d00,h99999)
1654  if(grib=='grib2') then
1655  cfld=cfld+1
1656  fld_info(cfld)%ifld=iavblfld(iget(779))
1657 !$omp parallel do private(i,j,jj)
1658  do j=1,jend-jsta+1
1659  jj = jsta+j-1
1660  do i=1,im
1661  datapd(i,j,cfld) = grid1(i,jj)
1662  enddo
1663  enddo
1664  endif
1665  END IF
1666 
1667 ! HIGHEST -20C ISOTHERM RELATIVE HUMIDITY
1668  IF (iget(780)>0)THEN
1669  grid1=spval
1670 !$omp parallel do private(i,j)
1671  DO j=jsta,jend
1672  DO i=1,im
1673  IF(rh1d(i,j) < spval) grid1(i,j)=rh1d(i,j)*100.
1674  ENDDO
1675  ENDDO
1676  CALL bound(grid1,h1,h100)
1677  if(grib=='grib2') then
1678  cfld=cfld+1
1679  fld_info(cfld)%ifld=iavblfld(iget(780))
1680 !$omp parallel do private(i,j,jj)
1681  do j=1,jend-jsta+1
1682  jj = jsta+j-1
1683  do i=1,im
1684  datapd(i,j,cfld) = grid1(i,jj)
1685  enddo
1686  enddo
1687  endif
1688  END IF
1689 
1690 ! HIGHEST -20C ISOTHERM LEVEL PRESSURE
1691  IF (iget(781)>0) THEN
1692 !$omp parallel do private(i,j)
1693  DO j=jsta,jend
1694  DO i=1,im
1695  grid1(i,j)=p1d(i,j)
1696  ENDDO
1697  ENDDO
1698  if(grib=='grib2') then
1699  cfld=cfld+1
1700  fld_info(cfld)%ifld=iavblfld(iget(781))
1701 !$omp parallel do private(i,j,jj)
1702  do j=1,jend-jsta+1
1703  jj = jsta+j-1
1704  do i=1,im
1705  datapd(i,j,cfld) = grid1(i,jj)
1706  enddo
1707  enddo
1708  endif
1709  ENDIF
1710 
1711  ENDIF
1712 !
1713  allocate(pbnd(im,jsta:jend,nbnd), tbnd(im,jsta:jend,nbnd), &
1714  qbnd(im,jsta:jend,nbnd), ubnd(im,jsta:jend,nbnd), &
1715  vbnd(im,jsta:jend,nbnd), rhbnd(im,jsta:jend,nbnd), &
1716  wbnd(im,jsta:jend,nbnd))
1717 
1718 !
1719 ! ***BLOCK 5: BOUNDARY LAYER FIELDS.
1720 !
1721  IF ( (iget(067)>0).OR.(iget(068)>0).OR. &
1722  (iget(069)>0).OR.(iget(070)>0).OR. &
1723  (iget(071)>0).OR.(iget(072)>0).OR. &
1724  (iget(073)>0).OR.(iget(074)>0).OR. &
1725  (iget(088)>0).OR.(iget(089)>0).OR. &
1726  (iget(090)>0).OR.(iget(075)>0).OR. &
1727  (iget(109)>0).OR.(iget(110)>0).OR. &
1728  (iget(031)>0).OR.(iget(032)>0).OR. &
1729  (iget(573)>0).OR. &
1730  (iget(107)>0).OR.(iget(091)>0).OR. &
1731  (iget(092)>0).OR.(iget(093)>0).OR. &
1732  (iget(094)>0).OR.(iget(095)>0).OR. &
1733  (iget(096)>0).OR.(iget(097)>0).OR. &
1734  (iget(098)>0).OR.(iget(221)>0) ) THEN
1735 !
1736  allocate(omgbnd(im,jsta:jend,nbnd),pwtbnd(im,jsta:jend,nbnd), &
1737  qcnvbnd(im,jsta:jend,nbnd),lvlbnd(im,jsta:jend,nbnd), &
1738  lb2(im,jsta:jend))
1739 
1740 ! COMPUTE ETA BOUNDARY LAYER FIELDS.
1741  CALL bndlyr(pbnd,tbnd,qbnd,rhbnd,ubnd,vbnd, &
1742  wbnd,omgbnd,pwtbnd,qcnvbnd,lvlbnd)
1743 
1744 !$omp parallel do private(i,j)
1745  DO j=jsta,jend
1746  DO i=1,im
1747  egrid2(i,j) = spval
1748  ENDDO
1749  ENDDO
1750 
1751 !
1752 ! LOOP OVER NBND BOUNDARY LAYERS.
1753  DO 20 lbnd = 1,nbnd
1754 !
1755 ! BOUNDARY LAYER PRESSURE.
1756  IF (iget(067)>0) THEN
1757  IF (lvls(lbnd,iget(067))>0) THEN
1758 !$omp parallel do private(i,j)
1759  DO j=jsta,jend
1760  DO i=1,im
1761  grid1(i,j) = pbnd(i,j,lbnd)
1762  ENDDO
1763  ENDDO
1764  if(grib=='grib2') then
1765  cfld=cfld+1
1766  fld_info(cfld)%ifld=iavblfld(iget(067))
1767  fld_info(cfld)%lvl=lvlsxml(lbnd,iget(067))
1768 !$omp parallel do private(i,j,jj)
1769  do j=1,jend-jsta+1
1770  jj = jsta+j-1
1771  do i=1,im
1772  datapd(i,j,cfld) = grid1(i,jj)
1773  enddo
1774  enddo
1775  endif
1776  ENDIF
1777  ENDIF
1778 !
1779 ! BOUNDARY LAYER TEMPERATURE.
1780  IF (iget(068)>0) THEN
1781  IF (lvls(lbnd,iget(068))>0) THEN
1782 !$omp parallel do private(i,j)
1783  DO j=jsta,jend
1784  DO i=1,im
1785  grid1(i,j)=tbnd(i,j,lbnd)
1786  ENDDO
1787  ENDDO
1788  if(grib=='grib2') then
1789  cfld=cfld+1
1790  fld_info(cfld)%ifld=iavblfld(iget(068))
1791  fld_info(cfld)%lvl=lvlsxml(lbnd,iget(068))
1792 !$omp parallel do private(i,j,jj)
1793  do j=1,jend-jsta+1
1794  jj = jsta+j-1
1795  do i=1,im
1796  datapd(i,j,cfld) = grid1(i,jj)
1797  enddo
1798  enddo
1799  endif
1800  ENDIF
1801  ENDIF
1802 !
1803 ! BOUNDARY LAYER POTENTIAL TEMPERATURE.
1804  IF (iget(069)>0) THEN
1805  IF (lvls(lbnd,iget(069))>0) THEN
1806  CALL calpot(pbnd(1,jsta,lbnd),tbnd(1,jsta,lbnd),grid1(1,jsta))
1807  if(grib=='grib2') then
1808  cfld=cfld+1
1809  fld_info(cfld)%ifld=iavblfld(iget(069))
1810  fld_info(cfld)%lvl=lvlsxml(ifd,iget(069))
1811 !$omp parallel do private(i,j,jj)
1812  do j=1,jend-jsta+1
1813  jj = jsta+j-1
1814  do i=1,im
1815  datapd(i,j,cfld) = grid1(i,jj)
1816  enddo
1817  enddo
1818  endif
1819  ENDIF
1820  ENDIF
1821 !
1822 ! BOUNDARY LAYER RELATIVE HUMIDITY.
1823  IF (iget(072)>0) THEN
1824  IF (lvls(lbnd,iget(072))>0) THEN
1825 !$omp parallel do private(i,j)
1826  DO j=jsta,jend
1827  DO i=1,im
1828  grid1(i,j)=rhbnd(i,j,lbnd)
1829  ENDDO
1830  ENDDO
1831  CALL sclfld(grid1,h100,im,jm)
1832  CALL bound(grid1,h1,h100)
1833  if(grib=='grib2') then
1834  cfld=cfld+1
1835  fld_info(cfld)%lvl=lvlsxml(lbnd,iget(072))
1836  fld_info(cfld)%ifld=iavblfld(iget(072))
1837 !$omp parallel do private(i,j,jj)
1838  do j=1,jend-jsta+1
1839  jj = jsta+j-1
1840  do i=1,im
1841  datapd(i,j,cfld) = grid1(i,jj)
1842  enddo
1843  enddo
1844  endif
1845  ENDIF
1846  ENDIF
1847 !
1848 ! BOUNDARY LAYER DEWPOINT TEMPERATURE.
1849  IF (iget(070)>0) THEN
1850  IF (lvls(lbnd,iget(070))>0) THEN
1851  CALL caldwp(pbnd(1,jsta,lbnd), qbnd(1,jsta,lbnd), &
1852  grid1(1,jsta), tbnd(1,jsta,lbnd))
1853  if(grib=='grib2') then
1854  cfld=cfld+1
1855  fld_info(cfld)%ifld=iavblfld(iget(070))
1856  fld_info(cfld)%lvl=lvlsxml(lbnd,iget(070))
1857 !$omp parallel do private(i,j,jj)
1858  do j=1,jend-jsta+1
1859  jj = jsta+j-1
1860  do i=1,im
1861  datapd(i,j,cfld) = grid1(i,jj)
1862  enddo
1863  enddo
1864  endif
1865  ENDIF
1866  ENDIF
1867 !
1868 ! BOUNDARY LAYER SPECIFIC HUMIDITY.
1869  IF (iget(071)>0) THEN
1870  IF (lvls(lbnd,iget(071))>0) THEN
1871 !$omp parallel do private(i,j)
1872  DO j=jsta,jend
1873  DO i=1,im
1874  grid1(i,j)=qbnd(i,j,lbnd)
1875  ENDDO
1876  ENDDO
1877  CALL bound(grid1,h1m12,h99999)
1878  if(grib=='grib2') then
1879  cfld=cfld+1
1880  fld_info(cfld)%ifld=iavblfld(iget(071))
1881  fld_info(cfld)%lvl=lvlsxml(lbnd,iget(071))
1882 !$omp parallel do private(i,j,jj)
1883  do j=1,jend-jsta+1
1884  jj = jsta+j-1
1885  do i=1,im
1886  datapd(i,j,cfld) = grid1(i,jj)
1887  enddo
1888  enddo
1889  endif
1890  ENDIF
1891  ENDIF
1892 !
1893 ! BOUNDARY LAYER MOISTURE CONVERGENCE.
1894  IF (iget(088)>0) THEN
1895  IF (lvls(lbnd,iget(088))>0) THEN
1896 !$omp parallel do private(i,j)
1897  DO j=jsta,jend
1898  DO i=1,im
1899  grid1(i,j) = qcnvbnd(i,j,lbnd)
1900  ENDDO
1901  ENDDO
1902  if(grib=='grib2') then
1903  cfld=cfld+1
1904  fld_info(cfld)%ifld=iavblfld(iget(088))
1905  fld_info(cfld)%lvl=lvlsxml(lbnd,iget(088))
1906 !$omp parallel do private(i,j,jj)
1907  do j=1,jend-jsta+1
1908  jj = jsta+j-1
1909  do i=1,im
1910  datapd(i,j,cfld) = grid1(i,jj)
1911  enddo
1912  enddo
1913  endif
1914  ENDIF
1915  ENDIF
1916 !
1917 ! BOUNDARY LAYER U WIND AND/OR V WIND.
1918 !
1919  field1=.false.
1920  field2=.false.
1921 !
1922  IF(iget(073)>0)THEN
1923  IF(lvls(lbnd,iget(073))>0)field1=.true.
1924  ENDIF
1925  IF(iget(074)>0)THEN
1926  IF(lvls(lbnd,iget(074))>0)field2=.true.
1927  ENDIF
1928 !
1929  IF(field1.OR.field2)THEN
1930 !$omp parallel do private(i,j)
1931  DO j=jsta,jend
1932  DO i=1,im
1933  grid1(i,j) = ubnd(i,j,lbnd)
1934  grid2(i,j) = vbnd(i,j,lbnd)
1935  ENDDO
1936  ENDDO
1937 !
1938  IF (iget(073)>0) THEN
1939  IF (lvls(lbnd,iget(073))>0) then
1940  if(grib=='grib2') then
1941  cfld=cfld+1
1942  fld_info(cfld)%ifld=iavblfld(iget(073))
1943  fld_info(cfld)%lvl=lvlsxml(lbnd,iget(073))
1944 !$omp parallel do private(i,j,jj)
1945  do j=1,jend-jsta+1
1946  jj = jsta+j-1
1947  do i=1,im
1948  datapd(i,j,cfld) = grid1(i,jj)
1949  enddo
1950  enddo
1951  endif
1952  ENDIF
1953  ENDIF
1954  IF (iget(074)>0) THEN
1955  IF (lvls(lbnd,iget(074))>0) THEN
1956  if(grib=='grib2') then
1957  cfld=cfld+1
1958  fld_info(cfld)%ifld=iavblfld(iget(074))
1959  fld_info(cfld)%lvl=lvlsxml(lbnd,iget(074))
1960 !$omp parallel do private(i,j,jj)
1961  do j=1,jend-jsta+1
1962  jj = jsta+j-1
1963  do i=1,im
1964  datapd(i,j,cfld) = grid2(i,jj)
1965  enddo
1966  enddo
1967  endif
1968  ENDIF
1969  ENDIF
1970  ENDIF
1971 !
1972 ! BOUNDARY LAYER OMEGA.
1973  IF (iget(090)>0) THEN
1974  IF (lvls(lbnd,iget(090))>0) THEN
1975 !$omp parallel do private(i,j)
1976  DO j=jsta,jend
1977  DO i=1,im
1978  grid1(i,j) = omgbnd(i,j,lbnd)
1979  ENDDO
1980  ENDDO
1981  if(grib=='grib2') then
1982  cfld=cfld+1
1983  fld_info(cfld)%ifld=iavblfld(iget(090))
1984  fld_info(cfld)%lvl=lvlsxml(lbnd,iget(090))
1985 !$omp parallel do private(i,j,jj)
1986  do j=1,jend-jsta+1
1987  jj = jsta+j-1
1988  do i=1,im
1989  datapd(i,j,cfld) = grid1(i,jj)
1990  enddo
1991  enddo
1992  endiF
1993  ENDIF
1994  ENDIF
1995 !
1996 ! BOUNDARY LAYER PRECIPITBLE WATER.
1997  IF (iget(089)>0) THEN
1998  IF (lvls(lbnd,iget(089))>0) THEN
1999 !$omp parallel do private(i,j)
2000  DO j=jsta,jend
2001  DO i=1,im
2002  grid1(i,j) = pwtbnd(i,j,lbnd)
2003  ENDDO
2004  ENDDO
2005  CALL bound(grid1,d00,h99999)
2006  if(grib=='grib2') then
2007  cfld=cfld+1
2008  fld_info(cfld)%ifld=iavblfld(iget(089))
2009  fld_info(cfld)%lvl=lvlsxml(lbnd,iget(089))
2010 !$omp parallel do private(i,j,jj)
2011  do j=1,jend-jsta+1
2012  jj = jsta+j-1
2013  do i=1,im
2014  datapd(i,j,cfld) = grid1(i,jj)
2015  enddo
2016  enddo
2017  endif
2018  ENDIF
2019  ENDIF
2020 !
2021 ! BOUNDARY LAYER LIFTED INDEX.
2022  IF (iget(075)>0 .OR. iget(031)>0 .OR. iget(573)>0) THEN
2023  CALL otlft(pbnd(1,jsta,lbnd),tbnd(1,jsta,lbnd), &
2024  qbnd(1,jsta,lbnd),grid1(1,jsta))
2025  IF(iget(075)>0)THEN
2026  IF (lvls(lbnd,iget(075))>0) THEN
2027  if(grib=='grib2') then
2028  cfld=cfld+1
2029  fld_info(cfld)%ifld=iavblfld(iget(075))
2030  fld_info(cfld)%lvl=lvlsxml(lbnd,iget(075))
2031 !$omp parallel do private(i,j,jj)
2032  do j=1,jend-jsta+1
2033  jj = jsta+j-1
2034  do i=1,im
2035  datapd(i,j,cfld) = grid1(i,jj)
2036  enddo
2037  enddo
2038  endif
2039  ENDIF
2040  END IF
2041  IF(iget(031)>0 .or. iget(573)>0)THEN
2042 !$omp parallel do private(i,j)
2043  DO j=jsta,jend
2044  DO i=1,im
2045  egrid2(i,j) = min(egrid2(i,j),grid1(i,j))
2046  END DO
2047  END DO
2048  END IF
2049  ENDIF
2050 !
2051 ! END OF ETA BOUNDARY LAYER LOOP.
2052  20 CONTINUE
2053  deallocate(omgbnd,pwtbnd,qcnvbnd)
2054 !
2055 ! BEST LIFTED INDEX FROM BOUNDARY LAYER FIELDS.
2056 !
2057  IF (iget(031)>0 .OR. iget(573)>0 ) THEN
2058 ! DO J=JSTA,JEND
2059 ! DO I=1,IM
2060 ! EGRID1(I,J) = H99999
2061 ! EGRID2(I,J) = H99999
2062 ! ENDDO
2063 ! ENDDO
2064 !
2065 ! DO 50 LBND = 1,NBND
2066 ! CALL OTLFT(PBND(1,1,LBND),TBND(1,1,LBND), &
2067 ! QBND(1,1,LBND),EGRID2)
2068 ! DO J=JSTA,JEND
2069 ! DO I=1,IM
2070 ! EGRID1(I,J)=AMIN1(EGRID1(I,J),EGRID2(I,J))
2071 ! ENDDO
2072 ! ENDDO
2073 ! 50 CONTINUE
2074 !$omp parallel do private(i,j)
2075  DO j=jsta,jend
2076  DO i=1,im
2077  grid1(i,j)=egrid2(i,j)
2078  ENDDO
2079  ENDDO
2080 ! print*,'writting out best lifted index'
2081 
2082  if (iget(031)>0) then
2083  if(grib=='grib2') then
2084  cfld=cfld+1
2085  fld_info(cfld)%ifld=iavblfld(iget(031))
2086  datapd(1:im,1:jend-jsta+1,cfld)=grid1(1:im,jsta:jend)
2087  endif
2088  endif
2089 
2090  if(iget(573)> 0 ) THEN
2091  if(grib=='grib2') then
2092  cfld=cfld+1
2093  fld_info(cfld)%ifld=iavblfld(iget(573))
2094 !$omp parallel do private(i,j,jj)
2095  do j=1,jend-jsta+1
2096  jj = jsta+j-1
2097  do i=1,im
2098  datapd(i,j,cfld) = grid1(i,jj)
2099  enddo
2100  enddo
2101  endif
2102  endif
2103 
2104  END IF
2105 !
2106 ! BEST BOUNDARY LAYER CAPE AND CINS.
2107 !
2108  field1=.false.
2109  field2=.false.
2110 !
2111  IF(iget(032)>0)THEN
2112  IF(lvls(2,iget(032))>0)field1=.true.
2113  ENDIF
2114  IF(iget(107)>0)THEN
2115  IF(lvls(2,iget(107))>0)field2=.true.
2116  ENDIF
2117 !
2118  IF(iget(566)>0)THEN
2119  field1=.true.
2120  ENDIF
2121  IF(iget(567)>0)THEN
2122  field2=.true.
2123  ENDIF
2124 !
2125  !if(grib=="grib2") print *,'in MISCLN.f,iget(566)=', &
2126  ! iget(566), 'iget(567)=',iget(567),'LVLSXML(1,IGET(566)=', &
2127  ! LVLSXML(1,IGET(566)),'LVLSXML(1,IGET(567)=', &
2128  ! LVLSXML(1,IGET(567)),'field1=',field1,'field2=',field2
2129 !
2130  IF(field1.OR.field2)THEN
2131  itype = 2
2132 !
2133 !$omp parallel do private(i,j)
2134  DO j=jsta,jend
2135  DO i=1,im
2136  egrid1(i,j) = -h99999
2137  egrid2(i,j) = -h99999
2138  ENDDO
2139  ENDDO
2140 !
2141  DO 80 lbnd = 1,nbnd
2142  CALL calthte(pbnd(1,jsta,lbnd),tbnd(1,jsta,lbnd), &
2143  qbnd(1,jsta,lbnd),egrid1)
2144 !$omp parallel do private(i,j)
2145  DO j=jsta,jend
2146  DO i=1,im
2147  IF (egrid1(i,j) > egrid2(i,j)) THEN
2148  egrid2(i,j) = egrid1(i,j)
2149  lb2(i,j) = lvlbnd(i,j,lbnd)
2150  p1d(i,j) = pbnd(i,j,lbnd)
2151  t1d(i,j) = tbnd(i,j,lbnd)
2152  q1d(i,j) = qbnd(i,j,lbnd)
2153  ENDIF
2154  ENDDO
2155  ENDDO
2156  80 CONTINUE
2157 !
2158  dpbnd = 0.
2159  CALL calcape(itype,dpbnd,p1d,t1d,q1d,lb2,egrid1, &
2160  egrid2,egrid3,egrid4,egrid5)
2161 !
2162  IF (iget(566)>0) THEN
2163 ! dong add missing value for cape
2164  grid1=spval
2165 !$omp parallel do private(i,j)
2166  DO j=jsta,jend
2167  DO i=1,im
2168  IF(t1d(i,j) < spval) grid1(i,j) = egrid1(i,j)
2169  ENDDO
2170  ENDDO
2171  CALL bound(grid1,d00,h99999)
2172  if(grib=='grib2') then
2173  cfld=cfld+1
2174  fld_info(cfld)%ifld=iavblfld(iget(566))
2175  fld_info(cfld)%lvl=lvlsxml(1,iget(566))
2176 !$omp parallel do private(i,j,jj)
2177  do j=1,jend-jsta+1
2178  jj = jsta+j-1
2179  do i=1,im
2180  datapd(i,j,cfld) = grid1(i,jj)
2181  enddo
2182  enddo
2183  endif
2184  ENDIF
2185 !
2186  IF (iget(567) > 0) THEN
2187 ! dong add missing value for cape
2188  grid1=spval
2189 !$omp parallel do private(i,j)
2190  DO j=jsta,jend
2191  DO i=1,im
2192  IF(t1d(i,j) < spval) grid1(i,j) = - egrid2(i,j)
2193  ENDDO
2194  ENDDO
2195 !
2196  CALL bound(grid1,d00,h99999)
2197 !
2198 !$omp parallel do private(i,j)
2199  DO j=jsta,jend
2200  DO i=1,im
2201  IF(t1d(i,j) < spval) grid1(i,j) = - grid1(i,j)
2202  ENDDO
2203  ENDDO
2204 !
2205  if(grib=='grib2') then
2206  cfld=cfld+1
2207  fld_info(cfld)%ifld=iavblfld(iget(567))
2208  fld_info(cfld)%lvl=lvlsxml(1,iget(567))
2209 !$omp parallel do private(i,j,jj)
2210  do j=1,jend-jsta+1
2211  jj = jsta+j-1
2212  do i=1,im
2213  datapd(i,j,cfld) = grid1(i,jj)
2214  enddo
2215  enddo
2216  endif
2217  ENDIF
2218  ENDIF
2219 !
2220 
2221 ! PBL HEIGHT
2222  IF(iget(221) > 0) THEN
2223 !$omp parallel do private(i,j)
2224  DO j=jsta,jend
2225  DO i=1,im
2226  grid1(i,j) = pblh(i,j)
2227  ENDDO
2228  ENDDO
2229  if(grib=='grib2') then
2230  cfld=cfld+1
2231  fld_info(cfld)%ifld=iavblfld(iget(221))
2232 !$omp parallel do private(i,j,jj)
2233  do j=1,jend-jsta+1
2234  jj = jsta+j-1
2235  do i=1,im
2236  datapd(i,j,cfld) = grid1(i,jj)
2237  enddo
2238  enddo
2239  endif
2240  END IF
2241 ! BOUNDARY LAYER LIFTING CONDENSATION PRESSURE AND HEIGHT.
2242 ! EGRID1 IS LCL PRESSURE. EGRID2 IS LCL HEIGHT.
2243 !
2244  IF ( (iget(109)>0).OR.(iget(110)>0) ) THEN
2245  CALL callcl(pbnd(1,jsta,1),tbnd(1,jsta,1), &
2246  qbnd(1,jsta,1),egrid1,egrid2)
2247  IF (iget(109)>0) THEN
2248  grid1=spval
2249 !$omp parallel do private(i,j)
2250  DO j=jsta,jend
2251  DO i=1,im
2252  IF(tbnd(i,j,1) < spval) grid1(i,j) = egrid2(i,j)
2253  ENDDO
2254  ENDDO
2255  if(grib=='grib2') then
2256  cfld=cfld+1
2257  fld_info(cfld)%ifld=iavblfld(iget(109))
2258 !$omp parallel do private(i,j,jj)
2259  do j=1,jend-jsta+1
2260  jj = jsta+j-1
2261  do i=1,im
2262  datapd(i,j,cfld) = grid1(i,jj)
2263  enddo
2264  enddo
2265  endif
2266  ENDIF
2267  IF (iget(110)>0) THEN
2268  grid1=spval
2269 !$omp parallel do private(i,j)
2270  DO j=jsta,jend
2271  DO i=1,im
2272  IF(tbnd(i,j,1) < spval) grid1(i,j) = egrid1(i,j)
2273  ENDDO
2274  ENDDO
2275  if(grib=='grib2') then
2276  cfld=cfld+1
2277  fld_info(cfld)%ifld=iavblfld(iget(110))
2278 !$omp parallel do private(i,j,jj)
2279  do j=1,jend-jsta+1
2280  jj = jsta+j-1
2281  do i=1,im
2282  datapd(i,j,cfld) = grid1(i,jj)
2283  enddo
2284  enddo
2285  endif
2286  ENDIF
2287  ENDIF
2288 !
2289 ! NGM BOUNDARY LAYER FIELDS.
2290 !
2291  IF ( (iget(091)>0).OR.(iget(092)>0).OR. &
2292  (iget(093)>0).OR.(iget(094)>0).OR. &
2293  (iget(095)>0).OR.(iget(095)>0).OR. &
2294  (iget(096)>0).OR.(iget(097)>0).OR. &
2295  (iget(098)>0) ) THEN
2296 
2297  allocate(t78483(im,jsta:jend), t89671(im,jsta:jend), &
2298  p78483(im,jsta:jend), p89671(im,jsta:jend))
2299 !
2300 ! COMPUTE SIGMA 0.89671 AND 0.78483 TEMPERATURES
2301 ! INTERPOLATE LINEAR IN LOG P
2302  IF (iget(097)>0.OR.iget(098)>0) THEN
2303 !$omp parallel do private(i,j)
2304  DO j=jsta,jend
2305  DO i=1,im
2306  p78483(i,j) = log(pint(i,j,nint(lmh(i,j)))*0.78483)
2307  p89671(i,j) = log(pint(i,j,nint(lmh(i,j)))*0.89671)
2308  ENDDO
2309  ENDDO
2310  done =.false.
2311  done1=.false.
2312 !!$omp parallel do private(fac1,fac2,pkl1,pku1,t78483,t89671)
2313  DO l=2,lm
2314  DO j=jsta,jend
2315  DO i=1,im
2316  pkl1=0.5*(alpint(i,j,l)+alpint(i,j,l+1))
2317  pku1=0.5*(alpint(i,j,l)+alpint(i,j,l-1))
2318 ! IF(I==1 .AND. J==1)PRINT*,'L,P89671,PKL1,PKU1= ', &
2319 ! L,P89671(I,J), PKL1, PKU1
2320  IF(p78483(i,j) < pkl1.AND.p78483(i,j) > pku1) THEN
2321  fac1 = (pkl1-p78483(i,j))/(pkl1-pku1)
2322  fac2 = (p78483(i,j)-pku1)/(pkl1-pku1)
2323  t78483(i,j) = t(i,j,l)*fac2 + t(i,j,l-1)*fac1
2324  done1(i,j)=.true.
2325  ENDIF
2326  IF(p89671(i,j) < pkl1.AND.p89671(i,j) > pku1)THEN
2327  fac1 = (pkl1-p89671(i,j))/(pkl1-pku1)
2328  fac2 = (p89671(i,j)-pku1)/(pkl1-pku1)
2329  t89671(i,j) = t(i,j,l)*fac2 + t(i,j,l-1)*fac1
2330  done(i,j) = .true.
2331 ! IF(I==1 .AND. J==1)PRINT*,'done(1,1)= ',done(1,1)
2332  ENDIF
2333  ENDDO
2334  ENDDO
2335  ENDDO
2336 ! print*,'done(1,1)= ',done(1,1)
2337 !$omp parallel do private(i,j,pl,tl,ql,qsat,rhl)
2338  DO j=jsta,jend
2339  DO i=1,im
2340  IF(.NOT. done(i,j)) THEN
2341  pl = pint(i,j,lm-1)
2342  tl = 0.5*(t(i,j,lm-2)+t(i,j,lm-1))
2343  ql = 0.5*(q(i,j,lm-2)+q(i,j,lm-1))
2344  qsat = pq0/pl *exp(a2*(tl-a3)/(tl-a4))
2345 !
2346  rhl = ql/qsat
2347 !
2348  IF(rhl > 1.)THEN
2349  rhl = 1.
2350  ql = rhl*qsat
2351  ENDIF
2352 !
2353  IF(rhl < rhmin)THEN
2354  rhl = rhmin
2355  ql = rhl*qsat
2356  ENDIF
2357 !
2358 ! TVRL = TL*(1.+0.608*QL)
2359 ! TVRBLO = TVRL*(P89671(I,J)/PL)**RGAMOG
2360 ! T89671(I,J) = TVRBLO/(1.+0.608*QL)
2361 
2362  t89671(i,j) = tl * (p89671(i,j)/pl)**rgamog
2363 !
2364 
2365 ! PKL1=0.5*(ALPINT(I,J,LM)+ALPINT(I,J,LM+1))
2366 ! PKU1=0.5*(ALPINT(I,J,LM-1)+ALPINT(I,J,LM))
2367 ! T89671(I,J)=T(I,J,LM)+(T(I,J,LM)-T(I,J,LM-1))*
2368 ! + (P89671(I,J)-PKL1)/(PKL1-PKU1)
2369 
2370 ! print*,'Debug T89671= ',i,j
2371 ! + ,P89671(I,J), PKL1, PKU1
2372 ! + ,T89671(I,J),T(I,J,LM-1),T(I,J,LM)
2373  END IF
2374 
2375  IF(.NOT. done1(i,j))THEN
2376  pl = pint(i,j,lm-1)
2377  tl = 0.5*(t(i,j,lm-2)+t(i,j,lm-1))
2378  ql = 0.5*(q(i,j,lm-2)+q(i,j,lm-1))
2379  qsat = pq0/pl *exp(a2*(tl-a3)/(tl-a4))
2380 !
2381  rhl = ql/qsat
2382 !
2383  IF(rhl>1.)THEN
2384  rhl = 1.
2385  ql = rhl*qsat
2386  ENDIF
2387 !
2388  IF(rhl<rhmin)THEN
2389  rhl = rhmin
2390  ql = rhl*qsat
2391  ENDIF
2392 !
2393 ! TVRL = TL*(1.+0.608*QL)
2394 ! TVRBLO = TVRL*(P78483(I,J)/PL)**RGAMOG
2395 ! T78483(I,J) =TVRBLO/(1.+0.608*QL)
2396 
2397  t78483(i,j) = tl * (p78483(i,j)/pl)**rgamog
2398 !
2399  END IF
2400 
2401  END DO
2402  END DO
2403 !
2404 ! SIGMA 0.89671 TEMPERATURE
2405  IF (iget(097) > 0) THEN
2406  grid1=spval
2407 !$omp parallel do private(i,j)
2408  DO j=jsta,jend
2409  DO i=1,im
2410  IF(t(i,j,lm) < spval) grid1(i,j) = t89671(i,j)
2411 ! IF(T89671(I,J)>350.)PRINT*,'LARGE T89671 ', &
2412 ! I,J,T89671(I,J)
2413  ENDDO
2414  ENDDO
2415  if(grib=='grib2') then
2416  cfld=cfld+1
2417  fld_info(cfld)%ifld=iavblfld(iget(097))
2418  fld_info(cfld)%lvl=lvlsxml(1,iget(097))
2419 !$omp parallel do private(i,j,jj)
2420  do j=1,jend-jsta+1
2421  jj = jsta+j-1
2422  do i=1,im
2423  datapd(i,j,cfld) = grid1(i,jj)
2424  enddo
2425  enddo
2426  endif
2427  ENDIF
2428 !
2429 ! SIGMA 0.78483 TEMPERATURE
2430  IF (iget(098)>0) THEN
2431  grid1=spval
2432 !$omp parallel do private(i,j)
2433  DO j=jsta,jend
2434  DO i=1,im
2435  IF(t(i,j,lm) < spval) grid1(i,j) = t78483(i,j)
2436  ENDDO
2437  ENDDO
2438  if(grib=='grib2') then
2439  cfld=cfld+1
2440  fld_info(cfld)%ifld=iavblfld(iget(098))
2441  fld_info(cfld)%lvl=lvlsxml(1,iget(098))
2442 !$omp parallel do private(i,j,jj)
2443  do j=1,jend-jsta+1
2444  jj = jsta+j-1
2445  do i=1,im
2446  datapd(i,j,cfld) = grid1(i,jj)
2447  enddo
2448  enddo
2449  endif
2450  ENDIF
2451  deallocate(t78483, t89671, p78483, p89671)
2452  ENDIF
2453 !
2454 ! NGM SIGMA LAYER 0.98230 FIELDS. THESE FIELDS ARE
2455 ! THE FIRST ETA LAYER BOUNDARY LAYER FIELDS.
2456 !
2457 !
2458  IF ( (iget(091)>0).OR.(iget(092)>0).OR. &
2459  (iget(093)>0).OR.(iget(094)>0).OR. &
2460  (iget(095)>0).OR.(iget(095)>0).OR. &
2461  (iget(096)>0) ) THEN
2462 !
2463 !
2464 ! PRESSURE.
2465  IF (iget(091)>0) THEN
2466 !$omp parallel do private(i,j)
2467  DO j=jsta,jend
2468  DO i=1,im
2469  grid1(i,j) = pbnd(i,j,1)
2470  ENDDO
2471  ENDDO
2472  if(grib=='grib2') then
2473  cfld=cfld+1
2474  fld_info(cfld)%ifld=iavblfld(iget(091))
2475 !$omp parallel do private(i,j,jj)
2476  do j=1,jend-jsta+1
2477  jj = jsta+j-1
2478  do i=1,im
2479  datapd(i,j,cfld) = grid1(i,jj)
2480  enddo
2481  enddo
2482  endif
2483  ENDIF
2484 !
2485 ! TEMPERATURE.
2486  IF (iget(092)>0) THEN
2487 !$omp parallel do private(i,j)
2488  DO j=jsta,jend
2489  DO i=1,im
2490  grid1(i,j) = tbnd(i,j,1)
2491  ENDDO
2492  ENDDO
2493  if(grib=='grib2') then
2494  cfld=cfld+1
2495  fld_info(cfld)%ifld=iavblfld(iget(092))
2496  fld_info(cfld)%lvl=lvlsxml(1,iget(092))
2497 !$omp parallel do private(i,j,jj)
2498  do j=1,jend-jsta+1
2499  jj = jsta+j-1
2500  do i=1,im
2501  datapd(i,j,cfld) = grid1(i,jj)
2502  enddo
2503  enddo
2504  endif
2505  ENDIF
2506 !
2507 ! SPECIFIC HUMIDITY.
2508  IF (iget(093)>0) THEN
2509 !$omp parallel do private(i,j)
2510  DO j=jsta,jend
2511  DO i=1,im
2512  grid1(i,j) = qbnd(i,j,1)
2513  ENDDO
2514  ENDDO
2515  CALL bound(grid1,h1m12,h99999)
2516  if(grib=='grib2') then
2517  cfld=cfld+1
2518  fld_info(cfld)%ifld=iavblfld(iget(093))
2519  fld_info(cfld)%lvl=lvlsxml(1,iget(093))
2520 !$omp parallel do private(i,j,jj)
2521  do j=1,jend-jsta+1
2522  jj = jsta+j-1
2523  do i=1,im
2524  datapd(i,j,cfld) = grid1(i,jj)
2525  enddo
2526  enddo
2527  endif
2528  ENDIF
2529 !
2530 ! RELATIVE HUMIDITY.
2531  IF (iget(094)>0) THEN
2532 !$omp parallel do private(i,j)
2533  DO j=jsta,jend
2534  DO i=1,im
2535  grid1(i,j) = rhbnd(i,j,1)
2536  ENDDO
2537  ENDDO
2538  CALL sclfld(grid1,h100,im,jm)
2539  CALL bound(grid1,h1,h100)
2540  if(grib=='grib2') then
2541  cfld=cfld+1
2542  fld_info(cfld)%ifld=iavblfld(iget(094))
2543  fld_info(cfld)%lvl=lvlsxml(1,iget(094))
2544 !$omp parallel do private(i,j,jj)
2545  do j=1,jend-jsta+1
2546  jj = jsta+j-1
2547  do i=1,im
2548  datapd(i,j,cfld) = grid1(i,jj)
2549  enddo
2550  enddo
2551  endif
2552  ENDIF
2553 !
2554 ! U AND/OR V WIND.
2555  IF ((iget(095)>0).OR.(iget(096)>0)) THEN
2556 !$omp parallel do private(i,j)
2557  DO j=jsta,jend
2558  DO i=1,im
2559  grid1(i,j) = ubnd(i,j,1)
2560  grid2(i,j) = vbnd(i,j,1)
2561  ENDDO
2562  ENDDO
2563  IF (iget(095)>0) then
2564  if(grib=='grib2') then
2565  cfld=cfld+1
2566  fld_info(cfld)%ifld=iavblfld(iget(095))
2567  fld_info(cfld)%lvl=lvlsxml(1,iget(095))
2568 !$omp parallel do private(i,j,jj)
2569  do j=1,jend-jsta+1
2570  jj = jsta+j-1
2571  do i=1,im
2572  datapd(i,j,cfld) = grid1(i,jj)
2573  enddo
2574  enddo
2575  endif
2576  ENDIF
2577  IF (iget(096)>0) then
2578  if(grib=='grib2') then
2579  cfld=cfld+1
2580  fld_info(cfld)%ifld=iavblfld(iget(096))
2581  fld_info(cfld)%lvl=lvlsxml(1,iget(096))
2582 !$omp parallel do private(i,j,jj)
2583  do j=1,jend-jsta+1
2584  jj = jsta+j-1
2585  do i=1,im
2586  datapd(i,j,cfld) = grid2(i,jj)
2587  enddo
2588  enddo
2589  endif
2590  ENDIF
2591  ENDIF
2592  ENDIF
2593  ENDIF
2594 !
2595 ! ENDIF FOR BOUNDARY LAYER BLOCK.
2596 !
2597  ENDIF
2598 !
2599 !
2600 ! ***BLOCK 6: MISCELLANEOUS LAYER MEAN LFM AND NGM FIELDS.
2601 !
2602  IF ( (iget(066)>0).OR.(iget(081)>0).OR. &
2603  (iget(082)>0).OR.(iget(104)>0).OR. &
2604  (iget(099)>0).OR.(iget(100)>0).OR. &
2605  (iget(101)>0).OR.(iget(102)>0).OR. &
2606  (iget(103)>0) ) THEN
2607 !
2608 ! LFM "MEAN" RELATIVE HUMIDITIES AND PRECIPITABLE WATER.
2609 !
2610  IF ( (iget(066)>0).OR.(iget(081)>0).OR. &
2611  (iget(082)>0).OR.(iget(104)>0) ) THEN
2612  allocate(rh3310(im,jsta:jend),rh6610(im,jsta:jend), &
2613  rh3366(im,jsta:jend),pw3310(im,jsta:jend))
2614  CALL lfmfld(rh3310,rh6610,rh3366,pw3310)
2615 !
2616 ! SIGMA 0.33-1.00 MEAN RELATIVE HUMIIDITY.
2617  IF (iget(066)>0) THEN
2618 !$omp parallel do private(i,j)
2619  DO j=jsta,jend
2620  DO i=1,im
2621  grid1(i,j) = rh3310(i,j)
2622  ENDDO
2623  ENDDO
2624  CALL sclfld(grid1,h100,im,jm)
2625  CALL bound(grid1,h1,h100)
2626  if(grib=='grib2') then
2627  cfld=cfld+1
2628  fld_info(cfld)%ifld=iavblfld(iget(066))
2629  fld_info(cfld)%lvl=lvlsxml(1,iget(066))
2630 !$omp parallel do private(i,j,jj)
2631  do j=1,jend-jsta+1
2632  jj = jsta+j-1
2633  do i=1,im
2634  datapd(i,j,cfld) = grid1(i,jj)
2635  enddo
2636  enddo
2637 ! print *,'in miscln,RH0.33-1.0,cfld=',cfld,'fld=', &
2638 ! IAVBLFLD(IGET(066)),'lvl=',fld_info(cfld)%lvl
2639  endif
2640  ENDIF
2641 !
2642 ! SIGMA 0.66-1.00 MEAN RELATIVE HUMIIDITY.
2643  IF (iget(081)>0) THEN
2644 !$omp parallel do private(i,j)
2645  DO j=jsta,jend
2646  DO i=1,im
2647  grid1(i,j) = rh6610(i,j)
2648  ENDDO
2649  ENDDO
2650  CALL sclfld(grid1,h100,im,jm)
2651  CALL bound(grid1,h1,h100)
2652  if(grib=='grib2') then
2653  cfld=cfld+1
2654  fld_info(cfld)%ifld=iavblfld(iget(081))
2655  fld_info(cfld)%lvl=lvlsxml(1,iget(081))
2656 !$omp parallel do private(i,j,jj)
2657  do j=1,jend-jsta+1
2658  jj = jsta+j-1
2659  do i=1,im
2660  datapd(i,j,cfld) = grid1(i,jj)
2661  enddo
2662  enddo
2663  endif
2664  ENDIF
2665 !
2666 ! SIGMA 0.33-0.66 MEAN RELATIVE HUMIIDITY.
2667  IF (iget(082)>0) THEN
2668 !$omp parallel do private(i,j)
2669  DO j=jsta,jend
2670  DO i=1,im
2671  grid1(i,j) = rh3366(i,j)
2672  ENDDO
2673  ENDDO
2674  CALL sclfld(grid1,h100,im,jm)
2675  CALL bound(grid1,h1,h100)
2676  if(grib=='grib2') then
2677  cfld=cfld+1
2678  fld_info(cfld)%ifld=iavblfld(iget(082))
2679  fld_info(cfld)%lvl=lvlsxml(1,iget(082))
2680 !$omp parallel do private(i,j,jj)
2681  do j=1,jend-jsta+1
2682  jj = jsta+j-1
2683  do i=1,im
2684  datapd(i,j,cfld) = grid1(i,jj)
2685  enddo
2686  enddo
2687  endif
2688  ENDIF
2689 !
2690 ! SIGMA 0.33-1.00 PRECIPITABLE WATER.
2691  IF (iget(104)>0) THEN
2692 !$omp parallel do private(i,j)
2693  DO j=jsta,jend
2694  DO i=1,im
2695  grid1(i,j) = pw3310(i,j)
2696  ENDDO
2697  ENDDO
2698  CALL bound(grid1,d00,h99999)
2699  if(grib=='grib2') then
2700  cfld=cfld+1
2701  fld_info(cfld)%ifld=iavblfld(iget(104))
2702  fld_info(cfld)%lvl=lvlsxml(1,iget(104))
2703 !$omp parallel do private(i,j,jj)
2704  do j=1,jend-jsta+1
2705  jj = jsta+j-1
2706  do i=1,im
2707  datapd(i,j,cfld) = grid1(i,jj)
2708  enddo
2709  enddo
2710  endif
2711  ENDIF
2712  deallocate(rh3310,rh6610,rh3366,pw3310)
2713  ENDIF
2714 !
2715 ! VARIOUS LAYER MEAN NGM SIGMA FIELDS.
2716 !
2717  IF ( (iget(099)>0).OR.(iget(100)>0).OR. &
2718  (iget(101)>0).OR.(iget(102)>0).OR. &
2719  (iget(103)>0) ) THEN
2720  allocate(rh4710(im,jsta_2l:jend_2u),rh4796(im,jsta_2l:jend_2u), &
2721  rh1847(im,jsta_2l:jend_2u))
2722  allocate(rh8498(im,jsta_2l:jend_2u),qm8510(im,jsta_2l:jend_2u))
2723 
2724  CALL ngmfld(rh4710,rh4796,rh1847,rh8498,qm8510)
2725 !
2726 ! SIGMA 0.47191-1.00000 RELATIVE HUMIDITY.
2727  IF (iget(099)>0) THEN
2728 !$omp parallel do private(i,j)
2729  DO j=jsta,jend
2730  DO i=1,im
2731  grid1(i,j) = rh4710(i,j)
2732  ENDDO
2733  ENDDO
2734  CALL sclfld(grid1,h100,im,jm)
2735  CALL bound(grid1,h1,h100)
2736  if(grib=='grib2') then
2737  cfld=cfld+1
2738  fld_info(cfld)%ifld=iavblfld(iget(099))
2739  fld_info(cfld)%lvl=lvlsxml(1,iget(099))
2740 !$omp parallel do private(i,j,jj)
2741  do j=1,jend-jsta+1
2742  jj = jsta+j-1
2743  do i=1,im
2744  datapd(i,j,cfld) = grid1(i,jj)
2745  enddo
2746  enddo
2747  endif
2748  ENDIF
2749 !
2750 ! SIGMA 0.47191-0.96470 RELATIVE HUMIDITY.
2751  IF (iget(100)>0) THEN
2752 !$omp parallel do private(i,j)
2753  DO j=jsta,jend
2754  DO i=1,im
2755  grid1(i,j) = rh4796(i,j)
2756  ENDDO
2757  ENDDO
2758  CALL sclfld(grid1,h100,im,jm)
2759  CALL bound(grid1,h1,h100)
2760  if(grib=='grib2') then
2761  cfld=cfld+1
2762  fld_info(cfld)%ifld=iavblfld(iget(100))
2763  fld_info(cfld)%lvl=lvlsxml(1,iget(100))
2764 !$omp parallel do private(i,j,jj)
2765  do j=1,jend-jsta+1
2766  jj = jsta+j-1
2767  do i=1,im
2768  datapd(i,j,cfld) = grid1(i,jj)
2769  enddo
2770  enddo
2771  endif
2772  ENDIF
2773 !
2774 ! SIGMA 0.18019-0.47191 RELATIVE HUMIDITY.
2775  IF (iget(101)>0) THEN
2776 !$omp parallel do private(i,j)
2777  DO j=jsta,jend
2778  DO i=1,im
2779  grid1(i,j) = rh1847(i,j)
2780  ENDDO
2781  ENDDO
2782  CALL sclfld(grid1,h100,im,jm)
2783  CALL bound(grid1,h1,h100)
2784  if(grib=='grib2') then
2785  cfld=cfld+1
2786  fld_info(cfld)%ifld=iavblfld(iget(101))
2787  fld_info(cfld)%lvl=lvlsxml(1,iget(101))
2788 !$omp parallel do private(i,j,jj)
2789  do j=1,jend-jsta+1
2790  jj = jsta+j-1
2791  do i=1,im
2792  datapd(i,j,cfld) = grid1(i,jj)
2793  enddo
2794  enddo
2795  endif
2796  ENDIF
2797 !
2798 ! SIGMA 0.84368-0.98230 RELATIVE HUMIDITY.
2799  IF (iget(102)>0) THEN
2800 !$omp parallel do private(i,j)
2801  DO j=jsta,jend
2802  DO i=1,im
2803  grid1(i,j) = rh8498(i,j)
2804  ENDDO
2805  ENDDO
2806  CALL sclfld(grid1,h100,im,jm)
2807  CALL bound(grid1,h1,h100)
2808  if(grib=='grib2') then
2809  cfld=cfld+1
2810  fld_info(cfld)%ifld=iavblfld(iget(102))
2811  fld_info(cfld)%lvl=lvlsxml(1,iget(102))
2812 !$omp parallel do private(i,j,jj)
2813  do j=1,jend-jsta+1
2814  jj = jsta+j-1
2815  do i=1,im
2816  datapd(i,j,cfld) = grid1(i,jj)
2817  enddo
2818  enddo
2819  endif
2820  ENDIF
2821 !
2822 ! SIGMA 0.85000-1.00000 MOISTURE CONVERGENCE.
2823  IF (iget(103)>0) THEN
2824  grid1=spval
2825 ! CONVERT TO DIVERGENCE FOR GRIB
2826 !$omp parallel do private(i,j)
2827  DO j=jsta,jend
2828  DO i=1,im
2829  IF(qm8510(i,j) < spval) grid1(i,j) = -1.0*qm8510(i,j)
2830  ENDDO
2831  ENDDO
2832  if(grib=='grib2') then
2833  cfld=cfld+1
2834  fld_info(cfld)%ifld=iavblfld(iget(103))
2835  fld_info(cfld)%lvl=lvlsxml(1,iget(103))
2836 !$omp parallel do private(i,j,jj)
2837  do j=1,jend-jsta+1
2838  jj = jsta+j-1
2839  do i=1,im
2840  datapd(i,j,cfld) = grid1(i,jj)
2841  enddo
2842  enddo
2843  endif
2844  ENDIF
2845  deallocate(rh4710,rh4796,rh1847)
2846  deallocate(rh8498,qm8510)
2847  ENDIF
2848  ENDIF
2849 
2850  IF ( (iget(318)>0).OR.(iget(319)>0).OR. &
2851  (iget(320)>0))THEN
2852  allocate(rh4410(im,jsta:jend),rh7294(im,jsta:jend), &
2853  rh4472(im,jsta:jend),rh3310(im,jsta:jend))
2854  CALL lfmfld_gfs(rh4410,rh7294,rh4472,rh3310)
2855 !
2856 ! SIGMA 0.44-1.00 MEAN RELATIVE HUMIIDITY.
2857  IF (iget(318)>0) THEN
2858  grid1=spval
2859 !$omp parallel do private(i,j)
2860  DO j=jsta,jend
2861  DO i=1,im
2862  IF(rh4410(i,j) < spval) grid1(i,j) = rh4410(i,j)*100.
2863  ENDDO
2864  ENDDO
2865  CALL bound(grid1,d00,h100)
2866  if(grib=='grib2') then
2867  cfld=cfld+1
2868  fld_info(cfld)%ifld=iavblfld(iget(318))
2869  fld_info(cfld)%lvl=lvlsxml(1,iget(318))
2870 !$omp parallel do private(i,j,jj)
2871  do j=1,jend-jsta+1
2872  jj = jsta+j-1
2873  do i=1,im
2874  datapd(i,j,cfld) = grid1(i,jj)
2875  enddo
2876  enddo
2877  endif
2878  ENDIF
2879 !
2880 ! SIGMA 0.72-0.94 MEAN RELATIVE HUMIIDITY.
2881  IF (iget(319)>0) THEN
2882  grid1=spval
2883 !$omp parallel do private(i,j)
2884  DO j=jsta,jend
2885  DO i=1,im
2886  IF(rh7294(i,j) < spval) grid1(i,j) = rh7294(i,j)*100.
2887  ENDDO
2888  ENDDO
2889  CALL bound(grid1,d00,h100)
2890  if(grib=='grib2') then
2891  cfld=cfld+1
2892  fld_info(cfld)%ifld=iavblfld(iget(319))
2893  fld_info(cfld)%lvl=lvlsxml(1,iget(319))
2894 !$omp parallel do private(i,j,jj)
2895  do j=1,jend-jsta+1
2896  jj = jsta+j-1
2897  do i=1,im
2898  datapd(i,j,cfld) = grid1(i,jj)
2899  enddo
2900  enddo
2901  endif
2902  ENDIF
2903 !
2904 ! SIGMA 0.44-0.72 MEAN RELATIVE HUMIIDITY.
2905  IF (iget(320)>0) THEN
2906  grid1=spval
2907 !$omp parallel do private(i,j)
2908  DO j=jsta,jend
2909  DO i=1,im
2910  IF(rh4472(i,j) < spval) grid1(i,j)=rh4472(i,j)*100.
2911  ENDDO
2912  ENDDO
2913  CALL bound(grid1,d00,h100)
2914  if(grib=='grib2') then
2915  cfld=cfld+1
2916  fld_info(cfld)%ifld=iavblfld(iget(320))
2917  fld_info(cfld)%lvl=lvlsxml(1,iget(320))
2918 !$omp parallel do private(i,j,jj)
2919  do j=1,jend-jsta+1
2920  jj = jsta+j-1
2921  do i=1,im
2922  datapd(i,j,cfld) = grid1(i,jj)
2923  enddo
2924  enddo
2925  endif
2926  ENDIF
2927  deallocate(rh4410,rh7294,rh4472,rh3310)
2928  END IF
2929 
2930 ! GFS computes sigma=0.9950 T, THETA, U, V from lowest two model level fields
2931  IF ( (iget(321)>0).OR.(iget(322)>0).OR. &
2932  (iget(323)>0).OR.(iget(324)>0).OR. &
2933  (iget(325)>0).OR.(iget(326)>0)) THEN
2934 !$omp parallel do private(i,j)
2935  DO j=jsta,jend
2936  DO i=1,im
2937  egrid2(i,j) = 0.995*pint(i,j,lm+1)
2938  egrid1(i,j) = log(pmid(i,j,lm)/egrid2(i,j)) &
2939  / log(pmid(i,j,lm)/pmid(i,j,lm-1))
2940 
2941  IF (modelname == 'GFS' .OR. modelname == 'FV3R') THEN
2942  egrid1(i,j) = log(pmid(i,j,lm)/egrid2(i,j)) &
2943  / max(1.e-6,log(pmid(i,j,lm)/pmid(i,j,lm-1)))
2944  egrid1(i,j) =max(-10.0,min(egrid1(i,j), 10.0))
2945  IF ( abs(pmid(i,j,lm)-pmid(i,j,lm-1)) < 0.5 ) THEN
2946  egrid1(i,j) = -1.
2947  ENDIF
2948  ENDIF
2949 
2950  END DO
2951  END DO
2952 ! Temperature
2953  IF (iget(321)>0) THEN
2954  grid1=spval
2955 !$omp parallel do private(i,j)
2956  DO j=jsta,jend
2957  DO i=1,im
2958  IF(t(i,j,lm)<spval.and.t(i,j,lm-1)<spval.and.egrid1(i,j)<spval)&
2959  grid1(i,j) = t(i,j,lm)+(t(i,j,lm-1)-t(i,j,lm)) &
2960  * egrid1(i,j)
2961  ENDDO
2962  ENDDO
2963  if(grib=='grib2') then
2964  cfld=cfld+1
2965  fld_info(cfld)%ifld=iavblfld(iget(321))
2966  fld_info(cfld)%lvl=lvlsxml(1,iget(321))
2967 !$omp parallel do private(i,j,jj)
2968  do j=1,jend-jsta+1
2969  jj = jsta+j-1
2970  do i=1,im
2971  datapd(i,j,cfld) = grid1(i,jj)
2972  enddo
2973  enddo
2974  endif
2975 ! print *,'in miscln,iget(321,temp,sigmadata=',maxval(GRID1(1:im,jsta:jend)), &
2976 ! minval(GRID1(1:im,jsta:jend)),'grib=',grib
2977  ENDIF
2978 ! Potential Temperature
2979  IF (iget(322)>0) THEN
2980  grid2=spval
2981 !$omp parallel do private(i,j)
2982  DO j=jsta,jend
2983  DO i=1,im
2984  IF(t(i,j,lm)<spval.and.t(i,j,lm-1)<spval.and.egrid1(i,j)<spval)&
2985  grid2(i,j) = t(i,j,lm)+(t(i,j,lm-1)-t(i,j,lm)) &
2986  * egrid1(i,j)
2987  ENDDO
2988  ENDDO
2989  CALL calpot(egrid2,grid2(1,jsta),grid1(1,jsta))
2990  if(grib=='grib2') then
2991  cfld=cfld+1
2992  fld_info(cfld)%ifld=iavblfld(iget(322))
2993  fld_info(cfld)%lvl=lvlsxml(1,iget(322))
2994 !$omp parallel do private(i,j,jj)
2995  do j=1,jend-jsta+1
2996  jj = jsta+j-1
2997  do i=1,im
2998  datapd(i,j,cfld) = grid1(i,jj)
2999  enddo
3000  enddo
3001  endif
3002  ENDIF
3003 ! RH
3004  IF (iget(323)>0) THEN
3005  grid1=spval
3006 !$omp parallel do private(i,j,es1,qs1,rh1,es2,qs2,rh2)
3007  DO j=jsta,jend
3008  DO i=1,im
3009  IF(pmid(i,j,lm)<spval.and.pmid(i,j,lm-1)<spval.and.&
3010  q(i,j,lm)<spval.and.q(i,j,lm-1)<spval) THEN
3011  es1 = min(pmid(i,j,lm),fpvsnew(t(i,j,lm)))
3012  qs1 = con_eps*es1/(pmid(i,j,lm)+con_epsm1*es1)
3013  rh1 = q(i,j,lm)/qs1
3014  es2 = min(pmid(i,j,lm-1),fpvsnew(t(i,j,lm-1)))
3015  qs2 = con_eps*es2/(pmid(i,j,lm-1)+con_epsm1*es2)
3016  rh2 = q(i,j,lm-1)/qs2
3017  grid1(i,j) = (rh1+(rh2-rh1)*egrid1(i,j))*100.
3018  ENDIF
3019  ENDDO
3020  ENDDO
3021  CALL bound(grid1,d00,h100)
3022  if(grib=='grib2') then
3023  cfld=cfld+1
3024  fld_info(cfld)%ifld=iavblfld(iget(323))
3025  fld_info(cfld)%lvl=lvlsxml(1,iget(323))
3026 !$omp parallel do private(i,j,jj)
3027  do j=1,jend-jsta+1
3028  jj = jsta+j-1
3029  do i=1,im
3030  datapd(i,j,cfld) = grid1(i,jj)
3031  enddo
3032  enddo
3033  endif
3034  ENDIF
3035 ! U
3036  IF (iget(324)>0) THEN
3037  grid1=spval
3038 !$omp parallel do private(i,j)
3039  DO j=jsta,jend
3040  DO i=1,im
3041  IF(uh(i,j,lm)<spval.and.uh(i,j,lm-1)<spval.and.egrid1(i,j)<spval)&
3042  grid1(i,j) = uh(i,j,lm)+(uh(i,j,lm-1)-uh(i,j,lm)) &
3043  * egrid1(i,j)
3044  ENDDO
3045  ENDDO
3046  if(grib=='grib2') then
3047  cfld=cfld+1
3048  fld_info(cfld)%ifld=iavblfld(iget(324))
3049  fld_info(cfld)%lvl=lvlsxml(1,iget(324))
3050 !$omp parallel do private(i,j,jj)
3051  do j=1,jend-jsta+1
3052  jj = jsta+j-1
3053  do i=1,im
3054  datapd(i,j,cfld) = grid1(i,jj)
3055  enddo
3056  enddo
3057  endif
3058  ENDIF
3059 ! V
3060  IF (iget(325)>0) THEN
3061  grid1=spval
3062 !$omp parallel do private(i,j)
3063  DO j=jsta,jend
3064  DO i=1,im
3065  IF(vh(i,j,lm)<spval.and.vh(i,j,lm-1)<spval.and.egrid1(i,j)<spval)&
3066  grid1(i,j) = vh(i,j,lm)+(vh(i,j,lm-1)-vh(i,j,lm)) &
3067  * egrid1(i,j)
3068  ENDDO
3069  ENDDO
3070  if(grib=='grib2') then
3071  cfld=cfld+1
3072  fld_info(cfld)%ifld=iavblfld(iget(325))
3073  fld_info(cfld)%lvl=lvlsxml(1,iget(325))
3074 !$omp parallel do private(i,j,jj)
3075  do j=1,jend-jsta+1
3076  jj = jsta+j-1
3077  do i=1,im
3078  datapd(i,j,cfld) = grid1(i,jj)
3079  enddo
3080  enddo
3081  endif
3082  ENDIF
3083 ! OMGA
3084  IF (iget(326)>0) THEN
3085  grid1=spval
3086 !$omp parallel do private(i,j)
3087  DO j=jsta,jend
3088  DO i=1,im
3089  IF(omga(i,j,lm)<spval.and.omga(i,j,lm-1)<spval.and.egrid1(i,j)<spval)&
3090  grid1(i,j) = omga(i,j,lm)+(omga(i,j,lm-1)-omga(i,j,lm))&
3091  * egrid1(i,j)
3092  ENDDO
3093  ENDDO
3094  if(grib=='grib2') then
3095  cfld=cfld+1
3096  fld_info(cfld)%ifld=iavblfld(iget(326))
3097  fld_info(cfld)%lvl=lvlsxml(1,iget(326))
3098 !$omp parallel do private(i,j,jj)
3099  do j=1,jend-jsta+1
3100  jj = jsta+j-1
3101  do i=1,im
3102  datapd(i,j,cfld) = grid1(i,jj)
3103  enddo
3104  enddo
3105  endif
3106  ENDIF
3107  END IF
3108 
3109 ! MIXED LAYER CAPE AND CINS
3110 !
3111  field1=.false.
3112  field2=.false.
3113 !
3114  IF(iget(032)>0)THEN
3115  IF(lvls(3,iget(032))>0)field1=.true.
3116  ENDIF
3117  IF(iget(107)>0)THEN
3118  IF(lvls(3,iget(107))>0)field2=.true.
3119  ENDIF
3120 !
3121  IF(iget(582)>0)THEN
3122  field1=.true.
3123  ENDIF
3124  IF(iget(583)>0)THEN
3125  field2=.true.
3126  ENDIF
3127 
3128  IF(field1.OR.field2)THEN
3129  itype = 2
3130 !
3131 !$omp parallel do private(i,j)
3132  DO j=jsta,jend
3133  DO i=1,im
3134  egrid1(i,j) = -h99999
3135  egrid2(i,j) = -h99999
3136  lb2(i,j) = (lvlbnd(i,j,1) + lvlbnd(i,j,2) + &
3137  lvlbnd(i,j,3))/3
3138  p1d(i,j) = (pbnd(i,j,1) + pbnd(i,j,2) + pbnd(i,j,3))/3
3139  t1d(i,j) = (tbnd(i,j,1) + tbnd(i,j,2) + tbnd(i,j,3))/3
3140  q1d(i,j) = (qbnd(i,j,1) + qbnd(i,j,2) + qbnd(i,j,3))/3
3141  ENDDO
3142  ENDDO
3143 !
3144  dpbnd = 0.
3145  CALL calcape(itype,dpbnd,p1d,t1d,q1d,lb2,egrid1, &
3146  egrid2,egrid3,egrid4,egrid5)
3147 
3148  IF (iget(582)>0) THEN
3149 ! dong add missing value for cape
3150  grid1=spval
3151 !$omp parallel do private(i,j)
3152  DO j=jsta,jend
3153  DO i=1,im
3154  IF(t1d(i,j) < spval) THEN
3155  grid1(i,j) = egrid1(i,j)
3156  IF (submodelname == 'RTMA')mlcape(i,j)=grid1(i,j)
3157  ENDIF
3158  ENDDO
3159  ENDDO
3160  CALL bound(grid1,d00,h99999)
3161  if(grib=='grib2') then
3162  cfld=cfld+1
3163  fld_info(cfld)%ifld=iavblfld(iget(582))
3164  fld_info(cfld)%lvl=lvlsxml(1,iget(582))
3165 !$omp parallel do private(i,j,jj)
3166  do j=1,jend-jsta+1
3167  jj = jsta+j-1
3168  do i=1,im
3169  datapd(i,j,cfld) = grid1(i,jj)
3170  enddo
3171  enddo
3172  endif
3173  ENDIF
3174  IF (iget(583)>0) THEN
3175 ! dong add missing value for cape
3176  grid1=spval
3177 !$omp parallel do private(i,j)
3178  DO j=jsta,jend
3179  DO i=1,im
3180  IF(t1d(i,j) < spval) grid1(i,j) = - egrid2(i,j)
3181  ENDDO
3182  ENDDO
3183 !
3184  CALL bound(grid1,d00,h99999)
3185 !
3186 !$omp parallel do private(i,j)
3187  DO j=jsta,jend
3188  DO i=1,im
3189  IF(t1d(i,j) < spval) THEN
3190  grid1(i,j) = - grid1(i,j)
3191  IF (submodelname == 'RTMA') mlcin(i,j)=grid1(i,j)
3192  ENDIF
3193  ENDDO
3194  ENDDO
3195 !
3196  if(grib=='grib2') then
3197  cfld=cfld+1
3198  fld_info(cfld)%ifld=iavblfld(iget(583))
3199  fld_info(cfld)%lvl=lvlsxml(1,iget(583))
3200 !$omp parallel do private(i,j,jj)
3201  do j=1,jend-jsta+1
3202  jj = jsta+j-1
3203  do i=1,im
3204  datapd(i,j,cfld) = grid1(i,jj)
3205  enddo
3206  enddo
3207  endif
3208  ENDIF
3209  ENDIF
3210 
3211 ! MIXED LAYER LIFTING CONDENSATION PRESSURE AND HEIGHT.
3212 ! EGRID1 IS LCL PRESSURE. EGRID2 IS LCL HEIGHT.
3213 !
3214  IF ( (iget(109)>0).OR.(iget(110)>0) ) THEN
3215  CALL callcl(p1d,t1d,q1d,egrid1,egrid2)
3216  IF (iget(109)>0) THEN
3217  grid1=spval
3218  DO j=jsta,jend
3219  DO i=1,im
3220  IF(t1d(i,j) < spval) grid1(i,j)=egrid2(i,j)
3221  IF (submodelname == 'RTMA') mllcl(i,j) = grid1(i,j)
3222  ENDDO
3223  ENDDO
3224 !
3225 ! ID(1:25) = 0
3226 !
3227 ! CALL GRIBIT(IGET(109),1,
3228 ! X GRID1,IM,JM)
3229  ENDIF
3230 !
3231 ! IF (IGET(110)>0) THEN
3232 ! DO J=JSTA,JEND
3233 ! DO I=1,IM
3234 ! GRID1(I,J)=EGRID1(I,J)
3235 ! ENDDO
3236 ! ENDDO
3237 !
3238 ! ID(1:25) = 0
3239 !
3240 ! CALL GRIBIT(IGET(110),1,
3241 ! X GRID1,IM,JM)
3242 ! ENDIF
3243  ENDIF
3244 !
3245 ! MOST UNSTABLE CAPE-LOWEST 300 MB
3246 !
3247  field1=.false.
3248  field2=.false.
3249 !
3250  IF(iget(032)>0)THEN
3251  IF(lvls(4,iget(032))>0)field1=.true.
3252  ENDIF
3253 
3254  IF(iget(107)>0)THEN
3255  IF(lvls(4,iget(107))>0)field2=.true.
3256  ENDIF
3257 !
3258  IF(iget(584)>0)THEN
3259  field1=.true.
3260  ENDIF
3261  IF(iget(585)>0)THEN
3262  field2=.true.
3263  ENDIF
3264 !
3265  IF(field1.OR.field2)THEN
3266  itype = 1
3267 !
3268 !$omp parallel do private(i,j)
3269  DO j=jsta,jend
3270  DO i=1,im
3271  egrid1(i,j) = -h99999
3272  egrid2(i,j) = -h99999
3273  ENDDO
3274  ENDDO
3275 
3276  dpbnd = 300.e2
3277  CALL calcape(itype,dpbnd,p1d,t1d,q1d,lb2,egrid1, &
3278  egrid2,egrid3,egrid4,egrid5)
3279  IF (submodelname == 'RTMA') mumixr(i,j) = q1d(i,j)
3280  IF (iget(584)>0) THEN
3281 ! dong add missing value to cin
3282  grid1 = spval
3283 !$omp parallel do private(i,j)
3284  DO j=jsta,jend
3285  DO i=1,im
3286  IF(t1d(i,j) < spval) THEN
3287  grid1(i,j) = egrid1(i,j)
3288  IF (submodelname == 'RTMA') mucape(i,j)=grid1(i,j)
3289  ENDIF
3290  ENDDO
3291  ENDDO
3292  CALL bound(grid1,d00,h99999)
3293 ! IF (SUBMODELNAME == 'RTMA') THEN
3294 ! CALL BOUND(MUCAPE,D00,H99999)
3295 ! ENDIF
3296  if(grib=='grib2') then
3297  cfld=cfld+1
3298  fld_info(cfld)%ifld=iavblfld(iget(584))
3299  fld_info(cfld)%lvl=lvlsxml(1,iget(584))
3300 !$omp parallel do private(i,j,jj)
3301  do j=1,jend-jsta+1
3302  jj = jsta+j-1
3303  do i=1,im
3304  datapd(i,j,cfld) = grid1(i,jj)
3305  enddo
3306  enddo
3307  endif
3308 
3309  ENDIF
3310 
3311  IF (iget(585)>0) THEN
3312 ! dong add missing value to cin
3313  grid1 = spval
3314 !$omp parallel do private(i,j)
3315  DO j=jsta,jend
3316  DO i=1,im
3317  IF(t1d(i,j) < spval) grid1(i,j) = - egrid2(i,j)
3318  ENDDO
3319  ENDDO
3320  CALL bound(grid1,d00,h99999)
3321  DO j=jsta,jend
3322  DO i=1,im
3323  IF(t1d(i,j) < spval) THEN
3324  grid1(i,j) = - grid1(i,j)
3325  IF (submodelname == 'RTMA')THEN
3326  mucape(i,j) = grid1(i,j)
3327  muq1d(i,j) = q1d(i,j)
3328  ENDIF
3329  ENDIF
3330  ENDDO
3331  ENDDO
3332  if(grib=='grib2') then
3333  cfld=cfld+1
3334  fld_info(cfld)%ifld=iavblfld(iget(585))
3335  fld_info(cfld)%lvl=lvlsxml(1,iget(585))
3336 !$omp parallel do private(i,j,jj)
3337  do j=1,jend-jsta+1
3338  jj = jsta+j-1
3339  do i=1,im
3340  datapd(i,j,cfld) = grid1(i,jj)
3341  enddo
3342  enddo
3343  endif
3344 
3345  ENDIF
3346 
3347 ! EQUILLIBRIUM HEIGHT
3348  IF (iget(443)>0) THEN
3349  grid1 = spval
3350 !$omp parallel do private(i,j)
3351  DO j=jsta,jend
3352  DO i=1,im
3353  IF(t1d(i,j) < spval) grid1(i,j) = egrid4(i,j)
3354  ENDDO
3355  ENDDO
3356  if(grib=='grib2') then
3357  cfld=cfld+1
3358  fld_info(cfld)%ifld=iavblfld(iget(443))
3359  fld_info(cfld)%lvl=lvlsxml(1,iget(443))
3360 !$omp parallel do private(i,j,jj)
3361  do j=1,jend-jsta+1
3362  jj = jsta+j-1
3363  do i=1,im
3364  datapd(i,j,cfld) = grid1(i,jj)
3365  enddo
3366  enddo
3367  endif
3368  ENDIF
3369 !Equilibrium Temperature
3370  IF (iget(982)>0) THEN
3371  DO j=jsta,jend
3372  DO i=1,im
3373  grid1(i,j) = teql(i,j)
3374  ENDDO
3375  ENDDO
3376  if(grib=='grib2') then
3377  cfld=cfld+1
3378  fld_info(cfld)%ifld=iavblfld(iget(982))
3379  fld_info(cfld)%lvl=lvlsxml(1,iget(982))
3380 !$omp parallel do private(i,j,jj)
3381  do j=1,jend-jsta+1
3382  jj = jsta+j-1
3383  do i=1,im
3384  datapd(i,j,cfld) = grid1(i,jj)
3385  enddo
3386  enddo
3387  endif
3388  ENDIF
3389 
3390 
3391 ! PRESSURE OF LEVEL FROM WHICH 300 MB MOST UNSTABLE CAPE
3392 ! PARCEL WAS LIFTED (eq. PRESSURE OF LEVEL OF HIGHEST THETA-E)
3393  IF (iget(246)>0) THEN
3394  grid1 = spval
3395 !$omp parallel do private(i,j)
3396  DO j=jsta,jend
3397  DO i=1,im
3398  IF(t1d(i,j) < spval) grid1(i,j) = egrid3(i,j)
3399  ENDDO
3400  ENDDO
3401  CALL bound(grid1,d00,h99999)
3402 ! print *,'in miscln,PLPL=',maxval(grid1(1:im,jsta:jend)), &
3403 ! minval(grid1(1:im,jsta:jend))
3404  if(grib=='grib2') then
3405  cfld=cfld+1
3406  fld_info(cfld)%ifld=iavblfld(iget(246))
3407  fld_info(cfld)%lvl=lvlsxml(1,iget(246))
3408 !$omp parallel do private(i,j,jj)
3409  do j=1,jend-jsta+1
3410  jj = jsta+j-1
3411  do i=1,im
3412  datapd(i,j,cfld) = grid1(i,jj)
3413  enddo
3414  enddo
3415  endif
3416  ENDIF ! 246
3417 
3418 ! GENERAL THUNDER PARAMETER ??? 458 ???
3419  IF (iget(444)>0) THEN
3420  grid1 = spval
3421 !$omp parallel do private(i,j)
3422  DO j=jsta,jend
3423  DO i=1,im
3424  IF(cprate(i,j) < spval) THEN
3425  IF (cprate(i,j) > pthresh) THEN
3426  grid1(i,j) = egrid5(i,j)
3427  ELSE
3428  grid1(i,j) = 0
3429  ENDIF
3430  ENDIF
3431  ENDDO
3432  ENDDO
3433  CALL bound(grid1,d00,h99999)
3434  if(grib=='grib2') then
3435  cfld=cfld+1
3436  fld_info(cfld)%ifld=iavblfld(iget(444))
3437  fld_info(cfld)%lvl=lvlsxml(1,iget(444))
3438 !$omp parallel do private(i,j,jj)
3439  do j=1,jend-jsta+1
3440  jj = jsta+j-1
3441  do i=1,im
3442  datapd(i,j,cfld) = grid1(i,jj)
3443  enddo
3444  enddo
3445  endif
3446  ENDIF
3447  ENDIF
3448 
3449 
3450  IF (submodelname == 'RTMA')THEN
3451 
3452 !
3453 ! --- Effective (inflow) Layer (EL)
3454 !
3455 
3456  ALLOCATE(el_base(im,jsta_2l:jend_2u))
3457  ALLOCATE(el_tops(im,jsta_2l:jend_2u))
3458  ALLOCATE(found_base(im,jsta_2l:jend_2u))
3459  ALLOCATE(found_tops(im,jsta_2l:jend_2u))
3460 !$omp parallel do private(i,j)
3461  DO j=jsta,jend
3462  DO i=1,im
3463  el_base(i,j) = lm
3464  el_tops(i,j) = lm
3465  found_base(i,j) = .false.
3466  found_tops(i,j) = .false.
3467  ENDDO
3468  ENDDO
3469 
3470 !
3471 
3472  itype = 2
3473  dpbnd = 0.
3474 
3475  DO l = lm, 1, -1
3476 
3477 ! SET AIR PARCELS FOR LEVEL L
3478 !$omp parallel do private(i,j)
3479  DO j=jsta,jend
3480  DO i=1,im
3481  egrid1(i,j) = -h99999
3482  egrid2(i,j) = -h99999
3483  idummy(i,j) = 0
3484  p1d(i,j) = pmid(i,j,l)
3485  t1d(i,j) = t(i,j,l)
3486  q1d(i,j) = q(i,j,l)
3487  ENDDO
3488  ENDDO
3489 
3490 !--- CALCULATE CAPE/CIN FOR ALL AIR PARCELS on LEVEL L
3491  IF (debugprint) WRITE(1000+me,'(A,I3)') &
3492  ' CALCULATING CAPE/CINS ON LEVEL:',l
3493  CALL calcape(itype,dpbnd,p1d,t1d,q1d,idummy,egrid1, &
3494  egrid2,egrid3,egrid4,egrid5)
3495 
3496 !--- CHECK CAPE/CIN OF EACH AIR PARCELS WITH EL CRITERIA
3497 !$omp parallel do private(i,j)
3498  DO j=jsta,jend
3499  DO i=1,im
3500  IF ( .NOT. found_base(i,j) ) THEN
3501  IF ( egrid1(i,j) >= 100. .AND. egrid2(i,j) >= -250. ) THEN
3502  el_base(i,j) = l
3503  found_base(i,j) = .true.
3504  ELSE
3505  el_base(i,j) = lm
3506  found_base(i,j) = .false.
3507  END IF
3508  ELSE
3509  IF ( .NOT. found_tops(i,j) ) THEN
3510  IF ( egrid1(i,j) < 100. .OR. egrid2(i,j) < -250. ) THEN
3511  el_tops(i,j) = l + 1
3512  found_tops(i,j) = .true.
3513  ELSE
3514  el_tops(i,j) = lm
3515  found_tops(i,j) = .false.
3516  END IF
3517  END IF
3518  END IF
3519  ENDDO
3520  ENDDO
3521 
3522  END DO ! L
3523 
3524 
3525  IF (ALLOCATED(found_base)) DEALLOCATE(found_base)
3526  IF (ALLOCATED(found_tops)) DEALLOCATE(found_tops)
3527 
3528  IF (debugprint) THEN
3529  WRITE(im_ch,'(I5.5)') im
3530  WRITE(jsta_ch,'(I5.5)') jsta
3531  WRITE(jend_ch,'(I5.5)') jend
3532  effl_fname="EFFL_NEW_"//im_ch//"_"//jsta_ch//"_"//jend_ch &
3533  //".dat"
3534  effl_fname2="EFFL_NEW_LVLS_"//im_ch//"_"//jsta_ch//"_"//jend_ch &
3535  //".dat"
3536  iunit=10000+jsta
3537  iunit2=20000+jsta
3538  irec=0
3539  irec2=0
3540  OPEN(iunit,file=trim(adjustl(effl_fname)),form='FORMATTED')
3541 
3542 
3543  DO j=jsta,jend
3544  DO i=1,im
3545  irec = irec + 1
3546  irec2 = irec2 + 1
3547  WRITE(iunit,'(1x,I6,2x,I6,2(2x,I6,2x,F12.3))') i, j, &
3548  el_base(i,j),pmid(i,j,el_base(i,j)), &
3549  el_tops(i,j),pmid(i,j,el_tops(i,j))
3550  END DO
3551  ENDDO
3552  CLOSE(iunit)
3553  ENDIF
3554 
3555  IF(ALLOCATED(tpar_base)) DEALLOCATE(tpar_base)
3556  IF(ALLOCATED(tpar_tops)) DEALLOCATE(tpar_tops)
3557 
3558  ENDIF
3559 !
3560 ! EXPAND HRRR CAPE/CIN RELATED VARIABLES
3561 !
3562 ! CAPE AND CINS 0-3KM, FOLLOW ML PROCEDURE WITH HEIGHT 0-3KM
3563 !
3564  field1=.false.
3565  field2=.false.
3566 !
3567  IF(iget(032)>0)THEN
3568  IF(lvls(3,iget(032))>0)field1=.true.
3569  ENDIF
3570  IF(iget(107)>0)THEN
3571  IF(lvls(3,iget(107))>0)field2=.true.
3572  ENDIF
3573 !
3574  IF(iget(950)>0)THEN
3575  field1=.true.
3576  ENDIF
3577  IF(iget(951)>0)THEN
3578  field2=.true.
3579  ENDIF
3580 !
3581 ! IF(FIELD1)ITYPE=2
3582 ! IF(FIELD2)ITYPE=2
3583 
3584  IF(field1.OR.field2)THEN
3585  itype = 2
3586 
3587 !
3588 !$omp parallel do private(i,j)
3589  DO j=jsta,jend
3590  DO i=1,im
3591  egrid1(i,j) = -h99999
3592  egrid2(i,j) = -h99999
3593  egrid3(i,j) = -h99999
3594  egrid4(i,j) = -h99999
3595  egrid5(i,j) = -h99999
3596  egrid6(i,j) = -h99999
3597  egrid7(i,j) = -h99999
3598  egrid8(i,j) = -h99999
3599 ! ENDDO
3600 ! ENDDO
3601 ! DO J=JSTA,JEND
3602 ! DO I=1,IM
3603  lb2(i,j) = (lvlbnd(i,j,1) + lvlbnd(i,j,2) + &
3604  lvlbnd(i,j,3))/3
3605  p1d(i,j) = (pbnd(i,j,1) + pbnd(i,j,2) + pbnd(i,j,3))/3
3606  t1d(i,j) = (tbnd(i,j,1) + tbnd(i,j,2) + tbnd(i,j,3))/3
3607  q1d(i,j) = (qbnd(i,j,1) + qbnd(i,j,2) + qbnd(i,j,3))/3
3608  ENDDO
3609  ENDDO
3610 
3611  dpbnd = 0.
3612  CALL calcape2(itype,dpbnd,p1d,t1d,q1d,lb2, &
3613  egrid1,egrid2,egrid3,egrid4,egrid5, &
3614  egrid6,egrid7,egrid8)
3615 
3616 ! CAPE1, CINS2, LFC3, ESRHL4,ESRHH5,
3617 ! DCAPE6,DGLD7, ESP8)
3618 !
3619  IF (iget(950)>0) THEN
3620 ! dong add missing value for cape
3621  grid1=spval
3622 !$omp parallel do private(i,j)
3623  DO j=jsta,jend
3624  DO i=1,im
3625  IF(t1d(i,j) < spval) grid1(i,j) = egrid1(i,j)
3626  ENDDO
3627  ENDDO
3628  CALL bound(grid1,d00,h99999)
3629  if(grib=='grib2') then
3630  cfld=cfld+1
3631  fld_info(cfld)%ifld=iavblfld(iget(950))
3632  fld_info(cfld)%lvl=lvlsxml(1,iget(950))
3633 !$omp parallel do private(i,j,jj)
3634  do j=1,jend-jsta+1
3635  jj = jsta+j-1
3636  do i=1,im
3637  datapd(i,j,cfld) = grid1(i,jj)
3638  enddo
3639  enddo
3640  endif
3641  ENDIF !950
3642 !
3643  IF (iget(951)>0) THEN
3644 ! dong add missing value for cape
3645  grid1=spval
3646 !$omp parallel do private(i,j)
3647  DO j=jsta,jend
3648  DO i=1,im
3649  IF(t1d(i,j) < spval) grid1(i,j) = - egrid2(i,j)
3650  ENDDO
3651  ENDDO
3652 !
3653  CALL bound(grid1,d00,h99999)
3654 !
3655 !$omp parallel do private(i,j)
3656  DO j=jsta,jend
3657  DO i=1,im
3658  IF(t1d(i,j) < spval) grid1(i,j) = - grid1(i,j)
3659  ENDDO
3660  ENDDO
3661 !
3662  if(grib=='grib2') then
3663  cfld=cfld+1
3664  fld_info(cfld)%ifld=iavblfld(iget(951))
3665  fld_info(cfld)%lvl=lvlsxml(1,iget(951))
3666 !$omp parallel do private(i,j,jj)
3667  do j=1,jend-jsta+1
3668  jj = jsta+j-1
3669  do i=1,im
3670  datapd(i,j,cfld) = grid1(i,jj)
3671  enddo
3672  enddo
3673  endif
3674  ENDIF !951
3675 
3676 ! LFC HEIGHT
3677 
3678  IF (iget(952)>0) THEN
3679  grid1=spval
3680 !$omp parallel do private(i,j)
3681  DO j=jsta,jend
3682  DO i=1,im
3683  IF(t1d(i,j) < spval) grid1(i,j) = egrid3(i,j)
3684  ENDDO
3685  ENDDO
3686  CALL bound(grid1,d00,h99999)
3687  if(grib=='grib2') then
3688  cfld=cfld+1
3689  fld_info(cfld)%ifld=iavblfld(iget(952))
3690  fld_info(cfld)%lvl=lvlsxml(1,iget(952))
3691 !$omp parallel do private(i,j,jj)
3692  do j=1,jend-jsta+1
3693  jj = jsta+j-1
3694  do i=1,im
3695  datapd(i,j,cfld) = grid1(i,jj)
3696  enddo
3697  enddo
3698  endif
3699  ENDIF !952
3700 
3701 
3702 ! EFFECTIVE STORM RELATIVE HELICITY AND STORM MOTION.
3703 
3704  allocate(ust(im,jsta_2l:jend_2u),vst(im,jsta_2l:jend_2u), &
3705  heli(im,jsta_2l:jend_2u,2))
3706  allocate(llow(im,jsta_2l:jend_2u),lupp(im,jsta_2l:jend_2u), &
3707  cangle(im,jsta_2l:jend_2u))
3708 
3709  iget1 = iget(953)
3710  iget2 = -1
3711  iget3 = -1
3712  if (iget1 > 0) then
3713  iget2 = lvls(1,iget1)
3714  iget3 = lvls(2,iget1)
3715  endif
3716  if(me==0) write(0,*) '953 ',iget1,iget2,iget3
3717  IF (iget1 > 0 .OR. iget(162) > 0 .OR. iget(953) > 0) THEN
3718  depth(1) = 3000.0
3719  depth(2) = 1000.0
3720  IF (submodelname == 'RTMA') THEN
3721 !--- IF USSING EL BASE & TOP COMPUTED BY NEW SCHEME FOR THE
3722 !RELATED VARIABLES
3723 !$omp parallel do private(i,j)
3724  DO j=jsta,jend
3725  DO i=1,im
3726  llow(i,j) = el_base(i,j)
3727  lupp(i,j) = el_tops(i,j)
3728  ENDDO
3729  ENDDO
3730  ELSE
3731 !$omp parallel do private(i,j)
3732  DO j=jsta,jend
3733  DO i=1,im
3734  llow(i,j) = int(egrid4(i,j))
3735  lupp(i,j) = int(egrid5(i,j))
3736  ENDDO
3737  ENDDO
3738  ENDIF
3739 !--- OUTPUT EL BASE & TOP COMPUTED BY OLD SCHEME
3740  IF (debugprint) THEN
3741  WRITE(im_ch,'(I5.5)') im
3742  WRITE(jsta_ch,'(I5.5)') jsta
3743  WRITE(jend_ch,'(I5.5)') jend
3744  effl_fname="EFFL_OLD_"//im_ch//"_"//jsta_ch//"_"//jend_ch &
3745  //".dat"
3746  iunit=10000+jsta
3747  irec=0
3748  OPEN(iunit,file=trim(adjustl(effl_fname)),form='FORMATTED')
3749  DO j=jsta,jend
3750  DO i=1,im
3751  irec = irec + 1
3752 ! WRITE(IUNIT,'(1x,I6,2x,I6,2x,I6,2x,I6)')I,J,LLOW(I,J),LUPP(I,J)
3753  WRITE(iunit,'(1x,I6,2x,I6,2(2x,I6,2x,F12.3))') i, j, &
3754  llow(i,j),pmid(i,j,llow(i,j)), &
3755  lupp(i,j),pmid(i,j,lupp(i,j))
3756  END DO
3757  ENDDO
3758  CLOSE(iunit)
3759  ENDIF
3760 
3761 
3762 ! CALL CALHEL(DEPTH,UST,VST,HELI,USHR1,VSHR1,USHR6,VSHR6)
3763  CALL calhel2(llow,lupp,depth,ust,vst,heli,cangle)
3764 !
3765  IF (iget2 > 0) then
3766 !$omp parallel do private(i,j)
3767  DO j=jsta,jend
3768  DO i=1,im
3769  grid1(i,j) = heli(i,j,1)
3770  ! GRID1(I,J) = HELI(I,J,2)
3771  ENDDO
3772  ENDDO
3773  if(grib=='grib2') then
3774  cfld=cfld+1
3775  fld_info(cfld)%ifld=iavblfld(iget1)
3776  fld_info(cfld)%lvl=lvlsxml(1,iget1)
3777 !$omp parallel do private(i,j,jj)
3778  do j=1,jend-jsta+1
3779  jj = jsta+j-1
3780  do i=1,im
3781  datapd(i,j,cfld) = grid1(i,jj)
3782  enddo
3783  enddo
3784  endif
3785  ENDIF
3786 
3787  ENDIF !953
3788 
3789 
3790  IF (submodelname == 'RTMA') THEN !Start RTMA block
3791 
3792 !EL field allocation
3793 
3794  allocate(eshr(im,jsta_2l:jend_2u),uvect(im,jsta_2l:jend_2u),&
3795  vvect(im,jsta_2l:jend_2u),htsfc(im,jsta_2l:jend_2u))
3796  allocate(effust(im,jsta_2l:jend_2u),effvst(im,jsta_2l:jend_2u),&
3797  esrh(im,jsta_2l:jend_2u))
3798 
3799 !
3800  DO j=jsta,jend
3801  DO i=1,im
3802  maxthe(i,j)=-h99999
3803  the(i,j)=-h99999
3804  maxthepos(i,j)=0
3805  muq1d(i,j) = 0.
3806  ENDDO
3807  ENDDO
3808  DO l=lm,1,-1
3809 
3810  DO j=jsta,jend
3811  DO i=1,im
3812  egrid1(i,j) = -h99999
3813  p1d(i,j)=pmid(i,j,l)
3814  t1d(i,j)=t(i,j,l)
3815  q1d(i,j)=q(i,j,l)
3816  ENDDO
3817  ENDDO
3818  CALL calthte(p1d,t1d,q1d,egrid1)
3819  DO j=jsta,jend
3820  DO i=1,im
3821  the(i,j)=egrid1(i,j)
3822  IF(the(i,j)>=maxthe(i,j))THEN
3823  maxthe(i,j)=the(i,j)
3824  maxthepos(i,j)=l
3825  muq1d(i,j) = q(i,j,l) ! save the Q of air parcel with max theta-e (MU Parcel)
3826  ENDIF
3827  ENDDO
3828  ENDDO
3829 
3830  ENDDO
3831 
3832 !
3833 !get surface height
3834  IF(gridtype == 'E')THEN
3835  jvn = 1
3836  jvs = -1
3837  do j=jsta,jend
3838  ive(j) = mod(j,2)
3839  ivw(j) = ive(j)-1
3840  enddo
3841  istart = 2
3842  istop = im-1
3843  jstart = jsta_m
3844  jstop = jend_m
3845  ELSE IF(gridtype == 'B')THEN
3846  jvn = 1
3847  jvs = 0
3848  do j=jsta,jend
3849  ive(j)=1
3850  ivw(j)=0
3851  enddo
3852  istart = 2
3853  istop = im-1
3854  jstart = jsta_m
3855  jstop = jend_m
3856  ELSE
3857  jvn = 0
3858  jvs = 0
3859  do j=jsta,jend
3860  ive(j) = 0
3861  ivw(j) = 0
3862  enddo
3863  istart = 1
3864  istop = im
3865  jstart = jsta
3866  jstop = jend
3867  END IF
3868 
3869  IF(gridtype /= 'A') CALL exch(fis(1:im,jsta:jend))
3870  DO j=jstart,jstop
3871  DO i=istart,istop
3872  ie = i+ive(j)
3873  iw = i+ivw(j)
3874  jn = j+jvn
3875  js = j+jvs
3876  IF (gridtype=='B')THEN
3877  htsfc(i,j)=(0.25/g)*(fis(iw,j)+fis(ie,j)+fis(i,jn)+fis(ie,jn))
3878  ELSE
3879  htsfc(i,j)=(0.25/g)*(fis(iw,j)+fis(ie,j)+fis(i,jn)+fis(i,js))
3880  ENDIF
3881  ENDDO
3882  ENDDO
3883 
3884 !Height of effbot
3885  IF (iget(979)>0) THEN
3886  grid1=spval
3887  DO j=jsta,jend
3888  DO i=1,im
3889  IF(zint(i,j,llow(i,j))<spval.and.htsfc(i,j)<spval)&
3890  grid1(i,j) = zint(i,j,llow(i,j)) - htsfc(i,j)
3891  ENDDO
3892  ENDDO
3893  if(grib=='grib2') then
3894  cfld=cfld+1
3895  fld_info(cfld)%ifld=iavblfld(iget(979))
3896  fld_info(cfld)%lvl=lvlsxml(1,iget(979))
3897 !$omp parallel do private(i,j,jj)
3898  do j=1,jend-jsta+1
3899  jj = jsta+j-1
3900  do i=1,im
3901  datapd(i,j,cfld) = grid1(i,jj)
3902  enddo
3903  enddo
3904  endif
3905  ENDIF
3906 !Height of effbot
3907  IF (iget(980)>0) THEN
3908  grid1=spval
3909  DO j=jsta,jend
3910  DO i=1,im
3911  IF(zint(i,j,lupp(i,j))<spval.and.htsfc(i,j)<spval)&
3912  grid1(i,j) = zint(i,j,lupp(i,j)) - htsfc(i,j)
3913  ENDDO
3914  ENDDO
3915  if(grib=='grib2') then
3916  cfld=cfld+1
3917  fld_info(cfld)%ifld=iavblfld(iget(980))
3918  fld_info(cfld)%lvl=lvlsxml(1,iget(980))
3919 !$omp parallel do private(i,j,jj)
3920  do j=1,jend-jsta+1
3921  jj = jsta+j-1
3922  do i=1,im
3923  datapd(i,j,cfld) = grid1(i,jj)
3924  enddo
3925  enddo
3926  endif
3927  ENDIF
3928 
3929 !U inflow based to 50% EL shear vector
3930 
3931  IF (iget(983)>0) THEN
3932  grid1=spval
3933  DO j=jsta,jend
3934  DO i=1,im
3935  IF(llow(i,j)<spval.and.lupp(i,j)<spval) THEN
3936  midcal=int(llow(i,j)+d50*(lupp(i,j)-llow(i,j)))
3937  !mid-layer
3938  !vertical
3939  !index
3940  uvect(i,j)=uh(i,j,midcal)-uh(i,j,llow(i,j))
3941  grid1(i,j)=uvect(i,j)
3942  ENDIF
3943  ENDDO
3944  ENDDO
3945  if(grib=='grib2') then
3946  cfld=cfld+1
3947  fld_info(cfld)%ifld=iavblfld(iget(983))
3948  fld_info(cfld)%lvl=lvlsxml(1,iget(983))
3949 !$omp parallel do private(i,j,jj)
3950  do j=1,jend-jsta+1
3951  jj = jsta+j-1
3952  do i=1,im
3953  datapd(i,j,cfld) = grid1(i,jj)
3954  enddo
3955  enddo
3956  endif
3957  ENDIF
3958 
3959 !V inflow based to 50% EL shear vector
3960  IF (iget(984)>0) THEN
3961  grid1=spval
3962  DO j=jsta,jend
3963  DO i=1,im
3964  IF(llow(i,j)<spval.and.lupp(i,j)<spval.and.&
3965  vh(i,j,midcal)<spval.and.vh(i,j,llow(i,j))<spval)THEN
3966  midcal=int(llow(i,j)+d50*(ieql(i,j)-llow(i,j)))
3967  !mid-layer
3968  !vertical
3969  !index
3970  vvect(i,j)=vh(i,j,midcal)-vh(i,j,llow(i,j))
3971  grid1(i,j)=vvect(i,j)
3972  ENDIF
3973  ENDDO
3974  ENDDO
3975  if(grib=='grib2') then
3976  cfld=cfld+1
3977  fld_info(cfld)%ifld=iavblfld(iget(984))
3978  fld_info(cfld)%lvl=lvlsxml(1,iget(984))
3979 !$omp parallel do private(i,j,jj)
3980  do j=1,jend-jsta+1
3981  jj = jsta+j-1
3982  do i=1,im
3983  datapd(i,j,cfld) = grid1(i,jj)
3984  enddo
3985  enddo
3986  endif
3987  ENDIF
3988 
3989 !Inflow based (ESFC) to (50%) EL shear magnitude
3990  IF (iget(985)>0) THEN
3991  grid1=spval
3992  DO j=jsta,jend
3993  DO i=1,im
3994  IF(uvect(i,j)<spval.and.vvect(i,j)<spval) THEN
3995  eshr(i,j)=sqrt((uvect(i,j)**2)+(vvect(i,j))**2)
3996  !effshear
3997  !calc
3998  grid1(i,j)=eshr(i,j) !Effective
3999  ENDIF
4000  ENDDO
4001  ENDDO
4002  if(grib=='grib2') then
4003  cfld=cfld+1
4004  fld_info(cfld)%ifld=iavblfld(iget(985))
4005  fld_info(cfld)%lvl=lvlsxml(1,iget(985))
4006 !$omp parallel do private(i,j,jj)
4007  do j=1,jend-jsta+1
4008  jj = jsta+j-1
4009  do i=1,im
4010  datapd(i,j,cfld) = grid1(i,jj)
4011  enddo
4012  enddo
4013  endif
4014  ENDIF
4015 
4016 ! Effective Helicity
4017 
4018  CALL calhel3(llow,lupp,effust,effvst,esrh)
4019 
4020 !U Bunkers Effective right motion
4021 !
4022 
4023  IF (iget(986)>0) THEN
4024  grid1=spval
4025  DO j=jsta,jend
4026  DO i=1,im
4027  IF(llow(i,j)<spval.and.lupp(i,j)<spval)&
4028  grid1(i,j)=effust(i,j)
4029  ENDDO
4030  ENDDO
4031  if(grib=='grib2') then
4032  cfld=cfld+1
4033  fld_info(cfld)%ifld=iavblfld(iget(986))
4034  fld_info(cfld)%lvl=lvlsxml(1,iget(986))
4035 ! $omp parallel do private(i,j,jj)
4036  do j=1,jend-jsta+1
4037  jj = jsta+j-1
4038  do i=1,im
4039  datapd(i,j,cfld) = grid1(i,jj)
4040  enddo
4041  enddo
4042  endif
4043  ENDIF
4044 
4045 !V Bunkers Effective right motion
4046  IF (iget(987)>0) THEN
4047  grid1=spval
4048  DO j=jsta,jend
4049  DO i=1,im
4050  IF(llow(i,j)<spval.and.lupp(i,j)<spval)&
4051  grid1(i,j)=effvst(i,j)
4052  ENDDO
4053  ENDDO
4054  if(grib=='grib2') then
4055  cfld=cfld+1
4056  fld_info(cfld)%ifld=iavblfld(iget(987))
4057  fld_info(cfld)%lvl=lvlsxml(1,iget(987))
4058 ! $omp parallel do private(i,j,jj)
4059  do j=1,jend-jsta+1
4060  jj = jsta+j-1
4061  do i=1,im
4062  datapd(i,j,cfld) = grid1(i,jj)
4063  enddo
4064  enddo
4065  endif
4066  ENDIF
4067 
4068 !Effective layer helicity
4069  IF (iget(988)>0) THEN
4070  grid1=spval
4071  DO j=jsta,jend
4072  DO i=1,im
4073  IF(llow(i,j)<spval.and.lupp(i,j)<spval)&
4074  grid1(i,j)=esrh(i,j)
4075  ENDDO
4076  ENDDO
4077  if(grib=='grib2') then
4078  cfld=cfld+1
4079  fld_info(cfld)%ifld=iavblfld(iget(988))
4080  fld_info(cfld)%lvl=lvlsxml(1,iget(988))
4081 ! $omp parallel do private(i,j,jj)
4082  do j=1,jend-jsta+1
4083  jj = jsta+j-1
4084  do i=1,im
4085  datapd(i,j,cfld) = grid1(i,jj)
4086  enddo
4087  enddo
4088  endif
4089  ENDIF
4090 
4091 !Effective Layer Tornado Parameter
4092  IF (iget(989)>0) THEN
4093  DO j=jsta,jend
4094  DO i=1,im
4095  IF (mllcl(i,j)>d2000) THEN
4096  mllcltmp=d00
4097  ELSEIF (mllcl(i,j)<d1000) THEN
4098  mllcltmp=1.0
4099  ELSE
4100  mllcltmp=((d2000-mllcl(i,j))/d1000)
4101  ENDIF
4102  IF (eshr(i,j)<12.5) THEN
4103  eshrtmp=d00
4104  ELSEIF (eshr(i,j)>30.0) THEN
4105  eshrtmp=1.5
4106  ELSE
4107  eshrtmp=(eshr(i,j)/20.)
4108  ENDIF
4109  IF (mlcin(i,j)>-50.) THEN
4110  mlcintmp=1.0
4111  ELSEIF (mlcin(i,j)<-200.) THEN
4112  mlcintmp=d00
4113  ELSE
4114  mlcintmp=(200.+mlcin(i,j))/150.
4115  ENDIF
4116  stp=(mlcape(i,j)/d1500)*mllcltmp*(esrh(i,j)/150.)*&
4117  eshrtmp*mlcintmp
4118  grid1(i,j) = spval
4119  IF(llow(i,j)<spval.and.lupp(i,j)<spval) THEN
4120  IF (stp>0) THEN
4121  grid1(i,j)=stp
4122  ELSE
4123  grid1(i,j)=d00
4124  ENDIF
4125  ENDIF
4126  ENDDO
4127  ENDDO
4128  if(grib=='grib2') then
4129  cfld=cfld+1
4130  fld_info(cfld)%ifld=iavblfld(iget(989))
4131  fld_info(cfld)%lvl=lvlsxml(1,iget(989))
4132 ! $omp parallel do private(i,j,jj)
4133  do j=1,jend-jsta+1
4134  jj = jsta+j-1
4135  do i=1,im
4136  datapd(i,j,cfld) = grid1(i,jj)
4137  enddo
4138  enddo
4139  endif
4140  ENDIF
4141 
4142 !Fixed Layer Tornado Parameter
4143  IF (iget(990)>0) THEN
4144  DO j=jsta,jend
4145  DO i=1,im
4146  llmh = nint(lmh(i,j))
4147  p1d(i,j) = pmid(i,j,llmh)
4148  t1d(i,j) = t(i,j,llmh)
4149  q1d(i,j) = q(i,j,llmh)
4150  ENDDO
4151  ENDDO
4152  CALL callcl(p1d,t1d,q1d,egrid1,egrid2)
4153  DO j=jsta,jend
4154  DO i=1,im
4155  slcl(i,j)=egrid2(i,j)
4156  ENDDO
4157  ENDDO
4158  itype = 1
4159  dpbnd = 10.e2
4160  dummy = 0.
4161  idummy = 0
4162  CALL calcape(itype,dpbnd,dummy,dummy,dummy,&
4163  idummy,egrid1,egrid2,&
4164  egrid3,dummy,dummy)
4165 
4166  DO j=jsta,jend
4167  DO i=1,im
4168  IF (slcl(i,j)>d2000) THEN
4169  slcltmp=d00
4170  ELSEIF (slcl(i,j)<=d1000) THEN
4171  slcltmp=1.0
4172  ELSE
4173  slcltmp=((d2000-slcl(i,j))/d1000)
4174  ENDIF
4175  IF (fshr(i,j)<12.5) THEN
4176  fshrtmp=d00
4177  ELSEIF (fshr(i,j)>30.0) THEN
4178  fshrtmp=1.5
4179  ELSE
4180  fshrtmp=(fshr(i,j)/20.)
4181  ENDIF
4182  IF (egrid2(i,j)>-50.) THEN
4183  scintmp=1.0
4184  ELSEIF (egrid2(i,j)<-200.) THEN
4185  scintmp=d00
4186  ELSE
4187  scintmp=((200.+egrid2(i,j)/150.))
4188  ENDIF
4189  stp=(egrid1(i,j)/d1500)*slcltmp*(heli(i,j,2)/150.)*&
4190  fshrtmp*scintmp
4191  grid1(i,j) = spval
4192  IF(t1d(i,j) < spval) THEN
4193  IF (stp>0) THEN
4194  grid1(i,j)=stp
4195  ELSE
4196  grid1(i,j)=d00
4197  ENDIF
4198  ENDIF
4199  ENDDO
4200  ENDDO
4201  if(grib=='grib2') then
4202  cfld=cfld+1
4203  fld_info(cfld)%ifld=iavblfld(iget(990))
4204  fld_info(cfld)%lvl=lvlsxml(1,iget(990))
4205 ! $omp parallel do private(i,j,jj)
4206  do j=1,jend-jsta+1
4207  jj = jsta+j-1
4208  do i=1,im
4209  datapd(i,j,cfld) = grid1(i,jj)
4210  enddo
4211  enddo
4212  endif
4213  ENDIF
4214 
4215 !Effective Layer Supercell Parameter
4216  IF (iget(991)>0) THEN
4217  DO j=jsta,jend
4218  DO i=1,im
4219  IF (eshr(i,j)<10.) THEN
4220  eshrtmp=d00
4221  ELSEIF (eshr(i,j)>20.0) THEN
4222  eshrtmp=1
4223  ELSE
4224  eshrtmp=(eshr(i,j)/20.)
4225  ENDIF
4226  IF (mucin(i,j)>-40.) THEN
4227  mucintmp=1.0
4228  ELSE
4229  mucintmp=(-40./mucin(i,j))
4230  ENDIF
4231  stp=(mucape(i,j)/d1000)*(esrh(i,j)/50.)*&
4232  eshrtmp*mucintmp
4233  grid1(i,j) = spval
4234  IF(t1d(i,j) < spval) THEN
4235  IF (stp>0) THEN
4236  grid1(i,j)=stp
4237  ELSE
4238  grid1(i,j)=d00
4239  ENDIF
4240  ENDIF
4241  ENDDO
4242  ENDDO
4243  if(grib=='grib2') then
4244  cfld=cfld+1
4245  fld_info(cfld)%ifld=iavblfld(iget(991))
4246  fld_info(cfld)%lvl=lvlsxml(1,iget(991))
4247 ! $omp parallel do private(i,j,jj)
4248  do j=1,jend-jsta+1
4249  jj = jsta+j-1
4250  do i=1,im
4251  datapd(i,j,cfld) = grid1(i,jj)
4252  enddo
4253  enddo
4254  endif
4255  ENDIF
4256 
4257 !Mixed Layer (100 mb) Virtual LFC
4258 
4259  IF (iget(992)>0) THEN
4260 !$omp parallel do private(i,j)
4261  DO j=jsta,jend
4262  DO i=1,im
4263  egrid1(i,j) = -h99999
4264  egrid2(i,j) = -h99999
4265  egrid3(i,j) = -h99999
4266  egrid4(i,j) = -h99999
4267  egrid5(i,j) = -h99999
4268  egrid6(i,j) = -h99999
4269  egrid7(i,j) = -h99999
4270  egrid8(i,j) = -h99999
4271  lb2(i,j) = (lvlbnd(i,j,1) + lvlbnd(i,j,2) + &
4272  lvlbnd(i,j,3))/3
4273  p1d(i,j) = (pbnd(i,j,1) + pbnd(i,j,2) + pbnd(i,j,3))/3
4274  t1d(i,j) = (tvirtual(tbnd(i,j,1),qbnd(i,j,1)) + &
4275  tvirtual(tbnd(i,j,2),qbnd(i,j,2)) + &
4276  tvirtual(tbnd(i,j,3),qbnd(i,j,3)))/3
4277  q1d(i,j) = (qbnd(i,j,1) + qbnd(i,j,2) + qbnd(i,j,3))/3
4278  ENDDO
4279  ENDDO
4280 
4281  dpbnd = 0.
4282  itype = 2
4283 ! EGRID3 is Virtual LFC
4284  CALL calcape2(itype,dpbnd,p1d,t1d,q1d,lb2, &
4285  egrid1,egrid2,egrid3,egrid4,egrid5, &
4286  egrid6,egrid7,egrid8)
4287 
4288  grid1=spval
4289  DO j=jsta,jend
4290  DO i=1,im
4291  IF(t1d(i,j) < spval) grid1(i,j) = egrid3(i,j)
4292  ENDDO
4293  ENDDO
4294  CALL bound(grid1,d00,h99999)
4295  if(grib=='grib2') then
4296  cfld=cfld+1
4297  fld_info(cfld)%ifld=iavblfld(iget(992))
4298  fld_info(cfld)%lvl=lvlsxml(1,iget(992))
4299 !$omp parallel do private(i,j,jj)
4300  do j=1,jend-jsta+1
4301  jj = jsta+j-1
4302  do i=1,im
4303  datapd(i,j,cfld) = grid1(i,jj)
4304  enddo
4305  enddo
4306  endif
4307  ENDIF !992
4308 
4309 
4310  IF (iget(763)>0) THEN
4311 !$omp parallel do private(i,j)
4312 ! EGRID3 is Virtual LFC
4313  DO j=jsta,jend
4314  DO i=1,im
4315  grid1(i,j) = q1d(i,j)
4316  ENDDO
4317  ENDDO
4318  if(grib=='grib2') then
4319  cfld=cfld+1
4320  fld_info(cfld)%ifld=iavblfld(iget(763))
4321  fld_info(cfld)%lvl=lvlsxml(1,iget(763))
4322 !$omp parallel do private(i,j,jj)
4323  do j=1,jend-jsta+1
4324  jj = jsta+j-1
4325  do i=1,im
4326  datapd(i,j,cfld) = grid1(i,jj)
4327  enddo
4328  enddo
4329  endif
4330  ENDIF
4331 
4332 !Hail parameter
4333  IF (iget(993)>0) THEN
4334  grid1=spval
4335  DO j=jsta,jend
4336  DO i=1,im
4337  lapse=-((t700(i,j)-t500(i,j))/((z700(i,j)-z500(i,j))))
4338  ship=(mucape(i,j)*d1000*muq1d(i,j)*lapse*(t500(i,j)-k2c)*fshr(i,j))/hconst
4339  IF (mucape(i,j)<1300.)THEN
4340  ship=ship*(mucape(i,j)/1300.)
4341  ENDIF
4342  IF (lapse < 5.8)THEN
4343  ship=ship*(lapse/5.8)
4344  ENDIF
4345  IF (freezelvl(i,j) < 2400.)THEN
4346  ship=ship*(freezelvl(i,j)/2400.)
4347  ENDIF
4348  grid1(i,j)=ship
4349  ENDDO
4350  ENDDO
4351  if(grib=='grib2') then
4352  cfld=cfld+1
4353  fld_info(cfld)%ifld=iavblfld(iget(993))
4354  fld_info(cfld)%lvl=lvlsxml(1,iget(993))
4355 ! $omp parallel do private(i,j,jj)
4356  do j=1,jend-jsta+1
4357  jj = jsta+j-1
4358  do i=1,im
4359  datapd(i,j,cfld) = grid1(i,jj)
4360  enddo
4361  enddo
4362  endif
4363  ENDIF
4364 
4365  ENDIF !END RTMA BLOCK
4366 
4367 
4368 ! Critical Angle
4369 
4370  IF (iget(957)>0) THEN
4371  grid1=spval
4372 !$omp parallel do private(i,j)
4373  DO j=jsta,jend
4374  DO i=1,im
4375  IF(t1d(i,j) < spval ) grid1(i,j) = cangle(i,j)
4376  ! IF(EGRID1(I,J)<100. .OR. EGRID2(I,J)>-250.) THEN
4377  ! GRID1(I,J) = 0.
4378  ! ENDIF
4379  ENDDO
4380  ENDDO
4381  if(grib=='grib2') then
4382  cfld=cfld+1
4383  fld_info(cfld)%ifld=iavblfld(iget(957))
4384  fld_info(cfld)%lvl=lvlsxml(1,iget(957))
4385 !$omp parallel do private(i,j,jj)
4386  do j=1,jend-jsta+1
4387  jj = jsta+j-1
4388  do i=1,im
4389  datapd(i,j,cfld) = grid1(i,jj)
4390  enddo
4391  enddo
4392  endif
4393  ENDIF !957
4394 
4395 ! Dendritic Layer Depth, -17C < T < -12C
4396 
4397  IF (iget(955)>0) THEN
4398  grid1=spval
4399 !$omp parallel do private(i,j)
4400  DO j=jsta,jend
4401  DO i=1,im
4402  IF(t1d(i,j) < spval ) grid1(i,j) = egrid7(i,j)
4403  ENDDO
4404  ENDDO
4405  CALL bound(grid1,d00,h99999)
4406  if(grib=='grib2') then
4407  cfld=cfld+1
4408  fld_info(cfld)%ifld=iavblfld(iget(955))
4409  fld_info(cfld)%lvl=lvlsxml(1,iget(955))
4410 !$omp parallel do private(i,j,jj)
4411  do j=1,jend-jsta+1
4412  jj = jsta+j-1
4413  do i=1,im
4414  datapd(i,j,cfld) = grid1(i,jj)
4415  enddo
4416  enddo
4417  endif
4418  ENDIF !955
4419 
4420 ! Enhanced Stretching Potential
4421 
4422  IF (iget(956)>0) THEN
4423  grid1=spval
4424 !$omp parallel do private(i,j)
4425  DO j=jsta,jend
4426  DO i=1,im
4427  IF(t1d(i,j) < spval ) grid1(i,j) = egrid8(i,j)
4428  ENDDO
4429  ENDDO
4430  CALL bound(grid1,d00,h99999)
4431  if(grib=='grib2') then
4432  cfld=cfld+1
4433  fld_info(cfld)%ifld=iavblfld(iget(956))
4434  fld_info(cfld)%lvl=lvlsxml(1,iget(956))
4435 !$omp parallel do private(i,j,jj)
4436  do j=1,jend-jsta+1
4437  jj = jsta+j-1
4438  do i=1,im
4439  datapd(i,j,cfld) = grid1(i,jj)
4440  enddo
4441  enddo
4442  endif
4443  ENDIF !956
4444 
4445 ! Downdraft CAPE
4446 
4447 ! ITYPE = 1
4448 ! DO J=JSTA,JEND
4449 ! DO I=1,IM
4450 ! LB2(I,J) = (LVLBND(I,J,1) + LVLBND(I,J,2) + &
4451 ! LVLBND(I,J,3))/3
4452 ! P1D(I,J) = (PBND(I,J,1) + PBND(I,J,2) + PBND(I,J,3))/3
4453 ! T1D(I,J) = (TBND(I,J,1) + TBND(I,J,2) + TBND(I,J,3))/3
4454 ! Q1D(I,J) = (QBND(I,J,1) + QBND(I,J,2) + QBND(I,J,3))/3
4455 ! ENDDO
4456 ! ENDDO
4457 
4458 ! DPBND = 400.E2
4459 ! CALL CALCAPE2(ITYPE,DPBND,P1D,T1D,Q1D,LB2, &
4460 ! EGRID1,EGRID2,EGRID3,EGRID4,EGRID5, &
4461 ! EGRID6,EGRID7,EGRID8)
4462 
4463  IF (iget(954)>0) THEN
4464  grid1 = spval
4465 !$omp parallel do private(i,j)
4466  DO j=jsta,jend
4467  DO i=1,im
4468  IF(t1d(i,j) < spval) grid1(i,j) = -egrid6(i,j)
4469  ENDDO
4470  ENDDO
4471  CALL bound(grid1,d00,h99999)
4472  if(grib=='grib2') then
4473  cfld=cfld+1
4474  fld_info(cfld)%ifld=iavblfld(iget(954))
4475  fld_info(cfld)%lvl=lvlsxml(1,iget(954))
4476 !$omp parallel do private(i,j,jj)
4477  do j=1,jend-jsta+1
4478  jj = jsta+j-1
4479  do i=1,im
4480  datapd(i,j,cfld) = grid1(i,jj)
4481  enddo
4482  enddo
4483  endif
4484 
4485  ENDIF !954
4486 
4487  if (allocated(ushr1)) deallocate(ushr1)
4488  if (allocated(vshr1)) deallocate(vshr1)
4489  if (allocated(ushr6)) deallocate(ushr6)
4490  if (allocated(vshr6)) deallocate(vshr6)
4491  if (allocated(ust)) deallocate(ust)
4492  if (allocated(vst)) deallocate(vst)
4493  if (allocated(heli)) deallocate(heli)
4494  if (allocated(llow)) deallocate(llow)
4495  if (allocated(lupp)) deallocate(lupp)
4496  if (allocated(cangle))deallocate(cangle)
4497  if (allocated(effust))deallocate(effust)
4498  if (allocated(effvst))deallocate(effvst)
4499  if (allocated(eshr)) deallocate(eshr)
4500  if (allocated(uvect)) deallocate(uvect)
4501  if (allocated(vvect)) deallocate(vvect)
4502  if (allocated(esrh)) deallocate(esrh)
4503  if (allocated(htsfc)) deallocate(htsfc)
4504  if (allocated(fshr)) deallocate(fshr)
4505  ENDIF
4506 
4507  if (allocated(pbnd)) deallocate(pbnd)
4508  if (allocated(tbnd)) deallocate(tbnd)
4509  if (allocated(qbnd)) deallocate(qbnd)
4510  if (allocated(ubnd)) deallocate(ubnd)
4511  if (allocated(vbnd)) deallocate(vbnd)
4512  if (allocated(rhbnd)) deallocate(rhbnd)
4513  if (allocated(wbnd)) deallocate(wbnd)
4514  if (allocated(lvlbnd)) deallocate(lvlbnd)
4515  if (allocated(lb2)) deallocate(lb2)
4516 !
4517 !
4518 ! RELATIVE HUMIDITY WITH RESPECT TO PRECIPITABLE WATER
4519  IF (iget(749)>0) THEN
4520  CALL calrh_pw(grid1(1,jsta))
4521  if(grib=='grib2') then
4522  cfld=cfld+1
4523  fld_info(cfld)%ifld=iavblfld(iget(749))
4524 !$omp parallel do private(i,j,jj)
4525  do j=1,jend-jsta+1
4526  jj = jsta+j-1
4527  do i=1,im
4528  datapd(i,j,cfld) = grid1(i,jj)
4529  enddo
4530  enddo
4531  endif
4532  ENDIF
4533 
4534 
4535 ! END OF ROUTINE.
4536 !
4537  RETURN
4538  END
subroutine mxwind(km, p, u, v, t, h, pmw, umw, vmw, tmw, hmw)
mxwind() computes maximum wind level fields.
Definition: GFSPOST.F:431
Definition: MASKS_mod.f:1
elemental real function, public fpvsnew(t)
Definition: UPP_PHYSICS.f:341
subroutine tpause(km, p, u, v, t, h, ptp, utp, vtp, ttp, htp, shrtp)
tpause() computes tropopause level fields.
Definition: GFSPOST.F:344
calcape() computes CAPE/CINS and other storm related variables.
Definition: UPP_PHYSICS.f:27