[e200b7a] | 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 | |
---|
| 22 | subroutine verttransform(n,uuh,vvh,wwh,pvh) |
---|
[4d45639] | 23 | ! i i i i i |
---|
[8a65cb0] | 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 | !***************************************************************************** |
---|
[e200b7a] | 65 | |
---|
| 66 | use par_mod |
---|
| 67 | use com_mod |
---|
| 68 | use cmapf_mod, only: cc2gll |
---|
[4d45639] | 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 |
---|
[e200b7a] | 74 | |
---|
| 75 | implicit none |
---|
| 76 | |
---|
| 77 | integer :: ix,jy,kz,iz,n,kmin,kl,klp,ix1,jy1,ixp,jyp,ixm,jym |
---|
[4d45639] | 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 |
---|
[e200b7a] | 85 | real :: xlon,ylat,xlonr,dzdx,dzdy |
---|
[4d45639] | 86 | real :: dzdx1,dzdx2,dzdy1,dzdy2 |
---|
[e200b7a] | 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) |
---|
[4d45639] | 92 | real :: wzlev(0:nxmax-1,0:nymax-1,nwzmax) |
---|
[e200b7a] | 93 | real,parameter :: const=r_air/ga |
---|
[4d45639] | 94 | ! integer:: ihr, imin, isec, i100th,ihr2, imin2, isec2, i100th2 |
---|
[8a65cb0] | 95 | parameter (precmin = 0.002) ! minimum prec in mm/h for cloud diagnostics |
---|
[e200b7a] | 96 | |
---|
| 97 | logical :: init = .true. |
---|
| 98 | |
---|
[d6a0977] | 99 | !ZHG SEP 2014 tests |
---|
[8a65cb0] | 100 | integer :: cloud_ver,cloud_min, cloud_max |
---|
[d6a0977] | 101 | integer ::teller(5), convpteller=0, lspteller=0 |
---|
[8a65cb0] | 102 | real :: cloud_col_wat, cloud_water |
---|
[d6a0977] | 103 | !ZHG 2015 temporary variables for testing |
---|
[8a65cb0] | 104 | real :: rcw(0:nxmax-1,0:nymax-1) |
---|
| 105 | real :: rpc(0:nxmax-1,0:nymax-1) |
---|
[d6a0977] | 106 | character(len=60) :: zhgpath='/xnilu_wrk/flex_wrk/zhg/' |
---|
| 107 | character(len=60) :: fnameA,fnameB,fnameC,fnameD,fnameE,fnameF,fnameG,fnameH |
---|
| 108 | CHARACTER(LEN=3) :: aspec |
---|
| 109 | integer :: virr=0 |
---|
| 110 | real :: tot_cloud_h |
---|
| 111 | !ZHG |
---|
[e200b7a] | 112 | |
---|
[8a65cb0] | 113 | !************************************************************************* |
---|
| 114 | ! If verttransform is called the first time, initialize heights of the * |
---|
| 115 | ! z levels in meter. The heights are the heights of model levels, where * |
---|
[4d45639] | 116 | ! u,v,T and qv are given, and of the interfaces, where w is given. So, * |
---|
| 117 | ! the vertical resolution in the z system is doubled. As reference point,* |
---|
| 118 | ! the lower left corner of the grid is used. * |
---|
| 119 | ! Unlike in the eta system, no difference between heights for u,v and * |
---|
| 120 | ! heights for w exists. * |
---|
[8a65cb0] | 121 | !************************************************************************* |
---|
[e200b7a] | 122 | |
---|
[4d45639] | 123 | |
---|
| 124 | !eso measure CPU time |
---|
| 125 | ! call mpif_mtime('verttransform',0) |
---|
| 126 | |
---|
[e200b7a] | 127 | if (init) then |
---|
| 128 | |
---|
[8a65cb0] | 129 | ! Search for a point with high surface pressure (i.e. not above significant topography) |
---|
| 130 | ! Then, use this point to construct a reference z profile, to be used at all times |
---|
[4d45639] | 131 | !***************************************************************************** |
---|
[e200b7a] | 132 | |
---|
| 133 | do jy=0,nymin1 |
---|
| 134 | do ix=0,nxmin1 |
---|
| 135 | if (ps(ix,jy,1,n).gt.100000.) then |
---|
| 136 | ixm=ix |
---|
| 137 | jym=jy |
---|
| 138 | goto 3 |
---|
| 139 | endif |
---|
| 140 | end do |
---|
| 141 | end do |
---|
| 142 | 3 continue |
---|
| 143 | |
---|
| 144 | |
---|
[4d45639] | 145 | tvold(ixm,jym)=tt2(ixm,jym,1,n)*(1.+0.378*ew(td2(ixm,jym,1,n))/ & |
---|
[8a65cb0] | 146 | ps(ixm,jym,1,n)) |
---|
[4d45639] | 147 | pold(ixm,jym)=ps(ixm,jym,1,n) |
---|
[e200b7a] | 148 | height(1)=0. |
---|
| 149 | |
---|
| 150 | do kz=2,nuvz |
---|
| 151 | pint=akz(kz)+bkz(kz)*ps(ixm,jym,1,n) |
---|
| 152 | tv=tth(ixm,jym,kz,n)*(1.+0.608*qvh(ixm,jym,kz,n)) |
---|
| 153 | |
---|
[4d45639] | 154 | if (abs(tv(ixm,jym)-tvold(ixm,jym)).gt.0.2) then |
---|
| 155 | height(kz)= & |
---|
| 156 | height(kz-1)+const*log(pold(ixm,jym)/pint(ixm,jym))* & |
---|
| 157 | (tv(ixm,jym)-tvold(ixm,jym))/log(tv(ixm,jym)/tvold(ixm,jym)) |
---|
[e200b7a] | 158 | else |
---|
[4d45639] | 159 | height(kz)=height(kz-1)+ & |
---|
| 160 | const*log(pold(ixm,jym)/pint(ixm,jym))*tv(ixm,jym) |
---|
[e200b7a] | 161 | endif |
---|
| 162 | |
---|
[4d45639] | 163 | tvold(ixm,jym)=tv(ixm,jym) |
---|
| 164 | pold(ixm,jym)=pint(ixm,jym) |
---|
[e200b7a] | 165 | end do |
---|
| 166 | |
---|
| 167 | |
---|
[8a65cb0] | 168 | ! Determine highest levels that can be within PBL |
---|
| 169 | !************************************************ |
---|
[e200b7a] | 170 | |
---|
| 171 | do kz=1,nz |
---|
| 172 | if (height(kz).gt.hmixmax) then |
---|
| 173 | nmixz=kz |
---|
| 174 | goto 9 |
---|
| 175 | endif |
---|
| 176 | end do |
---|
| 177 | 9 continue |
---|
| 178 | |
---|
[8a65cb0] | 179 | ! Do not repeat initialization of the Cartesian z grid |
---|
| 180 | !***************************************************** |
---|
[e200b7a] | 181 | |
---|
| 182 | init=.false. |
---|
| 183 | |
---|
| 184 | endif |
---|
| 185 | |
---|
| 186 | |
---|
[8a65cb0] | 187 | ! Loop over the whole grid |
---|
| 188 | !************************* |
---|
[e200b7a] | 189 | |
---|
[4d45639] | 190 | |
---|
[e200b7a] | 191 | do jy=0,nymin1 |
---|
| 192 | do ix=0,nxmin1 |
---|
[4d45639] | 193 | tvold(ix,jy)=tt2(ix,jy,1,n)*(1.+0.378*ew(td2(ix,jy,1,n))/ & |
---|
| 194 | ps(ix,jy,1,n)) |
---|
| 195 | enddo |
---|
| 196 | enddo |
---|
| 197 | pold=ps(:,:,1,n) |
---|
| 198 | uvzlev(:,:,1)=0. |
---|
| 199 | wzlev(:,:,1)=0. |
---|
| 200 | rhoh(:,:,1)=pold/(r_air*tvold) |
---|
[e200b7a] | 201 | |
---|
| 202 | |
---|
[8a65cb0] | 203 | ! Compute heights of eta levels |
---|
| 204 | !****************************** |
---|
[e200b7a] | 205 | |
---|
[4d45639] | 206 | do kz=2,nuvz |
---|
| 207 | pint=akz(kz)+bkz(kz)*ps(:,:,1,n) |
---|
| 208 | tv=tth(:,:,kz,n)*(1.+0.608*qvh(:,:,kz,n)) |
---|
| 209 | rhoh(:,:,kz)=pint(:,:)/(r_air*tv) |
---|
[e200b7a] | 210 | |
---|
[4d45639] | 211 | where (abs(tv-tvold).gt.0.2) |
---|
| 212 | uvzlev(:,:,kz)=uvzlev(:,:,kz-1)+const*log(pold/pint)* & |
---|
| 213 | (tv-tvold)/log(tv/tvold) |
---|
| 214 | elsewhere |
---|
| 215 | uvzlev(:,:,kz)=uvzlev(:,:,kz-1)+const*log(pold/pint)*tv |
---|
| 216 | endwhere |
---|
[e200b7a] | 217 | |
---|
[4d45639] | 218 | tvold=tv |
---|
| 219 | pold=pint |
---|
| 220 | end do |
---|
[e200b7a] | 221 | |
---|
| 222 | |
---|
[4d45639] | 223 | do kz=2,nwz-1 |
---|
| 224 | wzlev(:,:,kz)=(uvzlev(:,:,kz+1)+uvzlev(:,:,kz))/2. |
---|
| 225 | end do |
---|
| 226 | wzlev(:,:,nwz)=wzlev(:,:,nwz-1)+ & |
---|
| 227 | uvzlev(:,:,nuvz)-uvzlev(:,:,nuvz-1) |
---|
[e200b7a] | 228 | |
---|
[8a65cb0] | 229 | ! pinmconv=(h2-h1)/(p2-p1) |
---|
[e200b7a] | 230 | |
---|
[4d45639] | 231 | pinmconv(:,:,1)=(uvzlev(:,:,2))/ & |
---|
| 232 | ((aknew(2)+bknew(2)*ps(:,:,1,n))- & |
---|
| 233 | (aknew(1)+bknew(1)*ps(:,:,1,n))) |
---|
| 234 | do kz=2,nz-1 |
---|
| 235 | pinmconv(:,:,kz)=(uvzlev(:,:,kz+1)-uvzlev(:,:,kz-1))/ & |
---|
| 236 | ((aknew(kz+1)+bknew(kz+1)*ps(:,:,1,n))- & |
---|
| 237 | (aknew(kz-1)+bknew(kz-1)*ps(:,:,1,n))) |
---|
| 238 | end do |
---|
| 239 | pinmconv(:,:,nz)=(uvzlev(:,:,nz)-uvzlev(:,:,nz-1))/ & |
---|
| 240 | ((aknew(nz)+bknew(nz)*ps(:,:,1,n))- & |
---|
| 241 | (aknew(nz-1)+bknew(nz-1)*ps(:,:,1,n))) |
---|
[e200b7a] | 242 | |
---|
[8a65cb0] | 243 | ! Levels, where u,v,t and q are given |
---|
| 244 | !************************************ |
---|
[e200b7a] | 245 | |
---|
[4d45639] | 246 | |
---|
| 247 | uu(:,:,1,n)=uuh(:,:,1) |
---|
| 248 | vv(:,:,1,n)=vvh(:,:,1) |
---|
| 249 | tt(:,:,1,n)=tth(:,:,1,n) |
---|
| 250 | qv(:,:,1,n)=qvh(:,:,1,n) |
---|
[8a65cb0] | 251 | !hg adding the cloud water |
---|
[4d45639] | 252 | clwc(:,:,1,n)=clwch(:,:,1,n) |
---|
| 253 | ciwc(:,:,1,n)=ciwch(:,:,1,n) |
---|
[8a65cb0] | 254 | !hg |
---|
[4d45639] | 255 | pv(:,:,1,n)=pvh(:,:,1) |
---|
| 256 | rho(:,:,1,n)=rhoh(:,:,1) |
---|
| 257 | uu(:,:,nz,n)=uuh(:,:,nuvz) |
---|
| 258 | vv(:,:,nz,n)=vvh(:,:,nuvz) |
---|
| 259 | tt(:,:,nz,n)=tth(:,:,nuvz,n) |
---|
| 260 | qv(:,:,nz,n)=qvh(:,:,nuvz,n) |
---|
[8a65cb0] | 261 | |
---|
| 262 | !hg adding the cloud water |
---|
[4d45639] | 263 | clwc(:,:,nz,n)=clwch(:,:,nuvz,n) |
---|
| 264 | ciwc(:,:,nz,n)=ciwch(:,:,nuvz,n) |
---|
[8a65cb0] | 265 | !hg |
---|
[4d45639] | 266 | pv(:,:,nz,n)=pvh(:,:,nuvz) |
---|
| 267 | rho(:,:,nz,n)=rhoh(:,:,nuvz) |
---|
| 268 | |
---|
| 269 | |
---|
| 270 | kmin=2 |
---|
| 271 | idx=kmin |
---|
| 272 | do iz=2,nz-1 |
---|
| 273 | do jy=0,nymin1 |
---|
| 274 | do ix=0,nxmin1 |
---|
| 275 | if(height(iz).gt.uvzlev(ix,jy,nuvz)) then |
---|
| 276 | uu(ix,jy,iz,n)=uu(ix,jy,nz,n) |
---|
| 277 | vv(ix,jy,iz,n)=vv(ix,jy,nz,n) |
---|
| 278 | tt(ix,jy,iz,n)=tt(ix,jy,nz,n) |
---|
| 279 | qv(ix,jy,iz,n)=qv(ix,jy,nz,n) |
---|
[8a65cb0] | 280 | !hg adding the cloud water |
---|
[4d45639] | 281 | clwc(ix,jy,iz,n)=clwc(ix,jy,nz,n) |
---|
| 282 | ciwc(ix,jy,iz,n)=ciwc(ix,jy,nz,n) |
---|
[8a65cb0] | 283 | !hg |
---|
[4d45639] | 284 | pv(ix,jy,iz,n)=pv(ix,jy,nz,n) |
---|
| 285 | rho(ix,jy,iz,n)=rho(ix,jy,nz,n) |
---|
| 286 | else |
---|
| 287 | innuvz: do kz=idx(ix,jy),nuvz |
---|
| 288 | if (idx(ix,jy) .le. kz .and. (height(iz).gt.uvzlev(ix,jy,kz-1)).and. & |
---|
| 289 | (height(iz).le.uvzlev(ix,jy,kz))) then |
---|
| 290 | idx(ix,jy)=kz |
---|
| 291 | exit innuvz |
---|
| 292 | endif |
---|
| 293 | enddo innuvz |
---|
| 294 | endif |
---|
| 295 | enddo |
---|
| 296 | enddo |
---|
| 297 | do jy=0,nymin1 |
---|
| 298 | do ix=0,nxmin1 |
---|
| 299 | if(height(iz).le.uvzlev(ix,jy,nuvz)) then |
---|
| 300 | kz=idx(ix,jy) |
---|
| 301 | dz1=height(iz)-uvzlev(ix,jy,kz-1) |
---|
| 302 | dz2=uvzlev(ix,jy,kz)-height(iz) |
---|
| 303 | dz=dz1+dz2 |
---|
| 304 | uu(ix,jy,iz,n)=(uuh(ix,jy,kz-1)*dz2+uuh(ix,jy,kz)*dz1)/dz |
---|
| 305 | vv(ix,jy,iz,n)=(vvh(ix,jy,kz-1)*dz2+vvh(ix,jy,kz)*dz1)/dz |
---|
| 306 | tt(ix,jy,iz,n)=(tth(ix,jy,kz-1,n)*dz2 & |
---|
| 307 | +tth(ix,jy,kz,n)*dz1)/dz |
---|
| 308 | qv(ix,jy,iz,n)=(qvh(ix,jy,kz-1,n)*dz2 & |
---|
| 309 | +qvh(ix,jy,kz,n)*dz1)/dz |
---|
[8a65cb0] | 310 | !hg adding the cloud water |
---|
[4d45639] | 311 | clwc(ix,jy,iz,n)=(clwch(ix,jy,kz-1,n)*dz2+clwch(ix,jy,kz,n)*dz1)/dz |
---|
| 312 | ciwc(ix,jy,iz,n)=(ciwch(ix,jy,kz-1,n)*dz2+ciwch(ix,jy,kz,n)*dz1)/dz |
---|
[8a65cb0] | 313 | !hg |
---|
[4d45639] | 314 | pv(ix,jy,iz,n)=(pvh(ix,jy,kz-1)*dz2+pvh(ix,jy,kz)*dz1)/dz |
---|
| 315 | rho(ix,jy,iz,n)=(rhoh(ix,jy,kz-1)*dz2+rhoh(ix,jy,kz)*dz1)/dz |
---|
| 316 | endif |
---|
| 317 | enddo |
---|
| 318 | enddo |
---|
| 319 | enddo |
---|
[e200b7a] | 320 | |
---|
| 321 | |
---|
[8a65cb0] | 322 | ! Levels, where w is given |
---|
| 323 | !************************* |
---|
[e200b7a] | 324 | |
---|
[4d45639] | 325 | ww(:,:,1,n)=wwh(:,:,1)*pinmconv(:,:,1) |
---|
| 326 | ww(:,:,nz,n)=wwh(:,:,nwz)*pinmconv(:,:,nz) |
---|
| 327 | kmin=2 |
---|
| 328 | idx=kmin |
---|
| 329 | do iz=2,nz |
---|
| 330 | do jy=0,nymin1 |
---|
| 331 | do ix=0,nxmin1 |
---|
| 332 | inn: do kz=idx(ix,jy),nwz |
---|
| 333 | if(idx(ix,jy) .le. kz .and. height(iz).gt.wzlev(ix,jy,kz-1).and. & |
---|
| 334 | height(iz).le.wzlev(ix,jy,kz)) then |
---|
| 335 | idx(ix,jy)=kz |
---|
| 336 | exit inn |
---|
[e200b7a] | 337 | endif |
---|
[4d45639] | 338 | enddo inn |
---|
| 339 | enddo |
---|
| 340 | enddo |
---|
| 341 | do jy=0,nymin1 |
---|
| 342 | do ix=0,nxmin1 |
---|
| 343 | kz=idx(ix,jy) |
---|
| 344 | dz1=height(iz)-wzlev(ix,jy,kz-1) |
---|
| 345 | dz2=wzlev(ix,jy,kz)-height(iz) |
---|
| 346 | dz=dz1+dz2 |
---|
| 347 | ww(ix,jy,iz,n)=(wwh(ix,jy,kz-1)*pinmconv(ix,jy,kz-1)*dz2 & |
---|
| 348 | +wwh(ix,jy,kz)*pinmconv(ix,jy,kz)*dz1)/dz |
---|
| 349 | enddo |
---|
| 350 | enddo |
---|
| 351 | enddo |
---|
[e200b7a] | 352 | |
---|
[8a65cb0] | 353 | ! Compute density gradients at intermediate levels |
---|
| 354 | !************************************************* |
---|
[e200b7a] | 355 | |
---|
[4d45639] | 356 | drhodz(:,:,1,n)=(rho(:,:,2,n)-rho(:,:,1,n))/ & |
---|
| 357 | (height(2)-height(1)) |
---|
| 358 | do kz=2,nz-1 |
---|
| 359 | drhodz(:,:,kz,n)=(rho(:,:,kz+1,n)-rho(:,:,kz-1,n))/ & |
---|
| 360 | (height(kz+1)-height(kz-1)) |
---|
[e200b7a] | 361 | end do |
---|
[4d45639] | 362 | drhodz(:,:,nz,n)=drhodz(:,:,nz-1,n) |
---|
| 363 | |
---|
| 364 | ! end do |
---|
| 365 | ! end do |
---|
[e200b7a] | 366 | |
---|
| 367 | |
---|
[8a65cb0] | 368 | !**************************************************************** |
---|
| 369 | ! Compute slope of eta levels in windward direction and resulting |
---|
| 370 | ! vertical wind correction |
---|
| 371 | !**************************************************************** |
---|
[e200b7a] | 372 | |
---|
| 373 | do jy=1,ny-2 |
---|
[4d45639] | 374 | cosf(jy)=1./cos((real(jy)*dy+ylat0)*pi180) |
---|
| 375 | enddo |
---|
| 376 | |
---|
| 377 | kmin=2 |
---|
| 378 | idx=kmin |
---|
| 379 | do iz=2,nz-1 |
---|
| 380 | do jy=1,ny-2 |
---|
| 381 | do ix=1,nx-2 |
---|
| 382 | |
---|
| 383 | inneta: do kz=idx(ix,jy),nz |
---|
| 384 | if (idx(ix,jy) .le. kz .and. (height(iz).gt.uvzlev(ix,jy,kz-1)).and. & |
---|
| 385 | (height(iz).le.uvzlev(ix,jy,kz))) then |
---|
| 386 | idx(ix,jy)=kz |
---|
| 387 | exit inneta |
---|
[e200b7a] | 388 | endif |
---|
[4d45639] | 389 | enddo inneta |
---|
| 390 | enddo |
---|
| 391 | enddo |
---|
| 392 | |
---|
| 393 | do jy=1,ny-2 |
---|
| 394 | do ix=1,nx-2 |
---|
| 395 | kz=idx(ix,jy) |
---|
| 396 | dz1=height(iz)-uvzlev(ix,jy,kz-1) |
---|
| 397 | dz2=uvzlev(ix,jy,kz)-height(iz) |
---|
| 398 | dz=dz1+dz2 |
---|
| 399 | ix1=ix-1 |
---|
[e200b7a] | 400 | jy1=jy-1 |
---|
| 401 | ixp=ix+1 |
---|
| 402 | jyp=jy+1 |
---|
| 403 | |
---|
[4d45639] | 404 | dzdx1=(uvzlev(ixp,jy,kz-1)-uvzlev(ix1,jy,kz-1))/2. |
---|
| 405 | dzdx2=(uvzlev(ixp,jy,kz)-uvzlev(ix1,jy,kz))/2. |
---|
[e200b7a] | 406 | dzdx=(dzdx1*dz2+dzdx2*dz1)/dz |
---|
| 407 | |
---|
[4d45639] | 408 | dzdy1=(uvzlev(ix,jyp,kz-1)-uvzlev(ix,jy1,kz-1))/2. |
---|
| 409 | dzdy2=(uvzlev(ix,jyp,kz)-uvzlev(ix,jy1,kz))/2. |
---|
[e200b7a] | 410 | dzdy=(dzdy1*dz2+dzdy2*dz1)/dz |
---|
| 411 | |
---|
[4d45639] | 412 | 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) |
---|
[e200b7a] | 413 | |
---|
| 414 | end do |
---|
| 415 | |
---|
| 416 | end do |
---|
| 417 | end do |
---|
| 418 | |
---|
[8a65cb0] | 419 | ! If north pole is in the domain, calculate wind velocities in polar |
---|
| 420 | ! stereographic coordinates |
---|
| 421 | !******************************************************************* |
---|
[e200b7a] | 422 | |
---|
| 423 | if (nglobal) then |
---|
[4d45639] | 424 | do iz=1,nz |
---|
| 425 | do jy=int(switchnorthg)-2,nymin1 |
---|
| 426 | ylat=ylat0+real(jy)*dy |
---|
| 427 | do ix=0,nxmin1 |
---|
| 428 | xlon=xlon0+real(ix)*dx |
---|
[e200b7a] | 429 | call cc2gll(northpolemap,ylat,xlon,uu(ix,jy,iz,n), & |
---|
[4d45639] | 430 | vv(ix,jy,iz,n),uupol(ix,jy,iz,n), & |
---|
| 431 | vvpol(ix,jy,iz,n)) |
---|
[e200b7a] | 432 | end do |
---|
| 433 | end do |
---|
| 434 | end do |
---|
| 435 | |
---|
| 436 | |
---|
| 437 | do iz=1,nz |
---|
| 438 | |
---|
[8a65cb0] | 439 | ! CALCULATE FFPOL, DDPOL FOR CENTRAL GRID POINT |
---|
| 440 | ! |
---|
| 441 | ! AMSnauffer Nov 18 2004 Added check for case vv=0 |
---|
| 442 | ! |
---|
[e200b7a] | 443 | xlon=xlon0+real(nx/2-1)*dx |
---|
| 444 | xlonr=xlon*pi/180. |
---|
[4d45639] | 445 | ffpol=sqrt(uu(nx/2-1,nymin1,iz,n)**2+ & |
---|
| 446 | vv(nx/2-1,nymin1,iz,n)**2) |
---|
[e200b7a] | 447 | if (vv(nx/2-1,nymin1,iz,n).lt.0.) then |
---|
[4d45639] | 448 | ddpol=atan(uu(nx/2-1,nymin1,iz,n)/ & |
---|
| 449 | vv(nx/2-1,nymin1,iz,n))-xlonr |
---|
[e200b7a] | 450 | else if (vv(nx/2-1,nymin1,iz,n).gt.0.) then |
---|
| 451 | ddpol=pi+atan(uu(nx/2-1,nymin1,iz,n)/ & |
---|
[8a65cb0] | 452 | vv(nx/2-1,nymin1,iz,n))-xlonr |
---|
[e200b7a] | 453 | else |
---|
| 454 | ddpol=pi/2-xlonr |
---|
| 455 | endif |
---|
| 456 | if(ddpol.lt.0.) ddpol=2.0*pi+ddpol |
---|
| 457 | if(ddpol.gt.2.0*pi) ddpol=ddpol-2.0*pi |
---|
| 458 | |
---|
[8a65cb0] | 459 | ! CALCULATE U,V FOR 180 DEG, TRANSFORM TO POLAR STEREOGRAPHIC GRID |
---|
[e200b7a] | 460 | xlon=180.0 |
---|
| 461 | xlonr=xlon*pi/180. |
---|
| 462 | ylat=90.0 |
---|
| 463 | uuaux=-ffpol*sin(xlonr+ddpol) |
---|
| 464 | vvaux=-ffpol*cos(xlonr+ddpol) |
---|
[4d45639] | 465 | call cc2gll(northpolemap,ylat,xlon,uuaux,vvaux,uupolaux, & |
---|
| 466 | vvpolaux) |
---|
| 467 | |
---|
[e200b7a] | 468 | jy=nymin1 |
---|
| 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 | |
---|
[8a65cb0] | 476 | ! Fix: Set W at pole to the zonally averaged W of the next equator- |
---|
| 477 | ! ward parallel of latitude |
---|
[e200b7a] | 478 | |
---|
[8a65cb0] | 479 | do iz=1,nz |
---|
[e200b7a] | 480 | wdummy=0. |
---|
| 481 | jy=ny-2 |
---|
| 482 | do ix=0,nxmin1 |
---|
| 483 | wdummy=wdummy+ww(ix,jy,iz,n) |
---|
| 484 | end do |
---|
| 485 | wdummy=wdummy/real(nx) |
---|
| 486 | jy=nymin1 |
---|
| 487 | do ix=0,nxmin1 |
---|
| 488 | ww(ix,jy,iz,n)=wdummy |
---|
| 489 | end do |
---|
[8a65cb0] | 490 | end do |
---|
[e200b7a] | 491 | |
---|
| 492 | endif |
---|
| 493 | |
---|
| 494 | |
---|
[8a65cb0] | 495 | ! If south pole is in the domain, calculate wind velocities in polar |
---|
| 496 | ! stereographic coordinates |
---|
| 497 | !******************************************************************* |
---|
[e200b7a] | 498 | |
---|
| 499 | if (sglobal) then |
---|
[4d45639] | 500 | do iz=1,nz |
---|
| 501 | do jy=0,int(switchsouthg)+3 |
---|
| 502 | ylat=ylat0+real(jy)*dy |
---|
| 503 | do ix=0,nxmin1 |
---|
| 504 | xlon=xlon0+real(ix)*dx |
---|
[e200b7a] | 505 | call cc2gll(southpolemap,ylat,xlon,uu(ix,jy,iz,n), & |
---|
[4d45639] | 506 | vv(ix,jy,iz,n),uupol(ix,jy,iz,n), & |
---|
| 507 | vvpol(ix,jy,iz,n)) |
---|
[e200b7a] | 508 | end do |
---|
| 509 | end do |
---|
| 510 | end do |
---|
| 511 | |
---|
| 512 | do iz=1,nz |
---|
| 513 | |
---|
[8a65cb0] | 514 | ! CALCULATE FFPOL, DDPOL FOR CENTRAL GRID POINT |
---|
| 515 | ! |
---|
| 516 | ! AMSnauffer Nov 18 2004 Added check for case vv=0 |
---|
| 517 | ! |
---|
[e200b7a] | 518 | xlon=xlon0+real(nx/2-1)*dx |
---|
| 519 | xlonr=xlon*pi/180. |
---|
[4d45639] | 520 | ffpol=sqrt(uu(nx/2-1,0,iz,n)**2+ & |
---|
| 521 | vv(nx/2-1,0,iz,n)**2) |
---|
[e200b7a] | 522 | if (vv(nx/2-1,0,iz,n).lt.0.) then |
---|
[4d45639] | 523 | ddpol=atan(uu(nx/2-1,0,iz,n)/ & |
---|
| 524 | vv(nx/2-1,0,iz,n))+xlonr |
---|
[e200b7a] | 525 | else if (vv(nx/2-1,0,iz,n).gt.0.) then |
---|
[4d45639] | 526 | ddpol=pi+atan(uu(nx/2-1,0,iz,n)/ & |
---|
| 527 | vv(nx/2-1,0,iz,n))+xlonr |
---|
[e200b7a] | 528 | else |
---|
| 529 | ddpol=pi/2-xlonr |
---|
| 530 | endif |
---|
| 531 | if(ddpol.lt.0.) ddpol=2.0*pi+ddpol |
---|
| 532 | if(ddpol.gt.2.0*pi) ddpol=ddpol-2.0*pi |
---|
| 533 | |
---|
[8a65cb0] | 534 | ! CALCULATE U,V FOR 180 DEG, TRANSFORM TO POLAR STEREOGRAPHIC GRID |
---|
[e200b7a] | 535 | xlon=180.0 |
---|
| 536 | xlonr=xlon*pi/180. |
---|
| 537 | ylat=-90.0 |
---|
| 538 | uuaux=+ffpol*sin(xlonr-ddpol) |
---|
| 539 | vvaux=-ffpol*cos(xlonr-ddpol) |
---|
[4d45639] | 540 | call cc2gll(northpolemap,ylat,xlon,uuaux,vvaux,uupolaux, & |
---|
| 541 | vvpolaux) |
---|
[e200b7a] | 542 | |
---|
| 543 | jy=0 |
---|
| 544 | do ix=0,nxmin1 |
---|
| 545 | uupol(ix,jy,iz,n)=uupolaux |
---|
| 546 | vvpol(ix,jy,iz,n)=vvpolaux |
---|
| 547 | end do |
---|
| 548 | end do |
---|
| 549 | |
---|
| 550 | |
---|
[8a65cb0] | 551 | ! Fix: Set W at pole to the zonally averaged W of the next equator- |
---|
| 552 | ! ward parallel of latitude |
---|
[e200b7a] | 553 | |
---|
| 554 | do iz=1,nz |
---|
| 555 | wdummy=0. |
---|
| 556 | jy=1 |
---|
| 557 | do ix=0,nxmin1 |
---|
| 558 | wdummy=wdummy+ww(ix,jy,iz,n) |
---|
| 559 | end do |
---|
| 560 | wdummy=wdummy/real(nx) |
---|
| 561 | jy=0 |
---|
| 562 | do ix=0,nxmin1 |
---|
| 563 | ww(ix,jy,iz,n)=wdummy |
---|
| 564 | end do |
---|
| 565 | end do |
---|
| 566 | endif |
---|
| 567 | |
---|
| 568 | |
---|
[d6a0977] | 569 | !*********************************************************************************** |
---|
[f75967d] | 570 | if (readclouds) then !HG METHOD |
---|
[d6a0977] | 571 | ! The method is loops all grids vertically and constructs the 3D matrix for clouds |
---|
| 572 | ! Cloud top and cloud bottom gid cells are assigned as well as the total column |
---|
| 573 | ! cloud water. For precipitating grids, the type and whether it is in or below |
---|
| 574 | ! cloud scavenging are assigned with numbers 2-5 (following the old metod). |
---|
| 575 | ! Distinction is done for lsp and convp though they are treated the same in regards |
---|
| 576 | ! to scavenging. Also clouds that are not precipitating are defined which may be |
---|
| 577 | ! to include future cloud processing by non-precipitating-clouds. |
---|
| 578 | !*********************************************************************************** |
---|
[f75967d] | 579 | write(*,*) 'using cloud water from ECMWF' |
---|
| 580 | clw(:,:,:,n)=0 |
---|
| 581 | icloud_stats(:,:,:,n)=0 |
---|
| 582 | clouds(:,:,:,n)=0 |
---|
| 583 | do jy=0,nymin1 |
---|
| 584 | do ix=0,nxmin1 |
---|
| 585 | lsp=lsprec(ix,jy,1,n) |
---|
| 586 | convp=convprec(ix,jy,1,n) |
---|
| 587 | prec=lsp+convp |
---|
| 588 | tot_cloud_h=0 |
---|
| 589 | ! Find clouds in the vertical |
---|
| 590 | do kz=1, nz-1 !go from top to bottom |
---|
| 591 | if (clwc(ix,jy,kz,n).gt.0) then |
---|
| 592 | ! assuming rho is in kg/m3 and hz in m gives: kg/kg * kg/m3 *m3/kg /m = m2/m3 |
---|
| 593 | clw(ix,jy,kz,n)=(clwc(ix,jy,kz,n)*rho(ix,jy,kz,n))*(height(kz+1)-height(kz)) |
---|
| 594 | tot_cloud_h=tot_cloud_h+(height(kz+1)-height(kz)) |
---|
| 595 | icloud_stats(ix,jy,4,n)= icloud_stats(ix,jy,4,n)+clw(ix,jy,kz,n) ! Column cloud water [m3/m3] |
---|
| 596 | icloud_stats(ix,jy,3,n)= min(height(kz+1),height(kz)) ! Cloud BOT height stats [m] |
---|
[d6a0977] | 597 | !ZHG 2015 extra for testing |
---|
[f75967d] | 598 | clh(ix,jy,kz,n)=height(kz+1)-height(kz) |
---|
| 599 | icloud_stats(ix,jy,1,n)=icloud_stats(ix,jy,1,n)+(height(kz+1)-height(kz)) ! Cloud total vertical extent [m] |
---|
| 600 | icloud_stats(ix,jy,2,n)= max(icloud_stats(ix,jy,2,n),height(kz)) ! Cloud TOP height [m] |
---|
[d6a0977] | 601 | !ZHG |
---|
[f75967d] | 602 | endif |
---|
| 603 | end do |
---|
[d6a0977] | 604 | |
---|
[f75967d] | 605 | ! If Precipitation. Define removal type in the vertical |
---|
| 606 | if ((lsp.gt.0.01).or.(convp.gt.0.01)) then ! cloud and precipitation |
---|
| 607 | |
---|
| 608 | do kz=nz,1,-1 !go Bottom up! |
---|
| 609 | if (clw(ix,jy,kz,n).gt. 0) then ! is in cloud |
---|
| 610 | cloudsh(ix,jy,n)=cloudsh(ix,jy,n)+height(kz)-height(kz-1) |
---|
| 611 | clouds(ix,jy,kz,n)=1 ! is a cloud |
---|
| 612 | if (lsp.ge.convp) then |
---|
| 613 | clouds(ix,jy,kz,n)=3 ! lsp in-cloud |
---|
| 614 | else |
---|
| 615 | clouds(ix,jy,kz,n)=2 ! convp in-cloud |
---|
| 616 | endif ! convective or large scale |
---|
| 617 | elseif((clw(ix,jy,kz,n).le.0) .and. (icloud_stats(ix,jy,3,n).ge.height(kz)) ) then ! is below cloud |
---|
| 618 | if (lsp.ge.convp) then |
---|
| 619 | clouds(ix,jy,kz,n)=5 ! lsp dominated washout |
---|
| 620 | else |
---|
| 621 | clouds(ix,jy,kz,n)=4 ! convp dominated washout |
---|
| 622 | endif ! convective or large scale |
---|
| 623 | endif |
---|
[8a65cb0] | 624 | |
---|
[f75967d] | 625 | if (height(kz).ge. 19000) then ! set a max height for removal |
---|
| 626 | clouds(ix,jy,kz,n)=0 |
---|
| 627 | endif !clw>0 |
---|
| 628 | end do !nz |
---|
| 629 | endif ! precipitation |
---|
| 630 | end do |
---|
| 631 | end do |
---|
[d6a0977] | 632 | !************************************************************************** |
---|
[f75967d] | 633 | else ! use old definitions |
---|
[d6a0977] | 634 | !************************************************************************** |
---|
| 635 | ! create a cloud and rainout/washout field, clouds occur where rh>80% |
---|
| 636 | ! total cloudheight is stored at level 0 |
---|
[f75967d] | 637 | ! if (.not.readclouds) write(*,*) 'using cloud water from Parameterization' |
---|
| 638 | write(*,*) 'using cloud water from Parameterization' |
---|
| 639 | do jy=0,nymin1 |
---|
| 640 | do ix=0,nxmin1 |
---|
[8a65cb0] | 641 | ! OLD METHOD |
---|
[f75967d] | 642 | rain_cloud_above(ix,jy)=0 |
---|
| 643 | lsp=lsprec(ix,jy,1,n) |
---|
| 644 | convp=convprec(ix,jy,1,n) |
---|
| 645 | cloudsh(ix,jy,n)=0 |
---|
| 646 | do kz_inv=1,nz-1 |
---|
| 647 | kz=nz-kz_inv+1 |
---|
| 648 | pressure=rho(ix,jy,kz,n)*r_air*tt(ix,jy,kz,n) |
---|
| 649 | rh=qv(ix,jy,kz,n)/f_qvsat(pressure,tt(ix,jy,kz,n)) |
---|
| 650 | clouds(ix,jy,kz,n)=0 |
---|
| 651 | if (rh.gt.0.8) then ! in cloud |
---|
| 652 | if ((lsp.gt.0.01).or.(convp.gt.0.01)) then ! cloud and precipitation |
---|
| 653 | rain_cloud_above(ix,jy)=1 |
---|
| 654 | cloudsh(ix,jy,n)=cloudsh(ix,jy,n)+ & |
---|
| 655 | height(kz)-height(kz-1) |
---|
| 656 | if (lsp.ge.convp) then |
---|
| 657 | clouds(ix,jy,kz,n)=3 ! lsp dominated rainout |
---|
| 658 | else |
---|
| 659 | clouds(ix,jy,kz,n)=2 ! convp dominated rainout |
---|
| 660 | endif |
---|
| 661 | else ! no precipitation |
---|
| 662 | clouds(ix,jy,kz,n)=1 ! cloud |
---|
[e200b7a] | 663 | endif |
---|
[f75967d] | 664 | else ! no cloud |
---|
| 665 | if (rain_cloud_above(ix,jy).eq.1) then ! scavenging |
---|
| 666 | if (lsp.ge.convp) then |
---|
| 667 | clouds(ix,jy,kz,n)=5 ! lsp dominated washout |
---|
| 668 | else |
---|
| 669 | clouds(ix,jy,kz,n)=4 ! convp dominated washout |
---|
| 670 | endif |
---|
[e200b7a] | 671 | endif |
---|
[8a65cb0] | 672 | endif |
---|
[f75967d] | 673 | end do |
---|
[8a65cb0] | 674 | !END OLD METHOD |
---|
[d6a0977] | 675 | end do |
---|
[f75967d] | 676 | end do |
---|
| 677 | endif !readclouds |
---|
[d6a0977] | 678 | |
---|
| 679 | !********* TEST *************** |
---|
| 680 | ! WRITE OUT SOME TEST VARIABLES |
---|
| 681 | !********* TEST ************'** |
---|
| 682 | !teller(:)=0 |
---|
| 683 | !virr=virr+1 |
---|
| 684 | !WRITE(aspec, '(i3.3)'), virr |
---|
| 685 | |
---|
| 686 | !if (readclouds) then |
---|
| 687 | !fnameH=trim(zhgpath)//trim(aspec)//'Vertical_placement.txt' |
---|
| 688 | !else |
---|
| 689 | !fnameH=trim(zhgpath)//trim(aspec)//'Vertical_placement_old.txt' |
---|
| 690 | !endif |
---|
| 691 | ! |
---|
| 692 | !OPEN(UNIT=118, FILE=fnameH,FORM='FORMATTED',STATUS = 'UNKNOWN') |
---|
| 693 | !do kz_inv=1,nz-1 |
---|
| 694 | ! kz=nz-kz_inv+1 |
---|
| 695 | ! !kz=91 |
---|
| 696 | ! do jy=0,nymin1 |
---|
| 697 | ! do ix=0,nxmin1 |
---|
| 698 | ! if (clouds(ix,jy,kz,n).eq.1) teller(1)=teller(1)+1 ! no precipitation cloud |
---|
| 699 | ! if (clouds(ix,jy,kz,n).eq.2) teller(2)=teller(2)+1 ! convp dominated rainout |
---|
| 700 | ! if (clouds(ix,jy,kz,n).eq.3) teller(3)=teller(3)+1 ! lsp dominated rainout |
---|
| 701 | ! if (clouds(ix,jy,kz,n).eq.4) teller(4)=teller(4)+1 ! convp dominated washout |
---|
| 702 | ! if (clouds(ix,jy,kz,n).eq.5) teller(5)=teller(5)+1 ! lsp dominated washout |
---|
| 703 | ! |
---|
| 704 | ! ! write(*,*) height(kz),teller |
---|
| 705 | ! end do |
---|
| 706 | ! end do |
---|
| 707 | ! write(118,*) height(kz),teller |
---|
| 708 | ! teller(:)=0 |
---|
| 709 | !end do |
---|
| 710 | !teller(:)=0 |
---|
| 711 | !write(*,*) teller |
---|
| 712 | !write(*,*) aspec |
---|
| 713 | ! |
---|
| 714 | !fnameA=trim(zhgpath)//trim(aspec)//'cloudV.txt' |
---|
| 715 | !fnameB=trim(zhgpath)//trim(aspec)//'cloudT.txt' |
---|
| 716 | !fnameC=trim(zhgpath)//trim(aspec)//'cloudB.txt' |
---|
| 717 | !fnameD=trim(zhgpath)//trim(aspec)//'cloudW.txt' |
---|
| 718 | !fnameE=trim(zhgpath)//trim(aspec)//'old_cloudV.txt' |
---|
| 719 | !fnameF=trim(zhgpath)//trim(aspec)//'lsp.txt' |
---|
| 720 | !fnameG=trim(zhgpath)//trim(aspec)//'convp.txt' |
---|
| 721 | !if (readclouds) then |
---|
| 722 | !OPEN(UNIT=111, FILE=fnameA,FORM='FORMATTED',STATUS = 'UNKNOWN') |
---|
| 723 | !OPEN(UNIT=112, FILE=fnameB,FORM='FORMATTED',STATUS = 'UNKNOWN') |
---|
| 724 | !OPEN(UNIT=113, FILE=fnameC,FORM='FORMATTED',STATUS = 'UNKNOWN') |
---|
| 725 | !OPEN(UNIT=114, FILE=fnameD,FORM='FORMATTED',STATUS = 'UNKNOWN') |
---|
| 726 | !else |
---|
| 727 | !OPEN(UNIT=115, FILE=fnameE,FORM='FORMATTED',STATUS = 'UNKNOWN') |
---|
| 728 | !OPEN(UNIT=116, FILE=fnameF,FORM='FORMATTED',STATUS = 'UNKNOWN') |
---|
| 729 | !OPEN(UNIT=117, FILE=fnameG,FORM='FORMATTED',STATUS = 'UNKNOWN') |
---|
| 730 | !endif |
---|
| 731 | ! |
---|
| 732 | !do ix=0,nxmin1 |
---|
| 733 | !if (readclouds) then |
---|
| 734 | !write(111,*) (icloud_stats(ix,jy,1,n),jy=0,nymin1) |
---|
| 735 | !write(112,*) (icloud_stats(ix,jy,2,n),jy=0,nymin1) |
---|
| 736 | !write(113,*) (icloud_stats(ix,jy,3,n),jy=0,nymin1) |
---|
| 737 | !write(114,*) (icloud_stats(ix,jy,4,n),jy=0,nymin1) |
---|
| 738 | !else |
---|
| 739 | !write(115,*) (cloudsh(ix,jy,n),jy=0,nymin1) !integer |
---|
| 740 | !write(116,*) (lsprec(ix,jy,1,n),jy=0,nymin1) !7.83691406E-02 |
---|
| 741 | !write(117,*) (convprec(ix,jy,1,n),jy=0,nymin1) !5.38330078E-02 |
---|
| 742 | !endif |
---|
| 743 | !end do |
---|
| 744 | ! |
---|
| 745 | !if (readclouds) then |
---|
| 746 | !CLOSE(111) |
---|
| 747 | !CLOSE(112) |
---|
| 748 | !CLOSE(113) |
---|
| 749 | !CLOSE(114) |
---|
| 750 | !else |
---|
| 751 | !CLOSE(115) |
---|
| 752 | !CLOSE(116) |
---|
| 753 | !CLOSE(117) |
---|
| 754 | !endif |
---|
| 755 | ! |
---|
| 756 | !END ********* TEST *************** END |
---|
| 757 | ! WRITE OUT SOME TEST VARIABLES |
---|
| 758 | !END ********* TEST *************** END |
---|
| 759 | |
---|
[8a65cb0] | 760 | |
---|
| 761 | ! PS 2012 |
---|
[4d45639] | 762 | ! lsp=lsprec(ix,jy,1,n) |
---|
| 763 | ! convp=convprec(ix,jy,1,n) |
---|
[8a65cb0] | 764 | ! prec=lsp+convp |
---|
| 765 | ! if (lsp.gt.convp) then ! prectype='lsp' |
---|
| 766 | ! lconvectprec = .false. |
---|
| 767 | ! else ! prectype='cp ' |
---|
| 768 | ! lconvectprec = .true. |
---|
[4d45639] | 769 | ! endif |
---|
[8a65cb0] | 770 | ! else ! windfields does not contain cloud data |
---|
| 771 | ! rhmin = 0.90 ! standard condition for presence of clouds |
---|
| 772 | !PS note that original by Sabine Eckhart was 80% |
---|
| 773 | !PS however, for T<-20 C we consider saturation over ice |
---|
| 774 | !PS so I think 90% should be enough |
---|
| 775 | ! icloudbot(ix,jy,n)=icmv |
---|
| 776 | ! icloudtop=icmv ! this is just a local variable |
---|
| 777 | !98 do kz=1,nz |
---|
| 778 | ! pressure=rho(ix,jy,kz,n)*r_air*tt(ix,jy,kz,n) |
---|
| 779 | ! rh=qv(ix,jy,kz,n)/f_qvsat(pressure,tt(ix,jy,kz,n)) |
---|
| 780 | !ps if (prec.gt.0.01) print*,'relhum',prec,kz,rh,height(kz) |
---|
| 781 | ! if (rh .gt. rhmin) then |
---|
| 782 | ! if (icloudbot(ix,jy,n) .eq. icmv) then |
---|
| 783 | ! icloudbot(ix,jy,n)=nint(height(kz)) |
---|
| 784 | ! endif |
---|
| 785 | ! icloudtop=nint(height(kz)) ! use int to save memory |
---|
| 786 | ! endif |
---|
[4d45639] | 787 | ! end do |
---|
[8a65cb0] | 788 | !PS try to get a cloud thicker than 50 m |
---|
| 789 | !PS if there is at least .01 mm/h - changed to 0.002 and put into |
---|
| 790 | !PS parameter precpmin |
---|
| 791 | ! if ((icloudbot(ix,jy,n) .eq. icmv .or. & |
---|
| 792 | ! icloudtop-icloudbot(ix,jy,n) .lt. 50) .and. & |
---|
| 793 | ! prec .gt. precmin) then |
---|
| 794 | ! rhmin = rhmin - 0.05 |
---|
| 795 | ! if (rhmin .ge. 0.30) goto 98 ! give up for <= 25% rel.hum. |
---|
[4d45639] | 796 | ! end if |
---|
[8a65cb0] | 797 | |
---|
| 798 | !PS is based on looking at a limited set of comparison data |
---|
| 799 | ! if (lconvectprec .and. icloudtop .lt. 6000 .and. & |
---|
| 800 | ! prec .gt. precmin) then |
---|
| 801 | ! |
---|
| 802 | ! if (convp .lt. 0.1) then |
---|
| 803 | ! icloudbot(ix,jy,n) = 500 |
---|
| 804 | ! icloudtop = 8000 |
---|
| 805 | ! else |
---|
| 806 | ! icloudbot(ix,jy,n) = 0 |
---|
| 807 | ! icloudtop = 10000 |
---|
| 808 | ! endif |
---|
| 809 | ! endif |
---|
| 810 | ! if (icloudtop .ne. icmv) then |
---|
| 811 | ! icloudthck(ix,jy,n) = icloudtop-icloudbot(ix,jy,n) |
---|
| 812 | ! else |
---|
| 813 | ! icloudthck(ix,jy,n) = icmv |
---|
| 814 | ! endif |
---|
| 815 | !PS get rid of too thin clouds |
---|
| 816 | ! if (icloudthck(ix,jy,n) .lt. 50) then |
---|
| 817 | ! icloudbot(ix,jy,n)=icmv |
---|
| 818 | ! icloudthck(ix,jy,n)=icmv |
---|
| 819 | ! endif |
---|
| 820 | !hg__________________________________ |
---|
| 821 | ! rcw(ix,jy)=2E-7*prec**0.36 |
---|
| 822 | ! rpc(ix,jy)=prec |
---|
| 823 | !hg end______________________________ |
---|
| 824 | |
---|
| 825 | ! endif !hg read clouds |
---|
| 826 | |
---|
[e200b7a] | 827 | |
---|
[d6a0977] | 828 | |
---|
[4d45639] | 829 | |
---|
| 830 | !eso measure CPU time |
---|
| 831 | ! call mpif_mtime('verttransform',1) |
---|
| 832 | |
---|
| 833 | !eso print out the same measure as in Leo's routine |
---|
| 834 | ! write(*,*) 'verttransform: ', & |
---|
| 835 | ! sum(tt(:,:,:,n)*tt(:,:,:,n)), & |
---|
| 836 | ! sum(uu(:,:,:,n)*uu(:,:,:,n)),sum(vv(:,:,:,n)*vv(:,:,:,n)),& |
---|
| 837 | ! sum(qv(:,:,:,n)*qv(:,:,:,n)),sum(pv(:,:,:,n)*pv(:,:,:,n)),& |
---|
| 838 | ! sum(rho(:,:,:,n)*rho(:,:,:,n)),sum(drhodz(:,:,:,n)*drhodz(:,:,:,n)),& |
---|
| 839 | ! sum(ww(:,:,:,n)*ww(:,:,:,n)), & |
---|
| 840 | ! sum(clouds(:,:,:,n)), sum(cloudsh(:,:,n)),sum(idx),sum(pinmconv) |
---|
[f75967d] | 841 | end subroutine verttransform |
---|
[4d45639] | 842 | |
---|