/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2011-2017 OpenFOAM Foundation
Copyright (C) 2019-2022 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 "PBiCG.H"
#include "PrecisionAdaptor.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
defineTypeNameAndDebug(PBiCG, 0);
lduMatrix::solver::addasymMatrixConstructorToTable
addPBiCGAsymMatrixConstructorToTable_;
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::PBiCG::PBiCG
(
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::PBiCG::solve
(
scalarField& psi_s,
const scalarField& source,
const direction cmpt
) const
{
PrecisionAdaptor tpsi(psi_s);
solveScalarField& psi = tpsi.ref();
// --- 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 wA(nCells);
solveScalar* __restrict__ wAPtr = wA.begin();
// --- Calculate A.psi
matrix_.Amul(wA, psi, interfaceBouCoeffs_, interfaces_, cmpt);
// --- Calculate initial residual field
ConstPrecisionAdaptor tsource(source);
solveScalarField rA(tsource() - wA);
solveScalar* __restrict__ rAPtr = rA.begin();
matrix().setResidualField
(
ConstPrecisionAdaptor(rA)(),
fieldName_,
true
);
// --- Calculate normalisation factor
const solveScalar normFactor = this->normFactor(psi, tsource(), wA, 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 pT(nCells, 0);
solveScalar* __restrict__ pTPtr = pT.begin();
solveScalarField wT(nCells);
solveScalar* __restrict__ wTPtr = wT.begin();
// --- Calculate T.psi
matrix_.Tmul(wT, psi, interfaceIntCoeffs_, interfaces_, cmpt);
// --- Calculate initial transpose residual field
solveScalarField rT(tsource() - wT);
solveScalar* __restrict__ rTPtr = rT.begin();
// --- Initial value not used
solveScalar wArT = 0;
// --- Select and construct the preconditioner
if (!preconPtr_)
{
preconPtr_ = lduMatrix::preconditioner::New
(
*this,
controlDict_
);
}
// --- Solver iteration
do
{
// --- Store previous wArT
const solveScalar wArTold = wArT;
// --- Precondition residuals
preconPtr_->precondition(wA, rA, cmpt);
preconPtr_->preconditionT(wT, rT, cmpt);
// --- Update search directions:
wArT = gSumProd(wA, rT, matrix().mesh().comm());
if (solverPerf.nIterations() == 0)
{
for (label cell=0; cell upperMaxIters)
{
FatalErrorInFunction
<< "PBiCG has failed to converge within the maximum iterations: "
<< upperMaxIters << nl
<< " Please try the more robust PBiCGStab solver."
<< exit(FatalError);
}
if (preconPtr_)
{
preconPtr_->setFinished(solverPerf);
}
matrix().setResidualField
(
ConstPrecisionAdaptor(rA)(),
fieldName_,
false
);
return solverPerf;
}
// ************************************************************************* //