ThunderEgg  1.0.0
CycleBuilder.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) 2020-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_CYCLEBUILDER_H
22 #define THUNDEREGG_GMG_CYCLEBUILDER_H
23 
30 #include <ThunderEgg/GMG/Level.h>
31 #include <ThunderEgg/GMG/VCycle.h>
32 #include <ThunderEgg/GMG/WCycle.h>
34 namespace ThunderEgg::GMG {
54 template<int D>
56 {
57 private:
61  bool has_finest = false;
65  bool has_coarsest = false;
69  CycleOpts opts;
70 
74  std::shared_ptr<Level<D>> finest_level;
78  std::shared_ptr<Level<D>> prev_level;
79 
80 public:
86  explicit CycleBuilder(const CycleOpts& opts_in)
87  : opts(opts_in)
88  {}
97  void addFinestLevel(const Operator<D>& op,
98  const Smoother<D>& smoother,
99  const Restrictor<D>& restrictor)
100  {
101  if (has_finest) {
102  throw RuntimeError("addFinestLevel was already called");
103  }
104  has_finest = true;
105 
106  finest_level = std::make_shared<Level<D>>();
107  finest_level->setOperator(op);
108  finest_level->setSmoother(smoother);
109  finest_level->setRestrictor(restrictor);
110 
111  prev_level = finest_level;
112  }
123  const Smoother<D>& smoother,
124  const Restrictor<D>& restrictor,
125  const Interpolator<D>& interpolator)
126  {
127  if (!has_finest) {
128  throw RuntimeError("addFinestLevel has not been called yet");
129  }
130  if (has_coarsest) {
131  throw RuntimeError("addCoarsestLevel has been called");
132  }
133 
134  auto new_level = std::make_shared<Level<D>>();
135  new_level->setOperator(op);
136  new_level->setSmoother(smoother);
137  new_level->setInterpolator(interpolator);
138  new_level->setRestrictor(restrictor);
139 
140  prev_level->setCoarser(new_level);
141 
142  prev_level = new_level;
143  }
152  Smoother<D>& smoother,
153  const Interpolator<D>& interpolator)
154  {
155  if (!has_finest) {
156  throw RuntimeError("addFinestLevel has not been called yet");
157  }
158  if (has_coarsest) {
159  throw RuntimeError("addCoarsestLevel has already been called");
160  }
161  has_coarsest = true;
162 
163  auto new_level = std::make_shared<Level<D>>();
164  new_level->setOperator(op);
165  new_level->setSmoother(smoother);
166  new_level->setInterpolator(interpolator);
167 
168  prev_level->setCoarser(new_level);
169  }
176  std::shared_ptr<Cycle<D>> getCycle() const
177  {
178  if (!has_coarsest) {
179  throw RuntimeError("addCoarsestLevel has not been called");
180  }
181  std::shared_ptr<Cycle<D>> cycle;
182  if (opts.cycle_type == "V") {
183  cycle.reset(new VCycle<D>(*finest_level, opts));
184  } else if (opts.cycle_type == "W") {
185  cycle.reset(new WCycle<D>(*finest_level, opts));
186  } else if (opts.cycle_type == "F") {
187  cycle.reset(new FMGCycle<D>(*finest_level, opts));
188  } else {
189  throw RuntimeError("Unsupported Cycle type: " + opts.cycle_type);
190  }
191  return cycle;
192  }
193 };
194 extern template class CycleBuilder<2>;
195 extern template class CycleBuilder<3>;
196 } // namespace ThunderEgg::GMG
197 
198 #endif
ThunderEgg::GMG::CycleBuilder::addIntermediateLevel
void addIntermediateLevel(const Operator< D > &op, const Smoother< D > &smoother, const Restrictor< D > &restrictor, const Interpolator< D > &interpolator)
Add the next intermediate level to the Cycle.
Definition: CycleBuilder.h:122
RuntimeError.h
RuntimeError struct.
FMGCycle.h
FMGCycle class.
ThunderEgg::GMG::VCycle
Implementation of a V-cycle.
Definition: VCycle.h:35
ThunderEgg::GMG::Restrictor
Abstract class for restriction operators.
Definition: Restrictor.h:36
ThunderEgg::GMG::CycleBuilder
Builder for GMG cycles.
Definition: CycleBuilder.h:55
Level.h
Level class.
ThunderEgg::GMG::CycleBuilder::addFinestLevel
void addFinestLevel(const Operator< D > &op, const Smoother< D > &smoother, const Restrictor< D > &restrictor)
Add the finest level to the Cycle.
Definition: CycleBuilder.h:97
WCycle.h
WCycle class.
VCycle.h
VCycle class.
ThunderEgg::GMG::FMGCycle
Implementation of a full multigrid cycle.
Definition: FMGCycle.h:35
ThunderEgg::GMG::CycleBuilder::addCoarsestLevel
void addCoarsestLevel(const Operator< D > &op, Smoother< D > &smoother, const Interpolator< D > &interpolator)
Add the next intermediate level to the Cycle.
Definition: CycleBuilder.h:151
ThunderEgg::GMG::CycleOpts::cycle_type
std::string cycle_type
Cycle type.
Definition: CycleOpts.h:54
ThunderEgg::GMG::CycleBuilder::getCycle
std::shared_ptr< Cycle< D > > getCycle() const
Get the completed Cycle object.
Definition: CycleBuilder.h:176
ThunderEgg::GMG::Interpolator
Abstract class for interpolation operators.
Definition: Interpolator.h:35
ThunderEgg::GMG::WCycle
Implementation of a W-cycle.
Definition: WCycle.h:35
ThunderEgg::Operator
Base class for operators.
Definition: Operator.h:37
ThunderEgg::GMG::CycleBuilder::CycleBuilder
CycleBuilder(const CycleOpts &opts_in)
Construct a new CycleBuilder object.
Definition: CycleBuilder.h:86
ThunderEgg::GMG::CycleOpts
Options for Cycle classes.
Definition: CycleOpts.h:33
ThunderEgg::GMG::Smoother
Abstract class for smoothing operators.
Definition: Smoother.h:34
ThunderEgg::GMG
Geometric-Multigrid classes.
Definition: Cycle.h:33
ThunderEgg::RuntimeError
ThunderEgg runtime exception.
Definition: RuntimeError.h:36