ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_WaterVaporSaturation.H
Go to the documentation of this file.
1 //
2 // This module provides an interface to water vapor saturation methods, providing
3 // saturation vapor pressure and related calculations to CAM.
4 //
5 // The original wv_saturation codes were introduced by J. J. Hack,
6 // February 1990. The code has been extensively rewritten since then,
7 // including a total refactoring in Summer 2012.
8 //
9 // Methods:
10 //
11 // Pure water/ice saturation vapor pressures are calculated on the
12 // fly, with the specific method determined by a runtime option.
13 //
14 // The default method for calculating SVP is determined by a namelist
15 // option, and used whenever svp_water/ice or qsat are called.
16 //
17 #ifndef ERF_WATER_VAPOR_SATURATION_H_
18 #define ERF_WATER_VAPOR_SATURATION_H_
19 
20 #include <string>
21 #include <vector>
22 #include <memory>
23 
24 #include "ERF_Constants.H"
25 #include "ERF_Sat_methods.H"
26 
27 // Radiation code interface class
29 public:
30  // Make these public parameters in case another module wants to see the
31  // extent of the table.
32  static constexpr amrex::Real tmin = 127.16;
33  static constexpr amrex::Real tmax = 375.16;
34 
35  // Set coefficients for polynomial approximation of difference
36  // between saturation vapor press over water and saturation pressure
37  // over ice for -ttrice < t < 0 (degrees C). NOTE: polynomial is
38  // valid in the range -40 < t < 0 (degrees C).
39  static constexpr int npcf = 5;
40  static constexpr const amrex::Real pcf[npcf] = {5.04469588506e-01,
41  -5.47288442819e+00,
42  -3.67471858735e-01,
43  -8.95963532403e-03,
44  -7.78053686625e-05};
45 
46  // Compute saturation vapor pressure over water
47  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
48  static amrex::Real svp_water (const amrex::Real& t) {
50  }
51 
52  // Compute saturation vapor pressure over ice
53  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
54  static amrex::Real svp_ice (const amrex::Real& t) {
56  }
57 
58  // Compute saturation vapor pressure with an ice-water transition
59  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
60  static amrex::Real svp_trans (const amrex::Real& t) {
62  }
63 
64  // Get enthalpy based only on temperature
65  // and specific humidity.
66  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
67  static amrex::Real tq_enthalpy (const amrex::Real& t,
68  const amrex::Real& q,
69  const amrex::Real& hltalt) {
70  return Cp_d * t + hltalt * q;
71  }
72 
73  //------------------------------------------------
74  // LATENT HEAT OF VAPORIZATION CORRECTIONS
75  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
76  static void no_ip_hltalt (const amrex::Real& t,
77  amrex::Real& hltalt) {
78  hltalt = lat_vap;
79  // Account for change of lat_vap with t above freezing where
80  // constant slope is given by -2369 j/(kg c) = cpv - cw
81  if (t >= tmelt) hltalt = hltalt - 2369.0*(t-tmelt);
82  }
83 
84  // Calculate latent heat of vaporization of water at a given
85  // temperature, taking into account the ice phase if temperature
86  // is below freezing.
87  // Optional argument also calculates a term used to calculate
88  // d(es)/dT within the water-ice transition range.
89  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
90  static void calc_hltalt (const amrex::Real& t,
91  amrex::Real& hltalt,
92  amrex::Real* tterm = nullptr) {
93  amrex::Real tc, weight;
94  no_ip_hltalt(t, hltalt);
95  if (t < tmelt) {
96  // Weighting of hlat accounts for transition from water to ice.
97  tc = t - tmelt;
98  if (tc >= -ttrice) {
99  weight = -tc/ttrice;
100 
101  // polynomial expression approximates difference between es
102  // over water and es over ice from 0 to -ttrice (C) (max of
103  // ttrice is 40): required for accurate estimate of es
104  // derivative in transition range from ice to water
105 
106  if (tterm)
107  {
108  for(auto i = npcf-1; i > 0; --i) *tterm = pcf[i] + tc*(*tterm);
109  *tterm = *tterm/ttrice;
110  }
111  }
112  else {
113  weight = 1.0;
114  }
115 
116  hltalt = hltalt + weight*lat_ice;
117  }
118  }
119 
120  // Temperature derivative outputs, for qsat_*
121  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
122  static void deriv_outputs (const amrex::Real& t,
123  const amrex::Real& p,
124  const amrex::Real& es,
125  const amrex::Real& qs,
126  const amrex::Real& hltalt,
127  const amrex::Real& tterm,
128  amrex::Real& gam,
129  amrex::Real& dqsdt) {
130 
131  // Local variables
132  amrex::Real desdt; // d(es)/dt
133  amrex::Real dqsdt_loc; // local copy of dqsdt
134 
135  if (qs == 1.0) {
136  dqsdt_loc = 0.;
137  }
138  else {
139  desdt = hltalt*es/(R_v*t*t) + tterm;
140  dqsdt_loc = qs*p*desdt/(es*(p-omeps*es));
141  }
142 
143  dqsdt = dqsdt_loc;
144  gam = dqsdt_loc * (hltalt/Cp_d);
145  }
146 
147  // Look up and return saturation vapor pressure from precomputed
148  // table, then calculate and return saturation specific humidity.
149  // Optionally return various temperature derivatives or enthalpy
150  // at saturation.
151  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
152  static void qsat (const amrex::Real& t,
153  const amrex::Real& p,
154  amrex::Real& es,
155  amrex::Real& qs,
156  amrex::Real* gam = nullptr,
157  amrex::Real* dqsdt = nullptr,
158  amrex::Real* enthalpy = nullptr) {
159  // Local variables
160  amrex::Real hltalt; // Modified latent heat for T derivatives
161  amrex::Real tterm = 0.0; // Account for d(es)/dT in transition region
162 
163  es = svp_trans(t); //estblf(t);
164 
165  qs = SatMethods::wv_sat_svp_to_qsat(es, p);
166 
167  // Ensures returned es is consistent with limiters on qs.
168  es = std::min(es, p);
169 
170  // Calculate optional arguments.
171  // "generalized" analytic expression for t derivative of es
172  // accurate to within 1 percent for 173.16 < t < 373.16
173  calc_hltalt(t, hltalt, &tterm);
174 
175  if (enthalpy)
176  {
177  *enthalpy = tq_enthalpy(t, qs, hltalt);
178  }
179 
180  if (gam && dqsdt)
181  {
182  deriv_outputs(t, p, es, qs, hltalt, tterm, *gam, *dqsdt);
183  }
184  }
185 
186  // Calculate SVP over water at a given temperature, and then
187  // calculate and return saturation specific humidity.
188  // Optionally return various temperature derivatives or enthalpy
189  // at saturation.
190  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
191  static void qsat_water (const amrex::Real& t,
192  const amrex::Real& p,
193  amrex::Real& es,
194  amrex::Real& qs,
195  amrex::Real* gam = nullptr,
196  amrex::Real* dqsdt = nullptr,
197  amrex::Real* enthalpy = nullptr) {
198  // Local variables
199  amrex::Real hltalt; // Modified latent heat for T derivatives
200 
201  SatMethods::wv_sat_qsat_water(t, p, es, qs);
202 
203  // "generalized" analytic expression for t derivative of es
204  // accurate to within 1 percent for 173.16 < t < 373.16
205 
206  no_ip_hltalt(t, hltalt);
207 
208  if (enthalpy)
209  {
210  *enthalpy = tq_enthalpy(t, qs, hltalt);
211  }
212 
213  if (gam && dqsdt)
214  {
215  // For pure water/ice transition term is 0.
216  deriv_outputs(t, p, es, qs, hltalt, 0., *gam, *dqsdt);
217  }
218  }
219 
220  // Calculate SVP over ice at a given temperature, and then
221  // calculate and return saturation specific humidity.
222  // Optionally return various temperature derivatives or enthalpy
223  // at saturation.
224  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
225  static void qsat_ice (const amrex::Real& t,
226  const amrex::Real& p,
227  amrex::Real& es,
228  amrex::Real& qs,
229  amrex::Real& gam,
230  amrex::Real& dqsdt,
231  amrex::Real& enthalpy) {
232  // Local variables
233  amrex::Real hltalt; // Modified latent heat for T derivatives
234 
235  SatMethods::wv_sat_qsat_ice(t, p, es, qs);
236 
237  // For pure ice, just add latent heats.
238  hltalt = lat_vap + lat_ice;
239 
240  enthalpy = tq_enthalpy(t, qs, hltalt);
241 
242  // For pure water/ice transition term is 0.
243  deriv_outputs(t, p, es, qs, hltalt, 0., gam, dqsdt);
244  }
245 
246  // find the wet bulb temperature for a given t and q
247  // in a longitude height section
248  // wet bulb temp is the temperature and spec humidity that is
249  // just saturated and has the same enthalpy
250  // if q > qs(t) then tsp > t and qsp = qs(tsp) < q
251  // if q < qs(t) then tsp < t and qsp = qs(tsp) > q
252  //
253  // Method:
254  // a Newton method is used
255  // first guess uses an algorithm provided by John Petch from the UKMO
256  // we exclude points where the physical situation is unrealistic
257  // e.g. where the temperature is outside the range of validity for the
258  // saturation vapor pressure, or where the water vapor pressure
259  // exceeds the ambient pressure, or the saturation specific humidity is
260  // unrealistic
261  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
262  static void findsp (const amrex::Real& q,
263  const amrex::Real& t,
264  const amrex::Real& p,
265  const bool& use_ice,
266  amrex::Real& tsp,
267  amrex::Real& qsp,
268  int& status) {
269 
270  //
271  // local variables
272  //
273  int iter = 8; // max number of times to iterate the calculation
274 
275  amrex::Real es; // sat. vapor pressure
276  amrex::Real gam; // change in sat spec. hum. wrt temperature (times hltalt/cpair)
277  amrex::Real dgdt; // work variable
278  amrex::Real g; // work variable
279  amrex::Real hltalt; // lat. heat. of vap.
280  amrex::Real qs; // spec. hum. of water vapor
281 
282  // work variables
283  amrex::Real t1, q1, dt, dq;
284  amrex::Real qvd;
285  amrex::Real r1b, c1, c2, c3;
286  constexpr amrex::Real dttol = 1.e-4; // the relative temp error tolerance required to quit the iteration
287  constexpr amrex::Real dqtol = 1.e-4; // the relative moisture error tolerance required to quit the iteration
288  amrex::Real enin, enout;
289 
290  c3 = 287.04*(7.5*log(10.))/Cp_d;
291 
292  // Saturation specific humidity at this temperature
293  if (use_ice) {
294  qsat(t, p, es, qs);
295  }
296  else {
297  qsat_water(t, p, es, qs);
298  }
299 
300  // make sure a meaningful calculation is possible
301  if (p <= 5.*es || qs <= 0. || qs >= 0.5 || t < tmin || t > tmax) {
302  status = 1;
303  // Keep initial parameters when conditions aren't suitable
304  tsp = t;
305  qsp = q;
306  enin = 1.;
307  enout = 1.;
308  return;
309  }
310 
311  // Prepare to iterate
312  status = 2;
313 
314  // Get initial enthalpy
315  if (use_ice)
316  calc_hltalt(t,hltalt);
317  else
318  no_ip_hltalt(t,hltalt);
319 
320  enin = tq_enthalpy(t, q, hltalt);
321 
322  // make a guess at the wet bulb temp using a UKMO algorithm (from J. Petch)
323  c1 = hltalt*c3;
324  c2 = std::pow(t + 36., 2);
325  r1b = c2/(c2 + c1*qs);
326  qvd = r1b * (q - qs);
327  tsp = t + ((hltalt/Cp_d)*qvd);
328 
329  // Generate qsp, gam, and enout from tsp.
330  if (use_ice)
331  qsat(tsp, p, es, qsp, &gam, &enout);
332  else
333  qsat_water(tsp, p, es, qsp, &gam, &enout);
334 
335  // iterate on first guess
336  for(auto l = 1; l < iter; ++l) {
337 
338  g = enin - enout;
339  dgdt = -Cp_d * (1 + gam);
340 
341  // New tsp
342  t1 = tsp - g/dgdt;
343  dt = abs(t1 - tsp)/t1;
344  tsp = t1;
345 
346  // bail out if past end of temperature range
347  if ( tsp < tmin ) {
348  tsp = tmin;
349  // Get latent heat and set qsp to a value
350  // that preserves enthalpy.
351  if (use_ice)
352  calc_hltalt(tsp,hltalt);
353  else
354  no_ip_hltalt(tsp,hltalt);
355 
356  qsp = (enin - Cp_d*tsp)/hltalt;
357  enout = tq_enthalpy(tsp, qsp, hltalt);
358  status = 4;
359  break;
360  }
361 
362  // Re-generate qsp, gam, and enout from new tsp.
363  if (use_ice)
364  qsat(tsp, p, es, q1, &gam, &enout);
365  else
366  qsat_water(tsp, p, es, q1, &gam, &enout);
367 
368  dq = abs(q1 - qsp)/std::max(q1,1.e-12);
369  qsp = q1;
370 
371  // if converged at this point, exclude it from more iterations
372  if (dt < dttol && dq < dqtol) {
373  status = 0;
374  break;
375  }
376  }
377  // Test for enthalpy conservation
378  if (abs((enin-enout)/(enin+enout)) > 1.e-4) status = 8;
379  return;
380  }
381 };
382 #endif // ERF_WATER_VAPOR_SATURATION_H_
constexpr amrex::Real R_v
Definition: ERF_Constants.H:11
constexpr amrex::Real lat_vap
Definition: ERF_Constants.H:85
constexpr amrex::Real Cp_d
Definition: ERF_Constants.H:12
constexpr amrex::Real tmelt
Definition: ERF_Constants.H:88
constexpr amrex::Real ttrice
Definition: ERF_Constants.H:91
constexpr amrex::Real lat_ice
Definition: ERF_Constants.H:86
constexpr amrex::Real omeps
Definition: ERF_Constants.H:93
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real wv_sat_svp_water(const amrex::Real &t, const int idx=1)
Definition: ERF_Sat_methods.H:97
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real wv_sat_svp_to_qsat(const amrex::Real &es, const amrex::Real &p)
Definition: ERF_Sat_methods.H:42
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE void wv_sat_qsat_water(const amrex::Real &t, const amrex::Real &p, amrex::Real &es, amrex::Real &qs, const int idx=1)
Definition: ERF_Sat_methods.H:51
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE void wv_sat_qsat_ice(const amrex::Real &t, const amrex::Real &p, amrex::Real &es, amrex::Real &qs, const int idx=1)
Definition: ERF_Sat_methods.H:66
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real wv_sat_svp_trans(const amrex::Real &t, const int idx=1)
Definition: ERF_Sat_methods.H:131
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real wv_sat_svp_ice(const amrex::Real &t, const int idx=1)
Definition: ERF_Sat_methods.H:114
Definition: ERF_WaterVaporSaturation.H:28
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real svp_ice(const amrex::Real &t)
Definition: ERF_WaterVaporSaturation.H:54
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE void findsp(const amrex::Real &q, const amrex::Real &t, const amrex::Real &p, const bool &use_ice, amrex::Real &tsp, amrex::Real &qsp, int &status)
Definition: ERF_WaterVaporSaturation.H:262
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real svp_water(const amrex::Real &t)
Definition: ERF_WaterVaporSaturation.H:48
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE void qsat_ice(const amrex::Real &t, const amrex::Real &p, amrex::Real &es, amrex::Real &qs, amrex::Real &gam, amrex::Real &dqsdt, amrex::Real &enthalpy)
Definition: ERF_WaterVaporSaturation.H:225
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real svp_trans(const amrex::Real &t)
Definition: ERF_WaterVaporSaturation.H:60
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE void calc_hltalt(const amrex::Real &t, amrex::Real &hltalt, amrex::Real *tterm=nullptr)
Definition: ERF_WaterVaporSaturation.H:90
static constexpr amrex::Real tmax
Definition: ERF_WaterVaporSaturation.H:33
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE void qsat(const amrex::Real &t, const amrex::Real &p, amrex::Real &es, amrex::Real &qs, amrex::Real *gam=nullptr, amrex::Real *dqsdt=nullptr, amrex::Real *enthalpy=nullptr)
Definition: ERF_WaterVaporSaturation.H:152
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE void deriv_outputs(const amrex::Real &t, const amrex::Real &p, const amrex::Real &es, const amrex::Real &qs, const amrex::Real &hltalt, const amrex::Real &tterm, amrex::Real &gam, amrex::Real &dqsdt)
Definition: ERF_WaterVaporSaturation.H:122
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE void no_ip_hltalt(const amrex::Real &t, amrex::Real &hltalt)
Definition: ERF_WaterVaporSaturation.H:76
static constexpr int npcf
Definition: ERF_WaterVaporSaturation.H:39
static constexpr amrex::Real tmin
Definition: ERF_WaterVaporSaturation.H:32
static constexpr const amrex::Real pcf[npcf]
Definition: ERF_WaterVaporSaturation.H:40
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE void qsat_water(const amrex::Real &t, const amrex::Real &p, amrex::Real &es, amrex::Real &qs, amrex::Real *gam=nullptr, amrex::Real *dqsdt=nullptr, amrex::Real *enthalpy=nullptr)
Definition: ERF_WaterVaporSaturation.H:191
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real tq_enthalpy(const amrex::Real &t, const amrex::Real &q, const amrex::Real &hltalt)
Definition: ERF_WaterVaporSaturation.H:67