Loading...
Searching...
No Matches
FPCG.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-2021 OpenCFD Ltd.
10 Copyright (C) 2023 Huawei (Yu Ankun)
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
30#include "FPCG.H"
31#include "PrecisionAdaptor.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
38
39 lduMatrix::solver::addsymMatrixConstructorToTable<FPCG>
41}
42
43
44// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45
46Foam::FPCG::FPCG
47(
48 const word& fieldName,
49 const lduMatrix& matrix,
50 const FieldField<Field, scalar>& interfaceBouCoeffs,
51 const FieldField<Field, scalar>& interfaceIntCoeffs,
52 const lduInterfaceFieldPtrsList& interfaces,
53 const dictionary& solverControls
54)
55:
57 (
58 fieldName,
59 matrix,
60 interfaceBouCoeffs,
61 interfaceIntCoeffs,
62 interfaces,
63 solverControls
64 )
65{}
66
67
68// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
69
70void Foam::FPCG::gSumMagProd
71(
72 FixedList<solveScalar, 2>& globalSum,
73 const solveScalarField& a,
74 const solveScalarField& b,
75 const label comm
76)
77{
78 const label nCells = a.size();
79
80 globalSum = 0.0;
81 for (label cell=0; cell<nCells; ++cell)
82 {
83 globalSum[0] += a[cell]*b[cell]; // sumProd(a, b)
84 globalSum[1] += mag(b[cell]); // sumMag(b)
85 }
86
88 (
89 globalSum.data(),
90 globalSum.size(),
92 UPstream::msgType(), // (ignored): direct MPI call
93 comm
94 );
95}
96
97// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
98
100(
102 const solveScalarField& source,
103 const direction cmpt
104) const
105{
106 // --- Setup class containing solver performance data
107 solverPerformance solverPerf
108 (
110 fieldName_
111 );
112
113 label nCells = psi.size();
114
115 solveScalar* __restrict__ psiPtr = psi.begin();
116
117 solveScalarField pA(nCells);
118 solveScalar* __restrict__ pAPtr = pA.begin();
119
120 solveScalarField wA(nCells);
121 solveScalar* __restrict__ wAPtr = wA.begin();
122
123 solveScalar wArA = solverPerf.great_;
124 solveScalar wArAold = wArA;
125
126 // --- Calculate A.psi
127 matrix_.Amul(wA, psi, interfaceBouCoeffs_, interfaces_, cmpt);
128
129 // --- Calculate initial residual field
130 solveScalarField rA(source - wA);
131 solveScalar* __restrict__ rAPtr = rA.begin();
132
133 matrix().setResidualField
134 (
136 fieldName_,
137 true
138 );
139
140 // --- Calculate normalisation factor
141 solveScalar normFactor = this->normFactor(psi, source, wA, pA);
142
143 if ((log_ >= 2) || (lduMatrix::debug >= 2))
144 {
145 Info<< " Normalisation factor = " << normFactor << endl;
146 }
147
148 // --- Calculate normalised residual norm
149 solverPerf.initialResidual() =
150 gSumMag(rA, matrix().mesh().comm())
151 /normFactor;
152 solverPerf.finalResidual() = solverPerf.initialResidual();
153
154 // --- Check convergence, solve if not converged
155 if
156 (
157 minIter_ > 0
158 || !solverPerf.checkConvergence(tolerance_, relTol_, log_)
159 )
160 {
161 // --- Select and construct the preconditioner
162 if (!preconPtr_)
163 {
165 (
166 *this,
167 controlDict_
168 );
169 }
170
172
173 // --- Solver iteration
174 do
175 {
176 // --- Store previous wArA
177 wArAold = wArA;
178
179 // --- Precondition residual
180 preconPtr_->precondition(wA, rA, cmpt);
181
182 // --- Update search directions and calculate residual:
183 gSumMagProd(globalSum, wA, rA, matrix().mesh().comm());
184
185 wArA = globalSum[0];
186
187 solverPerf.finalResidual() = globalSum[1]/normFactor;
188
189 // Check convergence (bypass if not enough iterations yet)
190 if
191 (
192 (minIter_ <= 0 || solverPerf.nIterations() >= minIter_)
193 && solverPerf.checkConvergence(tolerance_, relTol_, log_)
194 )
195 {
196 break;
197 }
198
199 if (solverPerf.nIterations() == 0)
200 {
201 for (label cell=0; cell<nCells; cell++)
202 {
203 pAPtr[cell] = wAPtr[cell];
204 }
205 }
206 else
207 {
208 const solveScalar beta = wArA/wArAold;
209
210 for (label cell=0; cell<nCells; cell++)
211 {
212 pAPtr[cell] = wAPtr[cell] + beta*pAPtr[cell];
213 }
214 }
215
216
217 // --- Update preconditioned residual
218 matrix_.Amul(wA, pA, interfaceBouCoeffs_, interfaces_, cmpt);
219
220 solveScalar wApA = gSumProd(wA, pA, matrix().mesh().comm());
221
222 // --- Test for singularity
223 if (solverPerf.checkSingularity(mag(wApA)/normFactor)) break;
224
225
226 // --- Update solution and residual:
227
228 const solveScalar alpha = wArA/wApA;
229
230 for (label cell=0; cell<nCells; cell++)
231 {
232 psiPtr[cell] += alpha*pAPtr[cell];
233 rAPtr[cell] -= alpha*wAPtr[cell];
234 }
235
236 } while
237 (
238 ++solverPerf.nIterations() < maxIter_
239 );
240 }
241
242 if (preconPtr_)
243 {
244 preconPtr_->setFinished(solverPerf);
245 }
246
247 matrix().setResidualField
248 (
250 fieldName_,
251 false
252 );
253
254 return solverPerf;
255}
256
257
258
260(
261 scalarField& psi_s,
262 const scalarField& source,
263 const direction cmpt
264) const
265{
267 return scalarSolve
268 (
269 tpsi.ref(),
271 cmpt
272 );
273}
274
275
276// ************************************************************************* //
A const Field/List wrapper with possible data conversion.
label size() const noexcept
The number of elements in list.
Definition DLListBase.H:194
A 'faster' preconditioned conjugate gradient solver for symmetric lduMatrices using a run-time select...
Definition FPCG.H:58
virtual solverPerformance scalarSolve(solveScalarField &psi, const solveScalarField &source, const direction cmpt=0) const
Solve the matrix with this solver.
Definition FPCG.C:93
virtual solverPerformance solve(scalarField &psi, const scalarField &source, const direction cmpt=0) const
Solve the matrix with this solver.
Definition FPCG.C:253
A field of fields is a PtrList of fields with reference counting.
Definition FieldField.H:77
A 1D vector of objects of type <T> with a fixed length <N>.
Definition FixedList.H:73
static constexpr label size() noexcept
Return the number of elements in the FixedList.
Definition FixedList.H:619
T * data() noexcept
Return pointer to the underlying array serving as data storage.
Definition FixedListI.H:121
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.
static const scalar great_
Large Type for the use in solvers.
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
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static int & msgType() noexcept
Message tag of standard messages.
Definition UPstream.H:1926
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 > & 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
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
Namespace for OpenFOAM.
scalarProduct< Type, Type >::type gSumProd(const UList< Type > &f1, const UList< Type > &f2, const label comm)
lduMatrix::solver::addsymMatrixConstructorToTable< FPCG > addFPCGSymMatrixConstructorToTable_
Definition FPCG.C:33
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)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
Field< solveScalar > solveScalarField
uint8_t direction
Definition direction.H:49
typeOfMag< Type >::type gSumMag(const FieldField< Field, Type > &f)
UPtrList< const lduInterfaceField > lduInterfaceFieldPtrsList
List of coupled interface fields to be used in coupling.
SolverPerformance< scalar > solverPerformance
SolverPerformance instantiated for a scalar.
volScalarField & alpha
volScalarField & b
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)