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-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 "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::scalar Foam::lduMatrix::defaultTolerance = 1e-6;
45 
46 const Foam::Enum
47 <
49 >
51 ({
52  { normTypes::NO_NORM, "none" },
53  { normTypes::DEFAULT_NORM, "default" },
54  { normTypes::L1_SCALED_NORM, "L1_scaled" },
55 });
56 
57 
58 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59 
60 Foam::lduMatrix::lduMatrix(const lduMesh& mesh)
61 :
62  lduMesh_(mesh),
63  lowerPtr_(nullptr),
64  diagPtr_(nullptr),
65  upperPtr_(nullptr)
66 {}
67 
68 
70 :
71  lduMesh_(A.lduMesh_),
72  lowerPtr_(nullptr),
73  diagPtr_(nullptr),
74  upperPtr_(nullptr)
75 {
76  if (A.lowerPtr_)
77  {
78  lowerPtr_ = new scalarField(*(A.lowerPtr_));
79  }
80 
81  if (A.diagPtr_)
82  {
83  diagPtr_ = new scalarField(*(A.diagPtr_));
84  }
85 
86  if (A.upperPtr_)
87  {
88  upperPtr_ = new scalarField(*(A.upperPtr_));
89  }
90 }
91 
92 
93 Foam::lduMatrix::lduMatrix(lduMatrix& A, bool reuse)
94 :
95  lduMesh_(A.lduMesh_),
96  lowerPtr_(nullptr),
97  diagPtr_(nullptr),
98  upperPtr_(nullptr)
99 {
100  if (reuse)
101  {
102  if (A.lowerPtr_)
103  {
104  lowerPtr_ = A.lowerPtr_;
105  A.lowerPtr_ = nullptr;
106  }
107 
108  if (A.diagPtr_)
109  {
110  diagPtr_ = A.diagPtr_;
111  A.diagPtr_ = nullptr;
112  }
113 
114  if (A.upperPtr_)
115  {
116  upperPtr_ = A.upperPtr_;
117  A.upperPtr_ = nullptr;
118  }
119  }
120  else
121  {
122  if (A.lowerPtr_)
123  {
124  lowerPtr_ = new scalarField(*(A.lowerPtr_));
125  }
126 
127  if (A.diagPtr_)
128  {
129  diagPtr_ = new scalarField(*(A.diagPtr_));
130  }
131 
132  if (A.upperPtr_)
133  {
134  upperPtr_ = new scalarField(*(A.upperPtr_));
135  }
136  }
137 }
138 
139 
140 Foam::lduMatrix::lduMatrix(const lduMesh& mesh, Istream& is)
141 :
142  lduMesh_(mesh),
143  lowerPtr_(nullptr),
144  diagPtr_(nullptr),
145  upperPtr_(nullptr)
146 {
147  Switch hasLow(is);
148  Switch hasDiag(is);
149  Switch hasUp(is);
150 
151  if (hasLow)
152  {
153  lowerPtr_ = new scalarField(is);
154  }
155  if (hasDiag)
156  {
157  diagPtr_ = new scalarField(is);
158  }
159  if (hasUp)
160  {
161  upperPtr_ = new scalarField(is);
162  }
163 }
164 
165 
167 {
168  if (lowerPtr_)
169  {
170  delete lowerPtr_;
171  }
172 
173  if (diagPtr_)
174  {
175  delete diagPtr_;
176  }
177 
178  if (upperPtr_)
179  {
180  delete upperPtr_;
181  }
182 }
183 
184 
186 {
187  if (!lowerPtr_)
188  {
189  if (upperPtr_)
190  {
191  lowerPtr_ = new scalarField(*upperPtr_);
192  }
193  else
194  {
195  lowerPtr_ = new scalarField(lduAddr().lowerAddr().size(), Zero);
196  }
197  }
198 
199  return *lowerPtr_;
200 }
201 
202 
204 {
205  if (!diagPtr_)
206  {
207  diagPtr_ = new scalarField(lduAddr().size(), Zero);
208  }
209 
210  return *diagPtr_;
211 }
212 
213 
215 {
216  if (!upperPtr_)
217  {
218  if (lowerPtr_)
219  {
220  upperPtr_ = new scalarField(*lowerPtr_);
221  }
222  else
223  {
224  upperPtr_ = new scalarField(lduAddr().lowerAddr().size(), Zero);
225  }
226  }
227 
228  return *upperPtr_;
229 }
230 
231 
232 Foam::scalarField& Foam::lduMatrix::lower(const label nCoeffs)
233 {
234  if (!lowerPtr_)
235  {
236  if (upperPtr_)
237  {
238  lowerPtr_ = new scalarField(*upperPtr_);
239  }
240  else
241  {
242  lowerPtr_ = new scalarField(nCoeffs, Zero);
243  }
244  }
245 
246  return *lowerPtr_;
247 }
248 
249 
250 Foam::scalarField& Foam::lduMatrix::diag(const label size)
251 {
252  if (!diagPtr_)
253  {
254  diagPtr_ = new scalarField(size, Zero);
255  }
256 
257  return *diagPtr_;
258 }
259 
260 
261 Foam::scalarField& Foam::lduMatrix::upper(const label nCoeffs)
262 {
263  if (!upperPtr_)
264  {
265  if (lowerPtr_)
266  {
267  upperPtr_ = new scalarField(*lowerPtr_);
268  }
269  else
270  {
271  upperPtr_ = new scalarField(nCoeffs, Zero);
272  }
273  }
274 
275  return *upperPtr_;
276 }
277 
278 
280 {
281  if (!lowerPtr_ && !upperPtr_)
282  {
284  << "lowerPtr_ or upperPtr_ unallocated"
285  << abort(FatalError);
286  }
287 
288  if (lowerPtr_)
289  {
290  return *lowerPtr_;
291  }
292  else
293  {
294  return *upperPtr_;
295  }
296 }
297 
298 
300 {
301  if (!diagPtr_)
302  {
304  << "diagPtr_ unallocated"
306  }
307 
308  return *diagPtr_;
309 }
310 
311 
313 {
314  if (!lowerPtr_ && !upperPtr_)
315  {
317  << "lowerPtr_ or upperPtr_ unallocated"
318  << abort(FatalError);
319  }
320 
321  if (upperPtr_)
322  {
323  return *upperPtr_;
324  }
325  else
326  {
327  return *lowerPtr_;
328  }
329 }
330 
331 
333 (
334  const scalarField& residual,
335  const word& fieldName,
336  const bool initial
337 ) const
338 {
339  if (!mesh().hasDb())
340  {
341  return;
342  }
343 
344  scalarIOField* residualPtr =
346  (
347  initial
348  ? IOobject::scopedName("initialResidual", fieldName)
349  : IOobject::scopedName("residual", fieldName)
350  );
351 
352  if (residualPtr)
353  {
354  const IOdictionary* dataPtr =
355  mesh().thisDb().findObject<IOdictionary>("data");
356 
357  if (dataPtr)
358  {
359  if (initial && dataPtr->found("firstIteration"))
360  {
361  *residualPtr = residual;
362  DebugInfo
363  << "Setting residual field for first solver iteration "
364  << "for solver field: " << fieldName << endl;
365  }
366  }
367  else
368  {
369  *residualPtr = residual;
370  DebugInfo
371  << "Setting residual field for solver field "
372  << fieldName << endl;
373  }
374  }
375 }
376 
377 
378 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
379 
380 Foam::Ostream& Foam::operator<<(Ostream& os, const lduMatrix& ldum)
381 {
382  Switch hasLow = ldum.hasLower();
383  Switch hasDiag = ldum.hasDiag();
384  Switch hasUp = ldum.hasUpper();
385 
386  os << hasLow << token::SPACE << hasDiag << token::SPACE
387  << hasUp << token::SPACE;
388 
389  if (hasLow)
390  {
391  os << ldum.lower();
392  }
393 
394  if (hasDiag)
395  {
396  os << ldum.diag();
397  }
398 
399  if (hasUp)
400  {
401  os << ldum.upper();
402  }
405 
406  return os;
407 }
408 
409 
410 Foam::Ostream& Foam::operator<<(Ostream& os, const InfoProxy<lduMatrix>& ip)
411 {
412  const lduMatrix& ldum = ip.t_;
413 
414  Switch hasLow = ldum.hasLower();
415  Switch hasDiag = ldum.hasDiag();
416  Switch hasUp = ldum.hasUpper();
417 
418  os << "Lower:" << hasLow
419  << " Diag:" << hasDiag
420  << " Upper:" << hasUp << endl;
421 
422  if (hasLow)
423  {
424  os << "lower:" << ldum.lower().size() << endl;
425  }
426  if (hasDiag)
427  {
428  os << "diag :" << ldum.diag().size() << endl;
429  }
430  if (hasUp)
431  {
432  os << "upper:" << ldum.upper().size() << endl;
433  }
434 
435 
436  //if (hasLow)
437  //{
438  // os << "lower contents:" << endl;
439  // forAll(ldum.lower(), i)
440  // {
441  // os << "i:" << i << "\t" << ldum.lower()[i] << endl;
442  // }
443  // os << endl;
444  //}
445  //if (hasDiag)
446  //{
447  // os << "diag contents:" << endl;
448  // forAll(ldum.diag(), i)
449  // {
450  // os << "i:" << i << "\t" << ldum.diag()[i] << endl;
451  // }
452  // os << endl;
453  //}
454  //if (hasUp)
455  //{
456  // os << "upper contents:" << endl;
457  // forAll(ldum.upper(), i)
458  // {
459  // os << "i:" << i << "\t" << ldum.upper()[i] << endl;
460  // }
461  // os << endl;
462  //}
463 
465 
466  return os;
467 }
468 
469 
470 // ************************************************************************* //
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:45
lduMatrix(const lduMesh &)
Construct given an LDU addressed mesh.
Definition: lduMatrix.C:53
void setResidualField(const scalarField &residual, const word &fieldName, const bool initial) const
Set the residual field using an IOField on the object registry if it exists.
Definition: lduMatrix.C:326
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
scalarField & upper()
Definition: lduMatrix.C:207
static const Enum< normTypes > normTypesNames_
Names for the normTypes.
Definition: lduMatrix.H:113
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:377
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
Type * getObjectPtr(const word &name, const bool recursive=false) const
Return non-const pointer to the object of the given Type, using a const-cast to have it behave like a...
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Space [isspace].
Definition: token.H:128
IOField< scalar > scalarIOField
scalarField with IO.
Definition: scalarIOField.H:38
static const scalar defaultTolerance
Default (absolute) tolerance (1e-6)
Definition: lduMatrix.H:123
normTypes
Enumerated matrix normalisation types.
Definition: lduMatrix.H:103
errorManip< error > abort(error &err)
Definition: errorManip.H:139
#define DebugInfo
Report an information message using Foam::Info.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:55
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
defineTypeNameAndDebug(combustionModel, 0)
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:76
lduMatrix is a general matrix class in which the coefficients are stored as three arrays...
Definition: lduMatrix.H:79
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
scalarField & lower()
Definition: lduMatrix.C:178
bool hasDiag() const noexcept
Definition: lduMatrix.H:746
scalarField & diag()
Definition: lduMatrix.C:196
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
Namespace for OpenFOAM.
~lduMatrix()
Destructor.
Definition: lduMatrix.C:159
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157