FV3 Bundle
init_hydro_nlm.F90
Go to the documentation of this file.
1 !***********************************************************************
2 !* GNU General Public License *
3 !* This file is a part of fvGFS. *
4 !* *
5 !* fvGFS is free software; you can redistribute it and/or modify it *
6 !* and are expected to follow the terms of the GNU General Public *
7 !* License as published by the Free Software Foundation; either *
8 !* version 2 of the License, or (at your option) any later version. *
9 !* *
10 !* fvGFS is distributed in the hope that it will be useful, but *
11 !* WITHOUT ANY WARRANTY; without even the implied warranty of *
12 !* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU *
13 !* General Public License for more details. *
14 !* *
15 !* For the full text of the GNU General Public License, *
16 !* write to: Free Software Foundation, Inc., *
17 !* 675 Mass Ave, Cambridge, MA 02139, USA. *
18 !* or see: http://www.gnu.org/licenses/gpl.html *
19 !***********************************************************************
20 
22 
23  use constants_mod, only: grav, rdgas, rvgas
24  use fv_grid_utils_nlm_mod, only: g_sum
25  use fv_mp_nlm_mod, only: is_master
28  use mpp_domains_mod, only: domain2d
29  use fv_arrays_nlm_mod, only: r_grid
30 ! use fv_diagnostics_nlm_mod, only: prt_maxmin
31 !!! DEBUG CODE
32  use mpp_mod, only: mpp_pe
33 !!! END DEBUG CODE
34 
35  implicit none
36  private
37 
38  public :: p_var, hydro_eq
39 
40 contains
41 
42 !-------------------------------------------------------------------------------
43  subroutine p_var(km, ifirst, ilast, jfirst, jlast, ptop, ptop_min, &
44  delp, delz, pt, ps, pe, peln, pk, pkz, cappa, q, ng, nq, area, &
45  dry_mass, adjust_dry_mass, mountain, moist_phys, &
46  hydrostatic, nwat, domain, make_nh)
47 
48 ! Given (ptop, delp) computes (ps, pk, pe, peln, pkz)
49 ! Input:
50  integer, intent(in):: km
51  integer, intent(in):: ifirst, ilast ! Longitude strip
52  integer, intent(in):: jfirst, jlast ! Latitude strip
53  integer, intent(in):: nq, nwat
54  integer, intent(in):: ng
55  logical, intent(in):: adjust_dry_mass, mountain, moist_phys, hydrostatic
56  real, intent(in):: dry_mass, cappa, ptop, ptop_min
57  real, intent(in ):: pt(ifirst-ng:ilast+ng,jfirst-ng:jlast+ng, km)
58  real, intent(inout):: delz(ifirst-ng:ilast+ng,jfirst-ng:jlast+ng, km)
59  real, intent(inout):: delp(ifirst-ng:ilast+ng,jfirst-ng:jlast+ng, km)
60  real, intent(inout):: q(ifirst-ng:ilast+ng,jfirst-ng:jlast+ng, km, nq)
61  real(kind=R_GRID), intent(IN) :: area(ifirst-ng:ilast+ng,jfirst-ng:jlast+ng)
62  logical, optional:: make_nh
63 ! Output:
64  real, intent(out) :: ps(ifirst-ng:ilast+ng, jfirst-ng:jlast+ng)
65  real, intent(out) :: pk(ifirst:ilast, jfirst:jlast, km+1)
66  real, intent(out) :: pe(ifirst-1:ilast+1,km+1,jfirst-1:jlast+1) ! Ghosted Edge pressure
67  real, intent(out) :: peln(ifirst:ilast, km+1, jfirst:jlast) ! Edge pressure
68  real, intent(out) :: pkz(ifirst:ilast, jfirst:jlast, km)
69  type(domain2d), intent(IN) :: domain
70 
71 ! Local
72  integer sphum, liq_wat, ice_wat
73  integer rainwat, snowwat, graupel ! Lin Micro-physics
74  real ratio(ifirst:ilast)
75  real pek, lnp, ak1, rdg, dpd, zvir
76  integer i, j, k
77 
78 ! Check dry air mass & compute the adjustment amount:
79  if ( adjust_dry_mass ) &
80  call drymadj(km, ifirst, ilast, jfirst, jlast, ng, cappa, ptop, ps, &
81  delp, q, nq, area, nwat, dry_mass, adjust_dry_mass, moist_phys, dpd, domain)
82 
83  pek = ptop ** cappa
84 
85 !$OMP parallel do default(none) shared(ifirst,ilast,jfirst,jlast,km,ptop,pek,pe,pk, &
86 !$OMP ps,adjust_dry_mass,dpd,delp,peln,cappa, &
87 !$OMP ptop_min,hydrostatic,pkz ) &
88 !$OMP private(ratio, ak1, lnp)
89  do j=jfirst,jlast
90  do i=ifirst,ilast
91  pe(i,1,j) = ptop
92  pk(i,j,1) = pek
93  enddo
94 
95  if ( adjust_dry_mass ) then
96  do i=ifirst,ilast
97  ratio(i) = 1. + dpd/(ps(i,j)-ptop)
98  enddo
99  do k=1,km
100  do i=ifirst,ilast
101  delp(i,j,k) = delp(i,j,k) * ratio(i)
102  enddo
103  enddo
104  endif
105 
106  do k=2,km+1
107  do i=ifirst,ilast
108  pe(i,k,j) = pe(i,k-1,j) + delp(i,j,k-1)
109  peln(i,k,j) = log(pe(i,k,j))
110  pk(i,j,k) = exp( cappa*peln(i,k,j) )
111  enddo
112  enddo
113 
114  do i=ifirst,ilast
115  ps(i,j) = pe(i,km+1,j)
116  enddo
117 
118  if( ptop < ptop_min ) then
119 !---- small ptop modification -------------
120  ak1 = (cappa + 1.) / cappa
121  do i=ifirst,ilast
122  peln(i,1,j) = peln(i,2,j) - ak1
123  enddo
124  else
125  lnp = log( ptop )
126  do i=ifirst,ilast
127  peln(i,1,j) = lnp
128  enddo
129  endif
130 
131  if ( hydrostatic ) then
132  do k=1,km
133  do i=ifirst,ilast
134  pkz(i,j,k) = (pk(i,j,k+1)-pk(i,j,k))/(cappa*(peln(i,k+1,j)-peln(i,k,j)))
135  enddo
136  enddo
137  endif
138  enddo
139 
140 
141  if ( .not.hydrostatic ) then
142 
143  rdg = -rdgas / grav
144  if ( present(make_nh) ) then
145  if ( make_nh ) then
146  delz = 1.e25
147 !$OMP parallel do default(none) shared(ifirst,ilast,jfirst,jlast,km,delz,rdg,pt,peln)
148  do k=1,km
149  do j=jfirst,jlast
150  do i=ifirst,ilast
151  delz(i,j,k) = rdg*pt(i,j,k)*(peln(i,k+1,j)-peln(i,k,j))
152  enddo
153  enddo
154  enddo
155  if(is_master()) write(*,*) 'delz computed from hydrostatic state'
156  endif
157  endif
158 
159  if ( moist_phys ) then
160 !------------------------------------------------------------------
161 ! The following form is the same as in "fv_update_phys.F90"
162 !------------------------------------------------------------------
163  zvir = rvgas/rdgas - 1.
164 #ifdef MAPL_MODE
165  sphum = 1
166 #else
167  sphum = get_tracer_index(model_atmos, 'sphum')
168 #endif
169 !$OMP parallel do default(none) shared(ifirst,ilast,jfirst,jlast,km,pkz,cappa,rdg, &
170 !$OMP delp,pt,zvir,q,sphum,delz)
171  do k=1,km
172  do j=jfirst,jlast
173  do i=ifirst,ilast
174  pkz(i,j,k) = exp( cappa*log(rdg*delp(i,j,k)*pt(i,j,k)* &
175  (1.+zvir*q(i,j,k,sphum))/delz(i,j,k)) )
176  enddo
177  enddo
178  enddo
179  else
180 !$OMP parallel do default(none) shared(ifirst,ilast,jfirst,jlast,km,pkz,cappa,rdg, &
181 !$OMP delp,pt,delz)
182  do k=1,km
183  do j=jfirst,jlast
184  do i=ifirst,ilast
185  pkz(i,j,k) = exp( cappa*log(rdg*delp(i,j,k)*pt(i,j,k)/delz(i,j,k)) )
186  enddo
187  enddo
188  enddo
189  endif
190 
191  endif
192 
193  end subroutine p_var
194 
195 
196 
197  subroutine drymadj(km, ifirst, ilast, jfirst, jlast, ng, &
198  cappa, ptop, ps, delp, q, nq, area, nwat, &
199  dry_mass, adjust_dry_mass, moist_phys, dpd, domain)
201 ! !INPUT PARAMETERS:
202  integer km
203  integer ifirst, ilast ! Long strip
204  integer jfirst, jlast ! Latitude strip
205  integer nq, ng, nwat
206  real, intent(in):: dry_mass
207  real, intent(in):: ptop
208  real, intent(in):: cappa
209  logical, intent(in):: adjust_dry_mass
210  logical, intent(in):: moist_phys
211  real(kind=R_GRID), intent(IN) :: area(ifirst-ng:ilast+ng, jfirst-ng:jlast+ng)
212  type(domain2d), intent(IN) :: domain
213 
214 ! !INPUT/OUTPUT PARAMETERS:
215  real, intent(in):: q(ifirst-ng:ilast+ng,jfirst-ng:jlast+ng,km,nq)
216  real, intent(in)::delp(ifirst-ng:ilast+ng,jfirst-ng:jlast+ng,km) !
217  real, intent(inout):: ps(ifirst-ng:ilast+ng,jfirst-ng:jlast+ng) ! surface pressure
218  real, intent(out):: dpd
219 ! Local
220  real psd(ifirst:ilast,jfirst:jlast) ! surface pressure due to dry air mass
221  real psmo, psdry
222  integer i, j, k
223 
224 !$OMP parallel do default(none) shared(ifirst,ilast,jfirst,jlast,km,ps,ptop,psd,delp,nwat,q)
225  do j=jfirst,jlast
226 
227  do i=ifirst,ilast
228  ps(i,j) = ptop
229  psd(i,j) = ptop
230  enddo
231 
232  do k=1,km
233  do i=ifirst,ilast
234  ps(i,j) = ps(i,j) + delp(i,j,k)
235  enddo
236  enddo
237 
238  if ( nwat>=1 ) then
239  do k=1,km
240  do i=ifirst,ilast
241  psd(i,j) = psd(i,j) + delp(i,j,k)*(1. - sum(q(i,j,k,1:nwat)))
242  enddo
243  enddo
244  else
245  do i=ifirst,ilast
246  psd(i,j) = ps(i,j)
247  enddo
248  endif
249  enddo
250 
251 ! Check global maximum/minimum
252 #ifndef QUICK_SUM
253  psdry = g_sum(domain, psd, ifirst, ilast, jfirst, jlast, ng, area, 1, .true.)
254  psmo = g_sum(domain, ps(ifirst:ilast,jfirst:jlast), ifirst, ilast, jfirst, jlast, &
255  ng, area, 1, .true.)
256 #else
257  psdry = g_sum(domain, psd, ifirst, ilast, jfirst, jlast, ng, area, 1)
258  psmo = g_sum(domain, ps(ifirst:ilast,jfirst:jlast), ifirst, ilast, jfirst, jlast, &
259  ng, area, 1)
260 #endif
261 
262  if(is_master()) then
263  write(*,*) 'Total surface pressure (mb) = ', 0.01*psmo
264  if ( moist_phys ) then
265  write(*,*) 'mean dry surface pressure = ', 0.01*psdry
266  write(*,*) 'Total Water (kg/m**2) =', real(psmo-psdry,4)/GRAV
267  endif
268  endif
269 
270  if( adjust_dry_mass ) Then
271  dpd = real(dry_mass - psdry,4)
272  if(is_master()) write(*,*) 'dry mass to be added (pascals) =', dpd
273  endif
274 
275  end subroutine drymadj
276 
277 
278 
279  subroutine hydro_eq(km, is, ie, js, je, ps, hs, drym, delp, ak, bk, &
280  pt, delz, area, ng, mountain, hydrostatic, hybrid_z, domain)
281 ! Input:
282  integer, intent(in):: is, ie, js, je, km, ng
283  real, intent(in):: ak(km+1), bk(km+1)
284  real, intent(in):: hs(is-ng:ie+ng,js-ng:je+ng)
285  real, intent(in):: drym
286  logical, intent(in):: mountain
287  logical, intent(in):: hydrostatic
288  logical, intent(in):: hybrid_z
289  real(kind=R_GRID), intent(IN) :: area(is-ng:ie+ng,js-ng:je+ng)
290  type(domain2d), intent(IN) :: domain
291 ! Output
292  real, intent(out):: ps(is-ng:ie+ng,js-ng:je+ng)
293  real, intent(out):: pt(is-ng:ie+ng,js-ng:je+ng,km)
294  real, intent(out):: delp(is-ng:ie+ng,js-ng:je+ng,km)
295  real, intent(inout):: delz(is-ng:ie+ng,js-ng:je+ng,km)
296 ! Local
297  real gz(is:ie,km+1)
298  real ph(is:ie,km+1)
299  real mslp, z1, t1, p1, t0, a0, psm
300  real ztop, c0
301 #ifdef INIT_4BYTE
302  real(kind=4) :: dps
303 #else
304  real dps ! note that different PEs will get differt dps during initialization
305  ! this has no effect after cold start
306 #endif
307  real p0, gztop, ptop
308  integer i,j,k
309 
310  if ( is_master() ) write(*,*) 'Initializing ATM hydrostatically'
311 
312  if ( is_master() ) write(*,*) 'Initializing Earth'
313 ! Given p1 and z1 (250mb, 10km)
314  p1 = 25000.
315  z1 = 10.e3 * grav
316  t1 = 200.
317  t0 = 300. ! sea-level temp.
318  a0 = (t1-t0)/z1
319  c0 = t0/a0
320 
321  if ( hybrid_z ) then
322  ptop = 100. ! *** hardwired model top ***
323  else
324  ptop = ak(1)
325  endif
326 
327  ztop = z1 + (rdgas*t1)*log(p1/ptop)
328  if(is_master()) write(*,*) 'ZTOP is computed as', ztop/grav*1.e-3
329 
330  if ( mountain ) then
331  mslp = 100917.4
332  do j=js,je
333  do i=is,ie
334  ps(i,j) = mslp*( c0/(hs(i,j)+c0))**(1./(a0*rdgas))
335  enddo
336  enddo
337  psm = g_sum(domain, ps(is:ie,js:je), is, ie, js, je, ng, area, 1, .true.)
338  dps = drym - psm
339  if(is_master()) write(*,*) 'Computed mean ps=', psm
340  if(is_master()) write(*,*) 'Correction delta-ps=', dps
341  else
342  mslp = drym ! 1000.E2
343  do j=js,je
344  do i=is,ie
345  ps(i,j) = mslp
346  enddo
347  enddo
348  dps = 0.
349  endif
350 
351 
352  do j=js,je
353  do i=is,ie
354  ps(i,j) = ps(i,j) + dps
355  gz(i, 1) = ztop
356  gz(i,km+1) = hs(i,j)
357  ph(i, 1) = ptop
358  ph(i,km+1) = ps(i,j)
359  enddo
360 
361  if ( hybrid_z ) then
362 !---------------
363 ! Hybrid Z
364 !---------------
365  do k=km,2,-1
366  do i=is,ie
367  gz(i,k) = gz(i,k+1) - delz(i,j,k)*grav
368  enddo
369  enddo
370 ! Correct delz at the top:
371  do i=is,ie
372  delz(i,j,1) = (gz(i,2) - ztop) / grav
373  enddo
374 
375  do k=2,km
376  do i=is,ie
377  if ( gz(i,k) >= z1 ) then
378 ! Isothermal
379  ph(i,k) = ptop*exp( (gz(i,1)-gz(i,k))/(rdgas*t1) )
380  else
381 ! Constant lapse rate region (troposphere)
382  ph(i,k) = ps(i,j)*((hs(i,j)+c0)/(gz(i,k)+c0))**(1./(a0*rdgas))
383  endif
384  enddo
385  enddo
386  else
387 !---------------
388 ! Hybrid sigma-p
389 !---------------
390  do k=2,km+1
391  do i=is,ie
392  ph(i,k) = ak(k) + bk(k)*ps(i,j)
393  enddo
394  enddo
395 
396  do k=2,km
397  do i=is,ie
398  if ( ph(i,k) <= p1 ) then
399 ! Isothermal
400  gz(i,k) = ztop + (rdgas*t1)*log(ptop/ph(i,k))
401  else
402 ! Constant lapse rate region (troposphere)
403  gz(i,k) = (hs(i,j)+c0)/(ph(i,k)/ps(i,j))**(a0*rdgas) - c0
404  endif
405  enddo
406  enddo
407  if ( .not. hydrostatic ) then
408  do k=1,km
409  do i=is,ie
410  delz(i,j,k) = ( gz(i,k+1) - gz(i,k) ) / grav
411  enddo
412  enddo
413  endif
414  endif ! end hybrid_z
415 
416 ! Convert geopotential to Temperature
417  do k=1,km
418  do i=is,ie
419  pt(i,j,k) = (gz(i,k)-gz(i,k+1))/(rdgas*(log(ph(i,k+1)/ph(i,k))))
420  pt(i,j,k) = max(t1, pt(i,j,k))
421  delp(i,j,k) = ph(i,k+1) - ph(i,k)
422  enddo
423  enddo
424  enddo ! j-loop
425 
426 
427  end subroutine hydro_eq
428 
429 
430 end module init_hydro_nlm_mod
subroutine drymadj(km, ifirst, ilast, jfirst, jlast, ng, cappa, ptop, ps, delp, q, nq, area, nwat, dry_mass, adjust_dry_mass, moist_phys, dpd, domain)
integer, parameter, public model_atmos
subroutine, public p_var(km, ifirst, ilast, jfirst, jlast, ptop, ptop_min, delp, delz, pt, ps, pe, peln, pk, pkz, cappa, q, ng, nq, area, dry_mass, adjust_dry_mass, mountain, moist_phys, hydrostatic, nwat, domain, make_nh)
real, parameter, public rdgas
Gas constant for dry air [J/kg/deg].
Definition: constants.F90:77
Definition: mpp.F90:39
real, parameter, public rvgas
Gas constant for water vapor [J/kg/deg].
Definition: constants.F90:78
subroutine, public hydro_eq(km, is, ie, js, je, ps, hs, drym, delp, ak, bk, pt, delz, area, ng, mountain, hydrostatic, hybrid_z, domain)
integer, parameter, public r_grid
real, parameter, public grav
Acceleration due to gravity [m/s^2].
Definition: constants.F90:76
real function, public g_sum(domain, p, ifirst, ilast, jfirst, jlast, ngc, area, mode, reproduce)
#define max(a, b)
Definition: mosaic_util.h:33