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-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 #include "PrecisionAdaptor.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineRunTimeSelectionTable(lduMatrix::solver, symMatrix);
38  defineRunTimeSelectionTable(lduMatrix::solver, asymMatrix);
39 }
40 
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
45 (
46  const word& fieldName,
47  const lduMatrix& matrix,
48  const FieldField<Field, scalar>& interfaceBouCoeffs,
49  const FieldField<Field, scalar>& interfaceIntCoeffs,
50  const lduInterfaceFieldPtrsList& interfaces,
51  const dictionary& solverControls
52 )
53 {
54  const word name(solverControls.get<word>("solver"));
55 
56  if (matrix.diagonal())
57  {
59  (
60  new diagonalSolver
61  (
62  fieldName,
63  matrix,
64  interfaceBouCoeffs,
65  interfaceIntCoeffs,
66  interfaces,
67  solverControls
68  )
69  );
70  }
71  else if (matrix.symmetric())
72  {
73  auto* ctorPtr = symMatrixConstructorTable(name);
74 
75  if (!ctorPtr)
76  {
78  (
79  solverControls,
80  "symmetric matrix solver",
81  name,
82  *symMatrixConstructorTablePtr_
83  ) << exit(FatalIOError);
84  }
85 
87  (
88  ctorPtr
89  (
90  fieldName,
91  matrix,
92  interfaceBouCoeffs,
93  interfaceIntCoeffs,
94  interfaces,
95  solverControls
96  )
97  );
98  }
99  else if (matrix.asymmetric())
100  {
101  auto* ctorPtr = asymMatrixConstructorTable(name);
102 
103  if (!ctorPtr)
104  {
106  (
107  solverControls,
108  "asymmetric matrix solver",
109  name,
110  *asymMatrixConstructorTablePtr_
111  ) << exit(FatalIOError);
112  }
113 
115  (
116  ctorPtr
117  (
118  fieldName,
119  matrix,
120  interfaceBouCoeffs,
121  interfaceIntCoeffs,
122  interfaces,
123  solverControls
124  )
125  );
126  }
127 
128  FatalIOErrorInFunction(solverControls)
129  << "cannot solve incomplete matrix, "
130  "no diagonal or off-diagonal coefficient"
131  << exit(FatalIOError);
132 
133  return nullptr;
134 }
135 
136 
137 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
138 
140 (
141  const word& fieldName,
142  const lduMatrix& matrix,
143  const FieldField<Field, scalar>& interfaceBouCoeffs,
144  const FieldField<Field, scalar>& interfaceIntCoeffs,
145  const lduInterfaceFieldPtrsList& interfaces,
146  const dictionary& solverControls
147 )
148 :
149  fieldName_(fieldName),
150  matrix_(matrix),
151  interfaceBouCoeffs_(interfaceBouCoeffs),
152  interfaceIntCoeffs_(interfaceIntCoeffs),
153  interfaces_(interfaces),
154  controlDict_(solverControls),
155  profiling_("lduMatrix::solver." + fieldName)
156 {
157  readControls();
158 }
159 
160 
161 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
162 
164 {
165  log_ = controlDict_.getOrDefault<int>("log", 1);
166  minIter_ = controlDict_.getOrDefault<label>("minIter", 0);
167  maxIter_ = controlDict_.getOrDefault<label>("maxIter", defaultMaxIter_);
168  tolerance_ = controlDict_.getOrDefault<scalar>("tolerance", 1e-6);
169  relTol_ = controlDict_.getOrDefault<scalar>("relTol", 0);
170 }
171 
172 
173 void Foam::lduMatrix::solver::read(const dictionary& solverControls)
174 {
175  controlDict_ = solverControls;
176  readControls();
177 }
178 
179 
181 (
183  const solveScalarField& source,
184  const direction cmpt
185 ) const
186 {
188  return solve
189  (
190  tpsi_s.ref(),
192  cmpt
193  );
194 }
195 
196 
198 (
199  const solveScalarField& psi,
200  const solveScalarField& source,
201  const solveScalarField& Apsi,
202  solveScalarField& tmpField
203 ) const
204 {
205  // --- Calculate A dot reference value of psi
206  matrix_.sumA(tmpField, interfaceBouCoeffs_, interfaces_);
207 
208  tmpField *= gAverage(psi, matrix_.mesh().comm());
209 
210  return
211  gSum
212  (
213  (mag(Apsi - tmpField) + mag(source - tmpField))(),
214  matrix_.mesh().comm()
215  )
217 
218  // At convergence this simpler method is equivalent to the above
219  // return 2*gSumMag(source) + solverPerformance::small_;
220 }
221 
222 
223 // ************************************************************************* //
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::FieldField
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:53
Foam::diagonalSolver
Foam::diagonalSolver.
Definition: diagonalSolver.H:53
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:604
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
Foam::lduMatrix
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition: lduMatrix.H:83
Foam::lduMatrix::symmetric
bool symmetric() const
Definition: lduMatrix.H:628
Foam::lduMatrix::solver::scalarSolve
virtual solverPerformance scalarSolve(solveScalarField &psi, const solveScalarField &source, const direction cmpt=0) const
Solve with given field and rhs (in solveScalar precision).
Definition: lduMatrixSolver.C:181
PrecisionAdaptor.H
Foam::FatalIOError
IOerror FatalIOError
lduMatrix.H
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::lduMatrix::solver::normFactor
solveScalarField::cmptType normFactor(const solveScalarField &psi, const solveScalarField &source, const solveScalarField &Apsi, solveScalarField &tmpField) const
Definition: lduMatrixSolver.C:198
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:107
Foam::lduMatrix::solver::minIter_
label minIter_
Minimum number of iterations in the solver.
Definition: lduMatrix.H:120
Foam::Field< solveScalar >::cmptType
pTraits< solveScalar >::cmptType cmptType
Component type.
Definition: Field.H:86
FatalIOErrorInLookup
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
Foam::lduMatrix::solver::read
virtual void read(const dictionary &)
Read and reset the solver parameters from the given stream.
Definition: lduMatrixSolver.C:173
Foam::Field< solveScalar >
Foam::lduMatrix::solver::maxIter_
label maxIter_
Maximum number of iterations in the solver.
Definition: lduMatrix.H:123
diagonalSolver.H
Foam::UPtrList< const lduInterfaceField >
Foam::solve
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
Solve returning the solution statistics given convergence tolerance.
Foam::lduMatrix::diagonal
bool diagonal() const
Definition: lduMatrix.H:623
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::lduMatrix::solver::readControls
virtual void readControls()
Read the control parameters from the controlDict_.
Definition: lduMatrixSolver.C:163
Foam::lduMatrix::solver::New
static autoPtr< solver > New(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.
Definition: lduMatrixSolver.C:45
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::lduMatrix::solver::relTol_
scalar relTol_
Convergence tolerance relative to the initial.
Definition: lduMatrix.H:129
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::SolverPerformance::small_
static const scalar small_
Small Type for the use in solvers.
Definition: SolverPerformance.H:108
Foam::autoPtr< Foam::lduMatrix::solver >
Foam::PrecisionAdaptor
A non-const Field/List wrapper with possible data conversion.
Definition: PrecisionAdaptor.H:186
Foam::lduMatrix::solver::controlDict_
dictionary controlDict_
Dictionary of controls.
Definition: lduMatrix.H:114
Foam::lduMatrix::solver::tolerance_
scalar tolerance_
Final convergence tolerance.
Definition: lduMatrix.H:126
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::ConstPrecisionAdaptor
A const Field/List wrapper with possible data conversion.
Definition: PrecisionAdaptor.H:57
Foam::direction
uint8_t direction
Definition: direction.H:52
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::lduMatrix::asymmetric
bool asymmetric() const
Definition: lduMatrix.H:633
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
Foam::refPtr< Container< Type > >::ref
Container< Type > & ref() const
Definition: refPtrI.H:203
Foam::lduMatrix::solver::log_
int log_
Level of verbosity in the solver output statements.
Definition: lduMatrix.H:117
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
psi
const volScalarField & psi
Definition: createFieldRefs.H:1
Foam::lduMatrix::solver::defaultMaxIter_
static const label defaultMaxIter_
Default maximum number of iterations in the solver.
Definition: lduMatrix.H:105
Foam::SolverPerformance
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
Definition: SolverPerformance.H:52
Foam::lduMatrix::solver::solver
solver(const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces, const dictionary &solverControls)
Definition: lduMatrixSolver.C:140