PISM, A Parallel Ice Sheet Model 2.2.1-cd005eec8 committed by Constantine Khrulev on 2025-03-07
Loading...
Searching...
No Matches
Array3D.hh
Go to the documentation of this file.
1/* Copyright (C) 2022, 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 PISM_ARRAY3D_H
21#define PISM_ARRAY3D_H
22
23#include "pism/util/array/Array.hh"
24
25namespace pism {
26
27namespace array {
28
29class Scalar;
30
31//! \brief A virtual class collecting methods common to ice and bedrock 3D
32//! fields.
33class Array3D : public Array {
34public:
35
36 // Three-dimensional array with a number of vertical levels
37 Array3D(std::shared_ptr<const Grid> grid,
38 const std::string &name,
39 Kind ghostedp,
40 const std::vector<double> &levels,
41 unsigned int stencil_width = 1);
42
43 virtual ~Array3D() = default;
44
45 std::shared_ptr<Array3D> duplicate(Kind ghostedp = WITHOUT_GHOSTS) const;
46
47 void regrid_impl(const File &file, io::Default default_value);
48
49 void set_column(int i, int j, double c);
50 void set_column(int i, int j, const double *input);
51 double* get_column(int i, int j);
52 const double* get_column(int i, int j) const;
53
54 double interpolate(int i, int j, double z) const;
55
56 void copy_from(const Array3D &input);
57};
58
59void extract_surface(const Array3D &data, double z, Scalar &output);
60void extract_surface(const Array3D &data, const Scalar &z, Scalar &output);
61
62void sum_columns(const Array3D &data, double A, double B, Scalar &output);
63
64} // end of namespace array
65} // end of namespace pism
66
67#endif /* PISM_ARRAY3D_H */
High-level PISM I/O class.
Definition File.hh:55
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
virtual ~Array3D()=default
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
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
const std::vector< double > & levels() const
Definition Array.cc:139
std::shared_ptr< const Grid > grid() const
Definition Array.cc:131
unsigned int stencil_width() const
Get the stencil width of the current Array. Returns 0 if ghosts are not available.
Definition Array.cc:302
Abstract class for reading, writing, allocating, and accessing a DA-based PETSc Vec (2D and 3D fields...
Definition Array.hh:207
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
Kind
What "kind" of a vector to create: with or without ghosts.
Definition Array.hh:61
@ WITHOUT_GHOSTS
Definition Array.hh:61