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