Loading [MathJax]/extensions/tex2jax.js
PISM, A Parallel Ice Sheet Model 2.2.2-d6b3a29ca committed by Constantine Khrulev on 2025-03-28
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
LingleClark.cc
Go to the documentation of this file.
1// Copyright (C) 2010--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 "pism/earth/LingleClark.hh"
20
21#include "pism/util/io/File.hh"
22#include "pism/util/Grid.hh"
23#include "pism/util/ConfigInterface.hh"
24#include "pism/util/error_handling.hh"
25#include "pism/util/pism_utilities.hh"
26#include "pism/util/fftw_utilities.hh"
27#include "pism/earth/LingleClarkSerial.hh"
28#include "pism/util/Context.hh"
29#include <memory>
30
31namespace pism {
32namespace bed {
33
34LingleClark::LingleClark(std::shared_ptr<const Grid> grid)
35 : BedDef(grid, "Lingle-Clark"),
36 m_total_displacement(m_grid, "bed_displacement"),
37 m_relief(m_grid, "bed_relief"),
38 m_elastic_displacement(grid, "elastic_bed_displacement") {
39
40 // A work vector. This storage is used to put thickness change on rank 0 and to get the plate
41 // displacement change back.
43 .long_name(
44 "total (viscous and elastic) displacement in the Lingle-Clark bed deformation model")
45 .units("meters");
46
48
50 .long_name("bed relief relative to the modeled bed displacement")
51 .units("meters");
52
53 bool use_elastic_model = m_config->get_flag("bed_deformation.lc.elastic_model");
54
56 .long_name(
57 "elastic part of the displacement in the Lingle-Clark bed deformation model; see :cite:`BLKfastearth`")
58 .units("meters");
60
61 const int
62 Mx = m_grid->Mx(),
63 My = m_grid->My(),
64 Z = m_config->get_number("bed_deformation.lc.grid_size_factor"),
65 Nx = Z*(Mx - 1) + 1,
66 Ny = Z*(My - 1) + 1;
67
68 const double
69 Lx = Z * (m_grid->x0() - m_grid->x(0)),
70 Ly = Z * (m_grid->y0() - m_grid->y(0));
71
72 m_extended_grid = Grid::Shallow(m_grid->ctx(), Lx, Ly, m_grid->x0(), m_grid->y0(), Nx, Ny,
74
76 std::make_shared<array::Scalar>(m_extended_grid, "viscous_bed_displacement");
77 m_viscous_displacement->metadata(0)
78 .long_name(
79 "bed displacement in the viscous half-space bed deformation model; see BuelerLingleBrown")
80 .units("meters");
81
82 // coordinate variables of the extended grid should have different names
83 m_viscous_displacement->metadata().x().set_name("x_lc");
84 m_viscous_displacement->metadata().y().set_name("y_lc");
85
86 // do not point to auxiliary coordinates "lon" and "lat".
87 m_viscous_displacement->metadata()["coordinates"] = "";
88
89 m_viscous_displacement0 = m_viscous_displacement->allocate_proc0_copy();
90
91 ParallelSection rank0(m_grid->com);
92 try {
93 if (m_grid->rank() == 0) {
94 m_serial_model.reset(new LingleClarkSerial(m_log, *m_config, use_elastic_model,
95 Mx, My,
96 m_grid->dx(), m_grid->dy(),
97 Nx, Ny));
98 }
99 } catch (...) {
100 rank0.failed();
101 }
102 rank0.check();
103}
104
106 // empty, but implemented here instead of using "= default" in the header to be able to
107 // use the forward declaration of LingleClarkSerial in LingleClark.hh
108}
109
110/*!
111 * Initialize the model by computing the viscous bed displacement using uplift and the elastic
112 * response using ice thickness.
113 *
114 * Then compute the bed relief as the difference between bed elevation and total bed displacement.
115 *
116 * This method has to initialize m_viscous_displacement, m_elastic_displacement,
117 * m_total_displacement, and m_relief.
118 */
120 const array::Scalar &bed_uplift,
121 const array::Scalar &ice_thickness,
122 const array::Scalar &sea_level_elevation) {
123
124 auto load_proc0 = m_load.allocate_proc0_copy();
125
127
128 // initialize the plate displacement
129 {
130 auto &uplift_proc0 = *m_work0;
131 bed_uplift.put_on_proc0(uplift_proc0);
132
133 m_load.set(0.0);
134 accumulate_load(bed_elevation, ice_thickness, sea_level_elevation, 1.0, m_load);
135 m_load.put_on_proc0(*load_proc0);
136
137 ParallelSection rank0(m_grid->com);
138 try {
139 if (m_grid->rank() == 0) {
140 PetscErrorCode ierr = 0;
141
142 m_serial_model->bootstrap(*load_proc0, uplift_proc0);
143
144 ierr = VecCopy(m_serial_model->total_displacement(), total_displacement);
145 PISM_CHK(ierr, "VecCopy");
146
147 ierr = VecCopy(m_serial_model->viscous_displacement(), *m_viscous_displacement0);
148 PISM_CHK(ierr, "VecCopy");
149
150 ierr = VecCopy(m_serial_model->elastic_displacement(), *m_elastic_displacement0);
151 PISM_CHK(ierr, "VecCopy");
152 }
153 } catch (...) {
154 rank0.failed();
155 }
156 rank0.check();
157 }
158
160
162
164
165 // compute bed relief
167}
168
169/*!
170 * Return the load response matrix for the elastic response.
171 *
172 * This method is used for testing only.
173 */
174std::shared_ptr<array::Scalar> LingleClark::elastic_load_response_matrix() const {
175 std::shared_ptr<array::Scalar> result(new array::Scalar(m_extended_grid, "lrm"));
176
177 int
178 Nx = m_extended_grid->Mx(),
179 Ny = m_extended_grid->My();
180
181 auto lrm0 = result->allocate_proc0_copy();
182
183 {
184 ParallelSection rank0(m_grid->com);
185 try {
186 if (m_grid->rank() == 0) {
187 std::vector<std::complex<double> > array(Nx * Ny);
188
189 m_serial_model->compute_load_response_matrix((fftw_complex*)array.data());
190
191 get_real_part((fftw_complex*)array.data(), 1.0, Nx, Ny, Nx, Ny, 0, 0, *lrm0);
192 }
193 } catch (...) {
194 rank0.failed();
195 }
196 rank0.check();
197 }
198
199 result->get_from_proc0(*lrm0);
200
201 return result;
202}
203
204/*! Initialize the Lingle-Clark bed deformation model.
205 *
206 * Inputs:
207 *
208 * - bed topography,
209 * - ice thickness,
210 * - plate displacement (either read from a file or bootstrapped using uplift) and
211 * possibly re-gridded.
212 */
213void LingleClark::init_impl(const InputOptions &opts, const array::Scalar &ice_thickness,
214 const array::Scalar &sea_level_elevation) {
215 if (opts.type == INIT_RESTART) {
216 // Set viscous displacement by reading from the input file.
217 m_viscous_displacement->read(opts.filename, opts.record);
218 // Set elastic displacement by reading from the input file.
220 } else if (opts.type == INIT_BOOTSTRAP) {
221 this->bootstrap(m_topg, m_uplift, ice_thickness, sea_level_elevation);
222 } else {
223 // do nothing
224 }
225
226 // Try re-gridding plate_displacement.
227 regrid("Lingle-Clark bed deformation model",
229 regrid("Lingle-Clark bed deformation model",
231
232 // Now that viscous displacement and elastic displacement are finally initialized,
233 // put them on rank 0 and initialize the serial model itself.
234 {
237
238 ParallelSection rank0(m_grid->com);
239 try {
240 if (m_grid->rank() == 0) { // only processor zero does the work
241 PetscErrorCode ierr = 0;
242
244
245 ierr = VecCopy(m_serial_model->total_displacement(), *m_work0);
246 PISM_CHK(ierr, "VecCopy");
247 }
248 } catch (...) {
249 rank0.failed();
250 }
251 rank0.check();
252 }
253
255
256 // compute bed relief
258}
259
260/*!
261 * Get total bed displacement on the PISM grid.
262 */
266
270
274
276 return m_relief;
277}
278
279void LingleClark::step(const array::Scalar &load_thickness,
280 double dt) {
281
282 load_thickness.put_on_proc0(*m_work0);
283
284 ParallelSection rank0(m_grid->com);
285 try {
286 if (m_grid->rank() == 0) { // only processor zero does the step
287 PetscErrorCode ierr = 0;
288
289 m_serial_model->step(dt, *m_work0);
290
291 ierr = VecCopy(m_serial_model->total_displacement(), *m_work0);
292 PISM_CHK(ierr, "VecCopy");
293
294 ierr = VecCopy(m_serial_model->viscous_displacement(), *m_viscous_displacement0);
295 PISM_CHK(ierr, "VecCopy");
296
297 ierr = VecCopy(m_serial_model->elastic_displacement(), *m_elastic_displacement0);
298 PISM_CHK(ierr, "VecCopy");
299 }
300 } catch (...) {
301 rank0.failed();
302 }
303 rank0.check();
304
306
308
310
311 // Update bed elevation using bed displacement and relief.
313}
314
315//! Update the Lingle-Clark bed deformation model.
316void LingleClark::update_impl(const array::Scalar &load, double /*t*/, double dt) {
317 step(load, dt);
318 // mark m_topg as "modified"
320}
321
328
329void LingleClark::write_model_state_impl(const File &output) const {
331
332 m_viscous_displacement->write(output);
334}
335
337 DiagnosticList result = {
338 {"viscous_bed_displacement", Diagnostic::wrap(*m_viscous_displacement)},
339 {"elastic_bed_displacement", Diagnostic::wrap(m_elastic_displacement)},
340 };
341 return combine(result, BedDef::diagnostics_impl());
342}
343
344} // end of namespace bed
345} // end of namespace pism
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
static std::shared_ptr< Grid > Shallow(std::shared_ptr< const Context > ctx, double Lx, double Ly, double x0, double y0, unsigned int Mx, unsigned int My, grid::Registration r, grid::Periodicity p)
Initialize a uniform, shallow (3 z-levels) grid with half-widths (Lx,Ly) and Mx by My nodes.
Definition Grid.cc:134
void failed()
Indicates a failure of a parallel section.
VariableMetadata & long_name(const std::string &input)
VariableMetadata & units(const std::string &input)
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 get_from_proc0(petsc::Vec &onp0)
Gets a local Array2 from processor 0.
Definition Array.cc:1058
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 inc_state_counter()
Increment the object state counter.
Definition Array.cc:150
std::shared_ptr< petsc::Vec > allocate_proc0_copy() const
Definition Array.cc:926
void put_on_proc0(petsc::Vec &onp0) const
Puts a local array::Scalar on processor 0.
Definition Array.cc:1015
SpatialVariableMetadata & metadata(unsigned int N=0)
Returns a reference to the SpatialVariableMetadata object containing metadata for the compoment N.
Definition Array.cc:476
virtual void write_model_state_impl(const File &output) const
The default (empty implementation).
Definition BedDef.cc:95
virtual DiagnosticList diagnostics_impl() const
Definition BedDef.cc:103
virtual void define_model_state_impl(const File &output) const
The default (empty implementation).
Definition BedDef.cc:79
void bootstrap(const array::Scalar &bed_elevation, const array::Scalar &bed_uplift, const array::Scalar &ice_thickness, const array::Scalar &sea_level_elevation)
Initialize using provided bed elevation and uplift.
Definition BedDef.cc:179
array::Scalar2 m_topg
current bed elevation
Definition BedDef.hh:78
array::Scalar m_uplift
bed uplift rate
Definition BedDef.hh:87
const array::Scalar & bed_elevation() const
Definition BedDef.cc:71
array::Scalar m_load
Definition BedDef.hh:83
PISM bed deformation model (base class).
Definition BedDef.hh:37
Class implementing the bed deformation model described in [BLKfastearth].
std::shared_ptr< petsc::Vec > m_elastic_displacement0
rank 0 storage for the elastic displacement
void bootstrap_impl(const array::Scalar &bed_elevation, const array::Scalar &bed_uplift, const array::Scalar &ice_thickness, const array::Scalar &sea_level_elevation)
void init_impl(const InputOptions &opts, const array::Scalar &ice_thickness, const array::Scalar &sea_level_elevation)
std::shared_ptr< petsc::Vec > m_work0
array::Scalar m_elastic_displacement
Elastic bed displacement (part of the model state)
const array::Scalar & elastic_displacement() const
std::shared_ptr< Grid > m_extended_grid
extended grid for the viscous plate displacement
virtual void write_model_state_impl(const File &output) const
The default (empty implementation).
void update_impl(const array::Scalar &load, double t, double dt)
Update the Lingle-Clark bed deformation model.
std::shared_ptr< array::Scalar > m_viscous_displacement
Viscous displacement on the extended grid (part of the model state).
array::Scalar m_relief
Bed relief relative to the bed displacement.
virtual void define_model_state_impl(const File &output) const
The default (empty implementation).
std::shared_ptr< petsc::Vec > m_viscous_displacement0
rank 0 storage using the extended grid
const array::Scalar & relief() const
void step(const array::Scalar &load_thickness, double dt)
LingleClark(std::shared_ptr< const Grid > g)
DiagnosticList diagnostics_impl() const
std::unique_ptr< LingleClarkSerial > m_serial_model
Serial viscoelastic bed deformation model.
const array::Scalar & total_displacement() const
array::Scalar m_total_displacement
Total (viscous and elastic) bed displacement.
const array::Scalar & viscous_displacement() const
std::shared_ptr< array::Scalar > elastic_load_response_matrix() const
#define PISM_CHK(errcode, name)
void accumulate_load(const array::Scalar &bed_elevation, const array::Scalar &ice_thickness, const array::Scalar &sea_level_elevation, double C, array::Scalar &result)
Definition BedDef.cc:301
@ NOT_PERIODIC
Definition Grid.hh:54
@ CELL_CORNER
Definition Grid.hh:56
@ PISM_DOUBLE
Definition IO_Flags.hh:52
@ INIT_BOOTSTRAP
Definition Component.hh:56
@ INIT_RESTART
Definition Component.hh:56
std::map< std::string, Diagnostic::Ptr > DiagnosticList
T combine(const T &a, const T &b)
void get_real_part(fftw_complex *input, double normalization, int Mx, int My, int Nx, int Ny, int i0, int j0, petsc::Vec &output)
Get the real part of input and put it in output.
InitializationType type
initialization type
Definition Component.hh:61
std::string filename
name of the input file (if applicable)
Definition Component.hh:63
unsigned int record
index of the record to re-start from
Definition Component.hh:65