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,
66 Type* __restrict__ psiPtr = psi.
begin();
69 Type* __restrict__ pAPtr = pA.
begin();
72 Type* __restrict__ wAPtr = wA.
begin();
78 this->matrix_.Amul(wA, psi);
82 Type* __restrict__ rAPtr = rA.
begin();
85 Type normFactor = this->normFactor(psi, wA, pA);
89 Info<<
" Normalisation factor = " << normFactor <<
endl;
94 solverPerf.finalResidual() = solverPerf.initialResidual();
100 || !solverPerf.checkConvergence(this->tolerance_, this->relTol_)
118 preconPtr->precondition(wA, rA);
146 this->matrix_.Amul(wA, pA);
154 solverPerf.checkSingularity
178 solverPerf.finalResidual() =
184 nIter++ < this->maxIter_
185 && !solverPerf.checkConvergence(this->tolerance_, this->relTol_)
187 || nIter < this->minIter_
191 solverPerf.nIterations() =
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
A list of keyword definitions, which are a keyword followed by any number of values (e...
void size(const label)
Override size to be inconsistent with allocated storage.
PCICG(const word &fieldName, const LduMatrix< Type, DType, LUType > &matrix, const dictionary &solverDict)
Construct from matrix components and solver data dictionary.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Traits class for primitives.
Type gSumCmptMag(const UList< Type > &f, const label comm)
virtual SolverPerformance< Type > solve(Field< Type > &psi) const
Solve the matrix with this solver.
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
Type gSumCmptProd(const UList< Type > &f1, const UList< Type > &f2, const label comm)
Pre-declare SubField and related Field type.
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
A class for handling words, derived from string.
iterator begin()
Return an iterator to begin traversing the UList.
Abstract base-class for LduMatrix solvers.
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
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.
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...
const volScalarField & psi
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class representing the concept of 1 (scalar(1)) used to avoid unnecessary manipulations for objects...