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