Go to the documentation of this file.
34 lduMatrix::solver::addasymMatrixConstructorToTable<PBiCGStab>
43 const word& fieldName,
75 lduMatrix::preconditioner::getName(controlDict_) + typeName,
81 scalar* __restrict__ psiPtr =
psi.begin();
84 scalar* __restrict__ pAPtr = pA.begin();
87 scalar* __restrict__ yAPtr = yA.begin();
90 matrix_.Amul(yA,
psi, interfaceBouCoeffs_, interfaces_, cmpt);
94 scalar* __restrict__ rAPtr = rA.begin();
97 const scalar normFactor = this->normFactor(
psi, source, yA, pA);
99 if (lduMatrix::debug >= 2)
101 Info<<
" Normalisation factor = " << normFactor <<
endl;
118 scalar* __restrict__ AyAPtr = AyA.begin();
121 scalar* __restrict__ sAPtr = sA.begin();
124 scalar* __restrict__ zAPtr = zA.begin();
127 scalar* __restrict__ tAPtr = tA.begin();
149 const scalar rA0rAold = rA0rA;
175 const scalar
beta = (rA0rA/rA0rAold)*(
alpha/omega);
185 preconPtr->precondition(yA, pA, cmpt);
188 matrix_.Amul(AyA, yA, interfaceBouCoeffs_, interfaces_, cmpt);
190 const scalar rA0AyA =
gSumProd(rA0, AyA, matrix().
mesh().comm());
192 alpha = rA0rA/rA0AyA;
217 preconPtr->precondition(zA, sA, cmpt);
220 matrix_.Amul(tA, zA, interfaceBouCoeffs_, interfaces_, cmpt);
222 const scalar tAtA =
gSumSqr(tA, matrix().
mesh().comm());
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
A class for handling words, derived from string.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
scalar gSumSqr(const UList< Type > &f, const label comm)
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Abstract base-class for lduMatrix solvers.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< scalar > mag(const dimensioned< Type > &)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
A list of keyword definitions, which are a keyword followed by any number of values (e....
scalar gSumMag(const FieldField< Field, Type > &f)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
const volScalarField & psi
lduMatrix::solver::addasymMatrixConstructorToTable< PBiCGStab > addPBiCGStabAsymMatrixConstructorToTable_
PBiCGStab(const PBiCGStab &)
Disallow default bitwise copy construct.
virtual solverPerformance solve(scalarField &psi, const scalarField &source, const direction cmpt=0) const
Solve the matrix with this solver.
scalar gSumProd(const UList< Type > &f1, const UList< Type > &f2, const label comm)
defineTypeNameAndDebug(combustionModel, 0)
A cell is defined as a list of faces with extra functionality.