30 template<
class Type,
class DType,
class LUType>
33 const word& fieldName,
49 template<
class Type,
class DType,
class LUType>
53 word preconditionerName(this->controlDict_.lookup(
"preconditioner"));
58 preconditionerName + typeName,
64 Type* __restrict__ psiPtr = psi.
begin();
67 Type* __restrict__ pAPtr = pA.
begin();
70 Type* __restrict__ pTPtr = pT.
begin();
73 Type* __restrict__ wAPtr = wA.
begin();
76 Type* __restrict__ wTPtr = wT.
begin();
82 this->matrix_.Amul(wA, psi);
83 this->matrix_.Tmul(wT, psi);
88 Type* __restrict__ rAPtr = rA.
begin();
89 Type* __restrict__ rTPtr = rT.begin();
92 Type normFactor = this->normFactor(psi, wA, pA);
96 Info<<
" Normalisation factor = " << normFactor <<
endl;
101 solverPerf.finalResidual() = solverPerf.initialResidual();
104 if (!solverPerf.checkConvergence(this->tolerance_, this->relTol_))
121 preconPtr->precondition(wA, rA);
122 preconPtr->preconditionT(wT, rT);
127 if (solverPerf.nIterations() == 0)
146 pTPtr[cell] = wTPtr[cell] +
cmptMultiply(beta, pTPtr[cell]);
152 this->matrix_.Amul(wA, pA);
153 this->matrix_.Tmul(wT, pT);
160 solverPerf.checkSingularity
185 solverPerf.finalResidual() =
190 solverPerf.nIterations()++ < this->maxIter_
191 && !(solverPerf.checkConvergence(this->tolerance_, this->relTol_))
A cell is defined as a list of faces with extra functionality.
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics...
A class for handling words, derived from string.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
void size(const label)
Override size to be inconsistent with allocated storage.
A list of keyword definitions, which are a keyword followed by any number of values (e...
virtual SolverPerformance< Type > solve(Field< Type > &psi) const
Solve the matrix with this solver.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
iterator begin()
Return an iterator to begin traversing the UList.
Abstract base-class for LduMatrix solvers.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
const volScalarField & psi
Pre-declare SubField and related Field type.
tmp< DimensionedField< scalar, GeoMesh > > stabilise(const DimensionedField< scalar, GeoMesh > &dsf, const dimensioned< scalar > &ds)
Traits class for primitives.
Type gSumCmptMag(const UList< Type > &f, const label comm)
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
Type gSumCmptProd(const UList< Type > &f1, const UList< Type > &f2, const label comm)
Preconditioned bi-conjugate gradient solver for asymmetric lduMatrices using a run-time selectable pr...
LduMatrix is a general matrix class in which the coefficients are stored as three arrays...
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...