45 if (psi_.needReference())
49 internalCoeffs_[patchi][facei].component(cmpt) +=
50 diag()[psi_.mesh().boundary()[patchi].faceCells()[facei]];
52 boundaryCoeffs_[patchi][facei].component(cmpt) +=
53 diag()[psi_.mesh().boundary()[patchi].faceCells()[facei]]
67 if (psi_.mesh().name() != polyMesh::defaultRegion)
76 <<
"fvMatrix<Type>::solveSegregatedOrCoupled"
77 "(const dictionary& solverControls) : "
78 "solving fvMatrix<Type>"
83 if (solverControls.getOrDefault<label>(
"maxIter", -1) == 0)
85 return SolverPerformance<Type>();
88 word
type(solverControls.getOrDefault<word>(
"type",
"segregated"));
90 if (type ==
"segregated")
92 return solveSegregated(solverControls);
94 else if (type ==
"coupled")
96 return solveCoupled(solverControls);
101 <<
"Unknown type " <<
type
102 <<
"; currently supported solver types are segregated and coupled"
103 <<
exit(FatalIOError);
119 <<
"Implicit option is not allowed for type: " << Type::typeName
126 <<
"fvMatrix<Type>::solveSegregated"
127 "(const dictionary& solverControls) : "
128 "solving fvMatrix<Type>"
133 solverControls.getOrDefault<
int>
136 SolverPerformance<Type>::debug
139 auto&
psi = psi_.constCast();
141 SolverPerformance<Type> solverPerfVec
143 "fvMatrix<Type>::solveSegregated",
149 Field<Type> source(source_);
154 addBoundarySource(source);
156 typename Type::labelType validComponents
158 psi.mesh().template validComponents<Type>()
161 for (direction cmpt=0; cmpt<Type::nComponents; cmpt++)
163 if (validComponents[cmpt] == -1)
continue;
168 addBoundaryDiag(
diag(), cmpt);
172 FieldField<Field, scalar> bouCoeffsCmpt
174 boundaryCoeffs_.component(cmpt)
177 FieldField<Field, scalar> intCoeffsCmpt
179 internalCoeffs_.component(cmpt)
183 psi.boundaryField().scalarInterfaces();
189 PrecisionAdaptor<solveScalar, scalar> sourceCmpt_ss(sourceCmpt);
190 ConstPrecisionAdaptor<solveScalar, scalar> psiCmpt_ss(psiCmpt);
192 const label startRequest = UPstream::nRequests();
204 updateMatrixInterfaces
219 solverPerf = lduMatrix::solver::New
221 psi.name() + pTraits<Type>::componentNames[cmpt],
227 )->solve(psiCmpt, sourceCmpt, cmpt);
234 solverPerfVec.
replace(cmpt, solverPerf);
235 solverPerfVec.solverName() = solverPerf.solverName();
237 psi.primitiveFieldRef().replace(cmpt, psiCmpt);
241 psi.correctBoundaryConditions();
243 psi.mesh().data().setSolverPerformance(
psi.name(), solverPerfVec);
245 return solverPerfVec;
257 Info.masterStream(this->
mesh().comm())
258 <<
"fvMatrix<Type>::solveCoupled"
259 "(const dictionary& solverControls) : "
260 "solving fvMatrix<Type>"
268 SolverPerformance<Type>::debug
271 auto&
psi = psi_.constCast();
273 LduMatrix<Type, scalar, scalar> coupledMatrix(
psi.mesh());
274 coupledMatrix.diag() = diag();
275 coupledMatrix.upper() = upper();
276 coupledMatrix.lower() = lower();
277 coupledMatrix.source() = source();
279 addBoundaryDiag(coupledMatrix.diag(), 0);
280 addBoundarySource(coupledMatrix.source(),
false);
282 coupledMatrix.interfaces() =
psi.boundaryFieldRef().interfaces();
283 coupledMatrix.interfacesUpper() = boundaryCoeffs().component(0);
284 coupledMatrix.interfacesLower() = internalCoeffs().component(0);
286 autoPtr<typename LduMatrix<Type, scalar, scalar>::solver>
289 LduMatrix<Type, scalar, scalar>::solver::New
297 SolverPerformance<Type> solverPerf
299 coupledMatrixSolver->solve(
psi)
304 solverPerf.print(Info.
masterStream(this->mesh().comm()));
307 psi.correctBoundaryConditions();
309 psi.mesh().data().setSolverPerformance(
psi.name(), solverPerf);
328 return psi_.mesh().solve(*
this, solverControls);
343 return solve(fvMat_.solverDict());
365 auto& res = tres.ref();
370 for (
direction cmpt=0; cmpt<Type::nComponents; cmpt++)
379 boundaryCoeffs_.component(cmpt)
388 res.component(cmpt) - boundaryDiagCmpt*psiCmpt,
A const Field/List wrapper with possible data conversion.
A field of fields is a PtrList of fields with reference counting.
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
lduInterfaceFieldPtrsList scalarInterfaces() const
Return a list of pointers for each patch field with only those pointing to interfaces being set.
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
static autoPtr< solver > New(const word &fieldName, const LduMatrix< Type, DType, LUType > &matrix, const dictionary &solverDict)
Return a new solver.
LduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
const Field< LUType > & upper() const
const FieldField< Field, LUType > & interfacesUpper() const noexcept
const Field< Type > & source() const
const Field< DType > & diag() const
const LduInterfaceFieldPtrsList< Type > & interfaces() const noexcept
Const access to the interfaces.
const FieldField< Field, LUType > & interfacesLower() const noexcept
const Field< LUType > & lower() const
A non-const Field/List wrapper with possible data conversion.
void size(const label n)
Older name for setAddressableSize.
static label nRequests() noexcept
Number of outstanding requests (on the internal list of requests).
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
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 getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
SolverPerformance< Type > solve()
Solve returning the solution statistics.
autoPtr< fvSolver > solver()
Construct and return the solver.
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
const FieldField< Field, Type > & internalCoeffs() const noexcept
fvBoundary scalar field containing pseudo-matrix coeffs for internal cells
SolverPerformance< Type > solve()
Solve returning the solution statistics.
autoPtr< fvSolver > solver(const dictionary &)
Construct and return the solver.
const FieldField< Field, Type > & boundaryCoeffs() const noexcept
fvBoundary scalar field containing pseudo-matrix coeffs for boundary cells
SolverPerformance< Type > solveSegregatedOrCoupled()
Solve segregated or coupled returning the solution statistics.
void addBoundarySource(Field< Type > &source, const bool couples=true) const
void setComponentReference(const label patchi, const label facei, const direction cmpt, const scalar value)
Set reference level for a component of the solution on a given patch face.
const dictionary & solverDict(const word &name) const
Return the solver dictionary (from fvSolution) for name.
SolverPerformance< Type > solveSegregatedOrCoupled(const dictionary &)
Solve segregated or coupled returning the solution statistics.
SolverPerformance< Type > solveCoupled(const dictionary &)
Solve coupled returning the solution statistics.
tmp< Field< Type > > residual() const
Return the matrix residual.
void addBoundaryDiag(scalarField &diag, const direction cmpt) const
SolverPerformance< Type > solveSegregated(const dictionary &)
Solve segregated returning the solution statistics.
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Field< Type > & source() noexcept
const dictionary & solverDict() const
Return the solver dictionary for psi, taking into account finalIteration.
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.
void residual(solveScalarField &rA, const solveScalarField &psi, const scalarField &source, const FieldField< Field, scalar > &interfaceBouCoeffs, const lduInterfaceFieldPtrsList &interfaces, const direction cmpt) const
const scalarField & diag() const
void initMatrixInterfaces(const bool add, const FieldField< Field, scalar > &interfaceCoeffs, const lduInterfaceFieldPtrsList &interfaces, const solveScalarField &psiif, solveScalarField &result, const direction cmpt) const
Initialise the update of interfaced interfaces for matrix operations.
const scalarField & upper() const
const scalarField & lower() const
void updateMatrixInterfaces(const bool add, const FieldField< Field, scalar > &interfaceCoeffs, const lduInterfaceFieldPtrsList &interfaces, const solveScalarField &psiif, solveScalarField &result, const direction cmpt, const label startRequest) const
Update interfaced interfaces for matrix operations.
OSstream & masterStream(int communicator)
Return OSstream for output operations on the master process only, Snull on other processes.
A traits class, which is primarily used for primitives and vector-space.
static word defaultRegion
Return the default region name.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
const volScalarField & psi
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Namespace for handling debugging switches.
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.
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
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).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
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 addProfiling(Name,...)
Define profiling trigger with specified name and description string. The description is generated by ...