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 | Copyright (C) 2011-2016 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 "PBiCCCG.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  Field<Type>& psi
54 ) const
55 {
56  word preconditionerName(this->controlDict_.lookup("preconditioner"));
57 
58  // --- Setup class containing solver performance data
59  SolverPerformance<Type> solverPerf
60  (
61  preconditionerName + typeName,
62  this->fieldName_
63  );
64 
65  label nCells = psi.size();
66 
67  Type* __restrict__ psiPtr = psi.begin();
68 
69  Field<Type> pA(nCells);
70  Type* __restrict__ pAPtr = pA.begin();
71 
72  Field<Type> pT(nCells, Zero);
73  Type* __restrict__ pTPtr = pT.begin();
74 
75  Field<Type> wA(nCells);
76  Type* __restrict__ wAPtr = wA.begin();
77 
78  Field<Type> wT(nCells);
79  Type* __restrict__ wTPtr = wT.begin();
80 
81  scalar wArT = 1e15; //this->matrix_.great_;
82  scalar wArTold = wArT;
83 
84  // --- Calculate A.psi and T.psi
85  this->matrix_.Amul(wA, psi);
86  this->matrix_.Tmul(wT, psi);
87 
88  // --- Calculate initial residual and transpose residual fields
89  Field<Type> rA(this->matrix_.source() - wA);
90  Field<Type> rT(this->matrix_.source() - wT);
91  Type* __restrict__ rAPtr = rA.begin();
92  Type* __restrict__ rTPtr = rT.begin();
93 
94  // --- Calculate normalisation factor
95  Type normFactor = this->normFactor(psi, wA, pA);
96 
98  {
99  Info<< " Normalisation factor = " << normFactor << endl;
100  }
101 
102  // --- Calculate normalised residual norm
103  solverPerf.initialResidual() = cmptDivide(gSumCmptMag(rA), normFactor);
104  solverPerf.finalResidual() = solverPerf.initialResidual();
105 
106  // --- Check convergence, solve if not converged
107  if
108  (
109  this->minIter_ > 0
110  || !solverPerf.checkConvergence(this->tolerance_, this->relTol_)
111  )
112  {
113  // --- Select and construct the preconditioner
116  (
117  *this,
118  this->controlDict_
119  );
120 
121  // --- Solver iteration
122  do
123  {
124  // --- Store previous wArT
125  wArTold = wArT;
126 
127  // --- Precondition residuals
128  preconPtr->precondition(wA, rA);
129  preconPtr->preconditionT(wT, rT);
130 
131  // --- Update search directions:
132  wArT = gSumProd(wA, rT);
133 
134  if (solverPerf.nIterations() == 0)
135  {
136  for (label cell=0; cell<nCells; cell++)
137  {
138  pAPtr[cell] = wAPtr[cell];
139  pTPtr[cell] = wTPtr[cell];
140  }
141  }
142  else
143  {
144  scalar beta = wArT/wArTold;
145 
146  for (label cell=0; cell<nCells; cell++)
147  {
148  pAPtr[cell] = wAPtr[cell] + (beta* pAPtr[cell]);
149  pTPtr[cell] = wTPtr[cell] + (beta* pTPtr[cell]);
150  }
151  }
152 
153 
154  // --- Update preconditioned residuals
155  this->matrix_.Amul(wA, pA);
156  this->matrix_.Tmul(wT, pT);
157 
158  scalar wApT = gSumProd(wA, pT);
159 
160  // --- Test for singularity
161  if
162  (
163  solverPerf.checkSingularity
164  (
165  cmptDivide(pTraits<Type>::one*mag(wApT), normFactor)
166  )
167  )
168  {
169  break;
170  }
171 
172 
173  // --- Update solution and residual:
174 
175  scalar alpha = wArT/wApT;
176 
177  for (label cell=0; cell<nCells; cell++)
178  {
179  psiPtr[cell] += (alpha* pAPtr[cell]);
180  rAPtr[cell] -= (alpha* wAPtr[cell]);
181  rTPtr[cell] -= (alpha* wTPtr[cell]);
182  }
183 
184  solverPerf.finalResidual() =
185  cmptDivide(gSumCmptMag(rA), normFactor);
186 
187  } while
188  (
189  (
190  solverPerf.nIterations()++ < this->maxIter_
191  && !solverPerf.checkConvergence(this->tolerance_, this->relTol_)
192  )
193  || solverPerf.nIterations() < this->minIter_
194  );
195  }
196 
197  return solverPerf;
198 }
199 
200 
201 // ************************************************************************* //
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
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
void size(const label)
Override size to be inconsistent with allocated storage.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Traits class for primitives.
Definition: pTraits.H:50
Type gSumCmptMag(const UList< Type > &f, const label comm)
Preconditioned bi-conjugate gradient solver for asymmetric lduMatrices using a run-time selectable pr...
Definition: PBiCCCG.H:50
virtual SolverPerformance< Type > solve(Field< Type > &psi) const
Solve the matrix with this solver.
Definition: PBiCCCG.C:52
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Pre-declare SubField and related Field type.
Definition: Field.H:57
A class for handling words, derived from string.
Definition: word.H:59
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics...
iterator begin()
Return an iterator to begin traversing the UList.
static const zero Zero
Definition: zero.H:91
scalar gSumProd(const UList< Type > &f1, const UList< Type > &f2, const label comm)
Abstract base-class for LduMatrix solvers.
Definition: LduMatrix.H:112
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
LduMatrix is a general matrix class in which the coefficients are stored as three arrays...
Definition: LduMatrix.H:69
messageStream Info
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
dimensioned< scalar > mag(const dimensioned< Type > &)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].