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-2024 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  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 
41 SourceFiles
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 Foam_lduMatrix_H
54 #define Foam_lduMatrix_H
55 
56 #include "lduMesh.H"
57 #include "primitiveFieldsFwd.H"
58 #include "FieldField.H"
60 #include "solverPerformance.H"
61 #include "typeInfo.H"
62 #include "autoPtr.H"
63 #include "runTimeSelectionTables.H"
64 #include "InfoProxy.H"
65 #include "Enum.H"
66 #include "profilingTrigger.H"
67 #include <functional> // For reference_wrapper
68 
69 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
70 
71 namespace Foam
72 {
73 
74 // Forward Declarations
75 class lduMatrix;
76 
77 Ostream& operator<<(Ostream&, const lduMatrix&);
78 Ostream& operator<<(Ostream&, const InfoProxy<lduMatrix>&);
79 
80 
81 /*---------------------------------------------------------------------------*\
82  Class lduMatrix Declaration
83 \*---------------------------------------------------------------------------*/
84 
85 class lduMatrix
86 {
87  // Private Data
88 
89  //- LDU mesh reference
90  //const lduMesh& lduMesh_;
91  std::reference_wrapper<const lduMesh> lduMesh_;
92 
93  //- Diagonal coefficients
94  std::unique_ptr<scalarField> diagPtr_;
95 
96  //- Off-diagonal coefficients (not including interfaces)
97  std::unique_ptr<scalarField> lowerPtr_;
98 
99  //- Off-diagonal coefficients (not including interfaces)
100  std::unique_ptr<scalarField> upperPtr_;
101 
102 
103 public:
104 
105  // Public Types
106 
107  //- Enumerated matrix normalisation types
108  enum class normTypes : char
109  {
110  NO_NORM,
111  DEFAULT_NORM,
113  };
115  //- Names for the normTypes
116  static const Enum<normTypes> normTypesNames_;
117 
118  //- Default maximum number of iterations for solvers (1000)
119  static constexpr const label defaultMaxIter = 1000;
120 
121  //- Default (absolute) tolerance (1e-6)
122  static const scalar defaultTolerance;
123 
125  // -----------------------------------------------------------------------
126  //- Abstract base-class for lduMatrix solvers
127  class solver
128  {
129  protected:
130 
131  // Protected Data
132 
138 
139  //- Dictionary of solution controls
142  //- Verbosity level for solver output statements
143  int log_;
144 
145  //- Minimum number of iterations in the solver
146  label minIter_;
148  //- Maximum number of iterations in the solver
149  label maxIter_;
151  //- The normalisation type
153 
154  //- Final convergence tolerance
155  scalar tolerance_;
157  //- Convergence tolerance relative to the initial
158  scalar relTol_;
159 
160  //- Profiling instrumentation
162 
163 
164  // Protected Member Functions
165 
166  //- Read the control parameters from controlDict_
167  virtual void readControls();
168 
169 
170  public:
172  //- Runtime type information
173  virtual const word& type() const = 0;
174 
175 
176  // Declare run-time constructor selection tables
177 
179  (
180  autoPtr,
182  symMatrix,
183  (
184  const word& fieldName,
185  const lduMatrix& matrix,
189  const dictionary& solverControls
190  ),
191  (
192  fieldName,
193  matrix,
196  interfaces,
197  solverControls
198  )
199  );
200 
202  (
203  autoPtr,
204  solver,
205  asymMatrix,
206  (
207  const word& fieldName,
208  const lduMatrix& matrix,
212  const dictionary& solverControls
213  ),
214  (
215  fieldName,
216  matrix,
219  interfaces,
220  solverControls
221  )
222  );
223 
224 
225  // Constructors
226 
227  //- Construct solver for given field name, matrix etc
228  solver
229  (
230  const word& fieldName,
231  const lduMatrix& matrix,
235  const dictionary& solverControls
236  );
237 
238  // Selectors
239 
240  //- Return a new solver of given type
241  static autoPtr<solver> New
242  (
243  const word& solverName,
244  const word& fieldName,
245  const lduMatrix& matrix,
249  const dictionary& solverControls
250  );
251 
252  //- Return a new solver given dictionary
253  static autoPtr<solver> New
254  (
255  const word& fieldName,
256  const lduMatrix& matrix,
260  const dictionary& solverControls
261  );
262 
263 
264 
265  //- Destructor
266  virtual ~solver() = default;
267 
268 
269  // Member Functions
270 
271  const word& fieldName() const noexcept
272  {
273  return fieldName_;
274  }
275 
276  const lduMatrix& matrix() const noexcept
277  {
278  return matrix_;
279  }
280 
282  {
283  return interfaceBouCoeffs_;
284  }
285 
286  const FieldField<Field, scalar>& interfaceIntCoeffs() const noexcept
287  {
288  return interfaceIntCoeffs_;
289  }
290 
292  {
293  return interfaces_;
294  }
295 
296 
297  //- Read and reset the solver parameters from the given stream
298  virtual void read(const dictionary&);
299 
300  //- Solve with given field and rhs
301  virtual solverPerformance solve
302  (
303  scalarField& psi,
304  const scalarField& source,
305  const direction cmpt=0
306  ) const = 0;
307 
308  //- Solve with given field and rhs (in solveScalar precision).
309  // Default is to call solve routine
311  (
313  const solveScalarField& source,
314  const direction cmpt=0
315  ) const;
316 
317  //- Return the matrix norm using the specified norm method
319  (
320  const solveScalarField& psi,
321  const solveScalarField& source,
322  const solveScalarField& Apsi,
323  solveScalarField& tmpField,
324  const lduMatrix::normTypes normType
325  ) const;
326 
327  //- Return the matrix norm used to normalise the residual for the
328  //- stopping criterion
330  (
331  const solveScalarField& psi,
332  const solveScalarField& source,
333  const solveScalarField& Apsi,
334  solveScalarField& tmpField
335  ) const
336  {
337  return this->normFactor(psi, source, Apsi, tmpField, normType_);
338  }
339  };
340 
341 
342  // -----------------------------------------------------------------------
343  //- Abstract base-class for lduMatrix smoothers
344  class smoother
345  {
346  protected:
347 
348  // Protected Data
349 
351  const lduMatrix& matrix_;
355 
356 
357  public:
358 
359  //- Find the smoother name (directly or from a sub-dictionary)
360  static word getName(const dictionary&);
361 
362  //- Runtime type information
363  virtual const word& type() const = 0;
364 
365 
366  // Declare run-time constructor selection tables
367 
369  (
370  autoPtr,
371  smoother,
372  symMatrix,
373  (
374  const word& fieldName,
375  const lduMatrix& matrix,
379  ),
380  (
381  fieldName,
382  matrix,
385  interfaces
386  )
387  );
388 
390  (
391  autoPtr,
392  smoother,
393  asymMatrix,
394  (
395  const word& fieldName,
396  const lduMatrix& matrix,
400  ),
401  (
402  fieldName,
403  matrix,
407  )
408  );
410 
411  // Constructors
412 
413  //- Construct for given field name, matrix etc
414  smoother
415  (
416  const word& fieldName,
417  const lduMatrix& matrix,
421  );
422 
423 
424  // Selectors
425 
426  //- Return a new smoother
427  static autoPtr<smoother> New
428  (
429  const word& fieldName,
430  const lduMatrix& matrix,
434  const dictionary& solverControls
435  );
436 
437 
438  //- Destructor
439  virtual ~smoother() = default;
440 
441 
442  // Member Functions
443 
444  const word& fieldName() const noexcept
445  {
446  return fieldName_;
447  }
448 
449  const lduMatrix& matrix() const noexcept
450  {
451  return matrix_;
452  }
453 
455  {
456  return interfaceBouCoeffs_;
457  }
458 
459  const FieldField<Field, scalar>& interfaceIntCoeffs() const noexcept
460  {
461  return interfaceIntCoeffs_;
462  }
463 
465  {
466  return interfaces_;
467  }
468 
469 
470  //- Smooth the solution for a given number of sweeps
471  virtual void smooth
472  (
474  const scalarField& source,
475  const direction cmpt,
476  const label nSweeps
477  ) const = 0;
478 
479  //- Smooth the solution for a given number of sweeps
480  virtual void scalarSmooth
481  (
483  const solveScalarField& source,
484  const direction cmpt,
485  const label nSweeps
486  ) const = 0;
487  };
488 
489 
490  // -----------------------------------------------------------------------
491  //- Abstract base-class for lduMatrix preconditioners
492  class preconditioner
493  {
494  protected:
495 
496  // Protected Data
497 
498  //- Reference to the base-solver this preconditioner is used with
499  const solver& solver_;
500 
501 
502  public:
503 
504  //- Find the preconditioner name (directly or from a sub-dictionary)
505  static word getName(const dictionary&);
506 
507  //- Runtime type information
508  virtual const word& type() const = 0;
510 
511  // Declare run-time constructor selection tables
512 
514  (
515  autoPtr,
517  symMatrix,
518  (
519  const solver& sol,
520  const dictionary& solverControls
521  ),
522  (sol, solverControls)
523  );
526  (
527  autoPtr,
529  asymMatrix,
530  (
531  const solver& sol,
532  const dictionary& solverControls
533  ),
534  (sol, solverControls)
535  );
536 
537 
538  // Constructors
539 
540  //- Construct for given solver
541  explicit preconditioner(const solver& sol)
542  :
543  solver_(sol)
544  {}
545 
546 
547  // Selectors
548 
549  //- Return a new preconditioner
551  (
552  const solver& sol,
553  const dictionary& solverControls
554  );
555 
556 
557  //- Destructor
558  virtual ~preconditioner() = default;
559 
560 
561  // Member Functions
562 
563  //- Read and reset the preconditioner parameters
564  //- from the given stream
565  virtual void read(const dictionary&)
566  {}
567 
568  //- Return wA the preconditioned form of residual rA
569  virtual void precondition
570  (
571  solveScalarField& wA,
572  const solveScalarField& rA,
573  const direction cmpt=0
574  ) const = 0;
575 
576  //- Return wT the transpose-matrix preconditioned form of
577  //- residual rT.
578  // This is only required for preconditioning asymmetric matrices.
579  virtual void preconditionT
580  (
581  solveScalarField& wT,
582  const solveScalarField& rT,
583  const direction cmpt=0
584  ) const
585  {
587  }
588 
589  //- Signal end of solver
590  virtual void setFinished(const solverPerformance& perf) const
591  {}
592  };
593 
594 
595  // -----------------------------------------------------------------------
596 
597  // Static Data
598 
599  // Declare name of the class and its debug switch
600  ClassName("lduMatrix");
601 
602 
603  // Constructors
604 
605  //- Construct (without coefficients) for an LDU addressed mesh.
606  // Not yet 'explicit' (legacy code may rely on implicit construct)
607  lduMatrix(const lduMesh& mesh);
608 
609  //- Copy construct
610  lduMatrix(const lduMatrix&);
611 
612  //- Move construct
613  lduMatrix(lduMatrix&&);
614 
615  //- Construct as copy or re-use as specified.
616  lduMatrix(lduMatrix&, bool reuse);
617 
618  //- Construct given an LDU addressed mesh and an Istream
619  //- from which the coefficients are read
621 
622 
623  //- Destructor
624  ~lduMatrix() = default;
625 
626 
627  // Member Functions
628 
629  // Addressing
630 
631  //- Return the LDU mesh from which the addressing is obtained
632  const lduMesh& mesh() const noexcept
633  {
634  return lduMesh_;
635  }
636 
637  //- Set the LDU mesh containing the addressing
638  void setLduMesh(const lduMesh& m)
639  {
640  lduMesh_ = m;
641  }
642 
643  //- Return the LDU addressing
644  const lduAddressing& lduAddr() const
645  {
646  return mesh().lduAddr();
647  }
648 
649  //- Return the patch evaluation schedule
650  const lduSchedule& patchSchedule() const
651  {
652  return mesh().lduAddr().patchSchedule();
653  }
654 
655 
656  // Coefficients
657 
658  const scalarField& diag() const;
659  const scalarField& upper() const;
660  const scalarField& lower() const;
661 
662  scalarField& diag();
663  scalarField& upper();
664  scalarField& lower();
665 
666  // Size with externally provided sizes
667  // (for constructing with 'fake' mesh in GAMG)
668 
669  scalarField& diag(label size);
670  scalarField& upper(label nCoeffs);
671  scalarField& lower(label nCoeffs);
672 
673 
674  // Characteristics
675 
676  //- The matrix type (empty, diagonal, symmetric, ...)
677  word matrixTypeName() const;
678 
679  bool hasDiag() const noexcept { return bool(diagPtr_); }
680  bool hasUpper() const noexcept { return bool(upperPtr_); }
681  bool hasLower() const noexcept { return bool(lowerPtr_); }
683  //- Matrix has diagonal only
684  bool diagonal() const noexcept
685  {
686  return (diagPtr_ && !lowerPtr_ && !upperPtr_);
687  }
688 
689  //- Matrix is symmetric
690  bool symmetric() const noexcept
691  {
692  return (diagPtr_ && !lowerPtr_ && upperPtr_);
693  }
694 
695  //- Matrix is asymmetric (ie, full)
696  bool asymmetric() const noexcept
697  {
698  return (diagPtr_ && lowerPtr_ && upperPtr_);
699  }
700 
701 
702  // Operations
703 
704  void sumDiag();
705  void negSumDiag();
706 
707  void sumMagOffDiag(scalarField& sumOff) const;
708 
709  //- Matrix multiplication with updated interfaces.
710  void Amul
711  (
713  const tmp<solveScalarField>&,
714  const FieldField<Field, scalar>&,
716  const direction cmpt
717  ) const;
718 
719  //- Matrix transpose multiplication with updated interfaces.
720  void Tmul
721  (
723  const tmp<solveScalarField>&,
724  const FieldField<Field, scalar>&,
726  const direction cmpt
727  ) const;
728 
729 
730  //- Sum the coefficients on each row of the matrix
731  void sumA
732  (
734  const FieldField<Field, scalar>&,
736  ) const;
737 
738 
739  void residual
740  (
741  solveScalarField& rA,
742  const solveScalarField& psi,
743  const scalarField& source,
744  const FieldField<Field, scalar>& interfaceBouCoeffs,
745  const lduInterfaceFieldPtrsList& interfaces,
746  const direction cmpt
747  ) const;
748 
750  (
751  const solveScalarField& psi,
752  const scalarField& source,
753  const FieldField<Field, scalar>& interfaceBouCoeffs,
754  const lduInterfaceFieldPtrsList& interfaces,
755  const direction cmpt
756  ) const;
757 
758 
759  //- Initialise the update of interfaced interfaces
760  //- for matrix operations
762  (
763  const bool add,
764  const FieldField<Field, scalar>& interfaceCoeffs,
765  const lduInterfaceFieldPtrsList& interfaces,
766  const solveScalarField& psiif,
767  solveScalarField& result,
768  const direction cmpt
769  ) const;
770 
771  //- Update interfaced interfaces for matrix operations
773  (
774  const bool add,
775  const FieldField<Field, scalar>& interfaceCoeffs,
776  const lduInterfaceFieldPtrsList& interfaces,
777  const solveScalarField& psiif,
778  solveScalarField& result,
779  const direction cmpt,
780  const label startRequest // starting request (for non-blocking)
781  ) const;
782 
783  //- Set the residual field using an IOField on the object registry
784  //- if it exists
785  void setResidualField
786  (
787  const scalarField& residual,
788  const word& fieldName,
789  const bool initial
790  ) const;
791 
792  template<class Type>
793  tmp<Field<Type>> H(const Field<Type>&) const;
795  template<class Type>
796  tmp<Field<Type>> H(const tmp<Field<Type>>&) const;
797 
798  tmp<scalarField> H1() const;
799 
800  template<class Type>
802 
803  template<class Type>
804  tmp<Field<Type>> faceH(const tmp<Field<Type>>&) const;
805 
806 
807  // Info
808 
809  //- Return info proxy,
810  //- used to print matrix information to a stream
811  InfoProxy<lduMatrix> info() const noexcept { return *this; }
812 
813 
814  // Member Operators
815 
816  //- Copy assignment
817  void operator=(const lduMatrix&);
818 
819  //- Move assignment
820  void operator=(lduMatrix&&);
821 
822  void negate();
823 
824  void operator+=(const lduMatrix&);
825  void operator-=(const lduMatrix&);
826 
827  void operator*=(const scalarField&);
828  void operator*=(scalar);
829 
830 
831  // Ostream Operators
832 
833  friend Ostream& operator<<(Ostream&, const lduMatrix&);
834  friend Ostream& operator<<(Ostream&, const InfoProxy<lduMatrix>&);
835 };
836 
837 
838 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
839 
840 } // End namespace Foam
841 
842 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
843 
844 #ifdef NoRepository
845  #include "lduMatrixTemplates.C"
846 #endif
847 
848 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
849 
850 #endif
851 
852 // ************************************************************************* //
UPtrList< const lduInterfaceField > lduInterfaceFieldPtrsList
List of coupled interface fields to be used in coupling.
virtual void read(const dictionary &)
Read and reset the preconditioner parameters from the given stream.
Definition: lduMatrix.H:650
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.
const lduInterfaceFieldPtrsList & interfaces_
Definition: lduMatrix.H:409
const FieldField< Field, scalar > & interfaceBouCoeffs() const noexcept
Definition: lduMatrix.H:323
const scalarField & diag() const
Definition: lduMatrix.C:163
tmp< Field< Type > > faceH(const Field< Type > &) const
const lduInterfaceFieldPtrsList & interfaces() const noexcept
Definition: lduMatrix.H:529
void Tmul(solveScalarField &, const tmp< solveScalarField > &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &, const direction cmpt) const
Matrix transpose multiplication with updated interfaces.
bool symmetric() const noexcept
Matrix is symmetric.
Definition: lduMatrix.H:809
word matrixTypeName() const
The matrix type (empty, diagonal, symmetric, ...)
Definition: lduMatrix.C:146
virtual void read(const dictionary &)
Read and reset the solver parameters from the given stream.
const FieldField< Field, scalar > & interfaceIntCoeffs() const noexcept
Definition: lduMatrix.H:524
uint8_t direction
Definition: direction.H:46
const FieldField< Field, scalar > & interfaceIntCoeffs_
Definition: lduMatrix.H:150
const lduMatrix & matrix() const noexcept
Definition: lduMatrix.H:514
preconditioner(const solver &sol)
Construct for given solver.
Definition: lduMatrix.H:620
~lduMatrix()=default
Destructor.
Field< solveScalar > solveScalarField
virtual const word & type() const =0
Runtime type information.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
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:572
void residual(solveScalarField &rA, const solveScalarField &psi, const scalarField &source, const FieldField< Field, scalar > &interfaceBouCoeffs, const lduInterfaceFieldPtrsList &interfaces, const direction cmpt) const
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:335
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))
virtual solverPerformance solve(scalarField &psi, const scalarField &source, const direction cmpt=0) const =0
Solve with given field and rhs.
const FieldField< Field, scalar > & interfaceBouCoeffs_
Definition: lduMatrix.H:149
virtual void preconditionT(solveScalarField &wT, const solveScalarField &rT, const direction cmpt=0) const
Return wT the transpose-matrix preconditioned form of residual rT.
Definition: lduMatrix.H:670
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
pTraits< solveScalar >::cmptType cmptType
Component type.
Definition: Field.H:123
const lduMatrix & matrix() const noexcept
Definition: lduMatrix.H:318
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:53
InfoProxy< lduMatrix > info() const noexcept
Return info proxy, used to print matrix information to a stream.
Definition: lduMatrix.H:946
Triggers for starting/stopping code profiling.
ClassName("lduMatrix")
tmp< scalarField > H1() const
void sumMagOffDiag(scalarField &sumOff) const
static constexpr const label defaultMaxIter
Default maximum number of iterations for solvers (1000)
Definition: lduMatrix.H:129
int log_
Verbosity level for solver output statements.
Definition: lduMatrix.H:161
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:51
static const Enum< normTypes > normTypesNames_
Names for the normTypes.
Definition: lduMatrix.H:124
label minIter_
Minimum number of iterations in the solver.
Definition: lduMatrix.H:166
lduMatrix::normTypes normType_
The normalisation type.
Definition: lduMatrix.H:176
virtual ~preconditioner()=default
Destructor.
const lduMatrix & matrix_
Definition: lduMatrix.H:406
static word getName(const dictionary &)
Find the smoother name (directly or from a sub-dictionary)
const lduSchedule & patchSchedule() const
Return the patch evaluation schedule.
Definition: lduMatrix.H:763
Abstract base-class for lduMatrix smoothers.
Definition: lduMatrix.H:399
Forward declarations of the specialisations of Field<T> for scalar, vector and tensor.
scalar tolerance_
Final convergence tolerance.
Definition: lduMatrix.H:181
void Amul(solveScalarField &, const tmp< solveScalarField > &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &, const direction cmpt) const
Matrix multiplication with updated interfaces.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics...
label maxIter_
Maximum number of iterations in the solver.
Definition: lduMatrix.H:171
static autoPtr< preconditioner > New(const solver &sol, const dictionary &solverControls)
Return a new preconditioner.
Abstract base-class for lduMatrix solvers.
Definition: lduMatrix.H:141
const word & fieldName() const noexcept
Definition: lduMatrix.H:509
lduMatrix member H operations.
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.
void initMatrixInterfaces(const bool add, const FieldField< Field, scalar > &interfaceCoeffs, const lduInterfaceFieldPtrsList &interfaces, const solveScalarField &psiif, solveScalarField &result, const direction cmpt) const
Initialise the update of interfaced interfaces for matrix operations.
const scalarField & lower() const
Definition: lduMatrix.C:268
friend Ostream & operator<<(Ostream &, const lduMatrix &)
virtual const lduAddressing & lduAddr() const =0
Return ldu addressing.
void sumA(solveScalarField &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &) const
Sum the coefficients on each row of the matrix.
const lduInterfaceFieldPtrsList & interfaces() const noexcept
Definition: lduMatrix.H:333
const FieldField< Field, scalar > & interfaceBouCoeffs() const noexcept
Definition: lduMatrix.H:519
static const scalar defaultTolerance
Default (absolute) tolerance (1e-6)
Definition: lduMatrix.H:134
solver(const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces, const dictionary &solverControls)
Construct solver for given field name, matrix etc.
virtual const word & type() const =0
Runtime type information.
const FieldField< Field, scalar > & interfaceBouCoeffs_
Definition: lduMatrix.H:407
normTypes
Enumerated matrix normalisation types.
Definition: lduMatrix.H:114
void operator=(const lduMatrix &)
Copy assignment.
const lduMatrix & matrix_
Definition: lduMatrix.H:148
bool hasUpper() const noexcept
Definition: lduMatrix.H:795
bool asymmetric() const noexcept
Matrix is asymmetric (ie, full)
Definition: lduMatrix.H:817
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
smoother(const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces)
Construct for given field name, matrix etc.
const direction noexcept
Definition: Scalar.H:258
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
void operator*=(const scalarField &)
static autoPtr< solver > New(const word &solverName, 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 of given type.
lduInterfaceFieldPtrsList interfaces_
Definition: lduMatrix.H:151
solveScalarField::cmptType normFactor(const solveScalarField &psi, const solveScalarField &source, const solveScalarField &Apsi, solveScalarField &tmpField, const lduMatrix::normTypes normType) const
Return the matrix norm using the specified norm method.
bool diagonal() const noexcept
Matrix has diagonal only.
Definition: lduMatrix.H:801
virtual ~smoother()=default
Destructor.
virtual const word & type() const =0
Runtime type information.
dictionary controlDict_
Dictionary of solution controls.
Definition: lduMatrix.H:156
"default" norm (== L1_scaled)
Basic run-time type information using word as the type&#39;s name. Used to enhance the standard RTTI to c...
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:77
lduMatrix is a general matrix class in which the coefficients are stored as three arrays...
Definition: lduMatrix.H:80
"none" norm (returns 1)
const scalarField & upper() const
Definition: lduMatrix.C:203
Abstract base-class for lduMatrix preconditioners.
Definition: lduMatrix.H:563
A helper class for outputting values to Ostream.
Definition: ensightCells.H:44
const lduMesh & mesh() const noexcept
Return the LDU mesh from which the addressing is obtained.
Definition: lduMatrix.H:739
const lduAddressing & lduAddr() const
Return the LDU addressing.
Definition: lduMatrix.H:755
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.
virtual void precondition(solveScalarField &wA, const solveScalarField &rA, const direction cmpt=0) const =0
Return wA the preconditioned form of residual rA.
void setLduMesh(const lduMesh &m)
Set the LDU mesh containing the addressing.
Definition: lduMatrix.H:747
virtual ~solver()=default
Destructor.
const FieldField< Field, scalar > & interfaceIntCoeffs() const noexcept
Definition: lduMatrix.H:328
const FieldField< Field, scalar > & interfaceIntCoeffs_
Definition: lduMatrix.H:408
lduMatrix(const lduMesh &mesh)
Construct (without coefficients) for an LDU addressed mesh.
Definition: lduMatrix.C:54
virtual solverPerformance scalarSolve(solveScalarField &psi, const solveScalarField &source, const direction cmpt=0) const
Solve with given field and rhs (in solveScalar precision).
static word getName(const dictionary &)
Find the preconditioner name (directly or from a sub-dictionary)
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
bool hasDiag() const noexcept
Definition: lduMatrix.H:794
const volScalarField & psi
The class contains the addressing required by the lduMatrix: upper, lower and losort.
tmp< Field< Type > > H(const Field< Type > &) const
Macros to ease declaration of run-time selection tables.
void operator-=(const lduMatrix &)
void operator+=(const lduMatrix &)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
profilingTrigger profiling_
Profiling instrumentation.
Definition: lduMatrix.H:191
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:696
bool hasLower() const noexcept
Definition: lduMatrix.H:796
const word & fieldName() const noexcept
Definition: lduMatrix.H:313
virtual void setFinished(const solverPerformance &perf) const
Signal end of solver.
Definition: lduMatrix.H:682
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.
Namespace for OpenFOAM.
SolverPerformance< scalar > solverPerformance
SolverPerformance instantiated for a scalar.
virtual void readControls()
Read the control parameters from controlDict_.
virtual const lduSchedule & patchSchedule() const =0
Return patch field evaluation schedule.
scalar relTol_
Convergence tolerance relative to the initial.
Definition: lduMatrix.H:186
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))