FV3-JEDI
moisture_variables_mod.f90
Go to the documentation of this file.
1 ! (C) Copyright 2018-2019 UCAR
2 !
3 ! This software is licensed under the terms of the Apache Licence Version 2.0
4 ! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
5 
7 
11 
13 
14 implicit none
15 private
16 
17 public crtm_ade_efr
18 public crtm_mixratio
19 public crtm_mixratio_tl
20 public crtm_mixratio_ad
21 public rh_to_q
22 public rh_to_q_tl
23 public rh_to_q_ad
24 public q_to_rh
25 public q_to_rh_tl
26 public q_to_rh_ad
27 public get_qsat
28 public q4_to_q2
29 public q2_to_q4
30 
31 contains
32 
33 !>----------------------------------------------------------------------------
34 !> Compute cloud area density and effective radius for the crtm --------------
35 !>----------------------------------------------------------------------------
36 
37 subroutine crtm_ade_efr( geom,p,T,delp,sea_frac,q,ql,qi,ql_ade,qi_ade,ql_efr,qi_efr)
38 
39 implicit none
40 
41 !Arguments
42 type(fv3jedi_geom) , intent(in) :: geom
43 real(kind=kind_real), intent(in) :: p(geom%isc:geom%iec,geom%jsc:geom%jec, 1:geom%npz) !Pressure | Pa
44 real(kind=kind_real), intent(in) :: t(geom%isc:geom%iec,geom%jsc:geom%jec, 1:geom%npz) !Temperature | K
45 real(kind=kind_real), intent(in) :: delp(geom%isc:geom%iec,geom%jsc:geom%jec, 1:geom%npz) !Layer thickness | Pa
46 real(kind=kind_real), intent(in) :: sea_frac(geom%isc:geom%iec,geom%jsc:geom%jec) !Sea fraction | 1
47 real(kind=kind_real), intent(in) :: q(geom%isc:geom%iec,geom%jsc:geom%jec, 1:geom%npz) !Specific humidity | kg/kg
48 real(kind=kind_real), intent(in) :: ql(geom%isc:geom%iec,geom%jsc:geom%jec, 1:geom%npz) !Mixing ratio of cloud liquid water | kg/kg
49 real(kind=kind_real), intent(in) :: qi(geom%isc:geom%iec,geom%jsc:geom%jec, 1:geom%npz) !Mixing ratio of cloud ice water | kg/kg
50 
51 real(kind=kind_real), intent(out) :: ql_ade(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz) !area density for cloud liquid water | kg/m^2
52 real(kind=kind_real), intent(out) :: qi_ade(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz) !area density for cloud ice water | kg/m^2
53 real(kind=kind_real), intent(out) :: ql_efr(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz) !efr for cloud liquid water | microns
54 real(kind=kind_real), intent(out) :: qi_efr(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz) !efr for cloud ice | microns
55 
56 !Locals
57 integer :: isc,iec,jsc,jec,npz
58 integer :: i,j,k
59 logical, allocatable :: seamask(:,:)
60 real(kind=kind_real) :: tem1, tem2, tem3, kgkg_to_kgm2
61 
62 
63 ! Grid convenience
64 ! ----------------
65 isc = geom%isc
66 iec = geom%iec
67 jsc = geom%jsc
68 jec = geom%jec
69 npz = geom%npz
70 
71 
72 ! Set outputs to zero
73 ! -------------------
74 ql_ade = 0.0_kind_real
75 qi_ade = 0.0_kind_real
76 ql_efr = 0.0_kind_real
77 qi_efr = 0.0_kind_real
78 
79 
80 ! Sea mask
81 ! --------
82 allocate(seamask(isc:iec,jsc:jec))
83 seamask = .false.
84 do j = jsc,jec
85  do i = isc,iec
86  seamask(i,j) = min(max(0.0_kind_real,sea_frac(i,j)),1.0_kind_real) >= 0.99_kind_real
87  enddo
88 enddo
89 
90 
91 ! Convert clouds from kg/kg into kg/m^2
92 ! -------------------------------------
93 do k = 1,npz
94  do j = jsc,jec
95  do i = isc,iec
96  if (seamask(i,j)) then
97 
98  kgkg_to_kgm2 = delp(i,j,k) / grav
99  ql_ade(i,j,k) = max(ql(i,j,k),0.0_kind_real) * kgkg_to_kgm2
100  qi_ade(i,j,k) = max(qi(i,j,k),0.0_kind_real) * kgkg_to_kgm2
101 
102  if (t(i,j,k) - tice > -20.0_kind_real) then
103  ql_ade(i,j,k) = max(1.001_kind_real*1.0e-6_kind_real,ql_ade(i,j,k))
104  endif
105 
106  if (t(i,j,k) < tice) then
107  qi_ade(i,j,k) = max(1.001_kind_real*1.0e-6_kind_real,qi_ade(i,j,k))
108  endif
109 
110  endif
111  enddo
112  enddo
113 enddo
114 
115 ! Cloud liquid water effective radius
116 ! -----------------------------------
117 do k = 1,npz
118  do j = jsc,jec
119  do i = isc,iec
120  if (seamask(i,j)) then
121  tem1 = max(0.0_kind_real,(tice-t(i,j,k))*0.05_kind_real)
122  ql_efr(i,j,k) = 5.0_kind_real + 5.0_kind_real * min(1.0_kind_real, tem1)
123  endif
124  enddo
125  enddo
126 enddo
127 
128 
129 ! Cloud ice water effective radius
130 ! ---------------------------------
131 do k = 1,npz
132  do j = jsc,jec
133  do i = isc,iec
134 
135  if (seamask(i,j)) then
136 
137  tem2 = t(i,j,k) - tice
138  tem1 = grav/rdry
139  tem3 = tem1 * qi_ade(i,j,k) * (p(i,j,k)/delp(i,j,k))/t(i,j,k) * (1.0_kind_real + zvir * max(q(i,j,k),0.0_kind_real))
140 
141  if (tem2 < -50.0_kind_real ) then
142  qi_efr(i,j,k) = (1250._kind_real/9.917_kind_real)*tem3**0.109_kind_real
143  elseif (tem2 < -40.0_kind_real ) then
144  qi_efr(i,j,k) = (1250._kind_real/9.337_kind_real)*tem3**0.08_kind_real
145  elseif (tem2 < -30.0_kind_real ) then
146  qi_efr(i,j,k) = (1250._kind_real/9.208_kind_real)*tem3**0.055_kind_real
147  else
148  qi_efr(i,j,k) = (1250._kind_real/9.387_kind_real)*tem3**0.031_kind_real
149  endif
150 
151  endif
152 
153  enddo
154  enddo
155 enddo
156 
157 
158 ql_ade = max(0.0_kind_real,ql_ade)
159 qi_ade = max(0.0_kind_real,qi_ade)
160 ql_efr = max(0.0_kind_real,ql_efr)
161 qi_efr = max(0.0_kind_real,qi_efr)
162 
163 deallocate(seamask)
164 
165 end subroutine crtm_ade_efr
166 
167 !----------------------------------------------------------------------------
168 ! Compute mixing ratio from specific humidity -------------------------------
169 !----------------------------------------------------------------------------
170 
171 subroutine crtm_mixratio(geom,q,qmr)
172 
173 implicit none
174 
175 !Arguments
176 type(fv3jedi_geom) , intent(in ) :: geom
177 real(kind=kind_real), intent(in ) :: q(geom%isc:geom%iec,geom%jsc:geom%jec, 1:geom%npz) !Specific humidity | kg/kg
178 real(kind=kind_real), intent(out) :: qmr(geom%isc:geom%iec,geom%jsc:geom%jec, 1:geom%npz) !Mixing ratio | 1
179 
180 !Locals
181 integer :: isc,iec,jsc,jec,npz
182 integer :: i,j,k
183 real(kind=kind_real) :: q_pos(geom%isc:geom%iec,geom%jsc:geom%jec, 1:geom%npz)
184 real(kind=kind_real) :: c3(geom%isc:geom%iec,geom%jsc:geom%jec, 1:geom%npz)
185 
186 
187 ! Grid convenience
188 ! ----------------
189 isc = geom%isc
190 iec = geom%iec
191 jsc = geom%jsc
192 jec = geom%jec
193 npz = geom%npz
194 
195 
196 ! Remove negative values
197 ! ----------------------
198 q_pos = q
199 do k = 1,npz
200  do j = jsc,jec
201  do i = isc,iec
202  if (q_pos(i,j,k) < 0.0_kind_real) then
203  q_pos(i,j,k) = 0.0_kind_real
204  endif
205  enddo
206  enddo
207 enddo
208 
209 
210 ! Convert specific humidity
211 ! -------------------------
212 c3 = 1.0_kind_real / (1.0_kind_real - q_pos)
213 qmr = 1000.0_kind_real * q_pos * c3
214 
215 end subroutine crtm_mixratio
216 
217 !----------------------------------------------------------------------------
218 
219 subroutine crtm_mixratio_tl(geom,q,q_tl,qmr_tl)
220 
221 implicit none
222 
223 !Arguments
224 type(fv3jedi_geom) , intent(in ) :: geom
225 real(kind=kind_real), intent(in ) :: q(geom%isc:geom%iec,geom%jsc:geom%jec, 1:geom%npz) !Specific humidity | kg/kg
226 real(kind=kind_real), intent(in ) :: q_tl(geom%isc:geom%iec,geom%jsc:geom%jec, 1:geom%npz) !Specific humidity | kg/kg
227 real(kind=kind_real), intent(out) :: qmr_tl(geom%isc:geom%iec,geom%jsc:geom%jec, 1:geom%npz) !Mixing ratio | 1
228 
229 !Locals
230 integer :: isc,iec,jsc,jec,npz
231 integer :: i,j,k
232 real(kind=kind_real) :: q_pos(geom%isc:geom%iec,geom%jsc:geom%jec, 1:geom%npz)
233 real(kind=kind_real) :: q_tl_pos(geom%isc:geom%iec,geom%jsc:geom%jec, 1:geom%npz)
234 real(kind=kind_real) :: c3(geom%isc:geom%iec,geom%jsc:geom%jec, 1:geom%npz)
235 real(kind=kind_real) :: c3_tl(geom%isc:geom%iec,geom%jsc:geom%jec, 1:geom%npz)
236 
237 
238 ! Grid convenience
239 ! ----------------
240 isc = geom%isc
241 iec = geom%iec
242 jsc = geom%jsc
243 jec = geom%jec
244 npz = geom%npz
245 
246 
247 ! Remove negative values
248 ! ----------------------
249 q_pos = q
250 q_tl_pos = q_tl
251 do k = 1,npz
252  do j = jsc,jec
253  do i = isc,iec
254  if (q_pos(i,j,k) < 0.0_kind_real) then
255  q_pos(i,j,k) = 0.0_kind_real
256  q_tl_pos(i,j,k) = 0.0_kind_real
257  endif
258  enddo
259  enddo
260 enddo
261 
262 
263 ! Convert specific humidity
264 ! -------------------------
265 c3 = 1.0_kind_real / (1.0_kind_real - q_pos)
266 c3_tl = -((-q_tl_pos)/(1.0_kind_real-q_pos)**2)
267 qmr_tl = 1000.0_kind_real*(q_tl_pos*c3+q_pos*c3_tl)
268 
269 end subroutine crtm_mixratio_tl
270 
271 !----------------------------------------------------------------------------
272 
273 subroutine crtm_mixratio_ad(geom,q,q_ad,qmr_ad)
274 
275 implicit none
276 
277 !Arguments
278 type(fv3jedi_geom) , intent(in ) :: geom
279 real(kind=kind_real), intent(in ) :: q(geom%isc:geom%iec,geom%jsc:geom%jec, 1:geom%npz) !Specific humidity | kg/kg
280 real(kind=kind_real), intent(inout) :: q_ad(geom%isc:geom%iec,geom%jsc:geom%jec, 1:geom%npz) !Specific humidity | kg/kg
281 real(kind=kind_real), intent(inout) :: qmr_ad(geom%isc:geom%iec,geom%jsc:geom%jec, 1:geom%npz) !Mixing ratio | 1
282 
283 !Locals
284 integer :: isc,iec,jsc,jec,npz
285 integer :: i,j,k
286 real(kind=kind_real) :: q_pos(geom%isc:geom%iec,geom%jsc:geom%jec, 1:geom%npz)
287 real(kind=kind_real) :: q_ad_pos(geom%isc:geom%iec,geom%jsc:geom%jec, 1:geom%npz)
288 real(kind=kind_real) :: c3(geom%isc:geom%iec,geom%jsc:geom%jec, 1:geom%npz)
289 real(kind=kind_real) :: c3_ad(geom%isc:geom%iec,geom%jsc:geom%jec, 1:geom%npz)
290 
291 
292 ! Grid convenience
293 ! ----------------
294 isc = geom%isc
295 iec = geom%iec
296 jsc = geom%jsc
297 jec = geom%jec
298 npz = geom%npz
299 
300 
301 ! Remove negative values
302 ! ----------------------
303 q_pos = q
304 do k = 1,npz
305  do j = jsc,jec
306  do i = isc,iec
307  if (q_pos(i,j,k) < 0.0_kind_real) then
308  q_pos(i,j,k) = 0.0_kind_real
309  endif
310  enddo
311  enddo
312 enddo
313 
314 
315 ! Convert specific humidity
316 ! -------------------------
317 c3 = 1.0_kind_real/(1.0_kind_real-q_pos)
318 c3_ad = 1000.0_kind_real*q_pos*qmr_ad
319 q_ad_pos = c3_ad/(1.0_kind_real-q_pos)**2 + 1000.0_kind_real*c3*qmr_ad
320 qmr_ad = 0.0_8
321 
322 
323 ! Remove negative values adjoint
324 ! ------------------------------
325 do k = 1,npz
326  do j = jsc,jec
327  do i = isc,iec
328  if (q_pos(i,j,k) < 0.0_kind_real) then
329  q_ad_pos(i,j,k) = 0.0_kind_real
330  endif
331  enddo
332  enddo
333 enddo
334 
335 q_ad = q_ad + q_ad_pos
336 
337 end subroutine crtm_mixratio_ad
338 
339 
340 !----------------------------------------------------------------------------
341 ! Relative to specific humidity ---------------------------------------------
342 !----------------------------------------------------------------------------
343 
344 subroutine rh_to_q(geom,qsat,rh,q)
345 
346  implicit none
347  type(fv3jedi_geom), intent(in) :: geom
348  real(kind=kind_real), intent(in) :: qsat(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
349  real(kind=kind_real), intent(inout) :: q(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
350  real(kind=kind_real), intent(inout) :: rh(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
351 
352  q = rh * qsat
353 
354 end subroutine rh_to_q
355 
356 !----------------------------------------------------------------------------
357 
358 subroutine rh_to_q_tl(geom,qsat,rh,q)
359 
360  implicit none
361  type(fv3jedi_geom), intent(in) :: geom
362  real(kind=kind_real), intent(in) :: qsat(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
363  real(kind=kind_real), intent(inout) :: q(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
364  real(kind=kind_real), intent(in) :: rh(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
365 
366  q = rh * qsat
367 
368 end subroutine rh_to_q_tl
369 
370 !----------------------------------------------------------------------------
371 
372 subroutine rh_to_q_ad(geom,qsat,rh,q)
373 
374  implicit none
375  type(fv3jedi_geom), intent(in) :: geom
376  real(kind=kind_real), intent(in) :: qsat(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
377  real(kind=kind_real), intent(in) :: q(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
378  real(kind=kind_real), intent(inout) :: rh(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
379 
380  rh = rh + q * qsat
381 
382 end subroutine rh_to_q_ad
383 
384 !----------------------------------------------------------------------------
385 
386 subroutine q4_to_q2(geom,qils,qicn,qlls,qlcn,qi,ql,qilsf,qicnf)
387 
388  implicit none
389  type(fv3jedi_geom), intent(in) :: geom
390  real(kind=kind_real), intent(in) :: qils(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
391  real(kind=kind_real), intent(in) :: qicn(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
392  real(kind=kind_real), intent(in) :: qlls(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
393  real(kind=kind_real), intent(in) :: qlcn(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
394  real(kind=kind_real), intent(out) :: ql(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
395  real(kind=kind_real), intent(out) :: qi(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
396  real(kind=kind_real), optional, intent(out) :: qilsf(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
397  real(kind=kind_real), optional, intent(out) :: qicnf(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
398 
399  ! Convert cloud for large scale and convective regions into totals
400 
401  qi = qils + qicn
402  ql = qlls + qlcn
403 
404  ! Record fractions contained in each
405  if (present(qilsf)) then
406  qilsf = 0.0_kind_real
407  where (qi > 0.0_kind_real) qilsf = qils / qi
408  endif
409  if (present(qicnf)) then
410  qicnf = 0.0_kind_real
411  where (qi > 0.0_kind_real) qicnf = qicn / qi
412  endif
413 
414 end subroutine q4_to_q2
415 
416 !----------------------------------------------------------------------------
417 
418 subroutine q2_to_q4(geom,qi,ql,qilsf,qicnf,qils,qicn,qlls,qlcn)
419 
420  implicit none
421  type(fv3jedi_geom), intent(in) :: geom
422  real(kind=kind_real), intent(in) :: ql(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
423  real(kind=kind_real), intent(in) :: qi(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
424  real(kind=kind_real), intent(in) :: qilsf(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
425  real(kind=kind_real), intent(in) :: qicnf(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
426  real(kind=kind_real), intent(out) :: qils(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
427  real(kind=kind_real), intent(out) :: qicn(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
428  real(kind=kind_real), intent(out) :: qlls(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
429  real(kind=kind_real), intent(out) :: qlcn(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
430 
431  ! Convert totals into large scale and convective (requies knowledge of some desired fractions)
432  qils = qi * qilsf
433  qlls = qi * (1.0_kind_real - qilsf)
434  qicn = qi * qicnf
435  qlcn = qi * (1.0_kind_real - qicnf)
436 
437 end subroutine q2_to_q4
438 
439 !----------------------------------------------------------------------------
440 ! Specific to relative humidity ---------------------------------------------
441 !----------------------------------------------------------------------------
442 
443 subroutine q_to_rh(geom,qsat,q,rh)
444 
445  implicit none
446  type(fv3jedi_geom), intent(in) :: geom
447  real(kind=kind_real), intent(in) :: qsat(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
448  real(kind=kind_real), intent(in) :: q(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
449  real(kind=kind_real), intent(inout) :: rh(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
450 
451  rh = q / qsat
452 
453 end subroutine q_to_rh
454 
455 !----------------------------------------------------------------------------
456 
457 subroutine q_to_rh_tl(geom,qsat,q,rh)
458 
459  implicit none
460  type(fv3jedi_geom), intent(in) :: geom
461  real(kind=kind_real), intent(in) :: qsat(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
462  real(kind=kind_real), intent(in) :: q(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
463  real(kind=kind_real), intent(inout) :: rh(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
464 
465  rh = q / qsat
466 
467 end subroutine q_to_rh_tl
468 
469 !----------------------------------------------------------------------------
470 
471 subroutine q_to_rh_ad(geom,qsat,q,rh)
472 
473  implicit none
474  type(fv3jedi_geom), intent(in) :: geom
475  real(kind=kind_real), intent(in) :: qsat(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
476  real(kind=kind_real), intent(inout) :: q(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
477  real(kind=kind_real), intent(in) :: rh(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
478 
479  q = q + rh / qsat
480 
481 end subroutine q_to_rh_ad
482 
483 !----------------------------------------------------------------------------
484 
485 subroutine get_qsat(geom,delp,t,q,qsat)
486 
487  implicit none
488  type(fv3jedi_geom), intent(in) :: geom
489  real(kind=kind_real), intent(in) :: delp(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
490  real(kind=kind_real), intent(in) :: t(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
491  real(kind=kind_real), intent(in) :: q(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
492  real(kind=kind_real), intent(out) :: qsat(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz)
493 
494  integer :: i,j
495  real(kind=kind_real), allocatable :: pe(:,:,:)
496  real(kind=kind_real), allocatable :: pm(:,:,:)
497 
498  allocate(pe(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz+1))
499  allocate(pm(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz ))
500 
501  call delp_to_pe_p_logp(geom,delp,pe,pm)
502 
503  do j = geom%jsc,geom%jec
504  do i = geom%isc,geom%iec
505  call qsmith(geom%npz,t(i,j,:),q(i,j,:),pm(i,j,:),qsat(i,j,:))
506  enddo
507  enddo
508 
509  deallocate(pm)
510  deallocate(pe)
511 
512 end subroutine get_qsat
513 
514 !----------------------------------------------------------------------------
515 
516 subroutine qsmith(nlev,t,sphum,pl,qs)
517 
518  implicit none
519  integer, intent(in) :: nlev
520  real(kind_real), intent(in) :: t(nlev)
521  real(kind_real), intent(in) :: sphum(nlev)
522  real(kind_real), intent(in) :: pl(nlev)
523  real(kind_real), intent(out) :: qs(nlev)
524 
525  real(kind_real), parameter :: esl = 0.621971831
526  real(kind_real), allocatable :: table(:), des(:)
527 
528  real(kind_real) :: es
529  real(kind_real) :: ap1, eps10
530  integer :: k, it
531  integer, parameter :: length=2621
532 
533  eps10 = 10.0_kind_real*esl
534 
535  allocate ( table(length) )
536  call qs_table(length,table)
537 
538  allocate ( des(length) )
539  do k=1,length-1
540  des(k) = table(k+1) - table(k)
541  enddo
542  des(length) = des(length-1)
543 
544  do k=1,nlev
545  ap1 = 10.0_kind_real*dim(t(k), tice-160.0_kind_real) + 1.0_kind_real
546  ap1 = min(2621.0_kind_real, ap1)
547  it = int(ap1)
548  es = table(it) + (ap1-it)*des(it)
549  qs(k) = esl*es*(1.0_kind_real+zvir*sphum(k))/(pl(k))
550  enddo
551 
552  deallocate(table,des)
553 
554 end subroutine qsmith
555 
556 !----------------------------------------------------------------------------
557 
558 subroutine qs_table(n,table)
559 
560  implicit none
561 
562  integer, intent(in) :: n
563  real(kind=kind_real), intent(inout) :: table(n)
564 
565  real(kind=kind_real) :: tem, aa, b, c, d, e
566  integer :: i
567  real(kind=kind_real), parameter :: dt=0.1_kind_real
568  real(kind=kind_real), parameter :: esbasw = 1013246.0_kind_real
569  real(kind=kind_real), parameter :: tbasw = 373.16_kind_real
570  real(kind=kind_real), parameter :: tbasi = 273.16_kind_real
571  real(kind=kind_real), parameter :: tmin = tbasi - 160.0_kind_real
572 
573  ! Compute es over water, see smithsonian meteorological tables page 350.
574  do i=1,n
575  tem = tmin+dt*real(i-1)
576  aa = -7.90298_kind_real*(tbasw/tem-1)
577  b = 5.02808_kind_real*log10(tbasw/tem)
578  c = -1.3816e-07_kind_real*(10.0_kind_real**((1.0_kind_real-tem/tbasw)*11.344_kind_real)-1.0_kind_real)
579  d = 8.1328e-03_kind_real*(10.0_kind_real**((tbasw/tem-1.0_kind_real)*(-3.49149_kind_real))-1.0_kind_real)
580  e = log10(esbasw)
581  table(i) = 0.1_kind_real*10.0_kind_real**(aa+b+c+d+e)
582  enddo
583 
584 end subroutine qs_table
585 
586 !----------------------------------------------------------------------------
587 
588 end module moisture_vt_mod
moisture_vt_mod::qs_table
subroutine qs_table(n, table)
Definition: moisture_variables_mod.f90:559
moisture_vt_mod
Definition: moisture_variables_mod.f90:6
moisture_vt_mod::q_to_rh_ad
subroutine, public q_to_rh_ad(geom, qsat, q, rh)
Definition: moisture_variables_mod.f90:472
fv3jedi_constants_mod::tice
real(kind=kind_real), parameter, public tice
Definition: fv3jedi_constants_mod.f90:49
moisture_vt_mod::rh_to_q
subroutine, public rh_to_q(geom, qsat, rh, q)
Definition: moisture_variables_mod.f90:345
moisture_vt_mod::crtm_mixratio_tl
subroutine, public crtm_mixratio_tl(geom, q, q_tl, qmr_tl)
Definition: moisture_variables_mod.f90:220
moisture_vt_mod::q2_to_q4
subroutine, public q2_to_q4(geom, qi, ql, qilsf, qicnf, qils, qicn, qlls, qlcn)
Definition: moisture_variables_mod.f90:419
fv3jedi_geom_mod
Fortran module handling geometry for the FV3 model.
Definition: fv3jedi_geom_mod.f90:8
moisture_vt_mod::qsmith
subroutine qsmith(nlev, t, sphum, pl, qs)
Definition: moisture_variables_mod.f90:517
moisture_vt_mod::crtm_ade_efr
subroutine, public crtm_ade_efr(geom, p, T, delp, sea_frac, q, ql, qi, ql_ade, qi_ade, ql_efr, qi_efr)
Definition: moisture_variables_mod.f90:38
fv3jedi_geom_mod::fv3jedi_geom
Fortran derived type to hold geometry data for the FV3JEDI model.
Definition: fv3jedi_geom_mod.f90:46
moisture_vt_mod::q_to_rh
subroutine, public q_to_rh(geom, qsat, q, rh)
Definition: moisture_variables_mod.f90:444
fv3jedi_constants_mod
Definition: fv3jedi_constants_mod.f90:6
moisture_vt_mod::crtm_mixratio
subroutine, public crtm_mixratio(geom, q, qmr)
Definition: moisture_variables_mod.f90:172
moisture_vt_mod::crtm_mixratio_ad
subroutine, public crtm_mixratio_ad(geom, q, q_ad, qmr_ad)
Definition: moisture_variables_mod.f90:274
pressure_vt_mod
Definition: pressure_variables_mod.f90:6
moisture_vt_mod::get_qsat
subroutine, public get_qsat(geom, delp, t, q, qsat)
Definition: moisture_variables_mod.f90:486
moisture_vt_mod::q_to_rh_tl
subroutine, public q_to_rh_tl(geom, qsat, q, rh)
Definition: moisture_variables_mod.f90:458
fv3jedi_constants_mod::zvir
real(kind=kind_real), parameter, public zvir
Definition: fv3jedi_constants_mod.f90:41
moisture_vt_mod::rh_to_q_ad
subroutine, public rh_to_q_ad(geom, qsat, rh, q)
Definition: moisture_variables_mod.f90:373
fv3jedi_kinds_mod::kind_real
integer, parameter, public kind_real
Definition: fv3jedi_kinds_mod.f90:14
fv3jedi_constants_mod::rdry
real(kind=kind_real), parameter, public rdry
Definition: fv3jedi_constants_mod.f90:28
moisture_vt_mod::q4_to_q2
subroutine, public q4_to_q2(geom, qils, qicn, qlls, qlcn, qi, ql, qilsf, qicnf)
Definition: moisture_variables_mod.f90:387
pressure_vt_mod::delp_to_pe_p_logp
subroutine, public delp_to_pe_p_logp(geom, delp, pe, p, logp)
Definition: pressure_variables_mod.f90:33
fv3jedi_kinds_mod
Definition: fv3jedi_kinds_mod.f90:6
moisture_vt_mod::rh_to_q_tl
subroutine, public rh_to_q_tl(geom, qsat, rh, q)
Definition: moisture_variables_mod.f90:359
fv3jedi_constants_mod::grav
real(kind=kind_real), parameter, public grav
Definition: fv3jedi_constants_mod.f90:17