Loading...
Searching...
No Matches
lduMatrixSolver.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-2017 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"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
39}
40
41
42// * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
43
45(
46 const word& solverName,
47 const word& fieldName,
48 const lduMatrix& matrix,
49 const FieldField<Field, scalar>& interfaceBouCoeffs,
50 const FieldField<Field, scalar>& interfaceIntCoeffs,
51 const lduInterfaceFieldPtrsList& interfaces,
52 const dictionary& solverControls
53)
54{
55 if (matrix.diagonal())
56 {
58 (
60 (
62 matrix,
66 solverControls
67 )
68 );
69 }
70 else if (matrix.symmetric())
71 {
72 auto* ctorPtr = symMatrixConstructorTable(solverName);
73
74 if (!ctorPtr)
75 {
77 (
78 solverControls,
79 "symmetric matrix solver",
80 solverName,
81 *symMatrixConstructorTablePtr_
82 ) << exit(FatalIOError);
83 }
84
85 return autoPtr<lduMatrix::solver>
86 (
87 ctorPtr
88 (
90 matrix,
94 solverControls
95 )
96 );
97 }
98 else if (matrix.asymmetric())
99 {
100 auto* ctorPtr = asymMatrixConstructorTable(solverName);
101
102 if (!ctorPtr)
103 {
105 (
106 solverControls,
107 "asymmetric matrix solver",
108 solverName,
109 *asymMatrixConstructorTablePtr_
110 ) << exit(FatalIOError);
111 }
112
113 return autoPtr<lduMatrix::solver>
114 (
115 ctorPtr
116 (
117 fieldName,
118 matrix,
122 solverControls
123 )
124 );
125 }
126
127 FatalIOErrorInFunction(solverControls)
128 << "cannot solve incomplete matrix, "
129 "no diagonal or off-diagonal coefficient"
130 << exit(FatalIOError);
131
132 return nullptr;
133}
134
135
137(
138 const word& fieldName,
139 const lduMatrix& matrix,
140 const FieldField<Field, scalar>& interfaceBouCoeffs,
141 const FieldField<Field, scalar>& interfaceIntCoeffs,
142 const lduInterfaceFieldPtrsList& interfaces,
143 const dictionary& solverControls
144)
145{
146 return New
147 (
148 solverControls.get<word>("solver"),
149 fieldName,
150 matrix,
151 interfaceBouCoeffs,
152 interfaceIntCoeffs,
153 interfaces,
154 solverControls
155 );
156}
157
158
159// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
160
162(
163 const word& fieldName,
164 const lduMatrix& matrix,
165 const FieldField<Field, scalar>& interfaceBouCoeffs,
166 const FieldField<Field, scalar>& interfaceIntCoeffs,
167 const lduInterfaceFieldPtrsList& interfaces,
168 const dictionary& solverControls
169)
170:
171 fieldName_(fieldName),
172 matrix_(matrix),
173 interfaceBouCoeffs_(interfaceBouCoeffs),
174 interfaceIntCoeffs_(interfaceIntCoeffs),
175 interfaces_(interfaces),
176 controlDict_(solverControls),
177
178 log_(1),
179 minIter_(0),
180 maxIter_(lduMatrix::defaultMaxIter),
181 normType_(lduMatrix::normTypes::DEFAULT_NORM),
182 tolerance_(lduMatrix::defaultTolerance),
183 relTol_(Zero),
184
185 profiling_("lduMatrix::solver." + fieldName)
187 readControls();
188}
189
190
191// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
192
194{
195 log_ = 1;
196 minIter_ = 0;
197 maxIter_ = lduMatrix::defaultMaxIter;
199 tolerance_ = lduMatrix::defaultTolerance;
200 relTol_ = 0;
201
202 controlDict_.readIfPresent("log", log_);
203 lduMatrix::normTypesNames_.readIfPresent("norm", controlDict_, normType_);
204 controlDict_.readIfPresent("minIter", minIter_);
205 controlDict_.readIfPresent("maxIter", maxIter_);
206 controlDict_.readIfPresent("tolerance", tolerance_);
207 controlDict_.readIfPresent("relTol", relTol_);
208}
209
210
211void Foam::lduMatrix::solver::read(const dictionary& solverControls)
212{
213 controlDict_ = solverControls;
214 readControls();
215}
216
217
219(
221 const solveScalarField& source,
222 const direction cmpt
223) const
224{
226 return solve
227 (
228 tpsi_s.ref(),
230 cmpt
231 );
232}
233
234
236(
237 const solveScalarField& psi,
238 const solveScalarField& source,
239 const solveScalarField& Apsi,
240 solveScalarField& tmpField,
241 const lduMatrix::normTypes normType
242) const
243{
244 switch (normType)
245 {
247 {
248 break;
249 }
250
253 {
254 // --- Calculate A dot reference value of psi
255 matrix_.sumA(tmpField, interfaceBouCoeffs_, interfaces_);
256
257 tmpField *= gAverage(psi, matrix_.mesh().comm());
258
259 return
260 gSum
261 (
262 (mag(Apsi - tmpField) + mag(source - tmpField))(),
263 matrix_.mesh().comm()
265
266 // Equivalent at convergence:
267 // return 2*gSumMag(source) + solverPerformance::small_;
268 break;
269 }
270 }
271
272 // Fall-through: no norm
274}
275
276
277// ************************************************************************* //
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
pTraits< solveScalar >::cmptType cmptType
Definition Field.H:178
A non-const Field/List wrapper with possible data conversion.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
Foam::diagonalSolver.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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.
Definition lduMatrix.H:152
const FieldField< Field, scalar > & interfaceIntCoeffs() const noexcept
Definition lduMatrix.H:338
const lduMatrix & matrix_
Definition lduMatrix.H:158
label maxIter_
Maximum number of iterations in the solver.
Definition lduMatrix.H:181
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.
Definition lduMatrix.H:201
scalar tolerance_
Final convergence tolerance.
Definition lduMatrix.H:191
const lduInterfaceFieldPtrsList & interfaces() const noexcept
Definition lduMatrix.H:343
lduInterfaceFieldPtrsList interfaces_
Definition lduMatrix.H:161
const FieldField< Field, scalar > & interfaceBouCoeffs_
Definition lduMatrix.H:159
label minIter_
Minimum number of iterations in the solver.
Definition lduMatrix.H:176
const lduMatrix & matrix() const noexcept
Definition lduMatrix.H:328
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.
Definition lduMatrix.H:186
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.
Definition lduMatrix.H:171
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.
Definition lduMatrix.H:196
const FieldField< Field, scalar > & interfaceIntCoeffs_
Definition lduMatrix.H:160
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
Definition lduMatrix.H:333
dictionary controlDict_
Dictionary of solution controls.
Definition lduMatrix.H:166
const word & fieldName() const noexcept
Definition lduMatrix.H:323
lduMatrix(const lduMesh &mesh)
Construct (without coefficients) for an LDU addressed mesh.
Definition lduMatrix.C:54
normTypes
Enumerated matrix normalisation types.
Definition lduMatrix.H:125
@ DEFAULT_NORM
"default" norm (== L1_scaled)
Definition lduMatrix.H:127
@ L1_SCALED_NORM
"L1_scaled" norm
Definition lduMatrix.H:128
@ NO_NORM
"none" norm (returns 1)
Definition lduMatrix.H:126
bool symmetric() const noexcept
Matrix is symmetric.
Definition lduMatrix.H:842
bool diagonal() const noexcept
Matrix has diagonal only.
Definition lduMatrix.H:834
static constexpr const label defaultMaxIter
Default maximum number of iterations for solvers (1000).
Definition lduMatrix.H:139
static const Enum< normTypes > normTypesNames_
Names for the normTypes.
Definition lduMatrix.H:134
static const scalar defaultTolerance
Default (absolute) tolerance (1e-6).
Definition lduMatrix.H:144
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition refPtrI.H:230
A class for handling words, derived from Foam::string.
Definition word.H:66
const volScalarField & psi
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition error.H:637
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
Namespace for OpenFOAM.
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
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
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 defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
CEqn solve()