ThunderEgg  1.0.0
Level.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) 2018-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_GMG_LEVEL_H
22 #define THUNDEREGG_GMG_LEVEL_H
23 
31 #include <ThunderEgg/Operator.h>
32 #include <memory>
33 namespace ThunderEgg::GMG {
37 template<int D>
38 class Level
39 {
40 private:
44  std::shared_ptr<const Operator<D>> op;
48  std::shared_ptr<const Restrictor<D>> restrictor;
52  std::shared_ptr<const Interpolator<D>> interpolator;
56  std::shared_ptr<const Smoother<D>> smoother;
60  std::shared_ptr<const Level> coarser;
61 
62 public:
66  Level() {}
72  void setRestrictor(const Restrictor<D>& restrictor)
73  {
74  this->restrictor.reset(restrictor.clone());
75  }
82  {
83  if (restrictor == nullptr) {
84  throw RuntimeError("This level does not have a restrictor");
85  }
86  return *restrictor;
87  }
93  void setInterpolator(const Interpolator<D>& interpolator)
94  {
95  this->interpolator.reset(interpolator.clone());
96  }
103  {
104  if (interpolator == nullptr) {
105  throw RuntimeError("This level does not have an interpolator");
106  }
107  return *interpolator;
108  }
114  void setOperator(const Operator<D>& op) { this->op.reset(op.clone()); }
120  const Operator<D>& getOperator() const
121  {
122  if (op == nullptr) {
123  throw RuntimeError("This level does not have an Operator");
124  }
125  return *op;
126  }
132  void setSmoother(const Smoother<D>& smoother) { this->smoother.reset(smoother.clone()); }
138  const Smoother<D>& getSmoother() const
139  {
140  if (smoother == nullptr) {
141  throw RuntimeError("This level does not have a smoother");
142  }
143  return *smoother;
144  }
150  void setCoarser(std::shared_ptr<const Level> coarser) { this->coarser = coarser; }
156  const Level& getCoarser() const
157  {
158  if (coarser == nullptr) {
159  throw RuntimeError("This level does not have a coarser level.");
160  }
161  return *coarser;
162  }
168  bool finest() const { return interpolator == nullptr; }
174  bool coarsest() const { return coarser == nullptr; }
175 };
176 extern template class Level<2>;
177 extern template class Level<3>;
178 } // namespace ThunderEgg::GMG
179 #endif
ThunderEgg::GMG::Interpolator::clone
virtual Interpolator< D > * clone() const =0
Clone this interpolator.
ThunderEgg::GMG::Restrictor::clone
virtual Restrictor< D > * clone() const =0
Clone this interpolator.
ThunderEgg::GMG::Level::getSmoother
const Smoother< D > & getSmoother() const
Get smoother operator for this level.
Definition: Level.h:138
ThunderEgg::GMG::Level::setSmoother
void setSmoother(const Smoother< D > &smoother)
Set the smoother for this level.
Definition: Level.h:132
Smoother.h
Smoother class.
ThunderEgg::Operator::clone
virtual Operator< D > * clone() const =0
Clone this operator.
ThunderEgg::GMG::Level::setRestrictor
void setRestrictor(const Restrictor< D > &restrictor)
Set the restriction operator for restricting from this level to the coarser level.
Definition: Level.h:72
Operator.h
Operator class.
ThunderEgg::GMG::Level::setOperator
void setOperator(const Operator< D > &op)
Set the operator (matrix) for this level.
Definition: Level.h:114
ThunderEgg::GMG::Level::getInterpolator
const Interpolator< D > & getInterpolator() const
Get the interpolation operator for this level.
Definition: Level.h:102
ThunderEgg::GMG::Restrictor
Abstract class for restriction operators.
Definition: Restrictor.h:36
ThunderEgg::GMG::Level::getOperator
const Operator< D > & getOperator() const
Get the operator for this level.
Definition: Level.h:120
Interpolator.h
Interpolator class.
ThunderEgg::GMG::Smoother::clone
virtual Smoother< D > * clone() const =0
Clone this smoother.
ThunderEgg::GMG::Level::getRestrictor
const Restrictor< D > & getRestrictor() const
Get the restriction operator for this level.
Definition: Level.h:81
ThunderEgg::GMG::Level::getCoarser
const Level & getCoarser() const
get reference to the coarser level.
Definition: Level.h:156
ThunderEgg::GMG::Interpolator
Abstract class for interpolation operators.
Definition: Interpolator.h:35
ThunderEgg::GMG::Level::finest
bool finest() const
Check if this level is the finest level.
Definition: Level.h:168
ThunderEgg::Operator
Base class for operators.
Definition: Operator.h:37
ThunderEgg::GMG::Level::setCoarser
void setCoarser(std::shared_ptr< const Level > coarser)
Set pointer to the coarser level.
Definition: Level.h:150
ThunderEgg::GMG::Level
Represents a level in geometric multi-grid.
Definition: Level.h:38
ThunderEgg::GMG::Level::Level
Level()
Create a Level object.
Definition: Level.h:66
ThunderEgg::GMG::Smoother
Abstract class for smoothing operators.
Definition: Smoother.h:34
ThunderEgg::GMG
Geometric-Multigrid classes.
Definition: Cycle.h:33
Restrictor.h
Restrictor class.
ThunderEgg::GMG::Level::coarsest
bool coarsest() const
Check if this level is the coarsest level.
Definition: Level.h:174
ThunderEgg::GMG::Level::setInterpolator
void setInterpolator(const Interpolator< D > &interpolator)
Set the interpolation operator for interpolating from this level to the finer level.
Definition: Level.h:93
ThunderEgg::RuntimeError
ThunderEgg runtime exception.
Definition: RuntimeError.h:36