/*---------------------------------------------------------------------------*\ ========= | \\ / 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-2021 OpenCFD Ltd. Copyright (C) 2023 Huawei (Yu Ankun) ------------------------------------------------------------------------------- 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 "FPCG.H" #include "PrecisionAdaptor.H" // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * // namespace Foam { defineTypeNameAndDebug(FPCG, 0); lduMatrix::solver::addsymMatrixConstructorToTable addFPCGSymMatrixConstructorToTable_; } // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // Foam::FPCG::FPCG ( 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 ) {} // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * // void Foam::FPCG::gSumMagProd ( FixedList& globalSum, const solveScalarField& a, const solveScalarField& b, const label comm ) { const label nCells = a.size(); globalSum = 0.0; for (label cell=0; cell(), UPstream::msgType(), // (ignored): direct MPI call comm ); } // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // Foam::solverPerformance Foam::FPCG::scalarSolve ( solveScalarField& psi, const solveScalarField& source, const direction cmpt ) const { // --- Setup class containing solver performance data solverPerformance solverPerf ( lduMatrix::preconditioner::getName(controlDict_) + typeName, fieldName_ ); 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(); solveScalar wArA = solverPerf.great_; solveScalar wArAold = wArA; // --- Calculate A.psi matrix_.Amul(wA, psi, interfaceBouCoeffs_, interfaces_, cmpt); // --- Calculate initial residual field solveScalarField rA(source - wA); solveScalar* __restrict__ rAPtr = rA.begin(); matrix().setResidualField ( ConstPrecisionAdaptor(rA)(), fieldName_, true ); // --- Calculate normalisation factor solveScalar normFactor = this->normFactor(psi, source, 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_) ) { // --- Select and construct the preconditioner if (!preconPtr_) { preconPtr_ = lduMatrix::preconditioner::New ( *this, controlDict_ ); } FixedList globalSum; // --- Solver iteration do { // --- Store previous wArA wArAold = wArA; // --- Precondition residual preconPtr_->precondition(wA, rA, cmpt); // --- Update search directions and calculate residual: gSumMagProd(globalSum, wA, rA, matrix().mesh().comm()); wArA = globalSum[0]; solverPerf.finalResidual() = globalSum[1]/normFactor; // Check convergence (bypass if not enough iterations yet) if ( (minIter_ <= 0 || solverPerf.nIterations() >= minIter_) && solverPerf.checkConvergence(tolerance_, relTol_, log_) ) { break; } if (solverPerf.nIterations() == 0) { for (label cell=0; cellsetFinished(solverPerf); } matrix().setResidualField ( ConstPrecisionAdaptor(rA)(), fieldName_, false ); return solverPerf; } Foam::solverPerformance Foam::FPCG::solve ( scalarField& psi_s, const scalarField& source, const direction cmpt ) const { PrecisionAdaptor tpsi(psi_s); return scalarSolve ( tpsi.ref(), ConstPrecisionAdaptor(source)(), cmpt ); } // ************************************************************************* //