From 51732b04768831519fff4f9f5fe09d5f7381dd21 Mon Sep 17 00:00:00 2001
From: Kutalmis Bercin <kutalmis.bercin@esi-group.com>
Date: Mon, 19 Jul 2021 12:24:24 +0100
Subject: [PATCH] ENH: linear solvers: add variable-specific debug flags

Introduces a new optional keyword of label type 'log'
to linear-solver dictionaries to enable variable-specific
debug statements. For example, in fvOptions file:

    solvers
    {
        p
        {
            solver GAMG;
            ...
            log    2;
        }

        U
        {
            ...
            log    0;
        }
    }

The meanings of values of 'log' are:

    log    0;    <!--  no output
    log    1;    <!--  standard output
    log    2;    <!--  debug output
    // values higher than 2 are expected to have no effect

This keyword does not directly affect the operations of various
DebugSwitches and backward compatibility has been ensured in exchange
of code cleanness. The related DebugSwitches are:

    DebugSwitches
    {
        SolverPerformance    0;
        GAMG                 0;
        PCG                  0;
        PBiCG                0;
        smoothSolver         0;
    }
---
 .../matrices/LduMatrix/LduMatrix/LduMatrix.H  |  5 ++++-
 .../LduMatrix/LduMatrix/LduMatrixSolver.C     |  4 +++-
 .../LduMatrix/LduMatrix/SolverPerformance.C   |  6 ++++--
 .../LduMatrix/LduMatrix/SolverPerformance.H   |  4 +++-
 .../LduMatrix/Solvers/PBiCCCG/PBiCCCG.C       | 17 +++++++++++++---
 .../LduMatrix/Solvers/PBiCICG/PBiCICG.C       | 20 ++++++++++++++++---
 .../matrices/LduMatrix/Solvers/PCICG/PCICG.C  | 17 +++++++++++++---
 .../Solvers/SmoothSolver/SmoothSolver.C       | 17 +++++++++++++---
 .../matrices/lduMatrix/lduMatrix/lduMatrix.H  |  5 ++++-
 .../lduMatrix/lduMatrix/lduMatrixSolver.C     |  3 ++-
 .../lduMatrix/solvers/GAMG/GAMGSolver.C       |  6 +++---
 .../lduMatrix/solvers/GAMG/GAMGSolverSolve.C  | 18 ++++++++---------
 .../matrices/lduMatrix/solvers/PBiCG/PBiCG.C  |  8 ++++----
 .../lduMatrix/solvers/PBiCGStab/PBiCGStab.C   | 10 +++++-----
 .../matrices/lduMatrix/solvers/PCG/PCG.C      |  8 ++++----
 .../matrices/lduMatrix/solvers/PPCG/PPCG.C    |  4 ++--
 .../solvers/smoothSolver/smoothSolver.C       |  8 ++++----
 .../fvMatrices/fvMatrix/fvMatrixSolve.C       | 20 ++++++++++++++++---
 .../fvScalarMatrix/fvScalarMatrix.C           | 10 ++++++++--
 19 files changed, 135 insertions(+), 55 deletions(-)

diff --git a/src/OpenFOAM/matrices/LduMatrix/LduMatrix/LduMatrix.H b/src/OpenFOAM/matrices/LduMatrix/LduMatrix/LduMatrix.H
index bae2bd10d59..c0673eb34f7 100644
--- a/src/OpenFOAM/matrices/LduMatrix/LduMatrix/LduMatrix.H
+++ b/src/OpenFOAM/matrices/LduMatrix/LduMatrix/LduMatrix.H
@@ -6,7 +6,7 @@
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
     Copyright (C) 2011-2017 OpenFOAM Foundation
-    Copyright (C) 2019 OpenCFD Ltd.
+    Copyright (C) 2019-2021 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
@@ -134,6 +134,9 @@ public:
             //- Minimum number of iterations in the solver
             label minIter_;
 
+            //- Level of verbosity in the solver output statements
+            label log_;
+
             //- Final convergence tolerance
             Type tolerance_;
 
diff --git a/src/OpenFOAM/matrices/LduMatrix/LduMatrix/LduMatrixSolver.C b/src/OpenFOAM/matrices/LduMatrix/LduMatrix/LduMatrixSolver.C
index 26e0dab299c..54a5923603c 100644
--- a/src/OpenFOAM/matrices/LduMatrix/LduMatrix/LduMatrixSolver.C
+++ b/src/OpenFOAM/matrices/LduMatrix/LduMatrix/LduMatrixSolver.C
@@ -6,7 +6,7 @@
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
     Copyright (C) 2011-2017 OpenFOAM Foundation
-    Copyright (C) 2019 OpenCFD Ltd.
+    Copyright (C) 2019-2021 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
@@ -131,6 +131,7 @@ Foam::LduMatrix<Type, DType, LUType>::solver::solver
 
     maxIter_(defaultMaxIter_),
     minIter_(0),
+    log_(1),
     tolerance_(1e-6*pTraits<Type>::one),
     relTol_(Zero)
 {
@@ -145,6 +146,7 @@ void Foam::LduMatrix<Type, DType, LUType>::solver::readControls()
 {
     readControl(controlDict_, maxIter_, "maxIter");
     readControl(controlDict_, minIter_, "minIter");
+    readControl(controlDict_, log_, "log");
     readControl(controlDict_, tolerance_, "tolerance");
     readControl(controlDict_, relTol_, "relTol");
 }
diff --git a/src/OpenFOAM/matrices/LduMatrix/LduMatrix/SolverPerformance.C b/src/OpenFOAM/matrices/LduMatrix/LduMatrix/SolverPerformance.C
index e291b0289ab..2c53aedbee8 100644
--- a/src/OpenFOAM/matrices/LduMatrix/LduMatrix/SolverPerformance.C
+++ b/src/OpenFOAM/matrices/LduMatrix/LduMatrix/SolverPerformance.C
@@ -6,6 +6,7 @@
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
     Copyright (C) 2011-2016 OpenFOAM Foundation
+    Copyright (C) 2021 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
@@ -62,10 +63,11 @@ template<class Type>
 bool Foam::SolverPerformance<Type>::checkConvergence
 (
     const Type& Tolerance,
-    const Type& RelTolerance
+    const Type& RelTolerance,
+    const label log
 )
 {
-    if (debug >= 2)
+    if ((log >= 2) || (debug >= 2))
     {
         Info<< solverName_
             << ":  Iteration " << nIterations_
diff --git a/src/OpenFOAM/matrices/LduMatrix/LduMatrix/SolverPerformance.H b/src/OpenFOAM/matrices/LduMatrix/LduMatrix/SolverPerformance.H
index 9d835599ce5..476e842cb3c 100644
--- a/src/OpenFOAM/matrices/LduMatrix/LduMatrix/SolverPerformance.H
+++ b/src/OpenFOAM/matrices/LduMatrix/LduMatrix/SolverPerformance.H
@@ -6,6 +6,7 @@
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
     Copyright (C) 2012-2016 OpenFOAM Foundation
+    Copyright (C) 2021 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
@@ -218,7 +219,8 @@ public:
         bool checkConvergence
         (
             const Type& tolerance,
-            const Type& relTolerance
+            const Type& relTolerance,
+            const label log
         );
 
         //- Singularity test
diff --git a/src/OpenFOAM/matrices/LduMatrix/Solvers/PBiCCCG/PBiCCCG.C b/src/OpenFOAM/matrices/LduMatrix/Solvers/PBiCCCG/PBiCCCG.C
index f20a64a0260..c219c276a60 100644
--- a/src/OpenFOAM/matrices/LduMatrix/Solvers/PBiCCCG/PBiCCCG.C
+++ b/src/OpenFOAM/matrices/LduMatrix/Solvers/PBiCCCG/PBiCCCG.C
@@ -6,6 +6,7 @@
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
     Copyright (C) 2011-2016 OpenFOAM Foundation
+    Copyright (C) 2021 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
@@ -98,7 +99,7 @@ Foam::PBiCCCG<Type, DType, LUType>::solve
     // --- Calculate normalisation factor
     Type normFactor = this->normFactor(psi, wA, pA);
 
-    if (LduMatrix<Type, DType, LUType>::debug >= 2)
+    if ((this->log_ >= 2) || (LduMatrix<Type, DType, LUType>::debug >= 2))
     {
         Info<< "   Normalisation factor = " << normFactor << endl;
     }
@@ -111,7 +112,12 @@ Foam::PBiCCCG<Type, DType, LUType>::solve
     if
     (
         this->minIter_ > 0
-     || !solverPerf.checkConvergence(this->tolerance_, this->relTol_)
+     || !solverPerf.checkConvergence
+        (
+            this->tolerance_,
+            this->relTol_,
+            this->log_
+        )
     )
     {
         // --- Select and construct the preconditioner
@@ -192,7 +198,12 @@ Foam::PBiCCCG<Type, DType, LUType>::solve
         (
             (
                 nIter++ < this->maxIter_
-            && !solverPerf.checkConvergence(this->tolerance_, this->relTol_)
+            && !solverPerf.checkConvergence
+                (
+                    this->tolerance_,
+                    this->relTol_,
+                    this->log_
+                )
             )
          || nIter < this->minIter_
         );
diff --git a/src/OpenFOAM/matrices/LduMatrix/Solvers/PBiCICG/PBiCICG.C b/src/OpenFOAM/matrices/LduMatrix/Solvers/PBiCICG/PBiCICG.C
index 5796597fac9..6dee6ad57ee 100644
--- a/src/OpenFOAM/matrices/LduMatrix/Solvers/PBiCICG/PBiCICG.C
+++ b/src/OpenFOAM/matrices/LduMatrix/Solvers/PBiCICG/PBiCICG.C
@@ -6,6 +6,7 @@
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
     Copyright (C) 2011-2016 OpenFOAM Foundation
+    Copyright (C) 2021 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
@@ -95,7 +96,7 @@ Foam::PBiCICG<Type, DType, LUType>::solve(Field<Type>& psi) const
     // --- Calculate normalisation factor
     Type normFactor = this->normFactor(psi, wA, pA);
 
-    if (LduMatrix<Type, DType, LUType>::debug >= 2)
+    if ((this->log_ >= 2) || (LduMatrix<Type, DType, LUType>::debug >= 2))
     {
         Info<< "   Normalisation factor = " << normFactor << endl;
     }
@@ -105,7 +106,15 @@ Foam::PBiCICG<Type, DType, LUType>::solve(Field<Type>& psi) const
     solverPerf.finalResidual() = solverPerf.initialResidual();
 
     // --- Check convergence, solve if not converged
-    if (!solverPerf.checkConvergence(this->tolerance_, this->relTol_))
+    if
+    (
+        !solverPerf.checkConvergence
+        (
+            this->tolerance_,
+            this->relTol_,
+            this->log_
+        )
+    )
     {
         // --- Select and construct the preconditioner
         autoPtr<typename LduMatrix<Type, DType, LUType>::preconditioner>
@@ -192,7 +201,12 @@ Foam::PBiCICG<Type, DType, LUType>::solve(Field<Type>& psi) const
         } while
         (
             nIter++ < this->maxIter_
-        && !(solverPerf.checkConvergence(this->tolerance_, this->relTol_))
+        && !solverPerf.checkConvergence
+            (
+                this->tolerance_,
+                this->relTol_,
+                this->log_
+            )
         );
     }
 
diff --git a/src/OpenFOAM/matrices/LduMatrix/Solvers/PCICG/PCICG.C b/src/OpenFOAM/matrices/LduMatrix/Solvers/PCICG/PCICG.C
index 7adf70b392f..6dec5d953d3 100644
--- a/src/OpenFOAM/matrices/LduMatrix/Solvers/PCICG/PCICG.C
+++ b/src/OpenFOAM/matrices/LduMatrix/Solvers/PCICG/PCICG.C
@@ -6,6 +6,7 @@
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
     Copyright (C) 2011-2016 OpenFOAM Foundation
+    Copyright (C) 2021 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
@@ -86,7 +87,7 @@ Foam::PCICG<Type, DType, LUType>::solve(Field<Type>& psi) const
     // --- Calculate normalisation factor
     Type normFactor = this->normFactor(psi, wA, pA);
 
-    if (LduMatrix<Type, DType, LUType>::debug >= 2)
+    if ((this->log_ >= 2) || (LduMatrix<Type, DType, LUType>::debug >= 2))
     {
         Info<< "   Normalisation factor = " << normFactor << endl;
     }
@@ -99,7 +100,12 @@ Foam::PCICG<Type, DType, LUType>::solve(Field<Type>& psi) const
     if
     (
         this->minIter_ > 0
-     || !solverPerf.checkConvergence(this->tolerance_, this->relTol_)
+     || !solverPerf.checkConvergence
+        (
+            this->tolerance_,
+            this->relTol_,
+            this->log_
+        )
     )
     {
         // --- Select and construct the preconditioner
@@ -184,7 +190,12 @@ Foam::PCICG<Type, DType, LUType>::solve(Field<Type>& psi) const
         (
             (
                 nIter++ < this->maxIter_
-            && !solverPerf.checkConvergence(this->tolerance_, this->relTol_)
+            && !solverPerf.checkConvergence
+                (
+                    this->tolerance_,
+                    this->relTol_,
+                    this->log_
+                )
             )
          || nIter < this->minIter_
         );
diff --git a/src/OpenFOAM/matrices/LduMatrix/Solvers/SmoothSolver/SmoothSolver.C b/src/OpenFOAM/matrices/LduMatrix/Solvers/SmoothSolver/SmoothSolver.C
index 9c5a5c80fd7..ed10c2ae121 100644
--- a/src/OpenFOAM/matrices/LduMatrix/Solvers/SmoothSolver/SmoothSolver.C
+++ b/src/OpenFOAM/matrices/LduMatrix/Solvers/SmoothSolver/SmoothSolver.C
@@ -6,6 +6,7 @@
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
     Copyright (C) 2011-2016 OpenFOAM Foundation
+    Copyright (C) 2021 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
@@ -110,7 +111,7 @@ Foam::SmoothSolver<Type, DType, LUType>::solve(Field<Type>& psi) const
             solverPerf.finalResidual() = solverPerf.initialResidual();
         }
 
-        if (LduMatrix<Type, DType, LUType>::debug >= 2)
+        if ((this->log_ >= 2) || (LduMatrix<Type, DType, LUType>::debug >= 2))
         {
             Info<< "   Normalisation factor = " << normFactor << endl;
         }
@@ -120,7 +121,12 @@ Foam::SmoothSolver<Type, DType, LUType>::solve(Field<Type>& psi) const
         if
         (
             this->minIter_ > 0
-         || !solverPerf.checkConvergence(this->tolerance_, this->relTol_)
+         || !solverPerf.checkConvergence
+            (
+                this->tolerance_,
+                this->relTol_,
+                this->log_
+            )
         )
         {
             autoPtr<typename LduMatrix<Type, DType, LUType>::smoother>
@@ -150,7 +156,12 @@ Foam::SmoothSolver<Type, DType, LUType>::solve(Field<Type>& psi) const
             (
                 (
                     (nIter += nSweeps_) < this->maxIter_
-                && !solverPerf.checkConvergence(this->tolerance_, this->relTol_)
+                && !solverPerf.checkConvergence
+                    (
+                        this->tolerance_,
+                        this->relTol_,
+                        this->log_
+                    )
                 )
              || nIter < this->minIter_
             );
diff --git a/src/OpenFOAM/matrices/lduMatrix/lduMatrix/lduMatrix.H b/src/OpenFOAM/matrices/lduMatrix/lduMatrix/lduMatrix.H
index 644de9caa0d..f7642bad246 100644
--- a/src/OpenFOAM/matrices/lduMatrix/lduMatrix/lduMatrix.H
+++ b/src/OpenFOAM/matrices/lduMatrix/lduMatrix/lduMatrix.H
@@ -6,7 +6,7 @@
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
     Copyright (C) 2011-2017 OpenFOAM Foundation
-    Copyright (C) 2016-2019 OpenCFD Ltd.
+    Copyright (C) 2016-2021 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
@@ -120,6 +120,9 @@ public:
             //- Minimum number of iterations in the solver
             label minIter_;
 
+            //- Level of verbosity in the solver output statements
+            label log_;
+
             //- Final convergence tolerance
             scalar tolerance_;
 
diff --git a/src/OpenFOAM/matrices/lduMatrix/lduMatrix/lduMatrixSolver.C b/src/OpenFOAM/matrices/lduMatrix/lduMatrix/lduMatrixSolver.C
index b750b74a65b..8a3f5ad6a1e 100644
--- a/src/OpenFOAM/matrices/lduMatrix/lduMatrix/lduMatrixSolver.C
+++ b/src/OpenFOAM/matrices/lduMatrix/lduMatrix/lduMatrixSolver.C
@@ -6,7 +6,7 @@
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
     Copyright (C) 2011-2017 OpenFOAM Foundation
-    Copyright (C) 2016-2020 OpenCFD Ltd.
+    Copyright (C) 2016-2021 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
@@ -164,6 +164,7 @@ void Foam::lduMatrix::solver::readControls()
 {
     maxIter_ = controlDict_.getOrDefault<label>("maxIter", defaultMaxIter_);
     minIter_ = controlDict_.getOrDefault<label>("minIter", 0);
+    log_ = controlDict_.getOrDefault<label>("log", 1);
     tolerance_ = controlDict_.getOrDefault<scalar>("tolerance", 1e-6);
     relTol_ = controlDict_.getOrDefault<scalar>("relTol", 0);
 }
diff --git a/src/OpenFOAM/matrices/lduMatrix/solvers/GAMG/GAMGSolver.C b/src/OpenFOAM/matrices/lduMatrix/solvers/GAMG/GAMGSolver.C
index 946ae1d49f6..bd79b8c07fd 100644
--- a/src/OpenFOAM/matrices/lduMatrix/solvers/GAMG/GAMGSolver.C
+++ b/src/OpenFOAM/matrices/lduMatrix/solvers/GAMG/GAMGSolver.C
@@ -6,6 +6,7 @@
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
     Copyright (C) 2011-2017 OpenFOAM Foundation
+    Copyright (C) 2021 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
@@ -208,8 +209,7 @@ Foam::GAMGSolver::GAMGSolver
         }
     }
 
-
-    if (debug & 2)
+    if ((log_ >= 2) || (debug & 2))
     {
         for
         (
@@ -368,7 +368,7 @@ void Foam::GAMGSolver::readControls()
     controlDict_.readIfPresent("scaleCorrection", scaleCorrection_);
     controlDict_.readIfPresent("directSolveCoarsest", directSolveCoarsest_);
 
-    if (debug)
+    if ((log_ >= 2) || debug)
     {
         Info<< "GAMGSolver settings :"
             << " cacheAgglomeration:" << cacheAgglomeration_
diff --git a/src/OpenFOAM/matrices/lduMatrix/solvers/GAMG/GAMGSolverSolve.C b/src/OpenFOAM/matrices/lduMatrix/solvers/GAMG/GAMGSolverSolve.C
index a8d1dbf175e..58076216e33 100644
--- a/src/OpenFOAM/matrices/lduMatrix/solvers/GAMG/GAMGSolverSolve.C
+++ b/src/OpenFOAM/matrices/lduMatrix/solvers/GAMG/GAMGSolverSolve.C
@@ -6,7 +6,7 @@
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
     Copyright (C) 2011-2017 OpenFOAM Foundation
-    Copyright (C) 2016-2019 OpenCFD Ltd.
+    Copyright (C) 2016-2021 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
@@ -59,7 +59,7 @@ Foam::solverPerformance Foam::GAMGSolver::solve
     solveScalar normFactor =
         this->normFactor(psi, tsource(), Apsi, finestCorrection);
 
-    if (debug >= 2)
+    if ((log_ >= 2) || (debug >= 2))
     {
         Pout<< "   Normalisation factor = " << normFactor << endl;
     }
@@ -87,7 +87,7 @@ Foam::solverPerformance Foam::GAMGSolver::solve
     if
     (
         minIter_ > 0
-     || !solverPerf.checkConvergence(tolerance_, relTol_)
+     || !solverPerf.checkConvergence(tolerance_, relTol_, log_)
     )
     {
         // Create coarse grid correction fields
@@ -144,7 +144,7 @@ Foam::solverPerformance Foam::GAMGSolver::solve
                 matrix().mesh().comm()
             )/normFactor;
 
-            if (debug >= 2)
+            if ((log_ >= 2) || (debug >= 2))
             {
                 solverPerf.print(Info.masterStream(matrix().mesh().comm()));
             }
@@ -152,7 +152,7 @@ Foam::solverPerformance Foam::GAMGSolver::solve
         (
             (
               ++solverPerf.nIterations() < maxIter_
-            && !solverPerf.checkConvergence(tolerance_, relTol_)
+            && !solverPerf.checkConvergence(tolerance_, relTol_, log_)
             )
          || solverPerf.nIterations() < minIter_
         );
@@ -193,7 +193,7 @@ void Foam::GAMGSolver::Vcycle
     // Restrict finest grid residual for the next level up.
     agglomeration_.restrictField(coarseSources[0], finestResidual, 0, true);
 
-    if (debug >= 2 && nPreSweeps_)
+    if (nPreSweeps_ && ((log_ >= 2) || (debug >= 2)))
     {
         Pout<< "Pre-smoothing scaling factors: ";
     }
@@ -274,7 +274,7 @@ void Foam::GAMGSolver::Vcycle
         }
     }
 
-    if (debug >= 2 && nPreSweeps_)
+    if (nPreSweeps_ && ((log_ >= 2) || (debug >= 2)))
     {
         Pout<< endl;
     }
@@ -290,7 +290,7 @@ void Foam::GAMGSolver::Vcycle
         );
     }
 
-    if (debug >= 2)
+    if ((log_ >= 2) || (debug >= 2))
     {
         Pout<< "Post-smoothing scaling factors: ";
     }
@@ -703,7 +703,7 @@ void Foam::GAMGSolver::solveCoarsestLevel
             )
         );
 
-        if (debug)
+        if ((log_ >= 2) || debug)
         {
             coarseSolverPerf.print(Info.masterStream(coarseComm));
         }
diff --git a/src/OpenFOAM/matrices/lduMatrix/solvers/PBiCG/PBiCG.C b/src/OpenFOAM/matrices/lduMatrix/solvers/PBiCG/PBiCG.C
index cba20eba933..d033d292fc3 100644
--- a/src/OpenFOAM/matrices/lduMatrix/solvers/PBiCG/PBiCG.C
+++ b/src/OpenFOAM/matrices/lduMatrix/solvers/PBiCG/PBiCG.C
@@ -6,7 +6,7 @@
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
     Copyright (C) 2011-2017 OpenFOAM Foundation
-    Copyright (C) 2019 OpenCFD Ltd.
+    Copyright (C) 2019-2021 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
@@ -111,7 +111,7 @@ Foam::solverPerformance Foam::PBiCG::solve
     // --- Calculate normalisation factor
     const solveScalar normFactor = this->normFactor(psi, tsource(), wA, pA);
 
-    if (lduMatrix::debug >= 2)
+    if ((log_ >= 2) || (lduMatrix::debug >= 2))
     {
         Info<< "   Normalisation factor = " << normFactor << endl;
     }
@@ -126,7 +126,7 @@ Foam::solverPerformance Foam::PBiCG::solve
     if
     (
         minIter_ > 0
-     || !solverPerf.checkConvergence(tolerance_, relTol_)
+     || !solverPerf.checkConvergence(tolerance_, relTol_, log_)
     )
     {
         solveScalarField pT(nCells, 0);
@@ -217,7 +217,7 @@ Foam::solverPerformance Foam::PBiCG::solve
         (
             (
               ++solverPerf.nIterations() < maxIter_
-            && !solverPerf.checkConvergence(tolerance_, relTol_)
+            && !solverPerf.checkConvergence(tolerance_, relTol_, log_)
             )
          || solverPerf.nIterations() < minIter_
         );
diff --git a/src/OpenFOAM/matrices/lduMatrix/solvers/PBiCGStab/PBiCGStab.C b/src/OpenFOAM/matrices/lduMatrix/solvers/PBiCGStab/PBiCGStab.C
index 086fba5c07b..ccf535206fe 100644
--- a/src/OpenFOAM/matrices/lduMatrix/solvers/PBiCGStab/PBiCGStab.C
+++ b/src/OpenFOAM/matrices/lduMatrix/solvers/PBiCGStab/PBiCGStab.C
@@ -6,7 +6,7 @@
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
     Copyright (C) 2016-2017 OpenFOAM Foundation
-    Copyright (C) 2019 OpenCFD Ltd.
+    Copyright (C) 2019-2021 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
@@ -110,7 +110,7 @@ Foam::solverPerformance Foam::PBiCGStab::scalarSolve
     // --- Calculate normalisation factor
     const solveScalar normFactor = this->normFactor(psi, source, yA, pA);
 
-    if (lduMatrix::debug >= 2)
+    if ((log_ >= 2) || (lduMatrix::debug >= 2))
     {
         Info<< "   Normalisation factor = " << normFactor << endl;
     }
@@ -125,7 +125,7 @@ Foam::solverPerformance Foam::PBiCGStab::scalarSolve
     if
     (
         minIter_ > 0
-     || !solverPerf.checkConvergence(tolerance_, relTol_)
+     || !solverPerf.checkConvergence(tolerance_, relTol_, log_)
     )
     {
         solveScalarField AyA(nCells);
@@ -219,7 +219,7 @@ Foam::solverPerformance Foam::PBiCGStab::scalarSolve
             if
             (
                 solverPerf.nIterations() >= minIter_
-             && solverPerf.checkConvergence(tolerance_, relTol_)
+             && solverPerf.checkConvergence(tolerance_, relTol_, log_)
             )
             {
                 for (label cell=0; cell<nCells; cell++)
@@ -258,7 +258,7 @@ Foam::solverPerformance Foam::PBiCGStab::scalarSolve
         (
             (
               ++solverPerf.nIterations() < maxIter_
-            && !solverPerf.checkConvergence(tolerance_, relTol_)
+            && !solverPerf.checkConvergence(tolerance_, relTol_, log_)
             )
          || solverPerf.nIterations() < minIter_
         );
diff --git a/src/OpenFOAM/matrices/lduMatrix/solvers/PCG/PCG.C b/src/OpenFOAM/matrices/lduMatrix/solvers/PCG/PCG.C
index 2b34a6aad37..ce33493e17c 100644
--- a/src/OpenFOAM/matrices/lduMatrix/solvers/PCG/PCG.C
+++ b/src/OpenFOAM/matrices/lduMatrix/solvers/PCG/PCG.C
@@ -6,7 +6,7 @@
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
     Copyright (C) 2011-2017 OpenFOAM Foundation
-    Copyright (C) 2019 OpenCFD Ltd.
+    Copyright (C) 2019-2021 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
@@ -110,7 +110,7 @@ Foam::solverPerformance Foam::PCG::scalarSolve
     // --- Calculate normalisation factor
     solveScalar normFactor = this->normFactor(psi, source, wA, pA);
 
-    if (lduMatrix::debug >= 2)
+    if ((log_ >= 2) || (lduMatrix::debug >= 2))
     {
         Info<< "   Normalisation factor = " << normFactor << endl;
     }
@@ -125,7 +125,7 @@ Foam::solverPerformance Foam::PCG::scalarSolve
     if
     (
         minIter_ > 0
-     || !solverPerf.checkConvergence(tolerance_, relTol_)
+     || !solverPerf.checkConvergence(tolerance_, relTol_, log_)
     )
     {
         // --- Select and construct the preconditioner
@@ -193,7 +193,7 @@ Foam::solverPerformance Foam::PCG::scalarSolve
         (
             (
               ++solverPerf.nIterations() < maxIter_
-            && !solverPerf.checkConvergence(tolerance_, relTol_)
+            && !solverPerf.checkConvergence(tolerance_, relTol_, log_)
             )
          || solverPerf.nIterations() < minIter_
         );
diff --git a/src/OpenFOAM/matrices/lduMatrix/solvers/PPCG/PPCG.C b/src/OpenFOAM/matrices/lduMatrix/solvers/PPCG/PPCG.C
index d231a2daac6..6ad660f76f1 100644
--- a/src/OpenFOAM/matrices/lduMatrix/solvers/PPCG/PPCG.C
+++ b/src/OpenFOAM/matrices/lduMatrix/solvers/PPCG/PPCG.C
@@ -107,7 +107,7 @@ Foam::solverPerformance Foam::PPCG::scalarSolveCG
     solveScalarField p(nCells);
     const solveScalar normFactor = this->normFactor(psi, source, w, p);
 
-    if (lduMatrix::debug >= 2)
+    if ((log_ >= 2) || (lduMatrix::debug >= 2))
     {
         Info<< "   Normalisation factor = " << normFactor << endl;
     }
@@ -190,7 +190,7 @@ Foam::solverPerformance Foam::PPCG::scalarSolveCG
         if
         (
             (minIter_ <= 0 || solverPerf.nIterations() >= minIter_)
-         && solverPerf.checkConvergence(tolerance_, relTol_)
+         && solverPerf.checkConvergence(tolerance_, relTol_, log_)
         )
         {
             break;
diff --git a/src/OpenFOAM/matrices/lduMatrix/solvers/smoothSolver/smoothSolver.C b/src/OpenFOAM/matrices/lduMatrix/solvers/smoothSolver/smoothSolver.C
index 78e6a7c62d8..1c25d08ccba 100644
--- a/src/OpenFOAM/matrices/lduMatrix/solvers/smoothSolver/smoothSolver.C
+++ b/src/OpenFOAM/matrices/lduMatrix/solvers/smoothSolver/smoothSolver.C
@@ -6,7 +6,7 @@
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
     Copyright (C) 2011-2014 OpenFOAM Foundation
-    Copyright (C) 2016-2020 OpenCFD Ltd.
+    Copyright (C) 2016-2021 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
@@ -149,7 +149,7 @@ Foam::solverPerformance Foam::smoothSolver::solve
             solverPerf.finalResidual() = solverPerf.initialResidual();
         }
 
-        if (lduMatrix::debug >= 2)
+        if ((log_ >= 2) || (lduMatrix::debug >= 2))
         {
             Info.masterStream(matrix().mesh().comm())
                 << "   Normalisation factor = " << normFactor << endl;
@@ -160,7 +160,7 @@ Foam::solverPerformance Foam::smoothSolver::solve
         if
         (
             minIter_ > 0
-         || !solverPerf.checkConvergence(tolerance_, relTol_)
+         || !solverPerf.checkConvergence(tolerance_, relTol_, log_)
         )
         {
             addProfiling(solve, "lduMatrix::smoother." + fieldName_);
@@ -203,7 +203,7 @@ Foam::solverPerformance Foam::smoothSolver::solve
             (
                 (
                     (solverPerf.nIterations() += nSweeps_) < maxIter_
-                && !solverPerf.checkConvergence(tolerance_, relTol_)
+                && !solverPerf.checkConvergence(tolerance_, relTol_, log_)
                 )
              || solverPerf.nIterations() < minIter_
             );
diff --git a/src/finiteVolume/fvMatrices/fvMatrix/fvMatrixSolve.C b/src/finiteVolume/fvMatrices/fvMatrix/fvMatrixSolve.C
index 07c552b52de..111761404ed 100644
--- a/src/finiteVolume/fvMatrices/fvMatrix/fvMatrixSolve.C
+++ b/src/finiteVolume/fvMatrices/fvMatrix/fvMatrixSolve.C
@@ -6,7 +6,7 @@
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
     Copyright (C) 2011-2016 OpenFOAM Foundation
-    Copyright (C) 2016-2020 OpenCFD Ltd.
+    Copyright (C) 2016-2021 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
@@ -222,7 +222,14 @@ Foam::SolverPerformance<Type> Foam::fvMatrix<Type>::solveSegregated
             solverControls
         )->solve(psiCmpt, sourceCmpt, cmpt);
 
-        if (SolverPerformance<Type>::debug)
+        const label log =
+            solverControls.getOrDefault<label>
+            (
+                "log",
+                SolverPerformance<Type>::debug
+            );
+
+        if (log)
         {
             solverPerf.print(Info.masterStream(this->mesh().comm()));
         }
@@ -289,7 +296,14 @@ Foam::SolverPerformance<Type> Foam::fvMatrix<Type>::solveCoupled
         coupledMatrixSolver->solve(psi)
     );
 
-    if (SolverPerformance<Type>::debug)
+    const label log =
+        solverControls.getOrDefault<label>
+        (
+            "log",
+            SolverPerformance<Type>::debug
+        );
+
+    if (log)
     {
         solverPerf.print(Info.masterStream(this->mesh().comm()));
     }
diff --git a/src/finiteVolume/fvMatrices/fvScalarMatrix/fvScalarMatrix.C b/src/finiteVolume/fvMatrices/fvScalarMatrix/fvScalarMatrix.C
index 99f4baf92c2..029b0d004e7 100644
--- a/src/finiteVolume/fvMatrices/fvScalarMatrix/fvScalarMatrix.C
+++ b/src/finiteVolume/fvMatrices/fvScalarMatrix/fvScalarMatrix.C
@@ -137,7 +137,10 @@ Foam::solverPerformance Foam::fvMatrix<Foam::scalar>::fvSolver::solve
         totalSource
     );
 
-    if (solverPerformance::debug)
+    const label log =
+        solverControls.getOrDefault<label>("log", solverPerformance::debug);
+
+    if (log)
     {
         solverPerf.print(Info.masterStream(fvMat_.mesh().comm()));
     }
@@ -264,7 +267,10 @@ Foam::solverPerformance Foam::fvMatrix<Foam::scalar>::solveSegregated
         }
     }
 
-    if (solverPerformance::debug)
+    const label log =
+        solverControls.getOrDefault<label>("log", solverPerformance::debug);
+
+    if (log)
     {
         solverPerf.print(Info.masterStream(mesh().comm()));
     }
-- 
GitLab