PISM, A Parallel Ice Sheet Model  stable v2.1.1 committed by Constantine Khrulev on 2024-12-04 13:36:58 -0900
TemperatureModel.cc
Go to the documentation of this file.
1 /* Copyright (C) 2016, 2017, 2018, 2019, 2020, 2022, 2023 PISM Authors
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 
20 #include "pism/energy/TemperatureModel.hh"
21 #include "pism/energy/tempSystem.hh"
22 #include "pism/energy/utilities.hh"
23 #include "pism/util/Vars.hh"
24 #include "pism/util/array/CellType.hh"
25 #include "pism/util/io/File.hh"
26 #include "pism/util/pism_utilities.hh"
27 
28 namespace pism {
29 namespace energy {
30 
32  std::shared_ptr<const Grid> grid,
33  std::shared_ptr<const stressbalance::StressBalance> stress_balance)
34  : EnergyModel(grid, stress_balance),
35  m_ice_temperature(m_grid, "temp", array::WITH_GHOSTS, m_grid->z()) {
36 
38  .long_name("ice temperature")
39  .units("K")
40  .standard_name("land_ice_temperature");
41  m_ice_temperature.metadata()["valid_min"] = {0.0};
42 }
43 
45  return m_ice_temperature;
46 }
47 
48 void TemperatureModel::restart_impl(const File &input_file, int record) {
49 
50  m_log->message(2, "* Restarting the temperature-based energy balance model from %s...\n",
51  input_file.filename().c_str());
52 
53  m_basal_melt_rate.read(input_file, record);
54 
55  const array::Scalar &ice_thickness = *m_grid->variables().get_2d_scalar("land_ice_thickness");
56 
57  if (input_file.find_variable(m_ice_temperature.metadata().get_name())) {
58  m_ice_temperature.read(input_file, record);
59  } else {
60  init_enthalpy(input_file, false, record);
62  }
63 
64  regrid("Temperature-based energy balance model", m_basal_melt_rate, REGRID_WITHOUT_REGRID_VARS);
65  regrid("Temperature-based energy balance model", m_ice_temperature, REGRID_WITHOUT_REGRID_VARS);
66 
68 }
69 
70 void TemperatureModel::bootstrap_impl(const File &input_file,
71  const array::Scalar &ice_thickness,
72  const array::Scalar &surface_temperature,
73  const array::Scalar &climatic_mass_balance,
74  const array::Scalar &basal_heat_flux) {
75 
76  m_log->message(2, "* Bootstrapping the temperature-based energy balance model from %s...\n",
77  input_file.filename().c_str());
78 
79  m_basal_melt_rate.regrid(input_file,
80  io::Default(m_config->get_number("bootstrapping.defaults.bmelt")));
81  regrid("Temperature-based energy balance model", m_basal_melt_rate, REGRID_WITHOUT_REGRID_VARS);
82 
83  int temp_revision = m_ice_temperature.state_counter();
84  regrid("Temperature-based energy balance model", m_ice_temperature, REGRID_WITHOUT_REGRID_VARS);
85 
86  if (temp_revision == m_ice_temperature.state_counter()) {
87  bootstrap_ice_temperature(ice_thickness, surface_temperature, climatic_mass_balance,
88  basal_heat_flux, m_ice_temperature);
89  }
90 
92 }
93 
95  const array::Scalar &ice_thickness,
96  const array::Scalar &surface_temperature,
97  const array::Scalar &climatic_mass_balance,
98  const array::Scalar &basal_heat_flux) {
99 
100  m_log->message(2, "* Bootstrapping the temperature-based energy balance model...\n");
101 
103  regrid("Temperature-based energy balance model", m_basal_melt_rate, REGRID_WITHOUT_REGRID_VARS);
104 
105  int temp_revision = m_ice_temperature.state_counter();
106  regrid("Temperature-based energy balance model", m_ice_temperature, REGRID_WITHOUT_REGRID_VARS);
107 
108  if (temp_revision == m_ice_temperature.state_counter()) {
109  bootstrap_ice_temperature(ice_thickness, surface_temperature, climatic_mass_balance,
110  basal_heat_flux, m_ice_temperature);
111  }
112 
114 }
115 
116 //! Takes a semi-implicit time-step for the temperature equation.
117 /*!
118 This method should be kept because it is worth having alternative physics, and
119  so that older results can be reproduced. In particular, this method is
120  documented by papers [\ref BBL,\ref BBssasliding]. The new browser page
121  \ref bombproofenth essentially documents the cold-ice-BOMBPROOF method here, but
122  the newer enthalpy-based method is slightly different and (we hope) a superior
123  implementation of the conservation of energy principle.
124 
125  The conservation of energy equation written in terms of temperature is
126  \f[ \rho c_p(T) \frac{dT}{dt} = k \frac{\partial^2 T}{\partial z^2} + \Sigma,\f]
127  where \f$T(t,x,y,z)\f$ is the temperature of the ice. This equation is the shallow approximation
128  of the full 3D conservation of energy. Note \f$dT/dt\f$ stands for the material
129  derivative, so advection is included. Here \f$\rho\f$ is the density of ice,
130  \f$c_p\f$ is its specific heat, and \f$k\f$ is its conductivity. Also \f$\Sigma\f$ is the volume
131  strain heating (with SI units of \f$J/(\text{s} \text{m}^3) = \text{W}\,\text{m}^{-3}\f$).
132 
133  We handle horizontal advection explicitly by first-order upwinding. We handle
134  vertical advection implicitly by centered differencing when possible, and retreat to
135  implicit first-order upwinding when necessary. There is a CFL condition
136  for the horizontal explicit upwinding [\ref MortonMayers]. We report
137  any CFL violations, but they are designed to not occur.
138 
139  The vertical conduction term is always handled implicitly (%i.e. by backward Euler).
140 
141  We work from the bottom of the column upward in building the system to solve
142  (in the semi-implicit time-stepping scheme). The excess energy above pressure melting
143  is converted to melt-water, and that a fraction of this melt water is transported to
144  the ice base according to the scheme in excessToFromBasalMeltLayer().
145 
146  The method uses equally-spaced calculation but the columnSystemCtx
147  methods coarse_to_fine(), fine_to_coarse() interpolate
148  back-and-forth from this equally-spaced computational grid to the
149  (usually) non-equally spaced storage grid.
150 
151  An instance of tempSystemCtx is used to solve the tridiagonal system set-up here.
152 
153  In this procedure two scalar fields are modified: basal_melt_rate and m_work.
154  But basal_melt_rate will never need to communicate ghosted values (horizontal stencil
155  neighbors). The ghosted values for m_ice_temperature are updated from the values in m_work in the
156  communication done by energy_step().
157 
158  The (older) scheme cold-ice-BOMBPROOF, implemented here, is very reliable, but there is
159  still an extreme and rare fjord situation which causes trouble. For example, it
160  occurs in one column of ice in one fjord perhaps only once
161  in a 200ka simulation of the whole sheet, in my (ELB) experience modeling the Greenland
162  ice sheet. It causes the discretized advection bulge to give temperatures below that
163  of the coldest ice anywhere, a continuum impossibility. So as a final protection
164  there is a "bulge limiter" which sets the temperature to the surface temperature of
165  the column minus the bulge maximum (15 K) if it is below that level. The number of
166  times this occurs is reported as a "BPbulge" percentage.
167  */
168 void TemperatureModel::update_impl(double t, double dt, const Inputs &inputs) {
169  // current time does not matter here
170  (void) t;
171 
172  using mask::ocean;
173 
174  Logger log(MPI_COMM_SELF, m_log->get_threshold());
175 
176  const double
177  ice_density = m_config->get_number("constants.ice.density"),
178  ice_c = m_config->get_number("constants.ice.specific_heat_capacity"),
179  L = m_config->get_number("constants.fresh_water.latent_heat_of_fusion"),
180  melting_point_temp = m_config->get_number("constants.fresh_water.melting_point_temperature"),
181  beta_CC_grad = m_config->get_number("constants.ice.beta_Clausius_Clapeyron") * ice_density * m_config->get_number("constants.standard_gravity");
182 
183  const bool allow_above_melting = m_config->get_flag("energy.allow_temperature_above_melting");
184 
185  // this is bulge limit constant in K; is max amount by which ice
186  // or bedrock can be lower than surface temperature
187  const double bulge_max = m_config->get_number("energy.enthalpy.cold_bulge_max") / ice_c;
188 
189  inputs.check();
190  const array::Array3D
191  &strain_heating3 = *inputs.volumetric_heating_rate,
192  &u3 = *inputs.u3,
193  &v3 = *inputs.v3,
194  &w3 = *inputs.w3;
195 
196  const auto &cell_type = *inputs.cell_type;
197 
198  const array::Scalar
199  &basal_frictional_heating = *inputs.basal_frictional_heating,
200  &basal_heat_flux = *inputs.basal_heat_flux,
201  &shelf_base_temp = *inputs.shelf_base_temp,
202  &ice_surface_temp = *inputs.surface_temp,
203  &till_water_thickness = *inputs.till_water_thickness;
204 
205  const array::Scalar1 &ice_thickness = *inputs.ice_thickness;
206 
207  array::AccessScope list{&ice_surface_temp, &shelf_base_temp, &ice_thickness,
208  &cell_type, &basal_heat_flux, &till_water_thickness, &basal_frictional_heating,
209  &u3, &v3, &w3, &strain_heating3, &m_basal_melt_rate, &m_ice_temperature, &m_work};
210 
211  energy::tempSystemCtx system(m_grid->z(), "temperature",
212  m_grid->dx(), m_grid->dy(), dt,
213  *m_config,
214  m_ice_temperature, u3, v3, w3, strain_heating3);
215 
216  double dz = system.dz();
217  const std::vector<double>& z_fine = system.z();
218  size_t Mz_fine = z_fine.size();
219  std::vector<double> x(Mz_fine);// space for solution of system
220  std::vector<double> Tnew(Mz_fine); // post-processed solution
221 
222  // counts unreasonably low temperature values; deprecated?
223  unsigned int maxLowTempCount = m_config->get_number("energy.max_low_temperature_count");
224  const double T_minimum = m_config->get_number("energy.minimum_allowed_temperature");
225 
226  double margin_threshold = m_config->get_number("energy.margin_ice_thickness_limit");
227 
228  ParallelSection loop(m_grid->com);
229  try {
230  for (auto p = m_grid->points(); p; p.next()) {
231  const int i = p.i(), j = p.j();
232 
233  MaskValue mask = static_cast<MaskValue>(cell_type.as_int(i,j));
234 
235  const double H = ice_thickness(i, j);
236  const double T_surface = ice_surface_temp(i, j);
237 
238  system.initThisColumn(i, j,
239  marginal(ice_thickness, i, j, margin_threshold),
240  mask, H);
241 
242  const int ks = system.ks();
243 
244  if (ks > 0) { // if there are enough points in ice to bother ...
245 
246  if (system.lambda() < 1.0) {
247  m_stats.reduced_accuracy_counter += 1; // count columns with lambda < 1
248  }
249 
250  // set boundary values for tridiagonal system
251  system.setSurfaceBoundaryValuesThisColumn(T_surface);
252  system.setBasalBoundaryValuesThisColumn(basal_heat_flux(i,j),
253  shelf_base_temp(i,j),
254  basal_frictional_heating(i,j));
255 
256  // solve the system for this column; melting not addressed yet
257  system.solveThisColumn(x);
258  } // end of "if there are enough points in ice to bother ..."
259 
260  // prepare for melting/refreezing
261  double bwatnew = till_water_thickness(i,j);
262 
263  // insert solution for generic ice segments
264  for (int k=1; k <= ks; k++) {
265  if (allow_above_melting) { // in the ice
266  Tnew[k] = x[k];
267  } else {
268  const double
269  Tpmp = melting_point_temp - beta_CC_grad * (H - z_fine[k]); // FIXME issue #15
270  if (x[k] > Tpmp) {
271  Tnew[k] = Tpmp;
272  double Texcess = x[k] - Tpmp; // always positive
273  column_drainage(ice_density, ice_c, L, z_fine[k], dz, &Texcess, &bwatnew);
274  // Texcess will always come back zero here; ignore it
275  } else {
276  Tnew[k] = x[k];
277  }
278  }
279  if (Tnew[k] < T_minimum) {
280  log.message(1,
281  " [[too low (<200) ice segment temp T = %f at %d, %d, %d;"
282  " proc %d; mask=%d; w=%f m year-1]]\n",
283  Tnew[k], i, j, k, m_grid->rank(), mask,
284  units::convert(m_sys, system.w(k), "m second-1", "m year-1"));
285 
287  }
288  if (Tnew[k] < T_surface - bulge_max) {
289  Tnew[k] = T_surface - bulge_max;
290  m_stats.bulge_counter += 1;
291  }
292  }
293 
294  // insert solution for ice base segment
295  if (ks > 0) {
296  if (allow_above_melting == true) { // ice/rock interface
297  Tnew[0] = x[0];
298  } else { // compute diff between x[k0] and Tpmp; melt or refreeze as appropriate
299  const double Tpmp = melting_point_temp - beta_CC_grad * H; // FIXME issue #15
300  double Texcess = x[0] - Tpmp; // positive or negative
301  if (ocean(mask)) {
302  // when floating, only half a segment has had its temperature raised
303  // above Tpmp
304  column_drainage(ice_density, ice_c, L, 0.0, dz/2.0, &Texcess, &bwatnew);
305  } else {
306  column_drainage(ice_density, ice_c, L, 0.0, dz, &Texcess, &bwatnew);
307  }
308  Tnew[0] = Tpmp + Texcess;
309  if (Tnew[0] > (Tpmp + 0.00001)) {
310  throw RuntimeError(PISM_ERROR_LOCATION, "updated temperature came out above Tpmp");
311  }
312  }
313  if (Tnew[0] < T_minimum) {
314  log.message(1,
315  " [[too low (<200) ice/bedrock segment temp T = %f at %d,%d;"
316  " proc %d; mask=%d; w=%f]]\n",
317  Tnew[0],i,j,m_grid->rank(), mask,
318  units::convert(m_sys, system.w(0), "m second-1", "m year-1"));
319 
321  }
322  if (Tnew[0] < T_surface - bulge_max) {
323  Tnew[0] = T_surface - bulge_max;
324  m_stats.bulge_counter += 1;
325  }
326  }
327 
328  // set to air temp above ice
329  for (unsigned int k = ks; k < Mz_fine; k++) {
330  Tnew[k] = T_surface;
331  }
332 
333  // transfer column into m_work; communication later
334  system.fine_to_coarse(Tnew, i, j, m_work);
335 
336  // basal_melt_rate(i,j) is rate of mass loss at bottom of ice
337  if (ocean(mask)) {
338  m_basal_melt_rate(i,j) = 0.0;
339  } else {
340  // basalMeltRate is rate of change of bwat; can be negative
341  // (subglacial water freezes-on); note this rate is calculated
342  // *before* limiting or other nontrivial modelling of bwat,
343  // which is Hydrology's job
344  m_basal_melt_rate(i,j) = (bwatnew - till_water_thickness(i,j)) / dt;
345  } // end of the grounded case
346  }
347  } catch (...) {
348  loop.failed();
349  }
350  loop.check();
351 
353  if (m_stats.low_temperature_counter > maxLowTempCount) {
354  throw RuntimeError::formatted(PISM_ERROR_LOCATION, "too many low temps: %d",
356  }
357 
358  // copy to m_ice_temperature, updating ghosts
360 
361  // Set ice enthalpy in place. EnergyModel::update will scatter ghosts
362  compute_enthalpy_cold(m_work, ice_thickness, m_work);
363 }
364 
365 void TemperatureModel::define_model_state_impl(const File &output) const {
368  // ice enthalpy is not a part of the model state
369 }
370 
371 void TemperatureModel::write_model_state_impl(const File &output) const {
372  m_ice_temperature.write(output);
373  m_basal_melt_rate.write(output);
374  // ice enthalpy is not a part of the model state
375 }
376 
377 //! Compute the melt water which should go to the base if \f$T\f$ is above pressure-melting.
378 void TemperatureModel::column_drainage(const double rho, const double c, const double L,
379  const double z, const double dz,
380  double *Texcess, double *bwat) const {
381 
382  const double
383  darea = m_grid->cell_area(),
384  dvol = darea * dz,
385  dE = rho * c * (*Texcess) * dvol,
386  massmelted = dE / L;
387 
388  if (*Texcess >= 0.0) {
389  // T is at or above pressure-melting temp, so temp needs to be set to
390  // pressure-melting; a fraction of excess energy is turned into melt water at base
391  // note massmelted is POSITIVE!
392  const double FRACTION_TO_BASE
393  = (z < 100.0) ? 0.2 * (100.0 - z) / 100.0 : 0.0;
394  // note: ice-equiv thickness:
395  *bwat += (FRACTION_TO_BASE * massmelted) / (rho * darea);
396  *Texcess = 0.0;
397  } else { // neither Texcess nor bwat needs to change if Texcess < 0.0
398  // Texcess negative; only refreeze (i.e. reduce bwat) if at base and bwat > 0.0
399  // note ONLY CALLED IF AT BASE! note massmelted is NEGATIVE!
400  if (z > 0.00001) {
401  throw RuntimeError(PISM_ERROR_LOCATION, "excessToBasalMeltLayer() called with z not at base and negative Texcess");
402  }
403  if (*bwat > 0.0) {
404  const double thicknessToFreezeOn = - massmelted / (rho * darea);
405  if (thicknessToFreezeOn <= *bwat) { // the water *is* available to freeze on
406  *bwat -= thicknessToFreezeOn;
407  *Texcess = 0.0;
408  } else { // only refreeze bwat thickness of water; update Texcess
409  *bwat = 0.0;
410  const double dTemp = L * (*bwat) / (c * dz);
411  *Texcess += dTemp;
412  }
413  }
414  // note: if *bwat == 0 and Texcess < 0.0 then Texcess unmolested; temp will go down
415  }
416 }
417 
418 } // end of namespace energy
419 } // end of namespace
const units::System::Ptr m_sys
unit system used by this component
Definition: Component.hh:160
const Config::ConstPtr m_config
configuration database used by this component
Definition: Component.hh:158
const Logger::ConstPtr m_log
logger (for easy access)
Definition: Component.hh:162
@ REGRID_WITHOUT_REGRID_VARS
Definition: Component.hh:151
const std::shared_ptr< const Grid > m_grid
grid used by this component
Definition: Component.hh:156
void regrid(const std::string &module_name, array::Array &variable, RegriddingFlag flag=NO_REGRID_WITHOUT_REGRID_VARS)
Definition: Component.cc:159
VariableLookupData find_variable(const std::string &short_name, const std::string &std_name) const
Find a variable using its standard name and/or short name.
Definition: File.cc:361
std::string filename() const
Definition: File.cc:307
High-level PISM I/O class.
Definition: File.hh:56
void message(int threshold, const char format[],...) const __attribute__((format(printf
Print a message to the log.
Definition: Logger.cc:49
A basic logging class.
Definition: Logger.hh:40
void failed()
Indicates a failure of a parallel section.
static RuntimeError formatted(const ErrorLocation &location, const char format[],...) __attribute__((format(printf
build a RuntimeError with a formatted message
VariableMetadata & standard_name(const std::string &input)
VariableMetadata & long_name(const std::string &input)
VariableMetadata & units(const std::string &input)
std::string get_name() const
Makes sure that we call begin_access() and end_access() for all accessed array::Arrays.
Definition: Array.hh:65
void copy_from(const Array2D< T > &source)
Definition: Array2D.hh:73
void copy_from(const Array3D &input)
Definition: Array3D.cc:208
A virtual class collecting methods common to ice and bedrock 3D fields.
Definition: Array3D.hh:33
void read(const std::string &filename, unsigned int time)
Definition: Array.cc:809
void define(const File &file, io::Type default_type) const
Define variables corresponding to an Array in a file opened using file.
Definition: Array.cc:540
void write(const std::string &filename) const
Definition: Array.cc:800
int state_counter() const
Get the object state counter.
Definition: Array.cc:128
void regrid(const std::string &filename, io::Default default_value)
Definition: Array.cc:814
SpatialVariableMetadata & metadata(unsigned int N=0)
Returns a reference to the SpatialVariableMetadata object containing metadata for the compoment N.
Definition: Array.cc:553
const std::vector< double > & z() const
unsigned int ks() const
void fine_to_coarse(const std::vector< double > &fine, int i, int j, array::Array3D &coarse) const
unsigned int reduced_accuracy_counter
Definition: EnergyModel.hh:70
unsigned int low_temperature_counter
Definition: EnergyModel.hh:71
void init_enthalpy(const File &input_file, bool regrid, int record)
Initialize enthalpy by reading it from a file, or by reading temperature and liquid water fraction,...
Definition: EnergyModel.cc:135
EnergyModelStats m_stats
Definition: EnergyModel.hh:145
const array::Scalar & basal_melt_rate() const
Basal melt rate in grounded areas. (It is set to zero elsewhere.)
Definition: EnergyModel.cc:307
array::Array3D m_ice_enthalpy
Definition: EnergyModel.hh:141
array::Scalar m_basal_melt_rate
Definition: EnergyModel.hh:143
const array::Scalar1 * ice_thickness
Definition: EnergyModel.hh:46
const array::Array3D * w3
Definition: EnergyModel.hh:55
void check() const
Definition: EnergyModel.cc:58
const array::Array3D * v3
Definition: EnergyModel.hh:54
const array::Scalar * shelf_base_temp
Definition: EnergyModel.hh:48
const array::Scalar * basal_heat_flux
Definition: EnergyModel.hh:45
const array::Scalar * basal_frictional_heating
Definition: EnergyModel.hh:44
const array::Scalar * till_water_thickness
Definition: EnergyModel.hh:50
const array::Array3D * u3
Definition: EnergyModel.hh:53
const array::Scalar * surface_temp
Definition: EnergyModel.hh:49
const array::CellType * cell_type
Definition: EnergyModel.hh:43
const array::Array3D * volumetric_heating_rate
Definition: EnergyModel.hh:52
void write_model_state_impl(const File &output) const
The default (empty implementation).
const array::Array3D & temperature() const
void initialize_impl(const array::Scalar &basal_melt_rate, const array::Scalar &ice_thickness, const array::Scalar &surface_temperature, const array::Scalar &climatic_mass_balance, const array::Scalar &basal_heat_flux)
void define_model_state_impl(const File &output) const
The default (empty implementation).
void restart_impl(const File &input_file, int record)
virtual void update_impl(double t, double dt, const Inputs &inputs)=0
TemperatureModel(std::shared_ptr< const Grid > grid, std::shared_ptr< const stressbalance::StressBalance > stress_balance)
void column_drainage(const double rho, const double c, const double L, const double z, const double dz, double *Texcess, double *bwat) const
void bootstrap_impl(const File &input_file, const array::Scalar &ice_thickness, const array::Scalar &surface_temperature, const array::Scalar &climatic_mass_balance, const array::Scalar &basal_heat_flux)
void initThisColumn(int i, int j, bool is_marginal, MaskValue new_mask, double ice_thickness)
Definition: tempSystem.cc:68
void setSurfaceBoundaryValuesThisColumn(double my_Ts)
Definition: tempSystem.cc:94
void solveThisColumn(std::vector< double > &x)
Definition: tempSystem.cc:125
void setBasalBoundaryValuesThisColumn(double my_G0, double my_Tshelfbase, double my_Rb)
Definition: tempSystem.cc:103
Tridiagonal linear system for vertical column of temperature-based conservation of energy problem.
Definition: tempSystem.hh:49
#define PISM_ERROR_LOCATION
static const double L
Definition: exactTestL.cc:40
#define rho
Definition: exactTestM.c:35
@ WITH_GHOSTS
Definition: Array.hh:62
void bootstrap_ice_temperature(const array::Scalar &ice_thickness, const array::Scalar &ice_surface_temp, const array::Scalar &surface_mass_balance, const array::Scalar &basal_heat_flux, array::Array3D &result)
Create a temperature field within the ice from provided ice thickness, surface temperature,...
Definition: utilities.cc:327
bool marginal(const array::Scalar1 &thickness, int i, int j, double threshold)
Definition: EnergyModel.cc:90
void compute_enthalpy_cold(const array::Array3D &temperature, const array::Scalar &ice_thickness, array::Array3D &result)
Compute ice enthalpy from temperature temperature by assuming the ice has zero liquid fraction.
Definition: utilities.cc:48
void compute_temperature(const array::Array3D &enthalpy, const array::Scalar &ice_thickness, array::Array3D &result)
Definition: utilities.cc:76
@ PISM_DOUBLE
Definition: IO_Flags.hh:52
bool ocean(int M)
An ocean cell (floating ice or ice-free).
Definition: Mask.hh:40
double convert(System::Ptr system, double input, const std::string &spec1, const std::string &spec2)
Convert a quantity from unit1 to unit2.
Definition: Units.cc:70
static const double k
Definition: exactTestP.cc:42
MaskValue
Definition: Mask.hh:30
void GlobalSum(MPI_Comm comm, double *local, double *result, int count)