Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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#include "meshState.H"
36
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
40{
42}
43
44
45const Foam::scalar Foam::lduMatrix::defaultTolerance = 1e-6;
46
47const Foam::Enum
48<
50>
52({
53 { normTypes::NO_NORM, "none" },
55 { normTypes::L1_SCALED_NORM, "L1_scaled" },
56});
57
58
59// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
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 if (A.lowerCSRPtr_)
87 {
88 lowerCSRPtr_ = std::make_unique<scalarField>(*(A.lowerCSRPtr_));
89 }
90
91 // No need to copy work
92}
93
94
96:
97 lduMesh_(A.lduMesh_),
98 diagPtr_(std::move(A.diagPtr_)),
99 lowerPtr_(std::move(A.lowerPtr_)),
100 upperPtr_(std::move(A.upperPtr_)),
101 lowerCSRPtr_(std::move(A.lowerCSRPtr_)),
102 workPtr_(std::move(A.workPtr_))
103{}
104
105
107:
108 lduMesh_(A.lduMesh_)
109{
110 if (reuse)
111 {
112 diagPtr_ = std::move(A.diagPtr_);
113 upperPtr_ = std::move(A.upperPtr_);
114 lowerPtr_ = std::move(A.lowerPtr_);
115 lowerCSRPtr_ = std::move(A.lowerCSRPtr_);
116 workPtr_ = std::move(A.workPtr_);
117 }
118 else
119 {
120 // Copy assignment
121 if (A.diagPtr_)
122 {
123 diagPtr_ = std::make_unique<scalarField>(*(A.diagPtr_));
124 }
125
126 if (A.upperPtr_)
127 {
128 upperPtr_ = std::make_unique<scalarField>(*(A.upperPtr_));
129 }
130
131 if (A.lowerPtr_)
132 {
133 lowerPtr_ = std::make_unique<scalarField>(*(A.lowerPtr_));
134 }
135
136 // Note: no real need to keep lowerCSR except we use it (hasLowerCSR())
137 // to trigger certain actions
138 if (A.lowerCSRPtr_)
140 lowerCSRPtr_ = std::make_unique<scalarField>(*(A.lowerCSRPtr_));
141 }
142 }
143}
144
145
147:
148 lduMesh_(mesh)
149{
150 bool withLower, withDiag, withUpper, withLowerCSR;
151
152 is >> withLower >> withDiag >> withUpper >> withLowerCSR;
153
154 if (withLower)
155 {
156 lowerPtr_ = std::make_unique<scalarField>(is);
157 }
158 if (withDiag)
159 {
160 diagPtr_ = std::make_unique<scalarField>(is);
161 }
162 if (withUpper)
163 {
164 upperPtr_ = std::make_unique<scalarField>(is);
165 }
166 if (withLowerCSR)
167 {
168 // Check if it has been sent over or we need to construct it. We need
169 // to set the lowerCSRPtr since it determines behaviour.
170 if (withLower)
171 {
172 // Construct from lower
173 (void)lowerCSR();
174 }
175 else
176 {
177 lowerCSRPtr_ = std::make_unique<scalarField>(is);
179 }
180}
181
182
183// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
184
185Foam::word Foam::lduMatrix::matrixTypeName() const
186{
187 if (diagPtr_)
188 {
189 return
190 (
191 (!upperPtr_)
192 ? (!lowerPtr_ ? "diagonal" : "diagonal-lower")
193 : (!lowerPtr_ ? "symmetric" : "asymmetric")
194 );
196
197 // is empty (or just wrong)
198 return (!upperPtr_ && !lowerPtr_ ? "empty" : "ill-defined");
199}
200
201
203{
204 if (!diagPtr_)
205 {
207 << "diagPtr_ unallocated"
209 }
210
211 return *diagPtr_;
212}
213
214
216{
217 if (!diagPtr_)
218 {
219 diagPtr_ =
220 std::make_unique<scalarField>(lduAddr().size(), Foam::zero{});
221 }
222
223 return *diagPtr_;
224}
225
226
228{
229 if (!diagPtr_)
230 {
231 // if (size < 0)
232 // {
233 // size = lduAddr().size();
234 // }
235 diagPtr_ = std::make_unique<scalarField>(size, Foam::zero{});
236 }
237
238 return *diagPtr_;
239}
240
241
243{
244 if (upperPtr_)
245 {
246 return *upperPtr_;
247 }
248 else
249 {
250 if (!lowerPtr_)
251 {
253 << "lowerPtr_ and upperPtr_ unallocated"
254 << abort(FatalError);
256
257 return *lowerPtr_;
258 }
259}
260
261
263{
264 if (!upperPtr_)
265 {
266 if (lowerPtr_)
267 {
268 upperPtr_ = std::make_unique<scalarField>(*lowerPtr_);
269 }
270 else
271 {
272 // no lowerPtr so any lowerCSR was constructed from upper
273 lowerCSRPtr_.reset(nullptr);
274
275 upperPtr_ =
276 std::make_unique<scalarField>
277 (
278 lduAddr().lowerAddr().size(),
279 Foam::zero{}
280 );
282 }
283
284 return *upperPtr_;
285}
286
287
289{
290 if (!upperPtr_)
291 {
292 if (lowerPtr_)
293 {
294 upperPtr_ = std::make_unique<scalarField>(*lowerPtr_);
295 }
296 else
297 {
298 // no lowerPtr so any lowerCSR was constructed from upper
299 lowerCSRPtr_.reset(nullptr);
300
301 // if (nCoeffs < 0)
302 // {
303 // nCoeffs = lduAddr().lowerAddr().size();
304 // }
305 upperPtr_ = std::make_unique<scalarField>(nCoeffs, Foam::zero{});
307 }
308
309 return *upperPtr_;
310}
311
312
314{
315 if (lowerPtr_)
316 {
317 return *lowerPtr_;
318 }
319 else
320 {
321 if (!upperPtr_)
322 {
324 << "lowerPtr_ and upperPtr_ unallocated"
325 << abort(FatalError);
327
328 return *upperPtr_;
329 }
330}
331
332
334{
335 if (!lowerPtr_)
336 {
337 lowerCSRPtr_.reset(nullptr);
338
339 if (upperPtr_)
340 {
341 lowerPtr_ = std::make_unique<scalarField>(*upperPtr_);
342 }
343 else
344 {
345 lowerPtr_ =
346 std::make_unique<scalarField>
347 (
348 lduAddr().lowerAddr().size(),
349 Foam::zero{}
350 );
352 }
353
354 return *lowerPtr_;
355}
356
357
359{
360 if (!lowerPtr_)
361 {
362 lowerCSRPtr_.reset(nullptr);
363
364 if (upperPtr_)
365 {
366 lowerPtr_ = std::make_unique<scalarField>(*upperPtr_);
367 }
368 else
369 {
370 // if (nCoeffs < 0)
371 // {
372 // nCoeffs = lduAddr().lowerAddr().size();
373 // }
374 lowerPtr_ =
375 std::make_unique<scalarField>(nCoeffs, Foam::zero{});
377 }
378
379 return *lowerPtr_;
380}
381
382
384{
385 if (!lowerCSRPtr_)
386 {
387 const label nLower = lduAddr().losortAddr().size();
388
389 lowerCSRPtr_ = std::make_unique<scalarField>(nLower);
390
391 if (lowerPtr_)
392 {
393 lduAddr().map(*lowerPtr_, *lowerCSRPtr_);
394 }
395 else if (upperPtr_)
396 {
397 lduAddr().map(*upperPtr_, *lowerCSRPtr_);
398 }
399 else
400 {
402 << "lowerPtr_ and upperPtr_ unallocated"
403 << abort(FatalError);
405 }
406
407 return *lowerCSRPtr_;
408}
409
410
412{
413 if (!lowerCSRPtr_)
414 {
415 const label nLower = lduAddr().losortAddr().size();
416
417 lowerCSRPtr_ = std::make_unique<scalarField>(nLower);
418
419 if (lowerPtr_)
420 {
421 lduAddr().map(*lowerPtr_, *lowerCSRPtr_);
422 }
423 else if (upperPtr_)
424 {
425 lduAddr().map(*upperPtr_, *lowerCSRPtr_);
426 }
427 else
428 {
430 << "lowerPtr_ and upperPtr_ unallocated"
431 << abort(FatalError);
433 }
434
435 return *lowerCSRPtr_;
436}
437
438
440{
441 if (!workPtr_ || workPtr_->size() != size)
442 {
443 workPtr_ = std::make_unique<solveScalarField>(size, Foam::zero{});
444 }
445
446 return *workPtr_;
447}
448
449
451{
452 if (!workPtr_)
453 {
455 << "workPtr_ unallocated"
457 }
458
459 return *workPtr_;
460}
461
462
464(
465 const scalarField& residual,
466 const word& fieldName,
467 const bool initial
468) const
469{
470 if (!mesh().hasDb())
471 {
472 return;
473 }
474
475 scalarIOField* residualPtr =
477 (
478 initial
479 ? IOobject::scopedName("initialResidual", fieldName)
480 : IOobject::scopedName("residual", fieldName)
481 );
482
483 if (residualPtr)
484 {
485 const auto* dataPtr = mesh().thisDb().findObject<meshState>("data");
486
487 if (dataPtr)
488 {
489 if (initial && dataPtr->isFirstIteration())
490 {
491 *residualPtr = residual;
493 << "Setting residual field for first solver iteration "
494 << "for solver field: " << fieldName << endl;
495 }
496 }
497 else
498 {
499 *residualPtr = residual;
501 << "Setting residual field for solver field "
502 << fieldName << endl;
504 }
505}
506
507
508// * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
509
510Foam::Ostream& Foam::operator<<(Ostream& os, const lduMatrix& mat)
511{
512 os << mat.hasLower() << token::SPACE
513 << mat.hasDiag() << token::SPACE
514 << mat.hasUpper() << token::SPACE
515 << mat.hasLowerCSR() << token::SPACE;
516
517 if (mat.hasLower())
518 {
519 os << mat.lower();
520 }
521
522 if (mat.hasDiag())
523 {
524 os << mat.diag();
525 }
526
527 if (mat.hasUpper())
528 {
529 os << mat.upper();
530 }
531
532 if (mat.hasLowerCSR() && !mat.hasLower())
533 {
534 // Only send over if can not be reconstructed locally
535 os << mat.lowerCSR();
536 }
538 os.check(FUNCTION_NAME);
539
540 return os;
541}
542
543
544Foam::Ostream& Foam::operator<<
545(
546 Ostream& os,
547 const InfoProxy<lduMatrix>& iproxy
548)
549{
550 const auto& mat = *iproxy;
551
552 os << "Lower:" << Switch::name(mat.hasLower())
553 << " Diag:" << Switch::name(mat.hasDiag())
554 << " Upper:" << Switch::name(mat.hasUpper())
555 << " lowerCSR:" << Switch::name(mat.hasLowerCSR())
556 << endl;
557
558 if (mat.hasLower())
559 {
560 os << "lower:" << mat.lower().size() << endl;
561 }
562 if (mat.hasDiag())
563 {
564 os << "diag :" << mat.diag().size() << endl;
565 }
566 if (mat.hasUpper())
567 {
568 os << "upper:" << mat.upper().size() << endl;
569 }
570
571
572 //if (hasLower)
573 //{
574 // os << "lower contents:" << endl;
575 // forAll(mat.lower(), i)
576 // {
577 // os << "i:" << i << "\t" << mat.lower()[i] << endl;
578 // }
579 // os << endl;
580 //}
581 //if (hasDiag)
582 //{
583 // os << "diag contents:" << endl;
584 // forAll(mat.diag(), i)
585 // {
586 // os << "i:" << i << "\t" << mat.diag()[i] << endl;
587 // }
588 // os << endl;
589 //}
590 //if (hasUpper)
591 //{
592 // os << "upper contents:" << endl;
593 // forAll(mat.upper(), i)
594 // {
595 // os << "i:" << i << "\t" << mat.upper()[i] << endl;
596 // }
597 // os << endl;
598 //}
599
601
602 return os;
603}
604
605
606// ************************************************************************* //
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...
if(patchID !=-1)
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition IOstream.C:45
A helper class for outputting values to Ostream.
Definition InfoProxy.H:49
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
static const char * name(bool b) noexcept
A string representation of bool as "false" / "true".
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition fvMesh.H:376
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition lduMatrix.H:81
lduMatrix(const lduMesh &mesh)
Construct (without coefficients) for an LDU addressed mesh.
Definition lduMatrix.C:54
bool hasUpper() const noexcept
Definition lduMatrix.H:828
const lduAddressing & lduAddr() const
Return the LDU addressing.
Definition lduMatrix.H:769
normTypes
Enumerated matrix normalisation types.
Definition lduMatrix.H:125
@ DEFAULT_NORM
"default" norm (== L1_scaled)
Definition lduMatrix.H:127
void residual(solveScalarField &rA, const solveScalarField &psi, const scalarField &source, const FieldField< Field, scalar > &interfaceBouCoeffs, const lduInterfaceFieldPtrsList &interfaces, const direction cmpt) const
const scalarField & diag() const
Definition lduMatrix.C:195
const scalarField & upper() const
Definition lduMatrix.C:235
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:457
const scalarField & lowerCSR() const
Definition lduMatrix.C:376
static const Enum< normTypes > normTypesNames_
Names for the normTypes.
Definition lduMatrix.H:134
bool hasLower() const noexcept
Definition lduMatrix.H:829
const scalarField & lower() const
Definition lduMatrix.C:306
const lduMesh & mesh() const noexcept
Return the LDU mesh from which the addressing is obtained.
Definition lduMatrix.H:753
static const scalar defaultTolerance
Default (absolute) tolerance (1e-6).
Definition lduMatrix.H:144
word matrixTypeName() const
The matrix type (empty, diagonal, symmetric, ...).
Definition lduMatrix.C:178
const solveScalarField & work() const
Work array.
Definition lduMatrix.C:443
bool hasDiag() const noexcept
Definition lduMatrix.H:827
bool hasLowerCSR() const noexcept
Definition lduMatrix.H:803
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition lduMesh.H:54
Database for mesh data, solution data, solver performance and other reduced data.
Definition meshState.H:56
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
Return non-const pointer to the object of the given Type, using a const-cast to have it behave like a...
@ SPACE
Space [isspace].
Definition token.H:144
A class for handling words, derived from Foam::string.
Definition word.H:66
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
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).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< solveScalar > solveScalarField
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
IOField< scalar > scalarIOField
IO for a Field of scalar.
volScalarField & e