25#include "pism/util/array/Array3D.hh"
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"
44 const std::vector<double> &levels,
unsigned int stencil_width)
45 :
Array(grid, name, ghostedp, 1, stencil_width, levels) {
57 double ***arr = (
double ***)
m_array;
60 ierr = PetscMemzero(arr[j][i],
levels().
size() *
sizeof(
double));
63 unsigned int nlevels =
levels().size();
64 for (
unsigned int k = 0;
k < nlevels;
k++) {
74 double ***arr = (
double ***)
m_array;
75 PetscErrorCode ierr = PetscMemcpy(arr[j][i], input,
m_impl->
zlevels.size() *
sizeof(
double));
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);
97 if (not legal_level(zs, z)) {
99 "Array3D interpolate(): level %f is not legal; name = %s", z,
106 if (z >= zs[N - 1]) {
107 return column[N - 1];
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);
125 return ((
double ***)
m_array)[j][i];
132 return ((
double ***)
m_array)[j][i];
141 for (
auto p = output.
grid()->points(); p; p.next()) {
142 const int i = p.i(), j = p.j();
159 for (
auto p = output.
grid()->points(); p; p.next()) {
160 const int i = p.i(), j = p.j();
192 const unsigned int Mz = data.
grid()->Mz();
196 for (
auto p = data.
grid()->points(); p; p.next()) {
197 const int i = p.i(), j = p.j();
201 double scalar_sum = 0.0;
202 for (
unsigned int k = 0;
k < Mz; ++
k) {
203 scalar_sum += column[
k];
205 output(i, j) = A * output(i, j) + B * scalar_sum;
221 for (
auto p =
m_impl->
grid->points(); p; p.next()) {
222 const int i = p.i(), j = p.j();
224#if PETSC_VERSION_LT(3, 12, 0)
246 result->metadata() = this->
metadata();
253 auto log =
grid()->ctx()->log();
257 auto V = file.
find_variable(variable.get_name(), variable[
"standard_name"]);
264 int last_record = -1;
265 interp->regrid(variable, file, last_record, *
grid(), tmp);
273 PetscErrorCode ierr = VecCopy(tmp,
vec());
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.
High-level PISM I/O class.
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.
void copy_from(const Array3D &input)
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).
void set_column(int i, int j, double c)
Set all values of scalar quantity to given a single value in a particular column.
void regrid_impl(const File &file, io::Default default_value)
Gets an Array from a file file, interpolating onto the current grid.
Array3D(std::shared_ptr< const Grid > grid, const std::string &name, Kind ghostedp, const std::vector< double > &levels, unsigned int stencil_width=1)
std::shared_ptr< Array3D > duplicate(Kind ghostedp=WITHOUT_GHOSTS) const
double * get_column(int i, int j)
A virtual class collecting methods common to ice and bedrock 3D fields.
unsigned int ndof() const
Returns the number of degrees of freedom per grid point.
const std::vector< double > & levels() const
size_t size() const
Return the total number of elements in the owned part of an array.
const std::string & get_name() const
Get the name of an Array object.
std::shared_ptr< const Grid > grid() const
void inc_state_counter()
Increment the object state counter.
std::shared_ptr< petsc::DM > dm() const
void set_begin_access_use_dof(bool flag)
void update_ghosts()
Updates ghost points.
void check_array_indices(int i, int j, unsigned int k) const
Check array indices and warn if they are out of range.
SpatialVariableMetadata & metadata(unsigned int N=0)
Returns a reference to the SpatialVariableMetadata object containing metadata for the compoment N.
Abstract class for reading, writing, allocating, and accessing a DA-based PETSc Vec (2D and 3D fields...
#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.
void sum_columns(const Array3D &data, double A, double B, Scalar &output)
void global_to_local(petsc::DM &dm, Vec source, Vec destination)
void set_default_value_or_stop(const VariableMetadata &variable, io::Default default_value, const Logger &log, Vec output)
Kind
What "kind" of a vector to create: with or without ghosts.
bool ghosted
true if this Array is ghosted
std::shared_ptr< const Grid > grid
The computational grid.
InterpolationType interpolation_type
gsl_interp_accel * bsearch_accel
std::vector< double > zlevels
Vertical levels (for 3D fields)