PISM, A Parallel Ice Sheet Model  stable v2.1-1-g6902d5502 committed by Ed Bueler on 2023-12-20 08:38:27 -0800
Array3D.cc
Go to the documentation of this file.
1 // Copyright (C) 2008--2018, 2020, 2021, 2022, 2023 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 
37 namespace pism {
38 namespace array {
39 
40 // this file contains method for derived class array::Array3D
41 
42 Array3D::Array3D(std::shared_ptr<const Grid> grid, const std::string &name, Kind ghostedp,
43  const std::vector<double> &levels, unsigned int stencil_width)
44  : Array(grid, name, ghostedp, 1, stencil_width, levels) {
46 }
47 
48 //! Set all values of scalar quantity to given a single value in a particular column.
49 void Array3D::set_column(int i, int j, double c) {
50  PetscErrorCode ierr;
51 #if (Pism_DEBUG == 1)
52  assert(m_array != NULL);
53  check_array_indices(i, j, 0);
54 #endif
55 
56  double ***arr = (double ***)m_array;
57 
58  if (c == 0.0) {
59  ierr = PetscMemzero(arr[j][i], levels().size() * sizeof(double));
60  PISM_CHK(ierr, "PetscMemzero");
61  } else {
62  unsigned int nlevels = levels().size();
63  for (unsigned int k = 0; k < nlevels; k++) {
64  arr[j][i][k] = c;
65  }
66  }
67 }
68 
69 void Array3D::set_column(int i, int j, const double *input) {
70 #if (Pism_DEBUG == 1)
71  check_array_indices(i, j, 0);
72 #endif
73  double ***arr = (double ***)m_array;
74  PetscErrorCode ierr = PetscMemcpy(arr[j][i], input, m_impl->zlevels.size() * sizeof(double));
75  PISM_CHK(ierr, "PetscMemcpy");
76 }
77 
78 #if (Pism_DEBUG == 1)
79 static bool legal_level(const std::vector<double> &levels, double z) {
80  double z_min = levels.front();
81  double z_max = levels.back();
82  const double eps = 1.0e-6;
83  return not(z < z_min - eps || z > z_max + eps);
84 }
85 #endif
86 
87 //! Return value of scalar quantity at level z (m) above base of ice (by linear interpolation).
88 double Array3D::interpolate(int i, int j, double z) const {
89  const auto &zs = levels();
90  auto N = zs.size();
91 
92 #if (Pism_DEBUG == 1)
93  assert(m_array != NULL);
94  check_array_indices(i, j, 0);
95 
96  if (not legal_level(zs, z)) {
98  "Array3D interpolate(): level %f is not legal; name = %s", z,
99  m_impl->name.c_str());
100  }
101 #endif
102 
103  const auto *column = get_column(i, j);
104 
105  if (z >= zs[N - 1]) {
106  return column[N - 1];
107  }
108 
109  if (z <= zs[0]) {
110  return column[0];
111  }
112 
113  auto mcurr = gsl_interp_accel_find(m_impl->bsearch_accel, zs.data(), N, z);
114 
115  const double incr = (z - zs[mcurr]) / (zs[mcurr + 1] - zs[mcurr]);
116  const double valm = column[mcurr];
117  return valm + incr * (column[mcurr + 1] - valm);
118 }
119 
120 double *Array3D::get_column(int i, int j) {
121 #if (Pism_DEBUG == 1)
122  check_array_indices(i, j, 0);
123 #endif
124  return ((double ***)m_array)[j][i];
125 }
126 
127 const double *Array3D::get_column(int i, int j) const {
128 #if (Pism_DEBUG == 1)
129  check_array_indices(i, j, 0);
130 #endif
131  return ((double ***)m_array)[j][i];
132 }
133 
134 //! Copies a horizontal slice at level z of an Array3D into `output`.
135 void extract_surface(const Array3D &data, double z, Scalar &output) {
136  array::AccessScope list{ &data, &output };
137 
138  ParallelSection loop(output.grid()->com);
139  try {
140  for (auto p = output.grid()->points(); p; p.next()) {
141  const int i = p.i(), j = p.j();
142  output(i, j) = data.interpolate(i, j, z);
143  }
144  } catch (...) {
145  loop.failed();
146  }
147  loop.check();
148 }
149 
150 
151 //! Copies the values of an Array3D at the ice surface (specified by the level `z`) into
152 //! `output`.
153 void extract_surface(const Array3D &data, const Scalar &z, Scalar &output) {
154  array::AccessScope list{ &data, &output, &z };
155 
156  ParallelSection loop(output.grid()->com);
157  try {
158  for (auto p = output.grid()->points(); p; p.next()) {
159  const int i = p.i(), j = p.j();
160  output(i, j) = data.interpolate(i, j, z(i, j));
161  }
162  } catch (...) {
163  loop.failed();
164  }
165  loop.check();
166 }
167 
168 /** Sum a 3-D vector in the Z direction to create a 2-D vector.
169 
170 Note that this sums up all the values in a column, including ones
171 above the ice. This may or may not be what you need. Also, take a look
172 at IceModel::compute_ice_enthalpy(PetscScalar &result) in iMreport.cc.
173 
174 As for the difference between array::Array2D and array::Scalar, the
175 former can store fields with more than 1 "degree of freedom" per grid
176 point (such as 2D fields on the "staggered" grid, with the first
177 degree of freedom corresponding to the i-offset and second to
178 j-offset).
179 
180 array::Scalar is just array::Array2D with "dof == 1", and
181 array::Vector is array::Array2D with "dof == 2". (Plus some extra
182 methods, of course.)
183 
184 Either one of array::Array2D and array::Scalar would work in this
185 case.
186 
187 Computes output = A*output + B*sum_columns(input) + C
188 
189 @see https://github.com/pism/pism/issues/229 */
190 void sum_columns(const Array3D &data, double A, double B, Scalar &output) {
191  const unsigned int Mz = data.grid()->Mz();
192 
193  AccessScope access{ &data, &output };
194 
195  for (auto p = data.grid()->points(); p; p.next()) {
196  const int i = p.i(), j = p.j();
197 
198  const double *column = data.get_column(i, j);
199 
200  double scalar_sum = 0.0;
201  for (unsigned int k = 0; k < Mz; ++k) {
202  scalar_sum += column[k];
203  }
204  output(i, j) = A * output(i, j) + B * scalar_sum;
205  }
206 }
207 
208 void Array3D::copy_from(const Array3D &input) {
209  array::AccessScope list{ this, &input };
210 
211  assert(levels().size() == input.levels().size());
212  assert(ndof() == input.ndof());
213 
214  // 3D arrays have more than one level and ndof() of 1, collections of fields have one
215  // level and ndof() > 1
216  auto N = std::max((size_t)ndof(), levels().size());
217 
218  ParallelSection loop(m_impl->grid->com);
219  try {
220  for (auto p = m_impl->grid->points(); p; p.next()) {
221  const int i = p.i(), j = p.j();
222 
223 #if PETSC_VERSION_LT(3, 12, 0)
224  PetscMemmove(this->get_column(i, j), const_cast<double *>(input.get_column(i, j)),
225  N * sizeof(double));
226 #else
227  PetscArraymove(this->get_column(i, j), input.get_column(i, j), N);
228 #endif
229  }
230  } catch (...) {
231  loop.failed();
232  }
233  loop.check();
234 
235  update_ghosts();
236 
238 }
239 
240 std::shared_ptr<Array3D> Array3D::duplicate(Kind ghostedp) const {
241 
242  auto result =
243  std::make_shared<Array3D>(this->grid(), this->get_name(), ghostedp, this->levels());
244 
245  result->metadata() = this->metadata();
246 
247  return result;
248 }
249 
250 void Array3D::regrid_impl(const File &file, io::Default default_value) {
251 
252  auto log = grid()->ctx()->log();
253 
254  auto variable = metadata(0);
255 
256  auto V = file.find_variable(variable.get_name(), variable["standard_name"]);
257 
259 
260  if (not V.exists) {
261  set_default_value_or_stop(file.filename(), variable, default_value, *log, tmp);
262  } else {
263  grid::InputGridInfo input_grid(file, V.name, variable.unit_system(), grid()->registration());
264 
265  input_grid.report(*log, 4, variable.unit_system());
266 
267  io::check_input_grid(input_grid, *grid(), levels());
268 
269  LocalInterpCtx lic(input_grid, *grid(), levels(), m_impl->interpolation_type);
270 
271  // Note: this call will read the last time record (the index is set in `lic` based on
272  // info in `input_grid`).
273  petsc::VecArray tmp_array(tmp);
274  io::regrid_spatial_variable(variable, *grid(), lic, file, tmp_array.get());
275  }
276 
277  if (m_impl->ghosted) {
278  global_to_local(*dm(), tmp, vec());
279  } else {
280  PetscErrorCode ierr = VecCopy(tmp, vec());
281  PISM_CHK(ierr, "VecCopy");
282  }
283 }
284 
285 } // end of namespace array
286 } // 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:361
std::string filename() const
Definition: File.cc:307
High-level PISM I/O class.
Definition: File.hh:56
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:65
void copy_from(const Array3D &input)
Definition: Array3D.cc:208
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:88
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:49
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:250
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:42
std::shared_ptr< Array3D > duplicate(Kind ghostedp=WITHOUT_GHOSTS) const
Definition: Array3D.cc:240
double * get_column(int i, int j)
Definition: Array3D.cc:120
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:136
const std::vector< double > & levels() const
Definition: Array.cc:140
size_t size() const
Return the total number of elements in the owned part of an array.
Definition: Array.cc:946
void * m_array
Definition: Array.hh:278
const std::string & get_name() const
Get the name of an Array object.
Definition: Array.cc:383
petsc::Vec & vec() const
Definition: Array.cc:339
std::shared_ptr< const Grid > grid() const
Definition: Array.cc:132
void inc_state_counter()
Increment the object state counter.
Definition: Array.cc:151
std::shared_ptr< petsc::DM > dm() const
Definition: Array.cc:353
void set_begin_access_use_dof(bool flag)
Definition: Array.cc:175
void update_ghosts()
Updates ghost points.
Definition: Array.cc:693
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:714
SpatialVariableMetadata & metadata(unsigned int N=0)
Returns a reference to the SpatialVariableMetadata object containing metadata for the compoment N.
Definition: Array.cc:553
Abstract class for reading, writing, allocating, and accessing a DA-based PETSc Vec (2D and 3D fields...
Definition: Array.hh:208
void report(const Logger &log, int threshold, std::shared_ptr< units::System > s) const
Definition: Grid.cc:1025
Contains parameters of an input file grid.
Definition: Grid.hh:65
double * get()
Definition: Vec.cc:54
Wrapper around VecGetArray and VecRestoreArray.
Definition: Vec.hh:44
#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:135
void sum_columns(const Array3D &data, double A, double B, Scalar &output)
Definition: Array3D.cc:190
void set_default_value_or_stop(const std::string &filename, const VariableMetadata &variable, io::Default default_value, const Logger &log, Vec output)
Definition: Array.cc:387
void global_to_local(petsc::DM &dm, Vec source, Vec destination)
Definition: Array.cc:49
double max(const array::Scalar &input)
Finds maximum over all the values in an array::Scalar object. Ignores ghosts.
Definition: Scalar.cc:165
Kind
What "kind" of a vector to create: with or without ghosts.
Definition: Array.hh:62
static Vector3 column(const double A[3][3], size_t k)
Definition: Element.cc:51
void regrid_spatial_variable(SpatialVariableMetadata &variable, const Grid &internal_grid, const LocalInterpCtx &lic, const File &file, double *output)
Regrid from a NetCDF file into a distributed array output.
Definition: io_helpers.cc:793
void check_input_grid(const grid::InputGridInfo &input_grid, const Grid &internal_grid, const std::vector< double > &internal_z_levels)
Check that x, y, and z coordinates of the input grid are strictly increasing.
Definition: io_helpers.cc:757
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
Definition: Array_impl.hh:105
gsl_interp_accel * bsearch_accel
Definition: Array_impl.hh:111
std::vector< double > zlevels
Vertical levels (for 3D fields)
Definition: Array_impl.hh:108