/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2016-2017 OpenFOAM Foundation
Copyright (C) 2019-2021 OpenCFD Ltd.
-------------------------------------------------------------------------------
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 .
\*---------------------------------------------------------------------------*/
#include "PBiCGStab.H"
#include "PrecisionAdaptor.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
defineTypeNameAndDebug(PBiCGStab, 0);
lduMatrix::solver::addsymMatrixConstructorToTable
addPBiCGStabSymMatrixConstructorToTable_;
lduMatrix::solver::addasymMatrixConstructorToTable
addPBiCGStabAsymMatrixConstructorToTable_;
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::PBiCGStab::PBiCGStab
(
const word& fieldName,
const lduMatrix& matrix,
const FieldField& interfaceBouCoeffs,
const FieldField& interfaceIntCoeffs,
const lduInterfaceFieldPtrsList& interfaces,
const dictionary& solverControls
)
:
lduMatrix::solver
(
fieldName,
matrix,
interfaceBouCoeffs,
interfaceIntCoeffs,
interfaces,
solverControls
)
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
Foam::solverPerformance Foam::PBiCGStab::scalarSolve
(
solveScalarField& psi,
const solveScalarField& source,
const direction cmpt
) const
{
// --- Setup class containing solver performance data
solverPerformance solverPerf
(
lduMatrix::preconditioner::getName(controlDict_) + typeName,
fieldName_
);
const label nCells = psi.size();
solveScalar* __restrict__ psiPtr = psi.begin();
solveScalarField pA(nCells);
solveScalar* __restrict__ pAPtr = pA.begin();
solveScalarField yA(nCells);
solveScalar* __restrict__ yAPtr = yA.begin();
// --- Calculate A.psi
matrix_.Amul(yA, psi, interfaceBouCoeffs_, interfaces_, cmpt);
// --- Calculate initial residual field
solveScalarField rA(source - yA);
solveScalar* __restrict__ rAPtr = rA.begin();
matrix().setResidualField
(
ConstPrecisionAdaptor(rA)(),
fieldName_,
true
);
// --- Calculate normalisation factor
const solveScalar normFactor = this->normFactor(psi, source, yA, pA);
if ((log_ >= 2) || (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_, log_)
)
{
solveScalarField AyA(nCells);
solveScalar* __restrict__ AyAPtr = AyA.begin();
solveScalarField sA(nCells);
solveScalar* __restrict__ sAPtr = sA.begin();
solveScalarField zA(nCells);
solveScalar* __restrict__ zAPtr = zA.begin();
solveScalarField tA(nCells);
solveScalar* __restrict__ tAPtr = tA.begin();
// --- Store initial residual
const solveScalarField rA0(rA);
// --- Initial values not used
solveScalar rA0rA = 0;
solveScalar alpha = 0;
solveScalar omega = 0;
// --- Select and construct the preconditioner
if (!preconPtr_)
{
preconPtr_ = lduMatrix::preconditioner::New
(
*this,
controlDict_
);
}
// --- Solver iteration
do
{
// --- Store previous rA0rA
const solveScalar 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; cellprecondition(yA, pA, cmpt);
// --- Calculate AyA
matrix_.Amul(AyA, yA, interfaceBouCoeffs_, interfaces_, cmpt);
const solveScalar rA0AyA =
gSumProd(rA0, AyA, matrix().mesh().comm());
alpha = rA0rA/rA0AyA;
// --- Calculate sA
for (label cell=0; cell= minIter_
&& solverPerf.checkConvergence(tolerance_, relTol_, log_)
)
{
for (label cell=0; cellprecondition(zA, sA, cmpt);
// --- Calculate tA
matrix_.Amul(tA, zA, interfaceBouCoeffs_, interfaces_, cmpt);
const solveScalar tAtA = gSumSqr(tA, matrix().mesh().comm());
// --- Calculate omega from tA and sA
// (cheaper than using zA with preconditioned tA)
omega = gSumProd(tA, sA, matrix().mesh().comm())/tAtA;
// --- Update solution and residual
for (label cell=0; cellsetFinished(solverPerf);
}
matrix().setResidualField
(
ConstPrecisionAdaptor(rA)(),
fieldName_,
false
);
return solverPerf;
}
Foam::solverPerformance Foam::PBiCGStab::solve
(
scalarField& psi_s,
const scalarField& source,
const direction cmpt
) const
{
PrecisionAdaptor tpsi(psi_s);
return scalarSolve
(
tpsi.ref(),
ConstPrecisionAdaptor(source)(),
cmpt
);
}
// ************************************************************************* //