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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "PCICG.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class Type, class DType, class LUType>
34 (
35  const word& fieldName,
36  const LduMatrix<Type, DType, LUType>& matrix,
37  const dictionary& solverDict
38 )
39 :
41  (
42  fieldName,
43  matrix,
44  solverDict
45  )
46 {}
47 
48 
49 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
50 
51 template<class Type, class DType, class LUType>
54 {
55  const word preconditionerName(this->controlDict_.getWord("preconditioner"));
56 
57  // --- Setup class containing solver performance data
58  SolverPerformance<Type> solverPerf
59  (
60  preconditionerName + typeName,
61  this->fieldName_
62  );
63 
64  label nIter = 0;
65 
66  label nCells = psi.size();
67 
68  Type* __restrict__ psiPtr = psi.begin();
69 
70  Field<Type> pA(nCells);
71  Type* __restrict__ pAPtr = pA.begin();
72 
73  Field<Type> wA(nCells);
74  Type* __restrict__ wAPtr = wA.begin();
75 
76  Type wArA = solverPerf.great_*pTraits<Type>::one;
77  Type wArAold = wArA;
78 
79  // --- Calculate A.psi
80  this->matrix_.Amul(wA, psi);
81 
82  // --- Calculate initial residual field
83  Field<Type> rA(this->matrix_.source() - wA);
84  Type* __restrict__ rAPtr = rA.begin();
85 
86  // --- Calculate normalisation factor
87  Type normFactor = this->normFactor(psi, wA, pA);
88 
90  {
91  Info<< " Normalisation factor = " << normFactor << endl;
92  }
93 
94  // --- Calculate normalised residual norm
95  solverPerf.initialResidual() = cmptDivide(gSumCmptMag(rA), normFactor);
96  solverPerf.finalResidual() = solverPerf.initialResidual();
97 
98  // --- Check convergence, solve if not converged
99  if
100  (
101  this->minIter_ > 0
102  || !solverPerf.checkConvergence(this->tolerance_, this->relTol_)
103  )
104  {
105  // --- Select and construct the preconditioner
108  (
109  *this,
110  this->controlDict_
111  );
112 
113  // --- Solver iteration
114  do
115  {
116  // --- Store previous wArA
117  wArAold = wArA;
118 
119  // --- Precondition residual
120  preconPtr->precondition(wA, rA);
121 
122  // --- Update search directions:
123  wArA = gSumCmptProd(wA, rA);
124 
125  if (nIter == 0)
126  {
127  for (label cell=0; cell<nCells; cell++)
128  {
129  pAPtr[cell] = wAPtr[cell];
130  }
131  }
132  else
133  {
134  Type beta = cmptDivide
135  (
136  wArA,
137  stabilise(wArAold, solverPerf.vsmall_)
138  );
139 
140  for (label cell=0; cell<nCells; cell++)
141  {
142  pAPtr[cell] = wAPtr[cell] + cmptMultiply(beta, pAPtr[cell]);
143  }
144  }
145 
146 
147  // --- Update preconditioned residual
148  this->matrix_.Amul(wA, pA);
149 
150  Type wApA = gSumCmptProd(wA, pA);
151 
152 
153  // --- Test for singularity
154  if
155  (
156  solverPerf.checkSingularity
157  (
158  cmptDivide(cmptMag(wApA), normFactor)
159  )
160  )
161  {
162  break;
163  }
164 
165 
166  // --- Update solution and residual:
167 
168  Type alpha = cmptDivide
169  (
170  wArA,
171  stabilise(wApA, solverPerf.vsmall_)
172  );
173 
174  for (label cell=0; cell<nCells; cell++)
175  {
176  psiPtr[cell] += cmptMultiply(alpha, pAPtr[cell]);
177  rAPtr[cell] -= cmptMultiply(alpha, wAPtr[cell]);
178  }
179 
180  solverPerf.finalResidual() =
181  cmptDivide(gSumCmptMag(rA), normFactor);
182 
183  } while
184  (
185  (
186  nIter++ < this->maxIter_
187  && !solverPerf.checkConvergence(this->tolerance_, this->relTol_)
188  )
189  || nIter < this->minIter_
190  );
191  }
192 
193  solverPerf.nIterations() =
195 
196  return solverPerf;
197 }
198 
199 
200 // ************************************************************************* //
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::cmptMultiply
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::LduMatrix
LduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition: LduMatrix.H:72
Foam::one
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition: one.H:61
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::LduMatrix::solver
Abstract base-class for LduMatrix solvers.
Definition: LduMatrix.H:115
Foam::gSumCmptMag
Type gSumCmptMag(const UList< Type > &f, const label comm)
Definition: FieldFunctions.C:592
Foam::cmptMag
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:400
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::stabilise
tmp< DimensionedField< scalar, GeoMesh > > stabilise(const DimensionedField< scalar, GeoMesh > &dsf, const dimensioned< scalar > &ds)
Definition: DimensionedScalarField.C:43
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
beta
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::PCICG::solve
virtual SolverPerformance< Type > solve(Field< Type > &psi) const
Solve the matrix with this solver.
Definition: PCICG.C:53
Foam::PCICG
Preconditioned conjugate gradient solver for symmetric lduMatrices using a run-time selectable precon...
Definition: PCICG.H:53
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::cmptDivide
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
PCICG.H
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:54
psi
const volScalarField & psi
Definition: createFieldRefs.H:1
Foam::SolverPerformance
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
Definition: SolverPerformance.H:51
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:54
Foam::gSumCmptProd
Type gSumCmptProd(const UList< Type > &f1, const UList< Type > &f2, const label comm)
Definition: FieldFunctions.C:620