Go to the documentation of this file.
21 #ifndef THUNDEREGG_GMG_WCYCLE_H
22 #define THUNDEREGG_GMG_WCYCLE_H
38 int num_pre_sweeps = 1;
39 int num_post_sweeps = 1;
40 int num_coarse_sweeps = 1;
41 int num_mid_sweeps = 1;
53 for (
int i = 0; i < num_coarse_sweeps; i++) {
57 for (
int i = 0; i < num_pre_sweeps; i++) {
66 this->
visit(coarser_level, coarser_f, coarser_u);
70 for (
int i = 0; i < num_mid_sweeps; i++) {
74 coarser_f = this->
restrict(level, f, u);
76 this->
visit(coarser_level, coarser_f, coarser_u);
80 for (
int i = 0; i < num_post_sweeps; i++) {
93 :
Cycle<D>(finest_level)
107 extern template class WCycle<2>;
108 extern template class WCycle<3>;
const Smoother< D > & getSmoother() const
Get smoother operator for this level.
Definition: Level.h:138
const Interpolator< D > & getInterpolator() const
Get the interpolation operator for this level.
Definition: Level.h:102
WCycle< D > * clone() const override
Get a clone of this WCycle.
Definition: WCycle.h:105
int mid_sweeps
Number of sweeps inbetween up and down.
Definition: CycleOpts.h:46
void visit(const Level< D > &level, const Vector< D > &f, Vector< D > &u) const override
Implements W-cycle. Pre-smooth, visit coarser level, smooth, visit coarse level, and then post-smooth...
Definition: WCycle.h:50
int post_sweeps
Number of sweeps on up cycle.
Definition: CycleOpts.h:42
const Level & getCoarser() const
get reference to the coarser level.
Definition: Level.h:156
Base abstract class for cycles.
Definition: Cycle.h:40
int pre_sweeps
Number of sweeps on down cycle.
Definition: CycleOpts.h:38
Implementation of a W-cycle.
Definition: WCycle.h:35
Represents a level in geometric multi-grid.
Definition: Level.h:38
Vector class for use in thunderegg.
Definition: Vector.h:42
Options for Cycle classes.
Definition: CycleOpts.h:33
Vector< D > getZeroClone() const
Get a vector of the same length initialized to zero.
Definition: Vector.h:601
Vector< D > restrict(const Level< D > &level, const Vector< D > &f, const Vector< D > &u) const
Prepare vectors for coarser level.
Definition: Cycle.h:57
Geometric-Multigrid classes.
Definition: Cycle.h:33
WCycle(const Level< D > &finest_level, const CycleOpts &opts)
Create new W-cycle.
Definition: WCycle.h:92
int coarse_sweeps
Number of sweeps on coarse level.
Definition: CycleOpts.h:50
bool coarsest() const
Check if this level is the coarsest level.
Definition: Level.h:174