Loading...
Searching...
No Matches
PBiCG.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2019-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "PBiCG.H"
30#include "PrecisionAdaptor.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
37
38 lduMatrix::solver::addasymMatrixConstructorToTable<PBiCG>
40}
41
42
43// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44
45Foam::PBiCG::PBiCG
46(
47 const word& fieldName,
48 const lduMatrix& matrix,
49 const FieldField<Field, scalar>& interfaceBouCoeffs,
50 const FieldField<Field, scalar>& interfaceIntCoeffs,
51 const lduInterfaceFieldPtrsList& interfaces,
52 const dictionary& solverControls
53)
54:
56 (
57 fieldName,
58 matrix,
59 interfaceBouCoeffs,
60 interfaceIntCoeffs,
61 interfaces,
62 solverControls
63 )
64{}
65
66
67// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
68
70(
71 scalarField& psi_s,
72 const scalarField& source,
73 const direction cmpt
74) const
75{
77 solveScalarField& psi = tpsi.ref();
78
79 // --- Setup class containing solver performance data
80 solverPerformance solverPerf
81 (
83 fieldName_
84 );
85
86 const label nCells = psi.size();
87
88 solveScalar* __restrict__ psiPtr = psi.begin();
89
90 solveScalarField pA(nCells);
91 solveScalar* __restrict__ pAPtr = pA.begin();
92
93 solveScalarField wA(nCells);
94 solveScalar* __restrict__ wAPtr = wA.begin();
95
96 // --- Calculate A.psi
97 matrix_.Amul(wA, psi, interfaceBouCoeffs_, interfaces_, cmpt);
98
99 // --- Calculate initial residual field
101 solveScalarField rA(tsource() - wA);
102 solveScalar* __restrict__ rAPtr = rA.begin();
103
104 matrix().setResidualField
105 (
107 fieldName_,
108 true
109 );
110
111 // --- Calculate normalisation factor
112 const solveScalar normFactor = this->normFactor(psi, tsource(), wA, pA);
113
114 if ((log_ >= 2) || (lduMatrix::debug >= 2))
115 {
116 Info<< " Normalisation factor = " << normFactor << endl;
117 }
118
119 // --- Calculate normalised residual norm
120 solverPerf.initialResidual() =
121 gSumMag(rA, matrix().mesh().comm())
122 /normFactor;
123 solverPerf.finalResidual() = solverPerf.initialResidual();
124
125 // --- Check convergence, solve if not converged
126 if
127 (
128 minIter_ > 0
129 || !solverPerf.checkConvergence(tolerance_, relTol_, log_)
130 )
131 {
132 solveScalarField pT(nCells, 0);
133 solveScalar* __restrict__ pTPtr = pT.begin();
134
135 solveScalarField wT(nCells);
136 solveScalar* __restrict__ wTPtr = wT.begin();
137
138 // --- Calculate T.psi
139 matrix_.Tmul(wT, psi, interfaceIntCoeffs_, interfaces_, cmpt);
140
141 // --- Calculate initial transpose residual field
142 solveScalarField rT(tsource() - wT);
143 solveScalar* __restrict__ rTPtr = rT.begin();
144
145 // --- Initial value not used
146 solveScalar wArT = 0;
147
148 // --- Select and construct the preconditioner
149 if (!preconPtr_)
150 {
152 (
153 *this,
154 controlDict_
155 );
156 }
157
158 // --- Solver iteration
159 do
160 {
161 // --- Store previous wArT
162 const solveScalar wArTold = wArT;
163
164 // --- Precondition residuals
165 preconPtr_->precondition(wA, rA, cmpt);
166 preconPtr_->preconditionT(wT, rT, cmpt);
167
168 // --- Update search directions:
169 wArT = gSumProd(wA, rT, matrix().mesh().comm());
170
171 if (solverPerf.nIterations() == 0)
172 {
173 for (label cell=0; cell<nCells; cell++)
174 {
175 pAPtr[cell] = wAPtr[cell];
176 pTPtr[cell] = wTPtr[cell];
177 }
178 }
179 else
180 {
181 const solveScalar beta = wArT/wArTold;
182
183 for (label cell=0; cell<nCells; cell++)
184 {
185 pAPtr[cell] = wAPtr[cell] + beta*pAPtr[cell];
186 pTPtr[cell] = wTPtr[cell] + beta*pTPtr[cell];
187 }
188 }
189
190
191 // --- Update preconditioned residuals
192 matrix_.Amul(wA, pA, interfaceBouCoeffs_, interfaces_, cmpt);
193 matrix_.Tmul(wT, pT, interfaceIntCoeffs_, interfaces_, cmpt);
194
195 const solveScalar wApT = gSumProd(wA, pT, matrix().mesh().comm());
196
197 // --- Test for singularity
198 if (solverPerf.checkSingularity(mag(wApT)/normFactor))
199 {
200 break;
201 }
202
203
204 // --- Update solution and residual:
205
206 const solveScalar alpha = wArT/wApT;
207
208 for (label cell=0; cell<nCells; cell++)
209 {
210 psiPtr[cell] += alpha*pAPtr[cell];
211 rAPtr[cell] -= alpha*wAPtr[cell];
212 rTPtr[cell] -= alpha*wTPtr[cell];
213 }
214
215 solverPerf.finalResidual() =
216 gSumMag(rA, matrix().mesh().comm())
217 /normFactor;
218 } while
219 (
220 (
221 ++solverPerf.nIterations() < maxIter_
222 && !solverPerf.checkConvergence(tolerance_, relTol_, log_)
223 )
224 || solverPerf.nIterations() < minIter_
225 );
226 }
227
228 // Recommend PBiCGStab if PBiCG fails to converge
229 const label upperMaxIters = max(maxIter_, lduMatrix::defaultMaxIter);
230
231 if (solverPerf.nIterations() > upperMaxIters)
232 {
234 << "PBiCG has failed to converge within the maximum iterations: "
235 << upperMaxIters << nl
236 << " Please try the more robust PBiCGStab solver."
237 << exit(FatalError);
238 }
239
240 if (preconPtr_)
241 {
242 preconPtr_->setFinished(solverPerf);
243 }
244
245 matrix().setResidualField
246 (
248 fieldName_,
249 false
250 );
251
252 return solverPerf;
253}
254
255
256// ************************************************************************* //
A const Field/List wrapper with possible data conversion.
A field of fields is a PtrList of fields with reference counting.
Definition FieldField.H:77
Preconditioned bi-conjugate gradient solver for asymmetric lduMatrices using a run-time selectable pr...
Definition PBiCG.H:51
virtual solverPerformance solve(scalarField &psi, const scalarField &source, const direction cmpt=0) const
Solve the matrix with this solver.
Definition PBiCG.C:63
A non-const Field/List wrapper with possible data conversion.
const Type & finalResidual() const noexcept
Return final residual.
const labelType & nIterations() const noexcept
Return number of iterations.
bool checkSingularity(const Type &residual)
Singularity test.
const Type & initialResidual() const noexcept
Return initial residual.
bool checkConvergence(const Type &tolerance, const Type &relTolerance, const int logLevel=0)
Check, store and return convergence.
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition UListI.H:410
A cell is defined as a list of faces with extra functionality.
Definition cell.H:56
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
static autoPtr< preconditioner > New(const solver &sol, const dictionary &solverControls)
Return a new preconditioner.
static word getName(const dictionary &)
Find the preconditioner name (directly or from a sub-dictionary).
const FieldField< Field, scalar > & interfaceIntCoeffs() const noexcept
Definition lduMatrix.H:338
const lduMatrix & matrix_
Definition lduMatrix.H:158
label maxIter_
Maximum number of iterations in the solver.
Definition lduMatrix.H:181
solver(const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces, const dictionary &solverControls)
Construct solver for given field name, matrix etc.
scalar tolerance_
Final convergence tolerance.
Definition lduMatrix.H:191
const lduInterfaceFieldPtrsList & interfaces() const noexcept
Definition lduMatrix.H:343
lduInterfaceFieldPtrsList interfaces_
Definition lduMatrix.H:161
const FieldField< Field, scalar > & interfaceBouCoeffs_
Definition lduMatrix.H:159
label minIter_
Minimum number of iterations in the solver.
Definition lduMatrix.H:176
const lduMatrix & matrix() const noexcept
Definition lduMatrix.H:328
solveScalarField::cmptType normFactor(const solveScalarField &psi, const solveScalarField &source, const solveScalarField &Apsi, solveScalarField &tmpField, const lduMatrix::normTypes normType) const
Return the matrix norm using the specified norm method.
int log_
Verbosity level for solver output statements.
Definition lduMatrix.H:171
scalar relTol_
Convergence tolerance relative to the initial.
Definition lduMatrix.H:196
const FieldField< Field, scalar > & interfaceIntCoeffs_
Definition lduMatrix.H:160
const FieldField< Field, scalar > & interfaceBouCoeffs() const noexcept
Definition lduMatrix.H:333
dictionary controlDict_
Dictionary of solution controls.
Definition lduMatrix.H:166
const word & fieldName() const noexcept
Definition lduMatrix.H:323
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition lduMatrix.H:81
static constexpr const label defaultMaxIter
Default maximum number of iterations for solvers (1000).
Definition lduMatrix.H:139
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition refPtrI.H:230
Base solver class.
Definition solver.H:48
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
const volScalarField & psi
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Namespace for OpenFOAM.
scalarProduct< Type, Type >::type gSumProd(const UList< Type > &f1, const UList< Type > &f2, const label comm)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< solveScalar > solveScalarField
uint8_t direction
Definition direction.H:49
typeOfMag< Type >::type gSumMag(const FieldField< Field, Type > &f)
lduMatrix::solver::addasymMatrixConstructorToTable< PBiCG > addPBiCGAsymMatrixConstructorToTable_
Definition PBiCG.C:32
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UPtrList< const lduInterfaceField > lduInterfaceFieldPtrsList
List of coupled interface fields to be used in coupling.
SolverPerformance< scalar > solverPerformance
SolverPerformance instantiated for a scalar.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & alpha
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)