Class for multigrid V-cycle.
Stencil_T | The stencil used for the discrete operator |
OperatorCoarsening_T | The coarsening operator to use, defaults to direct coarsening |
Restrict_T | The restriction operator to use |
ProlongateAndCorrect_T | The prolongation and correction operator to use |
#include <VCycles.h>
Public Types | |
typedef GhostLayerField< real_t, 1 > | PdeField_T |
typedef std::vector< real_t > | Weight_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_t > | getSizeForLevel (const uint_t level, const shared_ptr< StructuredBlockStorage > &blocks, IBlock *const block) |
typedef GhostLayerField< real_t, 1 > walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::PdeField_T |
typedef GhostLayerField< real_t, Stencil_T::Size > walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::StencilField_T |
typedef std::vector< real_t > walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::Weight_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.
blocks | the block storage where the fields are stored |
uFieldId | the block data id of the solution field on the finest level |
fFieldId | the block data id of the right-hand side field on the finest level |
weights | vector of stencil weights for the discrete operator |
iterations | maximum number of V-cycles to perform |
numLvl | number of grid levels to use (including the finest level) |
preSmoothingIters | number of Gauss-Seidel iterations before restriction |
postSmoothingIters | number of Gauss-Seidel iterations after prolongation |
coarseIters | number of Conjugate Gradient iterations on coarsest grid |
residualNorm | function that returns the norm of the current residuum |
residualNormThreshold | norm threshold below which the iteration is terminated |
residualCheckFrequency | how often to check whether the threshold has been reached |
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.
blocks | the block storage where the fields are stored |
uFieldId | the block data id of the solution field on the finest level |
fFieldId | the block data id of the right-hand side field on the finest level |
stencilFieldId | the 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. |
operatorCoarsening | function that performs the stencil coarsening |
iterations | maximum number of V-cycles to perform |
numLvl | number of grid levels to use (including the finest level) |
preSmoothingIters | number of Gauss-Seidel iterations before restriction |
postSmoothingIters | number of Gauss-Seidel iterations after prolongation |
coarseIters | number of Conjugate Gradient iterations on coarsest grid |
residualNorm | function that returns the norm of the current residuum |
residualNormThreshold | norm threshold below which the iteration is terminated |
residualCheckFrequency | how often to check whether the threshold has been reached |
|
inline |
|
static |
|
inline |
void walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::operator() |
|
inline |
void walberla::pde::VCycles< Stencil_T, OperatorCoarsening_T, Restrict_T, ProlongateAndCorrect_T >::VCycle |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |