LduMatrix.H
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
27Class
28 Foam::LduMatrix
29
30Description
31 LduMatrix is a general matrix class in which the coefficients are
32 stored as three arrays, one for the upper triangle, one for the
33 lower triangle and a third for the diagonal.
34
35 Addressing arrays must be supplied for the upper and lower triangles.
36
37Note
38 It might be better if this class were organised as a hierarchy starting
39 from an empty matrix, then deriving diagonal, symmetric and asymmetric
40 matrices.
41
42SourceFiles
43 LduMatrixATmul.C
44 LduMatrix.C
45 LduMatrixOperations.C
46 LduMatrixSolver.C
47 LduMatrixPreconditioner.C
48 LduMatrixTests.C
49 LduMatrixUpdateMatrixInterfaces.C
50
51\*---------------------------------------------------------------------------*/
52
53#ifndef Foam_LduMatrix_H
54#define Foam_LduMatrix_H
55
56#include "lduMesh.H"
57#include "Field.H"
58#include "FieldField.H"
60#include "SolverPerformance.H"
61#include "typeInfo.H"
62#include "autoPtr.H"
64
65// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66
67namespace Foam
68{
69
70// Forward Declarations
71
72template<class Type, class DType, class LUType>
73class LduMatrix;
74
75template<class Type, class DType, class LUType>
76Ostream& operator<<
77(
78 Ostream&,
80);
81
82
83/*---------------------------------------------------------------------------*\
84 Class LduMatrix Declaration
85\*---------------------------------------------------------------------------*/
86
87template<class Type, class DType, class LUType>
88class LduMatrix
89{
90 // Private Data
91
92 //- LDU mesh reference
93 const lduMesh& lduMesh_;
94
95 //- Diagonal coefficients
96 Field<DType> *diagPtr_;
97
98 //- Off-diagonal coefficients
99 Field<LUType> *upperPtr_, *lowerPtr_;
100
101 //- Source
102 Field<Type> *sourcePtr_;
103
104 //- Field interfaces (processor patches etc.)
106
107 //- Off-diagonal coefficients for interfaces
108 FieldField<Field, LUType> interfacesUpper_, interfacesLower_;
109
110
111public:
112
113 friend class SolverPerformance<Type>;
114
115 //- Abstract base-class for LduMatrix solvers
116 class solver
117 {
118 protected:
119
120 // Protected Data
121
122 //- Default maximum number of iterations in the solver
123 static const label defaultMaxIter_ = 1000;
127
128 //- Dictionary of controls
130
131 //- Level of verbosity in the solver output statements
132 int log_;
133
134 //- Minimum number of iterations in the solver
135 label minIter_;
136
137 //- Maximum number of iterations in the solver
138 label maxIter_;
139
140 //- Final convergence tolerance
141 Type tolerance_;
142
143 //- Convergence tolerance relative to the initial
144 Type relTol_;
145
146
147 // Protected Member Functions
148
149 //- Read the control parameters from the controlDict_
150 virtual void readControls();
151
152
153 // Housekeeping
154
155 //- Deprecated(2021-09) Read control parameter from dictionary
156 // \deprecated(2021-09) - use dictionary methods directly
157 template<class T>
158 void readControl(const dictionary& dict, T& val, const word& key)
159 {
160 dict.readIfPresent(key, val);
161 }
162
163
164 public:
165
166 //- Runtime type information
167 virtual const word& type() const = 0;
168
169
170 // Declare run-time constructor selection tables
173 (
174 autoPtr,
175 solver,
176 symMatrix,
177 (
178 const word& fieldName,
180 const dictionary& solverDict
181 ),
182 (
183 fieldName,
184 matrix,
185 solverDict
186 )
187 );
190 (
191 autoPtr,
192 solver,
193 asymMatrix,
194 (
195 const word& fieldName,
197 const dictionary& solverDict
198 ),
199 (
200 fieldName,
201 matrix,
202 solverDict
203 )
204 );
205
206
207 // Constructors
208
209 solver
210 (
211 const word& fieldName,
213 const dictionary& solverDict
214 );
215
216
217 // Selectors
218
219 //- Return a new solver
220 static autoPtr<solver> New
221 (
222 const word& fieldName,
224 const dictionary& solverDict
225 );
226
227
228 // Destructor
230 virtual ~solver() = default;
231
232
233 // Member Functions
235 const word& fieldName() const noexcept
236 {
237 return fieldName_;
238 }
241 {
242 return matrix_;
243 }
244
245
246 //- Read and reset the solver parameters from the given dictionary
247 virtual void read(const dictionary& solverDict);
250 (
252 ) const = 0;
253
254 //- Return the matrix norm used to normalise the residual for the
255 // stopping criterion
256 Type normFactor
257 (
258 const Field<Type>& psi,
259 const Field<Type>& Apsi,
260 Field<Type>& tmpField
261 ) const;
262 };
263
264
265 //- Abstract base-class for LduMatrix smoothers
266 class smoother
267 {
268 protected:
269
270 // Protected Data
274
275
276 public:
277
278 //- Runtime type information
279 virtual const word& type() const = 0;
280
281
282 // Declare run-time constructor selection tables
285 (
286 autoPtr,
287 smoother,
288 symMatrix,
289 (
290 const word& fieldName,
292 ),
293 (
294 fieldName,
295 matrix
296 )
297 );
300 (
301 autoPtr,
302 smoother,
303 asymMatrix,
304 (
305 const word& fieldName,
307 ),
308 (
309 fieldName,
310 matrix
311 )
312 );
313
314
315 // Constructors
316
318 (
319 const word& fieldName,
321 );
322
323
324 // Selectors
325
326 //- Return a new smoother
328 (
329 const word& fieldName,
331 const dictionary& smootherDict
332 );
333
334
335 // Destructor
337 virtual ~smoother() = default;
338
339
340 // Member Functions
342 const word& fieldName() const noexcept
343 {
344 return fieldName_;
345 }
348 {
349 return matrix_;
350 }
351
352
353 //- Smooth the solution for a given number of sweeps
354 virtual void smooth
355 (
357 const label nSweeps
358 ) const = 0;
359 };
360
361
362 //- Abstract base-class for LduMatrix preconditioners
363 class preconditioner
364 {
365 protected:
366
367 // Protected Data
368
369 //- Reference to the base-solver this preconditioner is used with
370 const solver& solver_;
371
372
373 public:
374
375 //- Runtime type information
376 virtual const word& type() const = 0;
377
378
379 // Declare run-time constructor selection tables
382 (
383 autoPtr,
385 symMatrix,
386 (
387 const solver& sol,
388 const dictionary& preconditionerDict
389 ),
390 (sol, preconditionerDict)
391 );
394 (
395 autoPtr,
397 asymMatrix,
398 (
399 const solver& sol,
400 const dictionary& preconditionerDict
401 ),
402 (sol, preconditionerDict)
403 );
404
405
406 // Constructors
409 (
410 const solver& sol
411 )
412 :
413 solver_(sol)
414 {}
415
416
417 // Selectors
418
419 //- Return a new preconditioner
421 (
422 const solver& sol,
423 const dictionary& preconditionerDict
424 );
425
426
427 // Destructor
429 virtual ~preconditioner() = default;
430
431
432 // Member functions
433
434 //- Read and reset the preconditioner parameters
435 // from the given dictionary
436 virtual void read(const dictionary& preconditionerDict)
437 {}
438
439 //- Return wA the preconditioned form of residual rA
440 virtual void precondition
441 (
442 Field<Type>& wA,
443 const Field<Type>& rA
444 ) const = 0;
445
446 //- Return wT the transpose-matrix preconditioned form of
447 // residual rT.
448 // This is only required for preconditioning asymmetric matrices.
449 virtual void preconditionT
450 (
451 Field<Type>& wT,
452 const Field<Type>& rT
453 ) const
454 {
456 }
457 };
458
459
460 // Static Data
461
462 // Declare name of the class and its debug switch
463 ClassName("LduMatrix");
464
465
466 // Constructors
467
468 //- Construct given an LDU addressed mesh.
469 // The coefficients are initially empty for subsequent setting.
470 LduMatrix(const lduMesh&);
471
472 //- Construct as copy
474
475 //- Construct as copy or re-use as specified.
477
478 //- Construct given an LDU addressed mesh and an Istream
479 // from which the coefficients are read
480 LduMatrix(const lduMesh&, Istream&);
481
482
483 // Destructor
484
485 ~LduMatrix();
486
487
488 // Member functions
489
490 // Access to addressing
491
492 //- Return the LDU mesh from which the addressing is obtained
493 const lduMesh& mesh() const
494 {
495 return lduMesh_;
496 }
497
498 //- Return the LDU addressing
499 const lduAddressing& lduAddr() const
500 {
501 return lduMesh_.lduAddr();
502 }
503
504 //- Return the patch evaluation schedule
505 const lduSchedule& patchSchedule() const
506 {
507 return lduAddr().patchSchedule();
508 }
509
510 //- Return interfaces
512 {
513 return interfaces_;
514 }
515
516 //- Return interfaces
518 {
519 return interfaces_;
520 }
521
522
523 // Access to coefficients
524
531 {
532 return interfacesUpper_;
533 }
536 {
537 return interfacesLower_;
538 }
539
540
541 const Field<DType>& diag() const;
542 const Field<LUType>& upper() const;
543 const Field<LUType>& lower() const;
544 const Field<Type>& source() const;
547 {
548 return interfacesUpper_;
549 }
552 {
553 return interfacesLower_;
554 }
555
557 bool hasDiag() const
558 {
559 return (diagPtr_);
560 }
562 bool hasUpper() const
563 {
564 return (upperPtr_);
565 }
567 bool hasLower() const
568 {
569 return (lowerPtr_);
570 }
572 bool hasSource() const
573 {
574 return (sourcePtr_);
575 }
577 bool diagonal() const
578 {
579 return (diagPtr_ && !lowerPtr_ && !upperPtr_);
580 }
582 bool symmetric() const
583 {
584 return (diagPtr_ && (!lowerPtr_ && upperPtr_));
585 }
587 bool asymmetric() const
588 {
589 return (diagPtr_ && lowerPtr_ && upperPtr_);
590 }
591
592
593 // operations
594
595 void sumDiag();
596 void negSumDiag();
597
598 void sumMagOffDiag(Field<LUType>& sumOff) const;
599
600 //- Matrix multiplication
601 void Amul(Field<Type>&, const tmp<Field<Type>>&) const;
602
603 //- Matrix transpose multiplication
604 void Tmul(Field<Type>&, const tmp<Field<Type>>&) const;
605
606
607 //- Sum the coefficients on each row of the matrix
608 void sumA(Field<Type>&) const;
609
610
611 void residual(Field<Type>& rA, const Field<Type>& psi) const;
612
614
615
616 //- Initialise the update of interfaced interfaces
617 // for matrix operations
619 (
620 const bool add,
621 const FieldField<Field, LUType>& interfaceCoeffs,
622 const Field<Type>& psiif,
623 Field<Type>& result
624 ) const;
625
626 //- Update interfaced interfaces for matrix operations
628 (
629 const bool add,
630 const FieldField<Field, LUType>& interfaceCoeffs,
631 const Field<Type>& psiif,
632 Field<Type>& result
633 ) const;
634
635
636 tmp<Field<Type>> H(const Field<Type>&) const;
637 tmp<Field<Type>> H(const tmp<Field<Type>>&) const;
638
639 tmp<Field<Type>> faceH(const Field<Type>&) const;
640 tmp<Field<Type>> faceH(const tmp<Field<Type>>&) const;
641
642
643 // Member operators
644
646
647 void negate();
648
651
652 void operator*=(const scalarField&);
653 void operator*=(scalar);
654
655
656 // Ostream operator
658 friend Ostream& operator<< <Type, DType, LUType>
659 (
660 Ostream&,
662 );
663};
664
665
666// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
667
668} // End namespace Foam
669
670// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
672#define makeLduMatrix(Type, DType, LUType) \
673 \
674typedef Foam::LduMatrix<Type, DType, LUType> \
675 ldu##Type##DType##LUType##Matrix; \
676 \
677defineNamedTemplateTypeNameAndDebug(ldu##Type##DType##LUType##Matrix, 0); \
678 \
679 \
680typedef LduMatrix<Type, DType, LUType>::smoother \
681 ldu##Type##DType##LUType##Smoother; \
682 \
683defineTemplateRunTimeSelectionTable \
684( \
685 ldu##Type##DType##LUType##Smoother, \
686 symMatrix \
687); \
688 \
689defineTemplateRunTimeSelectionTable \
690( \
691 ldu##Type##DType##LUType##Smoother, \
692 asymMatrix \
693); \
694 \
695 \
696typedef LduMatrix<Type, DType, LUType>::preconditioner \
697 ldu##Type##DType##LUType##Preconditioner; \
698 \
699defineTemplateRunTimeSelectionTable \
700( \
701 ldu##Type##DType##LUType##Preconditioner, \
702 symMatrix \
703); \
704 \
705defineTemplateRunTimeSelectionTable \
706( \
707 ldu##Type##DType##LUType##Preconditioner, \
708 asymMatrix \
709); \
710 \
711 \
712typedef LduMatrix<Type, DType, LUType>::solver \
713 ldu##Type##DType##LUType##Solver; \
714 \
715defineTemplateRunTimeSelectionTable \
716( \
717 ldu##Type##DType##LUType##Solver, \
718 symMatrix \
719); \
720 \
721defineTemplateRunTimeSelectionTable \
722( \
723 ldu##Type##DType##LUType##Solver, \
724 asymMatrix \
725);
726
728#define makeLduPreconditioner(Precon, Type, DType, LUType) \
729 \
730typedef Precon<Type, DType, LUType> \
731 Precon##Type##DType##LUType##Preconditioner; \
732defineNamedTemplateTypeNameAndDebug \
733( \
734 Precon##Type##DType##LUType##Preconditioner, \
735 0 \
736);
738#define makeLduSymPreconditioner(Precon, Type, DType, LUType) \
739 \
740LduMatrix<Type, DType, LUType>::preconditioner:: \
741addsymMatrixConstructorToTable<Precon##Type##DType##LUType##Preconditioner> \
742add##Precon##Type##DType##LUType##PreconditionerSymMatrixConstructorToTable_;
744#define makeLduAsymPreconditioner(Precon, Type, DType, LUType) \
745 \
746LduMatrix<Type, DType, LUType>::preconditioner:: \
747addasymMatrixConstructorToTable<Precon##Type##DType##LUType##Preconditioner> \
748add##Precon##Type##DType##LUType##PreconditionerAsymMatrixConstructorToTable_;
749
751#define makeLduSmoother(Smoother, Type, DType, LUType) \
752 \
753typedef Smoother<Type, DType, LUType> \
754 Smoother##Type##DType##LUType##Smoother; \
755 \
756defineNamedTemplateTypeNameAndDebug \
757( \
758 Smoother##Type##DType##LUType##Smoother, \
759 0 \
760);
762#define makeLduSymSmoother(Smoother, Type, DType, LUType) \
763 \
764LduMatrix<Type, DType, LUType>::smoother:: \
765 addsymMatrixConstructorToTable<Smoother##Type##DType##LUType##Smoother> \
766 add##Smoother##Type##DType##LUType##SymMatrixConstructorToTable_;
768#define makeLduAsymSmoother(Smoother, Type, DType, LUType) \
769 \
770LduMatrix<Type, DType, LUType>::smoother:: \
771 addasymMatrixConstructorToTable<Smoother##Type##DType##LUType##Smoother> \
772 add##Smoother##Type##DType##LUType##AsymMatrixConstructorToTable_;
773
775#define makeLduSolver(Solver, Type, DType, LUType) \
776 \
777typedef Solver<Type, DType, LUType> \
778 Solver##Type##DType##LUType##Solver; \
779 \
780defineNamedTemplateTypeNameAndDebug \
781( \
782 Solver##Type##DType##LUType##Solver, \
783 0 \
784);
786#define makeLduSymSolver(Solver, Type, DType, LUType) \
787 \
788LduMatrix<Type, DType, LUType>::solver:: \
789 addsymMatrixConstructorToTable<Solver##Type##DType##LUType##Solver> \
790 add##Solver##Type##DType##LUType##SymMatrixConstructorToTable_;
792#define makeLduAsymSolver(Solver, Type, DType, LUType) \
793 \
794LduMatrix<Type, DType, LUType>::solver:: \
795 addasymMatrixConstructorToTable<Solver##Type##DType##LUType##Solver> \
796 add##Solver##Type##DType##LUType##AsymMatrixConstructorToTable_;
797
798
799// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
800
801#ifdef NoRepository
802 #include "LduMatrix.C"
803#endif
804
805// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
806
807#endif
808
809// ************************************************************************* //
List of coupled interface fields to be used in coupling.
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:80
Generic templated field type.
Definition: Field.H:82
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
Abstract base-class for LduMatrix preconditioners.
Definition: LduMatrix.H:363
const solver & solver_
Reference to the base-solver this preconditioner is used with.
Definition: LduMatrix.H:369
preconditioner(const solver &sol)
Definition: LduMatrix.H:408
declareRunTimeSelectionTable(autoPtr, preconditioner, symMatrix,(const solver &sol, const dictionary &preconditionerDict),(sol, preconditionerDict))
static autoPtr< preconditioner > New(const solver &sol, const dictionary &preconditionerDict)
Return a new preconditioner.
virtual void read(const dictionary &preconditionerDict)
Read and reset the preconditioner parameters.
Definition: LduMatrix.H:435
declareRunTimeSelectionTable(autoPtr, preconditioner, asymMatrix,(const solver &sol, const dictionary &preconditionerDict),(sol, preconditionerDict))
virtual void preconditionT(Field< Type > &wT, const Field< Type > &rT) const
Return wT the transpose-matrix preconditioned form of.
Definition: LduMatrix.H:449
virtual void precondition(Field< Type > &wA, const Field< Type > &rA) const =0
Return wA the preconditioned form of residual rA.
virtual const word & type() const =0
Runtime type information.
virtual ~preconditioner()=default
Abstract base-class for LduMatrix smoothers.
Definition: LduMatrix.H:266
declareRunTimeSelectionTable(autoPtr, smoother, symMatrix,(const word &fieldName, const LduMatrix< Type, DType, LUType > &matrix),(fieldName, matrix))
const LduMatrix< Type, DType, LUType > & matrix() const noexcept
Definition: LduMatrix.H:346
virtual ~smoother()=default
declareRunTimeSelectionTable(autoPtr, smoother, asymMatrix,(const word &fieldName, const LduMatrix< Type, DType, LUType > &matrix),(fieldName, matrix))
const LduMatrix< Type, DType, LUType > & matrix_
Definition: LduMatrix.H:272
virtual void smooth(Field< Type > &psi, const label nSweeps) const =0
Smooth the solution for a given number of sweeps.
virtual const word & type() const =0
Runtime type information.
const word & fieldName() const noexcept
Definition: LduMatrix.H:341
static autoPtr< smoother > New(const word &fieldName, const LduMatrix< Type, DType, LUType > &matrix, const dictionary &smootherDict)
Return a new smoother.
Abstract base-class for LduMatrix solvers.
Definition: LduMatrix.H:116
label maxIter_
Maximum number of iterations in the solver.
Definition: LduMatrix.H:137
static const label defaultMaxIter_
Default maximum number of iterations in the solver.
Definition: LduMatrix.H:122
const LduMatrix< Type, DType, LUType > & matrix() const noexcept
Definition: LduMatrix.H:239
Type relTol_
Convergence tolerance relative to the initial.
Definition: LduMatrix.H:143
void readControl(const dictionary &dict, T &val, const word &key)
Deprecated(2021-09) Read control parameter from dictionary.
Definition: LduMatrix.H:157
virtual void read(const dictionary &solverDict)
Read and reset the solver parameters from the given dictionary.
virtual SolverPerformance< Type > solve(Field< Type > &psi) const =0
Type normFactor(const Field< Type > &psi, const Field< Type > &Apsi, Field< Type > &tmpField) const
Return the matrix norm used to normalise the residual for the.
label minIter_
Minimum number of iterations in the solver.
Definition: LduMatrix.H:134
Type tolerance_
Final convergence tolerance.
Definition: LduMatrix.H:140
declareRunTimeSelectionTable(autoPtr, solver, asymMatrix,(const word &fieldName, const LduMatrix< Type, DType, LUType > &matrix, const dictionary &solverDict),(fieldName, matrix, solverDict))
int log_
Level of verbosity in the solver output statements.
Definition: LduMatrix.H:131
virtual ~solver()=default
virtual void readControls()
Read the control parameters from the controlDict_.
const LduMatrix< Type, DType, LUType > & matrix_
Definition: LduMatrix.H:125
declareRunTimeSelectionTable(autoPtr, solver, symMatrix,(const word &fieldName, const LduMatrix< Type, DType, LUType > &matrix, const dictionary &solverDict),(fieldName, matrix, solverDict))
static autoPtr< solver > New(const word &fieldName, const LduMatrix< Type, DType, LUType > &matrix, const dictionary &solverDict)
Return a new solver.
dictionary controlDict_
Dictionary of controls.
Definition: LduMatrix.H:128
virtual const word & type() const =0
Runtime type information.
const word & fieldName() const noexcept
Definition: LduMatrix.H:234
LduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition: LduMatrix.H:88
FieldField< Field, LUType > & interfacesLower()
Definition: LduMatrix.H:534
bool symmetric() const
Definition: LduMatrix.H:581
bool diagonal() const
Definition: LduMatrix.H:576
void Amul(Field< Type > &, const tmp< Field< Type > > &) const
Matrix multiplication.
LduMatrix(const LduMatrix< Type, DType, LUType > &)
Construct as copy.
FieldField< Field, LUType > & interfacesUpper()
Definition: LduMatrix.H:529
void sumA(Field< Type > &) const
Sum the coefficients on each row of the matrix.
const lduAddressing & lduAddr() const
Return the LDU addressing.
Definition: LduMatrix.H:498
const LduInterfaceFieldPtrsList< Type > & interfaces() const
Return interfaces.
Definition: LduMatrix.H:510
Field< Type > & source()
Definition: LduMatrix.C:250
bool hasUpper() const
Definition: LduMatrix.H:561
void updateMatrixInterfaces(const bool add, const FieldField< Field, LUType > &interfaceCoeffs, const Field< Type > &psiif, Field< Type > &result) const
Update interfaced interfaces for matrix operations.
tmp< Field< Type > > H(const Field< Type > &) const
bool hasLower() const
Definition: LduMatrix.H:566
Field< LUType > & upper()
Definition: LduMatrix.C:204
bool hasSource() const
Definition: LduMatrix.H:571
LduMatrix(LduMatrix< Type, DType, LUType > &, bool reuse)
Construct as copy or re-use as specified.
const FieldField< Field, LUType > & interfacesUpper() const
Definition: LduMatrix.H:545
void operator+=(const LduMatrix< Type, DType, LUType > &)
const lduMesh & mesh() const
Return the LDU mesh from which the addressing is obtained.
Definition: LduMatrix.H:492
LduInterfaceFieldPtrsList< Type > & interfaces()
Return interfaces.
Definition: LduMatrix.H:516
const lduSchedule & patchSchedule() const
Return the patch evaluation schedule.
Definition: LduMatrix.H:504
ClassName("LduMatrix")
void operator=(const LduMatrix< Type, DType, LUType > &)
void operator*=(const scalarField &)
void sumMagOffDiag(Field< LUType > &sumOff) const
void residual(Field< Type > &rA, const Field< Type > &psi) const
tmp< Field< Type > > faceH(const Field< Type > &) const
void initMatrixInterfaces(const bool add, const FieldField< Field, LUType > &interfaceCoeffs, const Field< Type > &psiif, Field< Type > &result) const
Initialise the update of interfaced interfaces.
const FieldField< Field, LUType > & interfacesLower() const
Definition: LduMatrix.H:550
void operator-=(const LduMatrix< Type, DType, LUType > &)
void Tmul(Field< Type > &, const tmp< Field< Type > > &) const
Matrix transpose multiplication.
bool hasDiag() const
Definition: LduMatrix.H:556
Field< LUType > & lower()
Definition: LduMatrix.C:227
Field< DType > & diag()
Definition: LduMatrix.C:192
bool asymmetric() const
Definition: LduMatrix.H:586
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: UPtrList.H:71
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
The class contains the addressing required by the lduMatrix: upper, lower and losort.
virtual const lduSchedule & patchSchedule() const =0
Return patch field evaluation schedule.
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:63
virtual const lduAddressing & lduAddr() const =0
Return ldu addressing.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define ClassName(TypeNameString)
Add typeName information from argument TypeNameString to a class.
Definition: className.H:67
const volScalarField & psi
const volScalarField & T
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:517
Namespace for OpenFOAM.
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
const direction noexcept
Definition: Scalar.H:223
Macros to ease declaration of run-time selection tables.
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes)
dictionary dict
CEqn solve()