PCICG.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-2016 OpenFOAM Foundation
9 Copyright (C) 2021 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 "PCICG.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class Type, class DType, class LUType>
35(
36 const word& fieldName,
38 const dictionary& solverDict
39)
40:
41 LduMatrix<Type, DType, LUType>::solver
42 (
43 fieldName,
44 matrix,
45 solverDict
46 )
47{}
48
49
50// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
51
52template<class Type, class DType, class LUType>
55{
56 const word preconditionerName(this->controlDict_.getWord("preconditioner"));
57
58 // --- Setup class containing solver performance data
60 (
61 preconditionerName + typeName,
62 this->fieldName_
63 );
64
65 label nIter = 0;
66
67 label nCells = psi.size();
68
69 Type* __restrict__ psiPtr = psi.begin();
70
71 Field<Type> pA(nCells);
72 Type* __restrict__ pAPtr = pA.begin();
73
74 Field<Type> wA(nCells);
75 Type* __restrict__ wAPtr = wA.begin();
76
77 Type wArA = solverPerf.great_*pTraits<Type>::one;
78 Type wArAold = wArA;
79
80 // --- Calculate A.psi
81 this->matrix_.Amul(wA, psi);
82
83 // --- Calculate initial residual field
84 Field<Type> rA(this->matrix_.source() - wA);
85 Type* __restrict__ rAPtr = rA.begin();
86
87 // --- Calculate normalisation factor
88 Type normFactor = this->normFactor(psi, wA, pA);
89
90 if ((this->log_ >= 2) || (LduMatrix<Type, DType, LUType>::debug >= 2))
91 {
92 Info<< " Normalisation factor = " << normFactor << endl;
93 }
94
95 // --- Calculate normalised residual norm
96 solverPerf.initialResidual() = cmptDivide(gSumCmptMag(rA), normFactor);
97 solverPerf.finalResidual() = solverPerf.initialResidual();
98
99 // --- Check convergence, solve if not converged
100 if
101 (
102 this->minIter_ > 0
103 || !solverPerf.checkConvergence
104 (
105 this->tolerance_,
106 this->relTol_,
107 this->log_
108 )
109 )
110 {
111 // --- Select and construct the preconditioner
114 (
115 *this,
116 this->controlDict_
117 );
118
119 // --- Solver iteration
120 do
121 {
122 // --- Store previous wArA
123 wArAold = wArA;
124
125 // --- Precondition residual
126 preconPtr->precondition(wA, rA);
127
128 // --- Update search directions:
129 wArA = gSumCmptProd(wA, rA);
130
131 if (nIter == 0)
132 {
133 for (label cell=0; cell<nCells; cell++)
134 {
135 pAPtr[cell] = wAPtr[cell];
136 }
137 }
138 else
139 {
140 Type beta = cmptDivide
141 (
142 wArA,
143 stabilise(wArAold, solverPerf.vsmall_)
144 );
145
146 for (label cell=0; cell<nCells; cell++)
147 {
148 pAPtr[cell] = wAPtr[cell] + cmptMultiply(beta, pAPtr[cell]);
149 }
150 }
151
152
153 // --- Update preconditioned residual
154 this->matrix_.Amul(wA, pA);
155
156 Type wApA = gSumCmptProd(wA, pA);
157
158
159 // --- Test for singularity
160 if
161 (
162 solverPerf.checkSingularity
163 (
164 cmptDivide(cmptMag(wApA), normFactor)
165 )
166 )
167 {
168 break;
169 }
170
171
172 // --- Update solution and residual:
173
174 Type alpha = cmptDivide
175 (
176 wArA,
177 stabilise(wApA, solverPerf.vsmall_)
178 );
179
180 for (label cell=0; cell<nCells; cell++)
181 {
182 psiPtr[cell] += cmptMultiply(alpha, pAPtr[cell]);
183 rAPtr[cell] -= cmptMultiply(alpha, wAPtr[cell]);
184 }
185
186 solverPerf.finalResidual() =
187 cmptDivide(gSumCmptMag(rA), normFactor);
188
189 } while
190 (
191 (
192 nIter++ < this->maxIter_
193 && !solverPerf.checkConvergence
194 (
195 this->tolerance_,
196 this->relTol_,
197 this->log_
198 )
199 )
200 || nIter < this->minIter_
201 );
202 }
203
204 solverPerf.nIterations() =
206
207 return solverPerf;
208}
209
210
211// ************************************************************************* //
Generic templated field type.
Definition: Field.H:82
LduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition: LduMatrix.H:88
Preconditioned conjugate gradient solver for symmetric lduMatrices using a run-time selectable precon...
Definition: PCICG.H:56
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
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 vsmall_
Very small Type for the use in solvers.
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.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition: UListI.H:329
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition: one.H:62
A traits class, which is primarily used for primitives.
Definition: pTraits.H:59
Base class for solution control classes.
Definition: solver.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
const volScalarField & psi
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensionedScalar stabilise(const dimensionedScalar &x, const dimensionedScalar &y)
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
Type gSumCmptMag(const UList< Type > &f, const label comm)
Type gSumCmptProd(const UList< Type > &f1, const UList< Type > &f2, const label comm)
volScalarField & alpha
CEqn solve()
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)