FV3DYCORE  Version1.0.0
multi_gases.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</td>
36 ! </tr>
37 ! </table>
38 
39  use constants_mod, only: rdgas, cp_air
40  use fv_mp_mod, only: is_master
41 
42 
43  implicit none
44  integer num_gas
45  integer ind_gas
46  integer num_wat
47  integer sphum, sphump1
48  real, allocatable :: ri(:)
49  real, allocatable :: cpi(:)
50  real, allocatable :: vir(:)
51  real, allocatable :: vicp(:)
52  real, allocatable :: vicv(:)
53 
54  private num_wat, sphum, sphump1
55  public vir, vicp, vicv, ind_gas, num_gas
56  public multi_gases_init
57  public virq
58  public virq_max
59  public virqd
60  public vicpqd
61  public virq_nodq
62  public virq_qpz
63  public vicpqd_qpz
64  public vicvqd
65  public vicvqd_qpz
66 
67  CONTAINS
68 ! --------------------------------------------------------
69  subroutine multi_gases_init(ngas, nwat)
70 !--------------------------------------------
71 ! !OUTPUT PARAMETERS
72 ! Ouput: vir(i): ri/rdgas - r0/rdgas
73 ! vir(0): r0/rdgas
74 ! vicp(i): cpi/cp_air - cp0/cp_air
75 ! vicp(0): cp0/cp_air
76 ! cv_air = cp_air - rdgas
77 ! vicv(i): cvi/cv_air - cv0/cv_air
78 ! vicv(0): cv0/cv_air
79 !--------------------------------------------
80  integer, intent(in):: ngas, nwat
81 ! Local:
82  integer n
83  real cvi(0:ngas)
84  real cv_air
85  logical :: default_gas=.false.
86 
87  sphum = 1
88  sphump1 = sphum+1
89  num_wat = nwat
90  ind_gas = num_wat+1
91  do n=0,ngas
92  if( ri(n).ne.0.0 .or. cpi(n).ne.0.0 ) num_gas=n
93  enddo
94  if ( num_gas.eq.1 ) default_gas=.true.
95  allocate( vir(0:num_gas) )
96  allocate( vicp(0:num_gas) )
97  allocate( vicv(0:num_gas) )
98 
99  cv_air = cp_air - rdgas
100  do n=0,num_gas
101  cvi(n) = cpi(n) - ri(n)
102  enddo
103 
104  vir(0) = ri(0)/rdgas
105  vicp(0) = cpi(0)/cp_air
106  vicv(0) = cvi(0)/cv_air
107  if( default_gas ) then
108  vir(0) = 1.0
109  vicp(0) = 1.0
110  vicv(0) = 1.0
111  endif
112  do n=1,num_gas
113  vir(n) = 0.0
114  if( ri(n).gt.0.0 ) vir(n) = ri(n)/rdgas - vir(0)
115  vicp(n) = 0.0
116  if( cpi(n).gt.0.0 ) vicp(n) = cpi(n)/cp_air - vicp(0)
117  vicv(n) = 0.0
118  if( cvi(n).gt.0.0 ) vicv(n) = cvi(n)/cv_air - vicv(0)
119  enddo
120 
121  if( is_master() ) then
122  write(*,*) ' multi_gases_init with ind_gas=',ind_gas
123  write(*,*) ' multi_gases_init with num_gas=',num_gas
124  write(*,*) ' multi_gases_init with vir =',vir
125  write(*,*) ' multi_gases_init with vicp=',vicp
126  write(*,*) ' multi_gases_init with vicv=',vicv
127  endif
128 
129  return
130  end subroutine multi_gases_init
131 
132 ! ----------------------------------------------------------------
133 
134 ! --------------------------------------------------------
135  pure real function virq(q)
136 !--------------------------------------------
137 ! !OUTPUT PARAMETERS
138 ! Ouput: variable gas 1+zvir/(1-qc)
139 !--------------------------------------------
140  real, intent(in) :: q(num_gas)
141 ! Local:
142  integer :: n
143 
144  virq = vir(sphum)*q(sphum)
145  do n=ind_gas,num_gas
146  virq = virq+vir(n)*q(sphum+n-1)
147  end do
148  virq = vir(0)+virq/(1.0-sum(q(sphump1:sphum+num_wat-1)))
149 
150  return
151  end function
152 !--------------------------------------------
153 
154 ! --------------------------------------------------------
155  pure real function virq_nodq(q)
156 !--------------------------------------------
157 ! !OUTPUT PARAMETERS
158 ! Ouput: variable gas 1+zvir without dividing by 1-qv or 1-qv-qc
159 !--------------------------------------------
160  real, intent(in) :: q(num_gas)
161 ! Local:
162  integer :: n
163 
164  virq_nodq = vir(0)+vir(sphum)*q(sphum)
165  do n=ind_gas,num_gas
166  virq_nodq = virq_nodq+vir(n)*q(sphum+n-1)
167  end do
168 
169  return
170  end function
171 !--------------------------------------------
172 
173 ! --------------------------------------------------------
174  pure real function virq_max(q, qmin)
175 !--------------------------------------------
176 ! !OUTPUT PARAMETERS
177 ! Ouput: variable gas 1+zvir using max(qmin,q(sphum))
178 !--------------------------------------------
179  real, intent(in) :: q(num_gas)
180  real, intent(in) :: qmin
181 ! Local:
182  integer :: n
183 
184  virq_max = vir(sphum)*max(qmin,q(sphum))
185  do n=ind_gas,num_gas
186  virq_max = virq_max+vir(n)*q(sphum+n-1)
187  end do
188  virq_max = vir(0)+virq_max/(1.0-sum(q(sphump1:sphum+num_wat-1)))
189 
190 
191  return
192  end function
193 !--------------------------------------------
194 
195 ! --------------------------------------------------------
196  pure real function virq_qpz(q, qpz)
197 !--------------------------------------------
198 ! !OUTPUT PARAMETERS
199 ! Ouput: variable gas 1+zvir/(1.-qpz): qpz in place of qv+qc from q
200 !--------------------------------------------
201  real, intent(in) :: q(num_gas)
202  real, intent(in) :: qpz
203 ! Local:
204  integer :: n
205 
206  virq_qpz = vir(sphum)*q(sphum)
207  do n=ind_gas,num_gas
208  virq_qpz = virq_qpz+vir(n)*q(sphum+n-1)
209  end do
210  virq_qpz = vir(0)+virq_qpz/(1.0-qpz)
211 
212 
213  return
214  end function
215 !--------------------------------------------
216 
217 ! --------------------------------------------------------
218  pure real function virqd(q)
219 !--------------------------------------------
220 ! !OUTPUT PARAMETERS
221 ! Ouput: variable gas 1+zvir/(1-(qv+qc)) (dry)
222 !--------------------------------------------
223  real, intent(in) :: q(num_gas)
224 ! Local:
225  integer :: n
226 
227  virqd = 0.0
228  do n=ind_gas,num_gas
229  virqd = virqd+vir(n)*q(sphum+n-1)
230  end do
231  virqd = vir(0)+virqd/(1.0-sum(q(sphum:sphum+num_wat-1)))
232 
233  return
234  end function
235 !--------------------------------------------
236 
237 ! --------------------------------------------------------
238  pure real function vicpqd(q)
239 !--------------------------------------------
240 ! !OUTPUT PARAMETERS
241 ! Ouput: variable gas cp (dry)
242 !--------------------------------------------
243  real, intent(in) :: q(num_gas)
244 ! Local:
245  integer :: n
246 
247  vicpqd = 0.0
248  do n=ind_gas,num_gas
249  vicpqd = vicpqd+vicp(n)*q(sphum+n-1)
250  end do
251  vicpqd = vicp(0)+vicpqd/(1.0-sum(q(sphum:sphum+num_wat-1)))
252 
253  return
254  end function
255 !--------------------------------------------
256 
257 ! --------------------------------------------------------
258  pure real function vicpqd_qpz(q, qpz)
259 !--------------------------------------------
260 ! !OUTPUT PARAMETERS
261 ! Ouput: variable gas cp (dry) with qpz in place of qv+qc from q
262 !--------------------------------------------
263  real, intent(in) :: q(num_gas)
264  real, intent(in) :: qpz
265 ! Local:
266  integer :: n
267 
268  vicpqd_qpz = 0.0
269  do n=ind_gas,num_gas
270  vicpqd_qpz = vicpqd_qpz+vicp(n)*q(sphum+n-1)
271  end do
272  vicpqd_qpz = vicp(0)+vicpqd_qpz/(1.0-qpz)
273 
274  return
275  end function
276 !--------------------------------------------
277 
278 ! --------------------------------------------------------
279  pure real function vicvqd(q)
280 !--------------------------------------------
281 ! !OUTPUT PARAMETERS
282 ! Ouput: variable gas cv (dry)
283 !--------------------------------------------
284  real, intent(in) :: q(num_gas)
285 ! Local:
286  integer :: n
287 
288  vicvqd = 0.0
289  do n=ind_gas,num_gas
290  vicvqd = vicvqd+vicv(n)*q(sphum+n-1)
291  end do
292  vicvqd = vicv(0)+vicvqd/(1.0-sum(q(sphum:sphum+num_wat-1)))
293 
294  return
295  end function
296 !--------------------------------------------
297 
298 ! --------------------------------------------------------
299  pure real function vicvqd_qpz(q,qpz)
300 !--------------------------------------------
301 ! !OUTPUT PARAMETERS
302 ! Ouput: variable gas cv (dry) with qpz in place of qv+qc from q
303 !--------------------------------------------
304  real, intent(in) :: q(num_gas)
305  real, intent(in) :: qpz
306 ! Local:
307  integer :: n
308 
309  vicvqd_qpz = 0.0
310  do n=ind_gas,num_gas
311  vicvqd_qpz = vicvqd_qpz+vicv(n)*q(sphum+n-1)
312  end do
313  vicvqd_qpz = vicv(0)+vicvqd_qpz/(1.0-qpz)
314 
315  return
316  end function
317 !--------------------------------------------
318 
319 end module multi_gases_mod
real, dimension(:), allocatable, public vir
Definition: multi_gases.F90:50
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
integer, public num_gas
Definition: multi_gases.F90:44
The module &#39;multi_gases&#39; peforms multi constitutents computations.
Definition: multi_gases.F90:25
pure real function, public vicvqd_qpz(q, qpz)
pure real function, public virq_max(q, qmin)
pure real function, public virq_qpz(q, qpz)
real, dimension(:), allocatable ri
Definition: multi_gases.F90:48
pure real function, public virqd(q)
pure real function, public virq_nodq(q)
pure real function, public virq(q)
integer, public ind_gas
Definition: multi_gases.F90:45
integer, private sphum
Definition: multi_gases.F90:47
integer, private num_wat
Definition: multi_gases.F90:46
subroutine, public multi_gases_init(ngas, nwat)
Definition: multi_gases.F90:70
real, dimension(:), allocatable cpi
Definition: multi_gases.F90:49
real, dimension(:), allocatable, public vicv
Definition: multi_gases.F90:52
real, dimension(:), allocatable, public vicp
Definition: multi_gases.F90:51
pure real function, public vicpqd_qpz(q, qpz)
pure real function, public vicvqd(q)
integer, private sphump1
Definition: multi_gases.F90:47