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-------------------------------------------------------------------------------
11License
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
38namespace Foam
39{
41}
42
43
44const 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 (!mesh().hasDb())
329 {
330 return;
331 }
332
333 scalarIOField* residualPtr =
335 (
336 initial
337 ? IOobject::scopedName("initialResidual", fieldName)
338 : IOobject::scopedName("residual", fieldName)
339 );
340
341 if (residualPtr)
342 {
343 const IOdictionary* dataPtr =
344 mesh().thisDb().findObject<IOdictionary>("data");
345
346 if (dataPtr)
347 {
348 if (initial && dataPtr->found("firstIteration"))
349 {
350 *residualPtr = residual;
352 << "Setting residual field for first solver iteration "
353 << "for solver field: " << fieldName << endl;
354 }
355 }
356 else
357 {
358 *residualPtr = residual;
360 << "Setting residual field for solver field "
361 << fieldName << endl;
362 }
363 }
364}
365
366
367// * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
368
370{
371 Switch hasLow = ldum.hasLower();
372 Switch hasDiag = ldum.hasDiag();
373 Switch hasUp = ldum.hasUpper();
374
375 os << hasLow << token::SPACE << hasDiag << token::SPACE
376 << hasUp << token::SPACE;
377
378 if (hasLow)
379 {
380 os << ldum.lower();
381 }
382
383 if (hasDiag)
384 {
385 os << ldum.diag();
386 }
387
388 if (hasUp)
389 {
390 os << ldum.upper();
391 }
392
394
395 return os;
396}
397
398
400{
401 const lduMatrix& ldum = ip.t_;
402
403 Switch hasLow = ldum.hasLower();
404 Switch hasDiag = ldum.hasDiag();
405 Switch hasUp = ldum.hasUpper();
406
407 os << "Lower:" << hasLow
408 << " Diag:" << hasDiag
409 << " Upper:" << hasUp << endl;
410
411 if (hasLow)
412 {
413 os << "lower:" << ldum.lower().size() << endl;
414 }
415 if (hasDiag)
416 {
417 os << "diag :" << ldum.diag().size() << endl;
418 }
419 if (hasUp)
420 {
421 os << "upper:" << ldum.upper().size() << endl;
422 }
423
424
425 //if (hasLow)
426 //{
427 // os << "lower contents:" << endl;
428 // forAll(ldum.lower(), i)
429 // {
430 // os << "i:" << i << "\t" << ldum.lower()[i] << endl;
431 // }
432 // os << endl;
433 //}
434 //if (hasDiag)
435 //{
436 // os << "diag contents:" << endl;
437 // forAll(ldum.diag(), i)
438 // {
439 // os << "i:" << i << "\t" << ldum.diag()[i] << endl;
440 // }
441 // os << endl;
442 //}
443 //if (hasUp)
444 //{
445 // os << "upper contents:" << endl;
446 // forAll(ldum.upper(), i)
447 // {
448 // os << "i:" << i << "\t" << ldum.upper()[i] << endl;
449 // }
450 // os << endl;
451 //}
452
454
455 return os;
456}
457
458
459// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:47
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:58
A helper class for outputting values to Ostream.
Definition: InfoProxy.H:52
const T & t_
Definition: InfoProxy.H:55
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:78
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:302
static const label defaultMaxIter_
Default maximum number of iterations in the solver.
Definition: lduMatrix.H:105
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition: lduMatrix.H:84
scalarField & upper()
Definition: lduMatrix.C:203
bool hasUpper() const
Definition: lduMatrix.H:613
bool hasLower() const
Definition: lduMatrix.H:618
void setResidualField(const scalarField &residual, const word &fieldName, const bool initial) const
Definition: lduMatrix.C:322
scalarField & lower()
Definition: lduMatrix.C:174
~lduMatrix()
Destructor.
Definition: lduMatrix.C:155
scalarField & diag()
Definition: lduMatrix.C:192
bool hasDiag() const
Definition: lduMatrix.H:608
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:63
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
Type * getObjectPtr(const word &name, const bool recursive=false) const
@ SPACE
Space [isspace].
Definition: token.H:125
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
#define DebugInfo
Report an information message using Foam::Info.
#define FUNCTION_NAME
Namespace for OpenFOAM.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
errorManip< error > abort(error &err)
Definition: errorManip.H:144
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError