PISM, A Parallel Ice Sheet Model  stable v2.1-1-g6902d5502 committed by Ed Bueler on 2023-12-20 08:38:27 -0800
EnthalpyConverter.cc
Go to the documentation of this file.
1 // Copyright (C) 2009-2017, 2021, 2023 Andreas Aschwanden, Ed Bueler and Constantine Khroulev
2 //
3 // This file is part of PISM.
4 //
5 // PISM is free software; you can redistribute it and/or modify it under the
6 // terms of the GNU General Public License as published by the Free Software
7 // Foundation; either version 3 of the License, or (at your option) any later
8 // version.
9 //
10 // PISM is distributed in the hope that it will be useful, but WITHOUT ANY
11 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
12 // FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
13 // details.
14 //
15 // You should have received a copy of the GNU General Public License
16 // along with PISM; if not, write to the Free Software
17 // Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
18 
19 #include "pism/util/EnthalpyConverter.hh"
20 #include "pism/util/ConfigInterface.hh"
21 
22 #include "pism/util/error_handling.hh"
23 
24 namespace pism {
25 
26 /*! @class EnthalpyConverter
27 
28  Maps from @f$ (H,p) @f$ to @f$ (T,\omega,p) @f$ and back.
29 
30  Requirements:
31 
32  1. A converter has to implement an invertible map @f$ (H,p) \to (T,
33  \omega, p) @f$ *and* its inverse. Both the forward map and the
34  inverse need to be defined for all permissible values of @f$
35  (H,p) @f$ and @f$ (T, \omega, p) @f$, respectively.
36 
37  2. A converter has to be consistent with laws and parameterizations
38  used elsewhere in the model. This includes models coupled to
39  PISM.
40 
41  3. For a fixed volume of liquid (or solid) water and given two
42  energy states @f$ (H_1, p_1) @f$ and @f$ (H_2, p_2) @f$ , let @f$
43  \Delta U_H @f$ be the difference in internal energy of this
44  volume between the two states *computed using enthalpy*. We require
45  that @f$ \Delta U_T = \Delta U_H @f$ , where @f$ \Delta U_T @f$
46  is the difference in internal energy *computed using corresponding*
47  @f$ (T_1, \omega_1, p_1) @f$ *and* @f$ (T_2, \omega_2, p_2) @f$.
48 
49  4. We assume that ice and water are incompressible, so a change in
50  pressure does no work, and @f$ \Diff{H}{p} = 0 @f$. In addition
51  to this, for cold ice and liquid water @f$ \Diff{T}{p} = 0 @f$.
52 */
53 
55  m_p_air = config.get_number("surface.pressure"); // Pa
56  m_g = config.get_number("constants.standard_gravity"); // m s-2
57  m_beta = config.get_number("constants.ice.beta_Clausius_Clapeyron"); // K Pa-1
58  m_rho_i = config.get_number("constants.ice.density"); // kg m-3
59  m_c_i = config.get_number("constants.ice.specific_heat_capacity"); // J kg-1 K-1
60  m_c_w = config.get_number("constants.fresh_water.specific_heat_capacity"); // J kg-1 K-1
61  m_L = config.get_number("constants.fresh_water.latent_heat_of_fusion"); // J kg-1
62  m_T_melting = config.get_number("constants.fresh_water.melting_point_temperature"); // K
63  m_T_tolerance = config.get_number("enthalpy_converter.relaxed_is_temperate_tolerance"); // K
64  m_T_0 = config.get_number("enthalpy_converter.T_reference"); // K
65 
66  m_do_cold_ice_methods = config.get_flag("energy.temperature_based");
67 }
68 
69 //! Return `true` if ice at `(E, P)` is temperate.
70 //! Determines if E >= E_s(p), that is, if the ice is at the pressure-melting point.
71 bool EnthalpyConverter::is_temperate(double E, double P) const {
73  return is_temperate_relaxed(E, P);
74  }
75 
76  return (E >= enthalpy_cts(P));
77 }
78 
79 //! A relaxed version of `is_temperate()`.
80 /*! Returns `true` if the pressure melting temperature corresponding to `(E, P)` is within
81  `enthalpy_converter.relaxed_is_temperate_tolerance` from `fresh_water.melting_point_temperature`.
82  */
83 bool EnthalpyConverter::is_temperate_relaxed(double E, double P) const {
85 }
86 
87 
88 void EnthalpyConverter::validate_T_omega_P(double T, double omega, double P) const {
89 #if (Pism_DEBUG==1)
90  const double T_melting = melting_temperature(P);
91  if (T <= 0.0) {
92  throw RuntimeError::formatted(PISM_ERROR_LOCATION, "T = %f <= 0 is not a valid absolute temperature",T);
93  }
94  const double eps = 1.0e-6;
95  if ((omega < 0.0 - eps) || (1.0 + eps < omega)) {
96  throw RuntimeError::formatted(PISM_ERROR_LOCATION, "water fraction omega=%f not in range [0,1]",omega);
97  }
98  if (T > T_melting + eps) {
99  throw RuntimeError::formatted(PISM_ERROR_LOCATION, "T=%f exceeds T_melting=%f; not allowed",T,T_melting);
100  }
101  if ((T < T_melting - eps) && (omega > 0.0 + eps)) {
102  throw RuntimeError::formatted(PISM_ERROR_LOCATION, "T < T_melting AND omega > 0 is contradictory;"
103  " got T=%f, T_melting=%f, omega=%f",
104  T, T_melting, omega);
105  }
106 #else
107  (void) T;
108  (void) omega;
109  (void) P;
110 #endif
111 }
112 
113 void EnthalpyConverter::validate_E_P(double E, double P) const {
114 #if (Pism_DEBUG==1)
115  if (E >= enthalpy_liquid(P)) {
116  throw RuntimeError::formatted(PISM_ERROR_LOCATION, "E=%f J/kg at P=%f Pa equals or exceeds that of liquid water (%f J/kg)",
117  E, P, enthalpy_liquid(P));
118  }
119 #else
120  (void) E;
121  (void) P;
122 #endif
123 }
124 
125 
126 //! Get pressure in ice from depth below surface using the hydrostatic assumption.
127 /*! If \f$d\f$ is the depth then
128  \f[ p = p_{\text{air}} + \rho_i g d. \f]
129 Frequently \f$d\f$ is computed from the thickess minus a level in the ice,
130 something like `ice_thickness(i, j) - z[k]`. The input depth to this routine is allowed to
131 be negative, representing a position above the surface of the ice.
132  */
133 double EnthalpyConverter::pressure(double depth) const {
134  if (depth >= 0.0) {
135  return m_p_air + m_rho_i * m_g * depth;
136  }
137 
138  return m_p_air; // at or above surface of ice
139 }
140 
141 //! Compute pressure in a column of ice. Does not check validity of `depth`.
142 void EnthalpyConverter::pressure(const std::vector<double> &depth,
143  unsigned int ks,
144  std::vector<double> &result) const {
145  for (unsigned int k = 0; k <= ks; ++k) {
146  result[k] = m_p_air + m_rho_i * m_g * depth[k];
147  }
148 }
149 
150 //! Get melting temperature from pressure p.
151 /*!
152  \f[ T_m(p) = T_{melting} - \beta p. \f]
153  */
155  return m_T_melting - m_beta * P;
156 }
157 
158 
159 //! @brief Compute the maximum allowed value of ice enthalpy
160 //! (corresponds to @f$ \omega = 1 @f$).
161 double EnthalpyConverter::enthalpy_liquid(double P) const {
162  return enthalpy_cts(P) + L(melting_temperature(P));
163 }
164 
165 
166 //! Get absolute (not pressure-adjusted) ice temperature (K) from enthalpy and pressure.
167 /*! From \ref AschwandenBuelerKhroulevBlatter,
168  \f[ T= T(E,p) = \begin{cases}
169  c_i^{-1} E + T_0, & E < E_s(p), \\
170  T_m(p), & E_s(p) \le E < E_l(p).
171  \end{cases} \f]
172 
173 We do not allow liquid water (%i.e. water fraction \f$\omega=1.0\f$) so we
174 throw an exception if \f$E \ge E_l(p)\f$.
175  */
176 double EnthalpyConverter::temperature(double E, double P) const {
177  validate_E_P(E, P);
178 
179  if (E < enthalpy_cts(P)) {
180  return temperature_cold(E);
181  }
182 
183  return melting_temperature(P);
184 }
185 
186 
187 //! Get pressure-adjusted ice temperature, in Kelvin, from enthalpy and pressure.
188 /*!
189 The pressure-adjusted temperature is:
190  \f[ T_{pa}(E,p) = T(E,p) - T_m(p) + T_{melting}. \f]
191  */
192 double EnthalpyConverter::pressure_adjusted_temperature(double E, double P) const {
193  return temperature(E, P) - melting_temperature(P) + m_T_melting;
194 }
195 
196 
197 //! Get liquid water fraction from enthalpy and pressure.
198 /*!
199  From [@ref AschwandenBuelerKhroulevBlatter],
200  @f[
201  \omega(E,p) =
202  \begin{cases}
203  0.0, & E \le E_s(p), \\
204  (E-E_s(p)) / L, & E_s(p) < E < E_l(p).
205  \end{cases}
206  @f]
207 
208  We do not allow liquid water (i.e. water fraction @f$ \omega=1.0 @f$).
209  */
210 double EnthalpyConverter::water_fraction(double E, double P) const {
211  validate_E_P(E, P);
212 
213  double E_s = enthalpy_cts(P);
214  if (E <= E_s) {
215  return 0.0;
216  }
217 
218  return (E - E_s) / L(melting_temperature(P));
219 }
220 
221 
222 //! Compute enthalpy from absolute temperature, liquid water fraction, and pressure.
223 /*! This is an inverse function to the functions \f$T(E,p)\f$ and
224 \f$\omega(E,p)\f$ [\ref AschwandenBuelerKhroulevBlatter]. It returns:
225  \f[E(T,\omega,p) =
226  \begin{cases}
227  c_i (T - T_0), & T < T_m(p) \quad\text{and}\quad \omega = 0, \\
228  E_s(p) + \omega L, & T = T_m(p) \quad\text{and}\quad \omega \ge 0.
229  \end{cases} \f]
230 Certain cases are not allowed and throw exceptions:
231 - \f$T<=0\f$ (error code 1)
232 - \f$\omega < 0\f$ or \f$\omega > 1\f$ (error code 2)
233 - \f$T>T_m(p)\f$ (error code 3)
234 - \f$T<T_m(p)\f$ and \f$\omega > 0\f$ (error code 4)
235 These inequalities may be violated in the sixth digit or so, however.
236  */
237 double EnthalpyConverter::enthalpy(double T, double omega, double P) const {
238  validate_T_omega_P(T, omega, P);
239 
240  const double T_melting = melting_temperature(P);
241 
242  if (T < T_melting) {
243  return enthalpy_cold(T);
244  }
245 
246  return enthalpy_cts(P) + omega * L(T_melting);
247 }
248 
249 
250 //! Compute enthalpy more permissively than enthalpy().
251 /*! Computes enthalpy from absolute temperature, liquid water fraction, and
252 pressure. Use this form of enthalpy() when outside sources (e.g. information
253 from a coupler) might generate a temperature above the pressure melting point or
254 generate cold ice with a positive water fraction.
255 
256 Treats temperatures above pressure-melting point as \e at the pressure-melting
257 point. Interprets contradictory case of \f$T < T_m(p)\f$ and \f$\omega > 0\f$
258 as cold ice, ignoring the water fraction (\f$\omega\f$) value.
259 
260 Calls enthalpy(), which validates its inputs.
261 
262 Computes:
263  @f[
264  E =
265  \begin{cases}
266  E(T,0.0,p), & T < T_m(p) \quad \text{and} \quad \omega \ge 0, \\
267  E(T_m(p),\omega,p), & T \ge T_m(p) \quad \text{and} \quad \omega \ge 0,
268  \end{cases}
269  @f]
270  but ensures @f$ 0 \le \omega \le 1 @f$ in second case. Calls enthalpy() for
271  @f$ E(T,\omega,p) @f$.
272  */
273 double EnthalpyConverter::enthalpy_permissive(double T, double omega, double P) const {
274  const double T_m = melting_temperature(P);
275 
276  if (T < T_m) {
277  return enthalpy(T, 0.0, P);
278  }
279 
280  // T >= T_m(P) replaced with T = T_m(P)
281  return enthalpy(T_m, std::max(0.0, std::min(omega, 1.0)), P);
282 }
283 
285  : EnthalpyConverter(config) {
286  // turn on the "cold" enthalpy converter mode
287  m_do_cold_ice_methods = true;
288  // set melting temperature to one million Kelvin so that all ice is cold
289  m_T_melting = 1e6;
290  // disable pressure-dependence of the melting temperature by setting Clausius-Clapeyron beta to
291  // zero
292  m_beta = 0.0;
293 }
294 
295 //! Latent heat of fusion of water as a function of pressure melting
296 //! temperature.
297 /*!
298 
299  Following a re-interpretation of [@ref
300  AschwandenBuelerKhroulevBlatter], we require that @f$ \Diff{H}{p} =
301  0 @f$:
302 
303  @f[
304  \Diff{H}{p} = \diff{H_w}{p} + \diff{H_w}{p}\Diff{T}{p}
305  @f]
306 
307  We assume that water is incompressible, so @f$ \Diff{T}{p} = 0 @f$
308  and the second term vanishes.
309 
310  As for the first term, equation (5) of [@ref
311  AschwandenBuelerKhroulevBlatter] defines @f$ H_w @f$ as follows:
312 
313  @f[
314  H_w = \int_{T_0}^{T_m(p)} C_i(t) dt + L + \int_{T_m(p)}^T C_w(t)dt
315  @f]
316 
317  Using the fundamental theorem of Calculus, we get
318  @f[
319  \diff{H_w}{p} = (C_i(T_m(p)) - C_w(T_m(p))) \diff{T_m(p)}{p} + \diff{L}{p}
320  @f]
321 
322  Assuming that @f$ C_i(T) = c_i @f$ and @f$ C_w(T) = c_w @f$ (i.e. specific heat
323  capacities of ice and water do not depend on temperature) and using
324  the Clausius-Clapeyron relation
325  @f[
326  T_m(p) = T_m(p_{\text{air}}) - \beta p,
327  @f]
328 
329  we get
330  @f{align}{
331  \Diff{H}{p} &= (c_i - c_w)\diff{T_m(p)}{p} + \diff{L}{p}\\
332  &= \beta(c_w - c_i) + \diff{L}{p}\\
333  @f}
334  Requiring @f$ \Diff{H}{p} = 0 @f$ implies
335  @f[
336  \diff{L}{p} = -\beta(c_w - c_i),
337  @f]
338  and so
339  @f{align}{
340  L(p) &= -\beta p (c_w - c_i) + C\\
341  &= (T_m(p) - T_m(p_{\text{air}})) (c_w - c_i) + C.
342  @f}
343 
344  Letting @f$ p = p_{\text{air}} @f$ we find @f$ C = L(p_\text{air}) = L_0 @f$, so
345  @f[
346  L(p) = (T_m(p) - T_m(p_{\text{air}})) (c_w - c_i) + L_0,
347  @f]
348  where @f$ L_0 @f$ is the latent heat of fusion of water at atmospheric
349  pressure.
350 
351  Therefore a consistent interpretation of [@ref
352  AschwandenBuelerKhroulevBlatter] requires the temperature-dependent
353  approximation of the latent heat of fusion of water given above.
354 
355  Note that this form of @f$ L(p) @f$ also follows from Kirchhoff's
356  law of thermochemistry.
357 */
358 double EnthalpyConverter::L(double T_pm) const {
359  return m_L + (m_c_w - m_c_i) * (T_pm - 273.15);
360 }
361 
362 //! Specific heat capacity of ice.
363 double EnthalpyConverter::c() const {
364  return m_c_i;
365 }
366 
367 //! Get enthalpy E_s(p) at cold-temperate transition point from pressure p.
368 /*! Returns
369  \f[ E_s(p) = c_i (T_m(p) - T_0), \f]
370  */
371 double EnthalpyConverter::enthalpy_cts(double P) const {
372  return m_c_i * (melting_temperature(P) - m_T_0);
373 }
374 
375 //! Convert temperature into enthalpy (cold case).
376 double EnthalpyConverter::enthalpy_cold(double T) const {
377  return m_c_i * (T - m_T_0);
378 }
379 
380 //! Convert enthalpy into temperature (cold case).
381 double EnthalpyConverter::temperature_cold(double E) const {
382  return (E / m_c_i) + m_T_0;
383 }
384 
385 } // end of namespace pism
ColdEnthalpyConverter(const Config &config)
double get_number(const std::string &name, UseFlag flag=REMEMBER_THIS_USE) const
bool get_flag(const std::string &name, UseFlag flag=REMEMBER_THIS_USE) const
A class for storing and accessing PISM configuration flags and parameters.
double enthalpy(double T, double omega, double P) const
Compute enthalpy from absolute temperature, liquid water fraction, and pressure.
EnthalpyConverter(const Config &config)
void validate_E_P(double E, double P) const
double m_rho_i
density of ice
bool is_temperate(double E, double P) const
double melting_temperature(double P) const
Get melting temperature from pressure p.
double m_T_tolerance
temperature tolerance used in is_temperate() in cold ice mode
double temperature(double E, double P) const
Get absolute (not pressure-adjusted) ice temperature (K) from enthalpy and pressure.
double m_beta
beta in the Clausius-Clapeyron relation ( ).
bool m_do_cold_ice_methods
if cold ice methods are selected, use is_temperate() check based on temperature, not enthalpy
double enthalpy_cold(double T) const
Convert temperature into enthalpy (cold case).
double enthalpy_liquid(double P) const
Compute the maximum allowed value of ice enthalpy (corresponds to ).
double m_c_i
specific heat capacity of ice
double c() const
Specific heat capacity of ice.
double m_g
acceleration due to gravity
double enthalpy_permissive(double T, double omega, double P) const
Compute enthalpy more permissively than enthalpy().
double m_L
latent heat of fusion of water at atmospheric pressure
double pressure_adjusted_temperature(double E, double P) const
Get pressure-adjusted ice temperature, in Kelvin, from enthalpy and pressure.
double enthalpy_cts(double P) const
Get enthalpy E_s(p) at cold-temperate transition point from pressure p.
bool is_temperate_relaxed(double E, double P) const
A relaxed version of is_temperate().
double temperature_cold(double E) const
Convert enthalpy into temperature (cold case).
double pressure(double depth) const
Get pressure in ice from depth below surface using the hydrostatic assumption.
double m_T_melting
melting temperature of pure water at atmospheric pressure
double water_fraction(double E, double P) const
Get liquid water fraction from enthalpy and pressure.
double m_c_w
specific heat capacity of pure water
double m_p_air
atmospheric pressure
double L(double T_pm) const
void validate_T_omega_P(double T, double omega, double P) const
double m_T_0
reference temperature in the definition of ice enthalpy
Converts between specific enthalpy and temperature or liquid content.
static RuntimeError formatted(const ErrorLocation &location, const char format[],...) __attribute__((format(printf
build a RuntimeError with a formatted message
#define PISM_ERROR_LOCATION
double max(const array::Scalar &input)
Finds maximum over all the values in an array::Scalar object. Ignores ghosts.
Definition: Scalar.cc:165
double min(const array::Scalar &input)
Finds minimum over all the values in an array::Scalar object. Ignores ghosts.
Definition: Scalar.cc:193
static const double k
Definition: exactTestP.cc:42