lduMatrixSolver.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "lduMatrix.H"
27 #include "diagonalSolver.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
35 }
36 
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
41 (
42  const word& fieldName,
43  const lduMatrix& matrix,
44  const FieldField<Field, scalar>& interfaceBouCoeffs,
45  const FieldField<Field, scalar>& interfaceIntCoeffs,
46  const lduInterfaceFieldPtrsList& interfaces,
47  const dictionary& solverControls
48 )
49 {
50  const word name(solverControls.lookup("solver"));
51 
52  if (matrix.diagonal())
53  {
55  (
56  new diagonalSolver
57  (
58  fieldName,
59  matrix,
62  interfaces,
63  solverControls
64  )
65  );
66  }
67  else if (matrix.symmetric())
68  {
69  symMatrixConstructorTable::iterator constructorIter =
70  symMatrixConstructorTablePtr_->find(name);
71 
72  if (constructorIter == symMatrixConstructorTablePtr_->end())
73  {
74  FatalIOErrorInFunction(solverControls)
75  << "Unknown symmetric matrix solver " << name << nl << nl
76  << "Valid symmetric matrix solvers are :" << endl
77  << symMatrixConstructorTablePtr_->sortedToc()
78  << exit(FatalIOError);
79  }
80 
82  (
83  constructorIter()
84  (
85  fieldName,
86  matrix,
89  interfaces,
90  solverControls
91  )
92  );
93  }
94  else if (matrix.asymmetric())
95  {
96  asymMatrixConstructorTable::iterator constructorIter =
97  asymMatrixConstructorTablePtr_->find(name);
98 
99  if (constructorIter == asymMatrixConstructorTablePtr_->end())
100  {
101  FatalIOErrorInFunction(solverControls)
102  << "Unknown asymmetric matrix solver " << name << nl << nl
103  << "Valid asymmetric matrix solvers are :" << endl
104  << asymMatrixConstructorTablePtr_->sortedToc()
105  << exit(FatalIOError);
106  }
107 
109  (
110  constructorIter()
111  (
112  fieldName,
113  matrix,
116  interfaces,
117  solverControls
118  )
119  );
120  }
121  else
122  {
123  FatalIOErrorInFunction(solverControls)
124  << "cannot solve incomplete matrix, "
125  "no diagonal or off-diagonal coefficient"
126  << exit(FatalIOError);
127 
128  return autoPtr<lduMatrix::solver>(nullptr);
129  }
130 }
131 
132 
133 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
134 
136 (
137  const word& fieldName,
138  const lduMatrix& matrix,
139  const FieldField<Field, scalar>& interfaceBouCoeffs,
140  const FieldField<Field, scalar>& interfaceIntCoeffs,
141  const lduInterfaceFieldPtrsList& interfaces,
142  const dictionary& solverControls
143 )
144 :
145  fieldName_(fieldName),
146  matrix_(matrix),
147  interfaceBouCoeffs_(interfaceBouCoeffs),
148  interfaceIntCoeffs_(interfaceIntCoeffs),
149  interfaces_(interfaces),
150  controlDict_(solverControls)
151 {
152  readControls();
153 }
154 
155 
156 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
157 
159 {
160  maxIter_ = controlDict_.lookupOrDefault<label>("maxIter", defaultMaxIter_);
161  minIter_ = controlDict_.lookupOrDefault<label>("minIter", 0);
162  tolerance_ = controlDict_.lookupOrDefault<scalar>("tolerance", 1e-6);
163  relTol_ = controlDict_.lookupOrDefault<scalar>("relTol", 0);
164 }
165 
166 
167 void Foam::lduMatrix::solver::read(const dictionary& solverControls)
168 {
169  controlDict_ = solverControls;
170  readControls();
171 }
172 
173 
175 (
176  const scalarField& psi,
177  const scalarField& source,
178  const scalarField& Apsi,
179  scalarField& tmpField
180 ) const
181 {
182  // --- Calculate A dot reference value of psi
183  matrix_.sumA(tmpField, interfaceBouCoeffs_, interfaces_);
184 
185  tmpField *= gAverage(psi, matrix_.lduMesh_.comm());
186 
187  return
188  gSum
189  (
190  (mag(Apsi - tmpField) + mag(source - tmpField))(),
191  matrix_.lduMesh_.comm()
192  )
194 
195  // At convergence this simpler method is equivalent to the above
196  // return 2*gSumMag(source) + solverPerformance::small_;
197 }
198 
199 
200 // ************************************************************************* //
Generic field type.
Definition: FieldField.H:77
static const scalar small_
Small Type for the use in solvers.
Foam::diagonalSolver.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:860
Abstract base-class for lduMatrix solvers.
Definition: lduMatrix.H:94
solver(const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces, const dictionary &solverControls)
const lduInterfaceFieldPtrsList & interfaces() const
Definition: lduMatrix.H:241
const FieldField< Field, scalar > & interfaceIntCoeffs() const
Definition: lduMatrix.H:236
scalar normFactor(const scalarField &psi, const scalarField &source, const scalarField &Apsi, scalarField &tmpField) const
Return the matrix norm used to normalise the residual for the.
const word & fieldName() const
Definition: lduMatrix.H:221
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.
virtual void readControls()
Read the control parameters from the controlDict_.
virtual void read(const dictionary &)
Read and reset the solver parameters from the given stream.
const lduMatrix & matrix() const
Definition: lduMatrix.H:226
const FieldField< Field, scalar > & interfaceBouCoeffs() const
Definition: lduMatrix.H:231
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition: lduMatrix.H:80
bool symmetric() const
Definition: lduMatrix.H:609
bool diagonal() const
Definition: lduMatrix.H:595
bool asymmetric() const
Definition: lduMatrix.H:623
A class for handling words, derived from string.
Definition: word.H:62
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
const volScalarField & psi
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Type gSum(const FieldField< Field, Type > &f)
const doubleScalar e
Definition: doubleScalar.H:105
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
dimensioned< scalar > mag(const dimensioned< Type > &)
IOerror FatalIOError
Type gAverage(const FieldField< Field, Type > &f)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static const char nl
Definition: Ostream.H:260