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