FV3 Bundle
sat_vapor_pres_k.F90
Go to the documentation of this file.
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
10 !*
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14 !* for more details.
15 !*
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 
21 
22 ! This module is what I (pjp) think a kernel should be.
23 ! There have been many proposals as to what a kernel should look like.
24 ! If fact, so many different ideas have been expressed that the lack
25 ! of agreement has greatly hampered progress.
26 ! The only way to move forward is to limit the requirments for a kernel
27 ! to only what is widely agreeded upon.
28 ! I believe that there are only two things widely agreeded upon.
29 
30 ! 1) A kernel should be independent of the rest of FMS so that it can
31 ! easily be ported into another programming system.
32 ! This requires that a kernel does not access anything by use association.
33 ! The one exception is this kernel, because it is not practical for physics
34 ! modules to avoid using a module that computes the saturation vapor
35 ! pressure of water vapor.
36 
37 ! 2) For the sake of thread safety, module globals should be written only at initialization.
38 ! In this case, the module globals are the tables and a handful of scalars.
39 
40 ! 3) A kernel should not read from an external file.
41 
42 ! One of the things that was not widely agreeded upon is that a kernel should
43 ! not be a fortran module. This complicates things greatly for questionable
44 ! benefit and could be done as a second step anyway, if necessary.
45 
46  implicit none
47  private
48 
49 ! Include variable "version" to be written to log file.
50 #include<file_version.h>
51 
52  public :: sat_vapor_pres_init_k
53  public :: lookup_es_k
54  public :: lookup_des_k
55  public :: lookup_es_des_k
56  public :: lookup_es2_k
57  public :: lookup_des2_k
58  public :: lookup_es2_des2_k
59  public :: lookup_es3_k
60  public :: lookup_des3_k
61  public :: lookup_es3_des3_k
62  public :: compute_qs_k
63  public :: compute_mrs_k
64 
65  interface lookup_es_k
66  module procedure lookup_es_k_0d
67  module procedure lookup_es_k_1d
68  module procedure lookup_es_k_2d
69  module procedure lookup_es_k_3d
70  end interface
71 
72  interface lookup_des_k
73  module procedure lookup_des_k_0d
74  module procedure lookup_des_k_1d
75  module procedure lookup_des_k_2d
76  module procedure lookup_des_k_3d
77  end interface
78 
79  interface lookup_es_des_k
80  module procedure lookup_es_des_k_0d
81  module procedure lookup_es_des_k_1d
82  module procedure lookup_es_des_k_2d
83  module procedure lookup_es_des_k_3d
84  end interface
85 
86  interface lookup_es2_k
87  module procedure lookup_es2_k_0d
88  module procedure lookup_es2_k_1d
89  module procedure lookup_es2_k_2d
90  module procedure lookup_es2_k_3d
91  end interface
92 
93  interface lookup_des2_k
94  module procedure lookup_des2_k_0d
95  module procedure lookup_des2_k_1d
96  module procedure lookup_des2_k_2d
97  module procedure lookup_des2_k_3d
98  end interface
99 
101  module procedure lookup_es2_des2_k_0d
102  module procedure lookup_es2_des2_k_1d
103  module procedure lookup_es2_des2_k_2d
104  module procedure lookup_es2_des2_k_3d
105  end interface
106 
107  interface lookup_es3_k
108  module procedure lookup_es3_k_0d
109  module procedure lookup_es3_k_1d
110  module procedure lookup_es3_k_2d
111  module procedure lookup_es3_k_3d
112  end interface
113 
114  interface lookup_des3_k
115  module procedure lookup_des3_k_0d
116  module procedure lookup_des3_k_1d
117  module procedure lookup_des3_k_2d
118  module procedure lookup_des3_k_3d
119  end interface
120 
122  module procedure lookup_es3_des3_k_0d
123  module procedure lookup_es3_des3_k_1d
124  module procedure lookup_es3_des3_k_2d
125  module procedure lookup_es3_des3_k_3d
126  end interface
127 
128  interface compute_qs_k
129  module procedure compute_qs_k_0d
130  module procedure compute_qs_k_1d
131  module procedure compute_qs_k_2d
132  module procedure compute_qs_k_3d
133  end interface
134 
135  interface compute_mrs_k
136  module procedure compute_mrs_k_0d
137  module procedure compute_mrs_k_1d
138  module procedure compute_mrs_k_2d
139  module procedure compute_mrs_k_3d
140  end interface
141 
142  real :: dtres, tepsl, tminl, dtinvl
143  integer :: table_siz
144  real, dimension(:), allocatable :: table ! sat vapor pres (es)
145  real, dimension(:), allocatable :: dtable ! first derivative of es
146  real, dimension(:), allocatable :: d2table ! second derivative of es
147  real, dimension(:), allocatable :: table2 ! sat vapor pres (es)
148  real, dimension(:), allocatable :: dtable2 ! first derivative of es
149  real, dimension(:), allocatable :: d2table2 ! second derivative of es
150  real, dimension(:), allocatable :: table3 ! sat vapor pres (es)
151  real, dimension(:), allocatable :: dtable3 ! first derivative of es
152  real, dimension(:), allocatable :: d2table3 ! second derivative of es
153 
154  logical :: use_exact_qs
155  logical :: module_is_initialized = .false.
156 
157  contains
158 
159  subroutine sat_vapor_pres_init_k(table_size, tcmin, tcmax, TFREEZE, HLV, RVGAS, ES0, err_msg, &
160  use_exact_qs_input, do_simple, &
161  construct_table_wrt_liq, &
162  construct_table_wrt_liq_and_ice, &
163  teps, tmin, dtinv)
165 ! This routine has been generalized to return tables for any temperature range and resolution
166 
167  integer, intent(in) :: table_size
168  real, intent(in) :: tcmin ! TABLE(1) = sat vapor pressure at temperature tcmin (deg C)
169  real, intent(in) :: tcmax ! TABLE(table_size) = sat vapor pressure at temperature tcmax (deg C)
170  real, intent(in) :: tfreeze, hlv, rvgas, es0
171  logical, intent(in) :: use_exact_qs_input, do_simple
172  logical, intent(in) :: construct_table_wrt_liq
173  logical, intent(in) :: construct_table_wrt_liq_and_ice
174  character(len=*), intent(out) :: err_msg
175  real, intent(out), optional :: teps, tmin, dtinv
176 
177 ! increment used to generate derivative table
178  real, dimension(3) :: tem(3), es(3)
179  real :: hdtinv, tinrc, tfact
180  integer :: i
181 
182  err_msg = ''
183 
184  if (module_is_initialized) return
185 
186  if(allocated(table) .or. allocated(dtable) .or. allocated(d2table)) then
187  err_msg = 'Attempt to allocate sat vapor pressure tables when already allocated'
188  return
189  else
190  allocate(table(table_size), dtable(table_size), d2table(table_size))
191  endif
192 
193  if (construct_table_wrt_liq) then
194  if(allocated(table2) .or. allocated(dtable2) .or. allocated(d2table2)) then
195  err_msg = 'Attempt to allocate sat vapor pressure table2s when already allocated'
196  return
197  else
198  allocate(table2(table_size), dtable2(table_size), d2table2(table_size))
199  endif
200  endif
201 
202  if (construct_table_wrt_liq_and_ice) then
203  if(allocated(table3) .or. allocated(dtable3) .or. allocated(d2table3)) then
204  err_msg = 'Attempt to allocate sat vapor pressure table2s when already allocated'
205  return
206  else
207  allocate(table3(table_size), dtable3(table_size), d2table3(table_size))
208  endif
209  endif
210 
211  table_siz = table_size
212  dtres = (tcmax - tcmin)/real(table_size-1)
213  tminl = real(tcmin)+tfreeze ! minimum valid temp in table
214  dtinvl = 1./dtres
215  tepsl = .5*dtres
216  tinrc = .1*dtres
217  if(present(teps )) teps =tepsl
218  if(present(tmin )) tmin =tminl
219  if(present(dtinv)) dtinv=dtinvl
220 
221 ! To be able to compute tables for any temperature range and resolution,
222 ! and at the same time exactly reproduce answers from memphis revision,
223 ! it is necessary to compute ftact differently than it is in memphis.
224  tfact = 5.0*dtinvl
225 
226  hdtinv = dtinvl*0.5
227 
228 ! compute es tables from tcmin to tcmax
229 ! estimate es derivative with small +/- difference
230 
231  if (do_simple) then
232 
233  do i = 1, table_size
234  tem(1) = tminl + dtres*real(i-1)
235  table(i) = es0*610.78*exp(-hlv/rvgas*(1./tem(1) - 1./tfreeze))
236  dtable(i) = hlv*table(i)/rvgas/tem(1)**2.
237  enddo
238 
239  else
240 
241  do i = 1, table_size
242  tem(1) = tminl + dtres*real(i-1)
243  tem(2) = tem(1)-tinrc
244  tem(3) = tem(1)+tinrc
245  es = compute_es_k(tem, tfreeze)
246  table(i) = es(1)
247  dtable(i) = (es(3)-es(2))*tfact
248  enddo
249 
250  endif !if (do_simple)
251 
252 ! compute one-half second derivative using centered differences
253 ! differencing des values in the table
254 
255  do i = 2, table_size-1
256  d2table(i) = 0.25*dtinvl*(dtable(i+1)-dtable(i-1))
257  enddo
258  ! one-sided derivatives at boundaries
259 
260  d2table(1) = 0.50*dtinvl*(dtable(2)-dtable(1))
261 
262  d2table(table_size) = 0.50*dtinvl*&
263  (dtable(table_size)-dtable(table_size-1))
264 
265  if (construct_table_wrt_liq) then
266 ! compute es tables from tcmin to tcmax
267 ! estimate es derivative with small +/- difference
268 
269  do i = 1, table_size
270  tem(1) = tminl + dtres*real(i-1)
271  tem(2) = tem(1)-tinrc
272  tem(3) = tem(1)+tinrc
273 ! pass in flag to force all values to be wrt liquid
274  es = compute_es_liq_k(tem, tfreeze)
275  table2(i) = es(1)
276  dtable2(i) = (es(3)-es(2))*tfact
277  enddo
278 
279 ! compute one-half second derivative using centered differences
280 ! differencing des values in the table
281 
282  do i = 2, table_size-1
283  d2table2(i) = 0.25*dtinvl*(dtable2(i+1)-dtable2(i-1))
284  enddo
285 ! one-sided derivatives at boundaries
286 
287  d2table2(1) = 0.50*dtinvl*(dtable2(2)-dtable2(1))
288 
289  d2table2(table_size) = 0.50*dtinvl*&
290  (dtable2(table_size)-dtable2(table_size-1))
291  endif
292 
293 
294  if (construct_table_wrt_liq_and_ice) then
295 ! compute es tables from tcmin to tcmax
296 ! estimate es derivative with small +/- difference
297 
298  do i = 1, table_size
299  tem(1) = tminl + dtres*real(i-1)
300  tem(2) = tem(1)-tinrc
301  tem(3) = tem(1)+tinrc
302 ! pass in flag to force all values to be wrt liquid
303  es = compute_es_liq_ice_k(tem, tfreeze)
304  table3(i) = es(1)
305  dtable3(i) = (es(3)-es(2))*tfact
306  enddo
307 
308 ! compute one-half second derivative using centered differences
309 ! differencing des values in the table
310 
311  do i = 2, table_size-1
312  d2table3(i) = 0.25*dtinvl*(dtable3(i+1)-dtable3(i-1))
313  enddo
314 ! one-sided derivatives at boundaries
315 
316  d2table3(1) = 0.50*dtinvl*(dtable3(2)-dtable3(1))
317 
318  d2table3(table_size) = 0.50*dtinvl*&
319  (dtable3(table_size)-dtable3(table_size-1))
320  endif
321 
322  use_exact_qs = use_exact_qs_input
323  module_is_initialized = .true.
324 
325  end subroutine sat_vapor_pres_init_k
326 
327 !#######################################################################
328 
329  function compute_es_k(tem, TFREEZE) result (es)
330  real, intent(in) :: tem(:), tfreeze
331  real :: es(size(tem,1))
332 
333  real :: x, esice, esh2o, tbasw, tbasi
334  integer :: i
335  real, parameter :: esbasw = 101324.60
336  real, parameter :: esbasi = 610.71
337 
338  tbasw = tfreeze+100.
339  tbasi = tfreeze
340 
341  do i = 1, size(tem)
342 
343 ! compute es over ice
344 
345  if (tem(i) < tbasi) then
346  x = -9.09718*(tbasi/tem(i)-1.0) - 3.56654*log10(tbasi/tem(i)) &
347  +0.876793*(1.0-tem(i)/tbasi) + log10(esbasi)
348  esice =10.**(x)
349  else
350  esice = 0.
351  endif
352 
353 ! compute es over water greater than -20 c.
354 ! values over 100 c may not be valid
355 ! see smithsonian meteorological tables page 350.
356 
357  if (tem(i) > -20.+tbasi) then
358  x = -7.90298*(tbasw/tem(i)-1.0) + 5.02808*log10(tbasw/tem(i)) &
359  -1.3816e-07*(10.0**((1.0-tem(i)/tbasw)*11.344)-1.0) &
360  +8.1328e-03*(10.0**((tbasw/tem(i)-1.0)*(-3.49149))-1.0) &
361  +log10(esbasw)
362  esh2o = 10.**(x)
363  else
364  esh2o = 0.
365  endif
366 
367 ! derive blended es over ice and supercooled water between -20c and 0c
368 
369  if (tem(i) <= -20.+tbasi) then
370  es(i) = esice
371  else if (tem(i) >= tbasi) then
372  es(i) = esh2o
373  else
374  es(i) = 0.05*((tbasi-tem(i))*esice + (tem(i)-tbasi+20.)*esh2o)
375  endif
376 
377  enddo
378 
379  end function compute_es_k
380 
381 !#######################################################################
382 
383  function compute_es_liq_k(tem, TFREEZE) result (es)
384  real, intent(in) :: tem(:), tfreeze
385  real :: es(size(tem,1))
386 
387  real :: x, esh2o, tbasw
388  integer :: i
389  real, parameter :: esbasw = 101324.60
390 
391  tbasw = tfreeze+100.
392 
393  do i = 1, size(tem)
394 
395 
396 ! compute es over water for all temps.
397 ! values over 100 c may not be valid
398 ! see smithsonian meteorological tables page 350.
399 
400  x = -7.90298*(tbasw/tem(i)-1.0) + 5.02808*log10(tbasw/tem(i)) &
401  -1.3816e-07*(10.0**((1.0-tem(i)/tbasw)*11.344)-1.0) &
402  +8.1328e-03*(10.0**((tbasw/tem(i)-1.0)*(-3.49149))-1.0) &
403  +log10(esbasw)
404  esh2o = 10.**(x)
405 
406 
407  es(i) = esh2o
408 
409  enddo
410 
411  end function compute_es_liq_k
412 
413 !#######################################################################
414 
415  function compute_es_liq_ice_k(tem, TFREEZE) result (es)
416  real, intent(in) :: tem(:), tfreeze
417  real :: es(size(tem,1))
418 
419  real :: x, tbasw, tbasi
420  integer :: i
421  real, parameter :: esbasw = 101324.60
422  real, parameter :: esbasi = 610.71
423 
424  tbasw = tfreeze+100.
425  tbasi = tfreeze
426 
427  do i = 1, size(tem)
428 
429  if (tem(i) < tbasi) then
430 
431 ! compute es over ice
432 
433  x = -9.09718*(tbasi/tem(i)-1.0) - 3.56654*log10(tbasi/tem(i)) &
434  +0.876793*(1.0-tem(i)/tbasi) + log10(esbasi)
435  es(i) =10.**(x)
436  else
437 
438 ! compute es over water
439 ! values over 100 c may not be valid
440 ! see smithsonian meteorological tables page 350.
441 
442  x = -7.90298*(tbasw/tem(i)-1.0) + 5.02808*log10(tbasw/tem(i)) &
443  -1.3816e-07*(10.0**((1.0-tem(i)/tbasw)*11.344)-1.0) &
444  +8.1328e-03*(10.0**((tbasw/tem(i)-1.0)*(-3.49149))-1.0) &
445  +log10(esbasw)
446  es(i) = 10.**(x)
447  endif
448 
449  enddo
450 
451  end function compute_es_liq_ice_k
452 
453 !#######################################################################
454 
455  subroutine compute_qs_k_3d (temp, press, eps, zvir, qs, nbad, q, hc, &
456  dqsdT, esat, es_over_liq, es_over_liq_and_ice)
458  real, intent(in), dimension(:,:,:) :: temp, press
459  real, intent(in) :: eps, zvir
460  real, intent(out), dimension(:,:,:) :: qs
461  integer, intent(out) :: nbad
462  real, intent(in), dimension(:,:,:), optional :: q
463  real, intent(in), optional :: hc
464  real, intent(out), dimension(:,:,:), optional :: dqsdT, esat
465  logical,intent(in), optional :: es_over_liq
466  logical,intent(in), optional :: es_over_liq_and_ice
467 
468  real, dimension(size(temp,1), size(temp,2), size(temp,3)) :: &
469  esloc, desat, denom
470  integer :: i, j, k
471  real :: hc_loc
472 
473  if (present(hc)) then
474  hc_loc = hc
475  else
476  hc_loc = 1.0
477  endif
478  if (present(es_over_liq)) then
479  if (present (dqsdt)) then
480  call lookup_es2_des2_k (temp, esloc, desat, nbad)
481  desat = desat*hc_loc
482  else
483  call lookup_es2_k (temp, esloc, nbad)
484  endif
485  else if (present(es_over_liq_and_ice)) then
486  if (present (dqsdt)) then
487  call lookup_es3_des3_k (temp, esloc, desat, nbad)
488  desat = desat*hc_loc
489  else
490  call lookup_es3_k (temp, esloc, nbad)
491  endif
492  else
493  if (present (dqsdt)) then
494  call lookup_es_des_k (temp, esloc, desat, nbad)
495  desat = desat*hc_loc
496  else
497  call lookup_es_k (temp, esloc, nbad)
498  endif
499  endif
500  esloc = esloc*hc_loc
501  if (present (esat)) then
502  esat = esloc
503  endif
504  if (nbad == 0) then
505  if (present (q) .and. use_exact_qs) then
506  qs = (1.0 + zvir*q)*eps*esloc/press
507  if (present (dqsdt)) then
508  dqsdt = (1.0 + zvir*q)*eps*desat/press
509  endif
510  else ! (present(q))
511  denom = press - (1.0 - eps)*esloc
512  do k=1,size(qs,3)
513  do j=1,size(qs,2)
514  do i=1,size(qs,1)
515  if (denom(i,j,k) > 0.0) then
516  qs(i,j,k) = eps*esloc(i,j,k)/denom(i,j,k)
517  else
518  qs(i,j,k) = eps
519  endif
520  end do
521  end do
522  end do
523  if (present (dqsdt)) then
524  dqsdt = eps*press*desat/denom**2
525  endif
526  endif ! (present(q))
527  else ! (nbad = 0)
528  qs = -999.
529  if (present (dqsdt)) then
530  dqsdt = -999.
531  endif
532  if (present (esat)) then
533  esat = -999.
534  endif
535  endif ! (nbad = 0)
536 
537 
538  end subroutine compute_qs_k_3d
539 
540 !#######################################################################
541 
542  subroutine compute_qs_k_2d (temp, press, eps, zvir, qs, nbad, q, hc, &
543  dqsdT, esat, es_over_liq, es_over_liq_and_ice)
545  real, intent(in), dimension(:,:) :: temp, press
546  real, intent(in) :: eps, zvir
547  real, intent(out), dimension(:,:) :: qs
548  integer, intent(out) :: nbad
549  real, intent(in), dimension(:,:), optional :: q
550  real, intent(in), optional :: hc
551  real, intent(out), dimension(:,:), optional :: dqsdT, esat
552  logical,intent(in), optional :: es_over_liq
553  logical,intent(in), optional :: es_over_liq_and_ice
554 
555  real, dimension(size(temp,1), size(temp,2)) :: esloc, desat, denom
556  integer :: i, j
557  real :: hc_loc
558 
559  if (present(hc)) then
560  hc_loc = hc
561  else
562  hc_loc = 1.0
563  endif
564 
565  if (present(es_over_liq)) then
566  if (present (dqsdt)) then
567  call lookup_es2_des2_k (temp, esloc, desat, nbad)
568  desat = desat*hc_loc
569  else
570  call lookup_es2_k (temp, esloc, nbad)
571  endif
572  else if (present(es_over_liq_and_ice)) then
573  if (present (dqsdt)) then
574  call lookup_es3_des3_k (temp, esloc, desat, nbad)
575  desat = desat*hc_loc
576  else
577  call lookup_es3_k (temp, esloc, nbad)
578  endif
579  else
580  if (present (dqsdt)) then
581  call lookup_es_des_k (temp, esloc, desat, nbad)
582  desat = desat*hc_loc
583  else
584  call lookup_es_k (temp, esloc, nbad)
585  endif
586  endif
587  esloc = esloc*hc_loc
588  if (present (esat)) then
589  esat = esloc
590  endif
591  if (nbad == 0) then
592  if (present (q) .and. use_exact_qs) then
593  qs = (1.0 + zvir*q)*eps*esloc/press
594  if (present (dqsdt)) then
595  dqsdt = (1.0 + zvir*q)*eps*desat/press
596  endif
597  else ! (present(q))
598  denom = press - (1.0 - eps)*esloc
599  do j=1,size(qs,2)
600  do i=1,size(qs,1)
601  if (denom(i,j) > 0.0) then
602  qs(i,j) = eps*esloc(i,j)/denom(i,j)
603  else
604  qs(i,j) = eps
605  endif
606  end do
607  end do
608  if (present (dqsdt)) then
609  dqsdt = eps*press*desat/denom**2
610  endif
611  endif ! (present(q))
612  else ! (nbad = 0)
613  qs = -999.
614  if (present (dqsdt)) then
615  dqsdt = -999.
616  endif
617  if (present (esat)) then
618  esat = -999.
619  endif
620  endif ! (nbad = 0)
621 
622 
623  end subroutine compute_qs_k_2d
624 
625 !#######################################################################
626 
627  subroutine compute_qs_k_1d (temp, press, eps, zvir, qs, nbad, q, hc, &
628  dqsdT, esat, es_over_liq, es_over_liq_and_ice)
630  real, intent(in), dimension(:) :: temp, press
631  real, intent(in) :: eps, zvir
632  real, intent(out), dimension(:) :: qs
633  integer, intent(out) :: nbad
634  real, intent(in), dimension(:), optional :: q
635  real, intent(in), optional :: hc
636  real, intent(out), dimension(:), optional :: dqsdT, esat
637  logical,intent(in), optional :: es_over_liq
638  logical,intent(in), optional :: es_over_liq_and_ice
639 
640  real, dimension(size(temp,1)) :: esloc, desat, denom
641  integer :: i
642  real :: hc_loc
643 
644  if (present(hc)) then
645  hc_loc = hc
646  else
647  hc_loc = 1.0
648  endif
649 
650  if (present(es_over_liq)) then
651  if (present (dqsdt)) then
652  call lookup_es2_des2_k (temp, esloc, desat, nbad)
653  desat = desat*hc_loc
654  else
655  call lookup_es2_k (temp, esloc, nbad)
656  endif
657  else if (present(es_over_liq_and_ice)) then
658  if (present (dqsdt)) then
659  call lookup_es3_des3_k (temp, esloc, desat, nbad)
660  desat = desat*hc_loc
661  else
662  call lookup_es3_k (temp, esloc, nbad)
663  endif
664  else
665  if (present (dqsdt)) then
666  call lookup_es_des_k (temp, esloc, desat, nbad)
667  desat = desat*hc_loc
668  else
669  call lookup_es_k (temp, esloc, nbad)
670  endif
671  endif
672  esloc = esloc*hc_loc
673  if (present (esat)) then
674  esat = esloc
675  endif
676  if (nbad == 0) then
677  if (present (q) .and. use_exact_qs) then
678  qs = (1.0 + zvir*q)*eps*esloc/press
679  if (present (dqsdt)) then
680  dqsdt = (1.0 + zvir*q)*eps*desat/press
681  endif
682  else ! (present(q))
683  denom = press - (1.0 - eps)*esloc
684  do i=1,size(qs,1)
685  if (denom(i) > 0.0) then
686  qs(i) = eps*esloc(i)/denom(i)
687  else
688  qs(i) = eps
689  endif
690  end do
691  if (present (dqsdt)) then
692  dqsdt = eps*press*desat/denom**2
693  endif
694  endif ! (present(q))
695  else ! (nbad = 0)
696  qs = -999.
697  if (present (dqsdt)) then
698  dqsdt = -999.
699  endif
700  if (present (esat)) then
701  esat = -999.
702  endif
703  endif ! (nbad = 0)
704 
705 
706  end subroutine compute_qs_k_1d
707 
708 !#######################################################################
709 
710  subroutine compute_qs_k_0d (temp, press, eps, zvir, qs, nbad, q, hc, &
711  dqsdT, esat, es_over_liq, es_over_liq_and_ice)
713  real, intent(in) :: temp, press
714  real, intent(in) :: eps, zvir
715  real, intent(out) :: qs
716  integer, intent(out) :: nbad
717  real, intent(in), optional :: q
718  real, intent(in), optional :: hc
719  real, intent(out), optional :: dqsdT, esat
720  logical,intent(in), optional :: es_over_liq
721  logical,intent(in), optional :: es_over_liq_and_ice
722 
723  real :: esloc, desat, denom
724  real :: hc_loc
725 
726  if (present(hc)) then
727  hc_loc = hc
728  else
729  hc_loc = 1.0
730  endif
731 
732  if (present(es_over_liq)) then
733  if (present (dqsdt)) then
734  call lookup_es2_des2_k (temp, esloc, desat, nbad)
735  desat = desat*hc_loc
736  else
737  call lookup_es2_k (temp, esloc, nbad)
738  endif
739  else if (present(es_over_liq_and_ice)) then
740  if (present (dqsdt)) then
741  call lookup_es3_des3_k (temp, esloc, desat, nbad)
742  desat = desat*hc_loc
743  else
744  call lookup_es3_k (temp, esloc, nbad)
745  endif
746  else
747  if (present (dqsdt)) then
748  call lookup_es_des_k (temp, esloc, desat, nbad)
749  desat = desat*hc_loc
750  else
751  call lookup_es_k (temp, esloc, nbad)
752  endif
753  endif
754  esloc = esloc*hc_loc
755  if (present (esat)) then
756  esat = esloc
757  endif
758  if (nbad == 0) then
759  if (present (q) .and. use_exact_qs) then
760  qs = (1.0 + zvir*q)*eps*esloc/press
761  if (present (dqsdt)) then
762  dqsdt = (1.0 + zvir*q)*eps*desat/press
763  endif
764  else ! (present(q))
765  denom = press - (1.0 - eps)*esloc
766  if (denom > 0.0) then
767  qs = eps*esloc/denom
768  else
769  qs = eps
770  endif
771  if (present (dqsdt)) then
772  dqsdt = eps*press*desat/denom**2
773  endif
774  endif ! (present(q))
775  else ! (nbad = 0)
776  qs = -999.
777  if (present (dqsdt)) then
778  dqsdt = -999.
779  endif
780  if (present (esat)) then
781  esat = -999.
782  endif
783  endif ! (nbad = 0)
784 
785 
786  end subroutine compute_qs_k_0d
787 
788 !#######################################################################
789 
790 !#######################################################################
791 
792  subroutine compute_mrs_k_3d (temp, press, eps, zvir, mrs, nbad, &
793  mr, hc, dmrsdT, esat,es_over_liq, es_over_liq_and_ice)
795  real, intent(in), dimension(:,:,:) :: temp, press
796  real, intent(in) :: eps, zvir
797  real, intent(out), dimension(:,:,:) :: mrs
798  integer, intent(out) :: nbad
799  real, intent(in), dimension(:,:,:), optional :: mr
800  real, intent(in), optional :: hc
801  real, intent(out), dimension(:,:,:), optional :: dmrsdT, esat
802  logical,intent(in), optional :: es_over_liq
803  logical,intent(in), optional :: es_over_liq_and_ice
804 
805  real, dimension(size(temp,1), size(temp,2), size(temp,3)) :: &
806  esloc, desat, denom
807  integer :: i, j, k
808  real :: hc_loc
809 
810  if (present(hc)) then
811  hc_loc = hc
812  else
813  hc_loc = 1.0
814  endif
815 
816  if (present (es_over_liq)) then
817  if (present (dmrsdt)) then
818  call lookup_es2_des2_k (temp, esloc, desat, nbad)
819  desat = desat*hc_loc
820  else
821  call lookup_es2_k (temp, esloc, nbad)
822  endif
823  else if (present(es_over_liq_and_ice)) then
824  if (present (dmrsdt)) then
825  call lookup_es3_des3_k (temp, esloc, desat, nbad)
826  desat = desat*hc_loc
827  else
828  call lookup_es3_k (temp, esloc, nbad)
829  endif
830  else
831  if (present (dmrsdt)) then
832  call lookup_es_des_k (temp, esloc, desat, nbad)
833  desat = desat*hc_loc
834  else
835  call lookup_es_k (temp, esloc, nbad)
836  endif
837  endif
838  esloc = esloc*hc_loc
839  if (present (esat)) then
840  esat = esloc
841  endif
842  if (nbad == 0) then
843  if (present (mr) .and. use_exact_qs) then
844  mrs = (eps + mr)*esloc/press
845  if (present (dmrsdt)) then
846  dmrsdt = (eps + mr)*desat/press
847  endif
848  else ! (present (mr))
849  denom = press - esloc
850  do k=1,size(mrs,3)
851  do j=1,size(mrs,2)
852  do i=1,size(mrs,1)
853  if (denom(i,j,k) > 0.0) then
854  mrs(i,j,k) = eps*esloc(i,j,k)/denom(i,j,k)
855  else
856  mrs(i,j,k) = eps
857  endif
858  end do
859  end do
860  end do
861  if (present (dmrsdt)) then
862  dmrsdt = eps*press*desat/denom**2
863  endif
864  endif !(present (mr))
865  else
866  mrs = -999.
867  if (present (dmrsdt)) then
868  dmrsdt = -999.
869  endif
870  if (present (esat)) then
871  esat = -999.
872  endif
873  endif
874 
875 
876  end subroutine compute_mrs_k_3d
877 
878 !#######################################################################
879 
880  subroutine compute_mrs_k_2d (temp, press, eps, zvir, mrs, nbad, &
881  mr, hc, dmrsdT, esat,es_over_liq, es_over_liq_and_ice)
883  real, intent(in), dimension(:,:) :: temp, press
884  real, intent(in) :: eps, zvir
885  real, intent(out), dimension(:,:) :: mrs
886  integer, intent(out) :: nbad
887  real, intent(in), dimension(:,:), optional :: mr
888  real, intent(in), optional :: hc
889  real, intent(out), dimension(:,:), optional :: dmrsdT, esat
890  logical,intent(in), optional :: es_over_liq
891  logical,intent(in), optional :: es_over_liq_and_ice
892 
893  real, dimension(size(temp,1), size(temp,2)) :: esloc, desat, denom
894  integer :: i, j
895  real :: hc_loc
896 
897  if (present(hc)) then
898  hc_loc = hc
899  else
900  hc_loc = 1.0
901  endif
902 
903  if (present (es_over_liq)) then
904  if (present (dmrsdt)) then
905  call lookup_es2_des2_k (temp, esloc, desat, nbad)
906  desat = desat*hc_loc
907  else
908  call lookup_es2_k (temp, esloc, nbad)
909  endif
910  else if (present(es_over_liq_and_ice)) then
911  if (present (dmrsdt)) then
912  call lookup_es3_des3_k (temp, esloc, desat, nbad)
913  desat = desat*hc_loc
914  else
915  call lookup_es3_k (temp, esloc, nbad)
916  endif
917  else
918  if (present (dmrsdt)) then
919  call lookup_es_des_k (temp, esloc, desat, nbad)
920  desat = desat*hc_loc
921  else
922  call lookup_es_k (temp, esloc, nbad)
923  endif
924  endif
925  esloc = esloc*hc_loc
926  if (present (esat)) then
927  esat = esloc
928  endif
929  if (nbad == 0) then
930  if (present (mr) .and. use_exact_qs) then
931  mrs = (eps + mr)*esloc/press
932  if (present (dmrsdt)) then
933  dmrsdt = (eps + mr)*desat/press
934  endif
935  else ! (present (mr))
936  denom = press - esloc
937  do j=1,size(mrs,2)
938  do i=1,size(mrs,1)
939  if (denom(i,j) > 0.0) then
940  mrs(i,j) = eps*esloc(i,j)/denom(i,j)
941  else
942  mrs(i,j) = eps
943  endif
944  end do
945  end do
946  if (present (dmrsdt)) then
947  dmrsdt = eps*press*desat/denom**2
948  endif
949  endif !(present (mr))
950  else
951  mrs = -999.
952  if (present (dmrsdt)) then
953  dmrsdt = -999.
954  endif
955  if (present (esat)) then
956  esat = -999.
957  endif
958  endif
959 
960 
961  end subroutine compute_mrs_k_2d
962 
963 !#######################################################################
964 
965  subroutine compute_mrs_k_1d (temp, press, eps, zvir, mrs, nbad, &
966  mr, hc, dmrsdT, esat,es_over_liq, es_over_liq_and_ice)
968  real, intent(in), dimension(:) :: temp, press
969  real, intent(in) :: eps, zvir
970  real, intent(out), dimension(:) :: mrs
971  integer, intent(out) :: nbad
972  real, intent(in), dimension(:), optional :: mr
973  real, intent(in), optional :: hc
974  real, intent(out), dimension(:), optional :: dmrsdT, esat
975  logical,intent(in), optional :: es_over_liq
976  logical,intent(in), optional :: es_over_liq_and_ice
977 
978  real, dimension(size(temp,1)) :: esloc, desat, denom
979  integer :: i
980  real :: hc_loc
981 
982  if (present(hc)) then
983  hc_loc = hc
984  else
985  hc_loc = 1.0
986  endif
987 
988  if (present (es_over_liq)) then
989  if (present (dmrsdt)) then
990  call lookup_es2_des2_k (temp, esloc, desat, nbad)
991  desat = desat*hc_loc
992  else
993  call lookup_es2_k (temp, esloc, nbad)
994  endif
995  else if (present(es_over_liq_and_ice)) then
996  if (present (dmrsdt)) then
997  call lookup_es3_des3_k (temp, esloc, desat, nbad)
998  desat = desat*hc_loc
999  else
1000  call lookup_es3_k (temp, esloc, nbad)
1001  endif
1002  else
1003  if (present (dmrsdt)) then
1004  call lookup_es_des_k (temp, esloc, desat, nbad)
1005  desat = desat*hc_loc
1006  else
1007  call lookup_es_k (temp, esloc, nbad)
1008  endif
1009  endif
1010  esloc = esloc*hc_loc
1011  if (present (esat)) then
1012  esat = esloc
1013  endif
1014  if (nbad == 0) then
1015  if (present (mr) .and. use_exact_qs) then
1016  mrs = (eps + mr)*esloc/press
1017  if (present (dmrsdt)) then
1018  dmrsdt = (eps + mr)*desat/press
1019  endif
1020  else ! (present (mr))
1021  denom = press - esloc
1022  do i=1,size(mrs,1)
1023  if (denom(i) > 0.0) then
1024  mrs(i) = eps*esloc(i)/denom(i)
1025  else
1026  mrs(i) = eps
1027  endif
1028  end do
1029  if (present (dmrsdt)) then
1030  dmrsdt = eps*press*desat/denom**2
1031  endif
1032  endif !(present (mr))
1033  else
1034  mrs = -999.
1035  if (present (dmrsdt)) then
1036  dmrsdt = -999.
1037  endif
1038  if (present (esat)) then
1039  esat = -999.
1040  endif
1041  endif
1042 
1043 
1044  end subroutine compute_mrs_k_1d
1045 
1046 !#######################################################################
1047 
1048  subroutine compute_mrs_k_0d (temp, press, eps, zvir, mrs, nbad, &
1049  mr, hc, dmrsdT, esat,es_over_liq, es_over_liq_and_ice)
1051  real, intent(in) :: temp, press
1052  real, intent(in) :: eps, zvir
1053  real, intent(out) :: mrs
1054  integer, intent(out) :: nbad
1055  real, intent(in), optional :: mr
1056  real, intent(in), optional :: hc
1057  real, intent(out), optional :: dmrsdT, esat
1058  logical,intent(in), optional :: es_over_liq
1059  logical,intent(in), optional :: es_over_liq_and_ice
1060 
1061  real :: esloc, desat, denom
1062  real :: hc_loc
1063 
1064  if (present(hc)) then
1065  hc_loc = hc
1066  else
1067  hc_loc = 1.0
1068  endif
1069 
1070  if (present (es_over_liq)) then
1071  if (present (dmrsdt)) then
1072  call lookup_es2_des2_k (temp, esloc, desat, nbad)
1073  desat = desat*hc_loc
1074  else
1075  call lookup_es2_k (temp, esloc, nbad)
1076  endif
1077  else if (present(es_over_liq_and_ice)) then
1078  if (present (dmrsdt)) then
1079  call lookup_es3_des3_k (temp, esloc, desat, nbad)
1080  desat = desat*hc_loc
1081  else
1082  call lookup_es3_k (temp, esloc, nbad)
1083  endif
1084  else
1085  if (present (dmrsdt)) then
1086  call lookup_es_des_k (temp, esloc, desat, nbad)
1087  desat = desat*hc_loc
1088  else
1089  call lookup_es_k (temp, esloc, nbad)
1090  endif
1091  endif
1092  esloc = esloc*hc_loc
1093  if (present (esat)) then
1094  esat = esloc
1095  endif
1096  if (nbad == 0) then
1097  if (present (mr) .and. use_exact_qs) then
1098  mrs = (eps + mr)*esloc/press
1099  if (present (dmrsdt)) then
1100  dmrsdt = (eps + mr)*desat/press
1101  endif
1102  else ! (present (mr))
1103  denom = press - esloc
1104  if (denom > 0.0) then
1105  mrs = eps*esloc/denom
1106  else
1107  mrs = eps
1108  endif
1109  if (present (dmrsdt)) then
1110  dmrsdt = eps*press*desat/denom**2
1111  endif
1112  endif !(present (mr))
1113  else
1114  mrs = -999.
1115  if (present (dmrsdt)) then
1116  dmrsdt = -999.
1117  endif
1118  if (present (esat)) then
1119  esat = -999.
1120  endif
1121  endif
1122 
1123 
1124  end subroutine compute_mrs_k_0d
1125 
1126 
1127 
1128 !#######################################################################
1129 
1130  subroutine lookup_es_des_k_3d (temp, esat, desat, nbad)
1131  real, intent(in), dimension(:,:,:) :: temp
1132  real, intent(out), dimension(:,:,:) :: esat, desat
1133  integer, intent(out) :: nbad
1134 
1135  real :: tmp, del
1136  integer :: ind, i, j, k
1137 
1138  nbad = 0
1139  do k = 1, size(temp,3)
1140  do j = 1, size(temp,2)
1141  do i = 1, size(temp,1)
1142  tmp = temp(i,j,k)-tminl
1143  ind = int(dtinvl*(tmp+tepsl))
1144  if (ind < 0 .or. ind >= table_siz) then
1145  nbad = nbad+1
1146  else
1147  del = tmp-dtres*real(ind)
1148  esat(i,j,k) = table(ind+1) + &
1149  del*(dtable(ind+1) + del*d2table(ind+1))
1150  desat(i,j,k) = dtable(ind+1) + 2.*del*d2table(ind+1)
1151  endif
1152  enddo
1153  enddo
1154  enddo
1155 
1156  end subroutine lookup_es_des_k_3d
1157 
1158 !#######################################################################
1159 
1160  subroutine lookup_es_des_k_2d (temp, esat, desat, nbad)
1161  real, intent(in), dimension(:,:) :: temp
1162  real, intent(out), dimension(:,:) :: esat, desat
1163  integer, intent(out) :: nbad
1164 
1165  real :: tmp, del
1166  integer :: ind, i, j
1167 
1168  nbad = 0
1169  do j = 1, size(temp,2)
1170  do i = 1, size(temp,1)
1171  tmp = temp(i,j)-tminl
1172  ind = int(dtinvl*(tmp+tepsl))
1173  if (ind < 0 .or. ind >= table_siz) then
1174  nbad = nbad+1
1175  else
1176  del = tmp-dtres*real(ind)
1177  esat(i,j) = table(ind+1) + &
1178  del*(dtable(ind+1) + del*d2table(ind+1))
1179  desat(i,j) = dtable(ind+1) + 2.*del*d2table(ind+1)
1180  endif
1181  enddo
1182  enddo
1183 
1184  end subroutine lookup_es_des_k_2d
1185 
1186 !#######################################################################
1187 
1188  subroutine lookup_es_des_k_1d (temp, esat, desat, nbad)
1189  real, intent(in), dimension(:) :: temp
1190  real, intent(out), dimension(:) :: esat, desat
1191  integer, intent(out) :: nbad
1192 
1193  real :: tmp, del
1194  integer :: ind, i
1195 
1196  nbad = 0
1197  do i = 1, size(temp,1)
1198  tmp = temp(i)-tminl
1199  ind = int(dtinvl*(tmp+tepsl))
1200  if (ind < 0 .or. ind >= table_siz) then
1201  nbad = nbad+1
1202  else
1203  del = tmp-dtres*real(ind)
1204  esat(i) = table(ind+1) + &
1205  del*(dtable(ind+1) + del*d2table(ind+1))
1206  desat(i) = dtable(ind+1) + 2.*del*d2table(ind+1)
1207  endif
1208  enddo
1209 
1210  end subroutine lookup_es_des_k_1d
1211 
1212 !#######################################################################
1213 
1214  subroutine lookup_es_des_k_0d (temp, esat, desat, nbad)
1215  real, intent(in) :: temp
1216  real, intent(out) :: esat, desat
1217  integer, intent(out) :: nbad
1218 
1219  real :: tmp, del
1220  integer :: ind
1221 
1222  nbad = 0
1223  tmp = temp-tminl
1224  ind = int(dtinvl*(tmp+tepsl))
1225  if (ind < 0 .or. ind >= table_siz) then
1226  nbad = nbad+1
1227  else
1228  del = tmp-dtres*real(ind)
1229  esat = table(ind+1) + &
1230  del*(dtable(ind+1) + del*d2table(ind+1))
1231  desat = dtable(ind+1) + 2.*del*d2table(ind+1)
1232  endif
1233 
1234  end subroutine lookup_es_des_k_0d
1235 
1236 !#######################################################################
1237 
1238  subroutine lookup_es_k_3d(temp, esat, nbad)
1239  real, intent(in), dimension(:,:,:) :: temp
1240  real, intent(out), dimension(:,:,:) :: esat
1241  integer, intent(out) :: nbad
1242  real :: tmp, del
1243  integer :: ind, i, j, k
1244 
1245  nbad = 0
1246  do k = 1, size(temp,3)
1247  do j = 1, size(temp,2)
1248  do i = 1, size(temp,1)
1249  tmp = temp(i,j,k)-tminl
1250  ind = int(dtinvl*(tmp+tepsl))
1251  if (ind < 0 .or. ind >= table_siz) then
1252  nbad = nbad+1
1253  else
1254  del = tmp-dtres*real(ind)
1255  esat(i,j,k) = table(ind+1) + &
1256  del*(dtable(ind+1) + del*d2table(ind+1))
1257  endif
1258  enddo
1259  enddo
1260  enddo
1261 
1262  end subroutine lookup_es_k_3d
1263 
1264 !#######################################################################
1265 
1266  subroutine lookup_des_k_3d(temp, desat, nbad)
1267  real, intent(in), dimension(:,:,:) :: temp
1268  real, intent(out), dimension(:,:,:) :: desat
1269  integer, intent(out) :: nbad
1270  real :: tmp, del
1271  integer :: ind, i, j, k
1272 
1273  nbad = 0
1274  do k = 1, size(temp,3)
1275  do j = 1, size(temp,2)
1276  do i = 1, size(temp,1)
1277  tmp = temp(i,j,k)-tminl
1278  ind = int(dtinvl*(tmp+tepsl))
1279  if (ind < 0 .or. ind >= table_siz) then
1280  nbad = nbad+1
1281  else
1282  del = tmp-dtres*real(ind)
1283  desat(i,j,k) = dtable(ind+1) + 2.*del*d2table(ind+1)
1284  endif
1285  enddo
1286  enddo
1287  enddo
1288 
1289  end subroutine lookup_des_k_3d
1290 
1291 !#######################################################################
1292  subroutine lookup_des_k_2d(temp, desat, nbad)
1293  real, intent(in), dimension(:,:) :: temp
1294  real, intent(out), dimension(:,:) :: desat
1295  integer, intent(out) :: nbad
1296  real :: tmp, del
1297  integer :: ind, i, j
1298 
1299  nbad = 0
1300  do j = 1, size(temp,2)
1301  do i = 1, size(temp,1)
1302  tmp = temp(i,j)-tminl
1303  ind = int(dtinvl*(tmp+tepsl))
1304  if (ind < 0 .or. ind >= table_siz) then
1305  nbad = nbad+1
1306  else
1307  del = tmp-dtres*real(ind)
1308  desat(i,j) = dtable(ind+1) + 2.*del*d2table(ind+1)
1309  endif
1310  enddo
1311  enddo
1312 
1313  end subroutine lookup_des_k_2d
1314 !#######################################################################
1315  subroutine lookup_es_k_2d(temp, esat, nbad)
1316  real, intent(in), dimension(:,:) :: temp
1317  real, intent(out), dimension(:,:) :: esat
1318  integer, intent(out) :: nbad
1319  real :: tmp, del
1320  integer :: ind, i, j
1321 
1322  nbad = 0
1323  do j = 1, size(temp,2)
1324  do i = 1, size(temp,1)
1325  tmp = temp(i,j)-tminl
1326  ind = int(dtinvl*(tmp+tepsl))
1327  if (ind < 0 .or. ind >= table_siz) then
1328  nbad = nbad+1
1329  else
1330  del = tmp-dtres*real(ind)
1331  esat(i,j) = table(ind+1) + del*(dtable(ind+1) + &
1332  del*d2table(ind+1))
1333  endif
1334  enddo
1335  enddo
1336 
1337  end subroutine lookup_es_k_2d
1338 !#######################################################################
1339  subroutine lookup_des_k_1d(temp, desat, nbad)
1340  real, intent(in), dimension(:) :: temp
1341  real, intent(out), dimension(:) :: desat
1342  integer, intent(out) :: nbad
1343  real :: tmp, del
1344  integer :: ind, i
1345 
1346  nbad = 0
1347  do i = 1, size(temp,1)
1348  tmp = temp(i)-tminl
1349  ind = int(dtinvl*(tmp+tepsl))
1350  if (ind < 0 .or. ind >= table_siz) then
1351  nbad = nbad+1
1352  else
1353  del = tmp-dtres*real(ind)
1354  desat(i) = dtable(ind+1) + 2.*del*d2table(ind+1)
1355  endif
1356  enddo
1357 
1358  end subroutine lookup_des_k_1d
1359 !#######################################################################
1360  subroutine lookup_es_k_1d(temp, esat, nbad)
1361  real, intent(in), dimension(:) :: temp
1362  real, intent(out), dimension(:) :: esat
1363  integer, intent(out) :: nbad
1364  real :: tmp, del
1365  integer :: ind, i
1366 
1367  nbad = 0
1368  do i = 1, size(temp,1)
1369  tmp = temp(i)-tminl
1370  ind = int(dtinvl*(tmp+tepsl))
1371  if (ind < 0 .or. ind >= table_siz) then
1372  nbad = nbad+1
1373  else
1374  del = tmp-dtres*real(ind)
1375  esat(i) = table(ind+1) + del*(dtable(ind+1) + del*d2table(ind+1))
1376  endif
1377  enddo
1378 
1379  end subroutine lookup_es_k_1d
1380 !#######################################################################
1381  subroutine lookup_des_k_0d(temp, desat, nbad)
1382  real, intent(in) :: temp
1383  real, intent(out) :: desat
1384  integer, intent(out) :: nbad
1385  real :: tmp, del
1386  integer :: ind
1387 
1388  nbad = 0
1389  tmp = temp-tminl
1390  ind = int(dtinvl*(tmp+tepsl))
1391  if (ind < 0 .or. ind >= table_siz) then
1392  nbad = nbad+1
1393  else
1394  del = tmp-dtres*real(ind)
1395  desat = dtable(ind+1) + 2.*del*d2table(ind+1)
1396  endif
1397 
1398  end subroutine lookup_des_k_0d
1399 !#######################################################################
1400  subroutine lookup_es_k_0d(temp, esat, nbad)
1401  real, intent(in) :: temp
1402  real, intent(out) :: esat
1403  integer, intent(out) :: nbad
1404  real :: tmp, del
1405  integer :: ind
1406 
1407  nbad = 0
1408  tmp = temp-tminl
1409  ind = int(dtinvl*(tmp+tepsl))
1410  if (ind < 0 .or. ind >= table_siz) then
1411  nbad = nbad+1
1412  else
1413  del = tmp-dtres*real(ind)
1414  esat = table(ind+1) + del*(dtable(ind+1) + del*d2table(ind+1))
1415  endif
1416 
1417  end subroutine lookup_es_k_0d
1418 !#######################################################################
1419 
1420  subroutine lookup_es2_des2_k_3d (temp, esat, desat, nbad)
1421  real, intent(in), dimension(:,:,:) :: temp
1422  real, intent(out), dimension(:,:,:) :: esat, desat
1423  integer, intent(out) :: nbad
1424 
1425  real :: tmp, del
1426  integer :: ind, i, j, k
1427 
1428  nbad = 0
1429  do k = 1, size(temp,3)
1430  do j = 1, size(temp,2)
1431  do i = 1, size(temp,1)
1432  tmp = temp(i,j,k)-tminl
1433  ind = int(dtinvl*(tmp+tepsl))
1434  if (ind < 0 .or. ind >= table_siz) then
1435  nbad = nbad+1
1436  else
1437  del = tmp-dtres*real(ind)
1438  esat(i,j,k) = table2(ind+1) + &
1439  del*(dtable2(ind+1) + del*d2table2(ind+1))
1440  desat(i,j,k) = dtable2(ind+1) + 2.*del*d2table2(ind+1)
1441  endif
1442  enddo
1443  enddo
1444  enddo
1445 
1446  end subroutine lookup_es2_des2_k_3d
1447 
1448 !#######################################################################
1449 
1450  subroutine lookup_es2_des2_k_2d (temp, esat, desat, nbad)
1451  real, intent(in), dimension(:,:) :: temp
1452  real, intent(out), dimension(:,:) :: esat, desat
1453  integer, intent(out) :: nbad
1454 
1455  real :: tmp, del
1456  integer :: ind, i, j
1457 
1458  nbad = 0
1459  do j = 1, size(temp,2)
1460  do i = 1, size(temp,1)
1461  tmp = temp(i,j)-tminl
1462  ind = int(dtinvl*(tmp+tepsl))
1463  if (ind < 0 .or. ind >= table_siz) then
1464  nbad = nbad+1
1465  else
1466  del = tmp-dtres*real(ind)
1467  esat(i,j) = table2(ind+1) + &
1468  del*(dtable2(ind+1) + del*d2table2(ind+1))
1469  desat(i,j) = dtable2(ind+1) + 2.*del*d2table2(ind+1)
1470  endif
1471  enddo
1472  enddo
1473 
1474  end subroutine lookup_es2_des2_k_2d
1475 
1476 !#######################################################################
1477 
1478  subroutine lookup_es2_des2_k_1d (temp, esat, desat, nbad)
1479  real, intent(in), dimension(:) :: temp
1480  real, intent(out), dimension(:) :: esat, desat
1481  integer, intent(out) :: nbad
1482 
1483  real :: tmp, del
1484  integer :: ind, i
1485 
1486  nbad = 0
1487  do i = 1, size(temp,1)
1488  tmp = temp(i)-tminl
1489  ind = int(dtinvl*(tmp+tepsl))
1490  if (ind < 0 .or. ind >= table_siz) then
1491  nbad = nbad+1
1492  else
1493  del = tmp-dtres*real(ind)
1494  esat(i) = table2(ind+1) + &
1495  del*(dtable2(ind+1) + del*d2table2(ind+1))
1496  desat(i) = dtable2(ind+1) + 2.*del*d2table2(ind+1)
1497  endif
1498  enddo
1499 
1500  end subroutine lookup_es2_des2_k_1d
1501 
1502 !#######################################################################
1503 
1504  subroutine lookup_es2_des2_k_0d (temp, esat, desat, nbad)
1505  real, intent(in) :: temp
1506  real, intent(out) :: esat, desat
1507  integer, intent(out) :: nbad
1508 
1509  real :: tmp, del
1510  integer :: ind
1511 
1512  nbad = 0
1513  tmp = temp-tminl
1514  ind = int(dtinvl*(tmp+tepsl))
1515  if (ind < 0 .or. ind >= table_siz) then
1516  nbad = nbad+1
1517  else
1518  del = tmp-dtres*real(ind)
1519  esat = table2(ind+1) + &
1520  del*(dtable2(ind+1) + del*d2table2(ind+1))
1521  desat = dtable2(ind+1) + 2.*del*d2table2(ind+1)
1522  endif
1523 
1524  end subroutine lookup_es2_des2_k_0d
1525 
1526 !#######################################################################
1527 
1528  subroutine lookup_es2_k_3d(temp, esat, nbad)
1529  real, intent(in), dimension(:,:,:) :: temp
1530  real, intent(out), dimension(:,:,:) :: esat
1531  integer, intent(out) :: nbad
1532  real :: tmp, del
1533  integer :: ind, i, j, k
1534 
1535  nbad = 0
1536  do k = 1, size(temp,3)
1537  do j = 1, size(temp,2)
1538  do i = 1, size(temp,1)
1539  tmp = temp(i,j,k)-tminl
1540  ind = int(dtinvl*(tmp+tepsl))
1541  if (ind < 0 .or. ind >= table_siz) then
1542  nbad = nbad+1
1543  else
1544  del = tmp-dtres*real(ind)
1545  esat(i,j,k) = table2(ind+1) + &
1546  del*(dtable2(ind+1) + del*d2table2(ind+1))
1547  endif
1548  enddo
1549  enddo
1550  enddo
1551 
1552  end subroutine lookup_es2_k_3d
1553 
1554 !#######################################################################
1555 
1556  subroutine lookup_des2_k_3d(temp, desat, nbad)
1557  real, intent(in), dimension(:,:,:) :: temp
1558  real, intent(out), dimension(:,:,:) :: desat
1559  integer, intent(out) :: nbad
1560  real :: tmp, del
1561  integer :: ind, i, j, k
1562 
1563  nbad = 0
1564  do k = 1, size(temp,3)
1565  do j = 1, size(temp,2)
1566  do i = 1, size(temp,1)
1567  tmp = temp(i,j,k)-tminl
1568  ind = int(dtinvl*(tmp+tepsl))
1569  if (ind < 0 .or. ind >= table_siz) then
1570  nbad = nbad+1
1571  else
1572  del = tmp-dtres*real(ind)
1573  desat(i,j,k) = dtable2(ind+1) + 2.*del*d2table2(ind+1)
1574  endif
1575  enddo
1576  enddo
1577  enddo
1578 
1579  end subroutine lookup_des2_k_3d
1580 
1581 !#######################################################################
1582  subroutine lookup_des2_k_2d(temp, desat, nbad)
1583  real, intent(in), dimension(:,:) :: temp
1584  real, intent(out), dimension(:,:) :: desat
1585  integer, intent(out) :: nbad
1586  real :: tmp, del
1587  integer :: ind, i, j
1588 
1589  nbad = 0
1590  do j = 1, size(temp,2)
1591  do i = 1, size(temp,1)
1592  tmp = temp(i,j)-tminl
1593  ind = int(dtinvl*(tmp+tepsl))
1594  if (ind < 0 .or. ind >= table_siz) then
1595  nbad = nbad+1
1596  else
1597  del = tmp-dtres*real(ind)
1598  desat(i,j) = dtable2(ind+1) + 2.*del*d2table2(ind+1)
1599  endif
1600  enddo
1601  enddo
1602 
1603  end subroutine lookup_des2_k_2d
1604 !#######################################################################
1605  subroutine lookup_es2_k_2d(temp, esat, nbad)
1606  real, intent(in), dimension(:,:) :: temp
1607  real, intent(out), dimension(:,:) :: esat
1608  integer, intent(out) :: nbad
1609  real :: tmp, del
1610  integer :: ind, i, j
1611 
1612  nbad = 0
1613  do j = 1, size(temp,2)
1614  do i = 1, size(temp,1)
1615  tmp = temp(i,j)-tminl
1616  ind = int(dtinvl*(tmp+tepsl))
1617  if (ind < 0 .or. ind >= table_siz) then
1618  nbad = nbad+1
1619  else
1620  del = tmp-dtres*real(ind)
1621  esat(i,j) = table2(ind+1) + del*(dtable2(ind+1) + &
1622  del*d2table2(ind+1))
1623  endif
1624  enddo
1625  enddo
1626 
1627  end subroutine lookup_es2_k_2d
1628 !#######################################################################
1629  subroutine lookup_des2_k_1d(temp, desat, nbad)
1630  real, intent(in), dimension(:) :: temp
1631  real, intent(out), dimension(:) :: desat
1632  integer, intent(out) :: nbad
1633  real :: tmp, del
1634  integer :: ind, i
1635 
1636  nbad = 0
1637  do i = 1, size(temp,1)
1638  tmp = temp(i)-tminl
1639  ind = int(dtinvl*(tmp+tepsl))
1640  if (ind < 0 .or. ind >= table_siz) then
1641  nbad = nbad+1
1642  else
1643  del = tmp-dtres*real(ind)
1644  desat(i) = dtable2(ind+1) + 2.*del*d2table2(ind+1)
1645  endif
1646  enddo
1647 
1648  end subroutine lookup_des2_k_1d
1649 !#######################################################################
1650  subroutine lookup_es2_k_1d(temp, esat, nbad)
1651  real, intent(in), dimension(:) :: temp
1652  real, intent(out), dimension(:) :: esat
1653  integer, intent(out) :: nbad
1654  real :: tmp, del
1655  integer :: ind, i
1656 
1657  nbad = 0
1658  do i = 1, size(temp,1)
1659  tmp = temp(i)-tminl
1660  ind = int(dtinvl*(tmp+tepsl))
1661  if (ind < 0 .or. ind >= table_siz) then
1662  nbad = nbad+1
1663  else
1664  del = tmp-dtres*real(ind)
1665  esat(i) = table2(ind+1) + del*(dtable2(ind+1) + del*d2table2(ind+1))
1666  endif
1667  enddo
1668 
1669  end subroutine lookup_es2_k_1d
1670 !#######################################################################
1671  subroutine lookup_des2_k_0d(temp, desat, nbad)
1672  real, intent(in) :: temp
1673  real, intent(out) :: desat
1674  integer, intent(out) :: nbad
1675  real :: tmp, del
1676  integer :: ind
1677 
1678  nbad = 0
1679  tmp = temp-tminl
1680  ind = int(dtinvl*(tmp+tepsl))
1681  if (ind < 0 .or. ind >= table_siz) then
1682  nbad = nbad+1
1683  else
1684  del = tmp-dtres*real(ind)
1685  desat = dtable2(ind+1) + 2.*del*d2table2(ind+1)
1686  endif
1687 
1688  end subroutine lookup_des2_k_0d
1689 !#######################################################################
1690  subroutine lookup_es2_k_0d(temp, esat, nbad)
1691  real, intent(in) :: temp
1692  real, intent(out) :: esat
1693  integer, intent(out) :: nbad
1694  real :: tmp, del
1695  integer :: ind
1696 
1697  nbad = 0
1698  tmp = temp-tminl
1699  ind = int(dtinvl*(tmp+tepsl))
1700  if (ind < 0 .or. ind >= table_siz) then
1701  nbad = nbad+1
1702  else
1703  del = tmp-dtres*real(ind)
1704  esat = table2(ind+1) + del*(dtable2(ind+1) + del*d2table2(ind+1))
1705  endif
1706 
1707  end subroutine lookup_es2_k_0d
1708 !#######################################################################
1709 
1710 !#######################################################################
1711 
1712  subroutine lookup_es3_des3_k_3d (temp, esat, desat, nbad)
1713  real, intent(in), dimension(:,:,:) :: temp
1714  real, intent(out), dimension(:,:,:) :: esat, desat
1715  integer, intent(out) :: nbad
1716 
1717  real :: tmp, del
1718  integer :: ind, i, j, k
1719 
1720  nbad = 0
1721  do k = 1, size(temp,3)
1722  do j = 1, size(temp,2)
1723  do i = 1, size(temp,1)
1724  tmp = temp(i,j,k)-tminl
1725  ind = int(dtinvl*(tmp+tepsl))
1726  if (ind < 0 .or. ind >= table_siz) then
1727  nbad = nbad+1
1728  else
1729  del = tmp-dtres*real(ind)
1730  esat(i,j,k) = table3(ind+1) + &
1731  del*(dtable3(ind+1) + del*d2table3(ind+1))
1732  desat(i,j,k) = dtable3(ind+1) + 2.*del*d2table3(ind+1)
1733  endif
1734  enddo
1735  enddo
1736  enddo
1737 
1738  end subroutine lookup_es3_des3_k_3d
1739 
1740 !#######################################################################
1741 
1742  subroutine lookup_es3_des3_k_2d (temp, esat, desat, nbad)
1743  real, intent(in), dimension(:,:) :: temp
1744  real, intent(out), dimension(:,:) :: esat, desat
1745  integer, intent(out) :: nbad
1746 
1747  real :: tmp, del
1748  integer :: ind, i, j
1749 
1750  nbad = 0
1751  do j = 1, size(temp,2)
1752  do i = 1, size(temp,1)
1753  tmp = temp(i,j)-tminl
1754  ind = int(dtinvl*(tmp+tepsl))
1755  if (ind < 0 .or. ind >= table_siz) then
1756  nbad = nbad+1
1757  else
1758  del = tmp-dtres*real(ind)
1759  esat(i,j) = table3(ind+1) + &
1760  del*(dtable3(ind+1) + del*d2table3(ind+1))
1761  desat(i,j) = dtable3(ind+1) + 2.*del*d2table3(ind+1)
1762  endif
1763  enddo
1764  enddo
1765 
1766  end subroutine lookup_es3_des3_k_2d
1767 
1768 !#######################################################################
1769 
1770  subroutine lookup_es3_des3_k_1d (temp, esat, desat, nbad)
1771  real, intent(in), dimension(:) :: temp
1772  real, intent(out), dimension(:) :: esat, desat
1773  integer, intent(out) :: nbad
1774 
1775  real :: tmp, del
1776  integer :: ind, i
1777 
1778  nbad = 0
1779  do i = 1, size(temp,1)
1780  tmp = temp(i)-tminl
1781  ind = int(dtinvl*(tmp+tepsl))
1782  if (ind < 0 .or. ind >= table_siz) then
1783  nbad = nbad+1
1784  else
1785  del = tmp-dtres*real(ind)
1786  esat(i) = table3(ind+1) + &
1787  del*(dtable3(ind+1) + del*d2table3(ind+1))
1788  desat(i) = dtable3(ind+1) + 2.*del*d2table3(ind+1)
1789  endif
1790  enddo
1791 
1792  end subroutine lookup_es3_des3_k_1d
1793 
1794 !#######################################################################
1795 
1796  subroutine lookup_es3_des3_k_0d (temp, esat, desat, nbad)
1797  real, intent(in) :: temp
1798  real, intent(out) :: esat, desat
1799  integer, intent(out) :: nbad
1800 
1801  real :: tmp, del
1802  integer :: ind
1803 
1804  nbad = 0
1805  tmp = temp-tminl
1806  ind = int(dtinvl*(tmp+tepsl))
1807  if (ind < 0 .or. ind >= table_siz) then
1808  nbad = nbad+1
1809  else
1810  del = tmp-dtres*real(ind)
1811  esat = table3(ind+1) + &
1812  del*(dtable3(ind+1) + del*d2table3(ind+1))
1813  desat = dtable3(ind+1) + 2.*del*d2table3(ind+1)
1814  endif
1815 
1816  end subroutine lookup_es3_des3_k_0d
1817 
1818 !#######################################################################
1819 
1820  subroutine lookup_es3_k_3d(temp, esat, nbad)
1821  real, intent(in), dimension(:,:,:) :: temp
1822  real, intent(out), dimension(:,:,:) :: esat
1823  integer, intent(out) :: nbad
1824  real :: tmp, del
1825  integer :: ind, i, j, k
1826 
1827  nbad = 0
1828  do k = 1, size(temp,3)
1829  do j = 1, size(temp,2)
1830  do i = 1, size(temp,1)
1831  tmp = temp(i,j,k)-tminl
1832  ind = int(dtinvl*(tmp+tepsl))
1833  if (ind < 0 .or. ind >= table_siz) then
1834  nbad = nbad+1
1835  else
1836  del = tmp-dtres*real(ind)
1837  esat(i,j,k) = table3(ind+1) + &
1838  del*(dtable3(ind+1) + del*d2table3(ind+1))
1839  endif
1840  enddo
1841  enddo
1842  enddo
1843 
1844  end subroutine lookup_es3_k_3d
1845 
1846 !#######################################################################
1847 
1848  subroutine lookup_des3_k_3d(temp, desat, nbad)
1849  real, intent(in), dimension(:,:,:) :: temp
1850  real, intent(out), dimension(:,:,:) :: desat
1851  integer, intent(out) :: nbad
1852  real :: tmp, del
1853  integer :: ind, i, j, k
1854 
1855  nbad = 0
1856  do k = 1, size(temp,3)
1857  do j = 1, size(temp,2)
1858  do i = 1, size(temp,1)
1859  tmp = temp(i,j,k)-tminl
1860  ind = int(dtinvl*(tmp+tepsl))
1861  if (ind < 0 .or. ind >= table_siz) then
1862  nbad = nbad+1
1863  else
1864  del = tmp-dtres*real(ind)
1865  desat(i,j,k) = dtable3(ind+1) + 2.*del*d2table3(ind+1)
1866  endif
1867  enddo
1868  enddo
1869  enddo
1870 
1871  end subroutine lookup_des3_k_3d
1872 
1873 !#######################################################################
1874  subroutine lookup_des3_k_2d(temp, desat, nbad)
1875  real, intent(in), dimension(:,:) :: temp
1876  real, intent(out), dimension(:,:) :: desat
1877  integer, intent(out) :: nbad
1878  real :: tmp, del
1879  integer :: ind, i, j
1880 
1881  nbad = 0
1882  do j = 1, size(temp,2)
1883  do i = 1, size(temp,1)
1884  tmp = temp(i,j)-tminl
1885  ind = int(dtinvl*(tmp+tepsl))
1886  if (ind < 0 .or. ind >= table_siz) then
1887  nbad = nbad+1
1888  else
1889  del = tmp-dtres*real(ind)
1890  desat(i,j) = dtable3(ind+1) + 2.*del*d2table3(ind+1)
1891  endif
1892  enddo
1893  enddo
1894 
1895  end subroutine lookup_des3_k_2d
1896 !#######################################################################
1897  subroutine lookup_es3_k_2d(temp, esat, nbad)
1898  real, intent(in), dimension(:,:) :: temp
1899  real, intent(out), dimension(:,:) :: esat
1900  integer, intent(out) :: nbad
1901  real :: tmp, del
1902  integer :: ind, i, j
1903 
1904  nbad = 0
1905  do j = 1, size(temp,2)
1906  do i = 1, size(temp,1)
1907  tmp = temp(i,j)-tminl
1908  ind = int(dtinvl*(tmp+tepsl))
1909  if (ind < 0 .or. ind >= table_siz) then
1910  nbad = nbad+1
1911  else
1912  del = tmp-dtres*real(ind)
1913  esat(i,j) = table3(ind+1) + del*(dtable3(ind+1) + &
1914  del*d2table3(ind+1))
1915  endif
1916  enddo
1917  enddo
1918 
1919  end subroutine lookup_es3_k_2d
1920 !#######################################################################
1921  subroutine lookup_des3_k_1d(temp, desat, nbad)
1922  real, intent(in), dimension(:) :: temp
1923  real, intent(out), dimension(:) :: desat
1924  integer, intent(out) :: nbad
1925  real :: tmp, del
1926  integer :: ind, i
1927 
1928  nbad = 0
1929  do i = 1, size(temp,1)
1930  tmp = temp(i)-tminl
1931  ind = int(dtinvl*(tmp+tepsl))
1932  if (ind < 0 .or. ind >= table_siz) then
1933  nbad = nbad+1
1934  else
1935  del = tmp-dtres*real(ind)
1936  desat(i) = dtable3(ind+1) + 2.*del*d2table3(ind+1)
1937  endif
1938  enddo
1939 
1940  end subroutine lookup_des3_k_1d
1941 !#######################################################################
1942  subroutine lookup_es3_k_1d(temp, esat, nbad)
1943  real, intent(in), dimension(:) :: temp
1944  real, intent(out), dimension(:) :: esat
1945  integer, intent(out) :: nbad
1946  real :: tmp, del
1947  integer :: ind, i
1948 
1949  nbad = 0
1950  do i = 1, size(temp,1)
1951  tmp = temp(i)-tminl
1952  ind = int(dtinvl*(tmp+tepsl))
1953  if (ind < 0 .or. ind >= table_siz) then
1954  nbad = nbad+1
1955  else
1956  del = tmp-dtres*real(ind)
1957  esat(i) = table3(ind+1) + del*(dtable3(ind+1) + del*d2table3(ind+1))
1958  endif
1959  enddo
1960 
1961  end subroutine lookup_es3_k_1d
1962 !#######################################################################
1963  subroutine lookup_des3_k_0d(temp, desat, nbad)
1964  real, intent(in) :: temp
1965  real, intent(out) :: desat
1966  integer, intent(out) :: nbad
1967  real :: tmp, del
1968  integer :: ind
1969 
1970  nbad = 0
1971  tmp = temp-tminl
1972  ind = int(dtinvl*(tmp+tepsl))
1973  if (ind < 0 .or. ind >= table_siz) then
1974  nbad = nbad+1
1975  else
1976  del = tmp-dtres*real(ind)
1977  desat = dtable3(ind+1) + 2.*del*d2table3(ind+1)
1978  endif
1979 
1980  end subroutine lookup_des3_k_0d
1981 !#######################################################################
1982  subroutine lookup_es3_k_0d(temp, esat, nbad)
1983  real, intent(in) :: temp
1984  real, intent(out) :: esat
1985  integer, intent(out) :: nbad
1986  real :: tmp, del
1987  integer :: ind
1988 
1989  nbad = 0
1990  tmp = temp-tminl
1991  ind = int(dtinvl*(tmp+tepsl))
1992  if (ind < 0 .or. ind >= table_siz) then
1993  nbad = nbad+1
1994  else
1995  del = tmp-dtres*real(ind)
1996  esat = table3(ind+1) + del*(dtable3(ind+1) + del*d2table3(ind+1))
1997  endif
1998 
1999  end subroutine lookup_es3_k_0d
2000 !#######################################################################
2001  end module sat_vapor_pres_k_mod
subroutine lookup_es3_k_2d(temp, esat, nbad)
real function, dimension(size(tem, 1)) compute_es_liq_k(tem, TFREEZE)
subroutine compute_qs_k_3d(temp, press, eps, zvir, qs, nbad, q, hc, dqsdT, esat, es_over_liq, es_over_liq_and_ice)
subroutine compute_mrs_k_2d(temp, press, eps, zvir, mrs, nbad, mr, hc, dmrsdT, esat, es_over_liq, es_over_liq_and_ice)
subroutine lookup_des2_k_1d(temp, desat, nbad)
subroutine lookup_es_des_k_3d(temp, esat, desat, nbad)
subroutine lookup_es3_k_3d(temp, esat, nbad)
subroutine compute_qs_k_0d(temp, press, eps, zvir, qs, nbad, q, hc, dqsdT, esat, es_over_liq, es_over_liq_and_ice)
real, dimension(:), allocatable d2table3
subroutine lookup_des2_k_3d(temp, desat, nbad)
subroutine lookup_es2_k_1d(temp, esat, nbad)
subroutine lookup_es_des_k_1d(temp, esat, desat, nbad)
subroutine lookup_es2_k_2d(temp, esat, nbad)
subroutine lookup_des_k_3d(temp, desat, nbad)
real, dimension(:), allocatable table
subroutine lookup_es_k_1d(temp, esat, nbad)
subroutine lookup_es_k_2d(temp, esat, 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, public 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)
subroutine lookup_es2_des2_k_0d(temp, esat, desat, nbad)
subroutine compute_mrs_k_3d(temp, press, eps, zvir, mrs, nbad, mr, hc, dmrsdT, esat, es_over_liq, es_over_liq_and_ice)
real, dimension(:), allocatable d2table
subroutine lookup_es2_k_0d(temp, esat, nbad)
subroutine lookup_es3_k_0d(temp, esat, nbad)
subroutine lookup_es2_des2_k_3d(temp, esat, desat, nbad)
real function, dimension(size(tem, 1)) compute_es_liq_ice_k(tem, TFREEZE)
subroutine lookup_es2_des2_k_2d(temp, esat, desat, nbad)
subroutine lookup_des2_k_2d(temp, desat, nbad)
subroutine lookup_des3_k_2d(temp, desat, nbad)
subroutine lookup_es3_des3_k_3d(temp, esat, desat, nbad)
subroutine lookup_es3_des3_k_1d(temp, esat, desat, nbad)
real, dimension(:), allocatable table2
subroutine lookup_es_des_k_2d(temp, esat, desat, nbad)
subroutine lookup_es2_des2_k_1d(temp, esat, desat, nbad)
subroutine lookup_des3_k_0d(temp, desat, nbad)
subroutine lookup_des_k_1d(temp, desat, nbad)
subroutine lookup_des_k_2d(temp, desat, nbad)
real, dimension(:), allocatable d2table2
subroutine lookup_es_k_0d(temp, esat, nbad)
real, dimension(:), allocatable dtable2
subroutine lookup_es2_k_3d(temp, esat, nbad)
subroutine lookup_es_k_3d(temp, esat, nbad)
real, dimension(:), allocatable table3
subroutine lookup_des_k_0d(temp, desat, nbad)
************************************************************************GNU Lesser General Public License **This file is part of the GFDL Flexible Modeling System(FMS). ! *! *FMS is free software without even the implied warranty of MERCHANTABILITY or *FITNESS FOR A PARTICULAR PURPOSE See the GNU General Public License *for more details **You should have received a copy of the GNU Lesser General Public *License along with FMS If see< http:! ***********************************************************************subroutine READ_RECORD_CORE_(unit, field, nwords, data, start, axsiz) integer, intent(in) ::unit type(fieldtype), intent(in) ::field integer, intent(in) ::nwords MPP_TYPE_, intent(inout) ::data(nwords) integer, intent(in) ::start(:), axsiz(:) integer(SHORT_KIND) ::i2vals(nwords)!rab used in conjunction with transfer intrinsic to determine size of a variable integer(KIND=1) ::one_byte(8) integer ::word_sz!#ifdef __sgi integer(INT_KIND) ::ivals(nwords) real(FLOAT_KIND) ::rvals(nwords)!#else! integer ::ivals(nwords)! real ::rvals(nwords)!#endif real(DOUBLE_KIND) ::r8vals(nwords) pointer(ptr1, i2vals) pointer(ptr2, ivals) pointer(ptr3, rvals) pointer(ptr4, r8vals) if(mpp_io_stack_size< nwords) call mpp_io_set_stack_size(nwords) call mpp_error(FATAL, 'MPP_READ currently requires use_netCDF option') end subroutine READ_RECORD_CORE_ subroutine READ_RECORD_(unit, field, nwords, data, time_level, domain, position, tile_count, start_in, axsiz_in)!routine that is finally called by all mpp_read routines to perform the read!a non-netCDF record contains:! field ID! a set of 4 coordinates(is:ie, js:je) giving the data subdomain! a timelevel and a timestamp(=NULLTIME if field is static)! 3D real data(stored as 1D)!if you are using direct access I/O, the RECL argument to OPEN must be large enough for the above!in a global direct access file, record position on PE is given by %record.!Treatment of timestamp:! We assume that static fields have been passed without a timestamp.! Here that is converted into a timestamp of NULLTIME.! For non-netCDF fields, field is treated no differently, but is written! with a timestamp of NULLTIME. There is no check in the code to prevent! the user from repeatedly writing a static field. integer, intent(in) ::unit, nwords type(fieldtype), intent(in) ::field MPP_TYPE_, intent(inout) ::data(nwords) integer, intent(in), optional ::time_level type(domain2D), intent(in), optional ::domain integer, intent(in), optional ::position, tile_count integer, intent(in), optional ::start_in(:), axsiz_in(:) integer, dimension(size(field%axes(:))) ::start, axsiz integer ::tlevel !, subdomain(4) integer ::i, error, is, ie, js, je, isg, ieg, jsg, jeg type(domain2d), pointer ::io_domain=> tlevel if(PRESENT(start_in) .AND. PRESENT(axsiz_in)) then if(size(start(! the data domain and compute domain must refer to the subdomain being passed ! In this ! since that attempts to gather all data on PE size(field%axes(:)) axsiz(i)
real function, dimension(size(tem, 1)) compute_es_k(tem, TFREEZE)
subroutine lookup_es3_des3_k_0d(temp, esat, desat, nbad)
real, dimension(:), allocatable dtable3
subroutine lookup_es3_des3_k_2d(temp, esat, desat, nbad)
subroutine compute_qs_k_1d(temp, press, eps, zvir, qs, nbad, q, hc, dqsdT, esat, es_over_liq, es_over_liq_and_ice)
subroutine lookup_es_des_k_0d(temp, esat, desat, nbad)
subroutine lookup_des3_k_1d(temp, desat, nbad)
subroutine lookup_des2_k_0d(temp, 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 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_3d(temp, desat, nbad)
real, dimension(:), allocatable dtable
subroutine lookup_es3_k_1d(temp, esat, nbad)