ThunderEgg  1.0.0
PatchArray.h
Go to the documentation of this file.
1 /***************************************************************************
2  * ThunderEgg, a library for solvers on adaptively refined block-structured
3  * Cartesian grids.
4  *
5  * Copyright (c) 2021 Scott Aiton
6  *
7  * This program is free software: you can redistribute it and/or modify
8  * it under the terms of the GNU General Public License as published by
9  * the Free Software Foundation, either version 3 of the License, or
10  * (at your option) any later version.
11  *
12  * This program is distributed in the hope that it will be useful,
13  * but WITHOUT ANY WARRANTY; without even the implied warranty of
14  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15  * GNU General Public License for more details.
16  *
17  * You should have received a copy of the GNU General Public License
18  * along with this program. If not, see <https://www.gnu.org/licenses/>.
19  ***************************************************************************/
20 
21 #ifndef THUNDEREGG_PATCHARRAY_H
22 #define THUNDEREGG_PATCHARRAY_H
23 
28 #include <ThunderEgg/PatchView.h>
29 namespace ThunderEgg {
35 template<int D>
37 {
38 private:
39  std::vector<double> vector;
41 
42 public:
47  PatchArray() = default;
55  PatchArray(const std::array<int, D>& lengths, int num_components, int num_ghost_cells)
56  {
57  std::array<int, D + 1> view_lengths;
58  for (int i = 0; i < D; i++) {
59  view_lengths[i] = lengths[i];
60  }
61  view_lengths[D] = num_components;
62  std::array<int, D + 1> strides;
63  strides[0] = 1;
64  for (int i = 1; i < D + 1; i++) {
65  strides[i] = strides[i - 1] * (lengths[i - 1] + 2 * num_ghost_cells);
66  }
67  int size = strides[D] * num_components;
68  vector.resize(size);
69  view = PatchView<double, D>(vector.data(), strides, view_lengths, num_ghost_cells);
70  }
71 
77  PatchArray(const PatchArray<D>& other)
78  : vector(other.vector)
79  , view(vector.data(),
80  other.getStrides(),
81  other.getGhostStart(),
82  other.getStart(),
83  other.getEnd(),
84  other.getGhostEnd())
85 
86  {}
87 
95  {
96  vector = other.vector;
97  view = PatchView<double, D>(vector.data(),
98  other.getStrides(),
99  other.getGhostStart(),
100  other.getStart(),
101  other.getEnd(),
102  other.getGhostEnd());
103  return *this;
104  }
105 
115  template<int M>
116  inline View<double, M + 1> getSliceOn(Face<D, M> f, const std::array<int, D - M>& offset)
117  {
118  return view.template getSliceOn<M>(f, offset);
119  }
120 
130  template<int M>
132  const std::array<int, D - M>& offset) const
133  {
134  return View<const double, M + 1>(view.template getSliceOn<M>(f, offset));
135  }
136 
146  template<int M>
148  const std::array<size_t, D - M>& offset) const
149  {
150  return view.template getGhostSliceOn<M>(f, offset);
151  }
152 
153  inline const double& operator[](const std::array<int, D + 1>& coord) const { return view[coord]; }
154  template<class... Types>
155  inline const double& operator()(Types... args) const
156  {
157  return view(args...);
158  }
159  inline void set(const std::array<int, D + 1>& coord, double value) const
160  {
161  view.set(coord, value);
162  }
163  inline double& operator[](const std::array<int, D + 1>& coord) { return view[coord]; }
164  template<class... Types>
165  inline double& operator()(Types... args)
166  {
167  return view(args...);
168  }
169  inline void set(const std::array<int, D + 1>& coord, double value) { view.set(coord, value); }
170 
174  inline const std::array<int, D + 1>& getStrides() const { return view.getStrides(); }
178  inline const std::array<int, D + 1>& getStart() const { return view.getStart(); }
182  inline const std::array<int, D + 1>& getEnd() const { return view.getEnd(); }
186  inline const std::array<int, D + 1>& getGhostStart() const { return view.getGhostStart(); }
190  inline const std::array<int, D + 1>& getGhostEnd() const { return view.getGhostEnd(); }
196  const PatchView<double, D>& getView() { return view; }
197 };
198 extern template class PatchArray<1>;
199 extern template class PatchArray<2>;
200 extern template class PatchArray<3>;
201 } // namespace ThunderEgg
202 #endif
ThunderEgg::View
Array for acessing data of a patch. It supports variable striding.
Definition: View.h:40
ThunderEgg::PatchArray::getGhostSliceOn
View< double, M+1 > getGhostSliceOn(Face< D, M > f, const std::array< size_t, D - M > &offset) const
Get the gosts slice on a given face.
Definition: PatchArray.h:147
ThunderEgg::View::getStrides
const std::array< int, D > & getStrides() const
Get the strides of the patch in each direction.
Definition: View.h:203
ThunderEgg::PatchArray::getStrides
const std::array< int, D+1 > & getStrides() const
Get the strides of the patch in each direction.
Definition: PatchArray.h:174
ThunderEgg::View::getGhostStart
const std::array< int, D > & getGhostStart() const
Get the coordinate of the first ghost cell element.
Definition: View.h:215
ThunderEgg::PatchArray::getSliceOn
View< double, M+1 > getSliceOn(Face< D, M > f, const std::array< int, D - M > &offset)
Get the slice on a given face.
Definition: PatchArray.h:116
ThunderEgg::PatchArray::getGhostStart
const std::array< int, D+1 > & getGhostStart() const
Get the coordinate of the first ghost cell element.
Definition: PatchArray.h:186
ThunderEgg::PatchArray::getView
const PatchView< double, D > & getView()
Get the View for the array.
Definition: PatchArray.h:196
PatchView.h
PatchView class.
ThunderEgg::PatchView< double, D >
ThunderEgg
The ThunderEgg namespace.
Definition: BiLinearGhostFiller.h:31
ThunderEgg::View::getEnd
const std::array< int, D > & getEnd() const
Get the coordinate of the last element.
Definition: View.h:211
ThunderEgg::PatchArray
Array for acessing data of a patch. It supports variable striding.
Definition: PatchArray.h:36
ThunderEgg::PatchArray::getEnd
const std::array< int, D+1 > & getEnd() const
Get the coordinate of the last element.
Definition: PatchArray.h:182
ThunderEgg::PatchArray::PatchArray
PatchArray(const std::array< int, D > &lengths, int num_components, int num_ghost_cells)
Construct a new PatchArray object.
Definition: PatchArray.h:55
ThunderEgg::View::set
void set(const std::array< int, D > &coord, T value) const
Set the value at a coordinate to the specified value.
Definition: View.h:174
ThunderEgg::PatchArray::PatchArray
PatchArray(const PatchArray< D > &other)
Copy constructor.
Definition: PatchArray.h:77
ThunderEgg::PatchArray::PatchArray
PatchArray()=default
Construct a new PatchArray object of size zero.
ThunderEgg::PatchArray::getGhostEnd
const std::array< int, D+1 > & getGhostEnd() const
Get the coordinate of the last ghost cell element.
Definition: PatchArray.h:190
ThunderEgg::PatchArray::operator=
PatchArray< D > & operator=(const PatchArray< D > &other)
Copy assignment.
Definition: PatchArray.h:94
ThunderEgg::PatchArray::getSliceOn
View< const double, M+1 > getSliceOn(Face< D, M > f, const std::array< int, D - M > &offset) const
Get the slice on a given face.
Definition: PatchArray.h:131
ThunderEgg::Face
Enum-style class for the faces of an n-dimensional cube.
Definition: Face.h:41
ThunderEgg::PatchArray::getStart
const std::array< int, D+1 > & getStart() const
Get the coordinate of the first element.
Definition: PatchArray.h:178
ThunderEgg::View::getGhostEnd
const std::array< int, D > & getGhostEnd() const
Get the coordinate of the last ghost cell element.
Definition: View.h:219
ThunderEgg::View::getStart
const std::array< int, D > & getStart() const
Get the coordinate of the first element.
Definition: View.h:207