source: flexpart.git/src/verttransform_gfs.f90 @ 3481cc1

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file since 3481cc1 was 3481cc1, checked in by Ignacio Pisso <ip@…>, 4 years ago

move license from headers to a different file

  • Property mode set to 100644
File size: 20.0 KB
Line 
1subroutine verttransform_gfs(n,uuh,vvh,wwh,pvh)
2  !                      i  i   i   i   i
3  !*****************************************************************************
4  !                                                                            *
5  !     This subroutine transforms temperature, dew point temperature and      *
6  !     wind components from eta to meter coordinates.                         *
7  !     The vertical wind component is transformed from Pa/s to m/s using      *
8  !     the conversion factor pinmconv.                                        *
9  !     In addition, this routine calculates vertical density gradients        *
10  !     needed for the parameterization of the turbulent velocities.           *
11  !                                                                            *
12  !     Author: A. Stohl, G. Wotawa                                            *
13  !                                                                            *
14  !     12 August 1996                                                         *
15  !     Update: 16 January 1998                                                *
16  !                                                                            *
17  !     Major update: 17 February 1999                                         *
18  !     by G. Wotawa                                                           *
19  !     CHANGE 17/11/2005 Caroline Forster, NCEP GFS version                   *
20  !                                                                            *
21  !   - Vertical levels for u, v and w are put together                        *
22  !   - Slope correction for vertical velocity: Modification of calculation    *
23  !     procedure                                                              *
24  !                                                                            *
25  !*****************************************************************************
26  !  Changes, Bernd C. Krueger, Feb. 2001:
27  !   Variables tth and qvh (on eta coordinates) from common block
28  !
29  !   Unified ECMWF and GFS builds                                     
30  !   Marian Harustak, 12.5.2017                                       
31  !     - Renamed routine from verttransform to verttransform_gfs
32  !
33  !*****************************************************************************
34  !                                                                            *
35  ! Variables:                                                                 *
36  ! nx,ny,nz                        field dimensions in x,y and z direction    *
37  ! uu(0:nxmax,0:nymax,nzmax,2)     wind components in x-direction [m/s]       *
38  ! vv(0:nxmax,0:nymax,nzmax,2)     wind components in y-direction [m/s]       *
39  ! ww(0:nxmax,0:nymax,nzmax,2)     wind components in z-direction [deltaeta/s]*
40  ! tt(0:nxmax,0:nymax,nzmax,2)     temperature [K]                            *
41  ! pv(0:nxmax,0:nymax,nzmax,2)     potential voriticity (pvu)                 *
42  ! ps(0:nxmax,0:nymax,2)           surface pressure [Pa]                      *
43  ! clouds(0:nxmax,0:nymax,0:nzmax,2) cloud field for wet deposition           *
44  !                                                                            *
45  !*****************************************************************************
46
47  use par_mod
48  use com_mod
49  use cmapf_mod
50
51  implicit none
52
53  integer :: ix,jy,kz,iz,n,kmin,kl,klp,ix1,jy1,ixp,jyp,ixm,jym
54  integer :: rain_cloud_above,kz_inv
55  real :: f_qvsat,pressure
56  real :: rh,lsp,cloudh_min,convp,prec
57  real :: rhoh(nuvzmax),pinmconv(nzmax)
58  real :: ew,pint,tv,tvold,pold,dz1,dz2,dz,ui,vi
59  real :: xlon,ylat,xlonr,dzdx,dzdy
60  real :: dzdx1,dzdx2,dzdy1,dzdy2,cosf
61  real :: uuaux,vvaux,uupolaux,vvpolaux,ddpol,ffpol,wdummy
62  real :: uuh(0:nxmax-1,0:nymax-1,nuvzmax)
63  real :: vvh(0:nxmax-1,0:nymax-1,nuvzmax)
64  real :: pvh(0:nxmax-1,0:nymax-1,nuvzmax)
65  real :: wwh(0:nxmax-1,0:nymax-1,nwzmax)
66  real :: wzlev(nwzmax),uvwzlev(0:nxmax-1,0:nymax-1,nzmax)
67  real,parameter :: const=r_air/ga
68
69  ! NCEP version
70  integer :: llev, i
71
72  logical :: init = .true.
73
74
75  !*************************************************************************
76  ! If verttransform is called the first time, initialize heights of the   *
77  ! z levels in meter. The heights are the heights of model levels, where  *
78  ! u,v,T and qv are given.                                                *
79  !*************************************************************************
80
81  if (init) then
82
83  ! Search for a point with high surface pressure (i.e. not above significant topography)
84  ! Then, use this point to construct a reference z profile, to be used at all times
85  !*****************************************************************************
86
87    do jy=0,nymin1
88      do ix=0,nxmin1
89        if (ps(ix,jy,1,n).gt.100000.) then
90          ixm=ix
91          jym=jy
92          goto 3
93        endif
94      end do
95    end do
963   continue
97
98
99    tvold=tt2(ixm,jym,1,n)*(1.+0.378*ew(td2(ixm,jym,1,n))/ &
100    ps(ixm,jym,1,n))
101    pold=ps(ixm,jym,1,n)
102    height(1)=0.
103
104    do kz=2,nuvz
105      pint=akz(kz)+bkz(kz)*ps(ixm,jym,1,n)
106      tv=tth(ixm,jym,kz,n)*(1.+0.608*qvh(ixm,jym,kz,n))
107
108      if (abs(tv-tvold).gt.0.2) then
109        height(kz)=height(kz-1)+const*log(pold/pint)* &
110        (tv-tvold)/log(tv/tvold)
111      else
112        height(kz)=height(kz-1)+const*log(pold/pint)*tv
113      endif
114
115      tvold=tv
116      pold=pint
117    end do
118
119
120  ! Determine highest levels that can be within PBL
121  !************************************************
122
123    do kz=1,nz
124      if (height(kz).gt.hmixmax) then
125        nmixz=kz
126        goto 9
127      endif
128    end do
1299   continue
130
131  ! Do not repeat initialization of the Cartesian z grid
132  !*****************************************************
133
134    init=.false.
135
136  endif
137
138
139  ! Loop over the whole grid
140  !*************************
141
142  do jy=0,nymin1
143    do ix=0,nxmin1
144
145  ! NCEP version: find first level above ground
146      llev = 0
147      do i=1,nuvz
148        if (ps(ix,jy,1,n).lt.akz(i)) llev=i
149      end do
150       llev = llev+1
151       if (llev.gt.nuvz-2) llev = nuvz-2
152  !     if (llev.eq.nuvz-2) write(*,*) 'verttransform
153  !    +WARNING: LLEV eq NUZV-2'
154  ! NCEP version
155
156
157  ! compute height of pressure levels above ground
158  !***********************************************
159
160      tvold=tth(ix,jy,llev,n)*(1.+0.608*qvh(ix,jy,llev,n))
161      pold=akz(llev)
162      wzlev(llev)=0.
163      uvwzlev(ix,jy,llev)=0.
164      rhoh(llev)=pold/(r_air*tvold)
165
166      do kz=llev+1,nuvz
167        pint=akz(kz)+bkz(kz)*ps(ix,jy,1,n)
168        tv=tth(ix,jy,kz,n)*(1.+0.608*qvh(ix,jy,kz,n))
169        rhoh(kz)=pint/(r_air*tv)
170
171        if (abs(tv-tvold).gt.0.2) then
172          uvwzlev(ix,jy,kz)=uvwzlev(ix,jy,kz-1)+const*log(pold/pint)* &
173          (tv-tvold)/log(tv/tvold)
174        else
175          uvwzlev(ix,jy,kz)=uvwzlev(ix,jy,kz-1)+const*log(pold/pint)*tv
176        endif
177        wzlev(kz)=uvwzlev(ix,jy,kz)
178
179        tvold=tv
180        pold=pint
181      end do
182
183  ! pinmconv=(h2-h1)/(p2-p1)
184
185      pinmconv(llev)=(uvwzlev(ix,jy,llev+1)-uvwzlev(ix,jy,llev))/ &
186           ((aknew(llev+1)+bknew(llev+1)*ps(ix,jy,1,n))- &
187           (aknew(llev)+bknew(llev)*ps(ix,jy,1,n)))
188      do kz=llev+1,nz-1
189        pinmconv(kz)=(uvwzlev(ix,jy,kz+1)-uvwzlev(ix,jy,kz-1))/ &
190             ((aknew(kz+1)+bknew(kz+1)*ps(ix,jy,1,n))- &
191             (aknew(kz-1)+bknew(kz-1)*ps(ix,jy,1,n)))
192      end do
193      pinmconv(nz)=(uvwzlev(ix,jy,nz)-uvwzlev(ix,jy,nz-1))/ &
194           ((aknew(nz)+bknew(nz)*ps(ix,jy,1,n))- &
195           (aknew(nz-1)+bknew(nz-1)*ps(ix,jy,1,n)))
196
197
198  ! Levels, where u,v,t and q are given
199  !************************************
200
201      uu(ix,jy,1,n)=uuh(ix,jy,llev)
202      vv(ix,jy,1,n)=vvh(ix,jy,llev)
203      tt(ix,jy,1,n)=tth(ix,jy,llev,n)
204      qv(ix,jy,1,n)=qvh(ix,jy,llev,n)
205! IP & SEC, 201812 add clouds
206      if (readclouds) then
207         clwc(ix,jy,1,n)=clwch(ix,jy,llev,n)
208      endif
209      pv(ix,jy,1,n)=pvh(ix,jy,llev)
210      rho(ix,jy,1,n)=rhoh(llev)
211      pplev(ix,jy,1,n)=akz(llev)
212      uu(ix,jy,nz,n)=uuh(ix,jy,nuvz)
213      vv(ix,jy,nz,n)=vvh(ix,jy,nuvz)
214      tt(ix,jy,nz,n)=tth(ix,jy,nuvz,n)
215      qv(ix,jy,nz,n)=qvh(ix,jy,nuvz,n)
216! IP & SEC, 201812 add clouds
217      if (readclouds) then
218         clwc(ix,jy,nz,n)=clwch(ix,jy,nuvz,n)
219      endif
220      pv(ix,jy,nz,n)=pvh(ix,jy,nuvz)
221      rho(ix,jy,nz,n)=rhoh(nuvz)
222      pplev(ix,jy,nz,n)=akz(nuvz)
223      kmin=llev+1
224      do iz=2,nz-1
225        do kz=kmin,nuvz
226          if(height(iz).gt.uvwzlev(ix,jy,nuvz)) then
227            uu(ix,jy,iz,n)=uu(ix,jy,nz,n)
228            vv(ix,jy,iz,n)=vv(ix,jy,nz,n)
229            tt(ix,jy,iz,n)=tt(ix,jy,nz,n)
230            qv(ix,jy,iz,n)=qv(ix,jy,nz,n)
231! IP & SEC, 201812 add clouds
232            if (readclouds) then
233               clwc(ix,jy,iz,n)=clwc(ix,jy,nz,n)
234            endif
235            pv(ix,jy,iz,n)=pv(ix,jy,nz,n)
236            rho(ix,jy,iz,n)=rho(ix,jy,nz,n)
237            pplev(ix,jy,iz,n)=pplev(ix,jy,nz,n)
238            goto 30
239          endif
240          if ((height(iz).gt.uvwzlev(ix,jy,kz-1)).and. &
241          (height(iz).le.uvwzlev(ix,jy,kz))) then
242            dz1=height(iz)-uvwzlev(ix,jy,kz-1)
243            dz2=uvwzlev(ix,jy,kz)-height(iz)
244            dz=dz1+dz2
245            uu(ix,jy,iz,n)=(uuh(ix,jy,kz-1)*dz2+uuh(ix,jy,kz)*dz1)/dz
246            vv(ix,jy,iz,n)=(vvh(ix,jy,kz-1)*dz2+vvh(ix,jy,kz)*dz1)/dz
247            tt(ix,jy,iz,n)=(tth(ix,jy,kz-1,n)*dz2 &
248            +tth(ix,jy,kz,n)*dz1)/dz
249            qv(ix,jy,iz,n)=(qvh(ix,jy,kz-1,n)*dz2 &
250            +qvh(ix,jy,kz,n)*dz1)/dz
251! IP & SEC, 201812 add clouds
252            if (readclouds) then
253               clwc(ix,jy,iz,n)=(clwch(ix,jy,kz-1,n)*dz2 &
254               +clwch(ix,jy,kz,n)*dz1)/dz
255            endif
256            pv(ix,jy,iz,n)=(pvh(ix,jy,kz-1)*dz2+pvh(ix,jy,kz)*dz1)/dz
257            rho(ix,jy,iz,n)=(rhoh(kz-1)*dz2+rhoh(kz)*dz1)/dz
258            pplev(ix,jy,iz,n)=(akz(kz-1)*dz2+akz(kz)*dz1)/dz
259          endif
260        end do
26130      continue
262      end do
263
264
265  ! Levels, where w is given
266  !*************************
267
268      ww(ix,jy,1,n)=wwh(ix,jy,llev)*pinmconv(llev)
269      ww(ix,jy,nz,n)=wwh(ix,jy,nwz)*pinmconv(nz)
270      kmin=llev+1
271      do iz=2,nz
272        do kz=kmin,nwz
273          if ((height(iz).gt.wzlev(kz-1)).and. &
274          (height(iz).le.wzlev(kz))) then
275            dz1=height(iz)-wzlev(kz-1)
276            dz2=wzlev(kz)-height(iz)
277            dz=dz1+dz2
278            ww(ix,jy,iz,n)=(wwh(ix,jy,kz-1)*pinmconv(kz-1)*dz2 &
279            +wwh(ix,jy,kz)*pinmconv(kz)*dz1)/dz
280          endif
281        end do
282      end do
283
284
285  ! Compute density gradients at intermediate levels
286  !*************************************************
287
288      drhodz(ix,jy,1,n)=(rho(ix,jy,2,n)-rho(ix,jy,1,n))/ &
289           (height(2)-height(1))
290      do kz=2,nz-1
291        drhodz(ix,jy,kz,n)=(rho(ix,jy,kz+1,n)-rho(ix,jy,kz-1,n))/ &
292        (height(kz+1)-height(kz-1))
293      end do
294      drhodz(ix,jy,nz,n)=drhodz(ix,jy,nz-1,n)
295
296    end do
297  end do
298
299
300  !****************************************************************
301  ! Compute slope of eta levels in windward direction and resulting
302  ! vertical wind correction
303  !****************************************************************
304
305  do jy=1,ny-2
306    cosf=cos((real(jy)*dy+ylat0)*pi180)
307    do ix=1,nx-2
308
309  ! NCEP version: find first level above ground
310      llev = 0
311      do i=1,nuvz
312       if (ps(ix,jy,1,n).lt.akz(i)) llev=i
313      end do
314       llev = llev+1
315       if (llev.gt.nuvz-2) llev = nuvz-2
316  !     if (llev.eq.nuvz-2) write(*,*) 'verttransform
317  !    +WARNING: LLEV eq NUZV-2'
318  ! NCEP version
319
320      kmin=llev+1
321      do iz=2,nz-1
322
323        ui=uu(ix,jy,iz,n)*dxconst/cosf
324        vi=vv(ix,jy,iz,n)*dyconst
325
326        do kz=kmin,nz
327          if ((height(iz).gt.uvwzlev(ix,jy,kz-1)).and. &
328          (height(iz).le.uvwzlev(ix,jy,kz))) then
329            dz1=height(iz)-uvwzlev(ix,jy,kz-1)
330            dz2=uvwzlev(ix,jy,kz)-height(iz)
331            dz=dz1+dz2
332            kl=kz-1
333            klp=kz
334            goto 47
335          endif
336        end do
337
33847      ix1=ix-1
339        jy1=jy-1
340        ixp=ix+1
341        jyp=jy+1
342
343        dzdx1=(uvwzlev(ixp,jy,kl)-uvwzlev(ix1,jy,kl))/2.
344        dzdx2=(uvwzlev(ixp,jy,klp)-uvwzlev(ix1,jy,klp))/2.
345        dzdx=(dzdx1*dz2+dzdx2*dz1)/dz
346
347        dzdy1=(uvwzlev(ix,jyp,kl)-uvwzlev(ix,jy1,kl))/2.
348        dzdy2=(uvwzlev(ix,jyp,klp)-uvwzlev(ix,jy1,klp))/2.
349        dzdy=(dzdy1*dz2+dzdy2*dz1)/dz
350
351        ww(ix,jy,iz,n)=ww(ix,jy,iz,n)+(dzdx*ui+dzdy*vi)
352
353      end do
354
355    end do
356  end do
357
358
359  ! If north pole is in the domain, calculate wind velocities in polar
360  ! stereographic coordinates
361  !*******************************************************************
362
363  if (nglobal) then
364    do jy=int(switchnorthg)-2,nymin1
365      ylat=ylat0+real(jy)*dy
366      do ix=0,nxmin1
367        xlon=xlon0+real(ix)*dx
368        do iz=1,nz
369          call cc2gll(northpolemap,ylat,xlon,uu(ix,jy,iz,n), &
370               vv(ix,jy,iz,n),uupol(ix,jy,iz,n), &
371               vvpol(ix,jy,iz,n))
372        end do
373      end do
374    end do
375
376
377    do iz=1,nz
378
379  ! CALCULATE FFPOL, DDPOL FOR CENTRAL GRID POINT
380      xlon=xlon0+real(nx/2-1)*dx
381      xlonr=xlon*pi/180.
382      ffpol=sqrt(uu(nx/2-1,nymin1,iz,n)**2+vv(nx/2-1,nymin1,iz,n)**2)
383      if (vv(nx/2-1,nymin1,iz,n).lt.0.) then
384        ddpol=atan(uu(nx/2-1,nymin1,iz,n)/vv(nx/2-1,nymin1,iz,n))-xlonr
385      elseif (vv(nx/2-1,nymin1,iz,n).gt.0.) then
386        ddpol=pi+atan(uu(nx/2-1,nymin1,iz,n)/ &
387        vv(nx/2-1,nymin1,iz,n))-xlonr
388      else
389        ddpol=pi/2-xlonr
390      endif
391      if(ddpol.lt.0.) ddpol=2.0*pi+ddpol
392      if(ddpol.gt.2.0*pi) ddpol=ddpol-2.0*pi
393
394  ! CALCULATE U,V FOR 180 DEG, TRANSFORM TO POLAR STEREOGRAPHIC GRID
395      xlon=180.0
396      xlonr=xlon*pi/180.
397      ylat=90.0
398      uuaux=-ffpol*sin(xlonr+ddpol)
399      vvaux=-ffpol*cos(xlonr+ddpol)
400      call cc2gll(northpolemap,ylat,xlon,uuaux,vvaux,uupolaux,vvpolaux)
401      jy=nymin1
402      do ix=0,nxmin1
403        uupol(ix,jy,iz,n)=uupolaux
404        vvpol(ix,jy,iz,n)=vvpolaux
405      end do
406    end do
407
408
409  ! Fix: Set W at pole to the zonally averaged W of the next equator-
410  ! ward parallel of latitude
411
412    do iz=1,nz
413      wdummy=0.
414      jy=ny-2
415      do ix=0,nxmin1
416        wdummy=wdummy+ww(ix,jy,iz,n)
417      end do
418      wdummy=wdummy/real(nx)
419      jy=nymin1
420      do ix=0,nxmin1
421        ww(ix,jy,iz,n)=wdummy
422      end do
423    end do
424
425  endif
426
427
428  ! If south pole is in the domain, calculate wind velocities in polar
429  ! stereographic coordinates
430  !*******************************************************************
431
432  if (sglobal) then
433    do jy=0,int(switchsouthg)+3
434      ylat=ylat0+real(jy)*dy
435      do ix=0,nxmin1
436        xlon=xlon0+real(ix)*dx
437        do iz=1,nz
438          call cc2gll(southpolemap,ylat,xlon,uu(ix,jy,iz,n), &
439          vv(ix,jy,iz,n),uupol(ix,jy,iz,n),vvpol(ix,jy,iz,n))
440        end do
441      end do
442    end do
443
444    do iz=1,nz
445
446  ! CALCULATE FFPOL, DDPOL FOR CENTRAL GRID POINT
447      xlon=xlon0+real(nx/2-1)*dx
448      xlonr=xlon*pi/180.
449      ffpol=sqrt(uu(nx/2-1,0,iz,n)**2+vv(nx/2-1,0,iz,n)**2)
450      if(vv(nx/2-1,0,iz,n).lt.0.) then
451        ddpol=atan(uu(nx/2-1,0,iz,n)/vv(nx/2-1,0,iz,n))+xlonr
452      elseif (vv(nx/2-1,0,iz,n).gt.0.) then
453        ddpol=pi+atan(uu(nx/2-1,0,iz,n)/vv(nx/2-1,0,iz,n))-xlonr
454      else
455        ddpol=pi/2-xlonr
456      endif
457      if(ddpol.lt.0.) ddpol=2.0*pi+ddpol
458      if(ddpol.gt.2.0*pi) ddpol=ddpol-2.0*pi
459
460  ! CALCULATE U,V FOR 180 DEG, TRANSFORM TO POLAR STEREOGRAPHIC GRID
461      xlon=180.0
462      xlonr=xlon*pi/180.
463      ylat=-90.0
464      uuaux=+ffpol*sin(xlonr-ddpol)
465      vvaux=-ffpol*cos(xlonr-ddpol)
466      call cc2gll(northpolemap,ylat,xlon,uuaux,vvaux,uupolaux,vvpolaux)
467
468      jy=0
469      do ix=0,nxmin1
470        uupol(ix,jy,iz,n)=uupolaux
471        vvpol(ix,jy,iz,n)=vvpolaux
472      end do
473    end do
474
475
476  ! Fix: Set W at pole to the zonally averaged W of the next equator-
477  ! ward parallel of latitude
478
479    do iz=1,nz
480      wdummy=0.
481      jy=1
482      do ix=0,nxmin1
483        wdummy=wdummy+ww(ix,jy,iz,n)
484      end do
485      wdummy=wdummy/real(nx)
486      jy=0
487      do ix=0,nxmin1
488        ww(ix,jy,iz,n)=wdummy
489      end do
490    end do
491  endif
492
493
494
495!***********************************************************************************
496! IP & SEC, 201812 GFS clouds read
497  if (readclouds) then
498! The method is loops all grids vertically and constructs the 3D matrix for clouds
499! Cloud top and cloud bottom gid cells are assigned as well as the total column
500! cloud water. For precipitating grids, the type and whether it is in or below
501! cloud scavenging are assigned with numbers 2-5 (following the old metod).
502! Distinction is done for lsp and convp though they are treated the same in regards
503! to scavenging. Also clouds that are not precipitating are defined which may be
504! to include future cloud processing by non-precipitating-clouds.
505!***********************************************************************************
506    write(*,*) 'Global NCEP fields: using cloud water'
507    clw(:,:,:,n)=0.0
508    ctwc(:,:,n)=0.0
509    clouds(:,:,:,n)=0
510! If water/ice are read separately into clwc and ciwc, store sum in clwc
511    do jy=0,nymin1
512      do ix=0,nxmin1
513        lsp=lsprec(ix,jy,1,n)
514        convp=convprec(ix,jy,1,n)
515        prec=lsp+convp
516! Find clouds in the vertical
517        do kz=1, nz-1 !go from top to bottom
518          if (clwc(ix,jy,kz,n).gt.0) then
519! assuming rho is in kg/m3 and hz in m gives: kg/kg * kg/m3 *m3/kg /m = m2/m3
520            clw(ix,jy,kz,n)=(clwc(ix,jy,kz,n)*rho(ix,jy,kz,n))*(height(kz+1)-height(kz))
521            ctwc(ix,jy,n) = ctwc(ix,jy,n)+clw(ix,jy,kz,n)
522            cloudh_min=min(height(kz+1),height(kz))
523          endif
524        end do
525
526! If Precipitation. Define removal type in the vertical
527        if ((lsp.gt.0.01).or.(convp.gt.0.01)) then ! cloud and precipitation
528
529          do kz=nz,2,-1 !go Bottom up!
530            if (clw(ix,jy,kz,n).gt. 0) then ! is in cloud
531              cloudsh(ix,jy,n)=cloudsh(ix,jy,n)+height(kz)-height(kz-1)
532              clouds(ix,jy,kz,n)=1                               ! is a cloud
533              if (lsp.ge.convp) then
534                clouds(ix,jy,kz,n)=3                             ! lsp in-cloud
535              else
536                clouds(ix,jy,kz,n)=2                             ! convp in-cloud
537              endif                                              ! convective or large scale
538            elseif((clw(ix,jy,kz,n).le.0) .and. (cloudh_min.ge.height(kz))) then ! is below cloud
539              if (lsp.ge.convp) then
540                clouds(ix,jy,kz,n)=5                             ! lsp dominated washout
541              else
542                clouds(ix,jy,kz,n)=4                             ! convp dominated washout
543              endif                                              ! convective or large scale
544            endif
545
546            if (height(kz).ge. 19000) then                        ! set a max height for removal
547              clouds(ix,jy,kz,n)=0
548            endif !clw>0
549          end do !nz
550        endif ! precipitation
551      end do
552    end do
553  else
554  write(*,*) 'Global NCEP fields: using cloud water from Parameterization'
555  !   write (*,*) 'initializing clouds, n:',n,nymin1,nxmin1,nz
556  !   create a cloud and rainout/washout field, clouds occur where rh>80%
557  !   total cloudheight is stored at level 0
558  do jy=0,nymin1
559    do ix=0,nxmin1
560      rain_cloud_above=0
561      lsp=lsprec(ix,jy,1,n)
562      convp=convprec(ix,jy,1,n)
563      cloudsh(ix,jy,n)=0
564      do kz_inv=1,nz-1
565         kz=nz-kz_inv+1
566         pressure=rho(ix,jy,kz,n)*r_air*tt(ix,jy,kz,n)
567         rh=qv(ix,jy,kz,n)/f_qvsat(pressure,tt(ix,jy,kz,n))
568         clouds(ix,jy,kz,n)=0
569         if (rh.gt.0.8) then ! in cloud
570           if ((lsp.gt.0.01).or.(convp.gt.0.01)) then ! cloud and precipitation
571              rain_cloud_above=1
572              cloudsh(ix,jy,n)=cloudsh(ix,jy,n)+height(kz)-height(kz-1)
573              if (lsp.ge.convp) then
574                 clouds(ix,jy,kz,n)=3 ! lsp dominated rainout
575              else
576                 clouds(ix,jy,kz,n)=2 ! convp dominated rainout
577              endif
578           else ! no precipitation
579             clouds(ix,jy,kz,n)=1 ! cloud
580           endif
581         else ! no cloud
582           if (rain_cloud_above.eq.1) then ! scavenging
583             if (lsp.ge.convp) then
584               clouds(ix,jy,kz,n)=5 ! lsp dominated washout
585             else
586               clouds(ix,jy,kz,n)=4 ! convp dominated washout
587             endif
588           endif
589         endif
590      end do
591    end do
592  end do
593  endif  ! IP & SEC 201812, GFS clouds read
594
595
596end subroutine verttransform_gfs
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG