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 | Copyright (C) 2011-2013 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 {
33  defineRunTimeSelectionTable(lduMatrix::solver, symMatrix);
34  defineRunTimeSelectionTable(lduMatrix::solver, asymMatrix);
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,
60  interfaceBouCoeffs,
61  interfaceIntCoeffs,
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  {
75  (
76  "lduMatrix::solver::New", solverControls
77  ) << "Unknown symmetric matrix solver " << name << nl << nl
78  << "Valid symmetric matrix solvers are :" << endl
79  << symMatrixConstructorTablePtr_->sortedToc()
80  << exit(FatalIOError);
81  }
82 
84  (
85  constructorIter()
86  (
87  fieldName,
88  matrix,
89  interfaceBouCoeffs,
90  interfaceIntCoeffs,
91  interfaces,
92  solverControls
93  )
94  );
95  }
96  else if (matrix.asymmetric())
97  {
98  asymMatrixConstructorTable::iterator constructorIter =
99  asymMatrixConstructorTablePtr_->find(name);
100 
101  if (constructorIter == asymMatrixConstructorTablePtr_->end())
102  {
104  (
105  "lduMatrix::solver::New", solverControls
106  ) << "Unknown asymmetric matrix solver " << name << nl << nl
107  << "Valid asymmetric matrix solvers are :" << endl
108  << asymMatrixConstructorTablePtr_->sortedToc()
109  << exit(FatalIOError);
110  }
111 
113  (
114  constructorIter()
115  (
116  fieldName,
117  matrix,
118  interfaceBouCoeffs,
119  interfaceIntCoeffs,
120  interfaces,
121  solverControls
122  )
123  );
124  }
125  else
126  {
128  (
129  "lduMatrix::solver::New", solverControls
130  ) << "cannot solve incomplete matrix, "
131  "no diagonal or off-diagonal coefficient"
132  << exit(FatalIOError);
133 
134  return autoPtr<lduMatrix::solver>(NULL);
135  }
136 }
137 
138 
139 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
140 
142 (
143  const word& fieldName,
144  const lduMatrix& matrix,
145  const FieldField<Field, scalar>& interfaceBouCoeffs,
146  const FieldField<Field, scalar>& interfaceIntCoeffs,
147  const lduInterfaceFieldPtrsList& interfaces,
148  const dictionary& solverControls
149 )
150 :
151  fieldName_(fieldName),
152  matrix_(matrix),
153  interfaceBouCoeffs_(interfaceBouCoeffs),
154  interfaceIntCoeffs_(interfaceIntCoeffs),
155  interfaces_(interfaces),
156  controlDict_(solverControls)
157 {
158  readControls();
159 }
160 
161 
162 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
163 
165 {
166  maxIter_ = controlDict_.lookupOrDefault<label>("maxIter", 1000);
167  minIter_ = controlDict_.lookupOrDefault<label>("minIter", 0);
168  tolerance_ = controlDict_.lookupOrDefault<scalar>("tolerance", 1e-6);
169  relTol_ = controlDict_.lookupOrDefault<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 (
182  const scalarField& psi,
183  const scalarField& source,
184  const scalarField& Apsi,
185  scalarField& tmpField
186 ) const
187 {
188  // --- Calculate A dot reference value of psi
189  matrix_.sumA(tmpField, interfaceBouCoeffs_, interfaces_);
190 
191  tmpField *= gAverage(psi, matrix_.lduMesh_.comm());
192 
193  return
194  gSum
195  (
196  (mag(Apsi - tmpField) + mag(source - tmpField))(),
197  matrix_.lduMesh_.comm()
198  )
200 
201  // At convergence this simpler method is equivalent to the above
202  // return 2*gSumMag(source) + solverPerformance::small_;
203 }
204 
205 
206 // ************************************************************************* //
virtual void read(const dictionary &)
Read and reset the solver parameters from the given stream.
virtual void readControls()
Read the control parameters from the controlDict_.
Generic field type.
Definition: FieldField.H:51
dimensioned< scalar > mag(const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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.
A class for handling words, derived from string.
Definition: word.H:59
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Namespace for OpenFOAM.
Type gSum(const FieldField< Field, Type > &f)
static const char nl
Definition: Ostream.H:260
const double e
Elementary charge.
Definition: doubleFloat.H:78
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
IOerror FatalIOError
bool asymmetric() const
Definition: lduMatrix.H:605
lduMatrix is a general matrix class in which the coefficients are stored as three arrays...
Definition: lduMatrix.H:77
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.
bool diagonal() const
Definition: lduMatrix.H:595
static const scalar small_
Small Type for the use in solvers.
Foam::diagonalSolver.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
Type gAverage(const FieldField< Field, Type > &f)
bool symmetric() const
Definition: lduMatrix.H:600
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
#define FatalIOErrorIn(functionName, ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:325
solver(const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces, const dictionary &solverControls)