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-2024 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 #include "meshState.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(lduMatrix, 1);
42 }
43 
44 
45 const Foam::scalar Foam::lduMatrix::defaultTolerance = 1e-6;
46 
47 const Foam::Enum
48 <
50 >
52 ({
53  { normTypes::NO_NORM, "none" },
54  { normTypes::DEFAULT_NORM, "default" },
55  { normTypes::L1_SCALED_NORM, "L1_scaled" },
56 });
57 
58 
59 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60 
62 :
63  lduMesh_(mesh)
64 {}
65 
66 
68 :
69  lduMesh_(A.lduMesh_)
70 {
71  if (A.diagPtr_)
72  {
73  diagPtr_ = std::make_unique<scalarField>(*(A.diagPtr_));
74  }
75 
76  if (A.upperPtr_)
77  {
78  upperPtr_ = std::make_unique<scalarField>(*(A.upperPtr_));
79  }
80 
81  if (A.lowerPtr_)
82  {
83  lowerPtr_ = std::make_unique<scalarField>(*(A.lowerPtr_));
84  }
85 }
86 
87 
88 Foam::lduMatrix::lduMatrix(lduMatrix&& A)
89 :
90  lduMesh_(A.lduMesh_),
91  diagPtr_(std::move(A.diagPtr_)),
92  lowerPtr_(std::move(A.lowerPtr_)),
93  upperPtr_(std::move(A.upperPtr_))
94 {}
95 
96 
98 :
99  lduMesh_(A.lduMesh_)
100 {
101  if (reuse)
102  {
103  diagPtr_ = std::move(A.diagPtr_);
104  upperPtr_ = std::move(A.upperPtr_);
105  lowerPtr_ = std::move(A.lowerPtr_);
106  }
107  else
108  {
109  // Copy assignment
110  if (A.diagPtr_)
111  {
112  diagPtr_ = std::make_unique<scalarField>(*(A.diagPtr_));
113  }
114 
115  if (A.upperPtr_)
116  {
117  upperPtr_ = std::make_unique<scalarField>(*(A.upperPtr_));
118  }
119 
120  if (A.lowerPtr_)
121  {
122  lowerPtr_ = std::make_unique<scalarField>(*(A.lowerPtr_));
123  }
124  }
125 }
126 
127 
128 Foam::lduMatrix::lduMatrix(const lduMesh& mesh, Istream& is)
129 :
130  lduMesh_(mesh)
131 {
132  bool withLower, withDiag, withUpper;
133 
134  is >> withLower >> withDiag >> withUpper;
135 
136  if (withLower)
137  {
138  lowerPtr_ = std::make_unique<scalarField>(is);
139  }
140  if (withDiag)
141  {
142  diagPtr_ = std::make_unique<scalarField>(is);
143  }
144  if (withUpper)
145  {
146  upperPtr_ = std::make_unique<scalarField>(is);
147  }
148 }
149 
150 
151 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
152 
154 {
155  if (diagPtr_)
156  {
157  return
158  (
159  (!upperPtr_)
160  ? (!lowerPtr_ ? "diagonal" : "diagonal-lower")
161  : (!lowerPtr_ ? "symmetric" : "asymmetric")
162  );
163  }
164 
165  // is empty (or just wrong)
166  return (!upperPtr_ && !lowerPtr_ ? "empty" : "ill-defined");
167 }
168 
169 
171 {
172  if (!diagPtr_)
173  {
175  << "diagPtr_ unallocated"
177  }
178 
179  return *diagPtr_;
180 }
181 
182 
184 {
185  if (!diagPtr_)
186  {
187  diagPtr_ =
188  std::make_unique<scalarField>(lduAddr().size(), Foam::zero{});
189  }
190 
191  return *diagPtr_;
192 }
193 
194 
196 {
197  if (!diagPtr_)
198  {
199  // if (size < 0)
200  // {
201  // size = lduAddr().size();
202  // }
203  diagPtr_ = std::make_unique<scalarField>(size, Foam::zero{});
204  }
205 
206  return *diagPtr_;
207 }
208 
209 
211 {
212  if (upperPtr_)
213  {
214  return *upperPtr_;
215  }
216  else
217  {
218  if (!lowerPtr_)
219  {
221  << "lowerPtr_ and upperPtr_ unallocated"
222  << abort(FatalError);
223  }
224 
225  return *lowerPtr_;
226  }
227 }
228 
229 
231 {
232  if (!upperPtr_)
233  {
234  if (lowerPtr_)
235  {
236  upperPtr_ = std::make_unique<scalarField>(*lowerPtr_);
237  }
238  else
239  {
240  upperPtr_ =
241  std::make_unique<scalarField>
242  (
243  lduAddr().lowerAddr().size(),
244  Foam::zero{}
245  );
246  }
247  }
248 
249  return *upperPtr_;
250 }
251 
252 
254 {
255  if (!upperPtr_)
256  {
257  if (lowerPtr_)
258  {
259  upperPtr_ = std::make_unique<scalarField>(*lowerPtr_);
260  }
261  else
262  {
263  // if (nCoeffs < 0)
264  // {
265  // nCoeffs = lduAddr().lowerAddr().size();
266  // }
267  upperPtr_ = std::make_unique<scalarField>(nCoeffs, Foam::zero{});
268  }
269  }
270 
271  return *upperPtr_;
272 }
273 
274 
276 {
277  if (lowerPtr_)
278  {
279  return *lowerPtr_;
280  }
281  else
282  {
283  if (!upperPtr_)
284  {
286  << "lowerPtr_ and upperPtr_ unallocated"
287  << abort(FatalError);
288  }
289 
290  return *upperPtr_;
291  }
292 }
293 
294 
296 {
297  if (!lowerPtr_)
298  {
299  if (upperPtr_)
300  {
301  lowerPtr_ = std::make_unique<scalarField>(*upperPtr_);
302  }
303  else
304  {
305  lowerPtr_ =
306  std::make_unique<scalarField>
307  (
308  lduAddr().lowerAddr().size(),
309  Foam::zero{}
310  );
311  }
312  }
313 
314  return *lowerPtr_;
315 }
316 
317 
319 {
320  if (!lowerPtr_)
321  {
322  if (upperPtr_)
323  {
324  lowerPtr_ = std::make_unique<scalarField>(*upperPtr_);
325  }
326  else
327  {
328  // if (nCoeffs < 0)
329  // {
330  // nCoeffs = lduAddr().lowerAddr().size();
331  // }
332  lowerPtr_ =
333  std::make_unique<scalarField>(nCoeffs, Foam::zero{});
334  }
335  }
336 
337  return *lowerPtr_;
338 }
339 
340 
342 (
343  const scalarField& residual,
344  const word& fieldName,
345  const bool initial
346 ) const
347 {
348  if (!mesh().hasDb())
349  {
350  return;
351  }
352 
353  scalarIOField* residualPtr =
355  (
356  initial
357  ? IOobject::scopedName("initialResidual", fieldName)
358  : IOobject::scopedName("residual", fieldName)
359  );
360 
361  if (residualPtr)
362  {
363  const auto* dataPtr = mesh().thisDb().findObject<meshState>("data");
364 
365  if (dataPtr)
366  {
367  if (initial && dataPtr->isFirstIteration())
368  {
369  *residualPtr = residual;
370  DebugInfo
371  << "Setting residual field for first solver iteration "
372  << "for solver field: " << fieldName << endl;
373  }
374  }
375  else
376  {
377  *residualPtr = residual;
378  DebugInfo
379  << "Setting residual field for solver field "
380  << fieldName << endl;
381  }
382  }
383 }
384 
385 
386 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
387 
388 Foam::Ostream& Foam::operator<<(Ostream& os, const lduMatrix& mat)
389 {
390  os << mat.hasLower() << token::SPACE
391  << mat.hasDiag() << token::SPACE
392  << mat.hasUpper() << token::SPACE;
393 
394  if (mat.hasLower())
395  {
396  os << mat.lower();
397  }
398 
399  if (mat.hasDiag())
400  {
401  os << mat.diag();
402  }
403 
404  if (mat.hasUpper())
405  {
406  os << mat.upper();
407  }
408 
410 
411  return os;
412 }
413 
414 
415 Foam::Ostream& Foam::operator<<
416 (
417  Ostream& os,
418  const InfoProxy<lduMatrix>& iproxy
419 )
420 {
421  const auto& mat = *iproxy;
422 
423  os << "Lower:" << Switch::name(mat.hasLower())
424  << " Diag:" << Switch::name(mat.hasDiag())
425  << " Upper:" << Switch::name(mat.hasUpper()) << endl;
426 
427  if (mat.hasLower())
428  {
429  os << "lower:" << mat.lower().size() << endl;
430  }
431  if (mat.hasDiag())
432  {
433  os << "diag :" << mat.diag().size() << endl;
434  }
435  if (mat.hasUpper())
436  {
437  os << "upper:" << mat.upper().size() << endl;
438  }
439 
440 
441  //if (hasLower)
442  //{
443  // os << "lower contents:" << endl;
444  // forAll(mat.lower(), i)
445  // {
446  // os << "i:" << i << "\t" << mat.lower()[i] << endl;
447  // }
448  // os << endl;
449  //}
450  //if (hasDiag)
451  //{
452  // os << "diag contents:" << endl;
453  // forAll(mat.diag(), i)
454  // {
455  // os << "i:" << i << "\t" << mat.diag()[i] << endl;
456  // }
457  // os << endl;
458  //}
459  //if (hasUpper)
460  //{
461  // os << "upper contents:" << endl;
462  // forAll(mat.upper(), i)
463  // {
464  // os << "i:" << i << "\t" << mat.upper()[i] << endl;
465  // }
466  // os << endl;
467  //}
468 
470 
471  return os;
472 }
473 
474 
475 // ************************************************************************* //
const scalarField & diag() const
Definition: lduMatrix.C:163
word matrixTypeName() const
The matrix type (empty, diagonal, symmetric, ...)
Definition: lduMatrix.C:146
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:608
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:45
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:335
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:53
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
IOField< scalar > scalarIOField
IO for a Field of scalar.
Definition: scalarIOField.H:32
static const Enum< normTypes > normTypesNames_
Names for the normTypes.
Definition: lduMatrix.H:124
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:376
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
static const char * name(const bool b) noexcept
A string representation of bool as "false" / "true".
Definition: Switch.C:141
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
A class for handling words, derived from Foam::string.
Definition: word.H:63
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:131
const scalarField & lower() const
Definition: lduMatrix.C:268
static const scalar defaultTolerance
Default (absolute) tolerance (1e-6)
Definition: lduMatrix.H:134
normTypes
Enumerated matrix normalisation types.
Definition: lduMatrix.H:114
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:56
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:77
lduMatrix is a general matrix class in which the coefficients are stored as three arrays...
Definition: lduMatrix.H:80
const scalarField & upper() const
Definition: lduMatrix.C:203
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
lduMatrix(const lduMesh &mesh)
Construct (without coefficients) for an LDU addressed mesh.
Definition: lduMatrix.C:54
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
Namespace for OpenFOAM.