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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "PBiCICG.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class Type, class DType, class LUType>
32 (
33  const word& fieldName,
34  const LduMatrix<Type, DType, LUType>& matrix,
35  const dictionary& solverDict
36 )
37 :
39  (
40  fieldName,
41  matrix,
42  solverDict
43  )
44 {}
45 
46 
47 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
48 
49 template<class Type, class DType, class LUType>
52 {
53  word preconditionerName(this->controlDict_.lookup("preconditioner"));
54 
55  // --- Setup class containing solver performance data
56  SolverPerformance<Type> solverPerf
57  (
58  preconditionerName + typeName,
59  this->fieldName_
60  );
61 
62  label nCells = psi.size();
63 
64  Type* __restrict__ psiPtr = psi.begin();
65 
66  Field<Type> pA(nCells);
67  Type* __restrict__ pAPtr = pA.begin();
68 
69  Field<Type> pT(nCells, pTraits<Type>::zero);
70  Type* __restrict__ pTPtr = pT.begin();
71 
72  Field<Type> wA(nCells);
73  Type* __restrict__ wAPtr = wA.begin();
74 
75  Field<Type> wT(nCells);
76  Type* __restrict__ wTPtr = wT.begin();
77 
78  Type wArT = solverPerf.great_*pTraits<Type>::one;
79  Type wArTold = wArT;
80 
81  // --- Calculate A.psi and T.psi
82  this->matrix_.Amul(wA, psi);
83  this->matrix_.Tmul(wT, psi);
84 
85  // --- Calculate initial residual and transpose residual fields
86  Field<Type> rA(this->matrix_.source() - wA);
87  Field<Type> rT(this->matrix_.source() - wT);
88  Type* __restrict__ rAPtr = rA.begin();
89  Type* __restrict__ rTPtr = rT.begin();
90 
91  // --- Calculate normalisation factor
92  Type normFactor = this->normFactor(psi, wA, pA);
93 
95  {
96  Info<< " Normalisation factor = " << normFactor << endl;
97  }
98 
99  // --- Calculate normalised residual norm
100  solverPerf.initialResidual() = cmptDivide(gSumCmptMag(rA), normFactor);
101  solverPerf.finalResidual() = solverPerf.initialResidual();
102 
103  // --- Check convergence, solve if not converged
104  if (!solverPerf.checkConvergence(this->tolerance_, this->relTol_))
105  {
106  // --- Select and construct the preconditioner
109  (
110  *this,
111  this->controlDict_
112  );
113 
114  // --- Solver iteration
115  do
116  {
117  // --- Store previous wArT
118  wArTold = wArT;
119 
120  // --- Precondition residuals
121  preconPtr->precondition(wA, rA);
122  preconPtr->preconditionT(wT, rT);
123 
124  // --- Update search directions:
125  wArT = gSumCmptProd(wA, rT);
126 
127  if (solverPerf.nIterations() == 0)
128  {
129  for (label cell=0; cell<nCells; cell++)
130  {
131  pAPtr[cell] = wAPtr[cell];
132  pTPtr[cell] = wTPtr[cell];
133  }
134  }
135  else
136  {
137  Type beta = cmptDivide
138  (
139  wArT,
140  stabilise(wArTold, solverPerf.vsmall_)
141  );
142 
143  for (label cell=0; cell<nCells; cell++)
144  {
145  pAPtr[cell] = wAPtr[cell] + cmptMultiply(beta, pAPtr[cell]);
146  pTPtr[cell] = wTPtr[cell] + cmptMultiply(beta, pTPtr[cell]);
147  }
148  }
149 
150 
151  // --- Update preconditioned residuals
152  this->matrix_.Amul(wA, pA);
153  this->matrix_.Tmul(wT, pT);
154 
155  Type wApT = gSumCmptProd(wA, pT);
156 
157  // --- Test for singularity
158  if
159  (
160  solverPerf.checkSingularity
161  (
162  cmptDivide(cmptMag(wApT), normFactor)
163  )
164  )
165  {
166  break;
167  }
168 
169 
170  // --- Update solution and residual:
171 
172  Type alpha = cmptDivide
173  (
174  wArT,
175  stabilise(wApT, solverPerf.vsmall_)
176  );
177 
178  for (label cell=0; cell<nCells; cell++)
179  {
180  psiPtr[cell] += cmptMultiply(alpha, pAPtr[cell]);
181  rAPtr[cell] -= cmptMultiply(alpha, wAPtr[cell]);
182  rTPtr[cell] -= cmptMultiply(alpha, wTPtr[cell]);
183  }
184 
185  solverPerf.finalResidual() =
186  cmptDivide(gSumCmptMag(rA), normFactor);
187 
188  } while
189  (
190  solverPerf.nIterations()++ < this->maxIter_
191  && !(solverPerf.checkConvergence(this->tolerance_, this->relTol_))
192  );
193  }
194 
195  return solverPerf;
196 }
197 
198 
199 // ************************************************************************* //
beta
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
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:216
Foam::SolverPerformance::great_
static const scalar great_
Large Type for the use in solvers.
Definition: SolverPerformance.H:99
Foam::SolverPerformance::finalResidual
const Type & finalResidual() const
Return final residual.
Definition: SolverPerformance.H:177
Foam::LduMatrix< Type, DType, LUType >
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::SolverPerformance::vsmall_
static const scalar vsmall_
Very small Type for the use in solvers.
Definition: SolverPerformance.H:105
Foam::gSumCmptMag
Type gSumCmptMag(const UList< Type > &f, const label comm)
Definition: FieldFunctions.C:548
Foam::cmptMag
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:359
Foam::SolverPerformance::checkConvergence
bool checkConvergence(const Type &tolerance, const Type &relTolerance)
Check, store and return convergence.
Definition: SolverPerformance.C:61
Foam::PBiCICG::PBiCICG
PBiCICG(const PBiCICG &)
Disallow default bitwise copy construct.
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field< Type >
Foam::stabilise
tmp< DimensionedField< scalar, GeoMesh > > stabilise(const DimensionedField< scalar, GeoMesh > &dsf, const dimensioned< scalar > &ds)
Definition: DimensionedScalarField.C:40
Foam::Info
messageStream Info
Foam::PBiCICG::solve
virtual SolverPerformance< Type > solve(Field< Type > &psi) const
Solve the matrix with this solver.
Definition: PBiCICG.C:51
Foam::SolverPerformance::initialResidual
const Type & initialResidual() const
Return initial residual.
Definition: SolverPerformance.H:164
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::cmptDivide
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
PBiCICG.H
psi
const volScalarField & psi
Definition: setRegionFluidFields.H:13
Foam::SolverPerformance::checkSingularity
bool checkSingularity(const Type &residual)
Singularity test.
Definition: SolverPerformance.C:33
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
Foam::SolverPerformance
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
Definition: SolverPerformance.H:49
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
Foam::gSumCmptProd
Type gSumCmptProd(const UList< Type > &f1, const UList< Type > &f2, const label comm)
Definition: FieldFunctions.C:567
Foam::SolverPerformance::nIterations
label nIterations() const
Return number of iterations.
Definition: SolverPerformance.H:190