Robustness improvement on top AMG level
This commit is contained in:
parent
c0003d35e6
commit
0ef0eb1e6a
2 changed files with 19 additions and 4 deletions
|
@ -525,6 +525,9 @@ Foam::clusterAmgPolicy::clusterAmgPolicy
|
||||||
"clusterAmgPolicy::clusterAmgPolicy\n"
|
"clusterAmgPolicy::clusterAmgPolicy\n"
|
||||||
"(\n"
|
"(\n"
|
||||||
" const lduMatrix& matrix,\n"
|
" const lduMatrix& matrix,\n"
|
||||||
|
" const FieldField<Field, scalar>& bouCoeffs,\n"
|
||||||
|
" const FieldField<Field, scalar>& intCoeffs,\n"
|
||||||
|
" const lduInterfaceFieldPtrsList& interfaceFields,\n"
|
||||||
" const label groupSize,\n"
|
" const label groupSize,\n"
|
||||||
" const label minCoarseEqns\n"
|
" const label minCoarseEqns\n"
|
||||||
")"
|
")"
|
||||||
|
|
|
@ -193,11 +193,9 @@ void Foam::coarseAmgLevel::solve
|
||||||
{
|
{
|
||||||
lduSolverPerformance coarseSolverPerf;
|
lduSolverPerformance coarseSolverPerf;
|
||||||
|
|
||||||
label maxIter = Foam::min(2*policyPtr_->minCoarseEqns(), 1000);
|
|
||||||
|
|
||||||
dictionary topLevelDict;
|
dictionary topLevelDict;
|
||||||
topLevelDict.add("minIter", 1);
|
topLevelDict.add("minIter", 1);
|
||||||
topLevelDict.add("maxIter", maxIter);
|
topLevelDict.add("maxIter", 500);
|
||||||
topLevelDict.add("tolerance", tolerance);
|
topLevelDict.add("tolerance", tolerance);
|
||||||
topLevelDict.add("relTol", relTol);
|
topLevelDict.add("relTol", relTol);
|
||||||
|
|
||||||
|
@ -220,7 +218,6 @@ void Foam::coarseAmgLevel::solve
|
||||||
if (matrixPtr_->matrix().symmetric())
|
if (matrixPtr_->matrix().symmetric())
|
||||||
{
|
{
|
||||||
// Note: must change preconditioner to C0. HJ. 10/Oct/2017
|
// Note: must change preconditioner to C0. HJ. 10/Oct/2017
|
||||||
// topLevelDict.add("preconditioner", "Cholesky");
|
|
||||||
topLevelDict.add("preconditioner", "ILUC0");
|
topLevelDict.add("preconditioner", "ILUC0");
|
||||||
|
|
||||||
coarseSolverPerf = cgSolver
|
coarseSolverPerf = cgSolver
|
||||||
|
@ -248,6 +245,21 @@ void Foam::coarseAmgLevel::solve
|
||||||
).solve(x, b, cmpt);
|
).solve(x, b, cmpt);
|
||||||
}
|
}
|
||||||
|
|
||||||
|
// Check for convergence
|
||||||
|
const scalar magInitialRes = mag(coarseSolverPerf.initialResidual());
|
||||||
|
const scalar magFinalRes = mag(coarseSolverPerf.finalResidual());
|
||||||
|
|
||||||
|
if (magFinalRes > magInitialRes && magInitialRes > SMALL)
|
||||||
|
{
|
||||||
|
if (blockLduMatrix::debug)
|
||||||
|
{
|
||||||
|
Info<< "Divergence in top AMG level" << endl;
|
||||||
|
coarseSolverPerf.print();
|
||||||
|
}
|
||||||
|
|
||||||
|
x = 0;
|
||||||
|
}
|
||||||
|
|
||||||
// Restore debug
|
// Restore debug
|
||||||
blockLduMatrix::debug = oldDebug;
|
blockLduMatrix::debug = oldDebug;
|
||||||
|
|
||||||
|
|
Reference in a new issue