PBiCCCG.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 "PBiCCCG.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 (
56 ) const
57 {
58  const word preconditionerName(this->controlDict_.getWord("preconditioner"));
59 
60  // --- Setup class containing solver performance data
61  SolverPerformance<Type> solverPerf
62  (
63  preconditionerName + typeName,
64  this->fieldName_
65  );
66 
67  label nIter = 0;
68 
69  label nCells = psi.size();
70 
71  Type* __restrict__ psiPtr = psi.begin();
72 
73  Field<Type> pA(nCells);
74  Type* __restrict__ pAPtr = pA.begin();
75 
76  Field<Type> pT(nCells, Zero);
77  Type* __restrict__ pTPtr = pT.begin();
78 
79  Field<Type> wA(nCells);
80  Type* __restrict__ wAPtr = wA.begin();
81 
82  Field<Type> wT(nCells);
83  Type* __restrict__ wTPtr = wT.begin();
84 
85  scalar wArT = 1e15; //this->matrix_.great_;
86  scalar wArTold = wArT;
87 
88  // --- Calculate A.psi and T.psi
89  this->matrix_.Amul(wA, psi);
90  this->matrix_.Tmul(wT, psi);
91 
92  // --- Calculate initial residual and transpose residual fields
93  Field<Type> rA(this->matrix_.source() - wA);
94  Field<Type> rT(this->matrix_.source() - wT);
95  Type* __restrict__ rAPtr = rA.begin();
96  Type* __restrict__ rTPtr = rT.begin();
97 
98  // --- Calculate normalisation factor
99  Type normFactor = this->normFactor(psi, wA, pA);
100 
102  {
103  Info<< " Normalisation factor = " << normFactor << endl;
104  }
105 
106  // --- Calculate normalised residual norm
107  solverPerf.initialResidual() = cmptDivide(gSumCmptMag(rA), normFactor);
108  solverPerf.finalResidual() = solverPerf.initialResidual();
109 
110  // --- Check convergence, solve if not converged
111  if
112  (
113  this->minIter_ > 0
114  || !solverPerf.checkConvergence(this->tolerance_, this->relTol_)
115  )
116  {
117  // --- Select and construct the preconditioner
120  (
121  *this,
122  this->controlDict_
123  );
124 
125  // --- Solver iteration
126  do
127  {
128  // --- Store previous wArT
129  wArTold = wArT;
130 
131  // --- Precondition residuals
132  preconPtr->precondition(wA, rA);
133  preconPtr->preconditionT(wT, rT);
134 
135  // --- Update search directions:
136  wArT = gSumProd(wA, rT);
137 
138  if (nIter == 0)
139  {
140  for (label cell=0; cell<nCells; cell++)
141  {
142  pAPtr[cell] = wAPtr[cell];
143  pTPtr[cell] = wTPtr[cell];
144  }
145  }
146  else
147  {
148  scalar beta = wArT/wArTold;
149 
150  for (label cell=0; cell<nCells; cell++)
151  {
152  pAPtr[cell] = wAPtr[cell] + (beta* pAPtr[cell]);
153  pTPtr[cell] = wTPtr[cell] + (beta* pTPtr[cell]);
154  }
155  }
156 
157 
158  // --- Update preconditioned residuals
159  this->matrix_.Amul(wA, pA);
160  this->matrix_.Tmul(wT, pT);
161 
162  scalar wApT = gSumProd(wA, pT);
163 
164  // --- Test for singularity
165  if
166  (
167  solverPerf.checkSingularity
168  (
169  cmptDivide(pTraits<Type>::one*mag(wApT), normFactor)
170  )
171  )
172  {
173  break;
174  }
175 
176 
177  // --- Update solution and residual:
178 
179  scalar alpha = wArT/wApT;
180 
181  for (label cell=0; cell<nCells; cell++)
182  {
183  psiPtr[cell] += (alpha* pAPtr[cell]);
184  rAPtr[cell] -= (alpha* wAPtr[cell]);
185  rTPtr[cell] -= (alpha* wTPtr[cell]);
186  }
187 
188  solverPerf.finalResidual() =
189  cmptDivide(gSumCmptMag(rA), normFactor);
190 
191  } while
192  (
193  (
194  nIter++ < this->maxIter_
195  && !solverPerf.checkConvergence(this->tolerance_, this->relTol_)
196  )
197  || nIter < this->minIter_
198  );
199  }
200 
201  solverPerf.nIterations() =
203 
204  return solverPerf;
205 }
206 
207 
208 // ************************************************************************* //
Foam::PBiCCCG::solve
virtual SolverPerformance< Type > solve(Field< Type > &psi) const
Solve the matrix with this solver.
Definition: PBiCCCG.C:54
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
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
PBiCCCG.H
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::Field
Generic templated field type.
Definition: Field.H:63
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
beta
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
Foam::gSumProd
scalarProduct< Type, Type >::type gSumProd(const UList< Type > &f1, const UList< Type > &f2, const label comm)
Definition: FieldFunctions.C:605
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
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:54
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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::PBiCCCG
Preconditioned bi-conjugate gradient solver for asymmetric lduMatrices using a run-time selectable pr...
Definition: PBiCCCG.H:53