walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T > Class Template Reference

Detailed Description

template<typename Stencil_T, typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
class walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >

Class for multigrid V-cycle.

Template Parameters
Stencil_TThe stencil used for the discrete operator
OperatorCoarsening_TThe coarsening operator to use, defaults to direct coarsening
Restrict_TThe restriction operator to use
ProlongateAndCorrect_TThe prolongation and correction operator to use

#include <VCycles.h>

Public Types

typedef GhostLayerField< real_t, 1 > PdeField_T
 
typedef std::vector< real_tWeight_T
 
typedef GhostLayerField< real_t, Stencil_T::Size > StencilField_T
 

Public Member Functions

 VCycles (shared_ptr< StructuredBlockForest > blocks, const BlockDataID &uFieldId, const BlockDataID &fFieldId, const Weight_T weights, const uint_t iterations, const uint_t numLvl, const uint_t preSmoothingIters, const uint_t postSmoothingIters, const uint_t coarseIters, const std::function< real_t() > &residualNorm, const real_t residualNormThreshold=real_t(0), const uint_t residualCheckFrequency=uint_t(1), const Set< SUID > &requiredSelectors=Set< SUID >::emptySet(), const Set< SUID > &incompatibleSelectors=Set< SUID >::emptySet())
 Creates a multigrid V-cycle with a fixed stencil. More...
 
 VCycles (shared_ptr< StructuredBlockForest > blocks, const BlockDataID &uFieldId, const BlockDataID &fFieldId, const BlockDataID &stencilFieldId, const OperatorCoarsening_T &operatorCoarsening, const uint_t iterations, const uint_t numLvl, const uint_t preSmoothingIters, const uint_t postSmoothingIters, const uint_t coarseIters, const std::function< real_t() > &residualNorm, const real_t residualNormThreshold=real_t(0), const uint_t residualCheckFrequency=uint_t(1), const Set< SUID > &requiredSelectors=Set< SUID >::emptySet(), const Set< SUID > &incompatibleSelectors=Set< SUID >::emptySet())
 Creates a multigrid V-cycle with a stencil field. More...
 
void operator() ()
 
void VCycle ()
 
uint_t iterationsPerformed () const
 
bool thresholdReached () const
 
const std::vector< real_t > & convergenceRate ()
 

Static Public Member Functions

static Vector3< uint_tgetSizeForLevel (const uint_t level, const shared_ptr< StructuredBlockStorage > &blocks, IBlock *const block)
 

Private Attributes

StructuredBlockForestblocks_
 
std::vector< Weight_Tweights_
 
uint_t iterations_
 
uint_t numLvl_
 
uint_t preSmoothingIters_
 
uint_t postSmoothingIters_
 
uint_t coarseIters_
 
real_t residualNormThreshold_
 
uint_t residualCheckFrequency_
 
uint_t iterationsPerformed_
 
bool thresholdReached_
 
std::function< real_t() > residualNorm_
 
std::vector< real_tconvergenceRate_
 
std::vector< BlockDataIDuId_
 
std::vector< BlockDataIDfId_
 
std::vector< BlockDataIDrId_
 
std::vector< BlockDataIDstencilId_
 
BlockDataID dId_
 
BlockDataID zId_
 
std::vector< shared_ptr< pde::RBGSFixedStencil< Stencil_T > > > RBGSFixedSweeps_
 
std::vector< shared_ptr< pde::RBGS< Stencil_T > > > RBGSSweeps_
 
std::vector< std::function< void() > > RBGSIteration_
 
std::function< void() > CGIteration_
 
std::vector< std::function< void(IBlock *) > > computeResidual_
 
std::vector< std::function< void(IBlock *) > > restrict_
 
std::vector< std::function< void(IBlock *) > > zeroize_
 
std::vector< std::function< void(IBlock *) > > prolongateAndCorrect_
 
std::vector< blockforest::communication::UniformBufferedScheme< Stencil_T > > communication_
 
