PISM, A Parallel Ice Sheet Model 2.3.0-79cae578d committed by Constantine Khrulev on 2026-03-22
Loading...
Searching...
No Matches
Pico.hh
Go to the documentation of this file.
1// Copyright (C) 2012-2016, 2018, 2020, 2021, 2022, 2023, 2025 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
27namespace pism {
28
29class VariableMetadata;
30
31namespace ocean {
32
33class PicoPhysics;
34
35//! Implements the PICO ocean model as submitted to The Cryosphere (March 2017).
36//!
37//! Generalizes the two dimensional ocean box model of [@ref OlbersHellmer2010] for
38//! use in PISM, i.e. three dimensions.
39//!
40class Pico : public CompleteOceanModel {
41public:
42 Pico(std::shared_ptr<const Grid> g);
43 virtual ~Pico() = default;
44
45protected:
46 void init_impl(const Geometry &geometry);
47 void update_impl(const Inputs &inputs, double t, double dt);
48 MaxTimestep max_timestep_impl(double t) const;
49
50 std::set<VariableMetadata> state_impl() const;
51 void write_state_impl(const OutputFile &output) const;
52
53 std::map<std::string, Diagnostic::Ptr> spatial_diagnostics_impl() const;
54
55private:
60
62
63 std::shared_ptr<array::Forcing> m_theta_ocean, m_salinity_ocean;
64
66 const array::Scalar &basin_mask,
67 const array::Scalar &continental_shelf_mask,
68 const array::Scalar &salinity_ocean,
69 const array::Scalar &theta_ocean,
70 std::vector<double> &temperature,
71 std::vector<double> &salinity) const;
72
73 void set_ocean_input_fields(const PicoPhysics &physics,
74 const array::Scalar &ice_thickness,
75 const array::CellType1 &mask,
76 const array::Scalar &basin_mask,
77 const array::Scalar &shelf_mask,
78 const std::vector<double> &basin_temperature,
79 const std::vector<double> &basin_salinity,
80 array::Scalar &Toc_box0,
81 array::Scalar &Soc_box0) const;
82
83 void process_box1(const PicoPhysics &physics,
84 const array::Scalar &ice_thickness,
85 const array::Scalar &shelf_mask,
86 const array::Scalar &box_mask,
87 const array::Scalar &Toc_box0,
88 const array::Scalar &Soc_box0,
89 array::Scalar &basal_melt_rate,
90 array::Scalar &basal_temperature,
91 array::Scalar &T_star,
92 array::Scalar &Toc,
93 array::Scalar &Soc,
94 array::Scalar &overturning);
95
96 void process_other_boxes(const PicoPhysics &physics,
97 const array::Scalar &ice_thickness,
98 const array::Scalar &shelf_mask,
99 const array::Scalar &box_mask,
100 array::Scalar &basal_melt_rate,
101 array::Scalar &basal_temperature,
102 array::Scalar &T_star,
103 array::Scalar &Toc,
104 array::Scalar &Soc) const;
105
106 void beckmann_goosse(const PicoPhysics &physics,
107 const array::Scalar &ice_thickness,
108 const array::Scalar &shelf_mask,
109 const array::CellType &cell_type,
110 const array::Scalar &Toc_box0,
111 const array::Scalar &Soc_box0,
112 array::Scalar &basal_melt_rate,
113 array::Scalar &basal_temperature,
114 array::Scalar &Toc,
115 array::Scalar &Soc);
116
117 void compute_box_average(int box_id,
118 const array::Scalar &field,
119 const array::Scalar &shelf_mask,
120 const array::Scalar &box_mask,
121 std::vector<double> &result) const;
122
123 void compute_box_area(int box_id,
124 const array::Scalar &shelf_mask,
125 const array::Scalar &box_mask,
126 std::vector<double> &result) const;
127
129};
130
131} // end of namespace ocean
132} // end of namespace pism
133
134#endif /* _POPICO_H_ */
Combines the max. time step with the flag indicating if a restriction is active. Makes is possible to...
"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:460
array::Scalar m_Toc_box0
Definition Pico.hh:57
virtual ~Pico()=default
void write_state_impl(const OutputFile &output) const
The default (empty implementation).
Definition Pico.cc:183
std::shared_ptr< array::Forcing > m_salinity_ocean
Definition Pico.hh:63
array::Scalar m_Soc
Definition Pico.hh:56
std::map< std::string, Diagnostic::Ptr > spatial_diagnostics_impl() const
Definition Pico.cc:784
MaxTimestep max_timestep_impl(double t) const
Definition Pico.cc:360
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:627
array::Scalar m_Toc
Definition Pico.hh:57
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:863
array::Scalar1 m_basal_melt_rate
Definition Pico.hh:59
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:373
array::Scalar m_Soc_box0
Definition Pico.hh:56
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:810
std::shared_ptr< array::Forcing > m_theta_ocean
Definition Pico.hh:63
std::set< VariableMetadata > state_impl() const
Definition Pico.cc:176
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:696
PicoGeometry m_geometry
Definition Pico.hh:61
void init_impl(const Geometry &geometry)
Definition Pico.cc:131
void update_impl(const Inputs &inputs, double t, double dt)
Definition Pico.cc:244
array::Scalar m_T_star
Definition Pico.hh:57
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:582
array::Scalar m_overturning
Definition Pico.hh:58
bool ocean(int M)
An ocean cell (floating ice or ice-free).
Definition Mask.hh:40
static const double g
Definition exactTestP.cc:36