1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
11 OpenFOAM is free software: you can redistribute it and/or modify it
12 under the terms of the GNU General Public License as published by
13 the Free Software Foundation, either version 3 of the License, or
14 (at your option) any later version.
16 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
28 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 defineTypeNameAndDebug(PCG, 0);
34 lduMatrix::solver::addsymMatrixConstructorToTable<PCG>
35 addPCGSymMatrixConstructorToTable_;
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 const word& fieldName,
44 const lduMatrix& matrix,
45 const FieldField<Field, scalar>& interfaceBouCoeffs,
46 const FieldField<Field, scalar>& interfaceIntCoeffs,
47 const lduInterfaceFieldPtrsList& interfaces,
48 const dictionary& solverControls
63 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
65 Foam::lduMatrix::solverPerformance Foam::PCG::solve
68 const scalarField& source,
72 // --- Setup class containing solver performance data
73 lduMatrix::solverPerformance solverPerf
75 lduMatrix::preconditioner::getName(controlDict_) + typeName,
79 register label nCells = psi.size();
81 scalar* __restrict__ psiPtr = psi.begin();
83 scalarField pA(nCells);
84 scalar* __restrict__ pAPtr = pA.begin();
86 scalarField wA(nCells);
87 scalar* __restrict__ wAPtr = wA.begin();
89 scalar wArA = matrix_.great_;
90 scalar wArAold = wArA;
92 // --- Calculate A.psi
93 matrix_.Amul(wA, psi, interfaceBouCoeffs_, interfaces_, cmpt);
95 // --- Calculate initial residual field
96 scalarField rA(source - wA);
97 scalar* __restrict__ rAPtr = rA.begin();
99 // --- Calculate normalisation factor
100 scalar normFactor = this->normFactor(psi, source, wA, pA);
102 if (lduMatrix::debug >= 2)
104 Info<< " Normalisation factor = " << normFactor << endl;
107 // --- Calculate normalised residual norm
108 solverPerf.initialResidual() = gSumMag(rA)/normFactor;
109 solverPerf.finalResidual() = solverPerf.initialResidual();
111 // --- Check convergence, solve if not converged
112 if (!solverPerf.checkConvergence(tolerance_, relTol_))
114 // --- Select and construct the preconditioner
115 autoPtr<lduMatrix::preconditioner> preconPtr =
116 lduMatrix::preconditioner::New
122 // --- Solver iteration
125 // --- Store previous wArA
128 // --- Precondition residual
129 preconPtr->precondition(wA, rA, cmpt);
131 // --- Update search directions:
132 wArA = gSumProd(wA, rA);
134 if (solverPerf.nIterations() == 0)
136 for (register label cell=0; cell<nCells; cell++)
138 pAPtr[cell] = wAPtr[cell];
143 scalar beta = wArA/wArAold;
145 for (register label cell=0; cell<nCells; cell++)
147 pAPtr[cell] = wAPtr[cell] + beta*pAPtr[cell];
152 // --- Update preconditioned residual
153 matrix_.Amul(wA, pA, interfaceBouCoeffs_, interfaces_, cmpt);
155 scalar wApA = gSumProd(wA, pA);
158 // --- Test for singularity
159 if (solverPerf.checkSingularity(mag(wApA)/normFactor)) break;
162 // --- Update solution and residual:
164 scalar alpha = wArA/wApA;
166 for (register label cell=0; cell<nCells; cell++)
168 psiPtr[cell] += alpha*pAPtr[cell];
169 rAPtr[cell] -= alpha*wAPtr[cell];
172 solverPerf.finalResidual() = gSumMag(rA)/normFactor;
176 solverPerf.nIterations()++ < maxIter_
177 && !(solverPerf.checkConvergence(tolerance_, relTol_))
185 // ************************************************************************* //