Set< SUID > requiredSelectors_
 
Set< SUID > incompatibleSelectors_
 

Member Typedef Documentation

◆ PdeField_T

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
typedef GhostLayerField< real_t, 1 > walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::PdeField_T

◆ StencilField_T

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
typedef GhostLayerField< real_t, Stencil_T::Size > walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::StencilField_T

◆ Weight_T

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
typedef std::vector< real_t > walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::Weight_T

Constructor & Destructor Documentation

◆ VCycles() [1/2]

template<typename Stencil_T , typename OperatorCoarsening_T , typename Restrict_T , typename ProlongateAndCorrect_T >
walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::VCycles ( shared_ptr< StructuredBlockForest blocks,
const BlockDataID uFieldId,
const BlockDataID fFieldId,
const Weight_T  weights,
const uint_t  iterations,
const uint_t  numLvl,
const uint_t  preSmoothingIters,
const uint_t  postSmoothingIters,
const uint_t  coarseIters,
const std::function< real_t() > &  residualNorm,
const real_t  residualNormThreshold = real_t(0),
const uint_t  residualCheckFrequency = uint_t(1),
const Set< SUID > &  requiredSelectors = Set<SUID>::emptySet(),
const Set< SUID > &  incompatibleSelectors = Set<SUID>::emptySet() 
)

Creates a multigrid V-cycle with a fixed stencil.

Parameters
blocksthe block storage where the fields are stored
uFieldIdthe block data id of the solution field on the finest level
fFieldIdthe block data id of the right-hand side field on the finest level
weightsvector of stencil weights for the discrete operator
iterationsmaximum number of V-cycles to perform
numLvlnumber of grid levels to use (including the finest level)
preSmoothingItersnumber of Gauss-Seidel iterations before restriction
postSmoothingItersnumber of Gauss-Seidel iterations after prolongation
coarseItersnumber of Conjugate Gradient iterations on coarsest grid
residualNormfunction that returns the norm of the current residuum
residualNormThresholdnorm threshold below which the iteration is terminated
residualCheckFrequencyhow often to check whether the threshold has been reached

◆ VCycles() [2/2]

template<typename Stencil_T , typename OperatorCoarsening_T , typename Restrict_T , typename ProlongateAndCorrect_T >
walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::VCycles ( shared_ptr< StructuredBlockForest blocks,
const BlockDataID uFieldId,
const BlockDataID fFieldId,
const BlockDataID stencilFieldId,
const OperatorCoarsening_T &  operatorCoarsening,
const uint_t  iterations,
const uint_t  numLvl,
const uint_t  preSmoothingIters,
const uint_t  postSmoothingIters,
const uint_t  coarseIters,
const std::function< real_t() > &  residualNorm,
const real_t  residualNormThreshold = real_t(0),
const uint_t  residualCheckFrequency = uint_t(1),
const Set< SUID > &  requiredSelectors = Set<SUID>::emptySet(),
const Set< SUID > &  incompatibleSelectors = Set<SUID>::emptySet() 
)

Creates a multigrid V-cycle with a stencil field.

Parameters
blocksthe block storage where the fields are stored
uFieldIdthe block data id of the solution field on the finest level
fFieldIdthe block data id of the right-hand side field on the finest level
stencilFieldIdthe block data id of the stencil field for the finest level. The values stored in the field must not change after this class has been constructed.
operatorCoarseningfunction that performs the stencil coarsening
iterationsmaximum number of V-cycles to perform
numLvlnumber of grid levels to use (including the finest level)
preSmoothingItersnumber of Gauss-Seidel iterations before restriction
postSmoothingItersnumber of Gauss-Seidel iterations after prolongation
coarseItersnumber of Conjugate Gradient iterations on coarsest grid
residualNormfunction that returns the norm of the current residuum
residualNormThresholdnorm threshold below which the iteration is terminated
residualCheckFrequencyhow often to check whether the threshold has been reached

Member Function Documentation

