PISM, A Parallel Ice Sheet Model  stable v2.1.1 committed by Constantine Khrulev on 2024-12-04 13:36:58 -0900
Pico.hh
Go to the documentation of this file.
1 // Copyright (C) 2012-2016, 2018, 2020, 2021, 2022, 2023 Ricarda Winkelmann, Ronja Reese, Torsten Albrecht
2 // 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 #ifndef _POPICO_H_
21 #define _POPICO_H_
22 
23 #include "pism/coupler/ocean/CompleteOceanModel.hh"
24 
25 #include "pism/coupler/ocean/PicoGeometry.hh"
26 
27 namespace pism {
28 
29 namespace ocean {
30 
31 class PicoPhysics;
32 
33 //! Implements the PICO ocean model as submitted to The Cryosphere (March 2017).
34 //!
35 //! Generalizes the two dimensional ocean box model of [@ref OlbersHellmer2010] for
36 //! use in PISM, i.e. three dimensions.
37 //!
38 class Pico : public CompleteOceanModel {
39 public:
40  Pico(std::shared_ptr<const Grid> g);
41  virtual ~Pico() = default;
42 
43 protected:
44  void init_impl(const Geometry &geometry);
45  void update_impl(const Geometry &geometry, double t, double dt);
46  MaxTimestep max_timestep_impl(double t) const;
47 
48  void define_model_state_impl(const File &output) const;
49  void write_model_state_impl(const File &output) const;
50 
51  std::map<std::string, Diagnostic::Ptr> diagnostics_impl() const;
52 
53 private:
58 
60 
61  std::shared_ptr<array::Forcing> m_theta_ocean, m_salinity_ocean;
62 
63  void compute_ocean_input_per_basin(const PicoPhysics &physics,
64  const array::Scalar &basin_mask,
65  const array::Scalar &continental_shelf_mask,
66  const array::Scalar &salinity_ocean,
67  const array::Scalar &theta_ocean,
68  std::vector<double> &temperature,
69  std::vector<double> &salinity) const;
70 
71  void set_ocean_input_fields(const PicoPhysics &physics,
72  const array::Scalar &ice_thickness,
73  const array::CellType1 &mask,
74  const array::Scalar &basin_mask,
75  const array::Scalar &shelf_mask,
76  const std::vector<double> &basin_temperature,
77  const std::vector<double> &basin_salinity,
78  array::Scalar &Toc_box0,
79  array::Scalar &Soc_box0) const;
80 
81  void process_box1(const PicoPhysics &physics,
82  const array::Scalar &ice_thickness,
83  const array::Scalar &shelf_mask,
84  const array::Scalar &box_mask,
85  const array::Scalar &Toc_box0,
86  const array::Scalar &Soc_box0,
87  array::Scalar &basal_melt_rate,
88  array::Scalar &basal_temperature,
89  array::Scalar &T_star,
90  array::Scalar &Toc,
91  array::Scalar &Soc,
92  array::Scalar &overturning);
93 
94  void process_other_boxes(const PicoPhysics &physics,
95  const array::Scalar &ice_thickness,
96  const array::Scalar &shelf_mask,
97  const array::Scalar &box_mask,
98  array::Scalar &basal_melt_rate,
99  array::Scalar &basal_temperature,
100  array::Scalar &T_star,
101  array::Scalar &Toc,
102  array::Scalar &Soc) const;
103 
104  void beckmann_goosse(const PicoPhysics &physics,
105  const array::Scalar &ice_thickness,
106  const array::Scalar &shelf_mask,
107  const array::CellType &cell_type,
108  const array::Scalar &Toc_box0,
109  const array::Scalar &Soc_box0,
110  array::Scalar &basal_melt_rate,
111  array::Scalar &basal_temperature,
112  array::Scalar &Toc,
113  array::Scalar &Soc);
114 
115  void compute_box_average(int box_id,
116  const array::Scalar &field,
117  const array::Scalar &shelf_mask,
118  const array::Scalar &box_mask,
119  std::vector<double> &result) const;
120 
121  void compute_box_area(int box_id,
122  const array::Scalar &shelf_mask,
123  const array::Scalar &box_mask,
124  std::vector<double> &result) const;
125 
127 };
128 
129 } // end of namespace ocean
130 } // end of namespace pism
131 
132 #endif /* _POPICO_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
"Cell type" mask. Adds convenience methods to array::Scalar.
Definition: CellType.hh:30
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
virtual ~Pico()=default
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
bool ocean(int M)
An ocean cell (floating ice or ice-free).
Definition: Mask.hh:40
static const double g
Definition: exactTestP.cc:36