Loading...
Searching...
No Matches
fvMatrixSolve.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2016-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "LduMatrix.H"
31#include "profiling.H"
32#include "PrecisionAdaptor.H"
33
34// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
35
36template<class Type>
38(
39 const label patchi,
40 const label facei,
41 const direction cmpt,
42 const scalar value
43)
44{
45 if (psi_.needReference())
46 {
47 if (Pstream::master())
48 {
49 internalCoeffs_[patchi][facei].component(cmpt) +=
50 diag()[psi_.mesh().boundary()[patchi].faceCells()[facei]];
51
52 boundaryCoeffs_[patchi][facei].component(cmpt) +=
53 diag()[psi_.mesh().boundary()[patchi].faceCells()[facei]]
54 *value;
55 }
56 }
57}
58
59
60template<class Type>
62(
63 const dictionary& solverControls
64)
65{
66 word regionName;
67 if (psi_.mesh().name() != polyMesh::defaultRegion)
68 {
69 regionName = psi_.mesh().name() + "::";
70 }
71 addProfiling(solve, "fvMatrix::solve.", regionName, psi_.name());
72
73 if (debug)
74 {
75 Info.masterStream(this->mesh().comm())
76 << "fvMatrix<Type>::solveSegregatedOrCoupled"
77 "(const dictionary& solverControls) : "
78 "solving fvMatrix<Type>"
79 << endl;
80 }
81
82 // Do not solve if maxIter == 0
83 if (solverControls.getOrDefault<label>("maxIter", -1) == 0)
84 {
85 return SolverPerformance<Type>();
86 }
87
88 word type(solverControls.getOrDefault<word>("type", "segregated"));
89
90 if (type == "segregated")
91 {
92 return solveSegregated(solverControls);
93 }
94 else if (type == "coupled")
95 {
96 return solveCoupled(solverControls);
97 }
98 else
99 {
100 FatalIOErrorInFunction(solverControls)
101 << "Unknown type " << type
102 << "; currently supported solver types are segregated and coupled"
103 << exit(FatalIOError);
106 }
107}
108
109
110template<class Type>
112(
113 const dictionary& solverControls
114)
115{
116 if (useImplicit_)
117 {
119 << "Implicit option is not allowed for type: " << Type::typeName
120 << exit(FatalError);
121 }
122
123 if (debug)
124 {
125 Info.masterStream(this->mesh().comm())
126 << "fvMatrix<Type>::solveSegregated"
127 "(const dictionary& solverControls) : "
128 "solving fvMatrix<Type>"
129 << endl;
130 }
131
132 const int logLevel =
133 solverControls.getOrDefault<int>
134 (
135 "log",
136 SolverPerformance<Type>::debug
137 );
138
139 auto& psi = psi_.constCast();
140
141 SolverPerformance<Type> solverPerfVec
142 (
143 "fvMatrix<Type>::solveSegregated",
144 psi.name()
145 );
146
147 scalarField saveDiag(diag());
148
149 Field<Type> source(source_);
150
151 // At this point include the boundary source from the coupled boundaries.
152 // This is corrected for the implicit part by updateMatrixInterfaces within
153 // the component loop.
154 addBoundarySource(source);
155
156 typename Type::labelType validComponents
157 (
158 psi.mesh().template validComponents<Type>()
159 );
160
161 for (direction cmpt=0; cmpt<Type::nComponents; cmpt++)
162 {
163 if (validComponents[cmpt] == -1) continue;
164
165 // copy field and source
166
167 scalarField psiCmpt(psi.primitiveField().component(cmpt));
168 addBoundaryDiag(diag(), cmpt);
169
170 scalarField sourceCmpt(source.component(cmpt));
171
172 FieldField<Field, scalar> bouCoeffsCmpt
173 (
174 boundaryCoeffs_.component(cmpt)
175 );
176
177 FieldField<Field, scalar> intCoeffsCmpt
178 (
179 internalCoeffs_.component(cmpt)
180 );
181
182 lduInterfaceFieldPtrsList interfaces =
183 psi.boundaryField().scalarInterfaces();
184
185 // Use the initMatrixInterfaces and updateMatrixInterfaces to correct
186 // bouCoeffsCmpt for the explicit part of the coupled boundary
187 // conditions
188 {
189 PrecisionAdaptor<solveScalar, scalar> sourceCmpt_ss(sourceCmpt);
190 ConstPrecisionAdaptor<solveScalar, scalar> psiCmpt_ss(psiCmpt);
191
192 const label startRequest = UPstream::nRequests();
193
194 initMatrixInterfaces
195 (
196 true,
197 bouCoeffsCmpt,
198 interfaces,
199 psiCmpt_ss(),
200 sourceCmpt_ss.ref(),
201 cmpt
202 );
203
204 updateMatrixInterfaces
205 (
206 true,
207 bouCoeffsCmpt,
208 interfaces,
209 psiCmpt_ss(),
210 sourceCmpt_ss.ref(),
211 cmpt,
212 startRequest
213 );
214 }
215
216 solverPerformance solverPerf;
217
218 // Solver call
219 solverPerf = lduMatrix::solver::New
220 (
221 psi.name() + pTraits<Type>::componentNames[cmpt],
222 *this,
223 bouCoeffsCmpt,
224 intCoeffsCmpt,
225 interfaces,
226 solverControls
227 )->solve(psiCmpt, sourceCmpt, cmpt);
228
229 if (logLevel)
230 {
231 solverPerf.print(Info.masterStream(this->mesh().comm()));
232 }
233
234 solverPerfVec.replace(cmpt, solverPerf);
235 solverPerfVec.solverName() = solverPerf.solverName();
236
237 psi.primitiveFieldRef().replace(cmpt, psiCmpt);
238 diag() = saveDiag;
239 }
240
241 psi.correctBoundaryConditions();
242
243 psi.mesh().data().setSolverPerformance(psi.name(), solverPerfVec);
244
245 return solverPerfVec;
246}
247
248
249template<class Type>
251(
252 const dictionary& solverControls
253)
254{
255 if (debug)
256 {
257 Info.masterStream(this->mesh().comm())
258 << "fvMatrix<Type>::solveCoupled"
259 "(const dictionary& solverControls) : "
260 "solving fvMatrix<Type>"
261 << endl;
262 }
263
264 const int logLevel =
265 solverControls.getOrDefault<int>
266 (
267 "log",
268 SolverPerformance<Type>::debug
269 );
270
271 auto& psi = psi_.constCast();
272
273 LduMatrix<Type, scalar, scalar> coupledMatrix(psi.mesh());
274 coupledMatrix.diag() = diag();
275 coupledMatrix.upper() = upper();
276 coupledMatrix.lower() = lower();
277 coupledMatrix.source() = source();
278
279 addBoundaryDiag(coupledMatrix.diag(), 0);
280 addBoundarySource(coupledMatrix.source(), false);
281
282 coupledMatrix.interfaces() = psi.boundaryFieldRef().interfaces();
283 coupledMatrix.interfacesUpper() = boundaryCoeffs().component(0);
284 coupledMatrix.interfacesLower() = internalCoeffs().component(0);
285
286 autoPtr<typename LduMatrix<Type, scalar, scalar>::solver>
287 coupledMatrixSolver
288 (
289 LduMatrix<Type, scalar, scalar>::solver::New
290 (
291 psi.name(),
292 coupledMatrix,
293 solverControls
294 )
295 );
296
297 SolverPerformance<Type> solverPerf
298 (
299 coupledMatrixSolver->solve(psi)
300 );
301
302 if (logLevel)
303 {
304 solverPerf.print(Info.masterStream(this->mesh().comm()));
305 }
306
307 psi.correctBoundaryConditions();
308
309 psi.mesh().data().setSolverPerformance(psi.name(), solverPerf);
310
311 return solverPerf;
312}
313
314
315template<class Type>
324(
325 const dictionary& solverControls
326)
328 return psi_.mesh().solve(*this, solverControls);
329}
330
331
332template<class Type>
338
339
340template<class Type>
345
346
347template<class Type>
352
353
354template<class Type>
356{
357 return this->solve(solverDict());
358}
359
360
361template<class Type>
363{
364 auto tres = tmp<Field<Type>>::New(source_);
365 auto& res = tres.ref();
366
368
369 // Loop over field components
370 for (direction cmpt=0; cmpt<Type::nComponents; cmpt++)
371 {
372 scalarField psiCmpt(psi_.primitiveField().component(cmpt));
373
374 scalarField boundaryDiagCmpt(psi_.size(), Zero);
375 addBoundaryDiag(boundaryDiagCmpt, cmpt);
376
377 FieldField<Field, scalar> bouCoeffsCmpt
378 (
379 boundaryCoeffs_.component(cmpt)
380 );
381
382 res.replace
383 (
384 cmpt,
386 (
387 psiCmpt,
388 res.component(cmpt) - boundaryDiagCmpt*psiCmpt,
389 bouCoeffsCmpt,
391 cmpt
392 )
393 );
394 }
395
396 return tres;
397}
398
399
400// ************************************************************************* //
A const Field/List wrapper with possible data conversion.
A field of fields is a PtrList of fields with reference counting.
Definition FieldField.H:77
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition Field.C:607
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,...
Definition LduMatrix.H:84
const Field< LUType > & upper() const
Definition LduMatrix.C:178
const FieldField< Field, LUType > & interfacesUpper() const noexcept
Definition LduMatrix.H:661
const Field< Type > & source() const
Definition LduMatrix.C:267
const Field< DType > & diag() const
Definition LduMatrix.C:151
const LduInterfaceFieldPtrsList< Type > & interfaces() const noexcept
Const access to the interfaces.
Definition LduMatrix.H:632
const FieldField< Field, LUType > & interfacesLower() const noexcept
Definition LduMatrix.H:666
const Field< LUType > & lower() const
Definition LduMatrix.C:223
A non-const Field/List wrapper with possible data conversion.
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
void replace(const label cmpt, const SolverPerformance< typename pTraits< Type >::cmptType > &sp)
Replace component based on the minimal SolverPerformance.
void print(Ostream &os) const
Print summary of solver performance to the given stream.
const word & solverName() const noexcept
Return solver name.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
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.
Definition UPstream.H:1714
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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
Definition fvMatrix.H:549
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
Definition fvMatrix.H:567
SolverPerformance< Type > solveSegregatedOrCoupled()
Solve segregated or coupled returning the solution statistics.
void addBoundarySource(Field< Type > &source, const bool couples=true) const
Definition fvMatrix.C:169
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.
Definition fvMatrix.C:1513
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
Definition fvMatrix.C:116
SolverPerformance< Type > solveSegregated(const dictionary &)
Solve segregated returning the solution statistics.
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition fvMatrix.H:487
Field< Type > & source() noexcept
Definition fvMatrix.H:535
const dictionary & solverDict() const
Return the solver dictionary for psi, taking into account finalIteration.
Definition fvMatrix.C:1522
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
Definition lduMatrix.C:195
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
Definition lduMatrix.C:235
const scalarField & lower() const
Definition lduMatrix.C:306
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.
Definition pTraits.H:64
static word defaultRegion
Return the default region name.
Definition polyMesh.H:406
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition refPtrI.H:230
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
const volScalarField & psi
dynamicFvMesh & mesh
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Namespace for handling debugging switches.
Definition debug.C:45
type
Types of root.
Definition Roots.H:53
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.
Definition Ostream.H:519
uint8_t direction
Definition direction.H:49
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).
Definition zero.H:127
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.
Definition exprTraits.C:127
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
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 ...
CEqn solve()