FMS 2025.01-dev
Flexible Modeling System
Loading...
Searching...
No Matches
sat_vapor_pres_k.inc
1!***********************************************************************
2!* GNU Lesser General Public License
3!*
4!* This file is part of the GFDL Flexible Modeling System (FMS).
5!*
6!* FMS is free software: you can redistribute it and/or modify it under
7!* the terms of the GNU Lesser General Public License as published by
8!* the Free Software Foundation, either version 3 of the License, or (at
9!* your option) any later version.
10!*
11!* FMS is distributed in the hope that it will be useful, but WITHOUT
12!* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13!* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14!* for more details.
15!*
16!* You should have received a copy of the GNU Lesser General Public
17!* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18!***********************************************************************
19!> @addtogroup sat_vapor_pres_k_mod
20!> @brief Kernel module to be used by @ref sat_vapor_pres_mod for
21!> @{
22
23!! table lookups and calculations
24
25!> This routine has been generalized to return tables for any temperature range and resolution
26!! The TABLEs for saturation vapor pressure are computed with r8_kind precision since
27!! these TABLES are module level variables that are decared as r8_kind.
28!! The subroutines compute_es_k, compute_es_k, compute_es_liq_k, and compute_es_liq_ice_k
29!! seem to be mostly used to compute the TABLE values (and thus all variables within them can be declared
30!! as r8_kind). However, these compute* subroutines have been modified to accept both r4_kind and r8_kind arguments
31!! for general usage and the math can be conducted in either r4_kind and r8_kind.
32!! In this initialization routine, r8_kind arguments are passed to these compute* subroutines.
33!! This routine does not assume the passed in arguments are in r8_precision.
34!! Thus all variables used for the computation of the TABLES (e.g. TABLE, DTABLE*, D2TABLE*) are promoted
35!! to r8_kind precision. All local variables used for computation are in r8_kind precision
36!! Thus the TABLEs are constructed as accurately as possible and are promoted down to r4_kind when users
37!! pass in r4_kind arguments to the LOOKUP* subroutines.
38 subroutine sat_vapor_pres_init_k_(table_size, tcmin, tcmax, TFREEZE, HLV, RVGAS, ES0, err_msg, &
39 use_exact_qs_input, do_simple, &
40 construct_table_wrt_liq, &
41 construct_table_wrt_liq_and_ice, &
42 teps, tmin, dtinv)
43
44 integer, intent(in) :: table_size
45 real(kind=fms_svp_kind_), intent(in) :: tcmin !< TABLE(1) = sat vapor pressure at temperature tcmin (deg C)
46 real(kind=fms_svp_kind_), intent(in) :: tcmax !< TABLE(table_size) = sat vapor pressure at temperature tcmax (deg C)
47 real(kind=fms_svp_kind_), intent(in) :: tfreeze !< Conversion to Kelvin
48 real(kind=fms_svp_kind_), intent(in) :: hlv !< Latent heat of evaporation [J/kg]
49 real(kind=fms_svp_kind_), intent(in) :: rvgas !< Gas constant for water vapor
50 real(kind=fms_svp_kind_), intent(in) :: es0 !< Humidity factor [dimensionless]
51 logical, intent(in) :: use_exact_qs_input
52 logical, intent(in) :: do_simple
53 logical, intent(in) :: construct_table_wrt_liq
54 logical, intent(in) :: construct_table_wrt_liq_and_ice
55 character(len=*), intent(out) :: err_msg
56 real(kind=fms_svp_kind_), intent(out), optional :: teps
57 real(kind=fms_svp_kind_), intent(out), optional :: tmin
58 real(kind=fms_svp_kind_), intent(out), optional :: dtinv
59
60
61!> increment used to generate derivative table
62!! the following variables are used in the computation
63!! of the *TABLES* (which is defined in r8_kind in sat_vapor_pres_mod)
64!! Thus these variables are declared with r8_kind
65 real(kind=r8_kind), dimension(3) :: tem(3), es(3)
66 real(kind=r8_kind) :: hdtinv, tinrc, tfact
67 integer :: i
68
69 err_msg = ''
70
71 if (module_is_initialized) return
72
73 if(allocated(table) .or. allocated(dtable) .or. allocated(d2table)) then
74 err_msg = 'Attempt to allocate sat vapor pressure tables when already allocated'
75 return
76 else
77 allocate(table(table_size), dtable(table_size), d2table(table_size))
78 endif
79
80 if (construct_table_wrt_liq) then
81 if(allocated(table2) .or. allocated(dtable2) .or. allocated(d2table2)) then
82 err_msg = 'Attempt to allocate sat vapor pressure table2s when already allocated'
83 return
84 else
85 allocate(table2(table_size), dtable2(table_size), d2table2(table_size))
86 endif
87 endif
88
89 if (construct_table_wrt_liq_and_ice) then
90 if(allocated(table3) .or. allocated(dtable3) .or. allocated(d2table3)) then
91 err_msg = 'Attempt to allocate sat vapor pressure table2s when already allocated'
92 return
93 else
94 allocate(table3(table_size), dtable3(table_size), d2table3(table_size))
95 endif
96 endif
97
98 table_siz = table_size
99 dtres = (real(tcmax,r8_kind)-real(tcmin,r8_kind))/real(table_size-1,r8_kind)
100 tminl = real(tcmin,r8_kind)+real(tfreeze,r8_kind) ! minimum valid temp in table
101 dtinvl = 1.0_r8_kind/dtres
102 tepsl = 0.5_r8_kind*dtres
103 tinrc = 0.1_r8_kind*dtres
104 if(present(teps )) teps =real(tepsl, fms_svp_kind_)
105 if(present(tmin )) tmin =real(tminl, fms_svp_kind_)
106 if(present(dtinv)) dtinv=real(dtinvl, fms_svp_kind_)
107
108!> To be able to compute tables for any temperature range and resolution,
109!! and at the same time exactly reproduce answers from memphis revision,
110!! it is necessary to compute ftact differently than it is in memphis.
111 tfact = 5.0_r8_kind*dtinvl
112
113 hdtinv = 0._r8_kind*dtinvl
114
115!> compute es tables from tcmin to tcmax
116!> estimate es derivative with small +/- difference
117
118 if (do_simple) then
119
120 !> TABLE = 610.78ES0*exp(-HLV/RGAS[1/tem - 1.TFREEZE])
121 !> DTABLE = HLV(TABLE/RVGAS)^2
122 do i = 1, table_size
123 tem(1) = tminl + dtres*real(i-1,r8_kind)
124 table(i) = real(es0,r8_kind)*610.78_r8_kind* &
125 exp( -real(hlv,r8_kind)/real(rvgas,r8_kind) * (1.0_r8_kind/tem(1) - 1._r8_kind/real(tfreeze,r8_kind)) )
126 dtable(i) = real(hlv,r8_kind)*table(i)/real(rvgas,r8_kind)/tem(1)**2._r8_kind
127 enddo
128
129 else
130
131 do i = 1, table_size
132 tem(1) = tminl + dtres*real(i-1,r8_kind)
133 tem(2) = tem(1)-tinrc
134 tem(3) = tem(1)+tinrc
135 es = compute_es_k(tem, real(tfreeze,r8_kind))
136 table(i) = es(1)
137 dtable(i) = (es(3)-es(2))*tfact
138 enddo
139
140 endif !if (do_simple)
141
142!> compute one-half second derivative using centered differences
143!! differencing des values in the table
144
145 do i = 2, table_size-1
146 d2table(i) = 0.25_r8_kind*dtinvl*(dtable(i+1)-dtable(i-1))
147 enddo
148!> one-sided derivatives at boundaries
149
150 d2table(1) = 0.50_r8_kind*dtinvl*(dtable(2)-dtable(1))
151
152 d2table(table_size) = 0.50_r8_kind*dtinvl*(dtable(table_size)-dtable(table_size-1))
153
154 if (construct_table_wrt_liq) then
155!> compute es tables from tcmin to tcmax
156!> estimate es derivative with small +/- difference
157
158 do i = 1, table_size
159 tem(1) = tminl + dtres*real(i-1,r8_kind)
160 tem(2) = tem(1)-tinrc
161 tem(3) = tem(1)+tinrc
162!> pass in flag to force all values to be wrt liquid
163 es = compute_es_liq_k(tem, real(tfreeze,r8_kind))
164 table2(i) = es(1)
165 dtable2(i) = (es(3)-es(2))*tfact
166 enddo
167
168!> compute one-half second derivative using centered differences
169!! differencing des values in the table
170
171 do i = 2, table_size-1
172 d2table2(i) = 0.25_r8_kind*dtinvl*(dtable2(i+1)-dtable2(i-1))
173 enddo
174!> one-sided derivatives at boundaries
175
176 d2table2(1) = 0.50_r8_kind*dtinvl*(dtable2(2)-dtable2(1))
177
178 d2table2(table_size) = 0.50_r8_kind*dtinvl*(dtable2(table_size)-dtable2(table_size-1))
179 endif
180
181
182 if (construct_table_wrt_liq_and_ice) then
183!> compute es tables from tcmin to tcmax
184!> estimate es derivative with small +/- difference
185
186 do i = 1, table_size
187 tem(1) = tminl + dtres*real(i-1,r8_kind)
188 tem(2) = tem(1)-tinrc
189 tem(3) = tem(1)+tinrc
190!> pass in flag to force all values to be wrt liquid
191 es = compute_es_liq_ice_k(tem, real(tfreeze,r8_kind))
192 table3(i) = es(1)
193 dtable3(i) = (es(3)-es(2))*tfact
194 enddo
195
196!> compute one-half second derivative using centered differences
197!! differencing des values in the table
198
199 do i = 2, table_size-1
200 d2table3(i) = 0.25_r8_kind*dtinvl*(dtable3(i+1)-dtable3(i-1))
201 enddo
202!> one-sided derivatives at boundaries
203
204 d2table3(1) = 0.50_r8_kind*dtinvl*(dtable3(2)-dtable3(1))
205
206 d2table3(table_size) = 0.50_r8_kind*dtinvl*(dtable3(table_size)-dtable3(table_size-1))
207 endif
208
209 use_exact_qs = use_exact_qs_input
210 module_is_initialized = .true.
211
212 end subroutine sat_vapor_pres_init_k_
213
214!#######################################################################
215
216 function compute_es_k_(tem, TFREEZE) result (es)
217 real(kind=fms_svp_kind_), intent(in) :: tem(:) !< temperature
218 real(kind=fms_svp_kind_), intent(in) :: tfreeze !< conversion to Kelvin
219 real(kind=fms_svp_kind_) :: es(size(tem,1)) !< saturation vapor pressure
220
221 real(kind=fms_svp_kind_) :: x
222 real(kind=fms_svp_kind_) :: esice
223 real(kind=fms_svp_kind_) :: esh2o
224 real(kind=fms_svp_kind_) :: tbasw
225 real(kind=fms_svp_kind_) :: tbasi
226 integer :: i
227
228 integer, parameter :: lkind=fms_svp_kind_ !< local kind parameter
229
230 real(kind=fms_svp_kind_), parameter :: esbasw = 101324.60_lkind
231 real(kind=fms_svp_kind_), parameter :: esbasi = 610.71_lkind
232
233 !> one and ten are declared for code readability. For example, one is easier to read
234 !! then 1.0_lkind where lkind=FMS_SVP_KIND_ throughout the code
235 real(fms_svp_kind_), parameter :: one=1.0_lkind
236 real(fms_svp_kind_), parameter :: ten=10.0_lkind
237
238 tbasw = tfreeze+100.0_lkind !to Kelvin
239 tbasi = tfreeze
240
241 do i = 1, size(tem)
242
243!> compute es over ice
244
245 !> x = -9.09718(TBASI/tem-1) - 3.56654log(TBASI/tem) + 0.876793(1-tem/TBASI) + log(ESBASI)
246 !! the coded equation below is the commented equation above
247 if (tem(i) < tbasi) then
248 x = -9.09718_lkind*(tbasi/tem(i)-one) &
249 -3.56654_lkind*log10(tbasi/tem(i)) &
250 +0.876793_lkind*(one-tem(i)/tbasi) + log10(esbasi)
251 esice =ten**(x)
252 else
253 esice = 0.0_lkind
254 endif
255
256!> compute es over water greater than -20 c.
257!! values over 100 c may not be valid
258!! see smithsonian meteorological tables page 350.
259
260 !> x = -7.90298(TBASW/tem-1) + 5.02808log(TBASW/tem)
261 !! -1.3816d-07*10^[11.344(1-tem/TBASW)-1]
262 !! +8.1328d-03*10^[-3.49149(TBASW/tem-1)-1] + log(ESBASW)
263 !! the coded equation below is the commented equation above
264 if (tem(i) > -20.0_lkind+tbasi) then
265 x = -7.90298_lkind*(tbasw/tem(i)-one) &
266 +5.02808_lkind*log10(tbasw/tem(i)) &
267 -1.3816e-07_lkind*(ten**((one-tem(i)/tbasw)*11.344_lkind)-one) &
268 +8.1328e-03_lkind*(ten**((tbasw/tem(i)-one)*(-3.49149_lkind))-one) &
269 +log10(esbasw)
270 esh2o = ten**(x)
271 else
272 esh2o = 0.0_lkind
273 endif
274
275!> derive blended es over ice and supercooled water between -20c and 0c
276
277 !> es = 0.05*[esice*(TBASI-10)+esh2o*(tem-TBASI+20)]
278 !! the coded equation below is the commented equation above
279 if (tem(i) <= -20.0_lkind+tbasi) then
280 es(i) = esice
281 else if (tem(i) >= tbasi) then
282 es(i) = esh2o
283 else
284 es(i) = 0.05_lkind*((tbasi-tem(i))*esice + (tem(i)-tbasi+20.0_lkind)*esh2o)
285 endif
286
287 enddo
288
289 end function compute_es_k_
290
291!#######################################################################
292
293 function compute_es_liq_k_(tem, TFREEZE) result (es)
294 real(kind=fms_svp_kind_), intent(in) :: tem(:) !< temperature
295 real(kind=fms_svp_kind_), intent(in) :: tfreeze !< conversion to Kelvin
296 real(kind=fms_svp_kind_) :: es(size(tem,1)) !< saturation vapor pressure
297
298 real(kind=fms_svp_kind_) :: x
299 real(kind=fms_svp_kind_) :: esh2o
300 real(kind=fms_svp_kind_) :: tbasw
301 integer :: i
302
303 !> local kind variable
304 !! one and ten are declared for code readability. For example, one is easier to read
305 !! then 1.0_lkind where lkind=FMS_SVP_KIND_ throughout the code
306 integer, parameter :: lkind=fms_svp_kind_
307 real(kind=fms_svp_kind_), parameter :: one=1.0_lkind
308 real(kind=fms_svp_kind_), parameter :: ten=10.0_lkind
309 real(kind=fms_svp_kind_), parameter :: esbasw = 101324.60_lkind
310
311 tbasw = tfreeze+100.0_lkind
312
313 do i = 1, size(tem)
314
315!> compute es over water for all temps.
316!! values over 100 c may not be valid
317!! see smithsonian meteorological tables page 350.
318
319 !> x = -7.90298(TBASW/tem-1) + 5.02808log(TBASW/tem)
320 !! -1.3816d-07*10^[11.344(1-tem/TBASW)-1]
321 !! +8.1328d-03*10^[-3.49149(TBASW/tem-1)-1] + log(ESBASW)
322 !! the coded equation below is the commented equation above
323 x = -7.90298_lkind*(tbasw/tem(i)-one) &
324 +5.02808_lkind*log10(tbasw/tem(i)) &
325 -1.3816e-07_lkind*(ten**((one-tem(i)/tbasw)*11.344_lkind)-one) &
326 +8.1328e-03_lkind*(ten**((tbasw/tem(i)-one)*(-3.49149_lkind))-one)&
327 +log10(esbasw)
328 esh2o = ten**(x)
329
330 es(i) = esh2o
331
332 enddo
333
334 end function compute_es_liq_k_
335
336!#######################################################################
337
338 function compute_es_liq_ice_k_(tem, TFREEZE) result (es)
339 real(kind=fms_svp_kind_), intent(in) :: tem(:) !< temperature
340 real(kind=fms_svp_kind_), intent(in) :: tfreeze !< conversion to Kelvin
341 real(kind=fms_svp_kind_) :: es(size(tem,1)) !< saturation vapor pressure
342
343 real(kind=fms_svp_kind_) :: x
344 real(kind=fms_svp_kind_) :: tbasw
345 real(kind=fms_svp_kind_) :: tbasi
346 integer :: i
347
348 integer, parameter :: lkind=fms_svp_kind_
349 real(kind=fms_svp_kind_), parameter :: esbasw = 101324.60_lkind
350 real(kind=fms_svp_kind_), parameter :: esbasi = 610.71_lkind
351 !> one and ten are declared for code readability. For example, one is easier to read
352 !! then 1.0_lkind where lkind=FMS_SVP_KIND_ throughout the code
353 real(kind=fms_svp_kind_), parameter :: one=1.0_lkind
354 real(kind=fms_svp_kind_), parameter :: ten=10.0_lkind
355
356 tbasw = tfreeze+100.0_lkind
357 tbasi = tfreeze
358
359 do i = 1, size(tem)
360
361 if (tem(i) < tbasi) then
362
363!> compute es over ice
364 !> x= -9.09718(TBASI/tem-1) -3.56654log(TBASI/tem) +0.87679(1-tem/TBASI)+log(EBASI)
365 !! the coded equation below is the commented equation above
366 x = -9.09718_lkind*(tbasi/tem(i)-one) &
367 -3.56654_lkind*log10(tbasi/tem(i)) &
368 +0.876793_lkind*(one-tem(i)/tbasi) + log10(esbasi)
369 es(i) =ten**(x)
370 else
371
372!> compute es over water
373!! values over 100 c may not be valid
374!! see smithsonian meteorological tables page 350.
375 !> x = -7.90298(TBASW/tem-1) + 5.02808log(TBASW/tem)
376 !! -1.3816d-07*10^[11.344(1-tem/TBASW)-1]
377 !! +8.1328d-03*10^[-3.49149(TBASW/tem-1)-1] + log(ESBASW)
378 !! the coded equation below is the commented equation above
379 x = -7.90298_lkind*(tbasw/tem(i)-one) &
380 +5.02808_lkind*log10(tbasw/tem(i)) &
381 -1.3816e-07_lkind*(ten**((one-tem(i)/tbasw)*11.344_lkind)-one) &
382 +8.1328e-03_lkind*(ten**((tbasw/tem(i)-one)*(-3.49149_lkind))-one) &
383 +log10(esbasw)
384 es(i) = ten**(x)
385 endif
386 enddo
387
388 end function compute_es_liq_ice_k_
389
390!#######################################################################
391
392 subroutine compute_qs_k_3d_ (temp, press, eps, zvir, qs, nbad, q, hc, &
393 dqsdT, esat, es_over_liq, es_over_liq_and_ice)
394
395 real(kind=fms_svp_kind_), intent(in), dimension(:,:,:) :: temp !< temperature
396 real(kind=fms_svp_kind_), intent(in), dimension(:,:,:) :: press !< pressure
397 real(kind=fms_svp_kind_), intent(in) :: eps !< EPSILO
398 real(kind=fms_svp_kind_), intent(in) :: zvir !< ZVIR
399 real(kind=fms_svp_kind_), intent(out), dimension(:,:,:) :: qs !< specific humidity
400 integer, intent(out) :: nbad !< if temperature is out of range
401 real(kind=fms_svp_kind_), intent(in), dimension(:,:,:), optional :: q !< vapor relative humidity
402 real(kind=fms_svp_kind_), intent(in), optional :: hc !< relative humidity
403 real(kind=fms_svp_kind_), intent(out), dimension(:,:,:), optional :: dqsdt !< d(qs)/dT
404 real(kind=fms_svp_kind_), intent(out), dimension(:,:,:), optional :: esat !< saturation vapor pressure
405 logical,intent(in), optional :: es_over_liq !< use es table wrt liquid only
406 logical,intent(in), optional :: es_over_liq_and_ice
407
408 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_, either r4_kind or r8_kind
409
410 real(kind=fms_svp_kind_), dimension(size(temp,1), size(temp,2), size(temp,3)) :: esloc
411 real(kind=fms_svp_kind_), dimension(size(temp,1), size(temp,2), size(temp,3)) :: desat
412 real(kind=fms_svp_kind_), dimension(size(temp,1), size(temp,2), size(temp,3)) :: denom
413 integer :: i, j, k
414 real(kind=fms_svp_kind_) :: hc_loc
415
416 if (present(hc)) then
417 hc_loc = hc
418 else
419 hc_loc = 1.0_lkind
420 endif
421 if (present(es_over_liq)) then
422 if (present (dqsdt)) then
423 call lookup_es2_des2_k (temp, esloc, desat, nbad)
424 desat = desat*hc_loc
425 else
426 call lookup_es2_k (temp, esloc, nbad)
427 endif
428 else if (present(es_over_liq_and_ice)) then
429 if (present (dqsdt)) then
430 call lookup_es3_des3_k (temp, esloc, desat, nbad)
431 desat = desat*hc_loc
432 else
433 call lookup_es3_k (temp, esloc, nbad)
434 endif
435 else
436 if (present (dqsdt)) then
437 call lookup_es_des_k (temp, esloc, desat, nbad)
438 desat = desat*hc_loc
439 else
440 call lookup_es_k (temp, esloc, nbad)
441 endif
442 endif
443 esloc = esloc*hc_loc
444 if (present (esat)) then
445 esat = esloc
446 endif
447 if (nbad == 0) then
448 if (present (q) .and. use_exact_qs) then
449 qs = (1.0_lkind + zvir*q)*eps*esloc/press
450 if (present (dqsdt)) then
451 dqsdt = (1.0_lkind + zvir*q)*eps*desat/press
452 endif
453 else ! (present(q))
454 denom = press - (1.0_lkind - eps)*esloc
455 do k=1,size(qs,3)
456 do j=1,size(qs,2)
457 do i=1,size(qs,1)
458 if (denom(i,j,k) > 0.0_lkind) then
459 qs(i,j,k) = eps*esloc(i,j,k)/denom(i,j,k)
460 else
461 qs(i,j,k) = eps
462 endif
463 end do
464 end do
465 end do
466 if (present (dqsdt)) then
467 dqsdt = eps*press*desat/denom**2
468 endif
469 endif ! (present(q))
470 else ! (nbad = 0)
471 qs = -999.0_lkind
472 if (present (dqsdt)) then
473 dqsdt = -999.0_lkind
474 endif
475 if (present (esat)) then
476 esat = -999.0_lkind
477 endif
478 endif ! (nbad = 0)
479
480
481 end subroutine compute_qs_k_3d_
482
483!#######################################################################
484
485 subroutine compute_qs_k_2d_ (temp, press, eps, zvir, qs, nbad, q, hc, &
486 dqsdT, esat, es_over_liq, es_over_liq_and_ice)
487
488 real(kind=fms_svp_kind_), intent(in), dimension(:,:) :: temp !< temperature
489 real(kind=fms_svp_kind_), intent(in), dimension(:,:) :: press !< pressure
490 real(kind=fms_svp_kind_), intent(in) :: eps !< EPSILO
491 real(kind=fms_svp_kind_), intent(in) :: zvir !< ZVIR
492 real(kind=fms_svp_kind_), intent(out), dimension(:,:) :: qs !< specific humidity
493 integer, intent(out) :: nbad !< if temperature is out of range
494 real(kind=fms_svp_kind_), intent(in), dimension(:,:), optional :: q !< vapor specific humidty
495 real(kind=fms_svp_kind_), intent(in), optional :: hc !< relative humidity
496 real(kind=fms_svp_kind_), intent(out), dimension(:,:), optional :: dqsdt !< d(qs)/dT
497 real(kind=fms_svp_kind_), intent(out), dimension(:,:), optional :: esat !<saturation vapor pressure
498 logical,intent(in), optional :: es_over_liq
499 logical,intent(in), optional :: es_over_liq_and_ice
500
501 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
502
503 real(kind=fms_svp_kind_), dimension(size(temp,1), size(temp,2)) :: esloc
504 real(kind=fms_svp_kind_), dimension(size(temp,1), size(temp,2)) :: desat
505 real(kind=fms_svp_kind_), dimension(size(temp,1), size(temp,2)) :: denom
506 integer :: i, j
507 real(kind=fms_svp_kind_) :: hc_loc
508
509 if (present(hc)) then
510 hc_loc = hc
511 else
512 hc_loc = 1.0_lkind
513 endif
514
515 if (present(es_over_liq)) then
516 if (present (dqsdt)) then
517 call lookup_es2_des2_k (temp, esloc, desat, nbad)
518 desat = desat*hc_loc
519 else
520 call lookup_es2_k (temp, esloc, nbad)
521 endif
522 else if (present(es_over_liq_and_ice)) then
523 if (present (dqsdt)) then
524 call lookup_es3_des3_k (temp, esloc, desat, nbad)
525 desat = desat*hc_loc
526 else
527 call lookup_es3_k (temp, esloc, nbad)
528 endif
529 else
530 if (present (dqsdt)) then
531 call lookup_es_des_k (temp, esloc, desat, nbad)
532 desat = desat*hc_loc
533 else
534 call lookup_es_k (temp, esloc, nbad)
535 endif
536 endif
537 esloc = esloc*hc_loc
538 if (present (esat)) then
539 esat = esloc
540 endif
541 if (nbad == 0) then
542 if (present (q) .and. use_exact_qs) then
543 qs = (1.0_lkind + zvir*q)*eps*esloc/press
544 if (present (dqsdt)) then
545 dqsdt = (1.0_lkind + zvir*q)*eps*desat/press
546 endif
547 else ! (present(q))
548 denom = press - (1.0_lkind - eps)*esloc
549 do j=1,size(qs,2)
550 do i=1,size(qs,1)
551 if (denom(i,j) > 0.0_lkind) then
552 qs(i,j) = eps*esloc(i,j)/denom(i,j)
553 else
554 qs(i,j) = eps
555 endif
556 end do
557 end do
558 if (present (dqsdt)) then
559 dqsdt = eps*press*desat/denom**2
560 endif
561 endif ! (present(q))
562 else ! (nbad = 0)
563 qs = -999.0_lkind
564 if (present (dqsdt)) then
565 dqsdt = -999.0_lkind
566 endif
567 if (present (esat)) then
568 esat = -999.0_lkind
569 endif
570 endif ! (nbad = 0)
571
572
573 end subroutine compute_qs_k_2d_
574
575!#######################################################################
576
577 subroutine compute_qs_k_1d_ (temp, press, eps, zvir, qs, nbad, q, hc, &
578 dqsdT, esat, es_over_liq, es_over_liq_and_ice)
579
580 real(kind=fms_svp_kind_), intent(in), dimension(:) :: temp !< temperature
581 real(kind=fms_svp_kind_), intent(in), dimension(:) :: press !< pressure
582 real(kind=fms_svp_kind_), intent(in) :: eps !< EPSILO
583 real(kind=fms_svp_kind_), intent(in) :: zvir !< ZVIR
584 real(kind=fms_svp_kind_), intent(out), dimension(:) :: qs !< specific humidity
585 integer, intent(out) :: nbad !< if temperature is out of range
586 real(kind=fms_svp_kind_), intent(in), dimension(:), optional :: q !< vapor specific humidity
587 real(kind=fms_svp_kind_), intent(in), optional :: hc !< relative humidity
588 real(kind=fms_svp_kind_), intent(out), dimension(:), optional :: dqsdt !< d(qs)/dt
589 real(kind=fms_svp_kind_), intent(out), dimension(:), optional :: esat !< saturation vapor pressure
590 logical,intent(in), optional :: es_over_liq
591 logical,intent(in), optional :: es_over_liq_and_ice
592
593 integer, parameter :: lkind=fms_svp_kind_
594
595 real(kind=fms_svp_kind_), dimension(size(temp,1)) :: esloc
596 real(kind=fms_svp_kind_), dimension(size(temp,1)) :: desat
597 real(kind=fms_svp_kind_), dimension(size(temp,1)) :: denom
598 integer :: i
599 real(kind=fms_svp_kind_) :: hc_loc
600
601 if (present(hc)) then
602 hc_loc = hc
603 else
604 hc_loc = 1.0_lkind
605 endif
606
607 if (present(es_over_liq)) then
608 if (present (dqsdt)) then
609 call lookup_es2_des2_k (temp, esloc, desat, nbad)
610 desat = desat*hc_loc
611 else
612 call lookup_es2_k (temp, esloc, nbad)
613 endif
614 else if (present(es_over_liq_and_ice)) then
615 if (present (dqsdt)) then
616 call lookup_es3_des3_k (temp, esloc, desat, nbad)
617 desat = desat*hc_loc
618 else
619 call lookup_es3_k (temp, esloc, nbad)
620 endif
621 else
622 if (present (dqsdt)) then
623 call lookup_es_des_k (temp, esloc, desat, nbad)
624 desat = desat*hc_loc
625 else
626 call lookup_es_k (temp, esloc, nbad)
627 endif
628 endif
629 esloc = esloc*hc_loc
630 if (present (esat)) then
631 esat = esloc
632 endif
633 if (nbad == 0) then
634 if (present (q) .and. use_exact_qs) then
635 qs = (1.0_lkind + zvir*q)*eps*esloc/press
636 if (present (dqsdt)) then
637 dqsdt = (1.0_lkind + zvir*q)*eps*desat/press
638 endif
639 else ! (present(q))
640 denom = press - (1.0_lkind - eps)*esloc
641 do i=1,size(qs,1)
642 if (denom(i) > 0.0_lkind) then
643 qs(i) = eps*esloc(i)/denom(i)
644 else
645 qs(i) = eps
646 endif
647 end do
648 if (present (dqsdt)) then
649 dqsdt = eps*press*desat/denom**2
650 endif
651 endif ! (present(q))
652 else ! (nbad = 0)
653 qs = -999.0_lkind
654 if (present (dqsdt)) then
655 dqsdt = -999.0_lkind
656 endif
657 if (present (esat)) then
658 esat = -999.0_lkind
659 endif
660 endif ! (nbad = 0)
661
662
663 end subroutine compute_qs_k_1d_
664
665!#######################################################################
666
667 subroutine compute_qs_k_0d_ (temp, press, eps, zvir, qs, nbad, q, hc, &
668 dqsdT, esat, es_over_liq, es_over_liq_and_ice)
669
670 real(kind=fms_svp_kind_), intent(in) :: temp !< temperature
671 real(kind=fms_svp_kind_), intent(in) :: press !< pressure
672 real(kind=fms_svp_kind_), intent(in) :: eps !< EPSILO
673 real(kind=fms_svp_kind_), intent(in) :: zvir !< ZVIR
674 real(kind=fms_svp_kind_), intent(out) :: qs !< specific humidity
675 integer, intent(out) :: nbad !< if temperature is out of range
676 real(kind=fms_svp_kind_), intent(in), optional :: q !< vapor specific humidity
677 real(kind=fms_svp_kind_), intent(in), optional :: hc !< relative humidity
678 real(kind=fms_svp_kind_), intent(out), optional :: dqsdt !< d(qs)/dT
679 real(kind=fms_svp_kind_), intent(out), optional :: esat !< saturation vapor pressure
680 logical,intent(in), optional :: es_over_liq
681 logical,intent(in), optional :: es_over_liq_and_ice
682
683 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
684
685 real(kind=fms_svp_kind_) :: esloc
686 real(kind=fms_svp_kind_) :: desat
687 real(kind=fms_svp_kind_) :: denom
688 real(kind=fms_svp_kind_) :: hc_loc
689
690 if (present(hc)) then
691 hc_loc = hc
692 else
693 hc_loc = 1.0_lkind
694 endif
695
696 if (present(es_over_liq)) then
697 if (present (dqsdt)) then
698 call lookup_es2_des2_k (temp, esloc, desat, nbad)
699 desat = desat*hc_loc
700 else
701 call lookup_es2_k (temp, esloc, nbad)
702 endif
703 else if (present(es_over_liq_and_ice)) then
704 if (present (dqsdt)) then
705 call lookup_es3_des3_k (temp, esloc, desat, nbad)
706 desat = desat*hc_loc
707 else
708 call lookup_es3_k (temp, esloc, nbad)
709 endif
710 else
711 if (present (dqsdt)) then
712 call lookup_es_des_k (temp, esloc, desat, nbad)
713 desat = desat*hc_loc
714 else
715 call lookup_es_k (temp, esloc, nbad)
716 endif
717 endif
718 esloc = esloc*hc_loc
719 if (present (esat)) then
720 esat = esloc
721 endif
722 if (nbad == 0) then
723 if (present (q) .and. use_exact_qs) then
724 qs = (1.0_lkind + zvir*q)*eps*esloc/press
725 if (present (dqsdt)) then
726 dqsdt = (1.0_lkind + zvir*q)*eps*desat/press
727 endif
728 else ! (present(q))
729 denom = press - (1.0_lkind - eps)*esloc
730 if (denom > 0.0_lkind) then
731 qs = eps*esloc/denom
732 else
733 qs = eps
734 endif
735 if (present (dqsdt)) then
736 dqsdt = eps*press*desat/denom**2
737 endif
738 endif ! (present(q))
739 else ! (nbad = 0)
740 qs = -999.0_lkind
741 if (present (dqsdt)) then
742 dqsdt = -999.0_lkind
743 endif
744 if (present (esat)) then
745 esat = -999.0_lkind
746 endif
747 endif ! (nbad = 0)
748
749
750 end subroutine compute_qs_k_0d_
751
752!#######################################################################
753
754 subroutine compute_mrs_k_3d_ (temp, press, eps, zvir, mrs, nbad, &
755 mr, hc, dmrsdT, esat,es_over_liq, es_over_liq_and_ice)
756
757 real(kind=fms_svp_kind_), intent(in), dimension(:,:,:) :: temp !< temperature
758 real(kind=fms_svp_kind_), intent(in), dimension(:,:,:) :: press !< pressure
759 real(kind=fms_svp_kind_), intent(in) :: eps !< EPSILO
760 real(kind=fms_svp_kind_), intent(in) :: zvir !< ZVIR
761 real(kind=fms_svp_kind_), intent(out), dimension(:,:,:) :: mrs !< mixing ratio at relative humidity
762 integer, intent(out) :: nbad !< if temperature is out of range
763 real(kind=fms_svp_kind_), intent(in), dimension(:,:,:), optional :: mr !< vapor mixing ratio
764 real(kind=fms_svp_kind_), intent(in), optional :: hc !< relative humidity
765 real(kind=fms_svp_kind_), intent(out), dimension(:,:,:), optional :: dmrsdt !< d(mrs)/dT
766 real(kind=fms_svp_kind_), intent(out), dimension(:,:,:), optional :: esat !< saturation vapor pressure
767 logical,intent(in), optional :: es_over_liq
768 logical,intent(in), optional :: es_over_liq_and_ice
769
770 real(FMS_SVP_KIND_), dimension(size(temp,1), size(temp,2), size(temp,3)) :: esloc
771 real(FMS_SVP_KIND_), dimension(size(temp,1), size(temp,2), size(temp,3)) :: desat
772 real(FMS_SVP_KIND_), dimension(size(temp,1), size(temp,2), size(temp,3)) :: denom
773
774 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
775
776 integer :: i, j, k
777 real(FMS_SVP_KIND_) :: hc_loc
778
779 if (present(hc)) then
780 hc_loc = hc
781 else
782 hc_loc = 1.0_lkind
783 endif
784
785 if (present (es_over_liq)) then
786 if (present (dmrsdt)) then
787 call lookup_es2_des2_k (temp, esloc, desat, nbad)
788 desat = desat*hc_loc
789 else
790 call lookup_es2_k (temp, esloc, nbad)
791 endif
792 else if (present(es_over_liq_and_ice)) then
793 if (present (dmrsdt)) then
794 call lookup_es3_des3_k (temp, esloc, desat, nbad)
795 desat = desat*hc_loc
796 else
797 call lookup_es3_k (temp, esloc, nbad)
798 endif
799 else
800 if (present (dmrsdt)) then
801 call lookup_es_des_k (temp, esloc, desat, nbad)
802 desat = desat*hc_loc
803 else
804 call lookup_es_k (temp, esloc, nbad)
805 endif
806 endif
807 esloc = esloc*hc_loc
808 if (present (esat)) then
809 esat = esloc
810 endif
811 if (nbad == 0) then
812 if (present (mr) .and. use_exact_qs) then
813 mrs = (eps + mr)*esloc/press
814 if (present (dmrsdt)) then
815 dmrsdt = (eps + mr)*desat/press
816 endif
817 else ! (present (mr))
818 denom = press - esloc
819 do k=1,size(mrs,3)
820 do j=1,size(mrs,2)
821 do i=1,size(mrs,1)
822 if (denom(i,j,k) > 0.0_lkind) then
823 mrs(i,j,k) = eps*esloc(i,j,k)/denom(i,j,k)
824 else
825 mrs(i,j,k) = eps
826 endif
827 end do
828 end do
829 end do
830 if (present (dmrsdt)) then
831 dmrsdt = eps*press*desat/denom**2
832 endif
833 endif !(present (mr))
834 else
835 mrs = -999.0_lkind
836 if (present (dmrsdt)) then
837 dmrsdt = -999.0_lkind
838 endif
839 if (present (esat)) then
840 esat = -999.0_lkind
841 endif
842 endif
843
844 end subroutine compute_mrs_k_3d_
845
846!#######################################################################
847
848 subroutine compute_mrs_k_2d_ (temp, press, eps, zvir, mrs, nbad, &
849 mr, hc, dmrsdT, esat,es_over_liq, es_over_liq_and_ice)
850
851 real(kind=fms_svp_kind_), intent(in), dimension(:,:) :: temp !< temperature
852 real(kind=fms_svp_kind_), intent(in), dimension(:,:) :: press !< pressure
853 real(kind=fms_svp_kind_), intent(in) :: eps !< EPSILO
854 real(kind=fms_svp_kind_), intent(in) :: zvir !< ZVIR
855 real(kind=fms_svp_kind_), intent(out), dimension(:,:) :: mrs !< mixing ratio at relative humidity
856 integer, intent(out) :: nbad !< if temperature is out of range
857 real(kind=fms_svp_kind_), intent(in), dimension(:,:), optional :: mr !< vapor mixing ratio
858 real(kind=fms_svp_kind_), intent(in), optional :: hc !< relative humidity
859 real(kind=fms_svp_kind_), intent(out), dimension(:,:), optional :: dmrsdt !< d(mrs)/dT
860 real(kind=fms_svp_kind_), intent(out), dimension(:,:), optional :: esat !< saturation vapor pressure
861 logical,intent(in), optional :: es_over_liq
862 logical,intent(in), optional :: es_over_liq_and_ice
863
864 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
865
866 real(kind=fms_svp_kind_), dimension(size(temp,1), size(temp,2)) :: esloc
867 real(kind=fms_svp_kind_), dimension(size(temp,1), size(temp,2)) :: desat
868 real(kind=fms_svp_kind_), dimension(size(temp,1), size(temp,2)) :: denom
869 integer :: i, j
870 real(kind=fms_svp_kind_) :: hc_loc
871
872 if (present(hc)) then
873 hc_loc = hc
874 else
875 hc_loc = 1.0_lkind
876 endif
877
878 if (present (es_over_liq)) then
879 if (present (dmrsdt)) then
880 call lookup_es2_des2_k (temp, esloc, desat, nbad)
881 desat = desat*hc_loc
882 else
883 call lookup_es2_k (temp, esloc, nbad)
884 endif
885 else if (present(es_over_liq_and_ice)) then
886 if (present (dmrsdt)) then
887 call lookup_es3_des3_k (temp, esloc, desat, nbad)
888 desat = desat*hc_loc
889 else
890 call lookup_es3_k (temp, esloc, nbad)
891 endif
892 else
893 if (present (dmrsdt)) then
894 call lookup_es_des_k (temp, esloc, desat, nbad)
895 desat = desat*hc_loc
896 else
897 call lookup_es_k (temp, esloc, nbad)
898 endif
899 endif
900 esloc = esloc*hc_loc
901 if (present (esat)) then
902 esat = esloc
903 endif
904 if (nbad == 0) then
905 if (present (mr) .and. use_exact_qs) then
906 mrs = (eps + mr)*esloc/press
907 if (present (dmrsdt)) then
908 dmrsdt = (eps + mr)*desat/press
909 endif
910 else ! (present (mr))
911 denom = press - esloc
912 do j=1,size(mrs,2)
913 do i=1,size(mrs,1)
914 if (denom(i,j) > 0.0_lkind) then
915 mrs(i,j) = eps*esloc(i,j)/denom(i,j)
916 else
917 mrs(i,j) = eps
918 endif
919 end do
920 end do
921 if (present (dmrsdt)) then
922 dmrsdt = eps*press*desat/denom**2
923 endif
924 endif !(present (mr))
925 else
926 mrs = -999.0_lkind
927 if (present (dmrsdt)) then
928 dmrsdt = -999.0_lkind
929 endif
930 if (present (esat)) then
931 esat = -999.0_lkind
932 endif
933 endif
934
935
936 end subroutine compute_mrs_k_2d_
937
938!#######################################################################
939
940 subroutine compute_mrs_k_1d_ (temp, press, eps, zvir, mrs, nbad, &
941 mr, hc, dmrsdT, esat,es_over_liq, es_over_liq_and_ice)
942
943 real(kind=fms_svp_kind_), intent(in), dimension(:) :: temp !< temperature
944 real(kind=fms_svp_kind_), intent(in), dimension(:) :: press !< pressure
945 real(kind=fms_svp_kind_), intent(in) :: eps !< EPSILO
946 real(kind=fms_svp_kind_), intent(in) :: zvir !< ZVIR
947 real(kind=fms_svp_kind_), intent(out), dimension(:) :: mrs !< mixing ratio at relative humidity
948 integer, intent(out) :: nbad !< if temperature is out of range
949 real(kind=fms_svp_kind_), intent(in), dimension(:), optional :: mr !< vapor mixing ratio
950 real(kind=fms_svp_kind_), intent(in), optional :: hc !< relative humidity
951 real(kind=fms_svp_kind_), intent(out), dimension(:), optional :: dmrsdt !< d(mrs)/dT
952 real(kind=fms_svp_kind_), intent(out), dimension(:), optional :: esat !< saturation vapor pressure
953 logical,intent(in), optional :: es_over_liq
954 logical,intent(in), optional :: es_over_liq_and_ice
955
956 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
957
958 real(kind=fms_svp_kind_), dimension(size(temp,1)) :: esloc
959 real(kind=fms_svp_kind_), dimension(size(temp,1)) :: desat
960 real(kind=fms_svp_kind_), dimension(size(temp,1)) :: denom
961 integer :: i
962 real(kind=fms_svp_kind_) :: hc_loc
963
964 if (present(hc)) then
965 hc_loc = hc
966 else
967 hc_loc = 1.0_lkind
968 endif
969
970 if (present (es_over_liq)) then
971 if (present (dmrsdt)) then
972 call lookup_es2_des2_k (temp, esloc, desat, nbad)
973 desat = desat*hc_loc
974 else
975 call lookup_es2_k (temp, esloc, nbad)
976 endif
977 else if (present(es_over_liq_and_ice)) then
978 if (present (dmrsdt)) then
979 call lookup_es3_des3_k (temp, esloc, desat, nbad)
980 desat = desat*hc_loc
981 else
982 call lookup_es3_k (temp, esloc, nbad)
983 endif
984 else
985 if (present (dmrsdt)) then
986 call lookup_es_des_k (temp, esloc, desat, nbad)
987 desat = desat*hc_loc
988 else
989 call lookup_es_k (temp, esloc, nbad)
990 endif
991 endif
992 esloc = esloc*hc_loc
993 if (present (esat)) then
994 esat = esloc
995 endif
996 if (nbad == 0) then
997 if (present (mr) .and. use_exact_qs) then
998 mrs = (eps + mr)*esloc/press
999 if (present (dmrsdt)) then
1000 dmrsdt = (eps + mr)*desat/press
1001 endif
1002 else ! (present (mr))
1003 denom = press - esloc
1004 do i=1,size(mrs,1)
1005 if (denom(i) > 0.0_lkind) then
1006 mrs(i) = eps*esloc(i)/denom(i)
1007 else
1008 mrs(i) = eps
1009 endif
1010 end do
1011 if (present (dmrsdt)) then
1012 dmrsdt = eps*press*desat/denom**2
1013 endif
1014 endif !(present (mr))
1015 else
1016 mrs = -999.0_lkind
1017 if (present (dmrsdt)) then
1018 dmrsdt = -999.0_lkind
1019 endif
1020 if (present (esat)) then
1021 esat = -999.0_lkind
1022 endif
1023 endif
1024
1025
1026 end subroutine compute_mrs_k_1d_
1027
1028!#######################################################################
1029
1030 subroutine compute_mrs_k_0d_ (temp, press, eps, zvir, mrs, nbad, &
1031 mr, hc, dmrsdT, esat,es_over_liq, es_over_liq_and_ice)
1032
1033 real(kind=fms_svp_kind_), intent(in) :: temp !< temperature
1034 real(kind=fms_svp_kind_), intent(in) :: press !< pressure
1035 real(kind=fms_svp_kind_), intent(in) :: eps !< EPSILO
1036 real(kind=fms_svp_kind_), intent(in) :: zvir !< ZVIR
1037 real(kind=fms_svp_kind_), intent(out) :: mrs !< mixing ratio at relative humidity
1038 integer, intent(out) :: nbad !< if temperature is out of range
1039 real(kind=fms_svp_kind_), intent(in), optional :: mr !< vapor mixing ratio
1040 real(kind=fms_svp_kind_), intent(in), optional :: hc !< relative humidity
1041 real(kind=fms_svp_kind_), intent(out), optional :: dmrsdt !< d(mrs)/dT
1042 real(kind=fms_svp_kind_), intent(out), optional :: esat !< saturation vapor pressure
1043 logical,intent(in), optional :: es_over_liq
1044 logical,intent(in), optional :: es_over_liq_and_ice
1045
1046 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
1047
1048 real(kind=fms_svp_kind_) :: esloc
1049 real(kind=fms_svp_kind_) :: desat
1050 real(kind=fms_svp_kind_) :: denom
1051 real(kind=fms_svp_kind_) :: hc_loc
1052
1053 if (present(hc)) then
1054 hc_loc = hc
1055 else
1056 hc_loc = 1.0_lkind
1057 endif
1058
1059 if (present (es_over_liq)) then
1060 if (present (dmrsdt)) then
1061 call lookup_es2_des2_k (temp, esloc, desat, nbad)
1062 desat = desat*hc_loc
1063 else
1064 call lookup_es2_k (temp, esloc, nbad)
1065 endif
1066 else if (present(es_over_liq_and_ice)) then
1067 if (present (dmrsdt)) then
1068 call lookup_es3_des3_k (temp, esloc, desat, nbad)
1069 desat = desat*hc_loc
1070 else
1071 call lookup_es3_k (temp, esloc, nbad)
1072 endif
1073 else
1074 if (present (dmrsdt)) then
1075 call lookup_es_des_k (temp, esloc, desat, nbad)
1076 desat = desat*hc_loc
1077 else
1078 call lookup_es_k (temp, esloc, nbad)
1079 endif
1080 endif
1081 esloc = esloc*hc_loc
1082 if (present (esat)) then
1083 esat = esloc
1084 endif
1085 if (nbad == 0) then
1086 if (present (mr) .and. use_exact_qs) then
1087 mrs = (eps + mr)*esloc/press
1088 if (present (dmrsdt)) then
1089 dmrsdt = (eps + mr)*desat/press
1090 endif
1091 else ! (present (mr))
1092 denom = press - esloc
1093 if (denom > 0.0_lkind) then
1094 mrs = eps*esloc/denom
1095 else
1096 mrs = eps
1097 endif
1098 if (present (dmrsdt)) then
1099 dmrsdt = eps*press*desat/denom**2
1100 endif
1101 endif !(present (mr))
1102 else
1103 mrs = -999.0_lkind
1104 if (present (dmrsdt)) then
1105 dmrsdt = -999.0_lkind
1106 endif
1107 if (present (esat)) then
1108 esat = -999.0_lkind
1109 endif
1110 endif
1111
1112
1113 end subroutine compute_mrs_k_0d_
1114
1115
1116!#######################################################################
1117
1118 subroutine lookup_es_des_k_3d_ (temp, esat, desat, nbad)
1119 real(kind=fms_svp_kind_), intent(in), dimension(:,:,:) :: temp !< temperature
1120 real(kind=fms_svp_kind_), intent(out), dimension(:,:,:) :: esat !<saturation vapor pressure
1121 real(kind=fms_svp_kind_), intent(out), dimension(:,:,:) :: desat !< derivative of saturation vapor pressure
1122 integer, intent(out) :: nbad !< if temperature is out of range
1123
1124 real(kind=fms_svp_kind_) :: tmp !< temp-TMIN
1125 real(kind=fms_svp_kind_) :: del !< delta T
1126 integer :: ind, i, j, k
1127 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
1128 !! for precision consistency and for code readability, the *ll variables are declared
1129 !! and used
1130 real(kind=fms_svp_kind_) :: dtresl
1131 real(kind=fms_svp_kind_) :: tepsll
1132 real(kind=fms_svp_kind_) :: tminll
1133 real(kind=fms_svp_kind_) :: dtinvll
1134 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
1135
1136 dtresl=real(dtres, fms_svp_kind_)
1137 tminll=real(tminl, fms_svp_kind_)
1138 tepsll=real(tepsl, fms_svp_kind_)
1139 dtinvll=real(dtinvl, fms_svp_kind_)
1140
1141 nbad = 0
1142 do k = 1, size(temp,3)
1143 do j = 1, size(temp,2)
1144 do i = 1, size(temp,1)
1145 tmp = temp(i,j,k)-tminll
1146 ind = int(dtinvll*(tmp+tepsll))
1147 if (ind < 0 .or. ind >= table_siz) then
1148 nbad = nbad+1
1149 else
1150 del = tmp-real(dtresl,fms_svp_kind_)*real(ind,fms_svp_kind_)
1151 !> esat = TABLE + del*(TABLE + del*D2TABLE)
1152 !! the coded equation below is the commented equation above
1153 esat(i,j,k) = real(table(ind+1),fms_svp_kind_) &
1154 + del*( real(dtable(ind+1),fms_svp_kind_) &
1155 + del*real(d2table(ind+1),fms_svp_kind_) )
1156 !> desat = DTABLE + 2del*D2TABLE
1157 !! the coded equation below is the commented equation above
1158 desat(i,j,k) = real(dtable(ind+1), fms_svp_kind_) &
1159 + 2.0_lkind*del*real(d2table(ind+1),fms_svp_kind_)
1160 endif
1161 enddo
1162 enddo
1163 enddo
1164
1165 end subroutine lookup_es_des_k_3d_
1166
1167!#######################################################################
1168
1169 subroutine lookup_es_des_k_2d_ (temp, esat, desat, nbad)
1170 real(kind=fms_svp_kind_), intent(in), dimension(:,:) :: temp !< temperature
1171 real(kind=fms_svp_kind_), intent(out), dimension(:,:) :: esat !< saturation vapor pressure
1172 real(kind=fms_svp_kind_), intent(out), dimension(:,:) :: desat !< derivative of the saturation vapor pressure
1173 integer, intent(out) :: nbad !< if temperature is out of range
1174 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
1175 real(kind=fms_svp_kind_) :: del !< delta T
1176 integer :: ind, i, j
1177 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
1178 !! for precision consistency and for code readability, the *ll variables are declared
1179 !! and used
1180 real(kind=fms_svp_kind_) :: dtresl
1181 real(kind=fms_svp_kind_) :: tepsll
1182 real(kind=fms_svp_kind_) :: tminll
1183 real(kind=fms_svp_kind_) :: dtinvll
1184 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
1185
1186 dtresl=real(dtres, fms_svp_kind_)
1187 tminll=real(tminl, fms_svp_kind_)
1188 tepsll=real(tepsl, fms_svp_kind_)
1189 dtinvll=real(dtinvl, fms_svp_kind_)
1190
1191
1192 nbad = 0
1193 do j = 1, size(temp,2)
1194 do i = 1, size(temp,1)
1195 tmp = temp(i,j)-tminll
1196 ind = int(dtinvll*(tmp+tepsll))
1197 if (ind < 0 .or. ind >= table_siz) then
1198 nbad = nbad+1
1199 else
1200 del = tmp-dtresl*real(ind,fms_svp_kind_)
1201 !> esat = TABLE + del*(TABLE + del*D2TABLE)
1202 !! the coded equation below is the commented equation above
1203 esat(i,j) = real(table(ind+1),fms_svp_kind_) &
1204 + del*( real(dtable(ind+1),fms_svp_kind_) &
1205 + del*real(d2table(ind+1),fms_svp_kind_) )
1206 !> desat = DTABLE + 2del*D2TABLE
1207 !! the coded equation below is the commented equation above
1208 desat(i,j) = real(dtable(ind+1),fms_svp_kind_) &
1209 + 2.0_lkind*del*real(d2table(ind+1),fms_svp_kind_)
1210 endif
1211 enddo
1212 enddo
1213
1214 end subroutine lookup_es_des_k_2d_
1215
1216!#######################################################################
1217
1218 subroutine lookup_es_des_k_1d_ (temp, esat, desat, nbad)
1219 real(kind=fms_svp_kind_), intent(in), dimension(:) :: temp !< temperature
1220 real(kind=fms_svp_kind_), intent(out), dimension(:) :: esat !< saturation vapor pressure
1221 real(kind=fms_svp_kind_), intent(out), dimension(:) :: desat !< derivative of the saturation vapor pressure
1222 integer, intent(out) :: nbad !< if temperature is out of range
1223
1224 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
1225 real(kind=fms_svp_kind_) :: del !< delta T
1226 integer :: ind, i
1227 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
1228 !! for precision consistency and for code readability, the *ll variables are declared
1229 !! and used
1230 real(kind=fms_svp_kind_) :: dtresl
1231 real(kind=fms_svp_kind_) :: tepsll
1232 real(kind=fms_svp_kind_) :: tminll
1233 real(kind=fms_svp_kind_) :: dtinvll
1234 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
1235
1236 dtresl=real(dtres, fms_svp_kind_)
1237 tminll=real(tminl, fms_svp_kind_)
1238 tepsll=real(tepsl, fms_svp_kind_)
1239 dtinvll=real(dtinvl, fms_svp_kind_)
1240
1241
1242 nbad = 0
1243 do i = 1, size(temp,1)
1244 tmp = temp(i)-tminll
1245 ind = int(dtinvll*(tmp+tepsll))
1246 if (ind < 0 .or. ind >= table_siz) then
1247 nbad = nbad+1
1248 else
1249 del = tmp-dtresl*real(ind,fms_svp_kind_)
1250 !> esat = TABLE + del*(TABLE + del*D2TABLE)
1251 !! the coded equation below is the commented equation above
1252 esat(i) = real(table(ind+1),fms_svp_kind_) &
1253 + del*( real(dtable(ind+1),fms_svp_kind_) &
1254 + del*real(d2table(ind+1),fms_svp_kind_) )
1255 !> desat = DTABLE + 2del*D2TABLE
1256 !! the coded equation below is the commented equation above
1257 desat(i) = real(dtable(ind+1),fms_svp_kind_) &
1258 + 2.0_lkind*del*real(d2table(ind+1),fms_svp_kind_)
1259 endif
1260 enddo
1261
1262 end subroutine lookup_es_des_k_1d_
1263
1264!#######################################################################
1265
1266 subroutine lookup_es_des_k_0d_ (temp, esat, desat, nbad)
1267 real(kind=fms_svp_kind_), intent(in) :: temp !< temperature
1268 real(kind=fms_svp_kind_), intent(out) :: esat !< saturation vapor pressure
1269 real(kind=fms_svp_kind_), intent(out) :: desat !< derivative of the saturation vapor pressure
1270 integer, intent(out) :: nbad !< if temperature is out of range
1271
1272 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
1273 real(kind=fms_svp_kind_) :: del !< delta T
1274 integer :: ind
1275 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
1276 !! for precision consistency and for code readability, the *ll variables are declared
1277 !! and used
1278 real(kind=fms_svp_kind_) :: dtresl
1279 real(kind=fms_svp_kind_) :: tepsll
1280 real(kind=fms_svp_kind_) :: tminll
1281 real(kind=fms_svp_kind_) :: dtinvll
1282 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
1283
1284 dtresl=real(dtres, fms_svp_kind_)
1285 tminll=real(tminl, fms_svp_kind_)
1286 tepsll=real(tepsl, fms_svp_kind_)
1287 dtinvll=real(dtinvl, fms_svp_kind_)
1288
1289 nbad = 0
1290 tmp = temp-tminll
1291 ind = int(dtinvll*(tmp+tepsll))
1292 if (ind < 0 .or. ind >= table_siz) then
1293 nbad = nbad+1
1294 else
1295 del = tmp-dtresl*real(ind,fms_svp_kind_)
1296 !> esat = TABLE + del*(TABLE + del*D2TABLE)
1297 !! the coded equation below is the commented equation above
1298 esat = real(table(ind+1),fms_svp_kind_) &
1299 + del*( real(dtable(ind+1),fms_svp_kind_) &
1300 + del*real(d2table(ind+1),fms_svp_kind_) )
1301 !> desat = DTABLE + 2del*D2TABLE
1302 !! the coded equation below is the commented equation above
1303 desat = real(dtable(ind+1),fms_svp_kind_) &
1304 + 2.0_lkind*del*real(d2table(ind+1),fms_svp_kind_)
1305 endif
1306
1307 end subroutine lookup_es_des_k_0d_
1308
1309!#######################################################################
1310
1311 subroutine lookup_es_k_3d_(temp, esat, nbad)
1312 real(kind=fms_svp_kind_), intent(in), dimension(:,:,:) :: temp !< temperature
1313 real(kind=fms_svp_kind_), intent(out), dimension(:,:,:) :: esat !< saturavation vapor pressure
1314 integer, intent(out) :: nbad !< if temperature is out of range
1315 real(kind=fms_svp_kind_) :: tmp !< temp - TMINLL
1316 real(kind=fms_svp_kind_) :: del !< delta T
1317 integer :: ind, i, j, k
1318 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
1319 !! for precision consistency and for code readability, the *ll variables are declared
1320 !! and used
1321 real(kind=fms_svp_kind_) :: dtresl
1322 real(kind=fms_svp_kind_) :: tepsll
1323 real(kind=fms_svp_kind_) :: tminll
1324 real(kind=fms_svp_kind_) :: dtinvll
1325
1326 dtresl=real(dtres, fms_svp_kind_)
1327 tminll=real(tminl, fms_svp_kind_)
1328 tepsll=real(tepsl, fms_svp_kind_)
1329 dtinvll=real(dtinvl, fms_svp_kind_)
1330
1331 nbad = 0
1332 do k = 1, size(temp,3)
1333 do j = 1, size(temp,2)
1334 do i = 1, size(temp,1)
1335 tmp = temp(i,j,k)-tminll
1336 ind = int(dtinvll*(tmp+tepsll))
1337 if (ind < 0 .or. ind >= table_siz) then
1338 nbad = nbad+1
1339 else
1340 del = tmp-dtresl*real(ind,fms_svp_kind_)
1341 !> esat = TABLE + del*(TABLE + del*D2TABLE)
1342 !! the coded equation below is the commented equation above
1343 esat(i,j,k) = real(table(ind+1),fms_svp_kind_) &
1344 + del*( real(dtable(ind+1),fms_svp_kind_) &
1345 + del*real(d2table(ind+1),fms_svp_kind_) )
1346 endif
1347 enddo
1348 enddo
1349 enddo
1350
1351 end subroutine lookup_es_k_3d_
1352
1353!#######################################################################
1354
1355 subroutine lookup_des_k_3d_(temp, desat, nbad)
1356 real(kind=fms_svp_kind_), intent(in), dimension(:,:,:) :: temp !< temperature
1357 real(kind=fms_svp_kind_), intent(out), dimension(:,:,:) :: desat !< derivative of sat vapor pressure
1358 integer, intent(out) :: nbad !< if temperature is out of range
1359 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
1360 real(kind=fms_svp_kind_) :: del !< delta T
1361 integer :: ind, i, j, k
1362 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
1363 !! for precision consistency and for code readability, the *ll variables are declared
1364 !! and used
1365 real(kind=fms_svp_kind_) :: dtresl
1366 real(kind=fms_svp_kind_) :: tepsll
1367 real(kind=fms_svp_kind_) :: tminll
1368 real(kind=fms_svp_kind_) :: dtinvll
1369 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
1370
1371 dtresl=real(dtres, fms_svp_kind_)
1372 tminll=real(tminl, fms_svp_kind_)
1373 tepsll=real(tepsl, fms_svp_kind_)
1374 dtinvll=real(dtinvl, fms_svp_kind_)
1375
1376 nbad = 0
1377 do k = 1, size(temp,3)
1378 do j = 1, size(temp,2)
1379 do i = 1, size(temp,1)
1380 tmp = temp(i,j,k)-tminll
1381 ind = int(dtinvll*(tmp+tepsll))
1382 if (ind < 0 .or. ind >= table_siz) then
1383 nbad = nbad+1
1384 else
1385 del = tmp-dtresl*real(ind,fms_svp_kind_)
1386 !> desat = DTABLE + 2del*D2TABLE
1387 !! the coded equation below is the commented equation above
1388 desat(i,j,k) = real(dtable(ind+1),fms_svp_kind_) &
1389 + 2.0_lkind*del*real(d2table(ind+1),fms_svp_kind_)
1390 endif
1391 enddo
1392 enddo
1393 enddo
1394
1395 end subroutine lookup_des_k_3d_
1396
1397!#######################################################################
1398 subroutine lookup_des_k_2d_(temp, desat, nbad)
1399 real(kind=fms_svp_kind_), intent(in), dimension(:,:) :: temp !< temperature
1400 real(kind=fms_svp_kind_), intent(out), dimension(:,:) :: desat !< derivative of sat vapor pressure
1401 integer, intent(out) :: nbad !< if temperature is out of range
1402 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
1403 real(kind=fms_svp_kind_) :: del !< delta T
1404 integer :: ind, i, j
1405
1406 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
1407 !! for precision consistency and for code readability, the *ll variables are declared
1408 !! and used
1409 real(kind=fms_svp_kind_) :: dtresl
1410 real(kind=fms_svp_kind_) :: tepsll
1411 real(kind=fms_svp_kind_) :: tminll
1412 real(kind=fms_svp_kind_) :: dtinvll
1413 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
1414
1415 dtresl=real(dtres, fms_svp_kind_)
1416 tminll=real(tminl, fms_svp_kind_)
1417 tepsll=real(tepsl, fms_svp_kind_)
1418 dtinvll=real(dtinvl, fms_svp_kind_)
1419
1420 nbad = 0
1421 do j = 1, size(temp,2)
1422 do i = 1, size(temp,1)
1423 tmp = temp(i,j)-tminll
1424 ind = int(dtinvll*(tmp+tepsll))
1425 if (ind < 0 .or. ind >= table_siz) then
1426 nbad = nbad+1
1427 else
1428 del = tmp-dtresl*real(ind,fms_svp_kind_)
1429 !> desat = DTABLE + 2del*D2TABLE
1430 !! the coded equation below is the commented equation above
1431 desat(i,j) = real(dtable(ind+1),fms_svp_kind_) &
1432 + 2.0_lkind*del*real(d2table(ind+1),fms_svp_kind_)
1433 endif
1434 enddo
1435 enddo
1436
1437 end subroutine lookup_des_k_2d_
1438!#######################################################################
1439 subroutine lookup_es_k_2d_(temp, esat, nbad)
1440 real(kind=fms_svp_kind_), intent(in), dimension(:,:) :: temp !< temperature
1441 real(kind=fms_svp_kind_), intent(out), dimension(:,:) :: esat !< saturation vapor pressure
1442 integer, intent(out) :: nbad !< if temperature is out of range
1443 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
1444 real(kind=fms_svp_kind_) :: del !< delta T
1445 integer :: ind, i, j
1446 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
1447 !! for precision consistency and for code readability, the *ll variables are declared
1448 !! and used
1449 real(kind=fms_svp_kind_) :: dtresl
1450 real(kind=fms_svp_kind_) :: tepsll
1451 real(kind=fms_svp_kind_) :: tminll
1452 real(kind=fms_svp_kind_) :: dtinvll
1453
1454 dtresl=real(dtres, fms_svp_kind_)
1455 tminll=real(tminl, fms_svp_kind_)
1456 tepsll=real(tepsl, fms_svp_kind_)
1457 dtinvll=real(dtinvl, fms_svp_kind_)
1458
1459
1460 nbad = 0
1461 do j = 1, size(temp,2)
1462 do i = 1, size(temp,1)
1463 tmp = temp(i,j)-tminll
1464 ind = int(dtinvll*(tmp+tepsll))
1465 if (ind < 0 .or. ind >= table_siz) then
1466 nbad = nbad+1
1467 else
1468 del = tmp-dtresl*real(ind,fms_svp_kind_)
1469 !> esat = TABLE + del*(TABLE + del*D2TABLE)
1470 !! the coded equation below is the commented equation above
1471 esat(i,j) = real(table(ind+1),fms_svp_kind_) &
1472 + del*(real(dtable(ind+1),fms_svp_kind_) &
1473 + del*real(d2table(ind+1),fms_svp_kind_) )
1474 endif
1475 enddo
1476 enddo
1477
1478 end subroutine lookup_es_k_2d_
1479!#######################################################################
1480 subroutine lookup_des_k_1d_(temp, desat, nbad)
1481 real(kind=fms_svp_kind_), intent(in), dimension(:) :: temp !< temperature
1482 real(kind=fms_svp_kind_), intent(out), dimension(:) :: desat !< derivative of sat vapor pressure
1483 integer, intent(out) :: nbad !< if temperature is out of range
1484 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
1485 real(kind=fms_svp_kind_) :: del !< delta T
1486 integer :: ind, i
1487 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
1488 !! for precision consistency and for code readability, the *ll variables are declared
1489 !! and used
1490 real(kind=fms_svp_kind_) :: dtresl
1491 real(kind=fms_svp_kind_) :: tepsll
1492 real(kind=fms_svp_kind_) :: tminll
1493 real(kind=fms_svp_kind_) :: dtinvll
1494 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
1495
1496 dtresl=real(dtres, fms_svp_kind_)
1497 tminll=real(tminl, fms_svp_kind_)
1498 tepsll=real(tepsl, fms_svp_kind_)
1499 dtinvll=real(dtinvl, fms_svp_kind_)
1500
1501 nbad = 0
1502 do i = 1, size(temp,1)
1503 tmp = temp(i)-tminll
1504 ind = int(dtinvll*(tmp+tepsll))
1505 if (ind < 0 .or. ind >= table_siz) then
1506 nbad = nbad+1
1507 else
1508 del = tmp-dtresl*real(ind,fms_svp_kind_)
1509 !> desat = DTABLE + 2del*D2TABLE
1510 !! the coded equation below is the commented equation above
1511 desat(i) = real(dtable(ind+1),fms_svp_kind_) &
1512 + 2.0_lkind*del*real(d2table(ind+1),fms_svp_kind_)
1513 endif
1514 enddo
1515
1516 end subroutine lookup_des_k_1d_
1517!#######################################################################
1518 subroutine lookup_es_k_1d_(temp, esat, nbad)
1519 real(kind=fms_svp_kind_), intent(in), dimension(:) :: temp !< temperature
1520 real(kind=fms_svp_kind_), intent(out), dimension(:) :: esat !< saturation vapor pressure
1521 integer, intent(out) :: nbad !< if temperature is out of range
1522 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
1523 real(kind=fms_svp_kind_) :: del !< delta T
1524 integer :: ind, i
1525 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
1526 !! for precision consistency and for code readability, the *ll variables are declared
1527 !! and used
1528 real(kind=fms_svp_kind_) :: dtresl
1529 real(kind=fms_svp_kind_) :: tepsll
1530 real(kind=fms_svp_kind_) :: tminll
1531 real(kind=fms_svp_kind_) :: dtinvll
1532
1533 dtresl=real(dtres, fms_svp_kind_)
1534 tminll=real(tminl, fms_svp_kind_)
1535 tepsll=real(tepsl, fms_svp_kind_)
1536 dtinvll=real(dtinvl, fms_svp_kind_)
1537
1538 nbad = 0
1539 do i = 1, size(temp,1)
1540 tmp = temp(i)-tminll
1541 ind = int(dtinvll*(tmp+tepsll))
1542 if (ind < 0 .or. ind >= table_siz) then
1543 nbad = nbad+1
1544 else
1545 del = tmp-dtresl*real(ind,fms_svp_kind_)
1546 !> esat = TABLE + del*(TABLE + del*D2TABLE)
1547 !! the coded equation below is the commented equation above
1548 esat(i) = real(table(ind+1),fms_svp_kind_) &
1549 + del*(real(dtable(ind+1),fms_svp_kind_) &
1550 + del*real(d2table(ind+1),fms_svp_kind_) )
1551 endif
1552 enddo
1553
1554 end subroutine lookup_es_k_1d_
1555!#######################################################################
1556 subroutine lookup_des_k_0d_(temp, desat, nbad)
1557 real(kind=fms_svp_kind_), intent(in) :: temp !< temperature
1558 real(kind=fms_svp_kind_), intent(out) :: desat !< derivative of sat vapor pressure
1559 integer, intent(out) :: nbad !< if temperature is out of range
1560 real(kind=fms_svp_kind_) :: tmp !< temp - TMINLL
1561 real(kind=fms_svp_kind_) :: del !< delta T
1562 integer :: ind
1563 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
1564 !! for precision consistency and for code readability, the *ll variables are declared
1565 !! and used
1566 real(kind=fms_svp_kind_) :: dtresl
1567 real(kind=fms_svp_kind_) :: tepsll
1568 real(kind=fms_svp_kind_) :: tminll
1569 real(kind=fms_svp_kind_) :: dtinvll
1570 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
1571
1572 dtresl=real(dtres, fms_svp_kind_)
1573 tminll=real(tminl, fms_svp_kind_)
1574 tepsll=real(tepsl, fms_svp_kind_)
1575 dtinvll=real(dtinvl, fms_svp_kind_)
1576
1577 nbad = 0
1578 tmp = temp-tminll
1579 ind = int(dtinvll*(tmp+tepsll))
1580 if (ind < 0 .or. ind >= table_siz) then
1581 nbad = nbad+1
1582 else
1583 del = tmp-dtresl*real(ind,fms_svp_kind_)
1584 !> desat = DTABLE + 2del*D2TABLE
1585 !! the coded equation below is the commented equation above
1586 desat = real(dtable(ind+1),fms_svp_kind_) &
1587 + 2.0_lkind*del*real(d2table(ind+1),fms_svp_kind_)
1588 endif
1589
1590 end subroutine lookup_des_k_0d_
1591!#######################################################################
1592 subroutine lookup_es_k_0d_(temp, esat, nbad)
1593 real(kind=fms_svp_kind_), intent(in) :: temp !< temperature
1594 real(kind=fms_svp_kind_), intent(out) :: esat !< saturation vapor pressure
1595 integer, intent(out) :: nbad !< if temperature is out of range
1596 real(kind=fms_svp_kind_) :: tmp !< temp - TMINLL
1597 real(kind=fms_svp_kind_) :: del !< delta T
1598 integer :: ind
1599 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
1600 !! for precision consistency and for code readability, the *ll variables are declared
1601 !! and used
1602 real(kind=fms_svp_kind_) :: dtresl
1603 real(kind=fms_svp_kind_) :: tepsll
1604 real(kind=fms_svp_kind_) :: tminll
1605 real(kind=fms_svp_kind_) :: dtinvll
1606
1607 dtresl=real(dtres, fms_svp_kind_)
1608 tminll=real(tminl, fms_svp_kind_)
1609 tepsll=real(tepsl, fms_svp_kind_)
1610 dtinvll=real(dtinvl, fms_svp_kind_)
1611
1612 nbad = 0
1613 tmp = temp-tminll
1614 ind = int(dtinvll*(tmp+tepsll))
1615 if (ind < 0 .or. ind >= table_siz) then
1616 nbad = nbad+1
1617 else
1618 del = tmp-dtresl*real(ind,fms_svp_kind_)
1619 !> esat = TABLE + del*(TABLE + del*D2TABLE)
1620 !! the coded equation below is the commented equation above
1621 esat = real(table(ind+1),fms_svp_kind_) &
1622 + del*( real(dtable(ind+1),fms_svp_kind_) &
1623 + del*real(d2table(ind+1),fms_svp_kind_) )
1624 endif
1625
1626 end subroutine lookup_es_k_0d_
1627!#######################################################################
1628
1629 subroutine lookup_es2_des2_k_3d_ (temp, esat, desat, nbad)
1630 real(kind=fms_svp_kind_), intent(in), dimension(:,:,:) :: temp !< temperature
1631 real(kind=fms_svp_kind_), intent(out), dimension(:,:,:) :: esat !< saturation vapor pressure
1632 real(kind=fms_svp_kind_), intent(out), dimension(:,:,:) :: desat !< derivative of sat vapor pressure
1633 integer, intent(out) :: nbad !< if temperature is out of range
1634
1635 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
1636 real(kind=fms_svp_kind_) :: del !< delta T
1637 integer :: ind, i, j, k
1638 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
1639 !! for precision consistency and for code readability, the *ll variables are declared
1640 !! and used
1641 real(kind=fms_svp_kind_) :: dtresl
1642 real(kind=fms_svp_kind_) :: tepsll
1643 real(kind=fms_svp_kind_) :: tminll
1644 real(kind=fms_svp_kind_) :: dtinvll
1645 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
1646
1647 dtresl=real(dtres, fms_svp_kind_)
1648 tminll=real(tminl, fms_svp_kind_)
1649 tepsll=real(tepsl, fms_svp_kind_)
1650 dtinvll=real(dtinvl, fms_svp_kind_)
1651
1652 nbad = 0
1653 do k = 1, size(temp,3)
1654 do j = 1, size(temp,2)
1655 do i = 1, size(temp,1)
1656 tmp = temp(i,j,k)-tminll
1657 ind = int(dtinvll*(tmp+tepsll))
1658 if (ind < 0 .or. ind >= table_siz) then
1659 nbad = nbad+1
1660 else
1661 del = tmp-dtresl*real(ind,fms_svp_kind_)
1662 !> esat = TABLE + del*(TABLE + del*D2TABLE)
1663 !! the coded equation below is the commented equation above
1664 esat(i,j,k) = real(table2(ind+1),fms_svp_kind_) &
1665 + del*( real(dtable2(ind+1),fms_svp_kind_) &
1666 + del*real(d2table2(ind+1),fms_svp_kind_) )
1667 !> desat = DTABLE + 2del*D2TABLE
1668 !! the coded equation below is the commented equation above
1669 desat(i,j,k) = real(dtable2(ind+1),fms_svp_kind_) &
1670 + 2.0_lkind*del*real(d2table2(ind+1),fms_svp_kind_)
1671 endif
1672 enddo
1673 enddo
1674 enddo
1675
1676 end subroutine lookup_es2_des2_k_3d_
1677
1678!#######################################################################
1679
1680 subroutine lookup_es2_des2_k_2d_ (temp, esat, desat, nbad)
1681 real(kind=fms_svp_kind_), intent(in), dimension(:,:) :: temp !< temperature
1682 real(kind=fms_svp_kind_), intent(out), dimension(:,:) :: esat !< saturation vapor pressure
1683 real(kind=fms_svp_kind_), intent(out), dimension(:,:) :: desat !< derivative of sat vapor pressure
1684 integer, intent(out) :: nbad !< if temperature is out of range
1685
1686 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
1687 real(kind=fms_svp_kind_) :: del !< delta T
1688 integer :: ind, i, j
1689 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
1690 !! for precision consistency and for code readability, the *ll variables are declared
1691 !! and used
1692 real(kind=fms_svp_kind_) :: dtresl
1693 real(kind=fms_svp_kind_) :: tepsll
1694 real(kind=fms_svp_kind_) :: tminll
1695 real(kind=fms_svp_kind_) :: dtinvll
1696 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
1697
1698 dtresl=real(dtres, fms_svp_kind_)
1699 tminll=real(tminl, fms_svp_kind_)
1700 tepsll=real(tepsl, fms_svp_kind_)
1701 dtinvll=real(dtinvl, fms_svp_kind_)
1702
1703 nbad = 0
1704 do j = 1, size(temp,2)
1705 do i = 1, size(temp,1)
1706 tmp = temp(i,j)-tminll
1707 ind = int(dtinvll*(tmp+tepsll))
1708 if (ind < 0 .or. ind >= table_siz) then
1709 nbad = nbad+1
1710 else
1711 del = tmp-dtresl*real(ind,fms_svp_kind_)
1712 !> esat = TABLE + del*(TABLE + del*D2TABLE)
1713 !! the coded equation below is the commented equation above
1714 esat(i,j) = real(table2(ind+1),fms_svp_kind_) &
1715 + del*( real(dtable2(ind+1),fms_svp_kind_) &
1716 + del*real(d2table2(ind+1),fms_svp_kind_) )
1717 !> desat = DTABLE + 2del*D2TABLE
1718 !! the coded equation below is the commented equation above
1719 desat(i,j) = real(dtable2(ind+1),fms_svp_kind_) &
1720 + 2.0_lkind*del*real(d2table2(ind+1),fms_svp_kind_)
1721 endif
1722 enddo
1723 enddo
1724
1725 end subroutine lookup_es2_des2_k_2d_
1726
1727!#######################################################################
1728
1729 subroutine lookup_es2_des2_k_1d_ (temp, esat, desat, nbad)
1730 real(kind=fms_svp_kind_), intent(in), dimension(:) :: temp !< temperature
1731 real(kind=fms_svp_kind_), intent(out), dimension(:) :: esat !< saturation vapor pressure
1732 real(kind=fms_svp_kind_), intent(out), dimension(:) :: desat !< derivative of sat vapor pressure
1733 integer, intent(out) :: nbad !< if temperature is out of range
1734
1735 real(kind=fms_svp_kind_) :: tmp !< temp - TMINLL
1736 real(kind=fms_svp_kind_) :: del !< delta T
1737 integer :: ind, i
1738 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
1739 !! for precision consistency and for code readability, the *ll variables are declared
1740 !! and used
1741 real(kind=fms_svp_kind_) :: dtresl
1742 real(kind=fms_svp_kind_) :: tepsll
1743 real(kind=fms_svp_kind_) :: tminll
1744 real(kind=fms_svp_kind_) :: dtinvll
1745 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
1746
1747 dtresl=real(dtres, fms_svp_kind_)
1748 tminll=real(tminl, fms_svp_kind_)
1749 tepsll=real(tepsl, fms_svp_kind_)
1750 dtinvll=real(dtinvl, fms_svp_kind_)
1751
1752 nbad = 0
1753 do i = 1, size(temp,1)
1754 tmp = temp(i)-tminll
1755 ind = int(dtinvll*(tmp+tepsll))
1756 if (ind < 0 .or. ind >= table_siz) then
1757 nbad = nbad+1
1758 else
1759 del = tmp-dtresl*real(ind,fms_svp_kind_)
1760 !> esat = TABLE + del*(TABLE + del*D2TABLE)
1761 !! the coded equation below is the commented equation above
1762 esat(i) = real(table2(ind+1),fms_svp_kind_) &
1763 + del*( real(dtable2(ind+1),fms_svp_kind_) &
1764 + del*real(d2table2(ind+1),fms_svp_kind_) )
1765 !> desat = DTABLE + 2del*D2TABLE
1766 desat(i) = real(dtable2(ind+1),fms_svp_kind_) &
1767 + 2.0_lkind*del*real(d2table2(ind+1),fms_svp_kind_)
1768 endif
1769 enddo
1770
1771 end subroutine lookup_es2_des2_k_1d_
1772
1773!#######################################################################
1774
1775 subroutine lookup_es2_des2_k_0d_ (temp, esat, desat, nbad)
1776 real(kind=fms_svp_kind_), intent(in) :: temp !< temperature
1777 real(kind=fms_svp_kind_), intent(out) :: esat !< saturation vapor pressure
1778 real(kind=fms_svp_kind_), intent(out) :: desat !< derivative of sat vapor pressure
1779 integer, intent(out) :: nbad !< if temperature is out of range
1780
1781 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
1782 real(kind=fms_svp_kind_) :: del !< delta T
1783 integer :: ind
1784 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
1785 !! for precision consistency and for code readability, the *ll variables are declared
1786 !! and used
1787 real(kind=fms_svp_kind_) :: dtresl
1788 real(kind=fms_svp_kind_) :: tepsll
1789 real(kind=fms_svp_kind_) :: tminll
1790 real(kind=fms_svp_kind_) :: dtinvll
1791 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
1792
1793 dtresl=real(dtres, fms_svp_kind_)
1794 tminll=real(tminl, fms_svp_kind_)
1795 tepsll=real(tepsl, fms_svp_kind_)
1796 dtinvll=real(dtinvl, fms_svp_kind_)
1797
1798
1799 nbad = 0
1800 tmp = temp-tminll
1801 ind = int(dtinvll*(tmp+tepsll))
1802 if (ind < 0 .or. ind >= table_siz) then
1803 nbad = nbad+1
1804 else
1805 del = tmp-dtresl*real(ind,fms_svp_kind_)
1806 !> esat = TABLE + del*(TABLE + del*D2TABLE)
1807 !! the coded equation below is the commented equation above
1808 esat = real(table2(ind+1),fms_svp_kind_) &
1809 + del*( real(dtable2(ind+1),fms_svp_kind_) &
1810 + del*real(d2table2(ind+1),fms_svp_kind_) )
1811 !> desat = DTABLE + 2del*D2TABLE
1812 !! the coded equation below is the commented equation above
1813 desat = real(dtable2(ind+1),fms_svp_kind_) &
1814 + 2.0_lkind*del*real(d2table2(ind+1),fms_svp_kind_)
1815 endif
1816
1817 end subroutine lookup_es2_des2_k_0d_
1818
1819!#######################################################################
1820
1821 subroutine lookup_es2_k_3d_(temp, esat, nbad)
1822 real(kind=fms_svp_kind_), intent(in), dimension(:,:,:) :: temp !< temperature
1823 real(kind=fms_svp_kind_), intent(out), dimension(:,:,:) :: esat !< saturation vapor pressure
1824 integer, intent(out) :: nbad !< if temperature is out of range
1825 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
1826 real(kind=fms_svp_kind_) :: del !< delta T
1827 integer :: ind, i, j, k
1828 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
1829 !! for precision consistency and for code readability, the *ll variables are declared
1830 !! and used
1831 real(kind=fms_svp_kind_) :: dtresl
1832 real(kind=fms_svp_kind_) :: tepsll
1833 real(kind=fms_svp_kind_) :: tminll
1834 real(kind=fms_svp_kind_) :: dtinvll
1835
1836 dtresl=real(dtres, fms_svp_kind_)
1837 tminll=real(tminl, fms_svp_kind_)
1838 tepsll=real(tepsl, fms_svp_kind_)
1839 dtinvll=real(dtinvl, fms_svp_kind_)
1840
1841 nbad = 0
1842 do k = 1, size(temp,3)
1843 do j = 1, size(temp,2)
1844 do i = 1, size(temp,1)
1845 tmp = temp(i,j,k)-tminll
1846 ind = int(dtinvll*(tmp+tepsll))
1847 if (ind < 0 .or. ind >= table_siz) then
1848 nbad = nbad+1
1849 else
1850 del = tmp-dtresl*real(ind,fms_svp_kind_)
1851 !> esat = TABLE + del*(TABLE + del*D2TABLE)
1852 !! the coded equation below is the commented equation above
1853 esat(i,j,k) = real(table2(ind+1),fms_svp_kind_) &
1854 + del*( real(dtable2(ind+1),fms_svp_kind_) &
1855 + del*real(d2table2(ind+1),fms_svp_kind_) )
1856 endif
1857 enddo
1858 enddo
1859 enddo
1860
1861 end subroutine lookup_es2_k_3d_
1862
1863!#######################################################################
1864
1865 subroutine lookup_des2_k_3d_(temp, desat, nbad)
1866 real(kind=fms_svp_kind_), intent(in), dimension(:,:,:) :: temp !< temperature
1867 real(kind=fms_svp_kind_), intent(out), dimension(:,:,:) :: desat !< derivative of sat vapor pressure
1868 integer, intent(out) :: nbad !< if temperature is out of range
1869 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
1870 real(kind=fms_svp_kind_) :: del !< delta T
1871 integer :: ind, i, j, k
1872 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
1873 !! for precision consistency and for code readability, the *ll variables are declared
1874 !! and used
1875 real(kind=fms_svp_kind_) :: dtresl
1876 real(kind=fms_svp_kind_) :: tepsll
1877 real(kind=fms_svp_kind_) :: tminll
1878 real(kind=fms_svp_kind_) :: dtinvll
1879 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
1880
1881 dtresl=real(dtres, fms_svp_kind_)
1882 tminll=real(tminl, fms_svp_kind_)
1883 tepsll=real(tepsl, fms_svp_kind_)
1884 dtinvll=real(dtinvl, fms_svp_kind_)
1885
1886 nbad = 0
1887 do k = 1, size(temp,3)
1888 do j = 1, size(temp,2)
1889 do i = 1, size(temp,1)
1890 tmp = temp(i,j,k)-tminll
1891 ind = int(dtinvll*(tmp+tepsll))
1892 if (ind < 0 .or. ind >= table_siz) then
1893 nbad = nbad+1
1894 else
1895 del = tmp-dtresl*real(ind,fms_svp_kind_)
1896 !> desat = DTABLE + 2del*D2TABLE
1897 !! the coded equation below is the commented equation above
1898 desat(i,j,k) = real(dtable2(ind+1),fms_svp_kind_) &
1899 + 2.0_lkind*del*real(d2table2(ind+1),fms_svp_kind_)
1900 endif
1901 enddo
1902 enddo
1903 enddo
1904
1905 end subroutine lookup_des2_k_3d_
1906
1907!#######################################################################
1908 subroutine lookup_des2_k_2d_(temp, desat, nbad)
1909 real(kind=fms_svp_kind_), intent(in), dimension(:,:) :: temp !< temperature
1910 real(kind=fms_svp_kind_), intent(out), dimension(:,:) :: desat !< derivative of sat vapor pressure
1911 integer, intent(out) :: nbad !< if temperature is out of range
1912 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
1913 real(kind=fms_svp_kind_) :: del !< delta T
1914 integer :: ind, i, j
1915 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
1916 !! for precision consistency and for code readability, the *ll variables are declared
1917 !! and used
1918 real(kind=fms_svp_kind_) :: dtresl
1919 real(kind=fms_svp_kind_) :: tepsll
1920 real(kind=fms_svp_kind_) :: tminll
1921 real(kind=fms_svp_kind_) :: dtinvll
1922 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
1923
1924 dtresl=real(dtres, fms_svp_kind_)
1925 tminll=real(tminl, fms_svp_kind_)
1926 tepsll=real(tepsl, fms_svp_kind_)
1927 dtinvll=real(dtinvl, fms_svp_kind_)
1928
1929 nbad = 0
1930 do j = 1, size(temp,2)
1931 do i = 1, size(temp,1)
1932 tmp = temp(i,j)-tminll
1933 ind = int(dtinvll*(tmp+tepsll))
1934 if (ind < 0 .or. ind >= table_siz) then
1935 nbad = nbad+1
1936 else
1937 del = tmp-dtresl*real(ind,fms_svp_kind_)
1938 !> desat = DTABLE + 2del*D2TABLE
1939 !! the coded equation below is the commented equation above
1940 desat(i,j) = real(dtable2(ind+1),fms_svp_kind_) &
1941 + 2.0_lkind*del*real(d2table2(ind+1),fms_svp_kind_)
1942 endif
1943 enddo
1944 enddo
1945
1946 end subroutine lookup_des2_k_2d_
1947!#######################################################################
1948 subroutine lookup_es2_k_2d_(temp, esat, nbad)
1949 real(kind=fms_svp_kind_), intent(in), dimension(:,:) :: temp !< temperature
1950 real(kind=fms_svp_kind_), intent(out), dimension(:,:) :: esat !< saturation vapor pressure
1951 integer, intent(out) :: nbad !< if temperature is out of range
1952 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
1953 real(kind=fms_svp_kind_) :: del !< delta T
1954 integer :: ind, i, j
1955 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
1956 !! for precision consistency and for code readability, the *ll variables are declared
1957 !! and used
1958 real(kind=fms_svp_kind_) :: dtresl
1959 real(kind=fms_svp_kind_) :: tepsll
1960 real(kind=fms_svp_kind_) :: tminll
1961 real(kind=fms_svp_kind_) :: dtinvll
1962
1963 dtresl=real(dtres, fms_svp_kind_)
1964 tminll=real(tminl, fms_svp_kind_)
1965 tepsll=real(tepsl, fms_svp_kind_)
1966 dtinvll=real(dtinvl, fms_svp_kind_)
1967
1968 nbad = 0
1969 do j = 1, size(temp,2)
1970 do i = 1, size(temp,1)
1971 tmp = temp(i,j)-tminll
1972 ind = int(dtinvll*(tmp+tepsll))
1973 if (ind < 0 .or. ind >= table_siz) then
1974 nbad = nbad+1
1975 else
1976 del = tmp-dtresl*real(ind,kind=fms_svp_kind_)
1977 !> esat = TABLE + del*(TABLE + del*D2TABLE)
1978 !! the coded equation below is the commented equation above
1979 esat(i,j) = real(table2(ind+1),fms_svp_kind_) &
1980 + del*( real(dtable2(ind+1),fms_svp_kind_) &
1981 + del*real(d2table2(ind+1),fms_svp_kind_) )
1982 endif
1983 enddo
1984 enddo
1985
1986 end subroutine lookup_es2_k_2d_
1987!#######################################################################
1988 subroutine lookup_des2_k_1d_(temp, desat, nbad)
1989 real(kind=fms_svp_kind_), intent(in), dimension(:) :: temp !< temperature
1990 real(kind=fms_svp_kind_), intent(out), dimension(:) :: desat !< derivative of sat vapor pressure
1991 integer, intent(out) :: nbad !< if temperature is out of range
1992 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
1993 real(kind=fms_svp_kind_) :: del !< delta T
1994 integer :: ind, i
1995 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
1996 !! for precision consistency and for code readability, the *ll variables are declared
1997 !! and used
1998 real(kind=fms_svp_kind_) :: dtresl
1999 real(kind=fms_svp_kind_) :: tepsll
2000 real(kind=fms_svp_kind_) :: tminll
2001 real(kind=fms_svp_kind_) :: dtinvll
2002 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
2003
2004 dtresl=real(dtres, fms_svp_kind_)
2005 tminll=real(tminl, fms_svp_kind_)
2006 tepsll=real(tepsl, fms_svp_kind_)
2007 dtinvll=real(dtinvl, fms_svp_kind_)
2008
2009 nbad = 0
2010 do i = 1, size(temp,1)
2011 tmp = temp(i)-tminll
2012 ind = int(dtinvll*(tmp+tepsll))
2013 if (ind < 0 .or. ind >= table_siz) then
2014 nbad = nbad+1
2015 else
2016 del = tmp-dtresl*real(ind,fms_svp_kind_)
2017 !> desat = DTABLE + 2del*D2TABLE
2018 !! the coded equation below is the commented equation above
2019 desat(i) = real(dtable2(ind+1),fms_svp_kind_) &
2020 + 2.0_lkind*del*real(d2table2(ind+1),fms_svp_kind_)
2021 endif
2022 enddo
2023
2024 end subroutine lookup_des2_k_1d_
2025!#######################################################################
2026 subroutine lookup_es2_k_1d_(temp, esat, nbad)
2027 real(kind=fms_svp_kind_), intent(in), dimension(:) :: temp !< temperature
2028 real(kind=fms_svp_kind_), intent(out), dimension(:) :: esat !< saturation vapor pressure
2029 integer, intent(out) :: nbad !< if temperature is out of range
2030 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
2031 real(kind=fms_svp_kind_) :: del !< delta T
2032 integer :: ind, i
2033 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
2034 !! for precision consistency and for code readability, the *ll variables are declared
2035 !! and used
2036 real(kind=fms_svp_kind_) :: dtresl
2037 real(kind=fms_svp_kind_) :: tepsll
2038 real(kind=fms_svp_kind_) :: tminll
2039 real(kind=fms_svp_kind_) :: dtinvll
2040
2041 dtresl=real(dtres, fms_svp_kind_)
2042 tminll=real(tminl, fms_svp_kind_)
2043 tepsll=real(tepsl, fms_svp_kind_)
2044 dtinvll=real(dtinvl, fms_svp_kind_)
2045
2046 nbad = 0
2047 do i = 1, size(temp,1)
2048 tmp = temp(i)-tminll
2049 ind = int(dtinvll*(tmp+tepsll))
2050 if (ind < 0 .or. ind >= table_siz) then
2051 nbad = nbad+1
2052 else
2053 del = tmp-dtresl*real(ind,fms_svp_kind_)
2054 !> esat = TABLE + del*(TABLE + del*D2TABLE)
2055 !! the coded equation below is the commented equation above
2056 esat(i) = real(table2(ind+1),fms_svp_kind_) &
2057 + del*( real(dtable2(ind+1),fms_svp_kind_) &
2058 + del*real(d2table2(ind+1),fms_svp_kind_) )
2059 endif
2060 enddo
2061
2062 end subroutine lookup_es2_k_1d_
2063!#######################################################################
2064 subroutine lookup_des2_k_0d_(temp, desat, nbad)
2065 real(kind=fms_svp_kind_), intent(in) :: temp !< temperature
2066 real(kind=fms_svp_kind_), intent(out) :: desat !< derivative of sat vapor pressure
2067 integer, intent(out) :: nbad !< if temperature is out of range
2068 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
2069 real(kind=fms_svp_kind_) :: del !< delta T
2070 integer :: ind
2071 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
2072 !! for precision consistency and for code readability, the *ll variables are declared
2073 !! and used
2074 real(kind=fms_svp_kind_) :: dtresl
2075 real(kind=fms_svp_kind_) :: tepsll
2076 real(kind=fms_svp_kind_) :: tminll
2077 real(kind=fms_svp_kind_) :: dtinvll
2078 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
2079
2080 dtresl=real(dtres, fms_svp_kind_)
2081 tminll=real(tminl, fms_svp_kind_)
2082 tepsll=real(tepsl, fms_svp_kind_)
2083 dtinvll=real(dtinvl, fms_svp_kind_)
2084
2085 nbad = 0
2086 tmp = temp-tminll
2087 ind = int(dtinvll*(tmp+tepsll))
2088 if (ind < 0 .or. ind >= table_siz) then
2089 nbad = nbad+1
2090 else
2091 del = tmp-dtresl*real(ind,fms_svp_kind_)
2092 !> desat = DTABLE + 2del*D2TABLE
2093 !! the coded equation below is the commented equation above
2094 desat = real(dtable2(ind+1),fms_svp_kind_) &
2095 + 2.0_lkind*del*real(d2table2(ind+1),fms_svp_kind_)
2096 endif
2097
2098 end subroutine lookup_des2_k_0d_
2099!#######################################################################
2100 subroutine lookup_es2_k_0d_(temp, esat, nbad)
2101 real(kind=fms_svp_kind_), intent(in) :: temp !< temperature
2102 real(kind=fms_svp_kind_), intent(out) :: esat !< saturation vapor pressure
2103 integer, intent(out) :: nbad !< if temperature is out of range
2104 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
2105 real(kind=fms_svp_kind_) :: del !< delta T
2106 integer :: ind
2107 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
2108 !! for precision consistency and for code readability, the *ll variables are declared
2109 !! and used
2110 real(kind=fms_svp_kind_) :: dtresl
2111 real(kind=fms_svp_kind_) :: tepsll
2112 real(kind=fms_svp_kind_) :: tminll
2113 real(kind=fms_svp_kind_) :: dtinvll
2114
2115 dtresl=real(dtres, fms_svp_kind_)
2116 tminll=real(tminl, fms_svp_kind_)
2117 tepsll=real(tepsl, fms_svp_kind_)
2118 dtinvll=real(dtinvl, fms_svp_kind_)
2119
2120 nbad = 0
2121 tmp = temp-tminll
2122 ind = int(dtinvll*(tmp+tepsll))
2123 if (ind < 0 .or. ind >= table_siz) then
2124 nbad = nbad+1
2125 else
2126 del = tmp-dtresl*real(ind,fms_svp_kind_)
2127 !> esat = TABLE + del*(TABLE + del*D2TABLE)
2128 !! the coded equation below is the commented equation above
2129 esat = real(table2(ind+1),fms_svp_kind_) &
2130 + del*( real(dtable2(ind+1),fms_svp_kind_) &
2131 + del*real(d2table2(ind+1),fms_svp_kind_) )
2132 endif
2133
2134 end subroutine lookup_es2_k_0d_
2135!#######################################################################
2136
2137!#######################################################################
2138
2139 subroutine lookup_es3_des3_k_3d_ (temp, esat, desat, nbad)
2140 real(kind=fms_svp_kind_), intent(in), dimension(:,:,:) :: temp !< temperature
2141 real(kind=fms_svp_kind_), intent(out), dimension(:,:,:) :: esat !< saturation vapor pressure
2142 real(kind=fms_svp_kind_), intent(out), dimension(:,:,:) :: desat !< derivative of esat
2143 integer, intent(out) :: nbad !< if temperature is out of range
2144
2145 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
2146 real(kind=fms_svp_kind_) :: del !< delta T
2147 integer :: ind, i, j, k
2148 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
2149 !! for precision consistency and for code readability, the *ll variables are declared
2150 !! and used
2151 real(kind=fms_svp_kind_) :: dtresl
2152 real(kind=fms_svp_kind_) :: tepsll
2153 real(kind=fms_svp_kind_) :: tminll
2154 real(kind=fms_svp_kind_) :: dtinvll
2155 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
2156
2157 dtresl=real(dtres, fms_svp_kind_)
2158 tminll=real(tminl, fms_svp_kind_)
2159 tepsll=real(tepsl, fms_svp_kind_)
2160 dtinvll=real(dtinvl, fms_svp_kind_)
2161
2162 nbad = 0
2163 do k = 1, size(temp,3)
2164 do j = 1, size(temp,2)
2165 do i = 1, size(temp,1)
2166 tmp = temp(i,j,k)-tminll
2167 ind = int(dtinvll*(tmp+tepsll))
2168 if (ind < 0 .or. ind >= table_siz) then
2169 nbad = nbad+1
2170 else
2171 del = tmp-dtresl*real(ind,fms_svp_kind_)
2172 !> esat = TABLE + del*(TABLE + del*D2TABLE)
2173 !! the coded equation below is the commented equation above
2174 esat(i,j,k) = real(table3(ind+1),fms_svp_kind_) &
2175 + del*( real(dtable3(ind+1),fms_svp_kind_) &
2176 + del*real(d2table3(ind+1),fms_svp_kind_) )
2177 !> desat = DTABLE + 2del*D2TABLE
2178 !! the coded equation below is the commented equation above
2179 desat(i,j,k) = real(dtable3(ind+1),fms_svp_kind_) &
2180 + 2.0_lkind*del*real(d2table3(ind+1),fms_svp_kind_)
2181 endif
2182 enddo
2183 enddo
2184 enddo
2185
2186 end subroutine lookup_es3_des3_k_3d_
2187
2188!#######################################################################
2189
2190 subroutine lookup_es3_des3_k_2d_ (temp, esat, desat, nbad)
2191 real(kind=fms_svp_kind_), intent(in), dimension(:,:) :: temp !< temperature
2192 real(kind=fms_svp_kind_), intent(out), dimension(:,:) :: esat !< saturation vapor pressure
2193 real(kind=fms_svp_kind_), intent(out), dimension(:,:) :: desat !< derivative of desat
2194 integer, intent(out) :: nbad !< if temperature is out of range
2195
2196 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
2197 real(kind=fms_svp_kind_) :: del !< delta T
2198 integer :: ind, i, j
2199 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
2200 !! for precision consistency and for code readability, the *ll variables are declared
2201 !! and used
2202 real(kind=fms_svp_kind_) :: dtresl
2203 real(kind=fms_svp_kind_) :: tepsll
2204 real(kind=fms_svp_kind_) :: tminll
2205 real(kind=fms_svp_kind_) :: dtinvll
2206 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
2207
2208 dtresl=real(dtres, fms_svp_kind_)
2209 tminll=real(tminl, fms_svp_kind_)
2210 tepsll=real(tepsl, fms_svp_kind_)
2211 dtinvll=real(dtinvl, fms_svp_kind_)
2212
2213 nbad = 0
2214 do j = 1, size(temp,2)
2215 do i = 1, size(temp,1)
2216 tmp = temp(i,j)-tminll
2217 ind = int(dtinvll*(tmp+tepsll))
2218 if (ind < 0 .or. ind >= table_siz) then
2219 nbad = nbad+1
2220 else
2221 del = tmp-dtresl*real(ind,fms_svp_kind_)
2222 !> esat = TABLE + del*(TABLE + del*D2TABLE)
2223 !! the coded equation below is the commented equation above
2224 esat(i,j) = real(table3(ind+1),fms_svp_kind_) &
2225 + del*( real(dtable3(ind+1),fms_svp_kind_) &
2226 + del*real(d2table3(ind+1),fms_svp_kind_) )
2227 !> desat = DTABLE + 2del*D2TABLE
2228 !! the coded equation below is the commented equation above
2229 desat(i,j) = real(dtable3(ind+1),fms_svp_kind_) &
2230 + 2.0_lkind*del*real(d2table3(ind+1),fms_svp_kind_)
2231 endif
2232 enddo
2233 enddo
2234
2235 end subroutine lookup_es3_des3_k_2d_
2236
2237!#######################################################################
2238
2239 subroutine lookup_es3_des3_k_1d_ (temp, esat, desat, nbad)
2240 real(kind=fms_svp_kind_), intent(in), dimension(:) :: temp !< temperature
2241 real(kind=fms_svp_kind_), intent(out), dimension(:) :: esat !< saturation vapor pressure
2242 real(kind=fms_svp_kind_), intent(out), dimension(:) :: desat !< derivative of esat
2243 integer, intent(out) :: nbad !< if temperature is out of range
2244
2245 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
2246 real(kind=fms_svp_kind_) :: del !< delta T
2247 integer :: ind, i
2248 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
2249 !! for precision consistency and for code readability, the *ll variables are declared
2250 !! and used
2251 real(kind=fms_svp_kind_) :: dtresl
2252 real(kind=fms_svp_kind_) :: tepsll
2253 real(kind=fms_svp_kind_) :: tminll
2254 real(kind=fms_svp_kind_) :: dtinvll
2255 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
2256
2257 dtresl=real(dtres, fms_svp_kind_)
2258 tminll=real(tminl, fms_svp_kind_)
2259 tepsll=real(tepsl, fms_svp_kind_)
2260 dtinvll=real(dtinvl, fms_svp_kind_)
2261
2262 nbad = 0
2263 do i = 1, size(temp,1)
2264 tmp = temp(i)-tminll
2265 ind = int(dtinvll*(tmp+tepsll))
2266 if (ind < 0 .or. ind >= table_siz) then
2267 nbad = nbad+1
2268 else
2269 del = tmp-dtresl*real(ind,fms_svp_kind_)
2270 !> esat = TABLE + del*(TABLE + del*D2TABLE)
2271 !! the coded equation below is the commented equation above
2272 esat(i) = real(table3(ind+1),fms_svp_kind_) &
2273 + del*( real(dtable3(ind+1),fms_svp_kind_) &
2274 + del*real(d2table3(ind+1),fms_svp_kind_) )
2275 !> desat = DTABLE + 2del*D2TABLE
2276 !! the coded equation below is the commented equation above
2277 desat(i) = real(dtable3(ind+1),fms_svp_kind_) &
2278 + 2.0_lkind*del*real(d2table3(ind+1),fms_svp_kind_)
2279 endif
2280 enddo
2281
2282 end subroutine lookup_es3_des3_k_1d_
2283
2284!#######################################################################
2285
2286 subroutine lookup_es3_des3_k_0d_ (temp, esat, desat, nbad)
2287 real(kind=fms_svp_kind_), intent(in) :: temp !< temperature
2288 real(kind=fms_svp_kind_), intent(out) :: esat !< saturation vapor pressure
2289 real(kind=fms_svp_kind_), intent(out) :: desat !< derivative of esat
2290 integer, intent(out) :: nbad !< if temperature is out of range
2291
2292 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
2293 real(kind=fms_svp_kind_) :: del !< delta T
2294 integer :: ind
2295 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
2296 !! for precision consistency and for code readability, the *ll variables are declared
2297 !! and used
2298 real(kind=fms_svp_kind_) :: dtresl
2299 real(kind=fms_svp_kind_) :: tepsll
2300 real(kind=fms_svp_kind_) :: tminll
2301 real(kind=fms_svp_kind_) :: dtinvll
2302 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
2303
2304 dtresl=real(dtres, fms_svp_kind_)
2305 tminll=real(tminl, fms_svp_kind_)
2306 tepsll=real(tepsl, fms_svp_kind_)
2307 dtinvll=real(dtinvl, fms_svp_kind_)
2308
2309 nbad = 0
2310 tmp = temp-tminll
2311 ind = int(dtinvll*(tmp+tepsll))
2312 if (ind < 0 .or. ind >= table_siz) then
2313 nbad = nbad+1
2314 else
2315 del = tmp-dtresl*real(ind,fms_svp_kind_)
2316 !> esat = TABLE + del*(TABLE + del*D2TABLE)
2317 !! the coded equation below is the commented equation above
2318 esat = real(table3(ind+1),fms_svp_kind_) &
2319 + del*( real(dtable3(ind+1),fms_svp_kind_) &
2320 + del*real(d2table3(ind+1),fms_svp_kind_) )
2321 !> desat = DTABLE + 2del*D2TABLE
2322 !! the coded equation below is the commented equation above
2323 desat = real(dtable3(ind+1),fms_svp_kind_) &
2324 + 2.0_lkind*del*real(d2table3(ind+1),fms_svp_kind_)
2325 endif
2326
2327 end subroutine lookup_es3_des3_k_0d_
2328
2329!#######################################################################
2330
2331 subroutine lookup_es3_k_3d_(temp, esat, nbad)
2332 real(kind=fms_svp_kind_), intent(in), dimension(:,:,:) :: temp !< temperature
2333 real(kind=fms_svp_kind_), intent(out), dimension(:,:,:) :: esat !< saturation vapor pressure
2334 integer, intent(out) :: nbad !< if temperature is out of range
2335 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
2336 real(kind=fms_svp_kind_) :: del !< delta T
2337 integer :: ind, i, j, k
2338 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
2339 !! for precision consistency and for code readability, the *ll variables are declared
2340 !! and used
2341 real(kind=fms_svp_kind_) :: dtresl
2342 real(kind=fms_svp_kind_) :: tepsll
2343 real(kind=fms_svp_kind_) :: tminll
2344 real(kind=fms_svp_kind_) :: dtinvll
2345
2346 dtresl=real(dtres, fms_svp_kind_)
2347 tminll=real(tminl, fms_svp_kind_)
2348 tepsll=real(tepsl, fms_svp_kind_)
2349 dtinvll=real(dtinvl, fms_svp_kind_)
2350
2351 nbad = 0
2352 do k = 1, size(temp,3)
2353 do j = 1, size(temp,2)
2354 do i = 1, size(temp,1)
2355 tmp = temp(i,j,k)-tminll
2356 ind = int(dtinvll*(tmp+tepsll))
2357 if (ind < 0 .or. ind >= table_siz) then
2358 nbad = nbad+1
2359 else
2360 del = tmp-dtresl*real(ind,fms_svp_kind_)
2361 !> esat = TABLE + del*(TABLE + del*D2TABLE)
2362 !! the coded equation below is the commented equation above
2363 esat(i,j,k) = real(table3(ind+1),fms_svp_kind_) &
2364 + del*( real(dtable3(ind+1),fms_svp_kind_) &
2365 + del*real(d2table3(ind+1),fms_svp_kind_) )
2366 endif
2367 enddo
2368 enddo
2369 enddo
2370
2371 end subroutine lookup_es3_k_3d_
2372
2373!#######################################################################
2374
2375 subroutine lookup_des3_k_3d_(temp, desat, nbad)
2376 real(kind=fms_svp_kind_), intent(in), dimension(:,:,:) :: temp !< temperature
2377 real(kind=fms_svp_kind_), intent(out), dimension(:,:,:) :: desat!< derivatove of saturation vap pressure
2378 integer, intent(out) :: nbad !< if temperature is out of range
2379 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
2380 real(kind=fms_svp_kind_) :: del !< delta T
2381 integer :: ind, i, j, k
2382 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
2383 !! for precision consistency and for code readability, the *ll variables are declared
2384 !! and used
2385 real(kind=fms_svp_kind_) :: dtresl
2386 real(kind=fms_svp_kind_) :: tepsll
2387 real(kind=fms_svp_kind_) :: tminll
2388 real(kind=fms_svp_kind_) :: dtinvll
2389 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
2390
2391 dtresl=real(dtres, fms_svp_kind_)
2392 tminll=real(tminl, fms_svp_kind_)
2393 tepsll=real(tepsl, fms_svp_kind_)
2394 dtinvll=real(dtinvl, fms_svp_kind_)
2395
2396 nbad = 0
2397 do k = 1, size(temp,3)
2398 do j = 1, size(temp,2)
2399 do i = 1, size(temp,1)
2400 tmp = temp(i,j,k)-tminll
2401 ind = int(dtinvll*(tmp+tepsll))
2402 if (ind < 0 .or. ind >= table_siz) then
2403 nbad = nbad+1
2404 else
2405 del = tmp-dtresl*real(ind,fms_svp_kind_)
2406 !> desat = DTABLE + 2del*D2TABLE
2407 !! the coded equation below is the commented equation above
2408 desat(i,j,k) = real(dtable3(ind+1),fms_svp_kind_) &
2409 + 2.0_lkind*del*real(d2table3(ind+1),fms_svp_kind_)
2410 endif
2411 enddo
2412 enddo
2413 enddo
2414
2415 end subroutine lookup_des3_k_3d_
2416
2417!#######################################################################
2418 subroutine lookup_des3_k_2d_(temp, desat, nbad)
2419 real(kind=fms_svp_kind_), intent(in), dimension(:,:) :: temp !< temperature
2420 real(kind=fms_svp_kind_), intent(out), dimension(:,:) :: desat !< derivative of sat vapor pressure
2421 integer, intent(out) :: nbad !< if temperature is out of range
2422 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
2423 real(kind=fms_svp_kind_) :: del !< delta T
2424 integer :: ind, i, j
2425 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
2426 !! for precision consistency and for code readability, the *ll variables are declared
2427 !! and used
2428 real(kind=fms_svp_kind_) :: dtresl
2429 real(kind=fms_svp_kind_) :: tepsll
2430 real(kind=fms_svp_kind_) :: tminll
2431 real(kind=fms_svp_kind_) :: dtinvll
2432 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
2433
2434 dtresl=real(dtres, fms_svp_kind_)
2435 tminll=real(tminl, fms_svp_kind_)
2436 tepsll=real(tepsl, fms_svp_kind_)
2437 dtinvll=real(dtinvl, fms_svp_kind_)
2438
2439 nbad = 0
2440 do j = 1, size(temp,2)
2441 do i = 1, size(temp,1)
2442 tmp = temp(i,j)-tminll
2443 ind = int(dtinvll*(tmp+tepsll))
2444 if (ind < 0 .or. ind >= table_siz) then
2445 nbad = nbad+1
2446 else
2447 del = tmp-dtresl*real(ind,fms_svp_kind_)
2448 !> desat = DTABLE + 2del*D2TABLE
2449 !! the coded equation below is the commented equation above
2450 desat(i,j) = real(dtable3(ind+1),fms_svp_kind_) &
2451 + 2.0_lkind*del*real(d2table3(ind+1),fms_svp_kind_)
2452 endif
2453 enddo
2454 enddo
2455
2456 end subroutine lookup_des3_k_2d_
2457!#######################################################################
2458 subroutine lookup_es3_k_2d_(temp, esat, nbad)
2459 real(kind=fms_svp_kind_), intent(in), dimension(:,:) :: temp !< temperature
2460 real(kind=fms_svp_kind_), intent(out), dimension(:,:) :: esat !< saturation vapor pressure
2461 integer, intent(out) :: nbad !< if temperature is out of range
2462 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
2463 real(kind=fms_svp_kind_) :: del !< delta T
2464 integer :: ind, i, j
2465 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
2466 !! for precision consistency and for code readability, the *ll variables are declared
2467 !! and used
2468 real(kind=fms_svp_kind_) :: dtresl
2469 real(kind=fms_svp_kind_) :: tepsll
2470 real(kind=fms_svp_kind_) :: tminll
2471 real(kind=fms_svp_kind_) :: dtinvll
2472
2473 dtresl=real(dtres, fms_svp_kind_)
2474 tminll=real(tminl, fms_svp_kind_)
2475 tepsll=real(tepsl, fms_svp_kind_)
2476 dtinvll=real(dtinvl, fms_svp_kind_)
2477
2478 nbad = 0
2479 do j = 1, size(temp,2)
2480 do i = 1, size(temp,1)
2481 tmp = temp(i,j)-tminll
2482 ind = int(dtinvll*(tmp+tepsll))
2483 if (ind < 0 .or. ind >= table_siz) then
2484 nbad = nbad+1
2485 else
2486 del = tmp-dtresl*real(ind,fms_svp_kind_)
2487 !> esat = TABLE + del*(TABLE + del*D2TABLE)
2488 !! the coded equation below is the commented equation above
2489 esat(i,j) = real(table3(ind+1),fms_svp_kind_) &
2490 + del*( real(dtable3(ind+1),fms_svp_kind_) &
2491 + del*real(d2table3(ind+1),fms_svp_kind_) )
2492 endif
2493 enddo
2494 enddo
2495
2496 end subroutine lookup_es3_k_2d_
2497!#######################################################################
2498 subroutine lookup_des3_k_1d_(temp, desat, nbad)
2499 real(kind=fms_svp_kind_), intent(in), dimension(:) :: temp !< temperature
2500 real(kind=fms_svp_kind_), intent(out), dimension(:) :: desat !< derivative of sat vapor pressure
2501 integer, intent(out) :: nbad !< if temperature is out of range
2502 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
2503 real(kind=fms_svp_kind_) :: del !< delta T
2504 integer :: ind, i
2505 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
2506 !! for precision consistency and for code readability, the *ll variables are declared
2507 !! and used
2508 real(kind=fms_svp_kind_) :: dtresl
2509 real(kind=fms_svp_kind_) :: tepsll
2510 real(kind=fms_svp_kind_) :: tminll
2511 real(kind=fms_svp_kind_) :: dtinvll
2512 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
2513
2514 dtresl=real(dtres, fms_svp_kind_)
2515 tminll=real(tminl, fms_svp_kind_)
2516 tepsll=real(tepsl, fms_svp_kind_)
2517 dtinvll=real(dtinvl, fms_svp_kind_)
2518
2519 nbad = 0
2520 do i = 1, size(temp,1)
2521 tmp = temp(i)-tminll
2522 ind = int(dtinvll*(tmp+tepsll))
2523 if (ind < 0 .or. ind >= table_siz) then
2524 nbad = nbad+1
2525 else
2526 del = tmp-dtresl*real(ind,fms_svp_kind_)
2527 !> desat = DTABLE + 2del*D2TABLE
2528 !! the coded equation below is the commented equation above
2529 desat(i) = real(dtable3(ind+1),fms_svp_kind_) &
2530 + 2.0_lkind*del*real(d2table3(ind+1),fms_svp_kind_)
2531 endif
2532 enddo
2533
2534 end subroutine lookup_des3_k_1d_
2535!#######################################################################
2536 subroutine lookup_es3_k_1d_(temp, esat, nbad)
2537 real(kind=fms_svp_kind_), intent(in), dimension(:) :: temp !< temperature
2538 real(kind=fms_svp_kind_), intent(out), dimension(:) :: esat !< saturation vapor pressure
2539 integer, intent(out) :: nbad !< if temperature is out of range
2540 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
2541 real(kind=fms_svp_kind_) :: del !< delta T
2542 integer :: ind, i
2543 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
2544 !! for precision consistency and for code readability, the *ll variables are declared
2545 !! and used
2546 real(kind=fms_svp_kind_) :: dtresl
2547 real(kind=fms_svp_kind_) :: tepsll
2548 real(kind=fms_svp_kind_) :: tminll
2549 real(kind=fms_svp_kind_) :: dtinvll
2550
2551 dtresl=real(dtres, fms_svp_kind_)
2552 tminll=real(tminl, fms_svp_kind_)
2553 tepsll=real(tepsl, fms_svp_kind_)
2554 dtinvll=real(dtinvl, fms_svp_kind_)
2555
2556 nbad = 0
2557 do i = 1, size(temp,1)
2558 tmp = temp(i)-tminll
2559 ind = int(dtinvll*(tmp+tepsll))
2560 if (ind < 0 .or. ind >= table_siz) then
2561 nbad = nbad+1
2562 else
2563 del = tmp-dtresl*real(ind,fms_svp_kind_)
2564 !> esat = TABLE + del*(TABLE + del*D2TABLE)
2565 !! the coded equation below is the commented equation above
2566 esat(i) = real(table3(ind+1),fms_svp_kind_) &
2567 + del*( real(dtable3(ind+1),fms_svp_kind_) &
2568 + del*real(d2table3(ind+1),fms_svp_kind_) )
2569 endif
2570 enddo
2571
2572 end subroutine lookup_es3_k_1d_
2573!#######################################################################
2574 subroutine lookup_des3_k_0d_(temp, desat, nbad)
2575 real(kind=fms_svp_kind_), intent(in) :: temp !< temperature
2576 real(kind=fms_svp_kind_), intent(out) :: desat !< derivative of sat vapor pressure
2577 integer, intent(out) :: nbad !< if temperature is out of range
2578 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
2579 real(kind=fms_svp_kind_) :: del !< delta T
2580 integer :: ind
2581 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
2582 !! for precision consistency and for code readability, the *ll variables are declared
2583 !! and used
2584 real(kind=fms_svp_kind_) :: dtresl
2585 real(kind=fms_svp_kind_) :: tepsll
2586 real(kind=fms_svp_kind_) :: tminll
2587 real(kind=fms_svp_kind_) :: dtinvll
2588 integer, parameter :: lkind=fms_svp_kind_ !< local FMS_SVP_KIND_
2589
2590 dtresl=real(dtres, fms_svp_kind_)
2591 tminll=real(tminl, fms_svp_kind_)
2592 tepsll=real(tepsl, fms_svp_kind_)
2593 dtinvll=real(dtinvl, fms_svp_kind_)
2594
2595 nbad = 0
2596 tmp = temp-tminll
2597 ind = int(dtinvll*(tmp+tepsll))
2598 if (ind < 0 .or. ind >= table_siz) then
2599 nbad = nbad+1
2600 else
2601 del = tmp-dtresl*real(ind,fms_svp_kind_)
2602 !> desat = DTABLE + 2del*D2TABLE
2603 !! the coded equation below is the commented equation above
2604 desat = real(dtable3(ind+1),fms_svp_kind_) &
2605 + 2.0_lkind*del*real(d2table3(ind+1),fms_svp_kind_)
2606 endif
2607
2608 end subroutine lookup_des3_k_0d_
2609!#######################################################################
2610 subroutine lookup_es3_k_0d_(temp, esat, nbad)
2611 real(kind=fms_svp_kind_), intent(in) :: temp !< temperature
2612 real(kind=fms_svp_kind_), intent(out) :: esat !< saturation vapor pressure
2613 integer, intent(out) :: nbad !< if temperature is out of range
2614 real(kind=fms_svp_kind_) :: tmp !< temp-TMINLL
2615 real(kind=fms_svp_kind_) :: del !< delta T
2616 integer :: ind
2617 !> dtres, tminl, tepsl, dtinvl are module level variables declared in r8_kind precision
2618 !! for precision consistency and for code readability, the *ll variables are declared
2619 !! and used
2620 real(kind=fms_svp_kind_) :: dtresl
2621 real(kind=fms_svp_kind_) :: tepsll
2622 real(kind=fms_svp_kind_) :: tminll
2623 real(kind=fms_svp_kind_) :: dtinvll
2624
2625 dtresl=real(dtres, fms_svp_kind_)
2626 tminll=real(tminl, fms_svp_kind_)
2627 tepsll=real(tepsl, fms_svp_kind_)
2628 dtinvll=real(dtinvl, fms_svp_kind_)
2629
2630 nbad = 0
2631 tmp = temp-tminll
2632 ind = int(dtinvll*(tmp+tepsll))
2633 if (ind < 0 .or. ind >= table_siz) then
2634 nbad = nbad+1
2635 else
2636 del = tmp-dtresl*real(ind,fms_svp_kind_)
2637 !> esat = TABLE + del*(TABLE + del*D2TABLE)
2638 !! the coded equation below is the commented equation above
2639 esat = real(table3(ind+1),fms_svp_kind_) &
2640 + del*( real(dtable3(ind+1),fms_svp_kind_) &
2641 + del*real(d2table3(ind+1),fms_svp_kind_) )
2642 endif
2643
2644 end subroutine lookup_es3_k_0d_
2645!#######################################################################
2646!> @}
2647! close documentation grouping
subroutine lookup_es2_k_3d_(temp, esat, nbad)
subroutine lookup_es3_des3_k_1d_(temp, esat, desat, nbad)
subroutine lookup_es2_k_2d_(temp, esat, nbad)
subroutine lookup_es2_des2_k_0d_(temp, esat, desat, nbad)
subroutine compute_mrs_k_3d_(temp, press, eps, zvir, mrs, nbad, mr, hc, dmrsdt, esat, es_over_liq, es_over_liq_and_ice)
subroutine lookup_es_des_k_3d_(temp, esat, desat, nbad)
subroutine lookup_es3_k_2d_(temp, esat, nbad)
subroutine compute_qs_k_1d_(temp, press, eps, zvir, qs, nbad, q, hc, dqsdt, esat, es_over_liq, es_over_liq_and_ice)
subroutine lookup_des_k_3d_(temp, desat, nbad)
subroutine sat_vapor_pres_init_k_(table_size, tcmin, tcmax, tfreeze, hlv, rvgas, es0, err_msg, use_exact_qs_input, do_simple, construct_table_wrt_liq, construct_table_wrt_liq_and_ice, teps, tmin, dtinv)
This routine has been generalized to return tables for any temperature range and resolution The TABLE...
real(kind=fms_svp_kind_) function, dimension(size(tem, 1)) compute_es_liq_k_(tem, tfreeze)
subroutine lookup_es3_k_0d_(temp, esat, nbad)
subroutine lookup_es3_des3_k_2d_(temp, esat, desat, nbad)
subroutine lookup_es_k_2d_(temp, esat, nbad)
subroutine lookup_es3_k_1d_(temp, esat, nbad)
subroutine lookup_es2_des2_k_3d_(temp, esat, desat, nbad)
subroutine lookup_es_k_3d_(temp, esat, nbad)
subroutine lookup_des2_k_2d_(temp, desat, nbad)
subroutine lookup_des2_k_0d_(temp, desat, nbad)
subroutine lookup_es_des_k_2d_(temp, esat, desat, nbad)
subroutine lookup_es_des_k_0d_(temp, esat, desat, nbad)
subroutine lookup_es_des_k_1d_(temp, esat, desat, nbad)
subroutine lookup_des_k_0d_(temp, desat, nbad)
subroutine lookup_es2_des2_k_2d_(temp, esat, desat, nbad)
subroutine lookup_des2_k_1d_(temp, desat, nbad)
subroutine lookup_es2_k_1d_(temp, esat, nbad)
real(kind=fms_svp_kind_) function, dimension(size(tem, 1)) compute_es_liq_ice_k_(tem, tfreeze)
subroutine lookup_es3_des3_k_0d_(temp, esat, desat, nbad)
subroutine lookup_es_k_0d_(temp, esat, nbad)
subroutine lookup_es3_k_3d_(temp, esat, nbad)
real(kind=fms_svp_kind_) function, dimension(size(tem, 1)) compute_es_k_(tem, tfreeze)
subroutine lookup_es3_des3_k_3d_(temp, esat, desat, nbad)
subroutine lookup_des_k_1d_(temp, desat, nbad)
subroutine lookup_es2_k_0d_(temp, esat, nbad)
subroutine lookup_des_k_2d_(temp, desat, nbad)
subroutine compute_qs_k_3d_(temp, press, eps, zvir, qs, nbad, q, hc, dqsdt, esat, es_over_liq, es_over_liq_and_ice)
subroutine lookup_des3_k_2d_(temp, desat, nbad)
subroutine lookup_es2_des2_k_1d_(temp, esat, desat, nbad)
subroutine lookup_des3_k_3d_(temp, desat, nbad)
subroutine lookup_es_k_1d_(temp, esat, nbad)
subroutine compute_mrs_k_0d_(temp, press, eps, zvir, mrs, nbad, mr, hc, dmrsdt, esat, es_over_liq, es_over_liq_and_ice)
subroutine lookup_des3_k_1d_(temp, desat, nbad)
subroutine compute_mrs_k_2d_(temp, press, eps, zvir, mrs, nbad, mr, hc, dmrsdt, esat, es_over_liq, es_over_liq_and_ice)
subroutine compute_qs_k_2d_(temp, press, eps, zvir, qs, nbad, q, hc, dqsdt, esat, es_over_liq, es_over_liq_and_ice)
subroutine compute_mrs_k_1d_(temp, press, eps, zvir, mrs, nbad, mr, hc, dmrsdt, esat, es_over_liq, es_over_liq_and_ice)
subroutine lookup_des2_k_3d_(temp, desat, nbad)
subroutine compute_qs_k_0d_(temp, press, eps, zvir, qs, nbad, q, hc, dqsdt, esat, es_over_liq, es_over_liq_and_ice)
subroutine lookup_des3_k_0d_(temp, desat, nbad)