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"
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"
40 std::shared_ptr<atmosphere::AtmosphereModel> 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") {
54 bool use_fausto_params =
m_config->get_flag(
"surface.pdd.fausto.enabled");
56 auto method =
m_config->get_string(
"surface.pdd.method");
58 if (method ==
"repeatable_random_process") {
60 }
else if (method ==
"random_process") {
66 if (use_fausto_params) {
71 auto sd_file =
m_config->get_string(
"surface.pdd.std_dev.file");
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");
90 .long_name(
"standard deviation of near-surface air temperature")
94 .
long_name(
"instantaneous surface mass balance (accumulation/ablation) rate")
96 .
standard_name(
"land_ice_surface_specific_mass_balance_flux");
101 .
long_name(
"snow cover depth (set to zero once a year)")
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");
132 " Computing number of positive degree-days by: %s.\n",
137 " Setting PDD parameters from [Faustoetal2009].\n");
140 " Using default PDD parameters.\n");
146 auto sd_file =
m_config->get_string(
"surface.pdd.std_dev.file");
147 if (sd_file.empty()) {
149 " Using constant standard deviation of near-surface temperature.\n");
158 " Reading standard deviation of near-surface temperature from '%s'...\n",
161 bool sd_periodic =
m_config->get_flag(
"surface.pdd.std_dev.periodic");
169 auto firn_file =
m_config->get_string(
"surface.pdd.firn_depth_file");
172 if (not firn_file.empty()) {
174 "surface.pdd.firn_depth_file is not allowed when"
175 " re-starting from a PISM output file.");
184 if (firn_file.empty()) {
193 if (firn_file.empty()) {
222 balance_year_start_day =
m_config->get_number(
"surface.mass_balance_year_start_day"),
225 balance_year_start = year_start + (balance_year_start_day - 1.0) * one_day;
227 if (balance_year_start >
time) {
228 return balance_year_start;
245 auto N =
static_cast<int>(
m_mbscheme->get_timeseries_length(dt));
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;
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");
271 if ((fausto_greve !=
nullptr) or sigmalapserate != 0.0) {
277 if (fausto_greve !=
nullptr) {
282 fausto_greve->
update_temp_mj(*surface_altitude, *latitude, *longitude);
291 const double ice_density =
m_config->get_number(
"constants.ice.density");
295 for (
auto p =
m_grid->points(); p; p.next()) {
296 const int i = p.i(), j = p.j();
301 if (mask.ice_free_ocean(i, j)) {
303 for (
int k = 0;
k < N; ++
k) {
313 for (
int k = 0;
k < N; ++
k) {
314 P[
k] = P[
k] / ice_density;
322 double tmp = (*m_air_temp_sd)(i, j);
323 for (
int k = 0;
k < N; ++
k) {
328 if (fausto_greve !=
nullptr) {
335 if (sigmalapserate != 0.0) {
336 double lat = (*latitude)(i, j);
337 for (
int k = 0;
k < N; ++
k) {
338 S[
k] += sigmalapserate * (lat - sigmabaselat);
340 (*m_air_temp_sd)(i, j) =
S[0];
345 for (
int k = 0;
k < N; ++
k) {
351 (*m_air_temp_sd)(i, j) =
S[0];
357 if (mask.ice_free_ocean(i, j)) {
358 for (
int k = 0;
k < N; ++
k) {
389 for (
int k = 0;
k < N; ++
k) {
390 if (ts[
k] >= next_snow_depth_reset) {
392 while (next_snow_depth_reset <= ts[
k]) {
431 (*m_accumulation)(i, j) = A * ice_density;
432 (*m_melt)(i, j) = M * ice_density;
433 (*m_runoff)(i, j) = R * ice_density;
440 if (mask.ice_free_ocean(i, j)) {
const units::System::Ptr m_sys
unit system used by this component
const Time & time() const
const Config::ConstPtr m_config
configuration database used by this component
const Logger::ConstPtr m_log
logger (for easy access)
const std::shared_ptr< const Grid > m_grid
grid used by this component
void regrid(const std::string &module_name, array::Array &variable, RegriddingFlag flag=NO_REGRID_WITHOUT_REGRID_VARS)
static Ptr wrap(const T &input)
High-level PISM I/O class.
array::Scalar2 ice_surface_elevation
array::CellType2 cell_type
array::Scalar2 ice_thickness
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
double calendar_year_start(double T) const
Returns the model time in seconds corresponding to the beginning of the year T falls into.
Makes sure that we call begin_access() and end_access() for all accessed array::Arrays.
void add(double alpha, const Array2D< T > &x)
void read(const std::string &filename, unsigned int time)
void define(const File &file, io::Type default_type) const
Define variables corresponding to an Array in a file opened using file.
void write(const std::string &filename) const
void set(double c)
Result: v[j] <- c for all j.
void regrid(const std::string &filename, io::Default default_value)
SpatialVariableMetadata & metadata(unsigned int N=0)
Returns a reference to the SpatialVariableMetadata object containing metadata for the compoment N.
static std::shared_ptr< Forcing > Constant(std::shared_ptr< const Grid > grid, const std::string &short_name, double value)
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
double m_next_balance_year_start
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_READONLY
open an existing file for reading only
double convert(System::Ptr system, double input, const std::string &spec1, const std::string &spec2)
Convert a quantity from unit1 to unit2.
InputOptions process_input_options(MPI_Comm com, Config::ConstPtr config)
std::map< std::string, Diagnostic::Ptr > DiagnosticList
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)