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