30 template<
class Type,
class DType,
class LUType>
33 const word& fieldName,
49 template<
class Type,
class DType,
class LUType>
56 word preconditionerName(this->controlDict_.lookup(
"preconditioner"));
61 preconditionerName + typeName,
67 Type* __restrict__ psiPtr = psi.
begin();
70 Type* __restrict__ pAPtr = pA.
begin();
73 Type* __restrict__ pTPtr = pT.
begin();
76 Type* __restrict__ wAPtr = wA.
begin();
79 Type* __restrict__ wTPtr = wT.
begin();
82 scalar wArTold = wArT;
85 this->matrix_.Amul(wA, psi);
86 this->matrix_.Tmul(wT, psi);
91 Type* __restrict__ rAPtr = rA.
begin();
92 Type* __restrict__ rTPtr = rT.begin();
95 Type normFactor = this->normFactor(psi, wA, pA);
99 Info<<
" Normalisation factor = " << normFactor <<
endl;
104 solverPerf.finalResidual() = solverPerf.initialResidual();
110 || !solverPerf.checkConvergence(this->tolerance_, this->relTol_)
128 preconPtr->precondition(wA, rA);
129 preconPtr->preconditionT(wT, rT);
134 if (solverPerf.nIterations() == 0)
144 scalar
beta = wArT/wArTold;
155 this->matrix_.Amul(wA, pA);
156 this->matrix_.Tmul(wT, pT);
163 solverPerf.checkSingularity
175 scalar
alpha = wArT/wApT;
179 psiPtr[
cell] += (alpha* pAPtr[
cell]);
180 rAPtr[
cell] -= (alpha* wAPtr[
cell]);
181 rTPtr[
cell] -= (alpha* wTPtr[
cell]);
184 solverPerf.finalResidual() =
190 solverPerf.nIterations()++ < this->maxIter_
191 && !solverPerf.checkConvergence(this->tolerance_, this->relTol_)
193 || solverPerf.nIterations() < this->minIter_
dimensioned< scalar > mag(const dimensioned< Type > &)
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...
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: [].
Preconditioned bi-conjugate gradient solver for asymmetric lduMatrices using a run-time selectable pr...
Pre-declare SubField and related Field type.
scalar gSumProd(const UList< Type > &f1, const UList< Type > &f2, const label comm)
Traits class for primitives.
Type gSumCmptMag(const UList< Type > &f, const label comm)
LduMatrix is a general matrix class in which the coefficients are stored as three arrays...
virtual SolverPerformance< Type > solve(Field< Type > &psi) const
Solve the matrix with this solver.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...