source: flexpart.git/src/verttransform.f90 @ 4d45639

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 4d45639 was 4d45639, checked in by Espen Sollum ATMOS <eso@…>, 8 years ago

Merged in optimized verttransform.f90 from author Leopold Haimberger

  • Property mode set to 100644
File size: 24.0 KB
Line 
1!**********************************************************************
2! Copyright 1998,1999,2000,2001,2002,2005,2007,2008,2009,2010         *
3! Andreas Stohl, Petra Seibert, A. Frank, Gerhard Wotawa,             *
4! Caroline Forster, Sabine Eckhardt, John Burkhart, Harald Sodemann   *
5!                                                                     *
6! This file is part of FLEXPART.                                      *
7!                                                                     *
8! FLEXPART is free software: you can redistribute it and/or modify    *
9! it under the terms of the GNU General Public License as published by*
10! the Free Software Foundation, either version 3 of the License, or   *
11! (at your option) any later version.                                 *
12!                                                                     *
13! FLEXPART is distributed in the hope that it will be useful,         *
14! but WITHOUT ANY WARRANTY; without even the implied warranty of      *
15! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the       *
16! GNU General Public License for more details.                        *
17!                                                                     *
18! You should have received a copy of the GNU General Public License   *
19! along with FLEXPART.  If not, see <http://www.gnu.org/licenses/>.   *
20!**********************************************************************
21
22subroutine verttransform(n,uuh,vvh,wwh,pvh)
23!                         i  i   i   i   i
24!*****************************************************************************
25!                                                                            *
26!     This subroutine transforms temperature, dew point temperature and      *
27!     wind components from eta to meter coordinates.                         *
28!     The vertical wind component is transformed from Pa/s to m/s using      *
29!     the conversion factor pinmconv.                                        *
30!     In addition, this routine calculates vertical density gradients        *
31!     needed for the parameterization of the turbulent velocities.           *
32!                                                                            *
33!     Author: A. Stohl, G. Wotawa                                            *
34!                                                                            *
35!     12 August 1996                                                         *
36!     Update: 16 January 1998                                                *
37!                                                                            *
38!     Major update: 17 February 1999                                         *
39!     by G. Wotawa                                                           *
40!                                                                            *
41!     - Vertical levels for u, v and w are put together                      *
42!     - Slope correction for vertical velocity: Modification of calculation  *
43!       procedure                                                            *
44!                                                                            *
45!*****************************************************************************
46!  Changes, Bernd C. Krueger, Feb. 2001:
47!   Variables tth and qvh (on eta coordinates) from common block
48!*****************************************************************************
49! Sabine Eckhardt, March 2007
50! added the variable cloud for use with scavenging - descr. in com_mod
51!*****************************************************************************
52!                                                                            *
53! Variables:                                                                 *
54! nx,ny,nz                        field dimensions in x,y and z direction    *
55! clouds(0:nxmax,0:nymax,0:nzmax,numwfmem) cloud field for wet deposition    *
56! uu(0:nxmax,0:nymax,nzmax,numwfmem)     wind components in x-direction [m/s]*
57! vv(0:nxmax,0:nymax,nzmax,numwfmem)     wind components in y-direction [m/s]*
58! ww(0:nxmax,0:nymax,nzmax,numwfmem)     wind components in z-direction      *
59!                                          [deltaeta/s]                      *
60! tt(0:nxmax,0:nymax,nzmax,numwfmem)     temperature [K]                     *
61! pv(0:nxmax,0:nymax,nzmax,numwfmem)     potential voriticity (pvu)          *
62! ps(0:nxmax,0:nymax,numwfmem)           surface pressure [Pa]               *
63!                                                                            *
64!*****************************************************************************
65
66  use par_mod
67  use com_mod
68  use cmapf_mod, only: cc2gll
69! eso TODO:
70! only used for timing of CPU measurement. Remove this (and calls to mpif_mtime below)
71! as this routine should not be dependent on MPI
72!  use mpi_mod
73! :TODO
74
75  implicit none
76
77  integer :: ix,jy,kz,iz,n,kmin,kl,klp,ix1,jy1,ixp,jyp,ixm,jym
78  integer :: rain_cloud_above(0:nxmax-1,0:nymax-1),kz_inv,idx(0:nxmax-1,0:nymax-1)
79  real :: f_qvsat,pressure,cosf(0:nymax-1)
80  real :: rh,lsp,convp,tim,tim2,rhmin,precmin,prec
81  real :: uvzlev(0:nxmax-1,0:nymax-1,nuvzmax),rhoh(0:nxmax-1,0:nymax-1,nuvzmax)
82  real :: pinmconv(0:nxmax-1,0:nymax-1,nzmax)
83  real :: ew,pint(0:nxmax-1,0:nymax-1),tv(0:nxmax-1,0:nymax-1)
84  real :: tvold(0:nxmax-1,0:nymax-1),pold(0:nxmax-1,0:nymax-1),dz1,dz2,dz,ui,vi
85  real :: xlon,ylat,xlonr,dzdx,dzdy
86  real :: dzdx1,dzdx2,dzdy1,dzdy2
87  real :: uuaux,vvaux,uupolaux,vvpolaux,ddpol,ffpol,wdummy
88  real :: uuh(0:nxmax-1,0:nymax-1,nuvzmax)
89  real :: vvh(0:nxmax-1,0:nymax-1,nuvzmax)
90  real :: pvh(0:nxmax-1,0:nymax-1,nuvzmax)
91  real :: wwh(0:nxmax-1,0:nymax-1,nwzmax)
92  real :: wzlev(0:nxmax-1,0:nymax-1,nwzmax)
93  real,parameter :: const=r_air/ga
94!  integer:: ihr, imin, isec, i100th,ihr2, imin2, isec2, i100th2
95  parameter (precmin = 0.002) ! minimum prec in mm/h for cloud diagnostics
96
97  logical :: init = .true.
98
99!hg
100  integer :: cloud_ver,cloud_min, cloud_max
101  real :: cloud_col_wat, cloud_water
102!hg temporary variables for testing
103  real :: rcw(0:nxmax-1,0:nymax-1)
104  real :: rpc(0:nxmax-1,0:nymax-1)
105!hg
106
107!*************************************************************************
108! If verttransform is called the first time, initialize heights of the   *
109! z levels in meter. The heights are the heights of model levels, where  *
110! u,v,T and qv are given, and of the interfaces, where w is given. So,   *
111! the vertical resolution in the z system is doubled. As reference point,*
112! the lower left corner of the grid is used.                             *
113! Unlike in the eta system, no difference between heights for u,v and    *
114! heights for w exists.                                                  *
115!*************************************************************************
116
117
118!eso measure CPU time
119!  call mpif_mtime('verttransform',0)
120
121  if (init) then
122
123! Search for a point with high surface pressure (i.e. not above significant topography)
124! Then, use this point to construct a reference z profile, to be used at all times
125!*****************************************************************************
126
127    do jy=0,nymin1
128      do ix=0,nxmin1
129        if (ps(ix,jy,1,n).gt.100000.) then
130          ixm=ix
131          jym=jy
132          goto 3
133        endif
134      end do
135    end do
1363   continue
137
138
139    tvold(ixm,jym)=tt2(ixm,jym,1,n)*(1.+0.378*ew(td2(ixm,jym,1,n))/ &
140         ps(ixm,jym,1,n))
141    pold(ixm,jym)=ps(ixm,jym,1,n)
142    height(1)=0.
143
144    do kz=2,nuvz
145      pint=akz(kz)+bkz(kz)*ps(ixm,jym,1,n)
146      tv=tth(ixm,jym,kz,n)*(1.+0.608*qvh(ixm,jym,kz,n))
147
148      if (abs(tv(ixm,jym)-tvold(ixm,jym)).gt.0.2) then
149        height(kz)= &
150             height(kz-1)+const*log(pold(ixm,jym)/pint(ixm,jym))* &
151             (tv(ixm,jym)-tvold(ixm,jym))/log(tv(ixm,jym)/tvold(ixm,jym))
152      else
153        height(kz)=height(kz-1)+ &
154             const*log(pold(ixm,jym)/pint(ixm,jym))*tv(ixm,jym)
155      endif
156
157      tvold(ixm,jym)=tv(ixm,jym)
158      pold(ixm,jym)=pint(ixm,jym)
159    end do
160
161
162! Determine highest levels that can be within PBL
163!************************************************
164
165    do kz=1,nz
166      if (height(kz).gt.hmixmax) then
167        nmixz=kz
168        goto 9
169      endif
170    end do
1719   continue
172
173! Do not repeat initialization of the Cartesian z grid
174!*****************************************************
175
176    init=.false.
177
178  endif
179
180
181! Loop over the whole grid
182!*************************
183
184
185  do jy=0,nymin1
186    do ix=0,nxmin1
187      tvold(ix,jy)=tt2(ix,jy,1,n)*(1.+0.378*ew(td2(ix,jy,1,n))/ &
188           ps(ix,jy,1,n))
189    enddo
190  enddo
191  pold=ps(:,:,1,n)
192  uvzlev(:,:,1)=0.
193  wzlev(:,:,1)=0.
194  rhoh(:,:,1)=pold/(r_air*tvold)
195
196
197! Compute heights of eta levels
198!******************************
199
200  do kz=2,nuvz
201    pint=akz(kz)+bkz(kz)*ps(:,:,1,n)
202    tv=tth(:,:,kz,n)*(1.+0.608*qvh(:,:,kz,n))
203    rhoh(:,:,kz)=pint(:,:)/(r_air*tv)
204
205    where (abs(tv-tvold).gt.0.2)
206      uvzlev(:,:,kz)=uvzlev(:,:,kz-1)+const*log(pold/pint)* &
207           (tv-tvold)/log(tv/tvold)
208    elsewhere
209      uvzlev(:,:,kz)=uvzlev(:,:,kz-1)+const*log(pold/pint)*tv
210    endwhere
211
212    tvold=tv
213    pold=pint
214  end do
215
216
217  do kz=2,nwz-1
218    wzlev(:,:,kz)=(uvzlev(:,:,kz+1)+uvzlev(:,:,kz))/2.
219  end do
220  wzlev(:,:,nwz)=wzlev(:,:,nwz-1)+ &
221       uvzlev(:,:,nuvz)-uvzlev(:,:,nuvz-1)
222
223! pinmconv=(h2-h1)/(p2-p1)
224
225  pinmconv(:,:,1)=(uvzlev(:,:,2))/ &
226       ((aknew(2)+bknew(2)*ps(:,:,1,n))- &
227       (aknew(1)+bknew(1)*ps(:,:,1,n)))
228  do kz=2,nz-1
229    pinmconv(:,:,kz)=(uvzlev(:,:,kz+1)-uvzlev(:,:,kz-1))/ &
230         ((aknew(kz+1)+bknew(kz+1)*ps(:,:,1,n))- &
231         (aknew(kz-1)+bknew(kz-1)*ps(:,:,1,n)))
232  end do
233  pinmconv(:,:,nz)=(uvzlev(:,:,nz)-uvzlev(:,:,nz-1))/ &
234       ((aknew(nz)+bknew(nz)*ps(:,:,1,n))- &
235       (aknew(nz-1)+bknew(nz-1)*ps(:,:,1,n)))
236
237! Levels, where u,v,t and q are given
238!************************************
239
240
241  uu(:,:,1,n)=uuh(:,:,1)
242  vv(:,:,1,n)=vvh(:,:,1)
243  tt(:,:,1,n)=tth(:,:,1,n)
244  qv(:,:,1,n)=qvh(:,:,1,n)
245!hg adding the cloud water
246  clwc(:,:,1,n)=clwch(:,:,1,n)
247  ciwc(:,:,1,n)=ciwch(:,:,1,n)   
248!hg
249  pv(:,:,1,n)=pvh(:,:,1)
250  rho(:,:,1,n)=rhoh(:,:,1)
251  uu(:,:,nz,n)=uuh(:,:,nuvz)
252  vv(:,:,nz,n)=vvh(:,:,nuvz)
253  tt(:,:,nz,n)=tth(:,:,nuvz,n)
254  qv(:,:,nz,n)=qvh(:,:,nuvz,n)
255
256!hg adding the cloud water
257  clwc(:,:,nz,n)=clwch(:,:,nuvz,n)
258  ciwc(:,:,nz,n)=ciwch(:,:,nuvz,n)
259!hg
260  pv(:,:,nz,n)=pvh(:,:,nuvz)
261  rho(:,:,nz,n)=rhoh(:,:,nuvz)
262
263
264  kmin=2
265  idx=kmin
266  do iz=2,nz-1
267    do jy=0,nymin1
268      do ix=0,nxmin1
269        if(height(iz).gt.uvzlev(ix,jy,nuvz)) then
270          uu(ix,jy,iz,n)=uu(ix,jy,nz,n)
271          vv(ix,jy,iz,n)=vv(ix,jy,nz,n)
272          tt(ix,jy,iz,n)=tt(ix,jy,nz,n)
273          qv(ix,jy,iz,n)=qv(ix,jy,nz,n)
274!hg adding the cloud water
275          clwc(ix,jy,iz,n)=clwc(ix,jy,nz,n)
276          ciwc(ix,jy,iz,n)=ciwc(ix,jy,nz,n)
277!hg
278          pv(ix,jy,iz,n)=pv(ix,jy,nz,n)
279          rho(ix,jy,iz,n)=rho(ix,jy,nz,n)
280        else
281          innuvz: do kz=idx(ix,jy),nuvz
282            if (idx(ix,jy) .le. kz .and. (height(iz).gt.uvzlev(ix,jy,kz-1)).and. &
283                 (height(iz).le.uvzlev(ix,jy,kz))) then
284              idx(ix,jy)=kz
285              exit innuvz
286            endif
287          enddo innuvz
288        endif
289      enddo
290    enddo
291    do jy=0,nymin1
292      do ix=0,nxmin1
293        if(height(iz).le.uvzlev(ix,jy,nuvz)) then
294          kz=idx(ix,jy)
295          dz1=height(iz)-uvzlev(ix,jy,kz-1)
296          dz2=uvzlev(ix,jy,kz)-height(iz)
297          dz=dz1+dz2
298          uu(ix,jy,iz,n)=(uuh(ix,jy,kz-1)*dz2+uuh(ix,jy,kz)*dz1)/dz
299          vv(ix,jy,iz,n)=(vvh(ix,jy,kz-1)*dz2+vvh(ix,jy,kz)*dz1)/dz
300          tt(ix,jy,iz,n)=(tth(ix,jy,kz-1,n)*dz2 &
301               +tth(ix,jy,kz,n)*dz1)/dz
302          qv(ix,jy,iz,n)=(qvh(ix,jy,kz-1,n)*dz2 &
303               +qvh(ix,jy,kz,n)*dz1)/dz
304!hg adding the cloud water
305          clwc(ix,jy,iz,n)=(clwch(ix,jy,kz-1,n)*dz2+clwch(ix,jy,kz,n)*dz1)/dz
306          ciwc(ix,jy,iz,n)=(ciwch(ix,jy,kz-1,n)*dz2+ciwch(ix,jy,kz,n)*dz1)/dz
307!hg
308          pv(ix,jy,iz,n)=(pvh(ix,jy,kz-1)*dz2+pvh(ix,jy,kz)*dz1)/dz
309          rho(ix,jy,iz,n)=(rhoh(ix,jy,kz-1)*dz2+rhoh(ix,jy,kz)*dz1)/dz
310        endif
311      enddo
312    enddo
313  enddo
314
315
316! Levels, where w is given
317!*************************
318
319  ww(:,:,1,n)=wwh(:,:,1)*pinmconv(:,:,1)
320  ww(:,:,nz,n)=wwh(:,:,nwz)*pinmconv(:,:,nz)
321  kmin=2
322  idx=kmin
323  do iz=2,nz
324    do jy=0,nymin1
325      do ix=0,nxmin1
326        inn:         do kz=idx(ix,jy),nwz
327          if(idx(ix,jy) .le. kz .and. height(iz).gt.wzlev(ix,jy,kz-1).and. &
328               height(iz).le.wzlev(ix,jy,kz)) then
329            idx(ix,jy)=kz
330            exit inn
331          endif
332        enddo inn
333      enddo
334    enddo
335    do jy=0,nymin1
336      do ix=0,nxmin1
337        kz=idx(ix,jy)
338        dz1=height(iz)-wzlev(ix,jy,kz-1)
339        dz2=wzlev(ix,jy,kz)-height(iz)
340        dz=dz1+dz2
341        ww(ix,jy,iz,n)=(wwh(ix,jy,kz-1)*pinmconv(ix,jy,kz-1)*dz2 &
342             +wwh(ix,jy,kz)*pinmconv(ix,jy,kz)*dz1)/dz
343      enddo
344    enddo
345  enddo
346
347! Compute density gradients at intermediate levels
348!*************************************************
349
350  drhodz(:,:,1,n)=(rho(:,:,2,n)-rho(:,:,1,n))/ &
351       (height(2)-height(1))
352  do kz=2,nz-1
353    drhodz(:,:,kz,n)=(rho(:,:,kz+1,n)-rho(:,:,kz-1,n))/ &
354         (height(kz+1)-height(kz-1))
355  end do
356  drhodz(:,:,nz,n)=drhodz(:,:,nz-1,n)
357
358!    end do
359!  end do
360
361
362!****************************************************************
363! Compute slope of eta levels in windward direction and resulting
364! vertical wind correction
365!****************************************************************
366
367  do jy=1,ny-2
368    cosf(jy)=1./cos((real(jy)*dy+ylat0)*pi180)
369  enddo
370
371  kmin=2
372  idx=kmin
373  do iz=2,nz-1
374    do jy=1,ny-2
375      do ix=1,nx-2
376
377        inneta: do kz=idx(ix,jy),nz
378          if (idx(ix,jy) .le. kz .and. (height(iz).gt.uvzlev(ix,jy,kz-1)).and. &
379               (height(iz).le.uvzlev(ix,jy,kz))) then
380            idx(ix,jy)=kz
381            exit inneta
382          endif
383        enddo inneta
384      enddo
385    enddo
386
387    do jy=1,ny-2
388      do ix=1,nx-2
389        kz=idx(ix,jy)
390        dz1=height(iz)-uvzlev(ix,jy,kz-1)
391        dz2=uvzlev(ix,jy,kz)-height(iz)
392        dz=dz1+dz2
393        ix1=ix-1
394        jy1=jy-1
395        ixp=ix+1
396        jyp=jy+1
397
398        dzdx1=(uvzlev(ixp,jy,kz-1)-uvzlev(ix1,jy,kz-1))/2.
399        dzdx2=(uvzlev(ixp,jy,kz)-uvzlev(ix1,jy,kz))/2.
400        dzdx=(dzdx1*dz2+dzdx2*dz1)/dz
401
402        dzdy1=(uvzlev(ix,jyp,kz-1)-uvzlev(ix,jy1,kz-1))/2.
403        dzdy2=(uvzlev(ix,jyp,kz)-uvzlev(ix,jy1,kz))/2.
404        dzdy=(dzdy1*dz2+dzdy2*dz1)/dz
405
406        ww(ix,jy,iz,n)=ww(ix,jy,iz,n)+(dzdx*uu(ix,jy,iz,n)*dxconst*cosf(jy)+dzdy*vv(ix,jy,iz,n)*dyconst)
407
408      end do
409
410    end do
411  end do
412
413! If north pole is in the domain, calculate wind velocities in polar
414! stereographic coordinates
415!*******************************************************************
416
417  if (nglobal) then
418    do iz=1,nz
419      do jy=int(switchnorthg)-2,nymin1
420        ylat=ylat0+real(jy)*dy
421        do ix=0,nxmin1
422          xlon=xlon0+real(ix)*dx
423          call cc2gll(northpolemap,ylat,xlon,uu(ix,jy,iz,n), &
424               vv(ix,jy,iz,n),uupol(ix,jy,iz,n), &
425               vvpol(ix,jy,iz,n))
426        end do
427      end do
428    end do
429
430
431    do iz=1,nz
432
433! CALCULATE FFPOL, DDPOL FOR CENTRAL GRID POINT
434!
435!   AMSnauffer Nov 18 2004 Added check for case vv=0
436!
437      xlon=xlon0+real(nx/2-1)*dx
438      xlonr=xlon*pi/180.
439      ffpol=sqrt(uu(nx/2-1,nymin1,iz,n)**2+ &
440           vv(nx/2-1,nymin1,iz,n)**2)
441      if (vv(nx/2-1,nymin1,iz,n).lt.0.) then
442        ddpol=atan(uu(nx/2-1,nymin1,iz,n)/ &
443             vv(nx/2-1,nymin1,iz,n))-xlonr
444      else if (vv(nx/2-1,nymin1,iz,n).gt.0.) then
445        ddpol=pi+atan(uu(nx/2-1,nymin1,iz,n)/ &
446             vv(nx/2-1,nymin1,iz,n))-xlonr
447      else
448        ddpol=pi/2-xlonr
449      endif
450      if(ddpol.lt.0.) ddpol=2.0*pi+ddpol
451      if(ddpol.gt.2.0*pi) ddpol=ddpol-2.0*pi
452
453! CALCULATE U,V FOR 180 DEG, TRANSFORM TO POLAR STEREOGRAPHIC GRID
454      xlon=180.0
455      xlonr=xlon*pi/180.
456      ylat=90.0
457      uuaux=-ffpol*sin(xlonr+ddpol)
458      vvaux=-ffpol*cos(xlonr+ddpol)
459      call cc2gll(northpolemap,ylat,xlon,uuaux,vvaux,uupolaux, &
460           vvpolaux)
461
462      jy=nymin1
463      do ix=0,nxmin1
464        uupol(ix,jy,iz,n)=uupolaux
465        vvpol(ix,jy,iz,n)=vvpolaux
466      end do
467    end do
468
469
470! Fix: Set W at pole to the zonally averaged W of the next equator-
471! ward parallel of latitude
472
473    do iz=1,nz
474      wdummy=0.
475      jy=ny-2
476      do ix=0,nxmin1
477        wdummy=wdummy+ww(ix,jy,iz,n)
478      end do
479      wdummy=wdummy/real(nx)
480      jy=nymin1
481      do ix=0,nxmin1
482        ww(ix,jy,iz,n)=wdummy
483      end do
484    end do
485
486  endif
487
488
489! If south pole is in the domain, calculate wind velocities in polar
490! stereographic coordinates
491!*******************************************************************
492
493  if (sglobal) then
494    do iz=1,nz
495      do jy=0,int(switchsouthg)+3
496        ylat=ylat0+real(jy)*dy
497        do ix=0,nxmin1
498          xlon=xlon0+real(ix)*dx
499          call cc2gll(southpolemap,ylat,xlon,uu(ix,jy,iz,n), &
500               vv(ix,jy,iz,n),uupol(ix,jy,iz,n), &
501               vvpol(ix,jy,iz,n))
502        end do
503      end do
504    end do
505
506    do iz=1,nz
507
508! CALCULATE FFPOL, DDPOL FOR CENTRAL GRID POINT
509!
510!   AMSnauffer Nov 18 2004 Added check for case vv=0
511!
512      xlon=xlon0+real(nx/2-1)*dx
513      xlonr=xlon*pi/180.
514      ffpol=sqrt(uu(nx/2-1,0,iz,n)**2+ &
515           vv(nx/2-1,0,iz,n)**2)
516      if (vv(nx/2-1,0,iz,n).lt.0.) then
517        ddpol=atan(uu(nx/2-1,0,iz,n)/ &
518             vv(nx/2-1,0,iz,n))+xlonr
519      else if (vv(nx/2-1,0,iz,n).gt.0.) then
520        ddpol=pi+atan(uu(nx/2-1,0,iz,n)/ &
521             vv(nx/2-1,0,iz,n))+xlonr
522      else
523        ddpol=pi/2-xlonr
524      endif
525      if(ddpol.lt.0.) ddpol=2.0*pi+ddpol
526      if(ddpol.gt.2.0*pi) ddpol=ddpol-2.0*pi
527
528! CALCULATE U,V FOR 180 DEG, TRANSFORM TO POLAR STEREOGRAPHIC GRID
529      xlon=180.0
530      xlonr=xlon*pi/180.
531      ylat=-90.0
532      uuaux=+ffpol*sin(xlonr-ddpol)
533      vvaux=-ffpol*cos(xlonr-ddpol)
534      call cc2gll(northpolemap,ylat,xlon,uuaux,vvaux,uupolaux, &
535           vvpolaux)
536
537      jy=0
538      do ix=0,nxmin1
539        uupol(ix,jy,iz,n)=uupolaux
540        vvpol(ix,jy,iz,n)=vvpolaux
541      end do
542    end do
543
544
545! Fix: Set W at pole to the zonally averaged W of the next equator-
546! ward parallel of latitude
547
548    do iz=1,nz
549      wdummy=0.
550      jy=1
551      do ix=0,nxmin1
552        wdummy=wdummy+ww(ix,jy,iz,n)
553      end do
554      wdummy=wdummy/real(nx)
555      jy=0
556      do ix=0,nxmin1
557        ww(ix,jy,iz,n)=wdummy
558      end do
559    end do
560  endif
561
562
563  if (readclouds) write(*,*) 'using cloud water from ECMWF'
564!  if (.not.readclouds) write(*,*) 'using cloud water from parameterization'
565
566  rcw(:,:)=0
567  rpc(:,:)=0
568
569  do jy=0,nymin1
570    do ix=0,nxmin1
571! OLD METHOD
572      rain_cloud_above(ix,jy)=0
573      lsp=lsprec(ix,jy,1,n)
574      convp=convprec(ix,jy,1,n)
575      cloudsh(ix,jy,n)=0
576      do kz_inv=1,nz-1
577        kz=nz-kz_inv+1
578        pressure=rho(ix,jy,kz,n)*r_air*tt(ix,jy,kz,n)
579        rh=qv(ix,jy,kz,n)/f_qvsat(pressure,tt(ix,jy,kz,n))
580        clouds(ix,jy,kz,n)=0
581        if (rh.gt.0.8) then ! in cloud
582          if ((lsp.gt.0.01).or.(convp.gt.0.01)) then ! cloud and precipitation
583            rain_cloud_above(ix,jy)=1
584            cloudsh(ix,jy,n)=cloudsh(ix,jy,n)+ &
585                 height(kz)-height(kz-1)
586            if (lsp.ge.convp) then
587              clouds(ix,jy,kz,n)=3 ! lsp dominated rainout
588            else
589              clouds(ix,jy,kz,n)=2 ! convp dominated rainout
590            endif
591          else ! no precipitation
592            clouds(ix,jy,kz,n)=1 ! cloud
593          endif
594        else ! no cloud
595          if (rain_cloud_above(ix,jy).eq.1) then ! scavenging
596            if (lsp.ge.convp) then
597              clouds(ix,jy,kz,n)=5 ! lsp dominated washout
598            else
599              clouds(ix,jy,kz,n)=4 ! convp dominated washout
600            endif
601          endif
602        endif
603      end do
604!END OLD METHOD
605
606! PS 2012
607!      lsp=lsprec(ix,jy,1,n)
608!      convp=convprec(ix,jy,1,n)
609!          prec=lsp+convp
610!          if (lsp.gt.convp) then !  prectype='lsp'
611!            lconvectprec = .false.
612!          else ! prectype='cp '
613!            lconvectprec = .true.
614!           endif
615!HG METHOD
616!readclouds =.true.
617!      if (readclouds) then
618!hg added APR 2014  Cloud Water=clwc(ix,jy,kz,n)  Cloud Ice=ciwc(ix,jy,kz,n)
619!hg Use the cloud water variables to determine existence of clouds. This makes the PS code obsolete
620!        cloud_min=99999
621!        cloud_max=-1
622!        cloud_col_wat=0
623
624!        do kz=1, nz
625!         !clw & ciw are given in kg/kg 
626          ! cloud_water=clwc(ix,jy,kz,n)+ciwc(ix,jy,kz,n)
627          ! if (cloud_water .gt. 0) then
628          !   cloud_min=min(nint(height(kz)),cloud_min) !hg needs reset each grid
629          !   cloud_max=max(nint(height(kz)),cloud_max) !hg needs reset each grid
630          !   cloud_col_wat=cloud_col_wat+cloud_water !hg needs reset each grid
631          ! endif
632          ! cloud_ver=max(0,cloud_max-cloud_min)
633
634          ! icloudbot(ix,jy,n)=cloud_min
635          ! icloudthck(ix,jy,n)=cloud_ver
636          ! rcw(ix,jy)=cloud_col_wat
637          ! rpc(ix,jy)=prec
638!write(*,*) 'Using clouds from ECMWF' !hg END Henrik Code
639!END HG METHOD
640
641
642
643!      else ! windfields does not contain cloud data
644!          rhmin = 0.90 ! standard condition for presence of clouds
645!PS       note that original by Sabine Eckhart was 80%
646!PS       however, for T<-20 C we consider saturation over ice
647!PS       so I think 90% should be enough         
648!          icloudbot(ix,jy,n)=icmv
649!          icloudtop=icmv ! this is just a local variable
650!98        do kz=1,nz
651!            pressure=rho(ix,jy,kz,n)*r_air*tt(ix,jy,kz,n)
652!            rh=qv(ix,jy,kz,n)/f_qvsat(pressure,tt(ix,jy,kz,n))
653!ps            if (prec.gt.0.01) print*,'relhum',prec,kz,rh,height(kz)
654!            if (rh .gt. rhmin) then
655!              if (icloudbot(ix,jy,n) .eq. icmv) then
656!                icloudbot(ix,jy,n)=nint(height(kz))
657!              endif
658!              icloudtop=nint(height(kz)) ! use int to save memory
659!            endif
660    ! end do
661!PS try to get a cloud thicker than 50 m
662!PS if there is at least .01 mm/h  - changed to 0.002 and put into
663!PS parameter precpmin       
664!          if ((icloudbot(ix,jy,n) .eq. icmv .or. &
665!              icloudtop-icloudbot(ix,jy,n) .lt. 50) .and. &
666!              prec .gt. precmin) then
667!            rhmin = rhmin - 0.05
668!            if (rhmin .ge. 0.30) goto 98 ! give up for <= 25% rel.hum.
669!   end if
670
671!PS is based on looking at a limited set of comparison data
672!          if (lconvectprec .and. icloudtop .lt. 6000 .and. &
673!             prec .gt. precmin) then
674!
675!            if (convp .lt. 0.1) then
676!              icloudbot(ix,jy,n) = 500
677!              icloudtop =         8000
678!            else
679!              icloudbot(ix,jy,n) = 0
680!              icloudtop =      10000
681!            endif
682!          endif
683!          if (icloudtop .ne. icmv) then
684!            icloudthck(ix,jy,n) = icloudtop-icloudbot(ix,jy,n)
685!          else
686!            icloudthck(ix,jy,n) = icmv
687!          endif
688!PS  get rid of too thin clouds     
689!          if (icloudthck(ix,jy,n) .lt. 50) then
690!            icloudbot(ix,jy,n)=icmv
691!            icloudthck(ix,jy,n)=icmv
692!          endif
693!hg__________________________________
694!           rcw(ix,jy)=2E-7*prec**0.36
695!           rpc(ix,jy)=prec
696!hg end______________________________
697
698!      endif !hg read clouds
699
700    end do
701  end do
702
703!do 102 kz=1,nuvz
704!write(an,'(i02)') kz+10
705!write(*,*) nuvz,nymin1,nxmin1,'--',an,'--'
706!open(4,file='/nilu_wrk2/sec/cloudtest/cloud'//an,form='formatted')
707!do 101 jy=0,nymin1
708!    write(4,*) (clouds(ix,jy,kz,n),ix=1,nxmin1)
709!101   continue
710! close(4)
711!102   continue
712
713! open(4,file='/nilu_wrk2/sec/cloudtest/height',form='formatted')
714! do 103 jy=0,nymin1
715!     write (4,*)
716!+       (height(kz),kz=1,nuvz)
717!103   continue
718! close(4)
719
720!open(4,file='/nilu_wrk2/sec/cloudtest/p',form='formatted')
721! do 104 jy=0,nymin1
722!     write (4,*)
723!+       (r_air*tt(ix,jy,1,n)*rho(ix,jy,1,n),ix=1,nxmin1)
724!104   continue
725! close(4)
726
727!eso measure CPU time
728!  call mpif_mtime('verttransform',1)
729
730!eso print out the same measure as in Leo's routine
731    ! write(*,*) 'verttransform: ', &
732    !      sum(tt(:,:,:,n)*tt(:,:,:,n)), &
733    !      sum(uu(:,:,:,n)*uu(:,:,:,n)),sum(vv(:,:,:,n)*vv(:,:,:,n)),&
734    !      sum(qv(:,:,:,n)*qv(:,:,:,n)),sum(pv(:,:,:,n)*pv(:,:,:,n)),&
735    !      sum(rho(:,:,:,n)*rho(:,:,:,n)),sum(drhodz(:,:,:,n)*drhodz(:,:,:,n)),&
736    !      sum(ww(:,:,:,n)*ww(:,:,:,n)), &
737    !      sum(clouds(:,:,:,n)), sum(cloudsh(:,:,n)),sum(idx),sum(pinmconv)
738  end subroutine verttransform
739
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG