FV3DYCORE  Version 2.0.0
init_hydro.F90
Go to the documentation of this file.
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the FV3 dynamical core.
5 !*
6 !* The FV3 dynamical core is free software: you can redistribute it
7 !* and/or modify it under the terms of the
8 !* GNU Lesser General Public License as published by the
9 !* Free Software Foundation, either version 3 of the License, or
10 !* (at your option) any later version.
11 !*
12 !* The FV3 dynamical core is distributed in the hope that it will be
13 !* useful, but WITHOUT ANYWARRANTY; without even the implied warranty
14 !* of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
15 !* See the GNU General Public License for more details.
16 !*
17 !* You should have received a copy of the GNU Lesser General Public
18 !* License along with the FV3 dynamical core.
19 !* If not, see <http://www.gnu.org/licenses/>.
20 !***********************************************************************
21 
23 
24 ! <table>
25 ! <tr>
26 ! <th>Module Name</th>
27 ! <th>Functions Included</th>
28 ! </tr>
29 ! <tr>
30 ! <td>constants_mod</td>
31 ! <td>grav, rdgas, rvgas</td>
32 ! </tr>
33 ! <tr>
34 ! <td>fv_grid_utils_mod</td>
35 ! <td>g_sum</td>
36 ! </tr>
37 ! <tr>
38 ! <td>fv_mp_mod</td>
39 ! <td>is_master</td>
40 ! </tr>
41 ! <tr>
42 ! <td>mpp_mod</td>
43 ! <td>mpp_chksum, stdout, mpp_error, FATAL, NOTE,get_unit, mpp_sum, mpp_broadcast,
44 ! mpp_get_current_pelist, mpp_npes, mpp_set_current_pelist, mpp_send, mpp_recv,
45 ! mpp_sync_self, mpp_npes, mpp_pe, mpp_sync</td>
46 ! </tr>
47 ! <tr>
48 ! <td>mpp_domains_mod</td>
49 ! <td> domain2d</td>
50 ! </tr>
51 ! <tr>
52 ! <td>tracer_manager_mod</td>
53 ! <td>get_tracer_index</td>
54 ! </tr>
55 ! </table>
56 
57 
58  use constants_mod, only: grav, rdgas, rvgas
59  use fv_grid_utils_mod, only: g_sum
60  use fv_mp_mod, only: is_master
61  use field_manager_mod, only: model_atmos
62  use tracer_manager_mod, only: get_tracer_index
63  use mpp_domains_mod, only: domain2d
64  use fv_arrays_mod, only: r_grid
65 ! use fv_diagnostics_mod, only: prt_maxmin
66 #ifdef MULTI_GASES
67  use multi_gases_mod, only: virq, virqd, vicpqd
68 #endif
69 
70  implicit none
71  private
72 
73  public :: p_var, hydro_eq
74 
75 contains
76 
77 !-------------------------------------------------------------------------------
82  subroutine p_var(km, ifirst, ilast, jfirst, jlast, ptop, ptop_min, &
83  delp, delz, pt, ps, pe, peln, pk, pkz, cappa, q, ng, nq, area, &
84  dry_mass, adjust_dry_mass, mountain, moist_phys, &
85  hydrostatic, nwat, domain, adiabatic, make_nh)
86 
87 ! Given (ptop, delp) computes (ps, pk, pe, peln, pkz)
88 ! Input:
89  integer, intent(in):: km
90  integer, intent(in):: ifirst, ilast
91  integer, intent(in):: jfirst, jlast
92  integer, intent(in):: nq, nwat
93  integer, intent(in):: ng
94  logical, intent(in):: adjust_dry_mass, mountain, moist_phys, hydrostatic, adiabatic
95  real, intent(in):: dry_mass, cappa, ptop, ptop_min
96  real, intent(in ):: pt(ifirst-ng:ilast+ng,jfirst-ng:jlast+ng, km)
97  real, intent(inout):: delz(ifirst:ilast,jfirst:jlast, km)
98  real, intent(inout):: delp(ifirst-ng:ilast+ng,jfirst-ng:jlast+ng, km)
99  real, intent(inout):: q(ifirst-ng:ilast+ng,jfirst-ng:jlast+ng, km, nq)
100  real(kind=R_GRID), intent(IN) :: area(ifirst-ng:ilast+ng,jfirst-ng:jlast+ng)
101  logical, optional:: make_nh
102 ! Output:
103  real, intent(out) :: ps(ifirst-ng:ilast+ng, jfirst-ng:jlast+ng)
104  real, intent(out) :: pk(ifirst:ilast, jfirst:jlast, km+1)
105  real, intent(out) :: pe(ifirst-1:ilast+1,km+1,jfirst-1:jlast+1)
106  real, intent(out) :: peln(ifirst:ilast, km+1, jfirst:jlast)
107  real, intent(out) :: pkz(ifirst:ilast, jfirst:jlast, km)
108  type(domain2d), intent(IN) :: domain
109 
110 ! Local
111  integer sphum, liq_wat, ice_wat
112  integer rainwat, snowwat, graupel ! GFDL Cloud Microphysics
113  real ratio(ifirst:ilast)
114  real pek, lnp, ak1, rdg, dpd, zvir
115  integer i, j, k, n
116 
117 ! Check dry air mass & compute the adjustment amount:
118  if ( adjust_dry_mass ) &
119  call drymadj(km, ifirst, ilast, jfirst, jlast, ng, cappa, ptop, ps, &
120  delp, q, nq, area, nwat, dry_mass, adjust_dry_mass, moist_phys, dpd, domain)
121 
122  pek = ptop ** cappa
123 
124 !$OMP parallel do default(none) shared(ifirst,ilast,jfirst,jlast,km,ptop,pek,pe,pk, &
125 !$OMP ps,adjust_dry_mass,dpd,delp,peln,cappa, &
126 !$OMP ptop_min,hydrostatic,pkz ) &
127 !$OMP private(ratio, ak1, lnp)
128  do j=jfirst,jlast
129  do i=ifirst,ilast
130  pe(i,1,j) = ptop
131  pk(i,j,1) = pek
132  enddo
133 
134  if ( adjust_dry_mass ) then
135  do i=ifirst,ilast
136  ratio(i) = 1. + dpd/(ps(i,j)-ptop)
137  enddo
138  do k=1,km
139  do i=ifirst,ilast
140  delp(i,j,k) = delp(i,j,k) * ratio(i)
141  enddo
142  enddo
143  endif
144 
145  do k=2,km+1
146  do i=ifirst,ilast
147  pe(i,k,j) = pe(i,k-1,j) + delp(i,j,k-1)
148  peln(i,k,j) = log(pe(i,k,j))
149  pk(i,j,k) = exp( cappa*peln(i,k,j) )
150  enddo
151  enddo
152 
153  do i=ifirst,ilast
154  ps(i,j) = pe(i,km+1,j)
155  enddo
156 
157  if( ptop < ptop_min ) then
158 !---- small ptop modification -------------
159  ak1 = (cappa + 1.) / cappa
160  do i=ifirst,ilast
161  peln(i,1,j) = peln(i,2,j) - ak1
162  enddo
163  else
164  lnp = log( ptop )
165  do i=ifirst,ilast
166  peln(i,1,j) = lnp
167  enddo
168  endif
169 
170  if ( hydrostatic ) then
171  do k=1,km
172  do i=ifirst,ilast
173  pkz(i,j,k) = (pk(i,j,k+1)-pk(i,j,k))/(cappa*(peln(i,k+1,j)-peln(i,k,j)))
174  enddo
175  enddo
176  endif
177  enddo
178 
179  if ( adiabatic ) then
180  zvir = 0.
181  else
182  zvir = rvgas/rdgas - 1.
183  endif
184  sphum = get_tracer_index(model_atmos, 'sphum')
185 
186  if ( .not.hydrostatic ) then
187 
188  rdg = -rdgas / grav
189  if ( present(make_nh) ) then
190  if ( make_nh ) then
191  delz = 1.e25
192 !$OMP parallel do default(none) shared(ifirst,ilast,jfirst,jlast,km,delz,rdg,pt,peln,zvir,sphum,q)
193  do k=1,km
194  do j=jfirst,jlast
195  do i=ifirst,ilast
196  delz(i,j,k) = rdg*pt(i,j,k)*(1.+zvir*q(i,j,k,sphum))*(peln(i,k+1,j)-peln(i,k,j))
197  enddo
198  enddo
199  enddo
200  if(is_master()) write(*,*) 'delz computed from hydrostatic state'
201  endif
202  endif
203 
204  if ( moist_phys ) then
205 !------------------------------------------------------------------
206 ! The following form is the same as in "fv_update_phys.F90"
207 !------------------------------------------------------------------
208 !$OMP parallel do default(none) shared(ifirst,ilast,jfirst,jlast,km,pkz,cappa,rdg, &
209 !$OMP delp,pt,zvir,q,sphum,delz)
210  do k=1,km
211  do j=jfirst,jlast
212  do i=ifirst,ilast
213 #ifdef MULTI_GASES
214  pkz(i,j,k) = exp( cappa*(virqd(q(i,j,k,:))/vicpqd(q(i,j,k,:))) * &
215  log(rdg*delp(i,j,k)*pt(i,j,k)* &
216  virq(q(i,j,k,:)) /delz(i,j,k)) )
217 #else
218  pkz(i,j,k) = exp( cappa*log(rdg*delp(i,j,k)*pt(i,j,k)* &
219  (1.+zvir*q(i,j,k,sphum))/delz(i,j,k)) )
220 #endif
221  enddo
222  enddo
223  enddo
224  else
225 !$OMP parallel do default(none) shared(ifirst,ilast,jfirst,jlast,km,pkz,cappa,rdg, &
226 #ifdef MULTI_GASES
227 !$OMP q, &
228 #endif
229 !$OMP delp,pt,delz)
230  do k=1,km
231  do j=jfirst,jlast
232  do i=ifirst,ilast
233 #ifdef MULTI_GASES
234  pkz(i,j,k) = exp( cappa * (virqd(q(i,j,k,:))/vicpqd(q(i,j,k,:))) * &
235  log(rdg*delp(i,j,k)*pt(i,j,k)/delz(i,j,k)) )
236 #else
237  pkz(i,j,k) = exp( cappa*log(rdg*delp(i,j,k)*pt(i,j,k)/delz(i,j,k)) )
238 #endif
239  enddo
240  enddo
241  enddo
242  endif
243 
244  endif
245 
246  end subroutine p_var
247 
248 
249 
250  subroutine drymadj(km, ifirst, ilast, jfirst, jlast, ng, &
251  cappa, ptop, ps, delp, q, nq, area, nwat, &
252  dry_mass, adjust_dry_mass, moist_phys, dpd, domain)
254 ! INPUT PARAMETERS:
255  integer km
256  integer ifirst, ilast
257  integer jfirst, jlast
258  integer nq, ng, nwat
259  real, intent(in):: dry_mass
260  real, intent(in):: ptop
261  real, intent(in):: cappa
262  logical, intent(in):: adjust_dry_mass
263  logical, intent(in):: moist_phys
264  real(kind=R_GRID), intent(IN) :: area(ifirst-ng:ilast+ng, jfirst-ng:jlast+ng)
265  type(domain2d), intent(IN) :: domain
266 
267 ! INPUT/OUTPUT PARAMETERS:
268  real, intent(in):: q(ifirst-ng:ilast+ng,jfirst-ng:jlast+ng,km,nq)
269  real, intent(in)::delp(ifirst-ng:ilast+ng,jfirst-ng:jlast+ng,km)
270  real, intent(inout):: ps(ifirst-ng:ilast+ng,jfirst-ng:jlast+ng) ! surface pressure
271  real, intent(out):: dpd
272 ! Local
273  real psd(ifirst:ilast,jfirst:jlast)
274  real psmo, psdry
275  integer i, j, k
276 
277 !$OMP parallel do default(none) shared(ifirst,ilast,jfirst,jlast,km,ps,ptop,psd,delp,nwat,q)
278  do j=jfirst,jlast
279 
280  do i=ifirst,ilast
281  ps(i,j) = ptop
282  psd(i,j) = ptop
283  enddo
284 
285  do k=1,km
286  do i=ifirst,ilast
287  ps(i,j) = ps(i,j) + delp(i,j,k)
288  enddo
289  enddo
290 
291  if ( nwat>=1 ) then
292  do k=1,km
293  do i=ifirst,ilast
294  psd(i,j) = psd(i,j) + delp(i,j,k)*(1. - sum(q(i,j,k,1:nwat)))
295  enddo
296  enddo
297  else
298  do i=ifirst,ilast
299  psd(i,j) = ps(i,j)
300  enddo
301  endif
302  enddo
303 
304 ! Check global maximum/minimum
305 #ifndef QUICK_SUM
306  psdry = g_sum(domain, psd, ifirst, ilast, jfirst, jlast, ng, area, 1, .true.)
307  psmo = g_sum(domain, ps(ifirst:ilast,jfirst:jlast), ifirst, ilast, jfirst, jlast, &
308  ng, area, 1, .true.)
309 #else
310  psdry = g_sum(domain, psd, ifirst, ilast, jfirst, jlast, ng, area, 1)
311  psmo = g_sum(domain, ps(ifirst:ilast,jfirst:jlast), ifirst, ilast, jfirst, jlast, &
312  ng, area, 1)
313 #endif
314 
315  if(is_master()) then
316  write(*,*) 'Total surface pressure (mb) = ', 0.01*psmo
317  if ( moist_phys ) then
318  write(*,*) 'mean dry surface pressure = ', 0.01*psdry
319  write(*,*) 'Total Water (kg/m**2) =', real(psmo-psdry,4)/GRAV
320  endif
321  endif
322 
323  if( adjust_dry_mass ) Then
324  dpd = real(dry_mass - psdry,4)
325  if(is_master()) write(*,*) 'dry mass to be added (pascals) =', dpd
326  endif
327 
328  end subroutine drymadj
329 
332  subroutine hydro_eq(km, is, ie, js, je, ps, hs, drym, delp, ak, bk, &
333  pt, delz, area, ng, mountain, hydrostatic, hybrid_z, domain)
334 ! Input:
335  integer, intent(in):: is, ie, js, je, km, ng
336  real, intent(in):: ak(km+1), bk(km+1)
337  real, intent(in):: hs(is-ng:ie+ng,js-ng:je+ng)
338  real, intent(in):: drym
339  logical, intent(in):: mountain
340  logical, intent(in):: hydrostatic
341  logical, intent(in):: hybrid_z
342  real(kind=R_GRID), intent(IN) :: area(is-ng:ie+ng,js-ng:je+ng)
343  type(domain2d), intent(IN) :: domain
344 ! Output
345  real, intent(out):: ps(is-ng:ie+ng,js-ng:je+ng)
346  real, intent(out):: pt(is-ng:ie+ng,js-ng:je+ng,km)
347  real, intent(out):: delp(is-ng:ie+ng,js-ng:je+ng,km)
348  real, intent(inout):: delz(is:,js:,1:)
349 ! Local
350  real gz(is:ie,km+1)
351  real ph(is:ie,km+1)
352  real mslp, z1, t1, p1, t0, a0, psm
353  real ztop, c0
354 #ifdef INIT_4BYTE
355  real(kind=4) :: dps
356 #else
357  real dps ! note that different PEs will get differt dps during initialization
358  ! this has no effect after cold start
359 #endif
360  real p0, gztop, ptop
361  integer i,j,k
362 
363  if ( is_master() ) write(*,*) 'Initializing ATM hydrostatically'
364 
365  if ( is_master() ) write(*,*) 'Initializing Earth'
366 ! Given p1 and z1 (250mb, 10km)
367  p1 = 25000.
368  z1 = 10.e3 * grav
369  t1 = 200.
370  t0 = 300. ! sea-level temp.
371  a0 = (t1-t0)/z1
372  c0 = t0/a0
373 
374  if ( hybrid_z ) then
375  ptop = 100. ! *** hardwired model top ***
376  else
377  ptop = ak(1)
378  endif
379 
380  ztop = z1 + (rdgas*t1)*log(p1/ptop)
381  if(is_master()) write(*,*) 'ZTOP is computed as', ztop/grav*1.e-3
382 
383  if ( mountain ) then
384  mslp = 100917.4
385  do j=js,je
386  do i=is,ie
387  ps(i,j) = mslp*( c0/(hs(i,j)+c0))**(1./(a0*rdgas))
388  enddo
389  enddo
390  psm = g_sum(domain, ps(is:ie,js:je), is, ie, js, je, ng, area, 1, .true.)
391  dps = drym - psm
392  if(is_master()) write(*,*) 'Computed mean ps=', psm
393  if(is_master()) write(*,*) 'Correction delta-ps=', dps
394  else
395  mslp = drym ! 1000.E2
396  do j=js,je
397  do i=is,ie
398  ps(i,j) = mslp
399  enddo
400  enddo
401  dps = 0.
402  endif
403 
404 
405  do j=js,je
406  do i=is,ie
407  ps(i,j) = ps(i,j) + dps
408  gz(i, 1) = ztop
409  gz(i,km+1) = hs(i,j)
410  ph(i, 1) = ptop
411  ph(i,km+1) = ps(i,j)
412  enddo
413 
414  if ( hybrid_z ) then
415 !---------------
416 ! Hybrid Z
417 !---------------
418  do k=km,2,-1
419  do i=is,ie
420  gz(i,k) = gz(i,k+1) - delz(i,j,k)*grav
421  enddo
422  enddo
423 ! Correct delz at the top:
424  do i=is,ie
425  delz(i,j,1) = (gz(i,2) - ztop) / grav
426  enddo
427 
428  do k=2,km
429  do i=is,ie
430  if ( gz(i,k) >= z1 ) then
431 ! Isothermal
432  ph(i,k) = ptop*exp( (gz(i,1)-gz(i,k))/(rdgas*t1) )
433  else
434 ! Constant lapse rate region (troposphere)
435  ph(i,k) = ps(i,j)*((hs(i,j)+c0)/(gz(i,k)+c0))**(1./(a0*rdgas))
436  endif
437  enddo
438  enddo
439  else
440 !---------------
441 ! Hybrid sigma-p
442 !---------------
443  do k=2,km+1
444  do i=is,ie
445  ph(i,k) = ak(k) + bk(k)*ps(i,j)
446  enddo
447  enddo
448 
449  do k=2,km
450  do i=is,ie
451  if ( ph(i,k) <= p1 ) then
452 ! Isothermal
453  gz(i,k) = ztop + (rdgas*t1)*log(ptop/ph(i,k))
454  else
455 ! Constant lapse rate region (troposphere)
456  gz(i,k) = (hs(i,j)+c0)/(ph(i,k)/ps(i,j))**(a0*rdgas) - c0
457  endif
458  enddo
459  enddo
460  if ( .not. hydrostatic ) then
461  do k=1,km
462  do i=is,ie
463  delz(i,j,k) = ( gz(i,k+1) - gz(i,k) ) / grav
464  enddo
465  enddo
466  endif
467  endif ! end hybrid_z
468 
469 ! Convert geopotential to Temperature
470  do k=1,km
471  do i=is,ie
472  pt(i,j,k) = (gz(i,k)-gz(i,k+1))/(rdgas*(log(ph(i,k+1)/ph(i,k))))
473  pt(i,j,k) = max(t1, pt(i,j,k))
474  delp(i,j,k) = ph(i,k+1) - ph(i,k)
475  enddo
476  enddo
477  enddo ! j-loop
478 
479 
480  end subroutine hydro_eq
481 
482 
483 end module init_hydro_mod
pure real function, public vicpqd(q)
The module &#39;fv_mp_mod&#39; is a single program multiple data (SPMD) parallel decompostion/communication m...
Definition: fv_mp_mod.F90:24
The module &#39;multi_gases&#39; peforms multi constitutents computations.
Definition: multi_gases.F90:25
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, adiabatic, make_nh)
the subroutine &#39;p_var&#39; computes auxiliary pressure variables for a hydrostatic state.
Definition: init_hydro.F90:86
real function, public g_sum(domain, p, ifirst, ilast, jfirst, jlast, ngc, area, mode, reproduce)
The function &#39;g_sum&#39; is the fast version of &#39;globalsum&#39;.
integer, parameter, public r_grid
Definition: fv_arrays.F90:34
pure real function, public virqd(q)
subroutine, public hydro_eq(km, is, ie, js, je, ps, hs, drym, delp, ak, bk, pt, delz, area, ng, mountain, hydrostatic, hybrid_z, domain)
The subroutine &#39;hydro_eq&#39; computes a hydrostatically balanced and isothermal basic state from input h...
Definition: init_hydro.F90:334
pure real function, public virq(q)
The module &#39;fv_arrays&#39; contains the &#39;fv_atmos_type&#39; and associated datatypes.
Definition: fv_arrays.F90:24
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)
Definition: init_hydro.F90:253
The module &#39;fv_grid_utils&#39; contains routines for setting up and computing grid-related quantities...