◆ convergenceRate()

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
const std::vector<real_t>& walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::convergenceRate ( )
inline

◆ getSizeForLevel()

template<typename Stencil_T , typename OperatorCoarsening_T , typename Restrict_T , typename ProlongateAndCorrect_T >
Vector3< uint_t > walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::getSizeForLevel ( const uint_t  level,
const shared_ptr< StructuredBlockStorage > &  blocks,
IBlock *const  block 
)
static

◆ iterationsPerformed()

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
uint_t walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::iterationsPerformed ( ) const
inline

◆ operator()()

template<typename Stencil_T , typename OperatorCoarsening_T , typename Restrict_T , typename ProlongateAndCorrect_T >
void walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::operator()

◆ thresholdReached()

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
bool walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::thresholdReached ( ) const
inline

◆ VCycle()

template<typename Stencil_T , typename OperatorCoarsening_T , typename Restrict_T , typename ProlongateAndCorrect_T >
void walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::VCycle

Member Data Documentation

◆ blocks_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
StructuredBlockForest& walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::blocks_
private

◆ CGIteration_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
std::function< void() > walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::CGIteration_
private

◆ coarseIters_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
uint_t walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::coarseIters_
private

◆ communication_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
std::vector< blockforest::communication::UniformBufferedScheme< Stencil_T > > walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::communication_
private

◆ computeResidual_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
std::vector<std::function< void(IBlock *) > > walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::computeResidual_
private

◆ convergenceRate_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
std::vector<real_t> walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::convergenceRate_
private

◆ dId_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
BlockDataID walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::dId_
private

◆ fId_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
std::vector<BlockDataID> walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::fId_
private

◆ incompatibleSelectors_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
Set< SUID > walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::incompatibleSelectors_
private

◆ iterations_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
uint_t walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::iterations_
private

◆ iterationsPerformed_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
uint_t walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::iterationsPerformed_
private

◆ numLvl_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
uint_t walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::numLvl_
private

◆ postSmoothingIters_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
uint_t walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::postSmoothingIters_
private

◆ preSmoothingIters_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
uint_t walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::preSmoothingIters_
private

◆ prolongateAndCorrect_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
std::vector<std::function< void(IBlock *) > > walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::prolongateAndCorrect_
private

◆ RBGSFixedSweeps_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
std::vector< shared_ptr<pde::RBGSFixedStencil< Stencil_T > > > walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::RBGSFixedSweeps_
private

◆ RBGSIteration_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
std::vector< std::function< void() > > walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::RBGSIteration_
private

◆ RBGSSweeps_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
std::vector< shared_ptr<pde::RBGS< Stencil_T > > > walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::RBGSSweeps_
private

◆ requiredSelectors_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
Set< SUID > walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::requiredSelectors_
private

◆ residualCheckFrequency_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
uint_t walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::residualCheckFrequency_
private

◆ residualNorm_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
std::function< real_t() > walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::residualNorm_
private

◆ residualNormThreshold_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
real_t walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::residualNormThreshold_
private

◆ restrict_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
std::vector<std::function< void(IBlock *) > > walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::restrict_
private

◆ rId_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
std::vector<BlockDataID> walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::rId_
private

◆ stencilId_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
std::vector<BlockDataID> walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::stencilId_
private

◆ thresholdReached_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
bool walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::thresholdReached_
private

◆ uId_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
std::vector<BlockDataID> walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::uId_
private

◆ weights_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
std::vector< Weight_T > walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::weights_
private

◆ zeroize_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
std::vector<std::function< void(IBlock *) > > walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::zeroize_
private

◆ zId_

template<typename Stencil_T , typename OperatorCoarsening_T = CoarsenStencilFieldsDCA<Stencil_T>, typename Restrict_T = Restrict<Stencil_T>, typename ProlongateAndCorrect_T = ProlongateAndCorrect<Stencil_T>>
BlockDataID walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::zId_
private

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