From e7b8b7d6ed30d2651e50bf7b099d49a2fa1014f6 Mon Sep 17 00:00:00 2001
From: mattijs <mattijs>
Date: Mon, 17 Jun 2019 10:15:24 +0100
Subject: [PATCH] ENH: GAMGSolver: specify coarsest level solver. Fixes #1342.

---
 .../lduMatrix/solvers/GAMG/GAMGSolver.C       | 60 +++++++++++++++++--
 .../lduMatrix/solvers/GAMG/GAMGSolver.H       |  3 +
 .../lduMatrix/solvers/GAMG/GAMGSolverSolve.C  | 40 ++-----------
 .../LES/motorBike/lesFiles/fvSolution         |  9 +++
 4 files changed, 74 insertions(+), 38 deletions(-)

diff --git a/src/OpenFOAM/matrices/lduMatrix/solvers/GAMG/GAMGSolver.C b/src/OpenFOAM/matrices/lduMatrix/solvers/GAMG/GAMGSolver.C
index 08414687e01..5a1edb952f2 100644
--- a/src/OpenFOAM/matrices/lduMatrix/solvers/GAMG/GAMGSolver.C
+++ b/src/OpenFOAM/matrices/lduMatrix/solvers/GAMG/GAMGSolver.C
@@ -27,6 +27,8 @@ License
 
 #include "GAMGSolver.H"
 #include "GAMGInterface.H"
+#include "PCG.H"
+#include "PBiCGStab.H"
 
 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
 
@@ -249,11 +251,11 @@ Foam::GAMGSolver::GAMGSolver
 
     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
                 (
@@ -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
diff --git a/src/OpenFOAM/matrices/lduMatrix/solvers/GAMG/GAMGSolver.H b/src/OpenFOAM/matrices/lduMatrix/solvers/GAMG/GAMGSolver.H
index 3736926650a..19ca1d0f4d8 100644
--- a/src/OpenFOAM/matrices/lduMatrix/solvers/GAMG/GAMGSolver.H
+++ b/src/OpenFOAM/matrices/lduMatrix/solvers/GAMG/GAMGSolver.H
@@ -135,6 +135,9 @@ class GAMGSolver
         //- LU decomposed coarsest matrix
         autoPtr<LUscalarMatrix> coarsestLUMatrixPtr_;
 
+        //- Sparse coarsest matrix solver
+        autoPtr<lduMatrix::solver> coarsestSolverPtr_;
+
 
     // Private Member Functions
 
diff --git a/src/OpenFOAM/matrices/lduMatrix/solvers/GAMG/GAMGSolverSolve.C b/src/OpenFOAM/matrices/lduMatrix/solvers/GAMG/GAMGSolverSolve.C
index 25fc7a84814..644de4eef9a 100644
--- a/src/OpenFOAM/matrices/lduMatrix/solvers/GAMG/GAMGSolverSolve.C
+++ b/src/OpenFOAM/matrices/lduMatrix/solvers/GAMG/GAMGSolverSolve.C
@@ -26,8 +26,6 @@ License
 \*---------------------------------------------------------------------------*/
 
 #include "GAMGSolver.H"
-#include "PCG.H"
-#include "PBiCGStab.H"
 #include "SubField.H"
 #include "PrecisionAdaptor.H"
 
@@ -695,42 +693,16 @@ void Foam::GAMGSolver::solveCoarsestLevel
     else
     {
         coarsestCorrField = 0;
-        solverPerformance coarseSolverPerf;
-
-        if (matrixLevels_[coarsestLevel].asymmetric())
-        {
-            coarseSolverPerf = PBiCGStab
-            (
-                "coarsestLevelCorr",
-                matrixLevels_[coarsestLevel],
-                interfaceLevelsBouCoeffs_[coarsestLevel],
-                interfaceLevelsIntCoeffs_[coarsestLevel],
-                interfaceLevels_[coarsestLevel],
-                PBiCGStabSolverDict(tolerance_, relTol_)
-            ).scalarSolve
-            (
-                coarsestCorrField,
-                coarsestSource
-            );
-        }
-        else
-        {
-            coarseSolverPerf = PCG
-            (
-                "coarsestLevelCorr",
-                matrixLevels_[coarsestLevel],
-                interfaceLevelsBouCoeffs_[coarsestLevel],
-                interfaceLevelsIntCoeffs_[coarsestLevel],
-                interfaceLevels_[coarsestLevel],
-                PCGsolverDict(tolerance_, relTol_)
-            ).scalarSolve
+        const solverPerformance coarseSolverPerf
+        (
+            coarsestSolverPtr_->solve
             (
                 coarsestCorrField,
                 coarsestSource
-            );
-        }
+            )
+        );
 
-        if (debug >= 2)
+        if (debug)
         {
             coarseSolverPerf.print(Info.masterStream(coarseComm));
         }
diff --git a/tutorials/incompressible/pisoFoam/LES/motorBike/lesFiles/fvSolution b/tutorials/incompressible/pisoFoam/LES/motorBike/lesFiles/fvSolution
index c94598b2732..a3f28e6951b 100644
--- a/tutorials/incompressible/pisoFoam/LES/motorBike/lesFiles/fvSolution
+++ b/tutorials/incompressible/pisoFoam/LES/motorBike/lesFiles/fvSolution
@@ -33,6 +33,15 @@ solvers
     {
         $p;
         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)"
-- 
GitLab