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
Array3D.cc
Go to the documentation of this file.
1// Copyright (C) 2008--2018, 2020, 2021, 2022, 2023, 2024 Ed Bueler and Constantine Khroulev
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 <cstring>
20#include <cstdlib>
21#include <cassert>
22#include <memory>
23#include <petscvec.h>
24
25#include "pism/util/array/Array3D.hh"
26
27#include "pism/util/ConfigInterface.hh"
28#include "pism/util/Context.hh"
29#include "pism/util/Grid.hh"
30#include "pism/util/VariableMetadata.hh"
31#include "pism/util/array/Array_impl.hh"
32#include "pism/util/array/Scalar.hh"
33#include "pism/util/error_handling.hh"
34#include "pism/util/io/IO_Flags.hh"
35#include "pism/util/io/io_helpers.hh"
36#include "pism/util/InputInterpolation.hh"
37
38namespace pism {
39namespace array {
40
41// this file contains method for derived class array::Array3D
42
43Array3D::Array3D(std::shared_ptr<const Grid> grid, const std::string &name, Kind ghostedp,
44 const std::vector<double> &levels, unsigned int stencil_width)
45 : Array(grid, name, ghostedp, 1, stencil_width, levels) {
47}
48
49//! Set all values of scalar quantity to given a single value in a particular column.
50void Array3D::set_column(int i, int j, double c) {
51 PetscErrorCode ierr;
52#if (Pism_DEBUG == 1)
53 assert(m_array != NULL);
54 check_array_indices(i, j, 0);
55#endif
56
57 double ***arr = (double ***)m_array;
58
59 if (c == 0.0) {
60 ierr = PetscMemzero(arr[j][i], levels().size() * sizeof(double));
61 PISM_CHK(ierr, "PetscMemzero");
62 } else {
63 unsigned int nlevels = levels().size();
64 for (unsigned int k = 0; k < nlevels; k++) {
65 arr[j][i][k] = c;
66 }
67 }
68}
69
70void Array3D::set_column(int i, int j, const double *input) {
71#if (Pism_DEBUG == 1)
72 check_array_indices(i, j, 0);
73#endif
74 double ***arr = (double ***)m_array;
75 PetscErrorCode ierr = PetscMemcpy(arr[j][i], input, m_impl->zlevels.size() * sizeof(double));
76 PISM_CHK(ierr, "PetscMemcpy");
77}
78
79#if (Pism_DEBUG == 1)
80static bool legal_level(const std::vector<double> &levels, double z) {
81 double z_min = levels.front();
82 double z_max = levels.back();
83 const double eps = 1.0e-6;
84 return not(z < z_min - eps || z > z_max + eps);
85}
86#endif
87
88//! Return value of scalar quantity at level z (m) above base of ice (by linear interpolation).
89double Array3D::interpolate(int i, int j, double z) const {
90 const auto &zs = levels();
91 auto N = zs.size();
92
93#if (Pism_DEBUG == 1)
94 assert(m_array != NULL);
95 check_array_indices(i, j, 0);
96
97 if (not legal_level(zs, z)) {
99 "Array3D interpolate(): level %f is not legal; name = %s", z,
100 m_impl->name.c_str());
101 }
102#endif
103
104 const auto *column = get_column(i, j);
105
106 if (z >= zs[N - 1]) {
107 return column[N - 1];
108 }
109
110 if (z <= zs[0]) {
111 return column[0];
112 }
113
114 auto mcurr = gsl_interp_accel_find(m_impl->bsearch_accel, zs.data(), N, z);
115
116 const double incr = (z - zs[mcurr]) / (zs[mcurr + 1] - zs[mcurr]);
117 const double valm = column[mcurr];
118 return valm + incr * (column[mcurr + 1] - valm);
119}
120
121double *Array3D::get_column(int i, int j) {
122#if (Pism_DEBUG == 1)
123 check_array_indices(i, j, 0);
124#endif
125 return ((double ***)m_array)[j][i];
126}
127
128const double *Array3D::get_column(int i, int j) const {
129#if (Pism_DEBUG == 1)
130 check_array_indices(i, j, 0);
131#endif
132 return ((double ***)m_array)[j][i];
133}
134
135//! Copies a horizontal slice at level z of an Array3D into `output`.
136void extract_surface(const Array3D &data, double z, Scalar &output) {
137 array::AccessScope list{ &data, &output };
138
139 ParallelSection loop(output.grid()->com);
140 try {
141 for (auto p = output.grid()->points(); p; p.next()) {
142 const int i = p.i(), j = p.j();
143 output(i, j) = data.interpolate(i, j, z);
144 }
145 } catch (...) {
146 loop.failed();
147 }
148 loop.check();
149}
150
151
152//! Copies the values of an Array3D at the ice surface (specified by the level `z`) into
153//! `output`.
154void extract_surface(const Array3D &data, const Scalar &z, Scalar &output) {
155 array::AccessScope list{ &data, &output, &z };
156
157 ParallelSection loop(output.grid()->com);
158 try {
159 for (auto p = output.grid()->points(); p; p.next()) {
160 const int i = p.i(), j = p.j();
161 output(i, j) = data.interpolate(i, j, z(i, j));
162 }
163 } catch (...) {
164 loop.failed();
165 }
166 loop.check();
167}
168
169/** Sum a 3-D vector in the Z direction to create a 2-D vector.
170
171Note that this sums up all the values in a column, including ones
172above the ice. This may or may not be what you need. Also, take a look
173at IceModel::compute_ice_enthalpy(PetscScalar &result) in iMreport.cc.
174
175As for the difference between array::Array2D and array::Scalar, the
176former can store fields with more than 1 "degree of freedom" per grid
177point (such as 2D fields on the "staggered" grid, with the first
178degree of freedom corresponding to the i-offset and second to
179j-offset).
180
181array::Scalar is just array::Array2D with "dof == 1", and
182array::Vector is array::Array2D with "dof == 2". (Plus some extra
183methods, of course.)
184
185Either one of array::Array2D and array::Scalar would work in this
186case.
187
188Computes output = A*output + B*sum_columns(input) + C
189
190@see https://github.com/pism/pism/issues/229 */
191void sum_columns(const Array3D &data, double A, double B, Scalar &output) {
192 const unsigned int Mz = data.grid()->Mz();
193
194 AccessScope access{ &data, &output };
195
196 for (auto p = data.grid()->points(); p; p.next()) {
197 const int i = p.i(), j = p.j();
198
199 const double *column = data.get_column(i, j);
200
201 double scalar_sum = 0.0;
202 for (unsigned int k = 0; k < Mz; ++k) {
203 scalar_sum += column[k];
204 }
205 output(i, j) = A * output(i, j) + B * scalar_sum;
206 }
207}
208
209void Array3D::copy_from(const Array3D &input) {
210 array::AccessScope list{ this, &input };
211
212 assert(levels().size() == input.levels().size());
213 assert(ndof() == input.ndof());
214
215 // 3D arrays have more than one level and ndof() of 1, collections of fields have one
216 // level and ndof() > 1
217 auto N = std::max((size_t)ndof(), levels().size());
218
219 ParallelSection loop(m_impl->grid->com);
220 try {
221 for (auto p = m_impl->grid->points(); p; p.next()) {
222 const int i = p.i(), j = p.j();
223
224#if PETSC_VERSION_LT(3, 12, 0)
225 PetscMemmove(this->get_column(i, j), const_cast<double *>(input.get_column(i, j)),
226 N * sizeof(double));
227#else
228 PetscArraymove(this->get_column(i, j), input.get_column(i, j), N);
229#endif
230 }
231 } catch (...) {
232 loop.failed();
233 }
234 loop.check();
235
237
239}
240
241std::shared_ptr<Array3D> Array3D::duplicate(Kind ghostedp) const {
242
243 auto result =
244 std::make_shared<Array3D>(this->grid(), this->get_name(), ghostedp, this->levels());
245
246 result->metadata() = this->metadata();
247
248 return result;
249}
250
251void Array3D::regrid_impl(const File &file, io::Default default_value) {
252
253 auto log = grid()->ctx()->log();
254
255 auto variable = metadata(0);
256
257 auto V = file.find_variable(variable.get_name(), variable["standard_name"]);
258
260
261 if (V.exists) {
262 auto interp = grid()->get_interpolation(levels(), file, V.name, m_impl->interpolation_type);
263
264 int last_record = -1;
265 interp->regrid(variable, file, last_record, *grid(), tmp);
266 } else {
267 set_default_value_or_stop(variable, default_value, *log, tmp);
268 }
269
270 if (m_impl->ghosted) {
271 global_to_local(*dm(), tmp, vec());
272 } else {
273 PetscErrorCode ierr = VecCopy(tmp, vec());
274 PISM_CHK(ierr, "VecCopy");
275 }
276}
277
278} // end of namespace array
279} // end of namespace pism
VariableLookupData find_variable(const std::string &short_name, const std::string &std_name) const
Find a variable using its standard name and/or short name.
Definition File.cc:328
std::string name() const
Definition File.cc:274
High-level PISM I/O class.
Definition File.hh:55
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
Makes sure that we call begin_access() and end_access() for all accessed array::Arrays.
Definition Array.hh:64
void copy_from(const Array3D &input)
Definition Array3D.cc:209
double interpolate(int i, int j, double z) const
Return value of scalar quantity at level z (m) above base of ice (by linear interpolation).
Definition Array3D.cc:89
void set_column(int i, int j, double c)
Set all values of scalar quantity to given a single value in a particular column.
Definition Array3D.cc:50
void regrid_impl(const File &file, io::Default default_value)
Gets an Array from a file file, interpolating onto the current grid.
Definition Array3D.cc:251
Array3D(std::shared_ptr< const Grid > grid, const std::string &name, Kind ghostedp, const std::vector< double > &levels, unsigned int stencil_width=1)
Definition Array3D.cc:43
std::shared_ptr< Array3D > duplicate(Kind ghostedp=WITHOUT_GHOSTS) const
Definition Array3D.cc:241
double * get_column(int i, int j)
Definition Array3D.cc:121
A virtual class collecting methods common to ice and bedrock 3D fields.
Definition Array3D.hh:33
unsigned int ndof() const
Returns the number of degrees of freedom per grid point.
Definition Array.cc:135
const std::vector< double > & levels() const
Definition Array.cc:139
size_t size() const
Return the total number of elements in the owned part of an array.
Definition Array.cc:909
const std::string & get_name() const
Get the name of an Array object.
Definition Array.cc:354
petsc::Vec & vec() const
Definition Array.cc:310
std::shared_ptr< const Grid > grid() const
Definition Array.cc:131
void inc_state_counter()
Increment the object state counter.
Definition Array.cc:150
std::shared_ptr< petsc::DM > dm() const
Definition Array.cc:324
void set_begin_access_use_dof(bool flag)
Definition Array.cc:174
void update_ghosts()
Updates ghost points.
Definition Array.cc:615
void check_array_indices(int i, int j, unsigned int k) const
Check array indices and warn if they are out of range.
Definition Array.cc:636
SpatialVariableMetadata & metadata(unsigned int N=0)
Returns a reference to the SpatialVariableMetadata object containing metadata for the compoment N.
Definition Array.cc:476
Abstract class for reading, writing, allocating, and accessing a DA-based PETSc Vec (2D and 3D fields...
Definition Array.hh:207
#define PISM_CHK(errcode, name)
#define PISM_ERROR_LOCATION
void extract_surface(const Array3D &data, double z, Scalar &output)
Copies a horizontal slice at level z of an Array3D into output.
Definition Array3D.cc:136
void sum_columns(const Array3D &data, double A, double B, Scalar &output)
Definition Array3D.cc:191
void global_to_local(petsc::DM &dm, Vec source, Vec destination)
Definition Array.cc:53
void set_default_value_or_stop(const VariableMetadata &variable, io::Default default_value, const Logger &log, Vec output)
Definition Array.cc:358
Kind
What "kind" of a vector to create: with or without ghosts.
Definition Array.hh:61
static const double k
Definition exactTestP.cc:42
bool ghosted
true if this Array is ghosted
Definition Array_impl.hh:86
std::shared_ptr< const Grid > grid
The computational grid.
Definition Array_impl.hh:77
InterpolationType interpolation_type
gsl_interp_accel * bsearch_accel
std::vector< double > zlevels
Vertical levels (for 3D fields)