lduMatrix.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 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 "IOstreams.H"
31 #include "Switch.H"
32 #include "objectRegistry.H"
33 #include "scalarIOField.H"
34 #include "Time.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(lduMatrix, 1);
41 }
42 
43 
44 const Foam::label Foam::lduMatrix::solver::defaultMaxIter_ = 1000;
45 
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
50 :
51  lduMesh_(mesh),
52  lowerPtr_(nullptr),
53  diagPtr_(nullptr),
54  upperPtr_(nullptr)
55 {}
56 
57 
59 :
60  lduMesh_(A.lduMesh_),
61  lowerPtr_(nullptr),
62  diagPtr_(nullptr),
63  upperPtr_(nullptr)
64 {
65  if (A.lowerPtr_)
66  {
67  lowerPtr_ = new scalarField(*(A.lowerPtr_));
68  }
69 
70  if (A.diagPtr_)
71  {
72  diagPtr_ = new scalarField(*(A.diagPtr_));
73  }
74 
75  if (A.upperPtr_)
76  {
77  upperPtr_ = new scalarField(*(A.upperPtr_));
78  }
79 }
80 
81 
83 :
84  lduMesh_(A.lduMesh_),
85  lowerPtr_(nullptr),
86  diagPtr_(nullptr),
87  upperPtr_(nullptr)
88 {
89  if (reuse)
90  {
91  if (A.lowerPtr_)
92  {
93  lowerPtr_ = A.lowerPtr_;
94  A.lowerPtr_ = nullptr;
95  }
96 
97  if (A.diagPtr_)
98  {
99  diagPtr_ = A.diagPtr_;
100  A.diagPtr_ = nullptr;
101  }
102 
103  if (A.upperPtr_)
104  {
105  upperPtr_ = A.upperPtr_;
106  A.upperPtr_ = nullptr;
107  }
108  }
109  else
110  {
111  if (A.lowerPtr_)
112  {
113  lowerPtr_ = new scalarField(*(A.lowerPtr_));
114  }
115 
116  if (A.diagPtr_)
117  {
118  diagPtr_ = new scalarField(*(A.diagPtr_));
119  }
120 
121  if (A.upperPtr_)
122  {
123  upperPtr_ = new scalarField(*(A.upperPtr_));
124  }
125  }
126 }
127 
128 
130 :
131  lduMesh_(mesh),
132  lowerPtr_(nullptr),
133  diagPtr_(nullptr),
134  upperPtr_(nullptr)
135 {
136  Switch hasLow(is);
137  Switch hasDiag(is);
138  Switch hasUp(is);
139 
140  if (hasLow)
141  {
142  lowerPtr_ = new scalarField(is);
143  }
144  if (hasDiag)
145  {
146  diagPtr_ = new scalarField(is);
147  }
148  if (hasUp)
149  {
150  upperPtr_ = new scalarField(is);
151  }
152 }
153 
154 
156 {
157  if (lowerPtr_)
158  {
159  delete lowerPtr_;
160  }
161 
162  if (diagPtr_)
163  {
164  delete diagPtr_;
165  }
166 
167  if (upperPtr_)
168  {
169  delete upperPtr_;
170  }
171 }
172 
173 
175 {
176  if (!lowerPtr_)
177  {
178  if (upperPtr_)
179  {
180  lowerPtr_ = new scalarField(*upperPtr_);
181  }
182  else
183  {
184  lowerPtr_ = new scalarField(lduAddr().lowerAddr().size(), Zero);
185  }
186  }
187 
188  return *lowerPtr_;
189 }
190 
191 
193 {
194  if (!diagPtr_)
195  {
196  diagPtr_ = new scalarField(lduAddr().size(), Zero);
197  }
198 
199  return *diagPtr_;
200 }
201 
202 
204 {
205  if (!upperPtr_)
206  {
207  if (lowerPtr_)
208  {
209  upperPtr_ = new scalarField(*lowerPtr_);
210  }
211  else
212  {
213  upperPtr_ = new scalarField(lduAddr().lowerAddr().size(), Zero);
214  }
215  }
216 
217  return *upperPtr_;
218 }
219 
220 
222 {
223  if (!lowerPtr_)
224  {
225  if (upperPtr_)
226  {
227  lowerPtr_ = new scalarField(*upperPtr_);
228  }
229  else
230  {
231  lowerPtr_ = new scalarField(nCoeffs, Zero);
232  }
233  }
234 
235  return *lowerPtr_;
236 }
237 
238 
240 {
241  if (!diagPtr_)
242  {
243  diagPtr_ = new scalarField(size, Zero);
244  }
245 
246  return *diagPtr_;
247 }
248 
249 
251 {
252  if (!upperPtr_)
253  {
254  if (lowerPtr_)
255  {
256  upperPtr_ = new scalarField(*lowerPtr_);
257  }
258  else
259  {
260  upperPtr_ = new scalarField(nCoeffs, Zero);
261  }
262  }
263 
264  return *upperPtr_;
265 }
266 
267 
269 {
270  if (!lowerPtr_ && !upperPtr_)
271  {
273  << "lowerPtr_ or upperPtr_ unallocated"
274  << abort(FatalError);
275  }
276 
277  if (lowerPtr_)
278  {
279  return *lowerPtr_;
280  }
281  else
282  {
283  return *upperPtr_;
284  }
285 }
286 
287 
289 {
290  if (!diagPtr_)
291  {
293  << "diagPtr_ unallocated"
294  << abort(FatalError);
295  }
296 
297  return *diagPtr_;
298 }
299 
300 
302 {
303  if (!lowerPtr_ && !upperPtr_)
304  {
306  << "lowerPtr_ or upperPtr_ unallocated"
307  << abort(FatalError);
308  }
309 
310  if (upperPtr_)
311  {
312  return *upperPtr_;
313  }
314  else
315  {
316  return *lowerPtr_;
317  }
318 }
319 
320 
322 (
323  const scalarField& residual,
324  const word& fieldName,
325  const bool initial
326 ) const
327 {
328  if (!lduMesh_.hasDb())
329  {
330  return;
331  }
332 
333  word lookupName;
334  if (initial)
335  {
336  lookupName = word("initialResidual:" + fieldName);
337  }
338  else
339  {
340  lookupName = word("residual:" + fieldName);
341  }
342 
343  scalarIOField* residualPtr =
344  lduMesh_.thisDb().getObjectPtr<scalarIOField>(lookupName);
345 
346  if (residualPtr)
347  {
348  const IOdictionary* dataPtr =
349  lduMesh_.thisDb().findObject<IOdictionary>("data");
350 
351  if (dataPtr)
352  {
353  if (initial && dataPtr->found("firstIteration"))
354  {
355  *residualPtr = residual;
356  DebugInfo
357  << "Setting residual field for first solver iteration "
358  << "for solver field: " << fieldName << endl;
359  }
360  }
361  else
362  {
363  *residualPtr = residual;
364  DebugInfo
365  << "Setting residual field for solver field "
366  << fieldName << endl;
367  }
368  }
369 }
370 
371 
372 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
373 
375 {
376  Switch hasLow = ldum.hasLower();
377  Switch hasDiag = ldum.hasDiag();
378  Switch hasUp = ldum.hasUpper();
379 
380  os << hasLow << token::SPACE << hasDiag << token::SPACE
381  << hasUp << token::SPACE;
382 
383  if (hasLow)
384  {
385  os << ldum.lower();
386  }
387 
388  if (hasDiag)
389  {
390  os << ldum.diag();
391  }
392 
393  if (hasUp)
394  {
395  os << ldum.upper();
396  }
397 
398  os.check(FUNCTION_NAME);
399 
400  return os;
401 }
402 
403 
405 {
406  const lduMatrix& ldum = ip.t_;
407 
408  Switch hasLow = ldum.hasLower();
409  Switch hasDiag = ldum.hasDiag();
410  Switch hasUp = ldum.hasUpper();
411 
412  os << "Lower:" << hasLow
413  << " Diag:" << hasDiag
414  << " Upper:" << hasUp << endl;
415 
416  if (hasLow)
417  {
418  os << "lower:" << ldum.lower().size() << endl;
419  }
420  if (hasDiag)
421  {
422  os << "diag :" << ldum.diag().size() << endl;
423  }
424  if (hasUp)
425  {
426  os << "upper:" << ldum.upper().size() << endl;
427  }
428 
429 
430  //if (hasLow)
431  //{
432  // os << "lower contents:" << endl;
433  // forAll(ldum.lower(), i)
434  // {
435  // os << "i:" << i << "\t" << ldum.lower()[i] << endl;
436  // }
437  // os << endl;
438  //}
439  //if (hasDiag)
440  //{
441  // os << "diag contents:" << endl;
442  // forAll(ldum.diag(), i)
443  // {
444  // os << "i:" << i << "\t" << ldum.diag()[i] << endl;
445  // }
446  // os << endl;
447  //}
448  //if (hasUp)
449  //{
450  // os << "upper contents:" << endl;
451  // forAll(ldum.upper(), i)
452  // {
453  // os << "i:" << i << "\t" << ldum.upper()[i] << endl;
454  // }
455  // os << endl;
456  //}
457 
458  os.check(FUNCTION_NAME);
459 
460  return os;
461 }
462 
463 
464 // ************************************************************************* //
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:77
IOstreams.H
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
Foam::lduMatrix::hasLower
bool hasLower() const
Definition: lduMatrix.H:612
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::InfoProxy
A helper class for outputting values to Ostream.
Definition: InfoProxy.H:47
Foam::IOField
A primitive field of type <T> with automated input and output.
Definition: foamVtkLagrangianWriter.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::lduMatrix::~lduMatrix
~lduMatrix()
Destructor.
Definition: lduMatrix.C:155
Foam::dictionary::found
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionary.C:364
Foam::lduMatrix::upper
scalarField & upper()
Definition: lduMatrix.C:203
Foam::lduMatrix
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition: lduMatrix.H:83
Foam::lduMatrix::hasDiag
bool hasDiag() const
Definition: lduMatrix.H:602
objectRegistry.H
A
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
lduMatrix.H
scalarIOField.H
Foam::lduMatrix::lduMatrix
lduMatrix(const lduMesh &)
Construct given an LDU addressed mesh.
Definition: lduMatrix.C:49
Foam::operator<<
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
Foam::lduMatrix::setResidualField
void setResidualField(const scalarField &residual, const word &fieldName, const bool initial) const
Definition: lduMatrix.C:322
Foam::Field< scalar >
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::InfoProxy::t_
const T & t_
Definition: InfoProxy.H:55
Switch.H
Foam::IOstream::check
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:51
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::lduMatrix::diag
scalarField & diag()
Definition: lduMatrix.C:192
Time.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:359
Foam::token::SPACE
Space [isspace].
Definition: token.H:117
FUNCTION_NAME
#define FUNCTION_NAME
Definition: messageStream.H:270
Foam::lduMatrix::hasUpper
bool hasUpper() const
Definition: lduMatrix.H:607
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::lduMatrix::solver::defaultMaxIter_
static const label defaultMaxIter_
Default maximum number of iterations in the solver.
Definition: lduMatrix.H:113
Foam::lduMatrix::lower
scalarField & lower()
Definition: lduMatrix.C:174
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::lduMesh
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:62