PISM, A Parallel Ice Sheet Model  stable v2.1-1-g6902d5502 committed by Ed Bueler on 2023-12-20 08:38:27 -0800
Pico.cc
Go to the documentation of this file.
1 // Copyright (C) 2012-2019, 2021, 2022, 2023 Constantine Khrulev, Ricarda Winkelmann, Ronja Reese, Torsten
2 // Albrecht, and Matthias Mengel
3 //
4 // This file is part of PISM.
5 //
6 // PISM is free software; you can redistribute it and/or modify it under the
7 // terms of the GNU General Public License as published by the Free Software
8 // Foundation; either version 2 of the License, or (at your option) any later
9 // version.
10 //
11 // PISM is distributed in the hope that it will be useful, but WITHOUT ANY
12 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
13 // FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
14 // details.
15 //
16 // You should have received a copy of the GNU General Public License
17 // along with PISM; if not, write to the Free Software
18 // Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
19 //
20 // Please cite this model as:
21 // 1.
22 // Antarctic sub-shelf melt rates via PICO
23 // R. Reese, T. Albrecht, M. Mengel, X. Asay-Davis and R. Winkelmann
24 // The Cryosphere, 12, 1969-1985, (2018)
25 // DOI: 10.5194/tc-12-1969-2018
26 //
27 // 2.
28 // A box model of circulation and melting in ice shelf caverns
29 // D. Olbers & H. Hellmer
30 // Ocean Dynamics (2010), Volume 60, Issue 1, pp 141–153
31 // DOI: 10.1007/s10236-009-0252-z
32 
33 #include <gsl/gsl_math.h> // GSL_NAN
34 
35 #include "pism/coupler/util/options.hh"
36 #include "pism/geometry/Geometry.hh"
37 #include "pism/util/ConfigInterface.hh"
38 #include "pism/util/Grid.hh"
39 #include "pism/util/Mask.hh"
40 #include "pism/util/Time.hh"
41 
42 #include "pism/coupler/ocean/Pico.hh"
43 #include "pism/coupler/ocean/PicoGeometry.hh"
44 #include "pism/coupler/ocean/PicoPhysics.hh"
45 #include "pism/util/array/Forcing.hh"
46 
47 namespace pism {
48 namespace ocean {
49 
50 Pico::Pico(std::shared_ptr<const Grid> grid)
51  : CompleteOceanModel(grid, std::shared_ptr<OceanModel>()),
52  m_Soc(m_grid, "pico_salinity"),
53  m_Soc_box0(m_grid, "pico_salinity_box0"),
54  m_Toc(m_grid, "pico_temperature"),
55  m_Toc_box0(m_grid, "pico_temperature_box0"),
56  m_T_star(m_grid, "pico_T_star"),
57  m_overturning(m_grid, "pico_overturning"),
58  m_basal_melt_rate(m_grid, "pico_basal_melt_rate"),
59  m_geometry(grid),
60  m_n_basins(0),
61  m_n_boxes(0),
62  m_n_shelves(0) {
63 
64  ForcingOptions opt(*m_grid->ctx(), "ocean.pico");
65 
66  {
67  auto buffer_size = static_cast<int>(m_config->get_number("input.forcing.buffer_size"));
68 
70 
71  m_theta_ocean = std::make_shared<array::Forcing>(m_grid,
72  file,
73  "theta_ocean",
74  "", // no standard name
75  buffer_size,
76  opt.periodic,
77  LINEAR);
78 
79  m_salinity_ocean = std::make_shared<array::Forcing>(m_grid,
80  file,
81  "salinity_ocean",
82  "", // no standard name
83  buffer_size,
84  opt.periodic,
85  LINEAR);
86  }
87 
88  m_theta_ocean->metadata(0)
89  .long_name("potential temperature of the adjacent ocean")
90  .units("Kelvin");
91 
92  m_salinity_ocean->metadata(0)
93  .long_name("salinity of the adjacent ocean")
94  .units("g/kg");
95 
96  // computed salinity in ocean boxes
97  m_Soc.metadata(0).long_name("ocean salinity field").units("g/kg");
98  m_Soc.metadata()["_FillValue"] = { 0.0 };
99 
100  // salinity input for box 1
101  m_Soc_box0.metadata(0).long_name("ocean base salinity field").units("g/kg");
102  m_Soc_box0.metadata()["_FillValue"] = { 0.0 };
103 
104  // computed temperature in ocean boxes
105  m_Toc.metadata(0).long_name("ocean temperature field").units("K");
106  m_Toc.metadata()["_FillValue"] = { 0.0 };
107 
108  // temperature input for box 1
109  m_Toc_box0.metadata(0).long_name("ocean base temperature").units("K");
110  m_Toc_box0.metadata()["_FillValue"] = { 0.0 };
111 
112  m_T_star.metadata(0).long_name("T_star field").units("Celsius");
113  m_T_star.metadata()["_FillValue"] = { 0.0 };
114 
115  m_overturning.metadata(0).long_name("cavity overturning").units("m^3 s-1");
116  m_overturning.metadata()["_FillValue"] = { 0.0 };
117 
119  .long_name("PICO sub-shelf melt rate")
120  .units("m s-1")
121  .output_units("m year-1");
122  m_basal_melt_rate.metadata()["_FillValue"] = {0.0};
123 
124  m_shelf_base_temperature->metadata()["_FillValue"] = {0.0};
125 
126  m_n_boxes = static_cast<int>(m_config->get_number("ocean.pico.number_of_boxes"));
127 }
128 
129 void Pico::init_impl(const Geometry &geometry) {
130  (void) geometry;
131 
132  m_log->message(2, "* Initializing the Potsdam Ice-shelf Cavity mOdel for the ocean ...\n");
133 
134  ForcingOptions opt(*m_grid->ctx(), "ocean.pico");
135 
136  m_theta_ocean->init(opt.filename, opt.periodic);
137  m_salinity_ocean->init(opt.filename, opt.periodic);
138 
139  // This initializes the basin_mask
140  m_geometry.init();
141 
142  // FIXME: m_n_basins is a misnomer
143  m_n_basins = static_cast<int>(max(m_geometry.basin_mask())) + 1;
144 
145  m_log->message(4, "PICO basin min=%f, max=%f\n",
148 
149  PicoPhysics physics(*m_config);
150 
151  m_log->message(2, " -Using %d drainage basins and values: \n"
152  " gamma_T= %.2e, overturning_coeff = %.2e... \n",
153  m_n_basins - 1, physics.gamma_T(), physics.overturning_coeff());
154 
155  m_log->message(2, " -Depth of continental shelf for computation of temperature and salinity input\n"
156  " is set for whole domain to continental_shelf_depth=%.0f meter\n",
157  physics.continental_shelf_depth());
158 
159  // read time-independent data right away:
160  if (m_theta_ocean->buffer_size() == 1 and m_salinity_ocean->buffer_size() == 1) {
161  m_theta_ocean->update(time().current(), 0.0);
162  m_salinity_ocean->update(time().current(), 0.0);
163  }
164 
165  double
166  ice_density = m_config->get_number("constants.ice.density"),
167  water_density = m_config->get_number("constants.sea_water.density"),
168  g = m_config->get_number("constants.standard_gravity");
169 
170  compute_average_water_column_pressure(geometry, ice_density, water_density, g,
172 }
173 
174 void Pico::define_model_state_impl(const File &output) const {
175 
180 
182 }
183 
184 void Pico::write_model_state_impl(const File &output) const {
185 
186  m_geometry.basin_mask().write(output);
187  m_Soc_box0.write(output);
188  m_Toc_box0.write(output);
189  m_overturning.write(output);
190 
192 }
193 
194 /*!
195 * Extend basal melt rates to grounded and ocean neighbors for consitency with subgl_melt.
196 * Note that melt rates are then simply interpolated into partially floating cells, they
197 * are not included in the calculations of PICO.
198 */
199 static void extend_basal_melt_rates(const array::CellType1 &cell_type,
200  array::Scalar1 &basal_melt_rate) {
201 
202  auto grid = basal_melt_rate.grid();
203 
204  // update ghosts of the basal melt rate so that we can use basal_melt_rate.box(i,j)
205  // below
206  basal_melt_rate.update_ghosts();
207 
208  array::AccessScope list{&cell_type, &basal_melt_rate};
209 
210  for (auto p = grid->points(); p; p.next()) {
211 
212  const int i = p.i(), j = p.j();
213 
214  auto M = cell_type.box(i, j);
215 
216  bool potential_partially_filled_cell =
217  ((M.c == MASK_GROUNDED or M.c == MASK_ICE_FREE_OCEAN) and
218  (M.w == MASK_FLOATING or M.e == MASK_FLOATING or
219  M.s == MASK_FLOATING or M.n == MASK_FLOATING or
220  M.sw == MASK_FLOATING or M.nw == MASK_FLOATING or
221  M.se == MASK_FLOATING or M.ne == MASK_FLOATING));
222 
223  if (potential_partially_filled_cell) {
224  auto BMR = basal_melt_rate.box(i, j);
225 
226  int N = 0;
227  double melt_sum = 0.0;
228 
229  melt_sum += M.nw == MASK_FLOATING ? (++N, BMR.nw) : 0.0;
230  melt_sum += M.n == MASK_FLOATING ? (++N, BMR.n) : 0.0;
231  melt_sum += M.ne == MASK_FLOATING ? (++N, BMR.ne) : 0.0;
232  melt_sum += M.e == MASK_FLOATING ? (++N, BMR.e) : 0.0;
233  melt_sum += M.se == MASK_FLOATING ? (++N, BMR.se) : 0.0;
234  melt_sum += M.s == MASK_FLOATING ? (++N, BMR.s) : 0.0;
235  melt_sum += M.sw == MASK_FLOATING ? (++N, BMR.sw) : 0.0;
236  melt_sum += M.w == MASK_FLOATING ? (++N, BMR.w) : 0.0;
237 
238  if (N != 0) { // If there are floating neigbors, return average melt rates
239  basal_melt_rate(i, j) = melt_sum / N;
240  }
241  }
242  } // end of the loop over grid points
243 }
244 
245 void Pico::update_impl(const Geometry &geometry, double t, double dt) {
246 
247  m_theta_ocean->update(t, dt);
248  m_salinity_ocean->update(t, dt);
249 
250  m_theta_ocean->average(t, dt);
251  m_salinity_ocean->average(t, dt);
252 
253  // set values that will be used outside of floating ice areas
254  {
255  double T_fill_value = m_config->get_number("constants.fresh_water.melting_point_temperature");
256  double Toc_fill_value = m_Toc.metadata().get_number("_FillValue");
257  double Soc_fill_value = m_Soc.metadata().get_number("_FillValue");
258  double M_fill_value = m_basal_melt_rate.metadata().get_number("_FillValue");
259  double O_fill_value = m_overturning.metadata().get_number("_FillValue");
260 
261  m_shelf_base_temperature->set(T_fill_value);
262  m_basal_melt_rate.set(M_fill_value);
263  m_Toc.set(Toc_fill_value);
264  m_Soc.set(Soc_fill_value);
265  m_overturning.set(O_fill_value);
266  m_T_star.set(Toc_fill_value);
267  }
268 
269  PicoPhysics physics(*m_config);
270 
271  const array::Scalar &ice_thickness = geometry.ice_thickness;
272  const auto &cell_type = geometry.cell_type;
273  const array::Scalar &bed_elevation = geometry.bed_elevation;
274 
275  // Geometric part of PICO
276  m_geometry.update(bed_elevation, cell_type);
277 
278  // FIXME: m_n_shelves is not really the number of shelves.
279  m_n_shelves = static_cast<int>(max(m_geometry.ice_shelf_mask())) + 1;
280 
281  // Physical part of PICO
282  {
283 
284  // prepare ocean input temperature and salinity
285  {
286  std::vector<double> basin_temperature(m_n_basins);
287  std::vector<double> basin_salinity(m_n_basins);
288 
293  *m_theta_ocean,
294  basin_temperature,
295  basin_salinity); // per basin
296 
297  set_ocean_input_fields(physics,
298  ice_thickness,
299  cell_type,
302  basin_temperature,
303  basin_salinity,
304  m_Toc_box0,
305  m_Soc_box0); // per shelf
306  }
307 
308  // Use the Beckmann-Goosse parameterization to set reasonable values throughout the
309  // domain.
310  beckmann_goosse(physics,
311  ice_thickness, // input
312  m_geometry.ice_shelf_mask(), // input
313  cell_type, // input
314  m_Toc_box0, // input
315  m_Soc_box0, // input
318  m_Toc,
319  m_Soc);
320 
321  // In ice shelves, replace Beckmann-Goosse values using the Olbers and Hellmer model.
322  process_box1(physics,
323  ice_thickness, // input
324  m_geometry.ice_shelf_mask(), // input
325  m_geometry.box_mask(), // input
326  m_Toc_box0, // input
327  m_Soc_box0, // input
330  m_T_star,
331  m_Toc,
332  m_Soc,
333  m_overturning);
334 
335  process_other_boxes(physics,
336  ice_thickness, // input
337  m_geometry.ice_shelf_mask(), // input
338  m_geometry.box_mask(), // input
341  m_T_star,
342  m_Toc,
343  m_Soc);
344  }
345 
347 
349  m_shelf_base_mass_flux->scale(physics.ice_density());
350 
351  double
352  ice_density = m_config->get_number("constants.ice.density"),
353  water_density = m_config->get_number("constants.sea_water.density"),
354  g = m_config->get_number("constants.standard_gravity");
355 
356  compute_average_water_column_pressure(geometry, ice_density, water_density, g,
358 }
359 
360 
362  (void) t;
363 
364  return MaxTimestep("ocean pico");
365 }
366 
367 
368 //! Compute temperature and salinity input from ocean data by averaging.
369 
370 //! We average the ocean data over the continental shelf reagion for each basin.
371 //! We use dummy ocean data if no such average can be calculated.
372 
373 
375  const array::Scalar &basin_mask,
376  const array::Scalar &continental_shelf_mask,
377  const array::Scalar &salinity_ocean,
378  const array::Scalar &theta_ocean,
379  std::vector<double> &temperature,
380  std::vector<double> &salinity) const {
381  std::vector<int> count(m_n_basins, 0);
382  // additional vectors to allreduce efficiently with IntelMPI
383  std::vector<int> countr(m_n_basins, 0);
384  std::vector<double> salinityr(m_n_basins);
385  std::vector<double> temperaturer(m_n_basins);
386 
387  temperature.resize(m_n_basins);
388  salinity.resize(m_n_basins);
389  for (int basin_id = 0; basin_id < m_n_basins; basin_id++) {
390  temperature[basin_id] = 0.0;
391  salinity[basin_id] = 0.0;
392  }
393 
394  array::AccessScope list{ &theta_ocean, &salinity_ocean, &basin_mask, &continental_shelf_mask };
395 
396  // compute the sum for each basin for region that intersects with the continental shelf
397  // area and is not covered by an ice shelf. (continental shelf mask excludes ice shelf
398  // areas)
399  for (auto p = m_grid->points(); p; p.next()) {
400  const int i = p.i(), j = p.j();
401 
402  if (continental_shelf_mask.as_int(i, j) == 2) {
403  int basin_id = basin_mask.as_int(i, j);
404 
405  count[basin_id] += 1;
406  salinity[basin_id] += salinity_ocean(i, j);
407  temperature[basin_id] += theta_ocean(i, j);
408  }
409  }
410 
411  // Divide by number of grid cells if more than zero cells belong to the basin. if no
412  // ocean_contshelf_mask values intersect with the basin, count is zero. In such case,
413  // use dummy temperature and salinity. This could happen, for example, if the ice shelf
414  // front advances beyond the continental shelf break.
415  GlobalSum(m_grid->com, count.data(), countr.data(), m_n_basins);
416  GlobalSum(m_grid->com, salinity.data(), salinityr.data(), m_n_basins);
417  GlobalSum(m_grid->com, temperature.data(), temperaturer.data(), m_n_basins);
418 
419  // copy values
420  count = countr;
421  salinity = salinityr;
422  temperature = temperaturer;
423 
424  // "dummy" basin
425  {
426  temperature[0] = physics.T_dummy();
427  salinity[0] = physics.S_dummy();
428  }
429 
430  for (int basin_id = 1; basin_id < m_n_basins; basin_id++) {
431 
432  if (count[basin_id] != 0) {
433  salinity[basin_id] /= count[basin_id];
434  temperature[basin_id] /= count[basin_id];
435 
436  m_log->message(5, " %d: temp =%.3f, salinity=%.3f\n", basin_id, temperature[basin_id], salinity[basin_id]);
437  } else {
438  m_log->message(
439  2,
440  "PICO WARNING: basin %d contains no cells with ocean data on continental shelf\n"
441  " (no values with ocean_contshelf_mask=2).\n"
442  " Using default temperature (%.3f K) and salinity (%.3f g/kg)\n"
443  " since mean salinity and temperature cannot be computed.\n"
444  " This may bias the basal melt rate estimate.\n"
445  " Please check your input data.\n",
446  basin_id, physics.T_dummy(), physics.S_dummy());
447 
448  temperature[basin_id] = physics.T_dummy();
449  salinity[basin_id] = physics.S_dummy();
450  }
451  }
452 }
453 
454 //! Set ocean ocean input from box 0 as boundary condition for box 1.
455 
456 //! Set ocean temperature and salinity (Toc_box0, Soc_box0)
457 //! from box 0 (in front of the ice shelf) as inputs for
458 //! box 1, which is the ocean box adjacent to the grounding line.
459 //!
460 //! We enforce that Toc_box0 is always at least the local pressure melting point.
462  const array::Scalar &ice_thickness,
463  const array::CellType1 &mask,
464  const array::Scalar &basin_mask,
465  const array::Scalar &shelf_mask,
466  const std::vector<double> &basin_temperature,
467  const std::vector<double> &basin_salinity,
468  array::Scalar &Toc_box0,
469  array::Scalar &Soc_box0) const {
470 
471  array::AccessScope list{ &ice_thickness, &basin_mask, &Soc_box0, &Toc_box0, &mask, &shelf_mask };
472 
473  std::vector<int> n_shelf_cells_per_basin(m_n_shelves * m_n_basins, 0);
474  std::vector<int> n_shelf_cells(m_n_shelves, 0);
475  std::vector<int> cfs_in_basins_per_shelf(m_n_shelves * m_n_basins, 0);
476  // additional vectors to allreduce efficiently with IntelMPI
477  std::vector<int> n_shelf_cells_per_basinr(m_n_shelves * m_n_basins, 0);
478  std::vector<int> n_shelf_cellsr(m_n_shelves, 0);
479  std::vector<int> cfs_in_basins_per_shelfr(m_n_shelves * m_n_basins, 0);
480 
481  // 1) count the number of cells in each shelf
482  // 2) count the number of cells in the intersection of each shelf with all the basins
483  {
484  for (auto p = m_grid->points(); p; p.next()) {
485  const int i = p.i(), j = p.j();
486  int s = shelf_mask.as_int(i, j);
487  int b = basin_mask.as_int(i, j);
488  n_shelf_cells_per_basin[s*m_n_basins+b]++;
489  n_shelf_cells[s]++;
490 
491  // find all basins b, in which the ice shelf s has a calving front with potential ocean water intrusion
492  if (mask.as_int(i, j) == MASK_FLOATING) {
493  auto M = mask.star(i, j);
494  if (M.n == MASK_ICE_FREE_OCEAN or
495  M.e == MASK_ICE_FREE_OCEAN or
496  M.s == MASK_ICE_FREE_OCEAN or
497  M.w == MASK_ICE_FREE_OCEAN) {
498  if (cfs_in_basins_per_shelf[s * m_n_basins + b] != b) {
499  cfs_in_basins_per_shelf[s * m_n_basins + b] = b;
500  }
501  }
502  }
503  }
504 
505  GlobalSum(m_grid->com, n_shelf_cells.data(),
506  n_shelf_cellsr.data(), m_n_shelves);
507  GlobalSum(m_grid->com, n_shelf_cells_per_basin.data(),
508  n_shelf_cells_per_basinr.data(), m_n_shelves*m_n_basins);
509  GlobalSum(m_grid->com, cfs_in_basins_per_shelf.data(),
510  cfs_in_basins_per_shelfr.data(), m_n_shelves*m_n_basins);
511  // copy data
512  n_shelf_cells = n_shelf_cellsr;
513  n_shelf_cells_per_basin = n_shelf_cells_per_basinr;
514  cfs_in_basins_per_shelf = cfs_in_basins_per_shelfr;
515 
516  for (int s = 0; s < m_n_shelves; s++) {
517  for (int b = 0; b < m_n_basins; b++) {
518  int sb = s * m_n_basins + b;
519  // remove ice shelf parts from the count that do not have a calving front in that basin
520  if (n_shelf_cells_per_basin[sb] > 0 and cfs_in_basins_per_shelf[sb] == 0) {
521  n_shelf_cells[s] -= n_shelf_cells_per_basin[sb];
522  n_shelf_cells_per_basin[sb] = 0;
523  }
524  }
525  }
526  }
527 
528  // now set potential temperature and salinity box 0:
529 
530  int low_temperature_counter = 0;
531  for (auto p = m_grid->points(); p; p.next()) {
532  const int i = p.i(), j = p.j();
533 
534  // make sure all temperatures are zero at the beginning of each time step
535  Toc_box0(i, j) = 0.0; // in K
536  Soc_box0(i, j) = 0.0; // in psu
537 
538  int s = shelf_mask.as_int(i, j);
539 
540  if (mask.as_int(i, j) == MASK_FLOATING and s > 0) {
541  // note: shelf_mask = 0 in lakes
542 
543  assert(n_shelf_cells[s] > 0);
544  double N = std::max(n_shelf_cells[s], 1); // protect from division by zero
545 
546  // weighted input depending on the number of shelf cells in each basin
547  for (int b = 1; b < m_n_basins; b++) { //Note: b=0 yields nan
548  int sb = s * m_n_basins + b;
549  Toc_box0(i, j) += basin_temperature[b] * n_shelf_cells_per_basin[sb] / N;
550  Soc_box0(i, j) += basin_salinity[b] * n_shelf_cells_per_basin[sb] / N;
551  }
552 
553  double theta_pm = physics.theta_pm(Soc_box0(i, j), physics.pressure(ice_thickness(i, j)));
554 
555  // temperature input for grounding line box should not be below pressure melting point
556  if (Toc_box0(i, j) < theta_pm) {
557  const double eps = 0.001;
558  // Setting Toc_box0 a little higher than theta_pm ensures that later equations are
559  // well solvable.
560  Toc_box0(i, j) = theta_pm + eps;
561  low_temperature_counter += 1;
562  }
563  }
564  }
565 
566  low_temperature_counter = GlobalSum(m_grid->com, low_temperature_counter);
567  if (low_temperature_counter > 0) {
568  m_log->message(
569  2,
570  "PICO WARNING: temperature has been below pressure melting temperature in %d cases,\n"
571  " setting it to pressure melting temperature\n",
572  low_temperature_counter);
573  }
574 }
575 
576 /*!
577  * Use the simpler parameterization due to [@ref BeckmannGoosse2003] to set default
578  * sub-shelf temperature and melt rate values.
579  *
580  * At grid points containing floating ice not connected to the ocean, set the basal melt
581  * rate to zero and set basal temperature to the pressure melting point.
582  */
583 void Pico::beckmann_goosse(const PicoPhysics &physics,
584  const array::Scalar &ice_thickness,
585  const array::Scalar &shelf_mask,
586  const array::CellType &cell_type,
587  const array::Scalar &Toc_box0,
588  const array::Scalar &Soc_box0,
589  array::Scalar &basal_melt_rate,
590  array::Scalar &basal_temperature,
591  array::Scalar &Toc,
592  array::Scalar &Soc) {
593 
594  const double T0 = m_config->get_number("constants.fresh_water.melting_point_temperature"),
595  beta_CC = m_config->get_number("constants.ice.beta_Clausius_Clapeyron"),
596  g = m_config->get_number("constants.standard_gravity"),
597  ice_density = m_config->get_number("constants.ice.density");
598 
599  array::AccessScope list{ &ice_thickness, &cell_type, &shelf_mask, &Toc_box0, &Soc_box0,
600  &Toc, &Soc, &basal_melt_rate, &basal_temperature };
601 
602  for (auto p = m_grid->points(); p; p.next()) {
603  const int i = p.i(), j = p.j();
604 
605  if (cell_type.floating_ice(i, j)) {
606  if (shelf_mask.as_int(i, j) > 0) {
607  double pressure = physics.pressure(ice_thickness(i, j));
608 
609  basal_melt_rate(i, j) =
610  physics.melt_rate_beckmann_goosse(physics.theta_pm(Soc_box0(i, j), pressure), Toc_box0(i, j));
611  basal_temperature(i, j) = physics.T_pm(Soc_box0(i, j), pressure);
612 
613  // diagnostic outputs
614  Toc(i, j) = Toc_box0(i, j); // in Kelvin
615  Soc(i, j) = Soc_box0(i, j); // in psu
616  } else {
617  // Floating ice cells not connected to the ocean.
618  const double pressure = ice_density * g * ice_thickness(i, j); // FIXME issue #15
619 
620  basal_temperature(i, j) = T0 - beta_CC * pressure;
621  basal_melt_rate(i, j) = 0.0;
622  }
623  }
624  }
625 }
626 
627 
628 void Pico::process_box1(const PicoPhysics &physics,
629  const array::Scalar &ice_thickness,
630  const array::Scalar &shelf_mask,
631  const array::Scalar &box_mask,
632  const array::Scalar &Toc_box0,
633  const array::Scalar &Soc_box0,
634  array::Scalar &basal_melt_rate,
635  array::Scalar &basal_temperature,
636  array::Scalar &T_star,
637  array::Scalar &Toc,
638  array::Scalar &Soc,
639  array::Scalar &overturning) {
640 
641  std::vector<double> box1_area(m_n_shelves);
642 
643  compute_box_area(1, shelf_mask, box_mask, box1_area);
644 
645  array::AccessScope list{ &ice_thickness, &shelf_mask, &box_mask, &T_star, &Toc_box0, &Toc,
646  &Soc_box0, &Soc, &overturning, &basal_melt_rate, &basal_temperature };
647 
648  int n_Toc_failures = 0;
649 
650  // basal melt rate, ambient temperature and salinity and overturning calculation
651  // for each box1 grid cell.
652  for (auto p = m_grid->points(); p; p.next()) {
653  const int i = p.i(), j = p.j();
654 
655  int shelf_id = shelf_mask.as_int(i, j);
656 
657  if (box_mask.as_int(i, j) == 1 and shelf_id > 0) {
658 
659  const double pressure = physics.pressure(ice_thickness(i, j));
660 
661  T_star(i, j) = physics.T_star(Soc_box0(i, j), Toc_box0(i, j), pressure);
662 
663  auto Toc_box1 = physics.Toc_box1(box1_area[shelf_id], T_star(i, j), Soc_box0(i, j), Toc_box0(i, j));
664 
665  // This can only happen if T_star > 0.25*p_coeff, in particular T_star > 0
666  // which can only happen for values of Toc_box0 close to the local pressure melting point
667  if (Toc_box1.failed) {
668  m_log->message(5,
669  "PICO WARNING: negative square root argument at %d, %d\n"
670  " probably because of positive T_star=%f \n"
671  " Not aborting, but setting square root to 0... \n",
672  i, j, T_star(i, j));
673 
674  n_Toc_failures += 1;
675  }
676 
677  Toc(i, j) = Toc_box1.value;
678  Soc(i, j) = physics.Soc_box1(Toc_box0(i, j), Soc_box0(i, j), Toc(i, j)); // in psu
679 
680  overturning(i, j) = physics.overturning(Soc_box0(i, j), Soc(i, j), Toc_box0(i, j), Toc(i, j));
681 
682  // main outputs
683  basal_melt_rate(i, j) = physics.melt_rate(physics.theta_pm(Soc(i, j), pressure), Toc(i, j));
684  basal_temperature(i, j) = physics.T_pm(Soc(i, j), pressure);
685  }
686  }
687 
688  n_Toc_failures = GlobalSum(m_grid->com, n_Toc_failures);
689  if (n_Toc_failures > 0) {
690  m_log->message(2,
691  "PICO WARNING: square-root argument for temperature calculation\n"
692  " has been negative in %d cases.\n",
693  n_Toc_failures);
694  }
695 }
696 
698  const array::Scalar &ice_thickness,
699  const array::Scalar &shelf_mask,
700  const array::Scalar &box_mask,
701  array::Scalar &basal_melt_rate,
702  array::Scalar &basal_temperature,
703  array::Scalar &T_star,
704  array::Scalar &Toc,
705  array::Scalar &Soc) const {
706 
707  std::vector<double> overturning(m_n_shelves, 0.0);
708  std::vector<double> salinity(m_n_shelves, 0.0);
709  std::vector<double> temperature(m_n_shelves, 0.0);
710 
711  // get average overturning from box 1 that is used as input later
712  compute_box_average(1, m_overturning, shelf_mask, box_mask, overturning);
713 
714  std::vector<bool> use_beckmann_goosse(m_n_shelves);
715 
716  array::AccessScope list{ &ice_thickness, &shelf_mask, &box_mask, &T_star, &Toc,
717  &Soc, &basal_melt_rate, &basal_temperature };
718 
719  // Iterate over all boxes i for i > 1
720  for (int box = 2; box <= m_n_boxes; ++box) {
721 
722  compute_box_average(box - 1, Toc, shelf_mask, box_mask, temperature);
723  compute_box_average(box - 1, Soc, shelf_mask, box_mask, salinity);
724 
725  // find all the shelves where we should fall back to the Beckmann-Goosse
726  // parameterization
727  for (int s = 1; s < m_n_shelves; ++s) {
728  use_beckmann_goosse[s] = (salinity[s] == 0.0 or
729  temperature[s] == 0.0 or
730  overturning[s] == 0.0);
731  }
732 
733  std::vector<double> box_area;
734  compute_box_area(box, shelf_mask, box_mask, box_area);
735 
736  int n_beckmann_goosse_cells = 0;
737 
738  for (auto p = m_grid->points(); p; p.next()) {
739  const int i = p.i(), j = p.j();
740 
741  int shelf_id = shelf_mask.as_int(i, j);
742 
743  if (box_mask.as_int(i, j) == box and shelf_id > 0) {
744 
745  if (use_beckmann_goosse[shelf_id]) {
746  n_beckmann_goosse_cells += 1;
747  continue;
748  }
749 
750  // Get the input from previous box
751  const double
752  S_previous = salinity[shelf_id],
753  T_previous = temperature[shelf_id],
754  overturning_box1 = overturning[shelf_id];
755 
756  {
757  double pressure = physics.pressure(ice_thickness(i, j));
758 
759  // diagnostic outputs
760  T_star(i, j) = physics.T_star(S_previous, T_previous, pressure);
761  Toc(i, j) = physics.Toc(box_area[shelf_id], T_previous, T_star(i, j), overturning_box1, S_previous);
762  Soc(i, j) = physics.Soc(S_previous, T_previous, Toc(i, j));
763 
764  // main outputs: basal melt rate and temperature
765  basal_melt_rate(i, j) = physics.melt_rate(physics.theta_pm(Soc(i, j), pressure), Toc(i, j));
766  basal_temperature(i, j) = physics.T_pm(Soc(i, j), pressure);
767  }
768  }
769  } // loop over grid points
770 
771  n_beckmann_goosse_cells = GlobalSum(m_grid->com, n_beckmann_goosse_cells);
772  if (n_beckmann_goosse_cells > 0) {
773  m_log->message(
774  2,
775  "PICO WARNING: [box %d]: switched to the Beckmann Goosse (2003) model at %d locations\n"
776  " (no boundary data from the previous box)\n",
777  box, n_beckmann_goosse_cells);
778  }
779 
780  } // loop over boxes
781 }
782 
783 
784 // Write diagnostic variables to extra files if requested
786 
787  DiagnosticList result = {
788  { "basins", Diagnostic::wrap(m_geometry.basin_mask()) },
789  { "pico_overturning", Diagnostic::wrap(m_overturning) },
790  { "pico_salinity_box0", Diagnostic::wrap(m_Soc_box0) },
791  { "pico_temperature_box0", Diagnostic::wrap(m_Toc_box0) },
792  { "pico_box_mask", Diagnostic::wrap(m_geometry.box_mask()) },
793  { "pico_shelf_mask", Diagnostic::wrap(m_geometry.ice_shelf_mask()) },
794  { "pico_ice_rise_mask", Diagnostic::wrap(m_geometry.ice_rise_mask()) },
795  { "pico_basal_melt_rate", Diagnostic::wrap(m_basal_melt_rate) },
796  { "pico_contshelf_mask", Diagnostic::wrap(m_geometry.continental_shelf_mask()) },
797  { "pico_salinity", Diagnostic::wrap(m_Soc) },
798  { "pico_temperature", Diagnostic::wrap(m_Toc) },
799  { "pico_T_star", Diagnostic::wrap(m_T_star) },
800  { "pico_basal_temperature", Diagnostic::wrap(*m_shelf_base_temperature) },
801  };
802 
803  return combine(result, OceanModel::diagnostics_impl());
804 }
805 
806 /*!
807  * For each shelf, compute average of a given field over the box with id `box_id`.
808  *
809  * This method is used to get inputs from a previous box for the next one.
810  */
812  const array::Scalar &field,
813  const array::Scalar &shelf_mask,
814  const array::Scalar &box_mask,
815  std::vector<double> &result) const {
816 
817  array::AccessScope list{ &field, &shelf_mask, &box_mask };
818 
819  // fill results with zeros
820  result.resize(m_n_shelves);
821  for (int s = 0; s < m_n_shelves; ++s) {
822  result[s] = 0.0;
823  }
824 
825  // compute the sum of field in each shelf's box box_id
826  std::vector<int> n_cells(m_n_shelves);
827  {
828  std::vector<int> n_cells_per_box(m_n_shelves, 0);
829  for (auto p = m_grid->points(); p; p.next()) {
830  const int i = p.i(), j = p.j();
831 
832  int shelf_id = shelf_mask.as_int(i, j);
833 
834  if (box_mask.as_int(i, j) == box_id) {
835  n_cells_per_box[shelf_id] += 1;
836  result[shelf_id] += field(i, j);
837  }
838  }
839  GlobalSum(m_grid->com, n_cells_per_box.data(), n_cells.data(), m_n_shelves);
840  }
841 
842  {
843  std::vector<double> tmp(m_n_shelves);
844  GlobalSum(m_grid->com, result.data(), tmp.data(), m_n_shelves);
845  // copy data
846  result = tmp;
847  }
848 
849  for (int s = 0; s < m_n_shelves; ++s) {
850  if (n_cells[s] > 0) {
851  result[s] /= static_cast<double>(n_cells[s]);
852  }
853  }
854 }
855 
856 /*!
857  * For all shelves compute areas of boxes with id `box_id`.
858  *
859  * @param[in] box_mask box index mask
860  * @param[out] result resulting box areas.
861  *
862  * Note: shelf and box indexes start from 1.
863  */
864 void Pico::compute_box_area(int box_id,
865  const array::Scalar &shelf_mask,
866  const array::Scalar &box_mask,
867  std::vector<double> &result) const {
868  result.resize(m_n_shelves);
869  array::AccessScope list{ &shelf_mask, &box_mask };
870 
871  auto cell_area = m_grid->cell_area();
872 
873  for (auto p = m_grid->points(); p; p.next()) {
874  const int i = p.i(), j = p.j();
875 
876  int shelf_id = shelf_mask.as_int(i, j);
877 
878  if (shelf_id > 0 and box_mask.as_int(i, j) == box_id) {
879  result[shelf_id] += cell_area;
880  }
881  }
882 
883  // compute GlobalSum from index 1 to index m_n_shelves-1
884  std::vector<double> result1(m_n_shelves);
885  GlobalSum(m_grid->com, &result[1], &result1[1], m_n_shelves-1);
886  // copy data
887  for (int i = 1; i < m_n_shelves; i++) {
888  result[i] = result1[i];
889  }
890 }
891 
892 } // end of namespace ocean
893 } // 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
static Ptr wrap(const T &input)
Definition: Diagnostic.hh:160
High-level PISM I/O class.
Definition: File.hh:56
array::CellType2 cell_type
Definition: Geometry.hh:55
array::Scalar2 ice_thickness
Definition: Geometry.hh:51
array::Scalar2 bed_elevation
Definition: Geometry.hh:47
Combines the max. time step with the flag indicating if a restriction is active. Makes is possible to...
Definition: MaxTimestep.hh:31
VariableMetadata & long_name(const std::string &input)
VariableMetadata & output_units(const std::string &input)
double get_number(const std::string &name) const
Get a single-valued scalar attribute.
VariableMetadata & units(const std::string &input)
Makes sure that we call begin_access() and end_access() for all accessed array::Arrays.
Definition: Array.hh:65
stencils::Box< T > box(int i, int j) const
Definition: Array2D.hh:93
stencils::Star< T > star(int i, int j) const
Definition: Array2D.hh:79
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
std::shared_ptr< const Grid > grid() const
Definition: Array.cc:132
void set(double c)
Result: v[j] <- c for all j.
Definition: Array.cc:707
void update_ghosts()
Updates ghost points.
Definition: Array.cc:693
SpatialVariableMetadata & metadata(unsigned int N=0)
Returns a reference to the SpatialVariableMetadata object containing metadata for the compoment N.
Definition: Array.cc:553
bool floating_ice(int i, int j) const
Definition: CellType.hh:50
"Cell type" mask. Adds convenience methods to array::Scalar.
Definition: CellType.hh:30
int as_int(int i, int j) const
Definition: Scalar.hh:45
std::shared_ptr< array::Scalar > m_shelf_base_mass_flux
std::shared_ptr< array::Scalar > m_shelf_base_temperature
std::shared_ptr< array::Scalar > m_water_column_pressure
Definition: OceanModel.hh:73
virtual void define_model_state_impl(const File &output) const
The default (empty implementation).
Definition: OceanModel.cc:122
virtual DiagnosticList diagnostics_impl() const
Definition: OceanModel.cc:198
virtual void write_model_state_impl(const File &output) const
The default (empty implementation).
Definition: OceanModel.cc:129
A very rudimentary PISM ocean model.
Definition: OceanModel.hh:33
const array::Scalar & continental_shelf_mask() const
Definition: PicoGeometry.cc:69
const array::Scalar & ice_shelf_mask() const
Definition: PicoGeometry.cc:77
const array::Scalar & ice_rise_mask() const
Definition: PicoGeometry.cc:81
void update(const array::Scalar &bed_elevation, const array::CellType1 &cell_type)
const array::Scalar & basin_mask() const
Definition: PicoGeometry.cc:85
const array::Scalar & box_mask() const
Definition: PicoGeometry.cc:73
double pressure(double ice_thickness) const
Definition: PicoPhysics.cc:90
double overturning_coeff() const
Definition: PicoPhysics.cc:195
double ice_density() const
Definition: PicoPhysics.cc:207
double S_dummy() const
Definition: PicoPhysics.cc:203
TocBox1 Toc_box1(double area, double T_star, double Soc_box0, double Toc_box0) const
Definition: PicoPhysics.cc:95
double Toc(double box_area, double temperature, double T_star, double overturning, double salinity) const
Definition: PicoPhysics.cc:120
double gamma_T() const
Definition: PicoPhysics.cc:191
double continental_shelf_depth() const
Definition: PicoPhysics.cc:211
double melt_rate(double pm_point, double Toc) const
equation 8 in the PICO paper.
Definition: PicoPhysics.cc:155
double overturning(double Soc_box0, double Soc, double Toc_box0, double Toc) const
equation 3 in the PICO paper. See also equation 4.
Definition: PicoPhysics.cc:169
double theta_pm(double salinity, double pressure) const
Definition: PicoPhysics.cc:142
double Soc_box1(double Toc_box0, double Soc_box0, double Toc) const
Definition: PicoPhysics.cc:129
double T_pm(double salinity, double pressure) const
Definition: PicoPhysics.cc:149
double melt_rate_beckmann_goosse(double pot_pm_point, double Toc) const
Beckmann & Goosse meltrate.
Definition: PicoPhysics.cc:161
double T_star(double salinity, double temperature, double pressure) const
See equation A6 and lines before in PICO paper.
Definition: PicoPhysics.cc:175
double T_dummy() const
Definition: PicoPhysics.cc:199
double Soc(double salinity, double temperature, double Toc) const
Definition: PicoPhysics.cc:134
void set_ocean_input_fields(const PicoPhysics &physics, const array::Scalar &ice_thickness, const array::CellType1 &mask, const array::Scalar &basin_mask, const array::Scalar &shelf_mask, const std::vector< double > &basin_temperature, const std::vector< double > &basin_salinity, array::Scalar &Toc_box0, array::Scalar &Soc_box0) const
Set ocean ocean input from box 0 as boundary condition for box 1.
Definition: Pico.cc:461
array::Scalar m_Toc_box0
Definition: Pico.hh:55
std::shared_ptr< array::Forcing > m_salinity_ocean
Definition: Pico.hh:61
array::Scalar m_Soc
Definition: Pico.hh:54
void write_model_state_impl(const File &output) const
The default (empty implementation).
Definition: Pico.cc:184
MaxTimestep max_timestep_impl(double t) const
Definition: Pico.cc:361
void process_box1(const PicoPhysics &physics, const array::Scalar &ice_thickness, const array::Scalar &shelf_mask, const array::Scalar &box_mask, const array::Scalar &Toc_box0, const array::Scalar &Soc_box0, array::Scalar &basal_melt_rate, array::Scalar &basal_temperature, array::Scalar &T_star, array::Scalar &Toc, array::Scalar &Soc, array::Scalar &overturning)
Definition: Pico.cc:628
array::Scalar m_Toc
Definition: Pico.hh:55
void compute_box_area(int box_id, const array::Scalar &shelf_mask, const array::Scalar &box_mask, std::vector< double > &result) const
Definition: Pico.cc:864
void define_model_state_impl(const File &output) const
The default (empty implementation).
Definition: Pico.cc:174
Pico(std::shared_ptr< const Grid > g)
Definition: Pico.cc:50
array::Scalar1 m_basal_melt_rate
Definition: Pico.hh:57
void update_impl(const Geometry &geometry, double t, double dt)
Definition: Pico.cc:245
void compute_ocean_input_per_basin(const PicoPhysics &physics, const array::Scalar &basin_mask, const array::Scalar &continental_shelf_mask, const array::Scalar &salinity_ocean, const array::Scalar &theta_ocean, std::vector< double > &temperature, std::vector< double > &salinity) const
Compute temperature and salinity input from ocean data by averaging.
Definition: Pico.cc:374
array::Scalar m_Soc_box0
Definition: Pico.hh:54
void compute_box_average(int box_id, const array::Scalar &field, const array::Scalar &shelf_mask, const array::Scalar &box_mask, std::vector< double > &result) const
Definition: Pico.cc:811
std::shared_ptr< array::Forcing > m_theta_ocean
Definition: Pico.hh:61
void process_other_boxes(const PicoPhysics &physics, const array::Scalar &ice_thickness, const array::Scalar &shelf_mask, const array::Scalar &box_mask, array::Scalar &basal_melt_rate, array::Scalar &basal_temperature, array::Scalar &T_star, array::Scalar &Toc, array::Scalar &Soc) const
Definition: Pico.cc:697
PicoGeometry m_geometry
Definition: Pico.hh:59
void init_impl(const Geometry &geometry)
Definition: Pico.cc:129
array::Scalar m_T_star
Definition: Pico.hh:55
void beckmann_goosse(const PicoPhysics &physics, const array::Scalar &ice_thickness, const array::Scalar &shelf_mask, const array::CellType &cell_type, const array::Scalar &Toc_box0, const array::Scalar &Soc_box0, array::Scalar &basal_melt_rate, array::Scalar &basal_temperature, array::Scalar &Toc, array::Scalar &Soc)
Definition: Pico.cc:583
std::map< std::string, Diagnostic::Ptr > diagnostics_impl() const
Definition: Pico.cc:785
array::Scalar m_overturning
Definition: Pico.hh:56
static const double beta_CC
Definition: exactTestO.c:23
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
@ PISM_NETCDF3
Definition: IO_Flags.hh:57
@ PISM_READONLY
open an existing file for reading only
Definition: IO_Flags.hh:72
@ PISM_DOUBLE
Definition: IO_Flags.hh:52
bool ocean(int M)
An ocean cell (floating ice or ice-free).
Definition: Mask.hh:40
static void extend_basal_melt_rates(const array::CellType1 &cell_type, array::Scalar1 &basal_melt_rate)
Definition: Pico.cc:199
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 const double g
Definition: exactTestP.cc:36
std::map< std::string, Diagnostic::Ptr > DiagnosticList
Definition: Diagnostic.hh:125
@ MASK_FLOATING
Definition: Mask.hh:34
@ MASK_ICE_FREE_OCEAN
Definition: Mask.hh:35
@ MASK_GROUNDED
Definition: Mask.hh:33
T combine(const T &a, const T &b)
void GlobalSum(MPI_Comm comm, double *local, double *result, int count)
std::string filename
Definition: options.hh:33
int count
Definition: test_cube.c:16