|
PISM, A Parallel Ice Sheet Model stable 0.4.1779
|
Converts between specific enthalpy and temperature or liquid content. More...
#include <enthalpyConverter.hh>

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 |
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.
| EnthalpyConverter | ( | const NCConfigVariable & | config | ) |
| virtual ~EnthalpyConverter | ( | ) | [inline, virtual] |
Definition at line 58 of file enthalpyConverter.hh.
| PetscErrorCode getAbsTemp | ( | double | E, |
| double | p, | ||
| double & | T | ||
| ) | const [virtual] |
Get absolute (not pressure-adjusted) ice temperature (K) from enthalpy and pressure.
From AschwandenBuelerKhroulevBlatter,
We do not allow liquid water (i.e. water fraction
) so we return an error code of 1 if
.
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
and
are the ice mixture enthalpy and the CTS enthalpy, respectively, then
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
and
[AschwandenBuelerKhroulevBlatter]. It returns:
Certain cases are not allowed and return errors:
(error code 1)
or
(error code 2)
(error code 3)
and
(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
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
is temperate if
:
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
Only the following case returns an error:
or
(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
and
as cold ice, ignoring the water fraction (
) value.
Checks if
and returns error code 1 if so.
Computes:
but ensures
in second case. Calls getEnth() for
.
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.
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:
Definition at line 183 of file enthalpyConverter.cc.
References getAbsTemp(), getMeltingTemp(), and T_melting.
Referenced by isTemperate().
| double getPressureFromDepth | ( | double | depth | ) | const [virtual] |
Get pressure in ice from depth below surface using the hydrostatic assumption.
If
is the depth then
Frequently
is computed from the thickess minus a level in the ice, something like "H[i][j] - z[k]". The input depth to this routine is allowed to be negative, representing a position above the surface of the ice.
Definition at line 86 of file enthalpyConverter.cc.
References g, p_air, and rho_i.
Referenced by IceFlowLaw::averagedHardness_from_enth(), SIA_Sliding::basalVelocitySIA(), IceModel::check_maximum_thickness(), IceModel_tempicethk_basal::compute(), IceModel::compute_enthalpy(), IceModel::compute_enthalpy_cold(), IceModel::compute_ice_area_cold(), IceModel::compute_ice_area_temperate(), IceModel::compute_ice_volume_cold(), IceModel::compute_ice_volume_temperate(), IceModel::compute_liquid_water_fraction(), SIAFD::compute_sigma(), IceModel::energyStats(), enthalpy_from_temperature_cold(), IceModel::enthalpyAndDrainageStep(), IceModel::get_bed_top_temp(), IceModel::getEnthalpyCTSColumn(), main(), IceModel::putTempAtDepth(), IceModel::setCTSFromEnthalpy(), and SIA_Sliding::update().
| PetscErrorCode getWaterFraction | ( | double | E, |
| double | p, | ||
| double & | omega | ||
| ) | const [virtual] |
Get liquid water fraction from enthalpy and pressure.
From AschwandenBuelerKhroulevBlatter,
We do not allow liquid water (i.e. water fraction
) so we return error code 1 if
, but we still compute
.
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] |
Determines if E >= E_s(p), that is, if the ice is at the pressure-melting point.
Reimplemented in ICMEnthalpyConverter.
Definition at line 142 of file enthalpyConverter.cc.
References do_cold_ice_methods, getEnthalpyCTS(), getPATemp(), T_melting, and T_tol.
Referenced by SIA_Sliding::basalVelocitySIA(), IceModel_tempicethk_basal::compute(), IceModel::compute_ice_area_cold(), IceModel::compute_ice_area_temperate(), IceModel::compute_ice_volume_cold(), IceModel::compute_ice_volume_temperate(), IceModel::energyStats(), and IceModel::enthalpyAndDrainageStep().
| 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().
double beta [protected] |
Definition at line 81 of file enthalpyConverter.hh.
Referenced by EnthalpyConverter(), getMeltingTemp(), and viewConstants().
double c_i [protected] |
Definition at line 81 of file enthalpyConverter.hh.
Referenced by EnthalpyConverter(), ICMEnthalpyConverter::getAbsTemp(), getAbsTemp(), ICMEnthalpyConverter::getEnth(), getEnth(), getEnthalpyCTS(), ICMEnthalpyConverter::getEnthPermissive(), and viewConstants().
bool do_cold_ice_methods [protected] |
Definition at line 83 of file enthalpyConverter.hh.
Referenced by EnthalpyConverter(), ICMEnthalpyConverter::ICMEnthalpyConverter(), isTemperate(), and viewConstants().
double g [protected] |
Definition at line 81 of file enthalpyConverter.hh.
Referenced by EnthalpyConverter(), getPressureFromDepth(), and viewConstants().
double L [protected] |
Definition at line 81 of file enthalpyConverter.hh.
Referenced by EnthalpyConverter(), getEnth(), getEnthalpyInterval(), getEnthAtWaterFraction(), getWaterFraction(), and viewConstants().
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] |
Definition at line 82 of file enthalpyConverter.hh.
Referenced by EnthalpyConverter(), ICMEnthalpyConverter::getAbsTemp(), getAbsTemp(), ICMEnthalpyConverter::getEnth(), getEnth(), getEnthalpyCTS(), ICMEnthalpyConverter::getEnthPermissive(), and viewConstants().
double T_melting [protected] |
Definition at line 81 of file enthalpyConverter.hh.
Referenced by EnthalpyConverter(), ICMEnthalpyConverter::getMeltingTemp(), getMeltingTemp(), getPATemp(), isTemperate(), and viewConstants().
double T_tol [protected] |
Definition at line 81 of file enthalpyConverter.hh.
Referenced by EnthalpyConverter(), isTemperate(), and viewConstants().
1.7.3