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-2023 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 Class
28  Foam::LduMatrix
29 
30 Description
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 Note
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 
42 SourceFiles
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 "lduMatrix.H"
58 #include "Field.H"
59 #include "FieldField.H"
61 #include "SolverPerformance.H"
62 #include "typeInfo.H"
63 #include "autoPtr.H"
64 #include "runTimeSelectionTables.H"
65 
66 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67 
68 namespace Foam
69 {
70 
71 // Forward Declarations
72 
73 template<class Type, class DType, class LUType> class LduMatrix;
74 
75 template<class Type, class DType, class LUType>
76 Ostream& operator<<
77 (
78  Ostream&,
80 );
81 
82 
83 /*---------------------------------------------------------------------------*\
84  Class LduMatrix Declaration
85 \*---------------------------------------------------------------------------*/
86 
87 template<class Type, class DType, class LUType>
88 class 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 
111 public:
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 
124 
125  //- Dictionary of solution controls
127 
128  //- Verbosity level for solver output statements
129  int log_;
130 
131  //- Minimum number of iterations in the solver
132  label minIter_;
133 
134  //- Maximum number of iterations in the solver
135  label maxIter_;
136 
137  //- The matrix normalisation type
139 
140  //- Final convergence tolerance
141  Type tolerance_;
143  //- Convergence tolerance relative to the initial
144  Type relTol_;
145 
146 
147  // Protected Member Functions
148 
149  //- Read the control parameters from controlDict_
150  virtual void readControls();
151 
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  }
163 
164  public:
165 
166  //- Runtime type information
167  virtual const word& type() const = 0;
168 
169 
170  // Declare run-time constructor selection tables
171 
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  );
188 
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  //- Construct for given field name, matrix and controls
210  solver
211  (
212  const word& fieldName,
214  const dictionary& solverDict
215  );
216 
217 
218  // Selectors
219 
220  //- Return a new solver
221  static autoPtr<solver> New
222  (
223  const word& fieldName,
225  const dictionary& solverDict
226  );
227 
228 
229  //- Destructor
230  virtual ~solver() = default;
231 
232 
233  // Member Functions
234 
235  const word& fieldName() const noexcept
236  {
237  return fieldName_;
238  }
239 
241  {
242  return matrix_;
243  }
244 
245 
246  //- Read and reset the solver parameters from the given dictionary
247  virtual void read(const dictionary&);
248 
250  (
252  ) const = 0;
253 
254  //- Return the matrix norm using the specified norm method
255  Type normFactor
256  (
257  const Field<Type>& psi,
258  const Field<Type>& Apsi,
259  Field<Type>& tmpField,
260  const lduMatrix::normTypes normType
261  ) const;
262 
263  //- Return the matrix norm used to normalise the residual for the
264  //- stopping criterion
265  Type normFactor
266  (
267  const Field<Type>& psi,
268  const Field<Type>& Apsi,
269  Field<Type>& tmpField
270  ) const
271  {
272  return this->normFactor(psi, Apsi, tmpField, normType_);
273  }
274  };
275 
277  //- Abstract base-class for LduMatrix smoothers
278  class smoother
279  {
280  protected:
281 
282  // Protected Data
283 
286 
287 
288  public:
289 
290  //- Runtime type information
291  virtual const word& type() const = 0;
292 
293 
294  // Declare run-time constructor selection tables
295 
297  (
298  autoPtr,
299  smoother,
300  symMatrix,
301  (
302  const word& fieldName,
304  ),
305  (
306  fieldName,
307  matrix
308  )
309  );
310 
312  (
313  autoPtr,
314  smoother,
315  asymMatrix,
316  (
317  const word& fieldName,
319  ),
320  (
321  fieldName,
323  )
324  );
325 
326 
327  // Constructors
329  //- Construct for given field name and matrix
330  smoother
331  (
332  const word& fieldName,
334  );
335 
336 
337  // Selectors
338 
339  //- Return a new smoother
340  static autoPtr<smoother> New
341  (
342  const word& fieldName,
344  const dictionary& smootherDict
345  );
346 
347 
348  //- Destructor
349  virtual ~smoother() = default;
350 
351 
352  // Member Functions
353 
354  const word& fieldName() const noexcept
355  {
356  return fieldName_;
357  }
358 
360  {
361  return matrix_;
362  }
363 
364 
365  //- Smooth the solution for a given number of sweeps
366  virtual void smooth
367  (
368  Field<Type>& psi,
369  const label nSweeps
370  ) const = 0;
371  };
372 
373 
374  //- Abstract base-class for LduMatrix preconditioners
375  class preconditioner
376  {
377  protected:
378 
379  // Protected Data
380 
381  //- Reference to the base-solver this preconditioner is used with
382  const solver& solver_;
383 
384 
385  public:
386 
387  //- Runtime type information
388  virtual const word& type() const = 0;
389 
390 
391  // Declare run-time constructor selection tables
392 
394  (
395  autoPtr,
397  symMatrix,
398  (
399  const solver& sol,
400  const dictionary& preconditionerDict
401  ),
402  (sol, preconditionerDict)
403  );
404 
406  (
407  autoPtr,
409  asymMatrix,
410  (
411  const solver& sol,
412  const dictionary& preconditionerDict
413  ),
414  (sol, preconditionerDict)
415  );
416 
417 
418  // Constructors
419 
420  //- Construct for given solver
421  preconditioner(const solver& sol)
422  :
423  solver_(sol)
424  {}
425 
426 
427  // Selectors
428 
429  //- Return a new preconditioner
431  (
432  const solver& sol,
433  const dictionary& preconditionerDict
434  );
435 
436 
437  //- Destructor
438  virtual ~preconditioner() = default;
439 
441  // Member functions
442 
443  //- Read and reset the preconditioner parameters
444  //- from the given dictionary
445  virtual void read(const dictionary&)
446  {}
447 
448  //- Return wA the preconditioned form of residual rA
449  virtual void precondition
450  (
451  Field<Type>& wA,
452  const Field<Type>& rA
453  ) const = 0;
454 
455  //- Return wT the transpose-matrix preconditioned form of
456  //- residual rT.
457  // This is only required for preconditioning asymmetric matrices.
458  virtual void preconditionT
459  (
460  Field<Type>& wT,
461  const Field<Type>& rT
462  ) const
463  {
465  }
466  };
467 
468 
469  // Static Data
470 
471  // Declare name of the class and its debug switch
472  ClassName("LduMatrix");
473 
474 
475  // Constructors
476 
477  //- Construct given an LDU addressed mesh.
478  // The coefficients are initially empty for subsequent setting.
479  LduMatrix(const lduMesh&);
480 
481  //- Construct as copy
482  LduMatrix(const LduMatrix<Type, DType, LUType>&);
484  //- Construct as copy or re-use as specified.
486 
487  //- Construct given an LDU addressed mesh and an Istream
488  // from which the coefficients are read
489  LduMatrix(const lduMesh&, Istream&);
490 
491 
492  //- Destructor
493  ~LduMatrix();
494 
495 
496  // Member Functions
497 
498  // Access to addressing
499 
500  //- Return the LDU mesh from which the addressing is obtained
501  const lduMesh& mesh() const noexcept
502  {
503  return lduMesh_;
504  }
505 
506  //- Return the LDU addressing
507  const lduAddressing& lduAddr() const
508  {
509  return lduMesh_.lduAddr();
510  }
511 
512  //- Return the patch evaluation schedule
513  const lduSchedule& patchSchedule() const
514  {
515  return lduAddr().patchSchedule();
516  }
517 
518  //- Return interfaces
520  {
521  return interfaces_;
522  }
523 
524  //- Return interfaces
526  {
527  return interfaces_;
528  }
529 
530 
531  // Access to coefficients
533  Field<DType>& diag();
534  Field<LUType>& upper();
535  Field<LUType>& lower();
536  Field<Type>& source();
537 
539  {
540  return interfacesUpper_;
541  }
542 
544  {
545  return interfacesLower_;
546  }
547 
548 
549  const Field<DType>& diag() const;
550  const Field<LUType>& upper() const;
551  const Field<LUType>& lower() const;
552  const Field<Type>& source() const;
553 
555  {
556  return interfacesUpper_;
557  }
558 
559  const FieldField<Field, LUType>& interfacesLower() const
560  {
561  return interfacesLower_;
562  }
563 
564 
565  bool hasDiag() const noexcept
566  {
567  return (diagPtr_);
568  }
569 
570  bool hasUpper() const noexcept
571  {
572  return (upperPtr_);
573  }
574 
575  bool hasLower() const noexcept
576  {
577  return (lowerPtr_);
578  }
579 
580  bool hasSource() const noexcept
581  {
582  return (sourcePtr_);
583  }
584 
585  bool diagonal() const noexcept
586  {
587  return (diagPtr_ && !lowerPtr_ && !upperPtr_);
588  }
589 
590  bool symmetric() const noexcept
591  {
592  return (diagPtr_ && (!lowerPtr_ && upperPtr_));
593  }
594 
595  bool asymmetric() const noexcept
596  {
597  return (diagPtr_ && lowerPtr_ && upperPtr_);
598  }
599 
600 
601  // Operations
602 
603  void sumDiag();
604  void negSumDiag();
605 
606  void sumMagOffDiag(Field<LUType>& sumOff) const;
607 
608  //- Matrix multiplication
609  void Amul(Field<Type>&, const tmp<Field<Type>>&) const;
610 
611  //- Matrix transpose multiplication
612  void Tmul(Field<Type>&, const tmp<Field<Type>>&) const;
613 
614 
615  //- Sum the coefficients on each row of the matrix
616  void sumA(Field<Type>&) const;
617 
618 
619  void residual(Field<Type>& rA, const Field<Type>& psi) const;
621  tmp<Field<Type>> residual(const Field<Type>& psi) const;
622 
623 
624  //- Initialise the update of interfaced interfaces
625  // for matrix operations
627  (
628  const bool add,
629  const FieldField<Field, LUType>& interfaceCoeffs,
630  const Field<Type>& psiif,
631  Field<Type>& result
632  ) const;
634  //- Update interfaced interfaces for matrix operations
636  (
637  const bool add,
638  const FieldField<Field, LUType>& interfaceCoeffs,
639  const Field<Type>& psiif,
640  Field<Type>& result,
641  const label startRequest // starting request (for non-blocking)
642  ) const;
643 
644 
645  tmp<Field<Type>> H(const Field<Type>&) const;
646  tmp<Field<Type>> H(const tmp<Field<Type>>&) const;
647 
648  tmp<Field<Type>> faceH(const Field<Type>&) const;
650 
651 
652  // Member Operators
653 
655 
656  void negate();
657 
661  void operator*=(const scalarField&);
662  void operator*=(scalar);
663 
664 
665  // Ostream Operator
666 
667  friend Ostream& operator<< <Type, DType, LUType>
668  (
669  Ostream&,
671  );
672 };
673 
674 
675 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
676 
677 } // End namespace Foam
678 
679 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
681 #define makeLduMatrix(Type, DType, LUType) \
682  \
683 typedef Foam::LduMatrix<Type, DType, LUType> \
684  ldu##Type##DType##LUType##Matrix; \
685  \
686 defineNamedTemplateTypeNameAndDebug(ldu##Type##DType##LUType##Matrix, 0); \
687  \
688  \
689 typedef LduMatrix<Type, DType, LUType>::smoother \
690  ldu##Type##DType##LUType##Smoother; \
691  \
692 defineTemplateRunTimeSelectionTable \
693 ( \
694  ldu##Type##DType##LUType##Smoother, \
695  symMatrix \
696 ); \
697  \
698 defineTemplateRunTimeSelectionTable \
699 ( \
700  ldu##Type##DType##LUType##Smoother, \
701  asymMatrix \
702 ); \
703  \
704  \
705 typedef LduMatrix<Type, DType, LUType>::preconditioner \
706  ldu##Type##DType##LUType##Preconditioner; \
707  \
708 defineTemplateRunTimeSelectionTable \
709 ( \
710  ldu##Type##DType##LUType##Preconditioner, \
711  symMatrix \
712 ); \
713  \
714 defineTemplateRunTimeSelectionTable \
715 ( \
716  ldu##Type##DType##LUType##Preconditioner, \
717  asymMatrix \
718 ); \
719  \
720  \
721 typedef LduMatrix<Type, DType, LUType>::solver \
722  ldu##Type##DType##LUType##Solver; \
723  \
724 defineTemplateRunTimeSelectionTable \
725 ( \
726  ldu##Type##DType##LUType##Solver, \
727  symMatrix \
728 ); \
729  \
730 defineTemplateRunTimeSelectionTable \
731 ( \
732  ldu##Type##DType##LUType##Solver, \
733  asymMatrix \
734 );
735 
736 
737 #define makeLduPreconditioner(Precon, Type, DType, LUType) \
738  \
739 typedef Precon<Type, DType, LUType> \
740  Precon##Type##DType##LUType##Preconditioner; \
741 defineNamedTemplateTypeNameAndDebug \
742 ( \
743  Precon##Type##DType##LUType##Preconditioner, \
744  0 \
745 );
746 
747 #define makeLduSymPreconditioner(Precon, Type, DType, LUType) \
748  \
749 LduMatrix<Type, DType, LUType>::preconditioner:: \
750 addsymMatrixConstructorToTable<Precon##Type##DType##LUType##Preconditioner> \
751 add##Precon##Type##DType##LUType##PreconditionerSymMatrixConstructorToTable_;
752 
753 #define makeLduAsymPreconditioner(Precon, Type, DType, LUType) \
754  \
755 LduMatrix<Type, DType, LUType>::preconditioner:: \
756 addasymMatrixConstructorToTable<Precon##Type##DType##LUType##Preconditioner> \
757 add##Precon##Type##DType##LUType##PreconditionerAsymMatrixConstructorToTable_;
758 
759 
760 #define makeLduSmoother(Smoother, Type, DType, LUType) \
761  \
762 typedef Smoother<Type, DType, LUType> \
763  Smoother##Type##DType##LUType##Smoother; \
764  \
765 defineNamedTemplateTypeNameAndDebug \
766 ( \
767  Smoother##Type##DType##LUType##Smoother, \
768  0 \
769 );
770 
771 #define makeLduSymSmoother(Smoother, Type, DType, LUType) \
772  \
773 LduMatrix<Type, DType, LUType>::smoother:: \
774  addsymMatrixConstructorToTable<Smoother##Type##DType##LUType##Smoother> \
775  add##Smoother##Type##DType##LUType##SymMatrixConstructorToTable_;
776 
777 #define makeLduAsymSmoother(Smoother, Type, DType, LUType) \
778  \
779 LduMatrix<Type, DType, LUType>::smoother:: \
780  addasymMatrixConstructorToTable<Smoother##Type##DType##LUType##Smoother> \
781  add##Smoother##Type##DType##LUType##AsymMatrixConstructorToTable_;
782 
783 
784 #define makeLduSolver(Solver, Type, DType, LUType) \
785  \
786 typedef Solver<Type, DType, LUType> \
787  Solver##Type##DType##LUType##Solver; \
788  \
789 defineNamedTemplateTypeNameAndDebug \
790 ( \
791  Solver##Type##DType##LUType##Solver, \
792  0 \
793 );
794 
795 #define makeLduSymSolver(Solver, Type, DType, LUType) \
796  \
797 LduMatrix<Type, DType, LUType>::solver:: \
798  addsymMatrixConstructorToTable<Solver##Type##DType##LUType##Solver> \
799  add##Solver##Type##DType##LUType##SymMatrixConstructorToTable_;
800 
801 #define makeLduAsymSolver(Solver, Type, DType, LUType) \
802  \
803 LduMatrix<Type, DType, LUType>::solver:: \
804  addasymMatrixConstructorToTable<Solver##Type##DType##LUType##Solver> \
805  add##Solver##Type##DType##LUType##AsymMatrixConstructorToTable_;
806 
807 
808 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
809 
810 #ifdef NoRepository
811  #include "LduMatrix.C"
812 #endif
813 
814 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
815 
816 #endif
817 
818 // ************************************************************************* //
smoother(const word &fieldName, const LduMatrix< Type, DType, LUType > &matrix)
Construct for given field name and matrix.
void initMatrixInterfaces(const bool add, const FieldField< Field, LUType > &interfaceCoeffs, const Field< Type > &psiif, Field< Type > &result) const
Initialise the update of interfaced interfaces.
virtual const word & type() const =0
Runtime type information.
void readControl(const dictionary &dict, T &val, const word &key)
Deprecated(2021-09) Read control parameter from dictionary.
Definition: LduMatrix.H:186
dictionary dict
virtual void preconditionT(Field< Type > &wT, const Field< Type > &rT) const
Return wT the transpose-matrix preconditioned form of residual rT.
Definition: LduMatrix.H:532
virtual const word & type() const =0
Runtime type information.
static autoPtr< preconditioner > New(const solver &sol, const dictionary &preconditionerDict)
Return a new preconditioner.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Field< Type > & source()
Definition: LduMatrix.C:243
FieldField< Field, LUType > & interfacesLower()
Definition: LduMatrix.H:638
const LduMatrix< Type, DType, LUType > & matrix() const noexcept
Definition: LduMatrix.H:411
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
LduMatrix(const lduMesh &)
Construct given an LDU addressed mesh.
Definition: LduMatrix.C:27
Type normFactor(const Field< Type > &psi, const Field< Type > &Apsi, Field< Type > &tmpField, const lduMatrix::normTypes normType) const
Return the matrix norm using the specified norm method.
void sumMagOffDiag(Field< LUType > &sumOff) const
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:53
static autoPtr< smoother > New(const word &fieldName, const LduMatrix< Type, DType, LUType > &matrix, const dictionary &smootherDict)
Return a new smoother.
virtual void precondition(Field< Type > &wA, const Field< Type > &rA) const =0
Return wA the preconditioned form of residual rA.
Field< LUType > & lower()
Definition: LduMatrix.C:220
label minIter_
Minimum number of iterations in the solver.
Definition: LduMatrix.H:147
const LduMatrix< Type, DType, LUType > & matrix() const noexcept
Definition: LduMatrix.H:276
void sumA(Field< Type > &) const
Sum the coefficients on each row of the matrix.
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:51
dictionary controlDict_
Dictionary of solution controls.
Definition: LduMatrix.H:137
virtual const word & type() const =0
Runtime type information.
ClassName("LduMatrix")
virtual void smooth(Field< Type > &psi, const label nSweeps) const =0
Smooth the solution for a given number of sweeps.
void Amul(Field< Type > &, const tmp< Field< Type >> &) const
Matrix multiplication.
int log_
Verbosity level for solver output statements.
Definition: LduMatrix.H:142
Abstract base-class for LduMatrix smoothers.
Definition: LduMatrix.H:322
A class for handling words, derived from Foam::string.
Definition: word.H:63
lduMatrix::normTypes normType_
The matrix normalisation type.
Definition: LduMatrix.H:157
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics...
Abstract base-class for LduMatrix preconditioners.
Definition: LduMatrix.H:431
bool hasSource() const noexcept
Definition: LduMatrix.H:675
tmp< Field< Type > > H(const Field< Type > &) const
bool asymmetric() const noexcept
Definition: LduMatrix.H:690
~LduMatrix()
Destructor.
Definition: LduMatrix.C:158
const lduSchedule & patchSchedule() const
Return the patch evaluation schedule.
Definition: LduMatrix.H:604
virtual void read(const dictionary &)
Read and reset the preconditioner parameters from the given dictionary.
Definition: LduMatrix.H:513
void operator*=(const scalarField &)
Field< LUType > & upper()
Definition: LduMatrix.C:197
void updateMatrixInterfaces(const bool add, const FieldField< Field, LUType > &interfaceCoeffs, const Field< Type > &psiif, Field< Type > &result, const label startRequest) const
Update interfaced interfaces for matrix operations.
normTypes
Enumerated matrix normalisation types.
Definition: lduMatrix.H:103
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: HashTable.H:106
virtual ~smoother()=default
Destructor.
void operator-=(const LduMatrix< Type, DType, LUType > &)
Abstract base-class for LduMatrix solvers.
Definition: LduMatrix.H:125
const word & fieldName() const noexcept
Definition: LduMatrix.H:406
void residual(Field< Type > &rA, const Field< Type > &psi) const
preconditioner(const solver &sol)
Construct for given solver.
Definition: LduMatrix.H:483
const LduMatrix< Type, DType, LUType > & matrix_
Definition: LduMatrix.H:329
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
const direction noexcept
Definition: Scalar.H:258
solver(const word &fieldName, const LduMatrix< Type, DType, LUType > &matrix, const dictionary &solverDict)
Construct for given field name, matrix and controls.
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
const word & fieldName() const noexcept
Definition: LduMatrix.H:271
void Tmul(Field< Type > &, const tmp< Field< Type >> &) const
Matrix transpose multiplication.
const lduAddressing & lduAddr() const
Return the LDU addressing.
Definition: LduMatrix.H:596
const solver & solver_
Reference to the base-solver this preconditioner is used with.
Definition: LduMatrix.H:440
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
bool hasUpper() const noexcept
Definition: LduMatrix.H:665
declareRunTimeSelectionTable(autoPtr, preconditioner, symMatrix,(const solver &sol, const dictionary &preconditionerDict),(sol, preconditionerDict))
virtual ~solver()=default
Destructor.
const LduInterfaceFieldPtrsList< Type > & interfaces() const
Return interfaces.
Definition: LduMatrix.H:612
Field< DType > & diag()
Definition: LduMatrix.C:185
declareRunTimeSelectionTable(autoPtr, smoother, symMatrix,(const word &fieldName, const LduMatrix< Type, DType, LUType > &matrix),(fieldName, matrix))
virtual void read(const dictionary &)
Read and reset the solver parameters from the given dictionary.
Type relTol_
Convergence tolerance relative to the initial.
Definition: LduMatrix.H:167
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
Definition: foamGltfBase.H:103
bool diagonal() const noexcept
Definition: LduMatrix.H:680
FieldField< Field, LUType > & interfacesUpper()
Definition: LduMatrix.H:633
LduMatrix is a general matrix class in which the coefficients are stored as three arrays...
Definition: LduMatrix.H:68
label maxIter_
Maximum number of iterations in the solver.
Definition: LduMatrix.H:152
bool hasLower() const noexcept
Definition: LduMatrix.H:670
void operator+=(const LduMatrix< Type, DType, LUType > &)
tmp< Field< Type > > faceH(const Field< Type > &) const
static autoPtr< solver > New(const word &fieldName, const LduMatrix< Type, DType, LUType > &matrix, const dictionary &solverDict)
Return a new solver.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
const lduMesh & mesh() const noexcept
Return the LDU mesh from which the addressing is obtained.
Definition: LduMatrix.H:588
const volScalarField & psi
The class contains the addressing required by the lduMatrix: upper, lower and losort.
virtual SolverPerformance< Type > solve(Field< Type > &psi) const =0
const LduMatrix< Type, DType, LUType > & matrix_
Definition: LduMatrix.H:132
Macros to ease declaration of run-time selection tables.
virtual void readControls()
Read the control parameters from controlDict_.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
virtual ~preconditioner()=default
Destructor.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:686
declareRunTimeSelectionTable(autoPtr, solver, symMatrix,(const word &fieldName, const LduMatrix< Type, DType, LUType > &matrix, const dictionary &solverDict),(fieldName, matrix, solverDict))
void operator=(const LduMatrix< Type, DType, LUType > &)
List of coupled interface fields to be used in coupling.
bool hasDiag() const noexcept
Definition: LduMatrix.H:660
Type tolerance_
Final convergence tolerance.
Definition: LduMatrix.H:162
bool symmetric() const noexcept
Definition: LduMatrix.H:685
Namespace for OpenFOAM.
virtual const lduSchedule & patchSchedule() const =0
Return patch field evaluation schedule.