Loading [MathJax]/extensions/tex2jax.js
PISM, A Parallel Ice Sheet Model 2.2.1-cd005eec8 committed by Constantine Khrulev on 2025-03-07
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
TemperatureIndex.cc
Go to the documentation of this file.
1// Copyright (C) 2011--2024 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#include <algorithm> // std::min
20
21#include "pism/coupler/surface/TemperatureIndex.hh"
22#include "pism/coupler/surface/localMassBalance.hh"
23#include "pism/util/Grid.hh"
24#include "pism/util/Time.hh"
25#include "pism/coupler/AtmosphereModel.hh"
26#include "pism/util/io/File.hh"
27
28#include "pism/util/error_handling.hh"
29#include "pism/util/pism_utilities.hh"
30#include "pism/util/array/CellType.hh"
31#include "pism/geometry/Geometry.hh"
32#include "pism/util/array/Forcing.hh"
33
34namespace pism {
35namespace surface {
36
37///// PISM surface model implementing a PDD scheme.
38
39TemperatureIndex::TemperatureIndex(std::shared_ptr<const Grid> g,
40 std::shared_ptr<atmosphere::AtmosphereModel> input)
41 : SurfaceModel(g, input),
42 m_mass_flux(m_grid, "climatic_mass_balance"),
43 m_firn_depth(m_grid, "firn_depth"),
44 m_snow_depth(m_grid, "snow_depth") {
45
46 m_base_ddf.snow = m_config->get_number("surface.pdd.factor_snow");
47 m_base_ddf.ice = m_config->get_number("surface.pdd.factor_ice");
48 m_base_ddf.refreeze_fraction = m_config->get_number("surface.pdd.refreeze");
49 m_base_pddStdDev = m_config->get_number("surface.pdd.std_dev.value");
50 m_sd_use_param = m_config->get_flag("surface.pdd.std_dev.use_param");
51 m_sd_param_a = m_config->get_number("surface.pdd.std_dev.param_a");
52 m_sd_param_b = m_config->get_number("surface.pdd.std_dev.param_b");
53
54 bool use_fausto_params = m_config->get_flag("surface.pdd.fausto.enabled");
55
56 auto method = m_config->get_string("surface.pdd.method");
57
58 if (method == "repeatable_random_process") {
60 } else if (method == "random_process") {
62 } else {
64 }
65
66 if (use_fausto_params) {
68 m_base_pddStdDev = 2.53;
69 }
70
71 auto sd_file = m_config->get_string("surface.pdd.std_dev.file");
72
73 if (not sd_file.empty()) {
74 bool sd_periodic = m_config->get_flag("surface.pdd.std_dev.periodic");
75 int max_buffer_size = (int) m_config->get_number("input.forcing.buffer_size");
76
77 File file(m_grid->com, sd_file, io::PISM_NETCDF3, io::PISM_READONLY);
78 m_air_temp_sd = std::make_shared<array::Forcing>(m_grid, file,
79 "air_temp_sd", "",
80 max_buffer_size,
81 sd_periodic,
82 LINEAR);
83 m_sd_file_set = true;
84 } else {
85 m_air_temp_sd = array::Forcing::Constant(m_grid, "air_temp_sd", 0.0);
86 m_sd_file_set = false;
87 }
88
89 m_air_temp_sd->metadata(0)
90 .long_name("standard deviation of near-surface air temperature")
91 .units("kelvin");
92
94 .long_name("instantaneous surface mass balance (accumulation/ablation) rate")
95 .units("kg m^-2 s^-1")
96 .standard_name("land_ice_surface_specific_mass_balance_flux");
97
98 m_mass_flux.metadata()["comment"] = "positive values correspond to ice gain";
99
101 .long_name("snow cover depth (set to zero once a year)")
102 .units("m");
103 m_snow_depth.set(0.0);
104
106 .long_name("firn cover depth")
107 .units("m");
108 m_firn_depth.metadata()["valid_min"] = {0.0};
109 m_firn_depth.set(0.0);
110
112
116}
117
119
120 // call the default implementation (not the interface method init())
121 SurfaceModel::init_impl(geometry);
122
123 // report user's modeling choices
124 {
125 m_log->message(2,
126 "* Initializing the default temperature-index, PDD-based surface processes scheme.\n"
127 " Precipitation and 2m air temperature provided by atmosphere are inputs.\n"
128 " Surface mass balance and ice upper surface temperature are outputs.\n"
129 " See PISM User's Manual for control of degree-day factors.\n");
130
131 m_log->message(2,
132 " Computing number of positive degree-days by: %s.\n",
133 m_mbscheme->method().c_str());
134
135 if (m_faustogreve) {
136 m_log->message(2,
137 " Setting PDD parameters from [Faustoetal2009].\n");
138 } else {
139 m_log->message(2,
140 " Using default PDD parameters.\n");
141 }
142 }
143
144 // initialize the spatially-variable air temperature standard deviation
145 {
146 auto sd_file = m_config->get_string("surface.pdd.std_dev.file");
147 if (sd_file.empty()) {
148 m_log->message(2,
149 " Using constant standard deviation of near-surface temperature.\n");
150
151 auto attributes = m_air_temp_sd->metadata();
152 // replace with a constant array::Forcing:
154 // restore metadata:
155 m_air_temp_sd->metadata() = attributes;
156 } else {
157 m_log->message(2,
158 " Reading standard deviation of near-surface temperature from '%s'...\n",
159 sd_file.c_str());
160
161 bool sd_periodic = m_config->get_flag("surface.pdd.std_dev.periodic");
162 m_air_temp_sd->init(sd_file, sd_periodic);
163 }
164 }
165
166 // initializing the model state
167 auto input = process_input_options(m_grid->com, m_config);
168
169 auto firn_file = m_config->get_string("surface.pdd.firn_depth_file");
170
171 if (input.type == INIT_RESTART) {
172 if (not firn_file.empty()) {
174 "surface.pdd.firn_depth_file is not allowed when"
175 " re-starting from a PISM output file.");
176 }
177
178 m_firn_depth.read(input.filename, input.record);
179 m_snow_depth.read(input.filename, input.record);
180 } else if (input.type == INIT_BOOTSTRAP) {
181
182 m_snow_depth.regrid(input.filename, io::Default(0.0));
183
184 if (firn_file.empty()) {
185 m_firn_depth.regrid(input.filename, io::Default(0.0));
186 } else {
188 }
189 } else {
190
191 m_snow_depth.set(0.0);
192
193 if (firn_file.empty()) {
194 m_firn_depth.set(0.0);
195 } else {
197 }
198 }
199
200 {
201 regrid("PDD surface model", m_snow_depth);
202 regrid("PDD surface model", m_firn_depth);
203 }
204
205 // finish up
206 {
208
209 m_accumulation->set(0.0);
210 m_melt->set(0.0);
211 m_runoff->set(0.0);
212 }
213}
214
216 return m_atmosphere->max_timestep(my_t);
217}
218
220 // compute the time corresponding to the beginning of the next balance year
221 double
222 balance_year_start_day = m_config->get_number("surface.mass_balance_year_start_day"),
223 one_day = units::convert(m_sys, 1.0, "days", "seconds"),
224 year_start = this->time().calendar_year_start(time),
225 balance_year_start = year_start + (balance_year_start_day - 1.0) * one_day;
226
227 if (balance_year_start > time) {
228 return balance_year_start;
229 }
230 return this->time().increment_date(balance_year_start, 1);
231}
232
233void TemperatureIndex::update_impl(const Geometry &geometry, double t, double dt) {
234
235 // make a copy of the pointer to convince clang static analyzer that its value does not
236 // change during the call
237 FaustoGrevePDDObject *fausto_greve = m_faustogreve.get();
238
239 // update to ensure that temperature and precipitation time series are correct:
240 m_atmosphere->update(geometry, t, dt);
241
242 m_temperature->copy_from(m_atmosphere->air_temperature());
243
244 // set up air temperature and precipitation time series
245 auto N = static_cast<int>(m_mbscheme->get_timeseries_length(dt));
246
247 const double dtseries = dt / N;
248 std::vector<double> ts(N), T(N), S(N), P(N), PDDs(N);
249 for (int k = 0; k < N; ++k) {
250 ts[k] = t + k * dtseries;
251 }
252
253 // update standard deviation time series
254 if (m_sd_file_set) {
255 m_air_temp_sd->update(t, dt);
256 m_air_temp_sd->init_interpolation(ts);
257 }
258
259 const auto &mask = geometry.cell_type;
260 const auto &H = geometry.ice_thickness;
261
262 array::AccessScope list{&mask, &H, m_air_temp_sd.get(), &m_mass_flux,
264 m_accumulation.get(), m_melt.get(), m_runoff.get()};
265
266 const double
267 sigmalapserate = m_config->get_number("surface.pdd.std_dev.lapse_lat_rate"),
268 sigmabaselat = m_config->get_number("surface.pdd.std_dev.lapse_lat_base");
269
270 const array::Scalar *latitude = &geometry.latitude;
271 if ((fausto_greve != nullptr) or sigmalapserate != 0.0) {
272 list.add(*latitude);
273 }
274
275 if (fausto_greve != nullptr) {
276 const array::Scalar
277 *longitude = &geometry.latitude,
278 *surface_altitude = &geometry.ice_surface_elevation;
279
280 fausto_greve->update_temp_mj(*surface_altitude, *latitude, *longitude);
281 }
282
284
285 m_atmosphere->init_timeseries(ts);
286
287 m_atmosphere->begin_pointwise_access();
288
289 const double ice_density = m_config->get_number("constants.ice.density");
290
291 ParallelSection loop(m_grid->com);
292 try {
293 for (auto p = m_grid->points(); p; p.next()) {
294 const int i = p.i(), j = p.j();
295
296 // the temperature time series from the AtmosphereModel and its modifiers
297 m_atmosphere->temp_time_series(i, j, T);
298
299 if (mask.ice_free_ocean(i, j)) {
300 // ignore precipitation over ice-free ocean
301 for (int k = 0; k < N; ++k) {
302 P[k] = 0.0;
303 }
304 } else {
305 // elsewhere, get precipitation from the atmosphere model
306 m_atmosphere->precip_time_series(i, j, P);
307 }
308
309 // convert precipitation from "kg m-2 second-1" to "m second-1" (PDDMassBalance expects
310 // accumulation in m/second ice equivalent)
311 for (int k = 0; k < N; ++k) {
312 P[k] = P[k] / ice_density;
313 // kg / (m^2 * second) / (kg / m^3) = m / second
314 }
315
316 // interpolate temperature standard deviation time series
317 if (m_sd_file_set) {
318 m_air_temp_sd->interp(i, j, S);
319 } else {
320 double tmp = (*m_air_temp_sd)(i, j);
321 for (int k = 0; k < N; ++k) {
322 S[k] = tmp;
323 }
324 }
325
326 if (fausto_greve != nullptr) {
327 // we have been asked to set mass balance parameters according to
328 // formula (6) in [\ref Faustoetal2009]; they overwrite ddf set above
329 ddf = fausto_greve->degree_day_factors(i, j, (*latitude)(i, j));
330 }
331
332 // apply standard deviation lapse rate on top of prescribed values
333 if (sigmalapserate != 0.0) {
334 double lat = (*latitude)(i, j);
335 for (int k = 0; k < N; ++k) {
336 S[k] += sigmalapserate * (lat - sigmabaselat);
337 }
338 (*m_air_temp_sd)(i, j) = S[0]; // ensure correct SD reporting
339 }
340
341 // apply standard deviation param over ice if in use
342 if (m_sd_use_param and mask.icy(i, j)) {
343 for (int k = 0; k < N; ++k) {
344 S[k] = m_sd_param_a * (T[k] - 273.15) + m_sd_param_b;
345 if (S[k] < 0.0) {
346 S[k] = 0.0 ;
347 }
348 }
349 (*m_air_temp_sd)(i, j) = S[0]; // ensure correct SD reporting
350 }
351
352 // Use temperature time series, the "positive" threshhold, and
353 // the standard deviation of the daily variability to get the
354 // number of positive degree days (PDDs)
355 if (mask.ice_free_ocean(i, j)) {
356 for (int k = 0; k < N; ++k) {
357 PDDs[k] = 0.0;
358 }
359 } else {
360 m_mbscheme->get_PDDs(dtseries, S, T, // inputs
361 PDDs); // output
362 }
363
364 // Use temperature time series to remove rainfall from precipitation
365 m_mbscheme->get_snow_accumulation(T, // air temperature (input)
366 P); // precipitation rate (input-output)
367
368 // Use degree-day factors, the number of PDDs, and the snow precipitation to get surface mass
369 // balance (and diagnostics: accumulation, melt, runoff)
370 {
371 double next_snow_depth_reset = m_next_balance_year_start;
372
373 // make copies of firn and snow depth values at this point to avoid accessing 2D
374 // fields in the inner loop
375 double
376 ice = H(i, j),
377 firn = m_firn_depth(i, j),
378 snow = m_snow_depth(i, j);
379
380 // accumulation, melt, runoff over this time-step
381 double
382 A = 0.0,
383 M = 0.0,
384 R = 0.0,
385 SMB = 0.0;
386
387 for (int k = 0; k < N; ++k) {
388 if (ts[k] >= next_snow_depth_reset) {
389 snow = 0.0;
390 while (next_snow_depth_reset <= ts[k]) {
391 next_snow_depth_reset = time().increment_date(next_snow_depth_reset, 1);
392 }
393 }
394
395 const double accumulation = P[k] * dtseries;
396
398 changes = m_mbscheme->step(ddf, PDDs[k],
399 ice, firn, snow, accumulation);
400
401 // update ice thickness
402 ice += changes.smb;
403 assert(ice >= 0);
404
405 // update firn depth
406 firn += changes.firn_depth;
407 assert(firn >= 0);
408
409 // update snow depth
410 snow += changes.snow_depth;
411 assert(snow >= 0);
412
413 // update total accumulation, melt, and runoff
414 {
415 A += accumulation;
416 M += changes.melt;
417 R += changes.runoff;
418 SMB += changes.smb;
419 }
420 } // end of the time-stepping loop
421
422 // set firn and snow depths
423 m_firn_depth(i, j) = firn;
424 m_snow_depth(i, j) = snow;
425
426 // set total accumulation, melt, and runoff, and SMB at this point, converting
427 // from "meters, ice equivalent" to "kg / m^2"
428 {
429 (*m_accumulation)(i, j) = A * ice_density;
430 (*m_melt)(i, j) = M * ice_density;
431 (*m_runoff)(i, j) = R * ice_density;
432 // m_mass_flux (unlike m_accumulation, m_melt, and m_runoff), is a
433 // rate. m * (kg / m^3) / second = kg / m^2 / second
434 m_mass_flux(i, j) = SMB * ice_density / dt;
435 }
436 }
437
438 if (mask.ice_free_ocean(i, j)) {
439 m_firn_depth(i, j) = 0.0; // no firn in the ocean
440 m_snow_depth(i, j) = 0.0; // snow over the ocean does not stick
441 }
442 }
443 } catch (...) {
444 loop.failed();
445 }
446 loop.check();
447
448 m_atmosphere->end_pointwise_access();
449
451}
452
456
460
464
466 return *m_melt;
467}
468
470 return *m_runoff;
471}
472
476
480
484
490
493 m_firn_depth.write(output);
494 m_snow_depth.write(output);
495}
496
498 DiagnosticList result = {
499 {"air_temp_sd", Diagnostic::wrap(*m_air_temp_sd)},
500 {"snow_depth", Diagnostic::wrap(m_snow_depth)},
501 {"firn_depth", Diagnostic::wrap(m_firn_depth)},
502 };
503
505
506 return result;
507}
508
509} // end of namespace surface
510} // end of namespace pism
const units::System::Ptr m_sys
unit system used by this component
Definition Component.hh:160
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
void regrid(const std::string &module_name, array::Array &variable, RegriddingFlag flag=NO_REGRID_WITHOUT_REGRID_VARS)
Definition Component.cc:159
static Ptr wrap(const T &input)
High-level PISM I/O class.
Definition File.hh:55
array::Scalar2 ice_surface_elevation
Definition Geometry.hh:57
array::CellType2 cell_type
Definition Geometry.hh:55
array::Scalar2 ice_thickness
Definition Geometry.hh:51
array::Scalar latitude
Definition Geometry.hh:41
Combines the max. time step with the flag indicating if a restriction is active. Makes is possible to...
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
double increment_date(double T, double years) const
Definition Time.cc:841
double calendar_year_start(double T) const
Returns the model time in seconds corresponding to the beginning of the year T falls into.
Definition Time.cc:832
VariableMetadata & long_name(const std::string &input)
VariableMetadata & units(const std::string &input)
VariableMetadata & standard_name(const std::string &input)
Makes sure that we call begin_access() and end_access() for all accessed array::Arrays.
Definition Array.hh:64
void add(double alpha, const Array2D< T > &x)
Definition Array2D.hh:65
void read(const std::string &filename, unsigned int time)
Definition Array.cc:731
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:463
void write(const std::string &filename) const
Definition Array.cc:722
void set(double c)
Result: v[j] <- c for all j.
Definition Array.cc:629
void regrid(const std::string &filename, io::Default default_value)
Definition Array.cc:736
SpatialVariableMetadata & metadata(unsigned int N=0)
Returns a reference to the SpatialVariableMetadata object containing metadata for the compoment N.
Definition Array.cc:476
static std::shared_ptr< Forcing > Constant(std::shared_ptr< const Grid > grid, const std::string &short_name, double value)
Definition Forcing.cc:148
static Default Nil()
Definition IO_Flags.hh:93
LocalMassBalance::DegreeDayFactors degree_day_factors(int i, int j, double latitude)
void update_temp_mj(const array::Scalar &surfelev, const array::Scalar &lat, const array::Scalar &lon)
Updates mean July near-surface air temperature.
A PDD implementation which computes the local mass balance based on an expectation integral.
An alternative PDD implementation which simulates a random process to get the number of PDDs.
static std::shared_ptr< array::Scalar > allocate_runoff(std::shared_ptr< const Grid > grid)
std::shared_ptr< atmosphere::AtmosphereModel > m_atmosphere
virtual DiagnosticList diagnostics_impl() const
virtual void write_model_state_impl(const File &output) const
The default (empty implementation).
virtual void define_model_state_impl(const File &output) const
The default (empty implementation).
static std::shared_ptr< array::Scalar > allocate_temperature(std::shared_ptr< const Grid > grid)
const array::Scalar & accumulation() const
Returns accumulation.
static std::shared_ptr< array::Scalar > allocate_accumulation(std::shared_ptr< const Grid > grid)
static std::shared_ptr< array::Scalar > allocate_melt(std::shared_ptr< const Grid > grid)
virtual void init_impl(const Geometry &geometry)
The interface of PISM's surface models.
std::shared_ptr< array::Scalar > m_accumulation
total accumulation during the last time step
virtual const array::Scalar & accumulation_impl() const
virtual void init_impl(const Geometry &geometry)
virtual const array::Scalar & runoff_impl() const
virtual void update_impl(const Geometry &geometry, double t, double dt)
virtual const array::Scalar & melt_impl() const
virtual const array::Scalar & mass_flux_impl() const
const array::Scalar & snow_depth() const
array::Scalar m_snow_depth
snow depth (reset once a year)
virtual const array::Scalar & temperature_impl() const
std::shared_ptr< array::Forcing > m_air_temp_sd
standard deviation of the daily variability of the air temperature
std::unique_ptr< FaustoGrevePDDObject > m_faustogreve
if not NULL then user wanted fausto PDD stuff
const array::Scalar & air_temp_sd() const
double compute_next_balance_year_start(double time)
virtual MaxTimestep max_timestep_impl(double t) const
double m_base_pddStdDev
K; daily amount of randomness.
array::Scalar m_mass_flux
cached surface mass balance rate
virtual void define_model_state_impl(const File &output) const
The default (empty implementation).
std::shared_ptr< array::Scalar > m_runoff
total runoff during the last time step
array::Scalar m_firn_depth
firn depth
std::shared_ptr< array::Scalar > m_melt
total melt during the last time step
virtual void write_model_state_impl(const File &output) const
The default (empty implementation).
std::shared_ptr< array::Scalar > m_temperature
virtual DiagnosticList diagnostics_impl() const
TemperatureIndex(std::shared_ptr< const Grid > g, std::shared_ptr< atmosphere::AtmosphereModel > input)
std::unique_ptr< LocalMassBalance > m_mbscheme
mass balance scheme to use
const array::Scalar & firn_depth() const
LocalMassBalance::DegreeDayFactors m_base_ddf
holds degree-day factors in location-independent case
#define PISM_ERROR_LOCATION
@ PISM_NETCDF3
Definition IO_Flags.hh:57
@ PISM_READONLY
open an existing file for reading only
Definition IO_Flags.hh:68
@ PISM_DOUBLE
Definition IO_Flags.hh:52
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
InputOptions process_input_options(MPI_Comm com, Config::ConstPtr config)
Definition Component.cc:43
static const double g
Definition exactTestP.cc:36
@ INIT_BOOTSTRAP
Definition Component.hh:56
@ INIT_RESTART
Definition Component.hh:56
std::map< std::string, Diagnostic::Ptr > DiagnosticList
static const double k
Definition exactTestP.cc:42
T combine(const T &a, const T &b)
double ice
m day^-1 K^-1; ice-equivalent amount of ice melted, per PDD
double refreeze_fraction
fraction of melted snow which refreezes as ice
double snow
m day^-1 K^-1; ice-equivalent amount of snow melted, per PDD
A struct which holds degree day factors.
static double S(unsigned n)
Definition test_cube.c:58