source: flexpart.git/src/verttransform_nests.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: 21.5 KB
Line 
1subroutine verttransform_nests(n,uuhn,vvhn,wwhn,pvhn)
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!     It is similar to verttransform, but makes the transformations for      *
12!     the nested grids.                                                      *
13!                                                                            *
14!     Author: A. Stohl, G. Wotawa                                            *
15!                                                                            *
16!     12 August 1996                                                         *
17!     Update: 16 January 1998                                                *
18!                                                                            *
19!     Major update: 17 February 1999                                         *
20!     by G. Wotawa                                                           *
21!                                                                            *
22!     - Vertical levels for u, v and w are put together                      *
23!     - Slope correction for vertical velocity: Modification of calculation  *
24!       procedure                                                            *
25!                                                                            *
26!*****************************************************************************
27!  Changes, Bernd C. Krueger, Feb. 2001:       (marked "C-cv")
28!   Variables tthn and qvhn (on eta coordinates) from common block
29!*****************************************************************************
30! Sabine Eckhardt, March 2007
31! add the variable cloud for use with scavenging - descr. in com_mod
32!*****************************************************************************
33! ESO, 2016
34! -note that divide-by-zero occurs when nxmaxn,nymaxn etc. are larger than
35!  the actual field dimensions
36!*****************************************************************************
37! Date: 2017-05-30 modification of a bug in ew. Don Morton (CTBTO project)   *
38!*****************************************************************************
39!                                                                            *
40! Variables:                                                                 *
41! nxn,nyn,nuvz,nwz                field dimensions in x,y and z direction    *
42! uun                             wind components in x-direction [m/s]       *
43! vvn                             wind components in y-direction [m/s]       *
44! wwn                             wind components in z-direction [deltaeta/s]*
45! ttn                             temperature [K]                            *
46! pvn                             potential vorticity (pvu)                  *
47! psn                             surface pressure [Pa]                      *
48!                                                                            *
49!*****************************************************************************
50
51  use par_mod
52  use com_mod
53
54  implicit none
55
56  real,intent(in),dimension(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests) :: uuhn,vvhn,pvhn
57  real,intent(in),dimension(0:nxmaxn-1,0:nymaxn-1,nwzmax,maxnests) :: wwhn
58
59  real,dimension(0:nxmaxn-1,0:nymaxn-1,nuvzmax) :: rhohn,uvzlev,wzlev
60  real,dimension(0:nxmaxn-1,0:nymaxn-1,nzmax) :: pinmconv
61  real,dimension(0:nxmaxn-1,0:nymaxn-1) :: tvold,pold,pint,tv
62  real,dimension(0:nymaxn-1) :: cosf
63
64  integer,dimension(0:nxmaxn-1,0:nymaxn-1) :: rain_cloud_above, idx
65
66  integer :: ix,jy,kz,iz,n,l,kmin,kl,klp,ix1,jy1,ixp,jyp,kz_inv
67  real :: f_qvsat,pressure,rh,lsp,convp,cloudh_min,prec
68
69  real :: ew,dz1,dz2,dz
70  real :: dzdx,dzdy
71  real :: dzdx1,dzdx2,dzdy1,dzdy2
72  real,parameter :: const=r_air/ga
73  real :: tot_cloud_h
74  integer :: nxm1, nym1
75
76!  real,parameter :: precmin = 0.002 ! minimum prec in mm/h for cloud diagnostics
77
78! Loop over all nests
79!********************
80
81  do l=1,numbnests
82    nxm1=nxn(l)-1
83    nym1=nyn(l)-1
84
85! Loop over the whole grid
86!*************************
87
88    do jy=0,nyn(l)-1
89      do ix=0,nxn(l)-1
90        tvold(ix,jy)=tt2n(ix,jy,1,n,l)*(1.+0.378*ew(td2n(ix,jy,1,n,l))/ &
91             psn(ix,jy,1,n,l))
92      end do
93    end do
94
95    pold(0:nxm1,0:nym1)=psn(0:nxm1,0:nym1,1,n,l)
96    uvzlev(0:nxm1,0:nym1,1)=0.
97    wzlev(0:nxm1,0:nym1,1)=0.
98    rhohn(0:nxm1,0:nym1,1)=pold(0:nxm1,0:nym1)/(r_air*tvold(0:nxm1,0:nym1))
99
100! Compute heights of eta levels
101!******************************
102
103    do kz=2,nuvz
104      pint(0:nxm1,0:nym1)=akz(kz)+bkz(kz)*psn(0:nxm1,0:nym1,1,n,l)
105      tv(0:nxm1,0:nym1)=tthn(0:nxm1,0:nym1,kz,n,l)*(1.+0.608*qvhn(0:nxm1,0:nym1,kz,n,l))
106      rhohn(0:nxm1,0:nym1,kz)=pint(0:nxm1,0:nym1)/(r_air*tv(0:nxm1,0:nym1))
107
108      where (abs(tv(0:nxm1,0:nym1)-tvold(0:nxm1,0:nym1)).gt.0.2)
109        uvzlev(0:nxm1,0:nym1,kz)=uvzlev(0:nxm1,0:nym1,kz-1)+const*&
110             &log(pold(0:nxm1,0:nym1)/pint(0:nxm1,0:nym1))* &
111             (tv(0:nxm1,0:nym1)-tvold(0:nxm1,0:nym1))/&
112             &log(tv(0:nxm1,0:nym1)/tvold(0:nxm1,0:nym1))
113      elsewhere
114        uvzlev(0:nxm1,0:nym1,kz)=uvzlev(0:nxm1,0:nym1,kz-1)+const*&
115             &log(pold(0:nxm1,0:nym1)/pint(0:nxm1,0:nym1))*tv(0:nxm1,0:nym1)
116      endwhere
117
118      tvold(0:nxm1,0:nym1)=tv(0:nxm1,0:nym1)
119      pold(0:nxm1,0:nym1)=pint(0:nxm1,0:nym1)
120
121    end do
122
123    do kz=2,nwz-1
124      wzlev(0:nxm1,0:nym1,kz)=(uvzlev(0:nxm1,0:nym1,kz+1)+uvzlev(0:nxm1,0:nym1,kz))/2.
125    end do
126    wzlev(0:nxm1,0:nym1,nwz)=wzlev(0:nxm1,0:nym1,nwz-1)+ &
127         uvzlev(0:nxm1,0:nym1,nuvz)-uvzlev(0:nxm1,0:nym1,nuvz-1)
128
129
130    pinmconv(0:nxm1,0:nym1,1)=(uvzlev(0:nxm1,0:nym1,2))/ &
131         ((aknew(2)+bknew(2)*psn(0:nxm1,0:nym1,1,n,l))- &
132         (aknew(1)+bknew(1)*psn(0:nxm1,0:nym1,1,n,l)))
133    do kz=2,nz-1
134      pinmconv(0:nxm1,0:nym1,kz)=(uvzlev(0:nxm1,0:nym1,kz+1)-uvzlev(0:nxm1,0:nym1,kz-1))/ &
135           ((aknew(kz+1)+bknew(kz+1)*psn(0:nxm1,0:nym1,1,n,l))- &
136           (aknew(kz-1)+bknew(kz-1)*psn(0:nxm1,0:nym1,1,n,l)))
137    end do
138    pinmconv(0:nxm1,0:nym1,nz)=(uvzlev(0:nxm1,0:nym1,nz)-uvzlev(0:nxm1,0:nym1,nz-1))/ &
139         ((aknew(nz)+bknew(nz)*psn(0:nxm1,0:nym1,1,n,l))- &
140         (aknew(nz-1)+bknew(nz-1)*psn(0:nxm1,0:nym1,1,n,l)))
141
142! Levels, where u,v,t and q are given
143!************************************
144
145    uun(0:nxm1,0:nym1,1,n,l)=uuhn(0:nxm1,0:nym1,1,l)
146    vvn(0:nxm1,0:nym1,1,n,l)=vvhn(0:nxm1,0:nym1,1,l)
147    ttn(0:nxm1,0:nym1,1,n,l)=tthn(0:nxm1,0:nym1,1,n,l)
148    qvn(0:nxm1,0:nym1,1,n,l)=qvhn(0:nxm1,0:nym1,1,n,l)
149    if (readclouds_nest(l)) then
150      clwcn(0:nxm1,0:nym1,1,n,l)=clwchn(0:nxm1,0:nym1,1,n,l)
151      if (.not.sumclouds_nest(l)) ciwcn(0:nxm1,0:nym1,1,n,l)=ciwchn(0:nxm1,0:nym1,1,n,l)
152    end if
153    pvn(0:nxm1,0:nym1,1,n,l)=pvhn(0:nxm1,0:nym1,1,l)
154    rhon(0:nxm1,0:nym1,1,n,l)=rhohn(0:nxm1,0:nym1,1)
155
156    uun(0:nxm1,0:nym1,nz,n,l)=uuhn(0:nxm1,0:nym1,nuvz,l)
157    vvn(0:nxm1,0:nym1,nz,n,l)=vvhn(0:nxm1,0:nym1,nuvz,l)
158    ttn(0:nxm1,0:nym1,nz,n,l)=tthn(0:nxm1,0:nym1,nuvz,n,l)
159    qvn(0:nxm1,0:nym1,nz,n,l)=qvhn(0:nxm1,0:nym1,nuvz,n,l)
160    if (readclouds_nest(l)) then
161      clwcn(0:nxm1,0:nym1,nz,n,l)=clwchn(0:nxm1,0:nym1,nuvz,n,l)
162      if (.not.sumclouds_nest(l)) ciwcn(0:nxm1,0:nym1,nz,n,l)=ciwchn(0:nxm1,0:nym1,nuvz,n,l)
163    end if
164    pvn(0:nxm1,0:nym1,nz,n,l)=pvhn(0:nxm1,0:nym1,nuvz,l)
165    rhon(0:nxm1,0:nym1,nz,n,l)=rhohn(0:nxm1,0:nym1,nuvz)
166
167
168    kmin=2
169    idx=kmin
170    do iz=2,nz-1
171      do jy=0,nyn(l)-1
172        do ix=0,nxn(l)-1
173          if(height(iz).gt.uvzlev(ix,jy,nuvz)) then
174            uun(ix,jy,iz,n,l)=uun(ix,jy,nz,n,l)
175            vvn(ix,jy,iz,n,l)=vvn(ix,jy,nz,n,l)
176            ttn(ix,jy,iz,n,l)=ttn(ix,jy,nz,n,l)
177            qvn(ix,jy,iz,n,l)=qvn(ix,jy,nz,n,l)
178!hg adding the cloud water
179            if (readclouds_nest(l)) then
180              clwcn(ix,jy,iz,n,l)=clwcn(ix,jy,nz,n,l)
181              if (.not.sumclouds_nest(l)) ciwcn(ix,jy,iz,n,l)=ciwcn(ix,jy,nz,n,l)
182            end if
183!hg
184            pvn(ix,jy,iz,n,l)=pvn(ix,jy,nz,n,l)
185            rhon(ix,jy,iz,n,l)=rhon(ix,jy,nz,n,l)
186          else
187            innuvz: do kz=idx(ix,jy),nuvz
188              if (idx(ix,jy) .le. kz .and. (height(iz).gt.uvzlev(ix,jy,kz-1)).and. &
189                   (height(iz).le.uvzlev(ix,jy,kz))) then
190                idx(ix,jy)=kz
191                exit innuvz
192              endif
193            enddo innuvz
194          endif
195        enddo
196      enddo
197      do jy=0,nyn(l)-1
198        do ix=0,nxn(l)-1
199          if(height(iz).le.uvzlev(ix,jy,nuvz)) then
200            kz=idx(ix,jy)
201            dz1=height(iz)-uvzlev(ix,jy,kz-1)
202            dz2=uvzlev(ix,jy,kz)-height(iz)
203            dz=dz1+dz2
204            uun(ix,jy,iz,n,l)=(uuhn(ix,jy,kz-1,l)*dz2+uuhn(ix,jy,kz,l)*dz1)/dz
205            vvn(ix,jy,iz,n,l)=(vvhn(ix,jy,kz-1,l)*dz2+vvhn(ix,jy,kz,l)*dz1)/dz
206            ttn(ix,jy,iz,n,l)=(tthn(ix,jy,kz-1,n,l)*dz2 &
207                 +tthn(ix,jy,kz,n,l)*dz1)/dz
208            qvn(ix,jy,iz,n,l)=(qvhn(ix,jy,kz-1,n,l)*dz2 &
209                 +qvhn(ix,jy,kz,n,l)*dz1)/dz
210!hg adding the cloud water
211            if (readclouds_nest(l)) then
212              clwcn(ix,jy,iz,n,l)=(clwchn(ix,jy,kz-1,n,l)*dz2+clwchn(ix,jy,kz,n,l)*dz1)/dz
213              if (.not.sumclouds_nest(l)) &
214                   &ciwcn(ix,jy,iz,n,l)=(ciwchn(ix,jy,kz-1,n,l)*dz2+ciwchn(ix,jy,kz,n,l)*dz1)/dz
215            end if
216!hg
217            pvn(ix,jy,iz,n,l)=(pvhn(ix,jy,kz-1,l)*dz2+pvhn(ix,jy,kz,l)*dz1)/dz
218            rhon(ix,jy,iz,n,l)=(rhohn(ix,jy,kz-1)*dz2+rhohn(ix,jy,kz)*dz1)/dz
219          endif
220        enddo
221      enddo
222    enddo
223
224
225
226!         do kz=kmin,nuvz
227!           if(height(iz).gt.uvwzlev(ix,jy,nuvz)) then
228!             uun(ix,jy,iz,n,l)=uun(ix,jy,nz,n,l)
229!             vvn(ix,jy,iz,n,l)=vvn(ix,jy,nz,n,l)
230!             ttn(ix,jy,iz,n,l)=ttn(ix,jy,nz,n,l)
231!             qvn(ix,jy,iz,n,l)=qvn(ix,jy,nz,n,l)
232!             pvn(ix,jy,iz,n,l)=pvn(ix,jy,nz,n,l)
233!             rhon(ix,jy,iz,n,l)=rhon(ix,jy,nz,n,l)
234!             goto 30
235!           endif
236!           if ((height(iz).gt.uvwzlev(ix,jy,kz-1)).and. &
237!           (height(iz).le.uvwzlev(ix,jy,kz))) then
238!            dz1=height(iz)-uvwzlev(ix,jy,kz-1)
239!            dz2=uvwzlev(ix,jy,kz)-height(iz)
240!            dz=dz1+dz2
241!            uun(ix,jy,iz,n,l)=(uuhn(ix,jy,kz-1,l)*dz2+ &
242!            uuhn(ix,jy,kz,l)*dz1)/dz
243!            vvn(ix,jy,iz,n,l)=(vvhn(ix,jy,kz-1,l)*dz2+ &
244!            vvhn(ix,jy,kz,l)*dz1)/dz
245!            ttn(ix,jy,iz,n,l)=(tthn(ix,jy,kz-1,n,l)*dz2+ &
246!            tthn(ix,jy,kz,n,l)*dz1)/dz
247!            qvn(ix,jy,iz,n,l)=(qvhn(ix,jy,kz-1,n,l)*dz2+ &
248!            qvhn(ix,jy,kz,n,l)*dz1)/dz
249!            pvn(ix,jy,iz,n,l)=(pvhn(ix,jy,kz-1,l)*dz2+ &
250!            pvhn(ix,jy,kz,l)*dz1)/dz
251!            rhon(ix,jy,iz,n,l)=(rhohn(kz-1)*dz2+rhohn(kz)*dz1)/dz
252!            kmin=kz
253!            goto 30
254!           endif
255!         end do
256! 30      continue
257!       end do
258
259
260! Levels, where w is given
261!*************************
262
263    wwn(0:nxm1,0:nym1,1,n,l)=wwhn(0:nxm1,0:nym1,1,l)*pinmconv(0:nxm1,0:nym1,1)
264    wwn(0:nxm1,0:nym1,nz,n,l)=wwhn(0:nxm1,0:nym1,nwz,l)*pinmconv(0:nxm1,0:nym1,nz)
265    kmin=2
266    idx=kmin
267    do iz=2,nz
268      do jy=0,nym1
269        do ix=0,nxm1
270          inn:         do kz=idx(ix,jy),nwz
271            if(idx(ix,jy) .le. kz .and. height(iz).gt.wzlev(ix,jy,kz-1).and. &
272                 height(iz).le.wzlev(ix,jy,kz)) then
273              idx(ix,jy)=kz
274              exit inn
275            endif
276          enddo inn
277        enddo
278      enddo
279      do jy=0,nym1
280        do ix=0,nxm1
281          kz=idx(ix,jy)
282          dz1=height(iz)-wzlev(ix,jy,kz-1)
283          dz2=wzlev(ix,jy,kz)-height(iz)
284          dz=dz1+dz2
285          wwn(ix,jy,iz,n,l)=(wwhn(ix,jy,kz-1,l)*pinmconv(ix,jy,kz-1)*dz2 &
286               +wwhn(ix,jy,kz,l)*pinmconv(ix,jy,kz)*dz1)/dz
287        enddo
288      enddo
289    enddo
290
291!       wwn(ix,jy,1,n,l)=wwhn(ix,jy,1,l)*pinmconv(1)
292!       wwn(ix,jy,nz,n,l)=wwhn(ix,jy,nwz,l)*pinmconv(nz)
293!       kmin=2
294!       do iz=2,nz
295!         do kz=kmin,nwz
296!           if ((height(iz).gt.wzlev(kz-1)).and. &
297!           (height(iz).le.wzlev(kz))) then
298!             dz1=height(iz)-wzlev(kz-1)
299!             dz2=wzlev(kz)-height(iz)
300!             dz=dz1+dz2
301!             wwn(ix,jy,iz,n,l)=(wwhn(ix,jy,kz-1,l)*pinmconv(kz-1)*dz2 &
302!             +wwhn(ix,jy,kz,l)*pinmconv(kz)*dz1)/dz
303!             kmin=kz
304!             goto 40
305!           endif
306!         end do
307! 40      continue
308!       end do
309
310! Compute density gradients at intermediate levels
311!*************************************************
312
313    drhodzn(0:nxm1,0:nym1,1,n,l)=(rhon(0:nxm1,0:nym1,2,n,l)-rhon(0:nxm1,0:nym1,1,n,l))/ &
314         (height(2)-height(1))
315    do kz=2,nz-1
316      drhodzn(0:nxm1,0:nym1,kz,n,l)=(rhon(0:nxm1,0:nym1,kz+1,n,l)-rhon(0:nxm1,0:nym1,kz-1,n,l))/ &
317           (height(kz+1)-height(kz-1))
318    end do
319    drhodzn(0:nxm1,0:nym1,nz,n,l)=drhodzn(0:nxm1,0:nym1,nz-1,n,l)
320
321
322! drhodzn(ix,jy,1,n,l)=(rhon(ix,jy,2,n,l)-rhon(ix,jy,1,n,l))/ &
323! (height(2)-height(1))
324! do kz=2,nz-1
325!   drhodzn(ix,jy,kz,n,l)=(rhon(ix,jy,kz+1,n,l)- &
326!   rhon(ix,jy,kz-1,n,l))/(height(kz+1)-height(kz-1))
327! end do
328! drhodzn(ix,jy,nz,n,l)=drhodzn(ix,jy,nz-1,n,l)
329
330
331
332!   end do
333! end do
334
335
336!****************************************************************
337! Compute slope of eta levels in windward direction and resulting
338! vertical wind correction
339!****************************************************************
340
341    do jy=1,nyn(l)-2
342!    cosf=cos((real(jy)*dyn(l)+ylat0n(l))*pi180)
343      cosf(jy)=1./cos((real(jy)*dyn(l)+ylat0n(l))*pi180)
344    end do
345
346    kmin=2
347    idx=kmin
348    do iz=2,nz-1
349      do jy=1,nyn(l)-2
350        do ix=1,nxn(l)-2
351
352          inneta: do kz=idx(ix,jy),nz
353            if (idx(ix,jy) .le. kz .and. (height(iz).gt.uvzlev(ix,jy,kz-1)).and. &
354                 (height(iz).le.uvzlev(ix,jy,kz))) then
355              idx(ix,jy)=kz
356              exit inneta
357            endif
358          enddo inneta
359        enddo
360      enddo
361
362      do jy=1,nyn(l)-2
363        do ix=1,nxn(l)-2
364          kz=idx(ix,jy)
365          dz1=height(iz)-uvzlev(ix,jy,kz-1)
366          dz2=uvzlev(ix,jy,kz)-height(iz)
367          dz=dz1+dz2
368          ix1=ix-1
369          jy1=jy-1
370          ixp=ix+1
371          jyp=jy+1
372
373          dzdx1=(uvzlev(ixp,jy,kz-1)-uvzlev(ix1,jy,kz-1))/2.
374          dzdx2=(uvzlev(ixp,jy,kz)-uvzlev(ix1,jy,kz))/2.
375          dzdx=(dzdx1*dz2+dzdx2*dz1)/dz
376
377          dzdy1=(uvzlev(ix,jyp,kz-1)-uvzlev(ix,jy1,kz-1))/2.
378          dzdy2=(uvzlev(ix,jyp,kz)-uvzlev(ix,jy1,kz))/2.
379          dzdy=(dzdy1*dz2+dzdy2*dz1)/dz
380
381          wwn(ix,jy,iz,n,l)=wwn(ix,jy,iz,n,l)+(dzdx*uun(ix,jy,iz,n,l)*dxconst*xresoln(l)*cosf(jy)+&
382               &dzdy*vvn(ix,jy,iz,n,l)*dyconst*yresoln(l))
383
384        end do
385
386      end do
387    end do
388
389
390!   do jy=1,nyn(l)-2
391!     do ix=1,nxn(l)-2
392!       kmin=2
393!       do iz=2,nz-1
394
395!         ui=uun(ix,jy,iz,n,l)*dxconst*xresoln(l)/cosf(jy)
396!         vi=vvn(ix,jy,iz,n,l)*dyconst*yresoln(l)
397
398!         do kz=kmin,nz
399!           if ((height(iz).gt.uvwzlev(ix,jy,kz-1)).and. &
400!                (height(iz).le.uvwzlev(ix,jy,kz))) then
401!             dz1=height(iz)-uvwzlev(ix,jy,kz-1)
402!             dz2=uvwzlev(ix,jy,kz)-height(iz)
403!             dz=dz1+dz2
404!             kl=kz-1
405!             klp=kz
406!             kmin=kz
407!             goto 47
408!           endif
409!         end do
410
411! 47      ix1=ix-1
412!         jy1=jy-1
413!         ixp=ix+1
414!         jyp=jy+1
415
416!         dzdx1=(uvwzlev(ixp,jy,kl)-uvwzlev(ix1,jy,kl))/2.
417!         dzdx2=(uvwzlev(ixp,jy,klp)-uvwzlev(ix1,jy,klp))/2.
418!         dzdx=(dzdx1*dz2+dzdx2*dz1)/dz
419
420!         dzdy1=(uvwzlev(ix,jyp,kl)-uvwzlev(ix,jy1,kl))/2.
421!         dzdy2=(uvwzlev(ix,jyp,klp)-uvwzlev(ix,jy1,klp))/2.
422!         dzdy=(dzdy1*dz2+dzdy2*dz1)/dz
423
424!         wwn(ix,jy,iz,n,l)=wwn(ix,jy,iz,n,l)+(dzdx*ui+dzdy*vi)
425
426!       end do
427
428!     end do
429!   end do
430
431
432!write (*,*) 'initializing nested cloudsn, n:',n
433!   create a cloud and rainout/washout field, cloudsn occur where rh>80%
434
435
436!*********************************************************************************** 
437    if (readclouds_nest(l)) then !HG METHOD
438! The method is loops all grids vertically and constructs the 3D matrix for clouds
439! Cloud top and cloud bottom gid cells are assigned as well as the total column
440! cloud water. For precipitating grids, the type and whether it is in or below
441! cloud scavenging are assigned with numbers 2-5 (following the old metod).
442! Distinction is done for lsp and convp though they are treated the same in regards
443! to scavenging. Also clouds that are not precipitating are defined which may be
444! to include future cloud processing by non-precipitating-clouds.
445!***********************************************************************************
446      write(*,*) 'Nested ECMWF fields: using cloud water'
447      clwn(0:nxm1,0:nym1,:,n,l)=0.0
448!    icloud_stats(0:nxm1,0:nym1,:,n)=0.0
449      ctwcn(0:nxm1,0:nym1,n,l)=0.0
450      cloudsn(0:nxm1,0:nym1,:,n,l)=0
451! If water/ice are read separately into clwc and ciwc, store sum in clwcn
452      if (.not.sumclouds_nest(l)) then
453        clwcn(0:nxm1,0:nym1,:,n,l) = clwcn(0:nxm1,0:nym1,:,n,l) + ciwcn(0:nxm1,0:nym1,:,n,l)
454      end if
455      do jy=0,nym1
456        do ix=0,nxm1
457          lsp=lsprecn(ix,jy,1,n,l)
458          convp=convprecn(ix,jy,1,n,l)
459          prec=lsp+convp
460          tot_cloud_h=0
461! Find clouds in the vertical
462          do kz=1, nz-1 !go from top to bottom
463            if (clwcn(ix,jy,kz,n,l).gt.0) then     
464! assuming rho is in kg/m3 and hz in m gives: kg/kg * kg/m3 *m3/kg /m = m2/m3
465              clwn(ix,jy,kz,n,l)=(clwcn(ix,jy,kz,n,l)*rhon(ix,jy,kz,n,l))*(height(kz+1)-height(kz))
466              tot_cloud_h=tot_cloud_h+(height(kz+1)-height(kz))
467              ctwcn(ix,jy,n,l) = ctwcn(ix,jy,n,l)+clwn(ix,jy,kz,n,l)
468!            icloud_stats(ix,jy,4,n)= icloud_stats(ix,jy,4,n)+clw(ix,jy,kz,n)          ! Column cloud water [m3/m3]
469!           icloud_stats(ix,jy,3,n)= min(height(kz+1),height(kz))                     ! Cloud BOT height stats      [m]
470              cloudh_min=min(height(kz+1),height(kz))
471!ZHG 2015 extra for testing
472!            clh(ix,jy,kz,n)=height(kz+1)-height(kz)
473!            icloud_stats(ix,jy,1,n)=icloud_stats(ix,jy,1,n)+(height(kz+1)-height(kz)) ! Cloud total vertical extent [m]
474!            icloud_stats(ix,jy,2,n)= max(icloud_stats(ix,jy,2,n),height(kz))          ! Cloud TOP height            [m]
475!ZHG
476            endif
477          end do
478
479! If Precipitation. Define removal type in the vertical
480          if ((lsp.gt.0.01).or.(convp.gt.0.01)) then ! cloud and precipitation
481
482            do kz=nz,1,-1 !go Bottom up!
483              if (clwn(ix,jy,kz,n,l).gt.0.0) then ! is in cloud
484                cloudshn(ix,jy,n,l)=cloudshn(ix,jy,n,l)+height(kz)-height(kz-1)
485                cloudsn(ix,jy,kz,n,l)=1                               ! is a cloud
486                if (lsp.ge.convp) then
487                  cloudsn(ix,jy,kz,n,l)=3                            ! lsp in-cloud
488                else
489                  cloudsn(ix,jy,kz,n,l)=2                             ! convp in-cloud
490                endif                                              ! convective or large scale
491              else if((clwn(ix,jy,kz,n,l).le.0.0) .and. (cloudh_min.ge.height(kz))) then ! is below cloud
492                if (lsp.ge.convp) then
493                  cloudsn(ix,jy,kz,n,l)=5                             ! lsp dominated washout
494                else
495                  cloudsn(ix,jy,kz,n,l)=4                             ! convp dominated washout
496                endif                                              ! convective or large scale
497              endif
498
499              if (height(kz).ge. 19000) then                        ! set a max height for removal
500                cloudsn(ix,jy,kz,n,l)=0
501              endif !clw>0
502            end do !nz
503          endif ! precipitation
504        end do
505      end do
506!********************************************************************
507    else ! old method:
508!********************************************************************
509      write(*,*) 'Nested fields: using cloud water from Parameterization'
510      do jy=0,nyn(l)-1
511        do ix=0,nxn(l)-1
512          rain_cloud_above(ix,jy)=0
513          lsp=lsprecn(ix,jy,1,n,l)
514          convp=convprecn(ix,jy,1,n,l)
515          cloudshn(ix,jy,n,l)=0
516          do kz_inv=1,nz-1
517            kz=nz-kz_inv+1
518            pressure=rhon(ix,jy,kz,n,l)*r_air*ttn(ix,jy,kz,n,l)
519            rh=qvn(ix,jy,kz,n,l)/f_qvsat(pressure,ttn(ix,jy,kz,n,l))
520            cloudsn(ix,jy,kz,n,l)=0
521            if (rh.gt.0.8) then ! in cloud
522              if ((lsp.gt.0.01).or.(convp.gt.0.01)) then
523                rain_cloud_above(ix,jy)=1
524                cloudshn(ix,jy,n,l)=cloudshn(ix,jy,n,l)+ &
525                     height(kz)-height(kz-1)
526                if (lsp.ge.convp) then
527                  cloudsn(ix,jy,kz,n,l)=3 ! lsp dominated rainout
528                else
529                  cloudsn(ix,jy,kz,n,l)=2 ! convp dominated rainout
530                endif
531              else ! no precipitation
532                cloudsn(ix,jy,kz,n,l)=1 ! cloud
533              endif
534            else ! no cloud
535              if (rain_cloud_above(ix,jy).eq.1) then ! scavenging
536                if (lsp.ge.convp) then
537                  cloudsn(ix,jy,kz,n,l)=5 ! lsp dominated washout
538                else
539                  cloudsn(ix,jy,kz,n,l)=4 ! convp dominated washout
540                endif
541              endif
542            endif
543          end do
544        end do
545      end do
546    end if
547
548  end do ! end loop over nests
549
550end subroutine verttransform_nests
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG