diff --git a/src/OpenFOAM/Make/files b/src/OpenFOAM/Make/files index e530cb2a00294b66f4ecb12171604b20375c9041..bcd75424a5572d983aeb2fb2de0b8de792df499a 100644 --- a/src/OpenFOAM/Make/files +++ b/src/OpenFOAM/Make/files @@ -271,6 +271,7 @@ $(lduMatrix)/solvers/diagonalSolver/diagonalSolver.C $(lduMatrix)/solvers/smoothSolver/smoothSolver.C $(lduMatrix)/solvers/PCG/PCG.C $(lduMatrix)/solvers/PBiCG/PBiCG.C +$(lduMatrix)/solvers/PBiCGStab/PBiCGStab.C $(lduMatrix)/smoothers/GaussSeidel/GaussSeidelSmoother.C $(lduMatrix)/smoothers/symGaussSeidel/symGaussSeidelSmoother.C diff --git a/src/OpenFOAM/matrices/lduMatrix/solvers/PBiCG/PBiCG.C b/src/OpenFOAM/matrices/lduMatrix/solvers/PBiCG/PBiCG.C index 3d43c64a0db691e2ab73c0b824435c18baaec7dd..092877e7a8a8be4c98300efd7b7efba06975f9e8 100644 --- a/src/OpenFOAM/matrices/lduMatrix/solvers/PBiCG/PBiCG.C +++ b/src/OpenFOAM/matrices/lduMatrix/solvers/PBiCG/PBiCG.C @@ -76,37 +76,25 @@ Foam::solverPerformance Foam::PBiCG::solve fieldName_ ); - label nCells = psi.size(); + const label nCells = psi.size(); scalar* __restrict__ psiPtr = psi.begin(); scalarField pA(nCells); scalar* __restrict__ pAPtr = pA.begin(); - scalarField pT(nCells, 0.0); - scalar* __restrict__ pTPtr = pT.begin(); - scalarField wA(nCells); scalar* __restrict__ wAPtr = wA.begin(); - scalarField wT(nCells); - scalar* __restrict__ wTPtr = wT.begin(); - - scalar wArT = solverPerf.great_; - scalar wArTold = wArT; - - // --- Calculate A.psi and T.psi + // --- Calculate A.psi matrix_.Amul(wA, psi, interfaceBouCoeffs_, interfaces_, cmpt); - matrix_.Tmul(wT, psi, interfaceIntCoeffs_, interfaces_, cmpt); - // --- Calculate initial residual and transpose residual fields + // --- Calculate initial residual field scalarField rA(source - wA); - scalarField rT(source - wT); scalar* __restrict__ rAPtr = rA.begin(); - scalar* __restrict__ rTPtr = rT.begin(); // --- Calculate normalisation factor - scalar normFactor = this->normFactor(psi, source, wA, pA); + const scalar normFactor = this->normFactor(psi, source, wA, pA); if (lduMatrix::debug >= 2) { @@ -126,6 +114,22 @@ Foam::solverPerformance Foam::PBiCG::solve || !solverPerf.checkConvergence(tolerance_, relTol_) ) { + scalarField pT(nCells, 0); + scalar* __restrict__ pTPtr = pT.begin(); + + scalarField wT(nCells); + scalar* __restrict__ wTPtr = wT.begin(); + + // --- Calculate T.psi + matrix_.Tmul(wT, psi, interfaceIntCoeffs_, interfaces_, cmpt); + + // --- Calculate initial transpose residual field + scalarField rT(source - wT); + scalar* __restrict__ rTPtr = rT.begin(); + + // --- Initial value not used + scalar wArT = 0; + // --- Select and construct the preconditioner autoPtr<lduMatrix::preconditioner> preconPtr = lduMatrix::preconditioner::New @@ -138,7 +142,7 @@ Foam::solverPerformance Foam::PBiCG::solve do { // --- Store previous wArT - wArTold = wArT; + const scalar wArTold = wArT; // --- Precondition residuals preconPtr->precondition(wA, rA, cmpt); @@ -157,7 +161,7 @@ Foam::solverPerformance Foam::PBiCG::solve } else { - scalar beta = wArT/wArTold; + const scalar beta = wArT/wArTold; for (label cell=0; cell<nCells; cell++) { @@ -171,7 +175,7 @@ Foam::solverPerformance Foam::PBiCG::solve matrix_.Amul(wA, pA, interfaceBouCoeffs_, interfaces_, cmpt); matrix_.Tmul(wT, pT, interfaceIntCoeffs_, interfaces_, cmpt); - scalar wApT = gSumProd(wA, pT, matrix().mesh().comm()); + const scalar wApT = gSumProd(wA, pT, matrix().mesh().comm()); // --- Test for singularity if (solverPerf.checkSingularity(mag(wApT)/normFactor)) @@ -182,7 +186,7 @@ Foam::solverPerformance Foam::PBiCG::solve // --- Update solution and residual: - scalar alpha = wArT/wApT; + const scalar alpha = wArT/wApT; for (label cell=0; cell<nCells; cell++) { diff --git a/src/OpenFOAM/matrices/lduMatrix/solvers/PBiCGStab/PBiCGStab.C b/src/OpenFOAM/matrices/lduMatrix/solvers/PBiCGStab/PBiCGStab.C new file mode 100644 index 0000000000000000000000000000000000000000..ff7c0cddd595ca97ad60cb23bb2967c561ba6222 --- /dev/null +++ b/src/OpenFOAM/matrices/lduMatrix/solvers/PBiCGStab/PBiCGStab.C @@ -0,0 +1,252 @@ +/*---------------------------------------------------------------------------*\ + ========= | + \\ / F ield | OpenFOAM: The Open Source CFD Toolbox + \\ / O peration | + \\ / A nd | Copyright (C) 2016 OpenFOAM Foundation + \\/ M anipulation | +------------------------------------------------------------------------------- +License + This file is part of OpenFOAM. + + OpenFOAM is free software: you can redistribute it and/or modify it + under the terms of the GNU General Public License as published by + the Free Software Foundation, either version 3 of the License, or + (at your option) any later version. + + OpenFOAM is distributed in the hope that it will be useful, but WITHOUT + ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or + FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License + for more details. + + You should have received a copy of the GNU General Public License + along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>. + +\*---------------------------------------------------------------------------*/ + +#include "PBiCGStab.H" + +// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * // + +namespace Foam +{ + defineTypeNameAndDebug(PBiCGStab, 0); + + lduMatrix::solver::addasymMatrixConstructorToTable<PBiCGStab> + addPBiCGStabAsymMatrixConstructorToTable_; +} + + +// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // + +Foam::PBiCGStab::PBiCGStab +( + const word& fieldName, + const lduMatrix& matrix, + const FieldField<Field, scalar>& interfaceBouCoeffs, + const FieldField<Field, scalar>& interfaceIntCoeffs, + const lduInterfaceFieldPtrsList& interfaces, + const dictionary& solverControls +) +: + lduMatrix::solver + ( + fieldName, + matrix, + interfaceBouCoeffs, + interfaceIntCoeffs, + interfaces, + solverControls + ) +{} + + +// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // + +Foam::solverPerformance Foam::PBiCGStab::solve +( + scalarField& psi, + const scalarField& source, + const direction cmpt +) const +{ + // --- Setup class containing solver performance data + solverPerformance solverPerf + ( + lduMatrix::preconditioner::getName(controlDict_) + typeName, + fieldName_ + ); + + const label nCells = psi.size(); + + scalar* __restrict__ psiPtr = psi.begin(); + + scalarField pA(nCells); + scalar* __restrict__ pAPtr = pA.begin(); + + scalarField yA(nCells); + scalar* __restrict__ yAPtr = yA.begin(); + + // --- Calculate A.psi + matrix_.Amul(yA, psi, interfaceBouCoeffs_, interfaces_, cmpt); + + // --- Calculate initial residual field + scalarField rA(source - yA); + scalar* __restrict__ rAPtr = rA.begin(); + + // --- Calculate normalisation factor + const scalar normFactor = this->normFactor(psi, source, yA, pA); + + if (lduMatrix::debug >= 2) + { + Info<< " Normalisation factor = " << normFactor << endl; + } + + // --- Calculate normalised residual norm + solverPerf.initialResidual() = + gSumMag(rA, matrix().mesh().comm()) + /normFactor; + solverPerf.finalResidual() = solverPerf.initialResidual(); + + // --- Check convergence, solve if not converged + if + ( + minIter_ > 0 + || !solverPerf.checkConvergence(tolerance_, relTol_) + ) + { + scalarField AyA(nCells); + scalar* __restrict__ AyAPtr = AyA.begin(); + + scalarField sA(nCells); + scalar* __restrict__ sAPtr = sA.begin(); + + scalarField zA(nCells); + scalar* __restrict__ zAPtr = zA.begin(); + + scalarField tA(nCells); + scalar* __restrict__ tAPtr = tA.begin(); + + // --- Store initial residual + const scalarField rA0(rA); + + // --- Initial values not used + scalar rA0rA = 0; + scalar alpha = 0; + scalar omega = 0; + + // --- Select and construct the preconditioner + autoPtr<lduMatrix::preconditioner> preconPtr = + lduMatrix::preconditioner::New + ( + *this, + controlDict_ + ); + + // --- Solver iteration + do + { + // --- Store previous rA0rA + const scalar rA0rAold = rA0rA; + + rA0rA = gSumProd(rA0, rA, matrix().mesh().comm()); + + // --- Test for singularity + if (solverPerf.checkSingularity(mag(rA0rA))) + { + break; + } + + // --- Update pA + if (solverPerf.nIterations() == 0) + { + for (label cell=0; cell<nCells; cell++) + { + pAPtr[cell] = rAPtr[cell]; + } + } + else + { + // --- Test for singularity + if (solverPerf.checkSingularity(mag(omega))) + { + break; + } + + const scalar beta = (rA0rA/rA0rAold)*(alpha/omega); + + for (label cell=0; cell<nCells; cell++) + { + pAPtr[cell] = + rAPtr[cell] + beta*(pAPtr[cell] - omega*AyAPtr[cell]); + } + } + + // --- Precondition pA + preconPtr->precondition(yA, pA, cmpt); + + // --- Calculate AyA + matrix_.Amul(AyA, yA, interfaceBouCoeffs_, interfaces_, cmpt); + + const scalar rA0AyA = gSumProd(rA0, AyA, matrix().mesh().comm()); + + alpha = rA0rA/rA0AyA; + + // --- Calculate sA + for (label cell=0; cell<nCells; cell++) + { + sAPtr[cell] = rAPtr[cell] - alpha*AyAPtr[cell]; + } + + // --- Test sA for convergence + solverPerf.finalResidual() = + gSumMag(sA, matrix().mesh().comm())/normFactor; + + if (solverPerf.checkConvergence(tolerance_, relTol_)) + { + for (label cell=0; cell<nCells; cell++) + { + psiPtr[cell] += alpha*yAPtr[cell]; + } + + solverPerf.nIterations()++; + + return solverPerf; + } + + // --- Precondition sA + preconPtr->precondition(zA, sA, cmpt); + + // --- Calculate tA + matrix_.Amul(tA, zA, interfaceBouCoeffs_, interfaces_, cmpt); + + const scalar tAtA = gSumSqr(tA, matrix().mesh().comm()); + + // --- Calculate omega from tA and sA + // (cheaper than using zA with preconditioned tA) + omega = gSumProd(tA, sA)/tAtA; + + // --- Update solution and residual + for (label cell=0; cell<nCells; cell++) + { + psiPtr[cell] += alpha*yAPtr[cell] + omega*zAPtr[cell]; + rAPtr[cell] = sAPtr[cell] - omega*tAPtr[cell]; + } + + solverPerf.finalResidual() = + gSumMag(rA, matrix().mesh().comm()) + /normFactor; + } while + ( + ( + solverPerf.nIterations()++ < maxIter_ + && !solverPerf.checkConvergence(tolerance_, relTol_) + ) + || solverPerf.nIterations() < minIter_ + ); + } + + return solverPerf; +} + + +// ************************************************************************* // diff --git a/src/OpenFOAM/matrices/lduMatrix/solvers/PBiCGStab/PBiCGStab.H b/src/OpenFOAM/matrices/lduMatrix/solvers/PBiCGStab/PBiCGStab.H new file mode 100644 index 0000000000000000000000000000000000000000..95dee23a500db506d448ef2c9da623fb4b03a486 --- /dev/null +++ b/src/OpenFOAM/matrices/lduMatrix/solvers/PBiCGStab/PBiCGStab.H @@ -0,0 +1,123 @@ +/*---------------------------------------------------------------------------*\ + ========= | + \\ / F ield | OpenFOAM: The Open Source CFD Toolbox + \\ / O peration | + \\ / A nd | Copyright (C) 2016 OpenFOAM Foundation + \\/ M anipulation | +------------------------------------------------------------------------------- +License + This file is part of OpenFOAM. + + OpenFOAM is free software: you can redistribute it and/or modify it + under the terms of the GNU General Public License as published by + the Free Software Foundation, either version 3 of the License, or + (at your option) any later version. + + OpenFOAM is distributed in the hope that it will be useful, but WITHOUT + ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or + FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License + for more details. + + You should have received a copy of the GNU General Public License + along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>. + +Class + Foam::PBiCGStab + +Description + Preconditioned bi-conjugate gradient stabilized solver for asymmetric + lduMatrices using a run-time selectable preconditiioner. + + References: + \verbatim + Van der Vorst, H. A. (1992). + Bi-CGSTAB: A fast and smoothly converging variant of Bi-CG + for the solution of nonsymmetric linear systems. + SIAM Journal on scientific and Statistical Computing, 13(2), 631-644. + + Barrett, R., Berry, M. W., Chan, T. F., Demmel, J., Donato, J., + Dongarra, J., Eijkhout, V., Pozo, R., Romine, C. & Van der Vorst, H. + (1994). + Templates for the solution of linear systems: + building blocks for iterative methods + (Vol. 43). Siam. + \endverbatim + +SourceFiles + PBiCGStab.C + +\*---------------------------------------------------------------------------*/ + +#ifndef PBiCGStab_H +#define PBiCGStab_H + +#include "lduMatrix.H" + +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +namespace Foam +{ + +/*---------------------------------------------------------------------------*\ + Class PBiCGStab Declaration +\*---------------------------------------------------------------------------*/ + +class PBiCGStab +: + public lduMatrix::solver +{ + // Private Member Functions + + //- Disallow default bitwise copy construct + PBiCGStab(const PBiCGStab&); + + //- Disallow default bitwise assignment + void operator=(const PBiCGStab&); + + +public: + + //- Runtime type information + TypeName("PBiCGStab"); + + + // Constructors + + //- Construct from matrix components and solver data stream + PBiCGStab + ( + const word& fieldName, + const lduMatrix& matrix, + const FieldField<Field, scalar>& interfaceBouCoeffs, + const FieldField<Field, scalar>& interfaceIntCoeffs, + const lduInterfaceFieldPtrsList& interfaces, + const dictionary& solverControls + ); + + + //- Destructor + virtual ~PBiCGStab() + {} + + + // Member Functions + + //- Solve the matrix with this solver + virtual solverPerformance solve + ( + scalarField& psi, + const scalarField& source, + const direction cmpt=0 + ) const; +}; + + +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +} // End namespace Foam + +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +#endif + +// ************************************************************************* // diff --git a/tutorials/basic/scalarTransportFoam/pitzDaily/system/fvSolution b/tutorials/basic/scalarTransportFoam/pitzDaily/system/fvSolution index 245dcb70fe23550e59cb2b5ff2390a01a9eeadaf..07d6974305ef13549f57ee0379813e942b833b84 100644 --- a/tutorials/basic/scalarTransportFoam/pitzDaily/system/fvSolution +++ b/tutorials/basic/scalarTransportFoam/pitzDaily/system/fvSolution @@ -19,7 +19,7 @@ solvers { T { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-06; relTol 0; diff --git a/tutorials/combustion/PDRFoam/flamePropagationWithObstacles/system/fvSolution b/tutorials/combustion/PDRFoam/flamePropagationWithObstacles/system/fvSolution index b76b24182ecfb7bcbb659b52e59a8cf2a662d5c7..a7a0765cd3eb5e3ebc9d4b60f576a85d3485f785 100644 --- a/tutorials/combustion/PDRFoam/flamePropagationWithObstacles/system/fvSolution +++ b/tutorials/combustion/PDRFoam/flamePropagationWithObstacles/system/fvSolution @@ -49,7 +49,7 @@ solvers "(b|Xi|ft|ha|hau|k|epsilon)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-06; relTol 0.1; @@ -57,7 +57,7 @@ solvers "(b|Xi|ft|ha|hau|k|epsilon)Final" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-06; relTol 0; diff --git a/tutorials/combustion/XiFoam/ras/moriyoshiHomogeneous/system/fvSolution b/tutorials/combustion/XiFoam/ras/moriyoshiHomogeneous/system/fvSolution index 706044c8b1cfd88a0aa68257d304f5a1eccd6751..800bf49f3529ef70635863b08619d9daff8293d8 100644 --- a/tutorials/combustion/XiFoam/ras/moriyoshiHomogeneous/system/fvSolution +++ b/tutorials/combustion/XiFoam/ras/moriyoshiHomogeneous/system/fvSolution @@ -34,7 +34,7 @@ solvers "(U|b|Su|Xi|ha|hau|k|epsilon)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-05; relTol 0.1; @@ -42,7 +42,7 @@ solvers "(U|b|Su|Xi|ha|hau|k|epsilon)Final" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-05; relTol 0; diff --git a/tutorials/combustion/chemFoam/gri/system/fvSolution b/tutorials/combustion/chemFoam/gri/system/fvSolution index c090558aa22c2a5bfd8c3a1f87a1444f7315c655..6b16dd4035184065e7d5c3f4b52f2c3c1fe89d30 100644 --- a/tutorials/combustion/chemFoam/gri/system/fvSolution +++ b/tutorials/combustion/chemFoam/gri/system/fvSolution @@ -19,7 +19,7 @@ solvers { Yi { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-12; relTol 0; diff --git a/tutorials/combustion/chemFoam/h2/system/fvSolution b/tutorials/combustion/chemFoam/h2/system/fvSolution index c090558aa22c2a5bfd8c3a1f87a1444f7315c655..6b16dd4035184065e7d5c3f4b52f2c3c1fe89d30 100644 --- a/tutorials/combustion/chemFoam/h2/system/fvSolution +++ b/tutorials/combustion/chemFoam/h2/system/fvSolution @@ -19,7 +19,7 @@ solvers { Yi { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-12; relTol 0; diff --git a/tutorials/combustion/chemFoam/ic8h18/system/fvSolution b/tutorials/combustion/chemFoam/ic8h18/system/fvSolution index c090558aa22c2a5bfd8c3a1f87a1444f7315c655..6b16dd4035184065e7d5c3f4b52f2c3c1fe89d30 100644 --- a/tutorials/combustion/chemFoam/ic8h18/system/fvSolution +++ b/tutorials/combustion/chemFoam/ic8h18/system/fvSolution @@ -19,7 +19,7 @@ solvers { Yi { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-12; relTol 0; diff --git a/tutorials/combustion/chemFoam/ic8h18_TDAC/system/fvSolution b/tutorials/combustion/chemFoam/ic8h18_TDAC/system/fvSolution index c090558aa22c2a5bfd8c3a1f87a1444f7315c655..6b16dd4035184065e7d5c3f4b52f2c3c1fe89d30 100644 --- a/tutorials/combustion/chemFoam/ic8h18_TDAC/system/fvSolution +++ b/tutorials/combustion/chemFoam/ic8h18_TDAC/system/fvSolution @@ -19,7 +19,7 @@ solvers { Yi { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-12; relTol 0; diff --git a/tutorials/combustion/chemFoam/nc7h16/system/fvSolution b/tutorials/combustion/chemFoam/nc7h16/system/fvSolution index c090558aa22c2a5bfd8c3a1f87a1444f7315c655..6b16dd4035184065e7d5c3f4b52f2c3c1fe89d30 100644 --- a/tutorials/combustion/chemFoam/nc7h16/system/fvSolution +++ b/tutorials/combustion/chemFoam/nc7h16/system/fvSolution @@ -19,7 +19,7 @@ solvers { Yi { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-12; relTol 0; diff --git a/tutorials/combustion/engineFoam/kivaTest/system/fvSolution b/tutorials/combustion/engineFoam/kivaTest/system/fvSolution index b164746874796c4bdd6576b211d66fde72d4571d..64ba5515cab21e96f72d4f569a6dde8882d10f76 100644 --- a/tutorials/combustion/engineFoam/kivaTest/system/fvSolution +++ b/tutorials/combustion/engineFoam/kivaTest/system/fvSolution @@ -49,7 +49,7 @@ solvers "(U|Xi|hau|eau|ft|b|ha|ea|k|epsilon)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-05; relTol 0.1; diff --git a/tutorials/combustion/fireFoam/les/flameSpreadWaterSuppressionPanel/system/filmRegion/fvSolution b/tutorials/combustion/fireFoam/les/flameSpreadWaterSuppressionPanel/system/filmRegion/fvSolution index e247ccbf42dc21cd8650ab503b283c1d749a6009..0dd063c8aa52389ab18a55b0aae72bb7457c5a06 100644 --- a/tutorials/combustion/fireFoam/les/flameSpreadWaterSuppressionPanel/system/filmRegion/fvSolution +++ b/tutorials/combustion/fireFoam/les/flameSpreadWaterSuppressionPanel/system/filmRegion/fvSolution @@ -19,21 +19,21 @@ solvers { hf { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 0; relTol 1e-3; } "(Uf|deltaf\*rhof)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-10; relTol 0; } deltaf { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-10; relTol 0; diff --git a/tutorials/combustion/fireFoam/les/flameSpreadWaterSuppressionPanel/system/fvSolution b/tutorials/combustion/fireFoam/les/flameSpreadWaterSuppressionPanel/system/fvSolution index 12f3b66435f7994e31a82213f9c9d077aee2882b..fc3a6c4ab9cc0d8cab86f8b92c7dca1d06296fd4 100644 --- a/tutorials/combustion/fireFoam/les/flameSpreadWaterSuppressionPanel/system/fvSolution +++ b/tutorials/combustion/fireFoam/les/flameSpreadWaterSuppressionPanel/system/fvSolution @@ -57,7 +57,7 @@ solvers "(U|Yi|h|k)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-8; relTol 0.1; diff --git a/tutorials/combustion/fireFoam/les/flameSpreadWaterSuppressionPanel/system/pyrolysisRegion/fvSolution b/tutorials/combustion/fireFoam/les/flameSpreadWaterSuppressionPanel/system/pyrolysisRegion/fvSolution index 45714f1a1d82034f501a098771d77ba2afad1d92..b400e253fb31f3ac7001120931346d67f0e0a9a0 100644 --- a/tutorials/combustion/fireFoam/les/flameSpreadWaterSuppressionPanel/system/pyrolysisRegion/fvSolution +++ b/tutorials/combustion/fireFoam/les/flameSpreadWaterSuppressionPanel/system/pyrolysisRegion/fvSolution @@ -26,7 +26,7 @@ solvers "Yi" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-06; relTol 0; diff --git a/tutorials/combustion/fireFoam/les/oppositeBurningPanels/system/fvSolution b/tutorials/combustion/fireFoam/les/oppositeBurningPanels/system/fvSolution index 1329cb1cab1ceb08796acec9a3bfeca96798b452..7bacd5c63a95696b350ae4419d36b9ff2e4422b9 100644 --- a/tutorials/combustion/fireFoam/les/oppositeBurningPanels/system/fvSolution +++ b/tutorials/combustion/fireFoam/les/oppositeBurningPanels/system/fvSolution @@ -47,7 +47,7 @@ solvers "(U|Yi|k|h|omega)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-6; relTol 0.1; diff --git a/tutorials/combustion/fireFoam/les/oppositeBurningPanels/system/panelRegion/fvSolution b/tutorials/combustion/fireFoam/les/oppositeBurningPanels/system/panelRegion/fvSolution index 88bd5c58e98a1b72616005ab16052fa21575b759..08c57ca2b0587dce5082df028e82853834653a2f 100644 --- a/tutorials/combustion/fireFoam/les/oppositeBurningPanels/system/panelRegion/fvSolution +++ b/tutorials/combustion/fireFoam/les/oppositeBurningPanels/system/panelRegion/fvSolution @@ -26,7 +26,7 @@ solvers "Yi" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-06; relTol 0; diff --git a/tutorials/combustion/fireFoam/les/smallPoolFire2D/system/fvSolution b/tutorials/combustion/fireFoam/les/smallPoolFire2D/system/fvSolution index dcf0353bf2e5acb8761930c76dc983f79308e3ea..4506a37939c095d03b1fe7d261c21241338193ef 100644 --- a/tutorials/combustion/fireFoam/les/smallPoolFire2D/system/fvSolution +++ b/tutorials/combustion/fireFoam/les/smallPoolFire2D/system/fvSolution @@ -56,7 +56,7 @@ solvers "(U|Yi|k|h|omega)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-6; relTol 0.1; diff --git a/tutorials/combustion/fireFoam/les/smallPoolFire3D/system/fvSolution b/tutorials/combustion/fireFoam/les/smallPoolFire3D/system/fvSolution index a26c2f30b91d5c6d828c7479cb0b4be37a03a47d..e89bec79c85fb56b2b5b5335395e2704fd857dfb 100644 --- a/tutorials/combustion/fireFoam/les/smallPoolFire3D/system/fvSolution +++ b/tutorials/combustion/fireFoam/les/smallPoolFire3D/system/fvSolution @@ -48,7 +48,7 @@ solvers "(U|Yi|k|h)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-6; relTol 0.1; diff --git a/tutorials/combustion/reactingFoam/laminar/counterFlowFlame2D/system/fvSolution b/tutorials/combustion/reactingFoam/laminar/counterFlowFlame2D/system/fvSolution index 7491893cc3ea0b75f0b6838e5d5ed6402b0f6cde..7ff37c8434d5caf54eb170a11543f915a259db9a 100644 --- a/tutorials/combustion/reactingFoam/laminar/counterFlowFlame2D/system/fvSolution +++ b/tutorials/combustion/reactingFoam/laminar/counterFlowFlame2D/system/fvSolution @@ -39,7 +39,7 @@ solvers "(U|h|k|epsilon)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-6; relTol 0.1; diff --git a/tutorials/combustion/reactingFoam/laminar/counterFlowFlame2DLTS/system/fvSolution b/tutorials/combustion/reactingFoam/laminar/counterFlowFlame2DLTS/system/fvSolution index 4d0db2d721e1377f65f469c4f5d4a9d4ba96017d..38d9937828bc8d5e99d92eaa8920ff4644a36e1a 100644 --- a/tutorials/combustion/reactingFoam/laminar/counterFlowFlame2DLTS/system/fvSolution +++ b/tutorials/combustion/reactingFoam/laminar/counterFlowFlame2DLTS/system/fvSolution @@ -38,7 +38,7 @@ solvers "(U|h|k|epsilon)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-6; relTol 0.1; diff --git a/tutorials/combustion/reactingFoam/laminar/counterFlowFlame2D_GRI/system/fvSolution b/tutorials/combustion/reactingFoam/laminar/counterFlowFlame2D_GRI/system/fvSolution index 7491893cc3ea0b75f0b6838e5d5ed6402b0f6cde..7ff37c8434d5caf54eb170a11543f915a259db9a 100644 --- a/tutorials/combustion/reactingFoam/laminar/counterFlowFlame2D_GRI/system/fvSolution +++ b/tutorials/combustion/reactingFoam/laminar/counterFlowFlame2D_GRI/system/fvSolution @@ -39,7 +39,7 @@ solvers "(U|h|k|epsilon)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-6; relTol 0.1; diff --git a/tutorials/combustion/reactingFoam/laminar/counterFlowFlame2D_GRI_TDAC/system/fvSolution b/tutorials/combustion/reactingFoam/laminar/counterFlowFlame2D_GRI_TDAC/system/fvSolution index 7491893cc3ea0b75f0b6838e5d5ed6402b0f6cde..7ff37c8434d5caf54eb170a11543f915a259db9a 100644 --- a/tutorials/combustion/reactingFoam/laminar/counterFlowFlame2D_GRI_TDAC/system/fvSolution +++ b/tutorials/combustion/reactingFoam/laminar/counterFlowFlame2D_GRI_TDAC/system/fvSolution @@ -39,7 +39,7 @@ solvers "(U|h|k|epsilon)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-6; relTol 0.1; diff --git a/tutorials/compressible/rhoCentralDyMFoam/movingCone/system/fvSolution b/tutorials/compressible/rhoCentralDyMFoam/movingCone/system/fvSolution index 5161002b216857457dd5c1b06362ba1b33354433..7a434bd220b44cef12a45e18160832ad30713ded 100644 --- a/tutorials/compressible/rhoCentralDyMFoam/movingCone/system/fvSolution +++ b/tutorials/compressible/rhoCentralDyMFoam/movingCone/system/fvSolution @@ -19,7 +19,7 @@ solvers { p { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-06; relTol 0.01; diff --git a/tutorials/compressible/sonicDyMFoam/movingCone/system/fvSolution b/tutorials/compressible/sonicDyMFoam/movingCone/system/fvSolution index 5161002b216857457dd5c1b06362ba1b33354433..7a434bd220b44cef12a45e18160832ad30713ded 100644 --- a/tutorials/compressible/sonicDyMFoam/movingCone/system/fvSolution +++ b/tutorials/compressible/sonicDyMFoam/movingCone/system/fvSolution @@ -19,7 +19,7 @@ solvers { p { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-06; relTol 0.01; diff --git a/tutorials/compressible/sonicFoam/ras/nacaAirfoil/system/fvSolution b/tutorials/compressible/sonicFoam/ras/nacaAirfoil/system/fvSolution index afc24d24a2e9106bb45a2a1013a763e503715c94..1a0146a6440835352dee48b8dfae77bc30122201 100644 --- a/tutorials/compressible/sonicFoam/ras/nacaAirfoil/system/fvSolution +++ b/tutorials/compressible/sonicFoam/ras/nacaAirfoil/system/fvSolution @@ -24,7 +24,7 @@ solvers "p.*" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-12; relTol 0; diff --git a/tutorials/financial/financialFoam/europeanCall/system/fvSolution b/tutorials/financial/financialFoam/europeanCall/system/fvSolution index 3b95a484798806305a349f9ce565a2f301298505..14d17e1f9771974a39577eea52c6b9c8f77e96e2 100644 --- a/tutorials/financial/financialFoam/europeanCall/system/fvSolution +++ b/tutorials/financial/financialFoam/europeanCall/system/fvSolution @@ -19,7 +19,7 @@ solvers { V { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-05; relTol 0; diff --git a/tutorials/heatTransfer/buoyantBoussinesqPimpleFoam/hotRoom/system/fvSolution b/tutorials/heatTransfer/buoyantBoussinesqPimpleFoam/hotRoom/system/fvSolution index 324df46ada9fcf1a29be08f2e58a2d71e7f4b3c3..dfaf067c8715f588ee7c6beda91bc79dcd496726 100644 --- a/tutorials/heatTransfer/buoyantBoussinesqPimpleFoam/hotRoom/system/fvSolution +++ b/tutorials/heatTransfer/buoyantBoussinesqPimpleFoam/hotRoom/system/fvSolution @@ -33,7 +33,7 @@ solvers "(U|T|k|epsilon|R)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-6; relTol 0.1; diff --git a/tutorials/heatTransfer/buoyantBoussinesqSimpleFoam/hotRoom/system/fvSolution b/tutorials/heatTransfer/buoyantBoussinesqSimpleFoam/hotRoom/system/fvSolution index a5df85e653117b7d6f318de1b11152b1583e0431..6965c865a5e01fe37978a72dbea9694b62b05748 100644 --- a/tutorials/heatTransfer/buoyantBoussinesqSimpleFoam/hotRoom/system/fvSolution +++ b/tutorials/heatTransfer/buoyantBoussinesqSimpleFoam/hotRoom/system/fvSolution @@ -27,7 +27,7 @@ solvers "(U|T|k|epsilon|R)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-05; relTol 0.1; diff --git a/tutorials/heatTransfer/buoyantBoussinesqSimpleFoam/iglooWithFridges/system/fvSolution b/tutorials/heatTransfer/buoyantBoussinesqSimpleFoam/iglooWithFridges/system/fvSolution index 627615c534bed3dab1b2f02f3a01f3cbef507916..5e871f1bb6914c10bf76214f526d9f2ae6e1d4a3 100644 --- a/tutorials/heatTransfer/buoyantBoussinesqSimpleFoam/iglooWithFridges/system/fvSolution +++ b/tutorials/heatTransfer/buoyantBoussinesqSimpleFoam/iglooWithFridges/system/fvSolution @@ -27,7 +27,7 @@ solvers "(U|T|k|epsilon)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-07; relTol 0.1; diff --git a/tutorials/heatTransfer/buoyantPimpleFoam/hotRoom/system/fvSolution b/tutorials/heatTransfer/buoyantPimpleFoam/hotRoom/system/fvSolution index 7127fdc8c6fa54fac6586cea7bf6ae1bcd6aa718..4e0487cd9e1fae91d127d7445bfbc65fe1041fd9 100644 --- a/tutorials/heatTransfer/buoyantPimpleFoam/hotRoom/system/fvSolution +++ b/tutorials/heatTransfer/buoyantPimpleFoam/hotRoom/system/fvSolution @@ -41,7 +41,7 @@ solvers "(U|h|e|k|epsilon|R)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-6; relTol 0.1; diff --git a/tutorials/heatTransfer/buoyantSimpleFoam/buoyantCavity/system/fvSolution b/tutorials/heatTransfer/buoyantSimpleFoam/buoyantCavity/system/fvSolution index 90ffdd49b0831663f15e784fec5dd5a05f4233fd..17183ab5c4d0efd855be5ded0f769f49f0bbd9b6 100644 --- a/tutorials/heatTransfer/buoyantSimpleFoam/buoyantCavity/system/fvSolution +++ b/tutorials/heatTransfer/buoyantSimpleFoam/buoyantCavity/system/fvSolution @@ -29,7 +29,7 @@ solvers "(U|h|k|epsilon|omega)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-8; relTol 0.1; diff --git a/tutorials/heatTransfer/buoyantSimpleFoam/circuitBoardCooling/system/fvSolution b/tutorials/heatTransfer/buoyantSimpleFoam/circuitBoardCooling/system/fvSolution index 66346036c90a7212e04d87d9b06cbf3db1640356..7f853270fc2b34ca5a0ac041cadbe75b8e9071fc 100644 --- a/tutorials/heatTransfer/buoyantSimpleFoam/circuitBoardCooling/system/fvSolution +++ b/tutorials/heatTransfer/buoyantSimpleFoam/circuitBoardCooling/system/fvSolution @@ -29,7 +29,7 @@ solvers "(U|h|k|epsilon|omega)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-7; relTol 0.01; diff --git a/tutorials/heatTransfer/buoyantSimpleFoam/externalCoupledCavity/system/fvSolution b/tutorials/heatTransfer/buoyantSimpleFoam/externalCoupledCavity/system/fvSolution index d9b7d98e9f075544b0d922e1d395a78ccec4dbd3..de519c1e1598be02fe539b75c12bfd2988d728fe 100644 --- a/tutorials/heatTransfer/buoyantSimpleFoam/externalCoupledCavity/system/fvSolution +++ b/tutorials/heatTransfer/buoyantSimpleFoam/externalCoupledCavity/system/fvSolution @@ -29,7 +29,7 @@ solvers "(U|h|k|epsilon|omega)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-8; relTol 0.1; diff --git a/tutorials/heatTransfer/buoyantSimpleFoam/hotRadiationRoom/system/fvSolution b/tutorials/heatTransfer/buoyantSimpleFoam/hotRadiationRoom/system/fvSolution index c4fd0cde06011fc0bbe94511b86954b0fb8effc2..24cdb30297c18f80a12bbde3ea2686bd48ea3619 100644 --- a/tutorials/heatTransfer/buoyantSimpleFoam/hotRadiationRoom/system/fvSolution +++ b/tutorials/heatTransfer/buoyantSimpleFoam/hotRadiationRoom/system/fvSolution @@ -27,7 +27,7 @@ solvers "(U|h|k|epsilon)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-05; relTol 0.1; diff --git a/tutorials/heatTransfer/buoyantSimpleFoam/hotRadiationRoomFvDOM/system/fvSolution b/tutorials/heatTransfer/buoyantSimpleFoam/hotRadiationRoomFvDOM/system/fvSolution index 42b3fef2852758f6ed38bd295e844d6ae94396a1..f610c939f53f24c26342989c0d247a1e9e11c759 100644 --- a/tutorials/heatTransfer/buoyantSimpleFoam/hotRadiationRoomFvDOM/system/fvSolution +++ b/tutorials/heatTransfer/buoyantSimpleFoam/hotRadiationRoomFvDOM/system/fvSolution @@ -37,7 +37,7 @@ solvers "(U|h|k|epsilon)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-05; relTol 0.1; diff --git a/tutorials/heatTransfer/chtMultiRegionFoam/multiRegionHeater/system/bottomWater/fvSolution b/tutorials/heatTransfer/chtMultiRegionFoam/multiRegionHeater/system/bottomWater/fvSolution index 07635c7d234d678dfe07090c67dd8c1ee94e51e8..0eb178deb7442d6482f7419353e96187dabf0148 100644 --- a/tutorials/heatTransfer/chtMultiRegionFoam/multiRegionHeater/system/bottomWater/fvSolution +++ b/tutorials/heatTransfer/chtMultiRegionFoam/multiRegionHeater/system/bottomWater/fvSolution @@ -50,7 +50,7 @@ solvers "(U|h|k|epsilon|R)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-7; relTol 0.1; diff --git a/tutorials/heatTransfer/chtMultiRegionFoam/multiRegionHeater/system/topAir/fvSolution b/tutorials/heatTransfer/chtMultiRegionFoam/multiRegionHeater/system/topAir/fvSolution index 07635c7d234d678dfe07090c67dd8c1ee94e51e8..0eb178deb7442d6482f7419353e96187dabf0148 100644 --- a/tutorials/heatTransfer/chtMultiRegionFoam/multiRegionHeater/system/topAir/fvSolution +++ b/tutorials/heatTransfer/chtMultiRegionFoam/multiRegionHeater/system/topAir/fvSolution @@ -50,7 +50,7 @@ solvers "(U|h|k|epsilon|R)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-7; relTol 0.1; diff --git a/tutorials/heatTransfer/chtMultiRegionFoam/snappyMultiRegionHeater/system/bottomAir/fvSolution b/tutorials/heatTransfer/chtMultiRegionFoam/snappyMultiRegionHeater/system/bottomAir/fvSolution index 1527d53d6d80319db5f9cd037f87f552acb68561..066616d5534bcb0645d1365cffa100d383aeb5e6 100644 --- a/tutorials/heatTransfer/chtMultiRegionFoam/snappyMultiRegionHeater/system/bottomAir/fvSolution +++ b/tutorials/heatTransfer/chtMultiRegionFoam/snappyMultiRegionHeater/system/bottomAir/fvSolution @@ -43,7 +43,7 @@ solvers "(U|h|k|epsilon|R)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-7; relTol 0.1; diff --git a/tutorials/heatTransfer/chtMultiRegionFoam/snappyMultiRegionHeater/system/topAir/fvSolution b/tutorials/heatTransfer/chtMultiRegionFoam/snappyMultiRegionHeater/system/topAir/fvSolution index 1527d53d6d80319db5f9cd037f87f552acb68561..066616d5534bcb0645d1365cffa100d383aeb5e6 100644 --- a/tutorials/heatTransfer/chtMultiRegionFoam/snappyMultiRegionHeater/system/topAir/fvSolution +++ b/tutorials/heatTransfer/chtMultiRegionFoam/snappyMultiRegionHeater/system/topAir/fvSolution @@ -43,7 +43,7 @@ solvers "(U|h|k|epsilon|R)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-7; relTol 0.1; diff --git a/tutorials/heatTransfer/chtMultiRegionSimpleFoam/heatExchanger/system/air/fvSolution b/tutorials/heatTransfer/chtMultiRegionSimpleFoam/heatExchanger/system/air/fvSolution index f8de0ffe7c8e45632d0e73ae52887b5bf5630d34..895988eace5b394caaad847f9d2f88e3b26708de 100644 --- a/tutorials/heatTransfer/chtMultiRegionSimpleFoam/heatExchanger/system/air/fvSolution +++ b/tutorials/heatTransfer/chtMultiRegionSimpleFoam/heatExchanger/system/air/fvSolution @@ -31,7 +31,7 @@ solvers "(U|h|e|k|epsilon)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-6; relTol 0.1; diff --git a/tutorials/heatTransfer/chtMultiRegionSimpleFoam/heatExchanger/system/porous/fvSolution b/tutorials/heatTransfer/chtMultiRegionSimpleFoam/heatExchanger/system/porous/fvSolution index f912142c5f32cf6bbca826bbe20be89159490eb1..945f8ea24715628718d49d19dccfc6b7a94baf96 100644 --- a/tutorials/heatTransfer/chtMultiRegionSimpleFoam/heatExchanger/system/porous/fvSolution +++ b/tutorials/heatTransfer/chtMultiRegionSimpleFoam/heatExchanger/system/porous/fvSolution @@ -30,7 +30,7 @@ solvers "(U|h|e|k|epsilon)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-7; relTol 0.1; diff --git a/tutorials/heatTransfer/chtMultiRegionSimpleFoam/multiRegionHeaterRadiation/system/bottomAir/fvSolution b/tutorials/heatTransfer/chtMultiRegionSimpleFoam/multiRegionHeaterRadiation/system/bottomAir/fvSolution index fc149dff1e764dbfee6c96f7c57b1099cc24d3e2..2f353374cdcebe05c339674695d1313b590f1ccf 100644 --- a/tutorials/heatTransfer/chtMultiRegionSimpleFoam/multiRegionHeaterRadiation/system/bottomAir/fvSolution +++ b/tutorials/heatTransfer/chtMultiRegionSimpleFoam/multiRegionHeaterRadiation/system/bottomAir/fvSolution @@ -36,7 +36,7 @@ solvers "(U|h|k|epsilon|G|Ii)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-7; relTol 0.1; diff --git a/tutorials/heatTransfer/chtMultiRegionSimpleFoam/multiRegionHeaterRadiation/system/topAir/fvSolution b/tutorials/heatTransfer/chtMultiRegionSimpleFoam/multiRegionHeaterRadiation/system/topAir/fvSolution index fc149dff1e764dbfee6c96f7c57b1099cc24d3e2..2f353374cdcebe05c339674695d1313b590f1ccf 100644 --- a/tutorials/heatTransfer/chtMultiRegionSimpleFoam/multiRegionHeaterRadiation/system/topAir/fvSolution +++ b/tutorials/heatTransfer/chtMultiRegionSimpleFoam/multiRegionHeaterRadiation/system/topAir/fvSolution @@ -36,7 +36,7 @@ solvers "(U|h|k|epsilon|G|Ii)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-7; relTol 0.1; diff --git a/tutorials/incompressible/simpleFoam/motorBike/system/fvSolution b/tutorials/incompressible/simpleFoam/motorBike/system/fvSolution index 1aa8af7b2726c5072aef90e9f6d073f8b5e1d801..e19623c35ddad2202c6a64c05521ca22d339764f 100644 --- a/tutorials/incompressible/simpleFoam/motorBike/system/fvSolution +++ b/tutorials/incompressible/simpleFoam/motorBike/system/fvSolution @@ -18,10 +18,10 @@ solvers { p { - solver GAMG; - tolerance 1e-7; - relTol 0.01; - smoother GaussSeidel; + solver GAMG; + smoother GaussSeidel; + tolerance 1e-7; + relTol 0.01; } Phi @@ -31,29 +31,29 @@ solvers U { - solver smoothSolver; - smoother GaussSeidel; - tolerance 1e-8; - relTol 0.1; - nSweeps 1; + solver smoothSolver; + smoother GaussSeidel; + tolerance 1e-8; + relTol 0.1; + nSweeps 1; } k { - solver smoothSolver; - smoother GaussSeidel; - tolerance 1e-8; - relTol 0.1; - nSweeps 1; + solver smoothSolver; + smoother GaussSeidel; + tolerance 1e-8; + relTol 0.1; + nSweeps 1; } omega { - solver smoothSolver; - smoother GaussSeidel; - tolerance 1e-8; - relTol 0.1; - nSweeps 1; + solver smoothSolver; + smoother GaussSeidel; + tolerance 1e-8; + relTol 0.1; + nSweeps 1; } } diff --git a/tutorials/lagrangian/reactingParcelFilmFoam/cylinder/system/fvSolution b/tutorials/lagrangian/reactingParcelFilmFoam/cylinder/system/fvSolution index 63c6e0ce0f8045650415279b5a3eb00a297cbc77..a1faf2da4e8194285e5cb7337f8201f7b20e847a 100644 --- a/tutorials/lagrangian/reactingParcelFilmFoam/cylinder/system/fvSolution +++ b/tutorials/lagrangian/reactingParcelFilmFoam/cylinder/system/fvSolution @@ -70,14 +70,14 @@ solvers "(h|Yi|O2|N2|H2O)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-6; relTol 0; } hFinal { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-6; relTol 0; diff --git a/tutorials/lagrangian/reactingParcelFilmFoam/hotBoxes/system/fvSolution b/tutorials/lagrangian/reactingParcelFilmFoam/hotBoxes/system/fvSolution index 59523edc89a5d1fa03175edde399d987c80aa132..cb05b065e84721778392d0ccea2643c937821076 100644 --- a/tutorials/lagrangian/reactingParcelFilmFoam/hotBoxes/system/fvSolution +++ b/tutorials/lagrangian/reactingParcelFilmFoam/hotBoxes/system/fvSolution @@ -34,7 +34,7 @@ solvers "(U|h|k|epsilon)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-5; relTol 0.1; diff --git a/tutorials/lagrangian/reactingParcelFilmFoam/splashPanel/system/fvSolution b/tutorials/lagrangian/reactingParcelFilmFoam/splashPanel/system/fvSolution index d10b1fda26f3cb22d7dbacab7fbea3cdc791532f..3bb49837a4f042e7a439653b8f0795611773ad72 100644 --- a/tutorials/lagrangian/reactingParcelFilmFoam/splashPanel/system/fvSolution +++ b/tutorials/lagrangian/reactingParcelFilmFoam/splashPanel/system/fvSolution @@ -66,14 +66,14 @@ solvers "(h|Yi|O2|N2|H2O)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-6; relTol 0; } hFinal { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-6; relTol 0; diff --git a/tutorials/lagrangian/reactingParcelFilmFoam/splashPanel/system/wallFilmRegion/fvSolution b/tutorials/lagrangian/reactingParcelFilmFoam/splashPanel/system/wallFilmRegion/fvSolution index ae07f4217a40b8a64f246a96602e83082bb664e7..0aab30d80a7ed3068b7e952f277ea2cacb279e59 100644 --- a/tutorials/lagrangian/reactingParcelFilmFoam/splashPanel/system/wallFilmRegion/fvSolution +++ b/tutorials/lagrangian/reactingParcelFilmFoam/splashPanel/system/wallFilmRegion/fvSolution @@ -19,14 +19,14 @@ solvers { "(Uf|hf|deltaf\*rhof)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-10; relTol 0; } deltaf { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-10; relTol 0; diff --git a/tutorials/lagrangian/reactingParcelFoam/filter/system/fvSolution b/tutorials/lagrangian/reactingParcelFoam/filter/system/fvSolution index 7248e662be8c6f22f0064da29d1e254ab59c2cbb..1914fbdcbb42be9e4d08bb9b5dd6aac17fbae639 100644 --- a/tutorials/lagrangian/reactingParcelFoam/filter/system/fvSolution +++ b/tutorials/lagrangian/reactingParcelFoam/filter/system/fvSolution @@ -64,7 +64,7 @@ solvers "(Yi|O2|N2|H2O)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-6; relTol 0; diff --git a/tutorials/lagrangian/reactingParcelFoam/parcelInBox/system/fvSolution b/tutorials/lagrangian/reactingParcelFoam/parcelInBox/system/fvSolution index 6d07e4ca22c052ddac5486f062a6606ddbcf2bee..b2d5ddbe448697fa1ea9e1f5a6814291b34f6c50 100644 --- a/tutorials/lagrangian/reactingParcelFoam/parcelInBox/system/fvSolution +++ b/tutorials/lagrangian/reactingParcelFoam/parcelInBox/system/fvSolution @@ -65,7 +65,7 @@ solvers "(Yi|O2|N2|H2O)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-6; relTol 0; diff --git a/tutorials/lagrangian/reactingParcelFoam/verticalChannel/system/fvSolution b/tutorials/lagrangian/reactingParcelFoam/verticalChannel/system/fvSolution index 404226dba78c1a8ac95ec873a41709a9c1bb8cd5..71738b6dc9bfb4d0cd0635402293a1a326f7c6be 100644 --- a/tutorials/lagrangian/reactingParcelFoam/verticalChannel/system/fvSolution +++ b/tutorials/lagrangian/reactingParcelFoam/verticalChannel/system/fvSolution @@ -69,7 +69,7 @@ solvers "(Yi|O2|N2|H2O)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-6; relTol 0; diff --git a/tutorials/lagrangian/sprayFoam/aachenBomb/system/fvSolution b/tutorials/lagrangian/sprayFoam/aachenBomb/system/fvSolution index a9aed26fba900d0731d4c731351070856b20e3a7..86340f5250969047d45dedbdad3a49715475e23e 100644 --- a/tutorials/lagrangian/sprayFoam/aachenBomb/system/fvSolution +++ b/tutorials/lagrangian/sprayFoam/aachenBomb/system/fvSolution @@ -64,7 +64,7 @@ solvers "(Yi|O2|N2|H2O)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-6; relTol 0; diff --git a/tutorials/mesh/moveDynamicMesh/SnakeRiverCanyon/system/fvSolution b/tutorials/mesh/moveDynamicMesh/SnakeRiverCanyon/system/fvSolution index 028c7d8818b348ca6bbca00eba798c1ce944661c..236f0297a11674af3e2002c5d71e15c14b525d7e 100644 --- a/tutorials/mesh/moveDynamicMesh/SnakeRiverCanyon/system/fvSolution +++ b/tutorials/mesh/moveDynamicMesh/SnakeRiverCanyon/system/fvSolution @@ -27,7 +27,7 @@ solvers U { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; tolerance 1e-05; relTol 0; diff --git a/tutorials/multiphase/driftFluxFoam/ras/dahl/system/fvSolution b/tutorials/multiphase/driftFluxFoam/ras/dahl/system/fvSolution index 5cedaa31cf6f2f4cce296e59311d696a6bb95f4e..820077fa967a8bd503ade3f0531ddf778f03cc13 100644 --- a/tutorials/multiphase/driftFluxFoam/ras/dahl/system/fvSolution +++ b/tutorials/multiphase/driftFluxFoam/ras/dahl/system/fvSolution @@ -59,9 +59,8 @@ solvers "(U|k|epsilon)" { - solver PBiCG; + solver PBiCGStab; preconditioner DILU; - smoother symGaussSeidel; tolerance 1e-7; relTol 0.1; minIter 1; diff --git a/tutorials/multiphase/driftFluxFoam/ras/mixerVessel2D/system/fvSolution b/tutorials/multiphase/driftFluxFoam/ras/mixerVessel2D/system/fvSolution index d2397321ec49d9c04401f497e4e58dda9b741f84..23a99325b0082178d93bb8ed36d23303b516980a 100644 --- a/tutorials/multiphase/driftFluxFoam/ras/mixerVessel2D/system/fvSolution +++ b/tutorials/multiphase/driftFluxFoam/ras/mixerVessel2D/system/fvSolution @@ -67,8 +67,8 @@ solvers "(U|k|epsilon)" { - solver smoothSolver; - smoother symGaussSeidel; + solver PBiCGStab; + preconditioner DILU; tolerance 1e-6; relTol 0.1; } diff --git a/tutorials/multiphase/driftFluxFoam/ras/tank3D/system/fvSolution b/tutorials/multiphase/driftFluxFoam/ras/tank3D/system/fvSolution index ed4bd47af88c83a9e16b5b618c5a7807702370a8..9c3d51f7c5bd45b224626efccdb496a5f09532be 100644 --- a/tutorials/multiphase/driftFluxFoam/ras/tank3D/system/fvSolution +++ b/tutorials/multiphase/driftFluxFoam/ras/tank3D/system/fvSolution @@ -59,8 +59,8 @@ solvers "(U|k|epsilon)" { - solver smoothSolver; - smoother symGaussSeidel; + solver PBiCGStab; + preconditioner DILU; tolerance 1e-7; relTol 0.1; minIter 1;