FV3DYCORE  Version 2.0.0
nh_core.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 
24 
26 
27 ! Modules Included:
28 ! <table>
29 ! <tr>
30 ! <th>Module Name</th>
31 ! <th>Functions Included</th>
32 ! </tr>
33 ! <tr>
34 ! <td>constants_mod</td>
35 ! <td>rdgas, cp_air, grav</td>
36 ! </tr>
37 ! <tr>
38 ! <td>nh_utils_mod</td>
39 ! <td>update_dz_c, update_dz_d, nest_halo_nh, sim3p0_solver, rim_2d,
40 ! sim_solver, sim1_solver, sim3_solver</td>
41 ! </tr>
42 ! <tr>
43 ! <td>tp_core_mod</td>
44 ! <td>fv_tp_2d</td>
45 ! </tr>
46 ! </table>
47 
48  use constants_mod, only: rdgas, cp_air, grav
49  use tp_core_mod, only: fv_tp_2d
53  use nh_utils_mod, only: riem_solver_c
54 
55  implicit none
56  private
57 
58  public riem_solver3, riem_solver_c, update_dz_c, update_dz_d, nh_bc
59  real, parameter:: r3 = 1./3.
60 
61 CONTAINS
62 
63  subroutine riem_solver3(ms, dt, is, ie, js, je, km, ng, &
64  isd, ied, jsd, jed, akap, cappa, cp, &
65 #ifdef MULTI_GASES
66  kapad, &
67 #endif
68  ptop, zs, q_con, w, delz, pt, &
69  delp, zh, pe, ppe, pk3, pk, peln, &
70  ws, scale_m, p_fac, a_imp, &
71  use_logp, last_call, fp_out)
72 !--------------------------------------------
73 ! !OUTPUT PARAMETERS
74 ! Ouput: gz: grav*height at edges
75 ! pe: full hydrostatic pressure
76 ! ppe: non-hydrostatic pressure perturbation
77 !--------------------------------------------
78  integer, intent(in):: ms, is, ie, js, je, km, ng
79  integer, intent(in):: isd, ied, jsd, jed
80  real, intent(in):: dt
81  real, intent(in):: akap, cp, ptop, p_fac, a_imp, scale_m
82  real, intent(in):: zs(isd:ied,jsd:jed)
83  logical, intent(in):: last_call, use_logp, fp_out
84  real, intent(in):: ws(is:ie,js:je)
85  real, intent(in), dimension(isd:,jsd:,1:):: q_con, cappa
86 #ifdef MULTI_GASES
87  real, intent(in), dimension(isd:ied,jsd:jed,km):: kapad
88 #endif
89  real, intent(in), dimension(isd:ied,jsd:jed,km):: delp, pt
90  real, intent(inout), dimension(isd:ied,jsd:jed,km+1):: zh
91  real, intent(inout), dimension(isd:ied,jsd:jed,km):: w
92  real, intent(inout):: pe(is-1:ie+1,km+1,js-1:je+1)
93  real, intent(out):: peln(is:ie,km+1,js:je)
94  real, intent(out), dimension(isd:ied,jsd:jed,km+1):: ppe
95  real, intent(out):: delz(is:ie,js:je,km)
96  real, intent(out):: pk(is:ie,js:je,km+1)
97  real, intent(out):: pk3(isd:ied,jsd:jed,km+1)
98 ! Local:
99  real, dimension(is:ie,km):: dm, dz2, pm2, w2, gm2, cp2
100  real, dimension(is:ie,km+1)::pem, pe2, peln2, peg, pelng
101 #ifdef MULTI_GASES
102  real, dimension(is:ie,km):: kapad2
103 #endif
104  real gama, rgrav, ptk, peln1
105  integer i, j, k
106 
107  gama = 1./(1.-akap)
108  rgrav = 1./grav
109  peln1 = log(ptop)
110  ptk = exp(akap*peln1)
111 
112 !$OMP parallel do default(none) shared(is,ie,js,je,km,delp,ptop,peln1,pk3,ptk,akap,rgrav,zh,pt, &
113 !$OMP w,a_imp,dt,gama,ws,p_fac,scale_m,ms,delz,last_call, &
114 #ifdef MULTI_GASES
115 !$OMP peln,pk,fp_out,ppe,use_logp,zs,pe,cappa,q_con,kapad ) &
116 !$OMP private(cp2, gm2, dm, dz2, pm2, pem, peg, pelng, pe2, peln2, w2,kapad2)
117 #else
118 !$OMP peln,pk,fp_out,ppe,use_logp,zs,pe,cappa,q_con ) &
119 !$OMP private(cp2, gm2, dm, dz2, pm2, pem, peg, pelng, pe2, peln2, w2)
120 #endif
121  do 2000 j=js, je
122 
123  do k=1,km
124  do i=is, ie
125  dm(i,k) = delp(i,j,k)
126 #ifdef MOIST_CAPPA
127  cp2(i,k) = cappa(i,j,k)
128 #endif
129 #ifdef MULTI_GASES
130  kapad2(i,k) = kapad(i,j,k)
131 #endif
132  enddo
133  enddo
134 
135  do i=is,ie
136  pem(i,1) = ptop
137  peln2(i,1) = peln1
138  pk3(i,j,1) = ptk
139 #ifdef USE_COND
140  peg(i,1) = ptop
141  pelng(i,1) = peln1
142 #endif
143  enddo
144  do k=2,km+1
145  do i=is,ie
146  pem(i,k) = pem(i,k-1) + dm(i,k-1)
147  peln2(i,k) = log(pem(i,k))
148 #ifdef USE_COND
149 ! Excluding contribution from condensates:
150 ! peln used during remap; pk3 used only for p_grad
151  peg(i,k) = peg(i,k-1) + dm(i,k-1)*(1.-q_con(i,j,k-1))
152  pelng(i,k) = log(peg(i,k))
153 #endif
154 !hmhj pk3 at interface , interface pk is using constant akap
155  pk3(i,j,k) = exp(akap*peln2(i,k))
156  enddo
157  enddo
158 
159  do k=1,km
160  do i=is, ie
161 #ifdef USE_COND
162  pm2(i,k) = (peg(i,k+1)-peg(i,k))/(pelng(i,k+1)-pelng(i,k))
163 
164 #ifdef MOIST_CAPPA
165  gm2(i,k) = 1. / (1.-cp2(i,k))
166 #endif
167 
168 #else
169  pm2(i,k) = dm(i,k)/(peln2(i,k+1)-peln2(i,k))
170 #endif
171  dm(i,k) = dm(i,k) * rgrav
172  dz2(i,k) = zh(i,j,k+1) - zh(i,j,k)
173  w2(i,k) = w(i,j,k)
174  enddo
175  enddo
176 
177 
178  if ( a_imp < -0.999 ) then
179  call sim3p0_solver(dt, is, ie, km, rdgas, gama, akap, &
180 #ifdef MULTI_GASES
181  kapad2, &
182 #endif
183  pe2, dm, &
184  pem, w2, dz2, pt(is:ie,j,1:km), ws(is,j), p_fac, scale_m )
185  elseif ( a_imp < -0.5 ) then
186  call sim3_solver(dt, is, ie, km, rdgas, gama, akap, &
187 #ifdef MULTI_GASES
188  kapad2, &
189 #endif
190  pe2, dm, &
191  pem, w2, dz2, pt(is:ie,j,1:km), ws(is,j), abs(a_imp), p_fac, scale_m)
192  elseif ( a_imp <= 0.5 ) then
193  call rim_2d(ms, dt, is, ie, km, rdgas, gama, gm2, &
194 #ifdef MULTI_GASES
195  kapad2, &
196 #endif
197  pe2, &
198  dm, pm2, w2, dz2, pt(is:ie,j,1:km), ws(is,j), .false.)
199  elseif ( a_imp > 0.999 ) then
200  call sim1_solver(dt, is, ie, km, rdgas, gama, gm2, cp2, akap, &
201 #ifdef MULTI_GASES
202  kapad2, &
203 #endif
204  pe2, dm, &
205  pm2, pem, w2, dz2, pt(is:ie,j,1:km), ws(is,j), p_fac)
206  else
207  call sim_solver(dt, is, ie, km, rdgas, gama, gm2, cp2, akap, &
208 #ifdef MULTI_GASES
209  kapad2, &
210 #endif
211  pe2, dm, &
212  pm2, pem, w2, dz2, pt(is:ie,j,1:km), ws(is,j), &
213  a_imp, p_fac, scale_m)
214  endif
215 
216 
217  do k=1, km
218  do i=is, ie
219  w(i,j,k) = w2(i,k)
220  delz(i,j,k) = dz2(i,k)
221  enddo
222  enddo
223 
224  if ( last_call ) then
225  do k=1,km+1
226  do i=is,ie
227  peln(i,k,j) = peln2(i,k)
228  pk(i,j,k) = pk3(i,j,k)
229  pe(i,k,j) = pem(i,k)
230  enddo
231  enddo
232  endif
233 
234  if( fp_out ) then
235  do k=1,km+1
236  do i=is, ie
237  ppe(i,j,k) = pe2(i,k) + pem(i,k)
238  enddo
239  enddo
240  else
241  do k=1,km+1
242  do i=is, ie
243  ppe(i,j,k) = pe2(i,k)
244  enddo
245  enddo
246  endif
247 
248  if ( use_logp ) then
249  do k=2,km+1
250  do i=is, ie
251  pk3(i,j,k) = peln2(i,k)
252  enddo
253  enddo
254  endif
255 
256  do i=is, ie
257  zh(i,j,km+1) = zs(i,j)
258  enddo
259  do k=km,1,-1
260  do i=is, ie
261  zh(i,j,k) = zh(i,j,k+1) - dz2(i,k)
262  enddo
263  enddo
264 
265 2000 continue
266 
267  end subroutine riem_solver3
268 
269 
270 end module nh_core_mod
subroutine, public sim_solver(dt, is, ie, km, rgas, gama, gm2, cp2, kappa, pe2, dm2, pm2, pem, w2, dz2, pt2, ws, alpha, p_fac, scale_m)
Definition: nh_utils.F90:1548
subroutine, public sim3_solver(dt, is, ie, km, rgas, gama, kappa, pe2, dm, pem, w2, dz2, pt2, ws, alpha, p_fac, scale_m)
Definition: nh_utils.F90:1010
subroutine, public fv_tp_2d(q, crx, cry, npx, npy, hord, fx, fy, xfx, yfx, gridstruct, bd, ra_x, ra_y, lim_fac, mfx, mfy, mass, nord, damp_c)
The subroutine &#39;fv_tp_2d&#39; contains the FV advection scheme .
Definition: tp_core.F90:109
real, parameter r3
Definition: nh_core.F90:59
subroutine, public nh_bc(ptop, grav, kappa, cp, delp, delzBC, pt, phis, pkc, gz, pk3, BC_step, BC_split, npx, npy, npz, bounded_domain, pkc_pertn, computepk3, fullhalo, bd)
Definition: nh_utils.F90:1910
The module &#39;nh_utils&#39; peforms non-hydrostatic computations.
Definition: nh_utils.F90:26
The module &#39;tp_core&#39; is a collection of routines to support FV transport.
Definition: tp_core.F90:24
subroutine, public update_dz_c(is, ie, js, je, km, ng, dt, dp0, zs, area, ut, vt, gz, ws, npx, npy, sw_corner, se_corner, ne_corner, nw_corner, bd, grid_type)
Definition: nh_utils.F90:76
subroutine, public riem_solver3(ms, dt, is, ie, js, je, km, ng, isd, ied, jsd, jed, akap, cappa, cp, ptop, zs, q_con, w, delz, pt, delp, zh, pe, ppe, pk3, pk, peln, ws, scale_m, p_fac, a_imp, use_logp, last_call, fp_out)
Definition: nh_core.F90:72
subroutine, public sim1_solver(dt, is, ie, km, rgas, gama, gm2, cp2, kappa, pe, dm2, pm2, pem, w2, dz2, pt2, ws, p_fac)
Definition: nh_utils.F90:1377
subroutine, public sim3p0_solver(dt, is, ie, km, rgas, gama, kappa, pe2, dm, pem, w2, dz2, pt2, ws, p_fac, scale_m)
Definition: nh_utils.F90:1197
subroutine, public rim_2d(ms, bdt, is, ie, km, rgas, gama, gm2, pe2, dm2, pm2, w2, dz2, pt2, ws, c_core)
Definition: nh_utils.F90:745
The module &#39;nh_core&#39; peforms non-hydrostatic computations include moisture effect in pt...
Definition: nh_core.F90:25
subroutine, public update_dz_d(ndif, damp, hord, is, ie, js, je, km, ng, npx, npy, area, rarea, dp0, zs, zh, crx, cry, xfx, yfx, ws, rdt, gridstruct, bd, lim_fac)
Definition: nh_utils.F90:216
subroutine, public riem_solver_c(ms, dt, is, ie, js, je, km, ng, akap, cappa, cp, ptop, hs, w3, pt, q_con, delp, gz, pef, ws, p_fac, a_imp, scale_m)
Definition: nh_utils.F90:336