ENH: GAMGSolver: specify coarsest level solver. Fixes #1342.

This commit is contained in:
mattijs
2019-06-17 10:15:24 +01:00
committed by Andrew Heather
parent ed10b19a2c
commit e7b8b7d6ed
4 changed files with 74 additions and 38 deletions

View File

@ -27,6 +27,8 @@ License
#include "GAMGSolver.H" #include "GAMGSolver.H"
#include "GAMGInterface.H" #include "GAMGInterface.H"
#include "PCG.H"
#include "PBiCGStab.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * // // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
@ -249,11 +251,11 @@ Foam::GAMGSolver::GAMGSolver
if (matrixLevels_.size()) if (matrixLevels_.size())
{ {
if (directSolveCoarsest_) const label coarsestLevel = matrixLevels_.size() - 1;
{
const label coarsestLevel = matrixLevels_.size() - 1;
if (matrixLevels_.set(coarsestLevel)) if (matrixLevels_.set(coarsestLevel))
{
if (directSolveCoarsest_)
{ {
coarsestLUMatrixPtr_.reset coarsestLUMatrixPtr_.reset
( (
@ -265,6 +267,56 @@ Foam::GAMGSolver::GAMGSolver
) )
); );
} }
else
{
entry* coarseEntry = controlDict_.findEntry
(
"coarsestLevelCorr",
keyType::LITERAL_RECURSIVE
);
if (coarseEntry && coarseEntry->isDict())
{
coarsestSolverPtr_ = lduMatrix::solver::New
(
"coarsestLevelCorr",
matrixLevels_[coarsestLevel],
interfaceLevelsBouCoeffs_[coarsestLevel],
interfaceLevelsIntCoeffs_[coarsestLevel],
interfaceLevels_[coarsestLevel],
coarseEntry->dict()
);
}
else if (matrixLevels_[coarsestLevel].asymmetric())
{
coarsestSolverPtr_.set
(
new PBiCGStab
(
"coarsestLevelCorr",
matrixLevels_[coarsestLevel],
interfaceLevelsBouCoeffs_[coarsestLevel],
interfaceLevelsIntCoeffs_[coarsestLevel],
interfaceLevels_[coarsestLevel],
PBiCGStabSolverDict(tolerance_, relTol_)
)
);
}
else
{
coarsestSolverPtr_.set
(
new PCG
(
"coarsestLevelCorr",
matrixLevels_[coarsestLevel],
interfaceLevelsBouCoeffs_[coarsestLevel],
interfaceLevelsIntCoeffs_[coarsestLevel],
interfaceLevels_[coarsestLevel],
PCGsolverDict(tolerance_, relTol_)
)
);
}
}
} }
} }
else else

View File

@ -135,6 +135,9 @@ class GAMGSolver
//- LU decomposed coarsest matrix //- LU decomposed coarsest matrix
autoPtr<LUscalarMatrix> coarsestLUMatrixPtr_; autoPtr<LUscalarMatrix> coarsestLUMatrixPtr_;
//- Sparse coarsest matrix solver
autoPtr<lduMatrix::solver> coarsestSolverPtr_;
// Private Member Functions // Private Member Functions

View File

@ -26,8 +26,6 @@ License
\*---------------------------------------------------------------------------*/ \*---------------------------------------------------------------------------*/
#include "GAMGSolver.H" #include "GAMGSolver.H"
#include "PCG.H"
#include "PBiCGStab.H"
#include "SubField.H" #include "SubField.H"
#include "PrecisionAdaptor.H" #include "PrecisionAdaptor.H"
@ -695,42 +693,16 @@ void Foam::GAMGSolver::solveCoarsestLevel
else else
{ {
coarsestCorrField = 0; coarsestCorrField = 0;
solverPerformance coarseSolverPerf; const solverPerformance coarseSolverPerf
(
if (matrixLevels_[coarsestLevel].asymmetric()) coarsestSolverPtr_->solve
{
coarseSolverPerf = PBiCGStab
(
"coarsestLevelCorr",
matrixLevels_[coarsestLevel],
interfaceLevelsBouCoeffs_[coarsestLevel],
interfaceLevelsIntCoeffs_[coarsestLevel],
interfaceLevels_[coarsestLevel],
PBiCGStabSolverDict(tolerance_, relTol_)
).scalarSolve
( (
coarsestCorrField, coarsestCorrField,
coarsestSource coarsestSource
); )
} );
else
{
coarseSolverPerf = PCG
(
"coarsestLevelCorr",
matrixLevels_[coarsestLevel],
interfaceLevelsBouCoeffs_[coarsestLevel],
interfaceLevelsIntCoeffs_[coarsestLevel],
interfaceLevels_[coarsestLevel],
PCGsolverDict(tolerance_, relTol_)
).scalarSolve
(
coarsestCorrField,
coarsestSource
);
}
if (debug >= 2) if (debug)
{ {
coarseSolverPerf.print(Info.masterStream(coarseComm)); coarseSolverPerf.print(Info.masterStream(coarseComm));
} }

View File

@ -33,6 +33,15 @@ solvers
{ {
$p; $p;
relTol 0; relTol 0;
// Explicit specify solver for coarse-level correction to override
// solution tolerance
coarsestLevelCorr
{
solver PCG;
preconditioner DIC;
relTol 0.05;
}
} }
"(U|k|B|nuTilda)" "(U|k|B|nuTilda)"