PISM, A Parallel Ice Sheet Model stable 0.4.1779
Public Member Functions | Protected Attributes

EnthalpyConverter Class Reference

Converts between specific enthalpy and temperature or liquid content. More...

#include <enthalpyConverter.hh>

Inheritance diagram for EnthalpyConverter:
Inheritance graph
[legend]

List of all members.

Public Member Functions

 EnthalpyConverter (const NCConfigVariable &config)
virtual ~EnthalpyConverter ()
virtual PetscErrorCode viewConstants (PetscViewer viewer) const
 Simple view of state of EnthalpyConverter. viewer==NULL sends to stdout.
virtual double getPressureFromDepth (double depth) const
 Get pressure in ice from depth below surface using the hydrostatic assumption.
virtual double getMeltingTemp (double p) const
 Get melting temperature from pressure p.
virtual double getEnthalpyCTS (double p) const
 Get enthalpy E_s(p) at cold-temperate transition point from pressure p.
virtual PetscErrorCode getEnthalpyInterval (double p, double &E_s, double &E_l) const
 Get enthalpies E_s(p) and E_l(p) (endpoints of temperate ice enthalpy range) from pressure p.
virtual double getCTS (double E, double p) const
 Computes the ratio CTS = E / E_s(p). The cold-temperate transition surface (CTS) is the level surface CTS = 1.
virtual bool isTemperate (double E, double p) const
 Determines if E >= E_s(p), that is, if the ice is at the pressure-melting point.
virtual bool isLiquified (double E, double p) const
 Is the ice mixture actually liquid water?
virtual PetscErrorCode getAbsTemp (double E, double p, double &T) const
 Get absolute (not pressure-adjusted) ice temperature (K) from enthalpy and pressure.
virtual PetscErrorCode getPATemp (double E, double p, double &T_pa) const
 Get pressure-adjusted ice temperature, in Kelvin, from enthalpy and pressure.
virtual PetscErrorCode getWaterFraction (double E, double p, double &omega) const
 Get liquid water fraction from enthalpy and pressure.
virtual PetscErrorCode getEnth (double T, double omega, double p, double &E) const
 Compute enthalpy from absolute temperature, liquid water fraction, and pressure.
virtual PetscErrorCode getEnthPermissive (double T, double omega, double p, double &E) const
 Compute enthalpy more permissively than getEnth().
virtual PetscErrorCode getEnthAtWaterFraction (double omega, double p, double &E) const
 Returns enthalpy for temperate ice with a given liquid fraction.

Protected Attributes

double T_melting
double L
double c_i
double rho_i
double g
double p_air
double beta
double T_tol
double T_0
bool do_cold_ice_methods

Detailed Description

Converts between specific enthalpy and temperature or liquid content.

Use this way, for example within IceModel with NCConfigVariable config member:

  #include "enthalpyConverter.hh"
  EnthalpyConverter EC(&config);  // runs constructor; do after initialization of NCConfigVariable config
  ...
  for (...) {
    ...
    E_s = EC.getEnthalpyCTS(p);
    ... etc ...
  }   

Some methods are functions which return the computed values or a boolean. These do no error checking. Others return PetscErrorCode and the modify pass-by-reference arguments. These check either that the enthalpy is below that of liquid water, or they check that the temperature (in K) is positive.

Specifically, getAbsTemp() gives return value 1 if the input enthalpy exceeded that of liquid water, but it puts the temperature of temperate ice into its computed value for T. Method getWaterFraction() gives return value of 1 under the same condition, but puts the maximum-liquid-content into its computed value for omega.

The three methods that get the enthalpy from temperatures and liquid fractions, namely getEnth(), getEnthPermissive(), getEnthAtWaterFraction(), are more strict about error checking. They call SETERRQ() if their arguments are invalid.

This class is documented by [AschwandenBuelerKhroulevBlatter].

Definition at line 55 of file enthalpyConverter.hh.


Constructor & Destructor Documentation

EnthalpyConverter ( const NCConfigVariable &  config)

Definition at line 24 of file enthalpyConverter.cc.

References beta, c_i, do_cold_ice_methods, g, L, p_air, rho_i, T_0, T_melting, and T_tol.

virtual ~EnthalpyConverter ( ) [inline, virtual]

Definition at line 58 of file enthalpyConverter.hh.


Member Function Documentation

PetscErrorCode getAbsTemp ( double  E,
double  p,
double &  T 
) const [virtual]

Get absolute (not pressure-adjusted) ice temperature (K) from enthalpy and pressure.

From AschwandenBuelerKhroulevBlatter,

