46 const word& solverName,
47 const word& fieldName,
72 auto* ctorPtr = symMatrixConstructorTable(solverName);
79 "symmetric matrix solver",
81 *symMatrixConstructorTablePtr_
85 return autoPtr<lduMatrix::solver>
98 else if (
matrix.asymmetric())
100 auto* ctorPtr = asymMatrixConstructorTable(solverName);
107 "asymmetric matrix solver",
109 *asymMatrixConstructorTablePtr_
113 return autoPtr<lduMatrix::solver>
128 <<
"cannot solve incomplete matrix, "
129 "no diagonal or off-diagonal coefficient"
138 const word& fieldName,
148 solverControls.
get<
word>(
"solver"),
163 const word& fieldName,
171 fieldName_(fieldName),
173 interfaceBouCoeffs_(interfaceBouCoeffs),
174 interfaceIntCoeffs_(interfaceIntCoeffs),
175 interfaces_(interfaces),
176 controlDict_(solverControls),
185 profiling_(
"lduMatrix::solver." + fieldName)
202 controlDict_.readIfPresent(
"log", log_);
255 matrix_.sumA(tmpField, interfaceBouCoeffs_, interfaces_);
262 (
mag(Apsi - tmpField) +
mag(source - tmpField))(),
263 matrix_.mesh().comm()
A const Field/List wrapper with possible data conversion.
A field of fields is a PtrList of fields with reference counting.
pTraits< solveScalar >::cmptType cmptType
A non-const Field/List wrapper with possible data conversion.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
Abstract base-class for lduMatrix solvers.
const FieldField< Field, scalar > & interfaceIntCoeffs() const noexcept
const lduMatrix & matrix_
label maxIter_
Maximum number of iterations in the solver.
solver(const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces, const dictionary &solverControls)
Construct solver for given field name, matrix etc.
profilingTrigger profiling_
Profiling instrumentation.
scalar tolerance_
Final convergence tolerance.
const lduInterfaceFieldPtrsList & interfaces() const noexcept
lduInterfaceFieldPtrsList interfaces_
const FieldField< Field, scalar > & interfaceBouCoeffs_
label minIter_
Minimum number of iterations in the solver.
const lduMatrix & matrix() const noexcept
solveScalarField::cmptType normFactor(const solveScalarField &psi, const solveScalarField &source, const solveScalarField &Apsi, solveScalarField &tmpField, const lduMatrix::normTypes normType) const
Return the matrix norm using the specified norm method.
lduMatrix::normTypes normType_
The normalisation type.
static autoPtr< solver > New(const word &solverName, const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces, const dictionary &solverControls)
Return a new solver of given type.
int log_
Verbosity level for solver output statements.
virtual solverPerformance scalarSolve(solveScalarField &psi, const solveScalarField &source, const direction cmpt=0) const
Solve with given field and rhs (in solveScalar precision).
scalar relTol_
Convergence tolerance relative to the initial.
const FieldField< Field, scalar > & interfaceIntCoeffs_
virtual void readControls()
Read the control parameters from controlDict_.
virtual void read(const dictionary &)
Read and reset the solver parameters from the given stream.
const FieldField< Field, scalar > & interfaceBouCoeffs() const noexcept
dictionary controlDict_
Dictionary of solution controls.
const word & fieldName() const noexcept
lduMatrix(const lduMesh &mesh)
Construct (without coefficients) for an LDU addressed mesh.
normTypes
Enumerated matrix normalisation types.
@ DEFAULT_NORM
"default" norm (== L1_scaled)
@ L1_SCALED_NORM
"L1_scaled" norm
@ NO_NORM
"none" norm (returns 1)
bool symmetric() const noexcept
Matrix is symmetric.
bool diagonal() const noexcept
Matrix has diagonal only.
static constexpr const label defaultMaxIter
Default maximum number of iterations for solvers (1000).
static const Enum< normTypes > normTypesNames_
Names for the normTypes.
static const scalar defaultTolerance
Default (absolute) tolerance (1e-6).
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
A class for handling words, derived from Foam::string.
const volScalarField & psi
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
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 gSum(const FieldField< Field, Type > &f)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< solveScalar > solveScalarField
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
static constexpr const zero Zero
Global zero (0).
errorManipArg< error, int > exit(error &err, const int errNo=1)
UPtrList< const lduInterfaceField > lduInterfaceFieldPtrsList
List of coupled interface fields to be used in coupling.
SolverPerformance< scalar > solverPerformance
SolverPerformance instantiated for a scalar.
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.