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 | 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 "PCICG.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> wA(nCells);
70  Type* __restrict__ wAPtr = wA.begin();
71 
72  Type wArA = solverPerf.great_*pTraits<Type>::one;
73  Type wArAold = wArA;
74 
75  // --- Calculate A.psi
76  this->matrix_.Amul(wA, psi);
77 
78  // --- Calculate initial residual field
79  Field<Type> rA(this->matrix_.source() - wA);
80  Type* __restrict__ rAPtr = rA.begin();
81 
82  // --- Calculate normalisation factor
83  Type normFactor = this->normFactor(psi, wA, pA);
84 
86  {
87  Info<< " Normalisation factor = " << normFactor << endl;
88  }
89 
90  // --- Calculate normalised residual norm
91  solverPerf.initialResidual() = cmptDivide(gSumCmptMag(rA), normFactor);
92  solverPerf.finalResidual() = solverPerf.initialResidual();
93 
94  // --- Check convergence, solve if not converged
95  if
96  (
97  this->minIter_ > 0
98  || !solverPerf.checkConvergence(this->tolerance_, this->relTol_)
99  )
100  {
101  // --- Select and construct the preconditioner
104  (
105  *this,
106  this->controlDict_
107  );
108 
109  // --- Solver iteration
110  do
111  {
112  // --- Store previous wArA
113  wArAold = wArA;
114 
115  // --- Precondition residual
116  preconPtr->precondition(wA, rA);
117 
118  // --- Update search directions:
119  wArA = gSumCmptProd(wA, rA);
120 
121  if (solverPerf.nIterations() == 0)
122  {
123  for (label cell=0; cell<nCells; cell++)
124  {
125  pAPtr[cell] = wAPtr[cell];
126  }
127  }
128  else
129  {
130  Type beta = cmptDivide
131  (
132  wArA,
133  stabilise(wArAold, solverPerf.vsmall_)
134  );
135 
136  for (label cell=0; cell<nCells; cell++)
137  {
138  pAPtr[cell] = wAPtr[cell] + cmptMultiply(beta, pAPtr[cell]);
139  }
140  }
141 
142 
143  // --- Update preconditioned residual
144  this->matrix_.Amul(wA, pA);
145 
146  Type wApA = gSumCmptProd(wA, pA);
147 
148 
149  // --- Test for singularity
150  if
151  (
152  solverPerf.checkSingularity
153  (
154  cmptDivide(cmptMag(wApA), normFactor)
155  )
156  )
157  {
158  break;
159  }
160 
161 
162  // --- Update solution and residual:
163 
164  Type alpha = cmptDivide
165  (
166  wArA,
167  stabilise(wApA, solverPerf.vsmall_)
168  );
169 
170  for (label cell=0; cell<nCells; cell++)
171  {
172  psiPtr[cell] += cmptMultiply(alpha, pAPtr[cell]);
173  rAPtr[cell] -= cmptMultiply(alpha, wAPtr[cell]);
174  }
175 
176  solverPerf.finalResidual() =
177  cmptDivide(gSumCmptMag(rA), normFactor);
178 
179  } while
180  (
181  (
182  solverPerf.nIterations()++ < this->maxIter_
183  && !solverPerf.checkConvergence(this->tolerance_, this->relTol_)
184  )
185  || solverPerf.nIterations() < this->minIter_
186  );
187  }
188 
189  return solverPerf;
190 }
191 
192 
193 // ************************************************************************* //
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)
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)
Type gSumCmptProd(const UList< Type > &f1, const UList< Type > &f2, const label comm)
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.
Preconditioned conjugate gradient solver for symmetric lduMatrices using a run-time selectable precon...
Definition: PCICG.H:50
Abstract base-class for LduMatrix solvers.
Definition: LduMatrix.H:112
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
virtual SolverPerformance< Type > solve(Field< Type > &psi) const
Solve the matrix with this solver.
Definition: PCICG.C:51
tmp< DimensionedField< scalar, GeoMesh > > stabilise(const DimensionedField< scalar, GeoMesh > &dsf, const dimensioned< scalar > &ds)
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)
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 volScalarField & psi
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].