\[ T= T(E,p) = \begin{cases} c_i^{-1} E + T_0, & E < E_s(p), \\ T_m(p), & E_s(p) \le E < E_l(p). \end{cases} \]

We do not allow liquid water (i.e. water fraction $\omega=1.0$) so we return an error code of 1 if $E \ge E_l(p)$.

Reimplemented in ICMEnthalpyConverter.

Definition at line 162 of file enthalpyConverter.cc.

References c_i, getEnthalpyInterval(), getMeltingTemp(), and T_0.

Referenced by IceModel::energyStats(), IceModel::get_bed_top_temp(), getPATemp(), and SIA_Sliding::update().

double getCTS ( double  E,
double  p 
) const [virtual]

Computes the ratio CTS = E / E_s(p). The cold-temperate transition surface (CTS) is the level surface CTS = 1.

If $E$ and $E_s$ are the ice mixture enthalpy and the CTS enthalpy, respectively, then

\[ CTS(E,p) = \frac{E}{E_s(p)}.\]

The level set CTS = 1 is the actual CTS. The value computed here is greater than one for temperate ice and less than one for cold ice. The output of this routine is normally only used for postprocessing.

Definition at line 135 of file enthalpyConverter.cc.

References getEnthalpyCTS().

Referenced by IceModel::setCTSFromEnthalpy().

PetscErrorCode getEnth ( double  T,
double  omega,
double  p,
double &  E 
) const [virtual]

Compute enthalpy from absolute temperature, liquid water fraction, and pressure.

This is an inverse function to the functions $T(E,p)$ and $\omega(E,p)$ [AschwandenBuelerKhroulevBlatter]. It returns:

\[E(T,\omega,p) = \begin{cases} c_i (T - T_0), & T < T_m(p) \quad\text{and}\quad \omega = 0, \\ E_s(p) + \omega L, & T = T_m(p) \quad\text{and}\quad \omega \ge 0. \end{cases} \]

Certain cases are not allowed and return errors:

  • $T<=0$ (error code 1)
  • $\omega < 0$ or $\omega > 1$ (error code 2)
  • $T>T_m(p)$ (error code 3)
  • $T<T_m(p)$ and $\omega > 0$ (error code 4) These inequalities may be violated in the sixth digit or so, however.

Reimplemented in ICMEnthalpyConverter.

Definition at line 241 of file enthalpyConverter.cc.

References c_i, e, getEnthalpyCTS(), getMeltingTemp(), L, and T_0.

Referenced by getEnthPermissive(), and main().

double getEnthalpyCTS ( double  p) const [virtual]

Get enthalpy E_s(p) at cold-temperate transition point from pressure p.

Returns

\[ E_s(p) = c_i (T_m(p) - T_0), \]

Definition at line 108 of file enthalpyConverter.cc.

References c_i, getMeltingTemp(), and T_0.

Referenced by IceModel_tempicethk_basal::compute(), getCTS(), getEnth(), IceModel::getEnthalpyCTSColumn(), getEnthalpyInterval(), ICMEnthalpyConverter::getEnthAtWaterFraction(), getEnthAtWaterFraction(), and isTemperate().

PetscErrorCode getEnthalpyInterval ( double  p,
double &  E_s,
double &  E_l 
) const [virtual]

Get enthalpies E_s(p) and E_l(p) (endpoints of temperate ice enthalpy range) from pressure p.

Ice at enthalpy $E$ is temperate if $E_s(p) < E < E_l(p)$:

\[ E_s(p) = c_i (T_m(p) - T_0), \]

\[ E_l(p) = E_s(p) + L. \]

Definition at line 118 of file enthalpyConverter.cc.

References getEnthalpyCTS(), and L.

Referenced by getAbsTemp(), getWaterFraction(), and isLiquified().

PetscErrorCode getEnthAtWaterFraction ( double  omega,
double  p,
double &  E 
) const [virtual]

Returns enthalpy for temperate ice with a given liquid fraction.

Computes

\[E = E_s(p) + \omega L.\]

Only the following case returns an error:

  • $\omega < 0$ or $\omega > 1$ (error code 2) These inequalities may be violated in the sixth digit or so, however.

Reimplemented in ICMEnthalpyConverter.

Definition at line 308 of file enthalpyConverter.cc.

References e, getEnthalpyCTS(), and L.

PetscErrorCode getEnthPermissive ( double  T,
double  omega,
double  p,
double &  E 
) const [virtual]

Compute enthalpy more permissively than getEnth().

Computes enthalpy from absolute temperature, liquid water fraction, and pressure. Use this form of getEnth() when outside sources (e.g. information from a coupler) might generate a temperature above the pressure melting point or generate cold ice with a positive water fraction.

