38 lduMatrix::solver::addsymMatrixConstructorToTable<PBiCGStab>
41 lduMatrix::solver::addasymMatrixConstructorToTable<PBiCGStab>
50 const word& fieldName,
86 const label nCells =
psi.
size();
88 solveScalar* __restrict__ psiPtr =
psi.
begin();
91 solveScalar* __restrict__ pAPtr = pA.
begin();
94 solveScalar* __restrict__ yAPtr = yA.
begin();
97 matrix_.Amul(yA,
psi, interfaceBouCoeffs_, interfaces_, cmpt);
101 solveScalar* __restrict__ rAPtr = rA.
begin();
103 matrix().setResidualField
111 const solveScalar normFactor = this->normFactor(
psi, source, yA, pA);
113 if ((log_ >= 2) || (lduMatrix::debug >= 2))
115 Info<<
" Normalisation factor = " << normFactor <<
endl;
132 solveScalar* __restrict__ AyAPtr = AyA.
begin();
135 solveScalar* __restrict__ sAPtr = sA.
begin();
138 solveScalar* __restrict__ zAPtr = zA.
begin();
141 solveScalar* __restrict__ tAPtr = tA.
begin();
147 solveScalar rA0rA = 0;
148 solveScalar
alpha = 0;
149 solveScalar omega = 0;
163 const solveScalar rA0rAold = rA0rA;
189 const solveScalar
beta = (rA0rA/rA0rAold)*(
alpha/omega);
199 preconPtr->precondition(yA, pA, cmpt);
202 matrix_.Amul(AyA, yA, interfaceBouCoeffs_, interfaces_, cmpt);
204 const solveScalar rA0AyA =
207 alpha = rA0rA/rA0AyA;
236 preconPtr->precondition(zA, sA, cmpt);
239 matrix_.Amul(tA, zA, interfaceBouCoeffs_, interfaces_, cmpt);
241 const solveScalar tAtA =
gSumSqr(tA, matrix().
mesh().comm());
267 matrix().setResidualField
A const Field/List wrapper with possible data conversion.
A field of fields is a PtrList of fields with reference counting.
const word & getName() const
Get name.
Preconditioned bi-conjugate gradient stabilized solver for asymmetric lduMatrices using a run-time se...
virtual solverPerformance scalarSolve(solveScalarField &psi, const solveScalarField &source, const direction cmpt=0) const
Solve the matrix with this solver.
A non-const Field/List wrapper with possible data conversion.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
iterator begin() noexcept
Return an iterator to begin traversing the UList.
void size(const label n)
Older name for setAddressableSize.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A cell is defined as a list of faces with extra functionality.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Base class for solution control classes.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
const volScalarField & psi
scalarProduct< Type, Type >::type gSumProd(const UList< Type > &f1, const UList< Type > &f2, const label comm)
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
lduMatrix::solver::addasymMatrixConstructorToTable< PBiCGStab > addPBiCGStabAsymMatrixConstructorToTable_
typeOfMag< Type >::type gSumMag(const FieldField< Field, Type > &f)
outerProduct1< Type >::type gSumSqr(const UList< Type > &f, const label comm)
lduMatrix::solver::addsymMatrixConstructorToTable< PBiCGStab > addPBiCGStabSymMatrixConstructorToTable_
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)