33 template<
class Type,
class DType,
class LUType>
36 const word& fieldName,
52 template<
class Type,
class DType,
class LUType>
56 const word preconditionerName(this->controlDict_.getWord(
"preconditioner"));
61 preconditionerName + typeName,
67 label nCells =
psi.size();
69 Type* __restrict__ psiPtr =
psi.begin();
72 Type* __restrict__ pAPtr = pA.begin();
75 Type* __restrict__ pTPtr = pT.begin();
78 Type* __restrict__ wAPtr = wA.begin();
81 Type* __restrict__ wTPtr = wT.begin();
87 this->matrix_.Amul(wA,
psi);
88 this->matrix_.Tmul(wT,
psi);
93 Type* __restrict__ rAPtr = rA.
begin();
94 Type* __restrict__ rTPtr = rT.begin();
97 Type normFactor = this->normFactor(
psi, wA, pA);
101 Info<<
" Normalisation factor = " << normFactor <<
endl;
106 solverPerf.finalResidual() = solverPerf.initialResidual();
111 !solverPerf.checkConvergence
120 autoPtr<typename LduMatrix<Type, DType, LUType>::preconditioner>
134 preconPtr->precondition(wA, rA);
135 preconPtr->preconditionT(wT, rT);
142 for (label cell=0; cell<nCells; cell++)
144 pAPtr[cell] = wAPtr[cell];
145 pTPtr[cell] = wTPtr[cell];
156 for (label cell=0; cell<nCells; cell++)
165 this->matrix_.Amul(wA, pA);
166 this->matrix_.Tmul(wT, pT);
173 solverPerf.checkSingularity
191 for (label cell=0; cell<nCells; cell++)
198 solverPerf.finalResidual() =
203 nIter++ < this->maxIter_
204 && !solverPerf.checkConvergence
213 solverPerf.nIterations() =
214 pTraits<typename pTraits<Type>::labelType>::one*nIter;
Base class for solution control classes.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Ostream & endl(Ostream &os)
Add newline and flush stream.
A traits class, which is primarily used for primitives.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Type gSumCmptMag(const UList< Type > &f, const label comm)
dimensionedScalar stabilise(const dimensionedScalar &x, const dimensionedScalar &y)
Preconditioned bi-conjugate gradient solver for asymmetric lduMatrices using a run-time selectable pr...
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)
Generic templated field type.
A class for handling words, derived from Foam::string.
iterator begin() noexcept
Return an iterator to begin traversing the UList.
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
LduMatrix is a general matrix class in which the coefficients are stored as three arrays...
messageStream Info
Information stream (stdout output on master, null elsewhere)
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
const volScalarField & psi
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
static constexpr const zero Zero
Global zero (0)