Treats temperatures above pressure-melting point as at the pressure-melting point. Interprets contradictory case of $T < T_m(p)$ and $\omega > 0$ as cold ice, ignoring the water fraction ( $\omega$) value.

Checks if $T <= 0$ and returns error code 1 if so.

Computes:

\[E = \begin{cases} E(T,0.0,p), & T < T_m(p) \quad \text{and} \quad \omega \ge 0, \\ E(T_m(p),\omega,p), & T \ge T_m(p) \quad \text{and} \quad \omega \ge 0, \end{cases} \]

but ensures $0\le \omega \le 1$ in second case. Calls getEnth() for $E(T,\omega,p)$.

Reimplemented in ICMEnthalpyConverter.

Definition at line 285 of file enthalpyConverter.cc.

References getEnth(), and getMeltingTemp().

Referenced by SIA_Sliding::basalVelocitySIA(), IceModel::check_maximum_thickness(), IceModel::compute_enthalpy(), IceModel::compute_enthalpy_cold(), enthalpy_from_temperature_cold(), IceModel::enthalpyAndDrainageStep(), and IceModel::putTempAtDepth().

double getMeltingTemp ( double  p) const [virtual]

Get melting temperature from pressure p.

\[ T_m(p) = T_{melting} - \beta p. \]

Reimplemented in ICMEnthalpyConverter.

Definition at line 99 of file enthalpyConverter.cc.

References beta, and T_melting.

Referenced by IceModel::enthalpyAndDrainageStep(), getAbsTemp(), getEnth(), getEnthalpyCTS(), getEnthPermissive(), getPATemp(), and main().

PetscErrorCode getPATemp ( double  E,
double  p,
double &  T_pa 
) const [virtual]

Get pressure-adjusted ice temperature, in Kelvin, from enthalpy and pressure.

The pressure-adjusted temperature is:

\[ T_{pa}(E,p) = T(E,p) - T_m(p) + T_{melting}. \]

Definition at line 183 of file enthalpyConverter.cc.

References getAbsTemp(), getMeltingTemp(), and T_melting.

Referenced by isTemperate().

double getPressureFromDepth ( double  depth) const [virtual]
PetscErrorCode getWaterFraction ( double  E,
double  p,
double &  omega 
) const [virtual]

Get liquid water fraction from enthalpy and pressure.

From AschwandenBuelerKhroulevBlatter,

\[ \omega(E,p) = \begin{cases} 0.0, & E \le E_s(p), \\ (E-E_s(p)) / L, & E_s(p) < E < E_l(p). \end{cases} \]

We do not allow liquid water (i.e. water fraction $\omega=1.0$) so we return error code 1 if $E \ge E_l(p)$, but we still compute $\omega=1.0$.

Reimplemented in ICMEnthalpyConverter.

Definition at line 202 of file enthalpyConverter.cc.

References getEnthalpyInterval(), and L.

Referenced by IceModel::compute_liquid_water_fraction(), and IceModel::enthalpyAndDrainageStep().

bool isLiquified ( double  E,
double  p 
) const [virtual]

Is the ice mixture actually liquid water?

Definition at line 219 of file enthalpyConverter.cc.

References getEnthalpyInterval().

bool isTemperate ( double  E,
double  p 
) const [virtual]
PetscErrorCode viewConstants ( PetscViewer  viewer) const [virtual]

Simple view of state of EnthalpyConverter. viewer==NULL sends to stdout.

Definition at line 40 of file enthalpyConverter.cc.

References beta, c_i, do_cold_ice_methods, g, L, p_air, rho_i, T_0, T_melting, and T_tol.

Referenced by IceModel::allocate_enthalpy_converter(), and IceModel::enthalpyAndDrainageStep().


Member Data Documentation

double beta [protected]

Definition at line 81 of file enthalpyConverter.hh.

Referenced by EnthalpyConverter(), getMeltingTemp(), and viewConstants().

double c_i [protected]
double g [protected]

Definition at line 81 of file enthalpyConverter.hh.

Referenced by EnthalpyConverter(), getPressureFromDepth(), and viewConstants().

double L [protected]
double p_air [protected]

Definition at line 81 of file enthalpyConverter.hh.

Referenced by EnthalpyConverter(), getPressureFromDepth(), and viewConstants().

double rho_i [protected]

Definition at line 81 of file enthalpyConverter.hh.

Referenced by EnthalpyConverter(), getPressureFromDepth(), and viewConstants().

double T_0 [protected]
double T_melting [protected]
double T_tol [protected]

Definition at line 81 of file enthalpyConverter.hh.

Referenced by EnthalpyConverter(), isTemperate(), and viewConstants().


The documentation for this class was generated from the following files:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines