Go to the documentation of this file.
33 template<
class Type,
class DType,
class LUType>
36 const word& fieldName,
52 template<
class Type,
class DType,
class LUType>
59 const word preconditionerName(this->controlDict_.getWord(
"preconditioner"));
64 preconditionerName + typeName,
70 label nCells =
psi.size();
72 Type* __restrict__ psiPtr =
psi.begin();
75 Type* __restrict__ pAPtr = pA.begin();
78 Type* __restrict__ pTPtr = pT.begin();
81 Type* __restrict__ wAPtr = wA.begin();
84 Type* __restrict__ wTPtr = wT.begin();
87 scalar wArTold = wArT;
90 this->matrix_.Amul(wA,
psi);
91 this->matrix_.Tmul(wT,
psi);
96 Type* __restrict__ rAPtr = rA.begin();
97 Type* __restrict__ rTPtr = rT.begin();
100 Type normFactor = this->normFactor(
psi, wA, pA);
104 Info<<
" Normalisation factor = " << normFactor <<
endl;
109 solverPerf.finalResidual() = solverPerf.initialResidual();
115 || !solverPerf.checkConvergence
138 preconPtr->precondition(wA, rA);
139 preconPtr->preconditionT(wT, rT);
154 scalar
beta = wArT/wArTold;
165 this->matrix_.Amul(wA, pA);
166 this->matrix_.Tmul(wT, pT);
173 solverPerf.checkSingularity
185 scalar
alpha = wArT/wApT;
194 solverPerf.finalResidual() =
200 nIter++ < this->maxIter_
201 && !solverPerf.checkConvergence
208 || nIter < this->minIter_
212 solverPerf.nIterations() =
virtual SolverPerformance< Type > solve(Field< Type > &psi) const
Solve the matrix with this solver.
A class for handling words, derived from Foam::string.
static constexpr const zero Zero
Global zero (0)
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
LduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Ostream & endl(Ostream &os)
Add newline and flush stream.
Abstract base-class for LduMatrix solvers.
Type gSumCmptMag(const UList< Type > &f, const label comm)
Generic templated field type.
messageStream Info
Information stream (stdout output on master, null elsewhere)
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
scalarProduct< Type, Type >::type gSumProd(const UList< Type > &f1, const UList< Type > &f2, const label comm)
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A traits class, which is primarily used for primitives.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const volScalarField & psi
A cell is defined as a list of faces with extra functionality.
Preconditioned bi-conjugate gradient solver for asymmetric lduMatrices using a run-time selectable pr...