PISM, A Parallel Ice Sheet Model  stable v2.1.1 committed by Constantine Khrulev on 2024-12-04 13:36:58 -0900
GivenTH.cc
Go to the documentation of this file.
1 // Copyright (C) 2011, 2012, 2013, 2014, 2015, 2016, 2017, 2018, 2019, 2020, 2021, 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 2 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 <gsl/gsl_poly.h>
20 #include <cassert>
21 
22 #include "pism/coupler/ocean/GivenTH.hh"
23 #include "pism/coupler/util/options.hh"
24 #include "pism/geometry/Geometry.hh"
25 #include "pism/util/ConfigInterface.hh"
26 #include "pism/util/Grid.hh"
27 #include "pism/util/Time.hh"
28 #include "pism/util/array/Forcing.hh"
29 
30 namespace pism {
31 namespace ocean {
32 
34  // coefficients of the in situ melting point temperature
35  // parameterization:
36  a[0] = -0.0575;
37  a[1] = 0.0901;
38  a[2] = -7.61e-4;
39  // coefficients of the in situ melting point potential temperature
40  // parameterization:
41  b[0] = -0.0575;
42  b[1] = 0.0921;
43  b[2] = -7.85e-4;
44 
45  // FIXME: this should not be hard-wired. Eventually we should be able
46  // to use the spatially-variable top-of-the-ice temperature.
47  shelf_top_surface_temperature = -20.0; // degrees Celsius
48 
49  gamma_T = config.get_number("ocean.th.gamma_T");
50  gamma_S = config.get_number("ocean.th.gamma_S");
51  water_latent_heat_fusion = config.get_number("constants.fresh_water.latent_heat_of_fusion");
52  sea_water_density = config.get_number("constants.sea_water.density");
53  sea_water_specific_heat_capacity = config.get_number("constants.sea_water.specific_heat_capacity");
54  ice_density = config.get_number("constants.ice.density");
55  ice_specific_heat_capacity = config.get_number("constants.ice.specific_heat_capacity");
56  ice_thermal_diffusivity = config.get_number("constants.ice.thermal_conductivity") / (ice_density * ice_specific_heat_capacity);
57  limit_salinity_range = config.get_flag("ocean.th.clip_salinity");
58 }
59 
60 GivenTH::GivenTH(std::shared_ptr<const Grid> g)
61  : CompleteOceanModel(g, std::shared_ptr<OceanModel>()) {
62 
63  ForcingOptions opt(*m_grid->ctx(), "ocean.th");
64 
65  {
66  unsigned int buffer_size = m_config->get_number("input.forcing.buffer_size");
67 
69 
70  m_theta_ocean = std::make_shared<array::Forcing>(m_grid,
71  file,
72  "theta_ocean",
73  "", // no standard name
74  buffer_size,
75  opt.periodic,
76  LINEAR);
77 
78  m_salinity_ocean = std::make_shared<array::Forcing>(m_grid,
79  file,
80  "salinity_ocean",
81  "", // no standard name
82  buffer_size,
83  opt.periodic,
84  LINEAR);
85  }
86 
87  m_theta_ocean->metadata(0)
88  .long_name("potential temperature of the adjacent ocean")
89  .units("Kelvin");
90 
91  m_salinity_ocean->metadata(0)
92  .long_name("salinity of the adjacent ocean")
93  .units("g/kg");
94 }
95 
96 void GivenTH::init_impl(const Geometry &geometry) {
97 
98  m_log->message(2,
99  "* Initializing the 3eqn melting parameterization ocean model\n"
100  " reading ocean temperature and salinity from a file...\n");
101 
102  ForcingOptions opt(*m_grid->ctx(), "ocean.th");
103 
104  // potential temperature is required
105  m_theta_ocean->init(opt.filename, opt.periodic);
106 
107  // read ocean salinity from a file if present, otherwise use a constant
108  {
110 
111  auto variable_name = m_salinity_ocean->metadata().get_name();
112 
113  if (input.find_variable(variable_name)) {
114  m_salinity_ocean->init(opt.filename, opt.periodic);
115  } else {
116  double salinity = m_config->get_number("constants.sea_water.salinity", "g / kg");
117 
118  m_salinity_ocean = array::Forcing::Constant(m_grid, variable_name, salinity);
119 
120  m_log->message(2, " Variable '%s' not found; using constant salinity: %f (g / kg).\n",
121  variable_name.c_str(), salinity);
122  }
123  }
124 
125  // read time-independent data right away:
126  if (m_theta_ocean->buffer_size() == 1 && m_salinity_ocean->buffer_size() == 1) {
127  update(geometry, time().current(), 0); // dt is irrelevant
128  }
129 
130  const double
131  ice_density = m_config->get_number("constants.ice.density"),
132  water_density = m_config->get_number("constants.sea_water.density"),
133  g = m_config->get_number("constants.standard_gravity");
134 
135  compute_average_water_column_pressure(geometry, ice_density, water_density, g,
137 }
138 
139 void GivenTH::update_impl(const Geometry &geometry, double t, double dt) {
140  m_theta_ocean->update(t, dt);
141  m_salinity_ocean->update(t, dt);
142 
143  m_theta_ocean->average(t, dt);
144  m_salinity_ocean->average(t, dt);
145 
146  Constants c(*m_config);
147 
148  const array::Scalar &ice_thickness = geometry.ice_thickness;
149 
150  array::Scalar &temperature = *m_shelf_base_temperature;
152 
153  array::AccessScope list{ &ice_thickness, m_theta_ocean.get(), m_salinity_ocean.get(),
154  &temperature, &mass_flux};
155 
156  for (auto p = m_grid->points(); p; p.next()) {
157  const int i = p.i(), j = p.j();
158 
159  double potential_temperature_celsius = (*m_theta_ocean)(i,j) - 273.15;
160 
161  double
162  shelf_base_temp_celsius = 0.0,
163  shelf_base_massflux = 0.0;
164 
166  (*m_salinity_ocean)(i,j),
167  potential_temperature_celsius,
168  ice_thickness(i,j),
169  &shelf_base_temp_celsius,
170  &shelf_base_massflux);
171 
172  // Convert from Celsius to Kelvin:
173  temperature(i,j) = shelf_base_temp_celsius + 273.15;
174  mass_flux(i,j) = shelf_base_massflux;
175  }
176 
177  // convert mass flux from [m s-1] to [kg m-2 s-1]:
178  m_shelf_base_mass_flux->scale(m_config->get_number("constants.ice.density"));
179 
180  const double
181  ice_density = m_config->get_number("constants.ice.density"),
182  water_density = m_config->get_number("constants.sea_water.density"),
183  g = m_config->get_number("constants.standard_gravity");
184 
185  compute_average_water_column_pressure(geometry, ice_density, water_density, g,
187 }
188 
190  (void) t;
191 
192  return MaxTimestep("ocean th");
193 }
194 
195 //* Evaluate the parameterization of the melting point temperature.
196 /** The value returned is in degrees Celsius.
197  */
199  double salinity, double ice_thickness) {
200  return c.a[0] * salinity + c.a[1] + c.a[2] * ice_thickness;
201 }
202 
203 /** Melt rate, obtained by solving the salt flux balance equation.
204  *
205  * @param c model constants
206  * @param sea_water_salinity sea water salinity
207  * @param basal_salinity shelf base salinity
208  *
209  * @return shelf base melt rate, in [m/s]
210  */
212  double sea_water_salinity, double basal_salinity) {
213 
214  return c.gamma_S * c.sea_water_density * (sea_water_salinity - basal_salinity) / (c.ice_density * basal_salinity);
215 }
216 
217 /** @brief Compute temperature and melt rate at the base of the shelf.
218  * Based on [@ref HellmerOlbers1989] and [@ref HollandJenkins1999].
219  *
220  * See the manual for details.
221  *
222  * @param[in] constants model constants
223  * @param[in] sea_water_salinity sea water salinity
224  * @param[in] sea_water_potential_temperature sea water potential temperature
225  * @param[in] thickness ice shelf thickness
226  * @param[out] shelf_base_temperature_out resulting basal temperature
227  * @param[out] shelf_base_melt_rate_out resulting basal melt rate
228  *
229  * @return 0 on success
230  */
231 void GivenTH::pointwise_update(const Constants &constants,
232  double sea_water_salinity,
233  double sea_water_potential_temperature,
234  double thickness,
235  double *shelf_base_temperature_out,
236  double *shelf_base_melt_rate_out) {
237 
238  assert(thickness >= 0.0);
239 
240  // This model works for sea water salinity in the range of [4, 40]
241  // psu. Ensure that input salinity is in this range.
242  const double
243  min_salinity = 4.0,
244  max_salinity = 40.0;
245 
246  if (constants.limit_salinity_range == true) {
247  if (sea_water_salinity < min_salinity) {
248  sea_water_salinity = min_salinity;
249  } else if (sea_water_salinity > max_salinity) {
250  sea_water_salinity = max_salinity;
251  }
252  }
253 
254  double basal_salinity = sea_water_salinity;
255  subshelf_salinity(constants, sea_water_salinity, sea_water_potential_temperature,
256  thickness, &basal_salinity);
257 
258  // Clip basal salinity so that we can use the freezing point
259  // temperature parameterization to recover shelf base temperature.
260  if (constants.limit_salinity_range == true) {
261  if (basal_salinity <= min_salinity) {
262  basal_salinity = min_salinity;
263  } else if (basal_salinity >= max_salinity) {
264  basal_salinity = max_salinity;
265  }
266  }
267 
268  *shelf_base_temperature_out = melting_point_temperature(constants, basal_salinity, thickness);
269 
270  *shelf_base_melt_rate_out = shelf_base_melt_rate(constants, sea_water_salinity, basal_salinity);
271 
272  // no melt if there is no ice
273  if (thickness == 0.0) {
274  *shelf_base_melt_rate_out = 0.0;
275  }
276 }
277 
278 
279 /** @brief Compute the basal salinity and make sure that it is
280  * consistent with the basal melt rate.
281  *
282  * @param[in] c constants
283  * @param[in] sea_water_salinity sea water salinity
284  * @param[in] sea_water_potential_temperature sea water potential temperature
285  * @param[in] thickness ice shelf thickness
286  * @param[out] shelf_base_salinity resulting shelf base salinity
287  *
288  * @return 0 on success
289  */
291  double sea_water_salinity,
292  double sea_water_potential_temperature,
293  double thickness,
294  double *shelf_base_salinity) {
295  double basal_salinity = sea_water_salinity;
296 
297  // first, assume that there is melt at the shelf base:
298  {
299  subshelf_salinity_melt(c, sea_water_salinity, sea_water_potential_temperature,
300  thickness, &basal_salinity);
301 
302  double basal_melt_rate = shelf_base_melt_rate(c, sea_water_salinity, basal_salinity);
303 
304  if (basal_melt_rate > 0.0) {
305  // computed basal melt rate is consistent with the assumption used
306  // to compute basal salinity
307  *shelf_base_salinity = basal_salinity;
308  return;
309  }
310  }
311 
312  // Assuming that there is melt resulted in an inconsistent
313  // (salinity, melt_rate) pair. Assume that there is freeze-on at the base.
314  {
315  subshelf_salinity_freeze_on(c, sea_water_salinity, sea_water_potential_temperature,
316  thickness, &basal_salinity);
317 
318  double basal_melt_rate = shelf_base_melt_rate(c, sea_water_salinity, basal_salinity);
319 
320  if (basal_melt_rate < 0.0) {
321  // computed basal melt rate is consistent with the assumption
322  // used to compute basal salinity
323  *shelf_base_salinity = basal_salinity;
324  return;
325  }
326  }
327 
328  // Both assumptions (above) resulted in inconsistencies. Revert to
329  // the "diffusion-only" case, which may be less accurate, but is
330  // generic and is always consistent.
331  {
332  subshelf_salinity_diffusion_only(c, sea_water_salinity, sea_water_potential_temperature,
333  thickness, &basal_salinity);
334 
335  *shelf_base_salinity = basal_salinity;
336  }
337 }
338 
339 /** Compute basal salinity in the basal melt case.
340  *
341  * We use the parameterization of the temperature gradient from [@ref
342  * Hellmeretal1998], equation 13:
343  *
344  * @f[ T_{\text{grad}} = -\Delta T\, \frac{\frac{\partial h}{\partial t}}{\kappa}, @f]
345  *
346  * where @f$ \Delta T @f$ is the difference between the ice
347  * temperature at the top of the ice column and its bottom:
348  * @f$ \Delta T = T^S - T^B. @f$ With this parameterization, we have
349  *
350  * @f[ Q_T^I = \rho_I\, c_{pI}\, {\frac{\partial h}{\partial t}}\, (T^S - T^B). @f]
351  *
352  * Then the coefficients of the quadratic equation for basal salinity
353  * (see pointwise_update()) are
354  *
355  * @f{align*}{
356  * A &= a_{0}\,\gamma_S\,c_{pI}-b_{0}\,\gamma_T\,c_{pW}\\
357  * B &= \gamma_S\,\left(L-c_{pI}\,\left(T^S+a_{0}\,S^W-a_{2}\,h-a_{1}\right)\right)+
358  * \gamma_T\,c_{pW}\,\left(\Theta^W-b_{2}\,h-b_{1}\right)\\
359  * C &= -\gamma_S\,S^W\,\left(L-c_{pI}\,\left(T^S-a_{2}\,h-a_{1}\right)\right)
360  * @f}
361  *
362  * @param[in] c physical constants, stored here to avoid looking them up in a double for loop
363  * @param[in] sea_water_salinity salinity of the ocean immediately adjacent to the shelf, [g/kg]
364  * @param[in] sea_water_potential_temperature potential temperature of the sea water, [degrees Celsius]
365  * @param[in] thickness thickness of the ice shelf, [meters]
366  * @param[out] shelf_base_salinity resulting shelf base salinity
367  *
368  * @return 0 on success
369  */
371  double sea_water_salinity,
372  double sea_water_potential_temperature,
373  double thickness,
374  double *shelf_base_salinity) {
375 
376  const double
381  S_W = sea_water_salinity,
382  Theta_W = sea_water_potential_temperature;
383 
384  // We solve a quadratic equation for Sb, the salinity at the shelf
385  // base.
386  //
387  // A*Sb^2 + B*Sb + C = 0
388  const double A = c.a[0] * c.gamma_S * c_pI - c.b[0] * c.gamma_T * c_pW;
389  const double B = (c.gamma_S * (L - c_pI * (T_S + c.a[0] * S_W - c.a[2] * thickness - c.a[1])) +
390  c.gamma_T * c_pW * (Theta_W - c.b[2] * thickness - c.b[1]));
391  const double C = -c.gamma_S * S_W * (L - c_pI * (T_S - c.a[2] * thickness - c.a[1]));
392 
393  double S1 = 0.0, S2 = 0.0;
394  const int n_roots = gsl_poly_solve_quadratic(A, B, C, &S1, &S2);
395 
396  assert(n_roots > 0);
397  assert(S2 > 0.0); // The bigger root should be positive.
398 
399  *shelf_base_salinity = S2;
400 }
401 
402 /** Compute basal salinity in the basal freeze-on case.
403  *
404  * In this case we assume that the temperature gradient at the shelf base is zero:
405  *
406  * @f[ T_{\text{grad}} = 0. @f]
407  *
408  * Please see pointwise_update() for details.
409  *
410  * In this case the coefficients of the quadratic equation for the
411  * basal salinity are:
412  *
413  * @f{align*}{
414  * A &= -b_{0}\,\gamma_T\,c_{pW} \\
415  * B &= \gamma_S\,L+\gamma_T\,c_{pW}\,\left(\Theta^W-b_{2}\,h-b_{1}\right) \\
416  * C &= -\gamma_S\,S^W\,L\\
417  * @f}
418  *
419  * @param[in] c model constants
420  * @param[in] sea_water_salinity sea water salinity
421  * @param[in] sea_water_potential_temperature sea water temperature
422  * @param[in] thickness ice shelf thickness
423  * @param[out] shelf_base_salinity resulting basal salinity
424  *
425  * @return 0 on success
426  */
428  double sea_water_salinity,
429  double sea_water_potential_temperature,
430  double thickness,
431  double *shelf_base_salinity) {
432 
433  const double
436  S_W = sea_water_salinity,
437  Theta_W = sea_water_potential_temperature,
438  h = thickness;
439 
440  // We solve a quadratic equation for Sb, the salinity at the shelf
441  // base.
442  //
443  // A*Sb^2 + B*Sb + C = 0
444  const double A = -c.b[0] * c.gamma_T * c_pW;
445  const double B = c.gamma_S * L + c.gamma_T * c_pW * (Theta_W - c.b[2] * h - c.b[1]);
446  const double C = -c.gamma_S * S_W * L;
447 
448  double S1 = 0.0, S2 = 0.0;
449  const int n_roots = gsl_poly_solve_quadratic(A, B, C, &S1, &S2);
450 
451  assert(n_roots > 0);
452  assert(S2 > 0.0); // The bigger root should be positive.
453 
454  *shelf_base_salinity = S2;
455 }
456 
457 /** @brief Compute basal salinity in the case of no basal melt and no
458  * freeze-on, with the diffusion-only temperature distribution in the
459  * ice column.
460  *
461  * In this case the temperature gradient at the base ([@ref
462  * HollandJenkins1999], equation 21) is
463  *
464  * @f[ T_{\text{grad}} = \frac{\Delta T}{h}, @f]
465  *
466  * where @f$ h @f$ is the ice shelf thickness and @f$ \Delta T = T^S -
467  * T^B @f$ is the difference between the temperature at the top and
468  * the bottom of the shelf.
469  *
470  * In this case the coefficients of the quadratic equation for the basal salinity are:
471  *
472  * @f{align*}{
473  * A &= - \frac{b_{0}\,\gamma_T\,h\,\rho_W\,c_{pW}-a_{0}\,\rho_I\,c_{pI}\,\kappa}{h\,\rho_W}\\
474  * B &= \frac{\rho_I\,c_{pI}\,\kappa\,\left(T^S-a_{2}\,h-a_{1}\right)}{h\,\rho_W}
475  +\gamma_S\,L+\gamma_T\,c_{pW}\,\left(\Theta^W-b_{2}\,h-b_{1}\right)\\
476  * C &= -\gamma_S\,S^W\,L\\
477  * @f}
478  *
479  * @param[in] c model constants
480  * @param[in] sea_water_salinity sea water salinity
481  * @param[in] sea_water_potential_temperature sea water potential temperature
482  * @param[in] thickness ice shelf thickness
483  * @param[out] shelf_base_salinity resulting basal salinity
484  *
485  * @return 0 on success
486  */
488  double sea_water_salinity,
489  double sea_water_potential_temperature,
490  double thickness,
491  double *shelf_base_salinity) {
492  const double
497  S_W = sea_water_salinity,
498  Theta_W = sea_water_potential_temperature,
499  h = thickness,
500  rho_W = c.sea_water_density,
501  rho_I = c.ice_density,
502  kappa = c.ice_thermal_diffusivity;
503 
504  // We solve a quadratic equation for Sb, the salinity at the shelf
505  // base.
506  //
507  // A*Sb^2 + B*Sb + C = 0
508  const double A = -(c.b[0] * c.gamma_T * h * rho_W * c_pW - c.a[0] * rho_I * c_pI * kappa) / (h * rho_W);
509  const double B = ((rho_I * c_pI * kappa * (T_S - c.a[2] * h - c.a[1])) / (h * rho_W) +
510  c.gamma_S * L + c.gamma_T * c_pW * (Theta_W - c.b[2] * h - c.b[1]));
511  const double C = -c.gamma_S * S_W * L;
512 
513  double S1 = 0.0, S2 = 0.0;
514  const int n_roots = gsl_poly_solve_quadratic(A, B, C, &S1, &S2);
515 
516  assert(n_roots > 0);
517  assert(S2 > 0.0); // The bigger root should be positive.
518 
519  *shelf_base_salinity = S2;
520 }
521 
522 } // end of namespace ocean
523 } // end of namespace pism
const Time & time() const
Definition: Component.cc:109
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
const std::shared_ptr< const Grid > m_grid
grid used by this component
Definition: Component.hh:156
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.
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
High-level PISM I/O class.
Definition: File.hh:56
array::Scalar2 ice_thickness
Definition: Geometry.hh:51
Combines the max. time step with the flag indicating if a restriction is active. Makes is possible to...
Definition: MaxTimestep.hh:31
Makes sure that we call begin_access() and end_access() for all accessed array::Arrays.
Definition: Array.hh:65
static std::shared_ptr< Forcing > Constant(std::shared_ptr< const Grid > grid, const std::string &short_name, double value)
Definition: Forcing.cc:147
std::shared_ptr< array::Scalar > m_shelf_base_mass_flux
std::shared_ptr< array::Scalar > m_shelf_base_temperature
Constants(const Config &config)
Definition: GivenTH.cc:33
double gamma_S
Turbulent salt transfer coefficient:
Definition: GivenTH.hh:49
double gamma_T
Turbulent heat transfer coefficient:
Definition: GivenTH.hh:47
GivenTH(std::shared_ptr< const Grid > g)
Definition: GivenTH.cc:60
void subshelf_salinity_diffusion_only(const Constants &constants, double sea_water_salinity, double sea_water_potential_temperature, double ice_thickness, double *shelf_base_salinity)
Compute basal salinity in the case of no basal melt and no freeze-on, with the diffusion-only tempera...
Definition: GivenTH.cc:487
void subshelf_salinity_melt(const Constants &constants, double sea_water_salinity, double sea_water_potential_temperature, double ice_thickness, double *shelf_base_salinity)
Definition: GivenTH.cc:370
std::shared_ptr< array::Forcing > m_salinity_ocean
Definition: GivenTH.hh:66
std::shared_ptr< array::Forcing > m_theta_ocean
Definition: GivenTH.hh:65
void subshelf_salinity_freeze_on(const Constants &constants, double sea_water_salinity, double sea_water_potential_temperature, double ice_thickness, double *shelf_base_salinity)
Definition: GivenTH.cc:427
MaxTimestep max_timestep_impl(double t) const
Definition: GivenTH.cc:189
void pointwise_update(const Constants &constants, double sea_water_salinity, double sea_water_potential_temperature, double ice_thickness, double *shelf_base_temperature_out, double *shelf_base_melt_rate_out)
Compute temperature and melt rate at the base of the shelf. Based on [HellmerOlbers1989] and [Holland...
Definition: GivenTH.cc:231
void subshelf_salinity(const Constants &constants, double sea_water_salinity, double sea_water_potential_temperature, double ice_thickness, double *shelf_base_salinity)
Compute the basal salinity and make sure that it is consistent with the basal melt rate.
Definition: GivenTH.cc:290
void update_impl(const Geometry &geometry, double t, double dt)
Definition: GivenTH.cc:139
void init_impl(const Geometry &geometry)
Definition: GivenTH.cc:96
std::shared_ptr< array::Scalar > m_water_column_pressure
Definition: OceanModel.hh:73
void update(const Geometry &geometry, double t, double dt)
Definition: OceanModel.cc:89
A very rudimentary PISM ocean model.
Definition: OceanModel.hh:33
static const double L
Definition: exactTestL.cc:40
@ PISM_GUESS
Definition: IO_Flags.hh:56
@ PISM_NETCDF3
Definition: IO_Flags.hh:57
@ PISM_READONLY
open an existing file for reading only
Definition: IO_Flags.hh:72
bool ocean(int M)
An ocean cell (floating ice or ice-free).
Definition: Mask.hh:40
static double melting_point_temperature(GivenTH::Constants c, double salinity, double ice_thickness)
Definition: GivenTH.cc:198
void compute_average_water_column_pressure(const Geometry &geometry, double ice_density, double water_density, double g, array::Scalar &result)
Definition: OceanModel.cc:217
static double shelf_base_melt_rate(GivenTH::Constants c, double sea_water_salinity, double basal_salinity)
Definition: GivenTH.cc:211
static const double g
Definition: exactTestP.cc:36
const int C[]
Definition: ssafd_code.cc:44
std::string filename
Definition: options.hh:33