PBiCICG.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 "PBiCICG.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> pT(nCells, Zero);
74  Type* __restrict__ pTPtr = pT.begin();
75 
76  Field<Type> wA(nCells);
77  Type* __restrict__ wAPtr = wA.begin();
78 
79  Field<Type> wT(nCells);
80  Type* __restrict__ wTPtr = wT.begin();
81 
82  Type wArT = solverPerf.great_*pTraits<Type>::one;
83  Type wArTold = wArT;
84 
85  // --- Calculate A.psi and T.psi
86  this->matrix_.Amul(wA, psi);
87  this->matrix_.Tmul(wT, psi);
88 
89  // --- Calculate initial residual and transpose residual fields
90  Field<Type> rA(this->matrix_.source() - wA);
91  Field<Type> rT(this->matrix_.source() - wT);
92  Type* __restrict__ rAPtr = rA.begin();
93  Type* __restrict__ rTPtr = rT.begin();
94 
95  // --- Calculate normalisation factor
96  Type normFactor = this->normFactor(psi, wA, pA);
97 
99  {
100  Info<< " Normalisation factor = " << normFactor << endl;
101  }
102 
103  // --- Calculate normalised residual norm
104  solverPerf.initialResidual() = cmptDivide(gSumCmptMag(rA), normFactor);
105  solverPerf.finalResidual() = solverPerf.initialResidual();
106 
107  // --- Check convergence, solve if not converged
108  if (!solverPerf.checkConvergence(this->tolerance_, this->relTol_))
109  {
110  // --- Select and construct the preconditioner
113  (
114  *this,
115  this->controlDict_
116  );
117 
118  // --- Solver iteration
119  do
120  {
121  // --- Store previous wArT
122  wArTold = wArT;
123 
124  // --- Precondition residuals
125  preconPtr->precondition(wA, rA);
126  preconPtr->preconditionT(wT, rT);
127 
128  // --- Update search directions:
129  wArT = gSumCmptProd(wA, rT);
130 
131  if (nIter == 0)
132  {
133  for (label cell=0; cell<nCells; cell++)
134  {
135  pAPtr[cell] = wAPtr[cell];
136  pTPtr[cell] = wTPtr[cell];
137  }
138  }
139  else
140  {
141  Type beta = cmptDivide
142  (
143  wArT,
144  stabilise(wArTold, solverPerf.vsmall_)
145  );
146 
147  for (label cell=0; cell<nCells; cell++)
148  {
149  pAPtr[cell] = wAPtr[cell] + cmptMultiply(beta, pAPtr[cell]);
150  pTPtr[cell] = wTPtr[cell] + cmptMultiply(beta, pTPtr[cell]);
151  }
152  }
153 
154 
155  // --- Update preconditioned residuals
156  this->matrix_.Amul(wA, pA);
157  this->matrix_.Tmul(wT, pT);
158 
159  Type wApT = gSumCmptProd(wA, pT);
160 
161  // --- Test for singularity
162  if
163  (
164  solverPerf.checkSingularity
165  (
166  cmptDivide(cmptMag(wApT), normFactor)
167  )
168  )
169  {
170  break;
171  }
172 
173 
174  // --- Update solution and residual:
175 
176  Type alpha = cmptDivide
177  (
178  wArT,
179  stabilise(wApT, solverPerf.vsmall_)
180  );
181 
182  for (label cell=0; cell<nCells; cell++)
183  {
184  psiPtr[cell] += cmptMultiply(alpha, pAPtr[cell]);
185  rAPtr[cell] -= cmptMultiply(alpha, wAPtr[cell]);
186  rTPtr[cell] -= cmptMultiply(alpha, wTPtr[cell]);
187  }
188 
189  solverPerf.finalResidual() =
190  cmptDivide(gSumCmptMag(rA), normFactor);
191 
192  } while
193  (
194  nIter++ < this->maxIter_
195  && !(solverPerf.checkConvergence(this->tolerance_, this->relTol_))
196  );
197  }
198 
199  solverPerf.nIterations() =
201 
202  return solverPerf;
203 }
204 
205 
206 // ************************************************************************* //
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::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
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::PBiCICG::solve
virtual SolverPerformance< Type > solve(Field< Type > &psi) const
Solve the matrix with this solver.
Definition: PBiCICG.C:53
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
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
PBiCICG.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
Foam::PBiCICG
Preconditioned bi-conjugate gradient solver for asymmetric lduMatrices using a run-time selectable pr...
Definition: PBiCICG.H:53