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 #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 
61 Foam::lduMatrix::lduMatrix(const lduMesh& mesh)
62 :
63  lduMesh_(mesh),
64  lowerPtr_(nullptr),
65  diagPtr_(nullptr),
66  upperPtr_(nullptr)
67 {}
68 
69 
71 :
72  lduMesh_(A.lduMesh_),
73  lowerPtr_(nullptr),
74  diagPtr_(nullptr),
75  upperPtr_(nullptr)
76 {
77  if (A.lowerPtr_)
78  {
79  lowerPtr_ = new scalarField(*(A.lowerPtr_));
80  }
81 
82  if (A.diagPtr_)
83  {
84  diagPtr_ = new scalarField(*(A.diagPtr_));
85  }
86 
87  if (A.upperPtr_)
88  {
89  upperPtr_ = new scalarField(*(A.upperPtr_));
90  }
91 }
92 
93 
94 Foam::lduMatrix::lduMatrix(lduMatrix& A, bool reuse)
95 :
96  lduMesh_(A.lduMesh_),
97  lowerPtr_(nullptr),
98  diagPtr_(nullptr),
99  upperPtr_(nullptr)
100 {
101  if (reuse)
102  {
103  if (A.lowerPtr_)
104  {
105  lowerPtr_ = A.lowerPtr_;
106  A.lowerPtr_ = nullptr;
107  }
108 
109  if (A.diagPtr_)
110  {
111  diagPtr_ = A.diagPtr_;
112  A.diagPtr_ = nullptr;
113  }
114 
115  if (A.upperPtr_)
116  {
117  upperPtr_ = A.upperPtr_;
118  A.upperPtr_ = nullptr;
119  }
120  }
121  else
122  {
123  if (A.lowerPtr_)
124  {
125  lowerPtr_ = new scalarField(*(A.lowerPtr_));
126  }
127 
128  if (A.diagPtr_)
129  {
130  diagPtr_ = new scalarField(*(A.diagPtr_));
131  }
132 
133  if (A.upperPtr_)
134  {
135  upperPtr_ = new scalarField(*(A.upperPtr_));
136  }
137  }
138 }
139 
140 
141 Foam::lduMatrix::lduMatrix(const lduMesh& mesh, Istream& is)
142 :
143  lduMesh_(mesh),
144  lowerPtr_(nullptr),
145  diagPtr_(nullptr),
146  upperPtr_(nullptr)
147 {
148  Switch hasLow(is);
149  Switch hasDiag(is);
150  Switch hasUp(is);
151 
152  if (hasLow)
153  {
154  lowerPtr_ = new scalarField(is);
155  }
156  if (hasDiag)
157  {
158  diagPtr_ = new scalarField(is);
159  }
160  if (hasUp)
161  {
162  upperPtr_ = new scalarField(is);
163  }
164 }
165 
166 
168 {
169  if (lowerPtr_)
170  {
171  delete lowerPtr_;
172  }
173 
174  if (diagPtr_)
175  {
176  delete diagPtr_;
177  }
178 
179  if (upperPtr_)
180  {
181  delete upperPtr_;
182  }
183 }
184 
185 
187 {
188  if (!lowerPtr_)
189  {
190  if (upperPtr_)
191  {
192  lowerPtr_ = new scalarField(*upperPtr_);
193  }
194  else
195  {
196  lowerPtr_ = new scalarField(lduAddr().lowerAddr().size(), Zero);
197  }
198  }
199 
200  return *lowerPtr_;
201 }
202 
203 
205 {
206  if (!diagPtr_)
207  {
208  diagPtr_ = new scalarField(lduAddr().size(), Zero);
209  }
210 
211  return *diagPtr_;
212 }
213 
214 
216 {
217  if (!upperPtr_)
218  {
219  if (lowerPtr_)
220  {
221  upperPtr_ = new scalarField(*lowerPtr_);
222  }
223  else
224  {
225  upperPtr_ = new scalarField(lduAddr().lowerAddr().size(), Zero);
226  }
227  }
228 
229  return *upperPtr_;
230 }
231 
232 
233 Foam::scalarField& Foam::lduMatrix::lower(const label nCoeffs)
234 {
235  if (!lowerPtr_)
236  {
237  if (upperPtr_)
238  {
239  lowerPtr_ = new scalarField(*upperPtr_);
240  }
241  else
242  {
243  lowerPtr_ = new scalarField(nCoeffs, Zero);
244  }
245  }
246 
247  return *lowerPtr_;
248 }
249 
250 
251 Foam::scalarField& Foam::lduMatrix::diag(const label size)
252 {
253  if (!diagPtr_)
254  {
255  diagPtr_ = new scalarField(size, Zero);
256  }
257 
258  return *diagPtr_;
259 }
260 
261 
262 Foam::scalarField& Foam::lduMatrix::upper(const label nCoeffs)
263 {
264  if (!upperPtr_)
265  {
266  if (lowerPtr_)
267  {
268  upperPtr_ = new scalarField(*lowerPtr_);
269  }
270  else
271  {
272  upperPtr_ = new scalarField(nCoeffs, Zero);
273  }
274  }
275 
276  return *upperPtr_;
277 }
278 
279 
281 {
282  if (!lowerPtr_ && !upperPtr_)
283  {
285  << "lowerPtr_ or upperPtr_ unallocated"
286  << abort(FatalError);
287  }
288 
289  if (lowerPtr_)
290  {
291  return *lowerPtr_;
292  }
293  else
294  {
295  return *upperPtr_;
296  }
297 }
298 
299 
301 {
302  if (!diagPtr_)
303  {
305  << "diagPtr_ unallocated"
307  }
308 
309  return *diagPtr_;
310 }
311 
312 
314 {
315  if (!lowerPtr_ && !upperPtr_)
316  {
318  << "lowerPtr_ or upperPtr_ unallocated"
319  << abort(FatalError);
320  }
321 
322  if (upperPtr_)
323  {
324  return *upperPtr_;
325  }
326  else
327  {
328  return *lowerPtr_;
329  }
330 }
331 
332 
334 (
335  const scalarField& residual,
336  const word& fieldName,
337  const bool initial
338 ) const
339 {
340  if (!mesh().hasDb())
341  {
342  return;
343  }
344 
345  scalarIOField* residualPtr =
347  (
348  initial
349  ? IOobject::scopedName("initialResidual", fieldName)
350  : IOobject::scopedName("residual", fieldName)
351  );
352 
353  if (residualPtr)
354  {
355  const auto* dataPtr = mesh().thisDb().findObject<meshState>("data");
356 
357  if (dataPtr)
358  {
359  if (initial && dataPtr->isFirstIteration())
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  }
403 
405 
406  return os;
407 }
408 
409 
410 Foam::Ostream& Foam::operator<<
411 (
412  Ostream& os,
413  const InfoProxy<lduMatrix>& iproxy
414 )
415 {
416  const auto& ldum = *iproxy;
417 
418  Switch hasLow = ldum.hasLower();
419  Switch hasDiag = ldum.hasDiag();
420  Switch hasUp = ldum.hasUpper();
421 
422  os << "Lower:" << hasLow
423  << " Diag:" << hasDiag
424  << " Upper:" << hasUp << endl;
425 
426  if (hasLow)
427  {
428  os << "lower:" << ldum.lower().size() << endl;
429  }
430  if (hasDiag)
431  {
432  os << "diag :" << ldum.diag().size() << endl;
433  }
434  if (hasUp)
435  {
436  os << "upper:" << ldum.upper().size() << endl;
437  }
438 
439 
440  //if (hasLow)
441  //{
442  // os << "lower contents:" << endl;
443  // forAll(ldum.lower(), i)
444  // {
445  // os << "i:" << i << "\t" << ldum.lower()[i] << endl;
446  // }
447  // os << endl;
448  //}
449  //if (hasDiag)
450  //{
451  // os << "diag contents:" << endl;
452  // forAll(ldum.diag(), i)
453  // {
454  // os << "i:" << i << "\t" << ldum.diag()[i] << endl;
455  // }
456  // os << endl;
457  //}
458  //if (hasUp)
459  //{
460  // os << "upper contents:" << endl;
461  // forAll(ldum.upper(), i)
462  // {
463  // os << "i:" << i << "\t" << ldum.upper()[i] << endl;
464  // }
465  // os << endl;
466  //}
467 
469 
470  return os;
471 }
472 
473 
474 // ************************************************************************* //
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:598
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:54
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:327
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
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
scalarField & upper()
Definition: lduMatrix.C:208
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: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
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:131
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: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:79
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
scalarField & lower()
Definition: lduMatrix.C:179
bool hasDiag() const noexcept
Definition: lduMatrix.H:766
scalarField & diag()
Definition: lduMatrix.C:197
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
Namespace for OpenFOAM.
~lduMatrix()
Destructor.
Definition: lduMatrix.C:160
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127