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) 2019-2022 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"
30#include "DiagonalSolver.H"
31
32// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33
34template<class Type, class DType, class LUType>
37(
38 const word& fieldName,
40 const dictionary& solverDict
41)
42{
43 const word solverName(solverDict.get<word>("solver"));
44
45 if (matrix.diagonal())
46 {
48 (
50 (
52 matrix,
53 solverDict
54 )
55 );
56 }
57 else if (matrix.symmetric())
58 {
59 auto* ctorPtr = symMatrixConstructorTable(solverName);
60
61 if (!ctorPtr)
62 {
64 (
65 solverDict,
66 "symmetric matrix solver",
67 solverName,
68 *symMatrixConstructorTablePtr_
69 ) << exit(FatalIOError);
70 }
71
73 (
74 ctorPtr
75 (
77 matrix,
78 solverDict
79 )
80 );
81 }
82 else if (matrix.asymmetric())
83 {
84 auto* ctorPtr = asymMatrixConstructorTable(solverName);
85
86 if (!ctorPtr)
87 {
89 (
90 solverDict,
91 "asymmetric matrix solver",
92 solverName,
93 *asymMatrixConstructorTablePtr_
94 ) << exit(FatalIOError);
95 }
96
98 (
99 ctorPtr
100 (
101 fieldName,
102 matrix,
103 solverDict
104 )
105 );
106 }
107
108 FatalIOErrorInFunction(solverDict)
109 << "cannot solve incomplete matrix, "
110 "no diagonal or off-diagonal coefficient"
111 << exit(FatalIOError);
112
113 return nullptr;
114}
115
116
117// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
118
119template<class Type, class DType, class LUType>
121(
122 const word& fieldName,
123 const LduMatrix<Type, DType, LUType>& matrix,
124 const dictionary& solverDict
125)
126:
127 fieldName_(fieldName),
128 matrix_(matrix),
129
130 controlDict_(solverDict),
131
132 log_(1),
133 minIter_(0),
134 maxIter_(lduMatrix::defaultMaxIter),
135 normType_(lduMatrix::normTypes::DEFAULT_NORM),
136 tolerance_(lduMatrix::defaultTolerance*pTraits<Type>::one),
137 relTol_(Zero)
138{
140}
141
142
143// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
144
145template<class Type, class DType, class LUType>
147{
148 controlDict_.readIfPresent("log", log_);
150 lduMatrix::normTypesNames_.readIfPresent("norm", controlDict_, normType_);
151 controlDict_.readIfPresent("minIter", minIter_);
152 controlDict_.readIfPresent("maxIter", maxIter_);
153 controlDict_.readIfPresent("tolerance", tolerance_);
154 controlDict_.readIfPresent("relTol", relTol_);
155}
156
157
158template<class Type, class DType, class LUType>
160(
161 const dictionary& solverDict
162)
164 controlDict_ = solverDict;
165 readControls();
166}
167
168
169template<class Type, class DType, class LUType>
171(
172 const Field<Type>& psi,
173 const Field<Type>& Apsi,
174 Field<Type>& tmpField,
175 const lduMatrix::normTypes normType
176) const
177{
178 switch (normType)
179 {
181 {
182 break;
183 }
184
187 {
188 // --- Calculate A dot reference value of psi
189 matrix_.sumA(tmpField);
190 cmptMultiply(tmpField, tmpField, gAverage(psi));
191
192 return stabilise
193 (
194 gSum
195 (
196 cmptMag(Apsi - tmpField)
197 + cmptMag(matrix_.source() - tmpField)
198 ),
200 );
201
202 // Equivalent at convergence:
203 // return stabilise
204 // (
205 // 2*gSumCmptMag(matrix_.source()), matrix_.small_
206 // );
207 break;
208 }
209 }
210
211 // Fall-through: no norm
212 return pTraits<Type>::one;
213}
214
215
216// ************************************************************************* //
Foam::DiagonalSolver.
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
label maxIter_
Maximum number of iterations in the solver.
Definition LduMatrix.H:163
solver(const word &fieldName, const LduMatrix< Type, DType, LUType > &matrix, const dictionary &solverDict)
Construct for given field name, matrix and controls.
const LduMatrix< Type, DType, LUType > & matrix() const noexcept
Definition LduMatrix.H:287
Type relTol_
Convergence tolerance relative to the initial.
Definition LduMatrix.H:178
Type normFactor(const Field< Type > &psi, const Field< Type > &Apsi, Field< Type > &tmpField, const lduMatrix::normTypes normType) const
Return the matrix norm using the specified norm method.
label minIter_
Minimum number of iterations in the solver.
Definition LduMatrix.H:158
Type tolerance_
Final convergence tolerance.
Definition LduMatrix.H:173
lduMatrix::normTypes normType_
The matrix normalisation type.
Definition LduMatrix.H:168
int log_
Verbosity level for solver output statements.
Definition LduMatrix.H:153
virtual void readControls()
Read the control parameters from controlDict_.
const LduMatrix< Type, DType, LUType > & matrix_
Definition LduMatrix.H:143
static autoPtr< solver > New(const word &fieldName, const LduMatrix< Type, DType, LUType > &matrix, const dictionary &solverDict)
Return a new solver.
virtual void read(const dictionary &)
Read and reset the solver parameters from the given dictionary.
dictionary controlDict_
Dictionary of solution controls.
Definition LduMatrix.H:148
const word & fieldName() const noexcept
Definition LduMatrix.H:282
LduMatrix(const lduMesh &mesh)
Construct given an LDU addressed mesh.
Definition LduMatrix.C:28
static const scalar small_
Small Type for the use in solvers.
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 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.
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition lduMatrix.H:81
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
static const Enum< normTypes > normTypesNames_
Names for the normTypes.
Definition lduMatrix.H:134
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition one.H:57
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
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
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
Type gSum(const FieldField< Field, Type > &f)
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
dimensionedScalar stabilise(const dimensionedScalar &x, const dimensionedScalar &y)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125