ThunderEgg  1.0.0
ThunderEgg::GMG::CycleBuilder< D > Class Template Reference

Builder for GMG cycles. More...

#include <CycleBuilder.h>

Public Member Functions

 CycleBuilder (const CycleOpts &opts_in)
 Construct a new CycleBuilder object. More...
 
void addFinestLevel (const Operator< D > &op, const Smoother< D > &smoother, const Restrictor< D > &restrictor)
 Add the finest level to the Cycle. More...
 
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. More...
 
void addCoarsestLevel (const Operator< D > &op, Smoother< D > &smoother, const Interpolator< D > &interpolator)
 Add the next intermediate level to the Cycle. More...
 
std::shared_ptr< Cycle< D > > getCycle () const
 Get the completed Cycle object. More...
 

Detailed Description

template<int D>
class ThunderEgg::GMG::CycleBuilder< D >

Builder for GMG cycles.

The user provides an Operator, Smoother, Restrictor, and Interpolator object for each level.

The user is expected to make these calls in the following order:

    builder.addFinestLevel()
    for each intermediate level (if any)
    {
        builder.addIntermediateLevel()
    }
    builder.addCoarsestLevel()
    M = builder.getCycle();

If these are called in the wrong order, an exception will be thrown.

Template Parameters
Dthe number of cartesian dimensions

Constructor & Destructor Documentation

◆ CycleBuilder()

template<int D>
ThunderEgg::GMG::CycleBuilder< D >::CycleBuilder ( const CycleOpts opts_in)
inlineexplicit

Construct a new CycleBuilder object.

Parameters
opts_inthe options to use for the GMG cycle

Member Function Documentation

◆ addCoarsestLevel()

template<int D>
void ThunderEgg::GMG::CycleBuilder< D >::addCoarsestLevel ( const Operator< D > &  op,
Smoother< D > &  smoother,
const Interpolator< D > &  interpolator 
)
inline

Add the next intermediate level to the Cycle.

Parameters
opthe Operator for the level
smootherthe Smoother for the level
interpolatorthe Interpolator that restricts from this level to the finer level

◆ addFinestLevel()

template<int D>
void ThunderEgg::GMG::CycleBuilder< D >::addFinestLevel ( const Operator< D > &  op,
const Smoother< D > &  smoother,
const Restrictor< D > &  restrictor 
)
inline

Add the finest level to the Cycle.

Parameters
opthe Operator for the level
smootherthe Smoother for the level
restrictorthe Restrictor that restricts from this level to the coarser level
vgthe VectorGenerator for the level

◆ addIntermediateLevel()

template<int D>
void ThunderEgg::GMG::CycleBuilder< D >::addIntermediateLevel ( const Operator< D > &  op,
const Smoother< D > &  smoother,
const Restrictor< D > &  restrictor,
const Interpolator< D > &  interpolator 
)
inline

Add the next intermediate level to the Cycle.

Parameters
opthe Operator for the level
smootherthe Smoother for the level
restrictorthe Restrictor that restricts from this level to the coarser level
interpolatorthe Interpolator that restricts from this level to the finer level
vgthe VectorGenerator for the level

◆ getCycle()

template<int D>
std::shared_ptr<Cycle<D> > ThunderEgg::GMG::CycleBuilder< D >::getCycle ( ) const
inline

Get the completed Cycle object.

Returns
std::shared_ptr<Cycle<D>> the completed Cycle, will throw an exception if it is incomplete

The documentation for this class was generated from the following file: