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) 2016-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
37 It might be better if this class were organised as a hierarchy starting
38 from an empty matrix, then deriving diagonal, symmetric and asymmetric
39 matrices.
40
41SourceFiles
42 lduMatrixATmul.C
43 lduMatrix.C
44 lduMatrixTemplates.C
45 lduMatrixOperations.C
46 lduMatrixSolver.C
47 lduMatrixPreconditioner.C
48 lduMatrixTests.C
49 lduMatrixUpdateMatrixInterfaces.C
50
51\*---------------------------------------------------------------------------*/
52
53#ifndef lduMatrix_H
54#define lduMatrix_H
55
56#include "lduMesh.H"
57#include "primitiveFieldsFwd.H"
58#include "FieldField.H"
60#include "typeInfo.H"
61#include "autoPtr.H"
63#include "solverPerformance.H"
64#include "InfoProxy.H"
65#include "profilingTrigger.H"
66
67// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68
69namespace Foam
70{
71
72// Forward Declarations
73
74class lduMatrix;
75
76Ostream& operator<<(Ostream&, const lduMatrix&);
77Ostream& operator<<(Ostream&, const InfoProxy<lduMatrix>&);
78
79
80/*---------------------------------------------------------------------------*\
81 Class lduMatrix Declaration
82\*---------------------------------------------------------------------------*/
84class lduMatrix
85{
86 // Private Data
87
88 //- LDU mesh reference
89 //const lduMesh& lduMesh_;
90 std::reference_wrapper<const lduMesh> lduMesh_;
91
92 //- Coefficients (not including interfaces)
93 scalarField *lowerPtr_, *diagPtr_, *upperPtr_;
94
95
96public:
97
98 //- Abstract base-class for lduMatrix solvers
99 class solver
100 {
101 protected:
102
103 // Protected Data
104
105 //- Default maximum number of iterations in the solver
106 static const label defaultMaxIter_;
113
114 //- Dictionary of controls
116
117 //- Level of verbosity in the solver output statements
118 int log_;
119
120 //- Minimum number of iterations in the solver
121 label minIter_;
122
123 //- Maximum number of iterations in the solver
124 label maxIter_;
125
126 //- Final convergence tolerance
127 scalar tolerance_;
128
129 //- Convergence tolerance relative to the initial
130 scalar relTol_;
133
134
135 // Protected Member Functions
136
137 //- Read the control parameters from the controlDict_
138 virtual void readControls();
139
140
141 public:
142
143 //- Runtime type information
144 virtual const word& type() const = 0;
145
146
147 // Declare run-time constructor selection tables
150 (
151 autoPtr,
152 solver,
153 symMatrix,
154 (
155 const word& fieldName,
156 const lduMatrix& matrix,
160 const dictionary& solverControls
161 ),
162 (
163 fieldName,
164 matrix,
168 solverControls
169 )
170 );
173 (
174 autoPtr,
175 solver,
176 asymMatrix,
177 (
178 const word& fieldName,
179 const lduMatrix& matrix,
183 const dictionary& solverControls
184 ),
185 (
186 fieldName,
187 matrix,
191 solverControls
192 )
193 );
194
195
196 // Constructors
197
198 solver
199 (
200 const word& fieldName,
201 const lduMatrix& matrix,
205 const dictionary& solverControls
206 );
207
208 // Selectors
209
210 //- Return a new solver
211 static autoPtr<solver> New
212 (
213 const word& fieldName,
214 const lduMatrix& matrix,
218 const dictionary& solverControls
219 );
220
221
222
223 //- Destructor
224 virtual ~solver() = default;
225
226
227 // Member Functions
229 const word& fieldName() const noexcept
230 {
231 return fieldName_;
232 }
234 const lduMatrix& matrix() const noexcept
235 {
236 return matrix_;
237 }
240 {
241 return interfaceBouCoeffs_;
242 }
245 {
246 return interfaceIntCoeffs_;
247 }
250 {
251 return interfaces_;
252 }
253
254
255 //- Read and reset the solver parameters from the given stream
256 virtual void read(const dictionary&);
257
258 //- Solve with given field and rhs
260 (
262 const scalarField& source,
263 const direction cmpt=0
264 ) const = 0;
265
266 //- Solve with given field and rhs (in solveScalar precision).
267 // Default is to call solve routine
269 (
271 const solveScalarField& source,
272 const direction cmpt=0
273 ) const;
274
275 //- Return the matrix norm used to normalise the residual for the
276 //- stopping criterion
278 (
279 const solveScalarField& psi,
280 const solveScalarField& source,
281 const solveScalarField& Apsi,
282 solveScalarField& tmpField
283 ) const;
284 };
285
286
287 //- Abstract base-class for lduMatrix smoothers
288 class smoother
289 {
290 protected:
291
292 // Protected Data
299
300
301 public:
302
303 //- Find the smoother name (directly or from a sub-dictionary)
304 static word getName(const dictionary&);
305
306 //- Runtime type information
307 virtual const word& type() const = 0;
308
309
310 // Declare run-time constructor selection tables
313 (
314 autoPtr,
315 smoother,
316 symMatrix,
317 (
318 const word& fieldName,
319 const lduMatrix& matrix,
323 ),
324 (
325 fieldName,
326 matrix,
330 )
331 );
334 (
335 autoPtr,
336 smoother,
337 asymMatrix,
338 (
339 const word& fieldName,
340 const lduMatrix& matrix,
344 ),
345 (
346 fieldName,
347 matrix,
351 )
352 );
353
354
355 // Constructors
356
358 (
359 const word& fieldName,
360 const lduMatrix& matrix,
364 );
365
366
367 // Selectors
368
369 //- Return a new smoother
371 (
372 const word& fieldName,
373 const lduMatrix& matrix,
377 const dictionary& solverControls
378 );
379
380
381 //- Destructor
382 virtual ~smoother() = default;
383
384
385 // Member Functions
387 const word& fieldName() const noexcept
388 {
389 return fieldName_;
390 }
392 const lduMatrix& matrix() const noexcept
393 {
394 return matrix_;
395 }
398 {
399 return interfaceBouCoeffs_;
400 }
403 {
404 return interfaceIntCoeffs_;
405 }
408 {
409 return interfaces_;
410 }
411
412
413 //- Smooth the solution for a given number of sweeps
414 virtual void smooth
415 (
417 const scalarField& source,
418 const direction cmpt,
419 const label nSweeps
420 ) const = 0;
421
422 //- Smooth the solution for a given number of sweeps
423 virtual void scalarSmooth
424 (
426 const solveScalarField& source,
427 const direction cmpt,
428 const label nSweeps
429 ) const = 0;
430 };
431
432
433 //- Abstract base-class for lduMatrix preconditioners
434 class preconditioner
435 {
436 protected:
437
438 // Protected Data
439
440 //- Reference to the base-solver this preconditioner is used with
441 const solver& solver_;
442
443
444 public:
445
446 //- Find the preconditioner name (directly or from a sub-dictionary)
447 static word getName(const dictionary&);
448
449 //- Runtime type information
450 virtual const word& type() const = 0;
451
452
453 // Declare run-time constructor selection tables
456 (
457 autoPtr,
459 symMatrix,
460 (
461 const solver& sol,
462 const dictionary& solverControls
463 ),
464 (sol, solverControls)
465 );
468 (
469 autoPtr,
471 asymMatrix,
472 (
473 const solver& sol,
474 const dictionary& solverControls
475 ),
476 (sol, solverControls)
477 );
478
479
480 // Constructors
483 (
484 const solver& sol
485 )
486 :
487 solver_(sol)
488 {}
489
490
491 // Selectors
492
493 //- Return a new preconditioner
495 (
496 const solver& sol,
497 const dictionary& solverControls
498 );
499
500
501 //- Destructor
502 virtual ~preconditioner() = default;
503
504
505 // Member Functions
506
507 //- Read and reset the preconditioner parameters
508 //- from the given stream
509 virtual void read(const dictionary&)
510 {}
511
512 //- Return wA the preconditioned form of residual rA
513 virtual void precondition
514 (
516 const solveScalarField& rA,
517 const direction cmpt=0
518 ) const = 0;
519
520 //- Return wT the transpose-matrix preconditioned form of
521 //- residual rT.
522 // This is only required for preconditioning asymmetric matrices.
523 virtual void preconditionT
524 (
526 const solveScalarField& rT,
527 const direction cmpt=0
528 ) const
529 {
531 }
532 };
533
534
535 // Static Data
536
537 // Declare name of the class and its debug switch
538 ClassName("lduMatrix");
539
540
541 // Constructors
542
543 //- Construct given an LDU addressed mesh.
544 // The coefficients are initially empty for subsequent setting.
545 lduMatrix(const lduMesh&);
546
547 //- Construct as copy
548 lduMatrix(const lduMatrix&);
549
550 //- Construct as copy or re-use as specified.
551 lduMatrix(lduMatrix&, bool reuse);
552
553 //- Construct given an LDU addressed mesh and an Istream
554 //- from which the coefficients are read
555 lduMatrix(const lduMesh&, Istream&);
556
557
558 //- Destructor
559 ~lduMatrix();
560
561
562 // Member Functions
563
564 // Access to addressing
565
566 //- Return the LDU mesh from which the addressing is obtained
567 const lduMesh& mesh() const
568 {
569 return lduMesh_;
570 }
571
572 //- Set the LDU mesh containing the addressing is obtained
573 void setLduMesh(const lduMesh& m)
574 {
575 lduMesh_ = m;
576 }
577
578 //- Return the LDU addressing
579 const lduAddressing& lduAddr() const
580 {
581 return mesh().lduAddr();
582 }
583
584 //- Return the patch evaluation schedule
585 const lduSchedule& patchSchedule() const
586 {
587 return lduAddr().patchSchedule();
588 }
589
590
591 // Access to coefficients
592
594 scalarField& diag();
596
597 // Size with externally provided sizes (for constructing with 'fake'
598 // mesh in GAMG)
599
600 scalarField& lower(const label size);
601 scalarField& diag(const label nCoeffs);
602 scalarField& upper(const label nCoeffs);
603
604
605 const scalarField& lower() const;
606 const scalarField& diag() const;
607 const scalarField& upper() const;
609 bool hasDiag() const
610 {
611 return (diagPtr_);
612 }
614 bool hasUpper() const
615 {
616 return (upperPtr_);
617 }
619 bool hasLower() const
620 {
621 return (lowerPtr_);
622 }
624 bool diagonal() const
625 {
626 return (diagPtr_ && !lowerPtr_ && !upperPtr_);
627 }
629 bool symmetric() const
630 {
631 return (diagPtr_ && (!lowerPtr_ && upperPtr_));
632 }
634 bool asymmetric() const
635 {
636 return (diagPtr_ && lowerPtr_ && upperPtr_);
637 }
638
639
640 // operations
641
642 void sumDiag();
643 void negSumDiag();
644
645 void sumMagOffDiag(scalarField& sumOff) const;
646
647 //- Matrix multiplication with updated interfaces.
648 void Amul
649 (
654 const direction cmpt
655 ) const;
656
657 //- Matrix transpose multiplication with updated interfaces.
658 void Tmul
659 (
664 const direction cmpt
665 ) const;
666
667
668 //- Sum the coefficients on each row of the matrix
669 void sumA
670 (
674 ) const;
675
676
677 void residual
678 (
680 const solveScalarField& psi,
681 const scalarField& source,
682 const FieldField<Field, scalar>& interfaceBouCoeffs,
683 const lduInterfaceFieldPtrsList& interfaces,
684 const direction cmpt
685 ) const;
686
688 (
689 const solveScalarField& psi,
690 const scalarField& source,
691 const FieldField<Field, scalar>& interfaceBouCoeffs,
692 const lduInterfaceFieldPtrsList& interfaces,
693 const direction cmpt
694 ) const;
695
696
697 //- Initialise the update of interfaced interfaces
698 //- for matrix operations
700 (
701 const bool add,
702 const FieldField<Field, scalar>& interfaceCoeffs,
703 const lduInterfaceFieldPtrsList& interfaces,
704 const solveScalarField& psiif,
705 solveScalarField& result,
706 const direction cmpt
707 ) const;
708
709 //- Update interfaced interfaces for matrix operations
711 (
712 const bool add,
713 const FieldField<Field, scalar>& interfaceCoeffs,
714 const lduInterfaceFieldPtrsList& interfaces,
715 const solveScalarField& psiif,
716 solveScalarField& result,
717 const direction cmpt,
718 const label startRequest // starting request (for non-blocking)
719 ) const;
720
721 //- Set the residual field using an IOField on the object registry
722 //- if it exists
724 (
725 const scalarField& residual,
726 const word& fieldName,
727 const bool initial
728 ) const;
729
730 template<class Type>
731 tmp<Field<Type>> H(const Field<Type>&) const;
732
733 template<class Type>
734 tmp<Field<Type>> H(const tmp<Field<Type>>&) const;
735
736 tmp<scalarField> H1() const;
737
738 template<class Type>
739 tmp<Field<Type>> faceH(const Field<Type>&) const;
740
741 template<class Type>
742 tmp<Field<Type>> faceH(const tmp<Field<Type>>&) const;
743
744
745 // Info
746
747 //- Return info proxy.
748 // Used to print matrix information to a stream
750 {
751 return *this;
752 }
753
754
755 // Member operators
756
757 void operator=(const lduMatrix&);
758
759 void negate();
760
761 void operator+=(const lduMatrix&);
762 void operator-=(const lduMatrix&);
763
764 void operator*=(const scalarField&);
765 void operator*=(scalar);
766
767
768 // Ostream operator
772};
773
774
775// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
776
777} // End namespace Foam
778
779// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
780
781#ifdef NoRepository
782 #include "lduMatrixTemplates.C"
783#endif
784
785// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
786
787#endif
788
789// ************************************************************************* //
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:80
A helper class for outputting values to Ostream.
Definition: InfoProxy.H:52
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
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
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
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
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 lduMatrix preconditioners.
Definition: lduMatrix.H:434
declareRunTimeSelectionTable(autoPtr, preconditioner, symMatrix,(const solver &sol, const dictionary &solverControls),(sol, solverControls))
const solver & solver_
Reference to the base-solver this preconditioner is used with.
Definition: lduMatrix.H:440
static autoPtr< preconditioner > New(const solver &sol, const dictionary &solverControls)
Return a new preconditioner.
preconditioner(const solver &sol)
Definition: lduMatrix.H:482
declareRunTimeSelectionTable(autoPtr, preconditioner, asymMatrix,(const solver &sol, const dictionary &solverControls),(sol, solverControls))
virtual void precondition(solveScalarField &wA, const solveScalarField &rA, const direction cmpt=0) const =0
Return wA the preconditioned form of residual rA.
virtual void preconditionT(solveScalarField &wT, const solveScalarField &rT, const direction cmpt=0) const
Definition: lduMatrix.H:523
virtual void read(const dictionary &)
Definition: lduMatrix.H:508
static word getName(const dictionary &)
Find the preconditioner name (directly or from a sub-dictionary)
virtual const word & type() const =0
Runtime type information.
virtual ~preconditioner()=default
Destructor.
Abstract base-class for lduMatrix smoothers.
Definition: lduMatrix.H:288
const FieldField< Field, scalar > & interfaceIntCoeffs() const noexcept
Definition: lduMatrix.H:401
const lduMatrix & matrix_
Definition: lduMatrix.H:294
const lduInterfaceFieldPtrsList & interfaces_
Definition: lduMatrix.H:297
const lduInterfaceFieldPtrsList & interfaces() const noexcept
Definition: lduMatrix.H:406
virtual ~smoother()=default
Destructor.
declareRunTimeSelectionTable(autoPtr, smoother, symMatrix,(const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces),(fieldName, matrix, interfaceBouCoeffs, interfaceIntCoeffs, interfaces))
const FieldField< Field, scalar > & interfaceBouCoeffs_
Definition: lduMatrix.H:295
virtual void scalarSmooth(solveScalarField &psi, const solveScalarField &source, const direction cmpt, const label nSweeps) const =0
Smooth the solution for a given number of sweeps.
const lduMatrix & matrix() const noexcept
Definition: lduMatrix.H:391
static autoPtr< smoother > New(const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces, const dictionary &solverControls)
Return a new smoother.
const FieldField< Field, scalar > & interfaceIntCoeffs_
Definition: lduMatrix.H:296
virtual void smooth(solveScalarField &psi, const scalarField &source, const direction cmpt, const label nSweeps) const =0
Smooth the solution for a given number of sweeps.
static word getName(const dictionary &)
Find the smoother name (directly or from a sub-dictionary)
const FieldField< Field, scalar > & interfaceBouCoeffs() const noexcept
Definition: lduMatrix.H:396
declareRunTimeSelectionTable(autoPtr, smoother, asymMatrix,(const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces),(fieldName, matrix, interfaceBouCoeffs, interfaceIntCoeffs, interfaces))
virtual const word & type() const =0
Runtime type information.
const word & fieldName() const noexcept
Definition: lduMatrix.H:386
Abstract base-class for lduMatrix solvers.
Definition: lduMatrix.H:99
const FieldField< Field, scalar > & interfaceIntCoeffs() const noexcept
Definition: lduMatrix.H:243
const lduMatrix & matrix_
Definition: lduMatrix.H:108
declareRunTimeSelectionTable(autoPtr, solver, symMatrix,(const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces, const dictionary &solverControls),(fieldName, matrix, interfaceBouCoeffs, interfaceIntCoeffs, interfaces, solverControls))
label maxIter_
Maximum number of iterations in the solver.
Definition: lduMatrix.H:123
profilingTrigger profiling_
Definition: lduMatrix.H:131
scalar tolerance_
Final convergence tolerance.
Definition: lduMatrix.H:126
const lduInterfaceFieldPtrsList & interfaces() const noexcept
Definition: lduMatrix.H:248
lduInterfaceFieldPtrsList interfaces_
Definition: lduMatrix.H:111
const FieldField< Field, scalar > & interfaceBouCoeffs_
Definition: lduMatrix.H:109
declareRunTimeSelectionTable(autoPtr, solver, asymMatrix,(const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces, const dictionary &solverControls),(fieldName, matrix, interfaceBouCoeffs, interfaceIntCoeffs, interfaces, solverControls))
label minIter_
Minimum number of iterations in the solver.
Definition: lduMatrix.H:120
const lduMatrix & matrix() const noexcept
Definition: lduMatrix.H:233
static autoPtr< solver > New(const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces, const dictionary &solverControls)
Return a new solver.
solveScalarField::cmptType normFactor(const solveScalarField &psi, const solveScalarField &source, const solveScalarField &Apsi, solveScalarField &tmpField) const
int log_
Level of verbosity in the solver output statements.
Definition: lduMatrix.H:117
virtual solverPerformance scalarSolve(solveScalarField &psi, const solveScalarField &source, const direction cmpt=0) const
Solve with given field and rhs (in solveScalar precision).
virtual ~solver()=default
Destructor.
scalar relTol_
Convergence tolerance relative to the initial.
Definition: lduMatrix.H:129
const FieldField< Field, scalar > & interfaceIntCoeffs_
Definition: lduMatrix.H:110
virtual void readControls()
Read the control parameters from the controlDict_.
virtual void read(const dictionary &)
Read and reset the solver parameters from the given stream.
static const label defaultMaxIter_
Default maximum number of iterations in the solver.
Definition: lduMatrix.H:105
const FieldField< Field, scalar > & interfaceBouCoeffs() const noexcept
Definition: lduMatrix.H:238
dictionary controlDict_
Dictionary of controls.
Definition: lduMatrix.H:114
virtual const word & type() const =0
Runtime type information.
const word & fieldName() const noexcept
Definition: lduMatrix.H:228
virtual solverPerformance solve(scalarField &psi, const scalarField &source, const direction cmpt=0) const =0
Solve with given field and rhs.
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition: lduMatrix.H:84
bool symmetric() const
Definition: lduMatrix.H:628
bool diagonal() const
Definition: lduMatrix.H:623
tmp< scalarField > H1() const
void operator=(const lduMatrix &)
scalarField & upper()
Definition: lduMatrix.C:203
tmp< Field< Type > > faceH(const Field< Type > &) const
const lduAddressing & lduAddr() const
Return the LDU addressing.
Definition: lduMatrix.H:578
friend Ostream & operator<<(Ostream &, const InfoProxy< lduMatrix > &)
ClassName("lduMatrix")
void residual(solveScalarField &rA, const solveScalarField &psi, const scalarField &source, const FieldField< Field, scalar > &interfaceBouCoeffs, const lduInterfaceFieldPtrsList &interfaces, const direction cmpt) const
bool hasUpper() const
Definition: lduMatrix.H:613
bool hasLower() const
Definition: lduMatrix.H:618
void Tmul(solveScalarField &, const tmp< solveScalarField > &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &, const direction cmpt) const
Matrix transpose multiplication with updated interfaces.
friend Ostream & operator<<(Ostream &, const lduMatrix &)
tmp< Field< Type > > faceH(const tmp< Field< Type > > &) const
void initMatrixInterfaces(const bool add, const FieldField< Field, scalar > &interfaceCoeffs, const lduInterfaceFieldPtrsList &interfaces, const solveScalarField &psiif, solveScalarField &result, const direction cmpt) const
void setLduMesh(const lduMesh &m)
Set the LDU mesh containing the addressing is obtained.
Definition: lduMatrix.H:572
const lduMesh & mesh() const
Return the LDU mesh from which the addressing is obtained.
Definition: lduMatrix.H:566
void sumA(solveScalarField &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &) const
Sum the coefficients on each row of the matrix.
void setResidualField(const scalarField &residual, const word &fieldName, const bool initial) const
Definition: lduMatrix.C:322
const lduSchedule & patchSchedule() const
Return the patch evaluation schedule.
Definition: lduMatrix.H:584
scalarField & lower()
Definition: lduMatrix.C:174
~lduMatrix()
Destructor.
Definition: lduMatrix.C:155
void operator*=(const scalarField &)
scalarField & diag()
Definition: lduMatrix.C:192
InfoProxy< lduMatrix > info() const
Return info proxy.
Definition: lduMatrix.H:748
tmp< Field< Type > > H(const Field< Type > &) const
void operator+=(const lduMatrix &)
void sumMagOffDiag(scalarField &sumOff) const
void Amul(solveScalarField &, const tmp< solveScalarField > &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &, const direction cmpt) const
Matrix multiplication with updated interfaces.
void operator-=(const lduMatrix &)
tmp< Field< Type > > H(const tmp< Field< Type > > &) const
bool hasDiag() const
Definition: lduMatrix.H:608
void updateMatrixInterfaces(const bool add, const FieldField< Field, scalar > &interfaceCoeffs, const lduInterfaceFieldPtrsList &interfaces, const solveScalarField &psiif, solveScalarField &result, const direction cmpt, const label startRequest) const
Update interfaced interfaces for matrix operations.
bool asymmetric() const
Definition: lduMatrix.H:633
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.
Triggers for starting/stopping code profiling.
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
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:517
lduMatrix member H operations.
Namespace for OpenFOAM.
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
uint8_t direction
Definition: direction.H:56
const direction noexcept
Definition: Scalar.H:223
Forward declarations of the specialisations of Field<T> for scalar, vector and tensor.
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)
CEqn solve()