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