PISM, A Parallel Ice Sheet Model  stable v2.1.1 committed by Constantine Khrulev on 2024-12-04 13:36:58 -0900
LingleClark.hh
Go to the documentation of this file.
1 /* Copyright (C) 2013, 2014, 2015, 2016, 2017, 2018, 2019, 2020, 2021, 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 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 
20 #ifndef _PBLINGLECLARK_H_
21 #define _PBLINGLECLARK_H_
22 
23 #include <memory> // std::unique_ptr
24 
25 #include "pism/earth/BedDef.hh"
26 
27 namespace pism {
28 namespace bed {
29 
30 class LingleClarkSerial;
31 
32 //! A wrapper class around LingleClarkSerial.
33 class LingleClark : public BedDef {
34 public:
35  LingleClark(std::shared_ptr<const Grid> g);
36  virtual ~LingleClark();
37 
38  const array::Scalar& total_displacement() const;
39 
40  const array::Scalar& viscous_displacement() const;
41 
42  const array::Scalar& elastic_displacement() const;
43 
44  const array::Scalar& relief() const;
45 
46  void step(const array::Scalar &ice_thickness,
47  const array::Scalar &sea_level_elevation,
48  double dt);
49 
50  std::shared_ptr<array::Scalar> elastic_load_response_matrix() const;
51 protected:
52  virtual void define_model_state_impl(const File &output) const;
53  virtual void write_model_state_impl(const File &output) const;
54 
56 
57  MaxTimestep max_timestep_impl(double t) const;
58  void init_impl(const InputOptions &opts, const array::Scalar &ice_thickness,
59  const array::Scalar &sea_level_elevation);
61  const array::Scalar &bed_uplift,
62  const array::Scalar &ice_thickness,
63  const array::Scalar &sea_level_elevation);
64  void update_impl(const array::Scalar &ice_thickness,
65  const array::Scalar &sea_level_elevation,
66  double t, double dt);
67 
68  //! Total (viscous and elastic) bed displacement.
70 
71  //! Storage on rank zero. Used to pass the load to the serial deformation model and get
72  //! bed displacement back.
73  std::shared_ptr<petsc::Vec> m_work0;
74 
75  //! Bed relief relative to the bed displacement.
77 
78  //! Ice-equivalent load thickness.
80 
81  //! Serial viscoelastic bed deformation model.
82  std::unique_ptr<LingleClarkSerial> m_serial_model;
83 
84  //! extended grid for the viscous plate displacement
85  std::shared_ptr<Grid> m_extended_grid;
86 
87  //! Viscous displacement on the extended grid (part of the model state).
88  std::shared_ptr<array::Scalar> m_viscous_displacement;
89  //! rank 0 storage using the extended grid
90  std::shared_ptr<petsc::Vec> m_viscous_displacement0;
91 
92  //! Elastic bed displacement (part of the model state)
94  //! rank 0 storage for the elastic displacement
95  std::shared_ptr<petsc::Vec> m_elastic_displacement0;
96 
97  //! time of the last bed deformation update
98  double m_t_last;
99  //! Update interval in seconds
101  //! Temporal resolution to use when checking whether it's time to update
102  double m_t_eps;
103  //! Name of the variable used to store the last update time.
104  std::string m_time_name;
105 };
106 
107 } // end of namespace bed
108 } // end of namespace pism
109 
110 #endif /* _PBLINGLECLARK_H_ */
High-level PISM I/O class.
Definition: File.hh:56
Combines the max. time step with the flag indicating if a restriction is active. Makes is possible to...
Definition: MaxTimestep.hh:31
const array::Scalar & bed_elevation() const
Definition: BedDef.cc:51
PISM bed deformation model (base class).
Definition: BedDef.hh:39
void step(const array::Scalar &ice_thickness, const array::Scalar &sea_level_elevation, double dt)
Definition: LingleClark.cc:338
std::shared_ptr< petsc::Vec > m_elastic_displacement0
rank 0 storage for the elastic displacement
Definition: LingleClark.hh:95
MaxTimestep max_timestep_impl(double t) const
Definition: LingleClark.cc:296
void bootstrap_impl(const array::Scalar &bed_elevation, const array::Scalar &bed_uplift, const array::Scalar &ice_thickness, const array::Scalar &sea_level_elevation)
Definition: LingleClark.cc:133
void init_impl(const InputOptions &opts, const array::Scalar &ice_thickness, const array::Scalar &sea_level_elevation)
Definition: LingleClark.cc:227
double m_t_eps
Temporal resolution to use when checking whether it's time to update.
Definition: LingleClark.hh:102
std::shared_ptr< petsc::Vec > m_work0
Definition: LingleClark.hh:73
array::Scalar m_elastic_displacement
Elastic bed displacement (part of the model state)
Definition: LingleClark.hh:93
const array::Scalar & elastic_displacement() const
Definition: LingleClark.cc:330
std::shared_ptr< Grid > m_extended_grid
extended grid for the viscous plate displacement
Definition: LingleClark.hh:85
std::string m_time_name
Name of the variable used to store the last update time.
Definition: LingleClark.hh:104
virtual void write_model_state_impl(const File &output) const
The default (empty implementation).
Definition: LingleClark.cc:421
std::shared_ptr< array::Scalar > m_viscous_displacement
Viscous displacement on the extended grid (part of the model state).
Definition: LingleClark.hh:88
double m_t_last
time of the last bed deformation update
Definition: LingleClark.hh:98
array::Scalar m_relief
Bed relief relative to the bed displacement.
Definition: LingleClark.hh:76
virtual void define_model_state_impl(const File &output) const
The default (empty implementation).
Definition: LingleClark.cc:406
std::shared_ptr< petsc::Vec > m_viscous_displacement0
rank 0 storage using the extended grid
Definition: LingleClark.hh:90
const array::Scalar & relief() const
Definition: LingleClark.cc:334
LingleClark(std::shared_ptr< const Grid > g)
Definition: LingleClark.cc:37
DiagnosticList diagnostics_impl() const
Definition: LingleClark.cc:430
std::unique_ptr< LingleClarkSerial > m_serial_model
Serial viscoelastic bed deformation model.
Definition: LingleClark.hh:82
const array::Scalar & total_displacement() const
Definition: LingleClark.cc:322
array::Scalar m_total_displacement
Total (viscous and elastic) bed displacement.
Definition: LingleClark.hh:69
const array::Scalar & viscous_displacement() const
Definition: LingleClark.cc:326
void update_impl(const array::Scalar &ice_thickness, const array::Scalar &sea_level_elevation, double t, double dt)
Update the Lingle-Clark bed deformation model.
Definition: LingleClark.cc:387
double m_update_interval
Update interval in seconds.
Definition: LingleClark.hh:100
std::shared_ptr< array::Scalar > elastic_load_response_matrix() const
Definition: LingleClark.cc:188
array::Scalar m_load_thickness
Ice-equivalent load thickness.
Definition: LingleClark.hh:79
A wrapper class around LingleClarkSerial.
Definition: LingleClark.hh:33
static const double g
Definition: exactTestP.cc:36
std::map< std::string, Diagnostic::Ptr > DiagnosticList
Definition: Diagnostic.hh:125