Go to the documentation of this file.
21 #ifndef THUNDEREGG_GMG_CYCLEBUILDER_H
22 #define THUNDEREGG_GMG_CYCLEBUILDER_H
61 bool has_finest =
false;
65 bool has_coarsest =
false;
74 std::shared_ptr<Level<D>> finest_level;
78 std::shared_ptr<Level<D>> prev_level;
102 throw RuntimeError(
"addFinestLevel was already called");
106 finest_level = std::make_shared<Level<D>>();
107 finest_level->setOperator(op);
108 finest_level->setSmoother(smoother);
109 finest_level->setRestrictor(restrictor);
111 prev_level = finest_level;
128 throw RuntimeError(
"addFinestLevel has not been called yet");
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);
140 prev_level->setCoarser(new_level);
142 prev_level = new_level;
156 throw RuntimeError(
"addFinestLevel has not been called yet");
159 throw RuntimeError(
"addCoarsestLevel has already been called");
163 auto new_level = std::make_shared<Level<D>>();
164 new_level->setOperator(op);
165 new_level->setSmoother(smoother);
166 new_level->setInterpolator(interpolator);
168 prev_level->setCoarser(new_level);
179 throw RuntimeError(
"addCoarsestLevel has not been called");
181 std::shared_ptr<Cycle<D>> cycle;
183 cycle.reset(
new VCycle<D>(*finest_level, opts));
185 cycle.reset(
new WCycle<D>(*finest_level, opts));
194 extern template class CycleBuilder<2>;
195 extern template class CycleBuilder<3>;
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
Implementation of a V-cycle.
Definition: VCycle.h:35
Abstract class for restriction operators.
Definition: Restrictor.h:36
Builder for GMG cycles.
Definition: CycleBuilder.h:55
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
Implementation of a full multigrid cycle.
Definition: FMGCycle.h:35
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
std::string cycle_type
Cycle type.
Definition: CycleOpts.h:54
std::shared_ptr< Cycle< D > > getCycle() const
Get the completed Cycle object.
Definition: CycleBuilder.h:176
Abstract class for interpolation operators.
Definition: Interpolator.h:35
Implementation of a W-cycle.
Definition: WCycle.h:35
Base class for operators.
Definition: Operator.h:37
CycleBuilder(const CycleOpts &opts_in)
Construct a new CycleBuilder object.
Definition: CycleBuilder.h:86
Options for Cycle classes.
Definition: CycleOpts.h:33
Abstract class for smoothing operators.
Definition: Smoother.h:34
Geometric-Multigrid classes.
Definition: Cycle.h:33
ThunderEgg runtime exception.
Definition: RuntimeError.h:36