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-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  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 "typeInfo.H"
61 #include "autoPtr.H"
62 #include "runTimeSelectionTables.H"
63 #include "solverPerformance.H"
64 #include "InfoProxy.H"
65 #include "Enum.H"
66 #include "profilingTrigger.H"
67 
68 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
69 
70 namespace Foam
71 {
72 
73 // Forward Declarations
74 class lduMatrix;
75 
76 Ostream& operator<<(Ostream&, const lduMatrix&);
77 Ostream& operator<<(Ostream&, const InfoProxy<lduMatrix>&);
78 
79 
80 /*---------------------------------------------------------------------------*\
81  Class lduMatrix Declaration
82 \*---------------------------------------------------------------------------*/
83 
84 class 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 
96 public:
97 
98  // Public Types
99 
100  //- Enumerated matrix normalisation types
101  enum class normTypes : char
102  {
104  DEFAULT_NORM,
106  };
107 
108  //- Names for the normTypes
109  static const Enum<normTypes> normTypesNames_;
110 
111  //- Default maximum number of iterations for solvers (1000)
112  static constexpr const label defaultMaxIter = 1000;
114  //- Default (absolute) tolerance (1e-6)
115  static const scalar defaultTolerance;
116 
117 
118  //- Abstract base-class for lduMatrix solvers
119  class solver
120  {
121  protected:
122 
123  // Protected Data
124 
126  const lduMatrix& matrix_;
130 
131  //- Dictionary of solution controls
133 
134  //- Verbosity level for solver output statements
135  int log_;
137  //- Minimum number of iterations in the solver
138  label minIter_;
140  //- Maximum number of iterations in the solver
141  label maxIter_;
142 
143  //- The normalisation type
145 
146  //- Final convergence tolerance
147  scalar tolerance_;
148 
149  //- Convergence tolerance relative to the initial
150  scalar relTol_;
151 
152  //- Profiling instrumentation
155 
156  // Protected Member Functions
157 
158  //- Read the control parameters from controlDict_
159  virtual void readControls();
160 
161 
162  public:
163 
164  //- Runtime type information
165  virtual const word& type() const = 0;
166 
167 
168  // Declare run-time constructor selection tables
171  (
172  autoPtr,
173  solver,
174  symMatrix,
175  (
176  const word& fieldName,
177  const lduMatrix& matrix,
181  const dictionary& solverControls
182  ),
183  (
184  fieldName,
185  matrix,
188  interfaces,
189  solverControls
190  )
191  );
192 
194  (
195  autoPtr,
196  solver,
197  asymMatrix,
198  (
199  const word& fieldName,
200  const lduMatrix& matrix,
204  const dictionary& solverControls
205  ),
206  (
207  fieldName,
208  matrix,
211  interfaces,
212  solverControls
213  )
214  );
215 
216 
217  // Constructors
218 
219  //- Construct solver for given field name, matrix etc
220  solver
221  (
222  const word& fieldName,
223  const lduMatrix& matrix,
227  const dictionary& solverControls
228  );
229 
230  // Selectors
231 
232  //- Return a new solver of given type
233  static autoPtr<solver> New
234  (
235  const word& solverName,
236  const word& fieldName,
237  const lduMatrix& matrix,
241  const dictionary& solverControls
242  );
243 
244  //- Return a new solver given dictionary
245  static autoPtr<solver> New
246  (
247  const word& fieldName,
248  const lduMatrix& matrix,
252  const dictionary& solverControls
253  );
254 
255 
256 
257  //- Destructor
258  virtual ~solver() = default;
259 
260 
261  // Member Functions
262 
263  const word& fieldName() const noexcept
264  {
265  return fieldName_;
266  }
267 
268  const lduMatrix& matrix() const noexcept
269  {
270  return matrix_;
271  }
272 
274  {
275  return interfaceBouCoeffs_;
276  }
277 
278  const FieldField<Field, scalar>& interfaceIntCoeffs() const noexcept
279  {
280  return interfaceIntCoeffs_;
281  }
282 
284  {
285  return interfaces_;
286  }
287 
288 
289  //- Read and reset the solver parameters from the given stream
290  virtual void read(const dictionary&);
291 
292  //- Solve with given field and rhs
293  virtual solverPerformance solve
294  (
295  scalarField& psi,
296  const scalarField& source,
297  const direction cmpt=0
298  ) const = 0;
299 
300  //- Solve with given field and rhs (in solveScalar precision).
301  // Default is to call solve routine
303  (
305  const solveScalarField& source,
306  const direction cmpt=0
307  ) const;
308 
309  //- Return the matrix norm using the specified norm method
311  (
312  const solveScalarField& psi,
313  const solveScalarField& source,
314  const solveScalarField& Apsi,
315  solveScalarField& tmpField,
316  const lduMatrix::normTypes normType
317  ) const;
318 
319  //- Return the matrix norm used to normalise the residual for the
320  //- stopping criterion
322  (
323  const solveScalarField& psi,
324  const solveScalarField& source,
325  const solveScalarField& Apsi,
326  solveScalarField& tmpField
327  ) const
328  {
329  return this->normFactor(psi, source, Apsi, tmpField, normType_);
330  }
331  };
332 
333 
334  //- Abstract base-class for lduMatrix smoothers
335  class smoother
336  {
337  protected:
338 
339  // Protected Data
340 
342  const lduMatrix& matrix_;
346 
347 
348  public:
349 
350  //- Find the smoother name (directly or from a sub-dictionary)
351  static word getName(const dictionary&);
352 
353  //- Runtime type information
354  virtual const word& type() const = 0;
355 
356 
357  // Declare run-time constructor selection tables
358 
360  (
361  autoPtr,
362  smoother,
363  symMatrix,
364  (
365  const word& fieldName,
366  const lduMatrix& matrix,
370  ),
371  (
372  fieldName,
373  matrix,
376  interfaces
377  )
378  );
379 
381  (
382  autoPtr,
383  smoother,
384  asymMatrix,
385  (
386  const word& fieldName,
387  const lduMatrix& matrix,
391  ),
392  (
397  interfaces
398  )
399  );
400 
401 
402  // Constructors
403 
404  //- Construct for given field name, matrix etc
405  smoother
406  (
407  const word& fieldName,
408  const lduMatrix& matrix,
412  );
413 
414 
415  // Selectors
416 
417  //- Return a new smoother
418  static autoPtr<smoother> New
419  (
420  const word& fieldName,
421  const lduMatrix& matrix,
425  const dictionary& solverControls
426  );
427 
428 
429  //- Destructor
430  virtual ~smoother() = default;
431 
432 
433  // Member Functions
434 
435  const word& fieldName() const noexcept
436  {
437  return fieldName_;
438  }
439 
440  const lduMatrix& matrix() const noexcept
441  {
442  return matrix_;
443  }
444 
446  {
447  return interfaceBouCoeffs_;
448  }
449 
450  const FieldField<Field, scalar>& interfaceIntCoeffs() const noexcept
451  {
452  return interfaceIntCoeffs_;
453  }
454 
456  {
457  return interfaces_;
458  }
459 
460 
461  //- Smooth the solution for a given number of sweeps
462  virtual void smooth
463  (
465  const scalarField& source,
466  const direction cmpt,
467  const label nSweeps
468  ) const = 0;
469 
470  //- Smooth the solution for a given number of sweeps
471  virtual void scalarSmooth
472  (
474  const solveScalarField& source,
475  const direction cmpt,
476  const label nSweeps
477  ) const = 0;
478  };
479 
480 
481  //- Abstract base-class for lduMatrix preconditioners
482  class preconditioner
483  {
484  protected:
485 
486  // Protected Data
487 
488  //- Reference to the base-solver this preconditioner is used with
489  const solver& solver_;
490 
491 
492  public:
493 
494  //- Find the preconditioner name (directly or from a sub-dictionary)
495  static word getName(const dictionary&);
497  //- Runtime type information
498  virtual const word& type() const = 0;
499 
500 
501  // Declare run-time constructor selection tables
502 
504  (
505  autoPtr,
507  symMatrix,
508  (
509  const solver& sol,
510  const dictionary& solverControls
511  ),
512  (sol, solverControls)
513  );
514 
516  (
517  autoPtr,
519  asymMatrix,
520  (
521  const solver& sol,
522  const dictionary& solverControls
523  ),
524  (sol, solverControls)
525  );
526 
527 
528  // Constructors
529 
530  //- Construct for given solver
531  explicit preconditioner(const solver& sol)
532  :
533  solver_(sol)
534  {}
535 
536 
537  // Selectors
538 
539  //- Return a new preconditioner
541  (
542  const solver& sol,
543  const dictionary& solverControls
544  );
545 
546 
547  //- Destructor
548  virtual ~preconditioner() = default;
550 
551  // Member Functions
552 
553  //- Read and reset the preconditioner parameters
554  //- from the given stream
555  virtual void read(const dictionary&)
556  {}
557 
558  //- Return wA the preconditioned form of residual rA
559  virtual void precondition
560  (
561  solveScalarField& wA,
562  const solveScalarField& rA,
563  const direction cmpt=0
564  ) const = 0;
565 
566  //- Return wT the transpose-matrix preconditioned form of
567  //- residual rT.
568  // This is only required for preconditioning asymmetric matrices.
569  virtual void preconditionT
570  (
571  solveScalarField& wT,
572  const solveScalarField& rT,
573  const direction cmpt=0
574  ) const
575  {
577  }
578 
579  //- Signal end of solver
580  virtual void setFinished(const solverPerformance& perf) const
581  {}
582  };
583 
584 
585  // Static Data
586 
587  // Declare name of the class and its debug switch
588  ClassName("lduMatrix");
589 
590 
591  // Constructors
592 
593  //- Construct given an LDU addressed mesh.
594  // The coefficients are initially empty for subsequent setting.
595  lduMatrix(const lduMesh&);
596 
597  //- Construct as copy
598  lduMatrix(const lduMatrix&);
599 
600  //- Construct as copy or re-use as specified.
601  lduMatrix(lduMatrix&, bool reuse);
602 
603  //- Construct given an LDU addressed mesh and an Istream
604  //- from which the coefficients are read
605  lduMatrix(const lduMesh&, Istream&);
607 
608  //- Destructor
609  ~lduMatrix();
610 
611 
612  // Member Functions
613 
614  // Access to addressing
615 
616  //- Return the LDU mesh from which the addressing is obtained
617  const lduMesh& mesh() const noexcept
618  {
619  return lduMesh_;
620  }
621 
622  //- Set the LDU mesh containing the addressing is obtained
623  void setLduMesh(const lduMesh& m)
624  {
625  lduMesh_ = m;
626  }
627 
628  //- Return the LDU addressing
629  const lduAddressing& lduAddr() const
630  {
631  return mesh().lduAddr();
632  }
633 
634  //- Return the patch evaluation schedule
635  const lduSchedule& patchSchedule() const
636  {
637  return lduAddr().patchSchedule();
638  }
639 
640 
641  // Access to coefficients
642 
643  scalarField& lower();
644  scalarField& diag();
645  scalarField& upper();
646 
647  // Size with externally provided sizes (for constructing with 'fake'
648  // mesh in GAMG)
649 
650  scalarField& lower(const label size);
651  scalarField& diag(const label nCoeffs);
652  scalarField& upper(const label nCoeffs);
653 
654 
655  const scalarField& lower() const;
656  const scalarField& diag() const;
657  const scalarField& upper() const;
658 
659  bool hasDiag() const noexcept
660  {
661  return (diagPtr_);
662  }
663 
664  bool hasUpper() const noexcept
665  {
666  return (upperPtr_);
667  }
669  bool hasLower() const noexcept
670  {
671  return (lowerPtr_);
672  }
673 
674  bool diagonal() const noexcept
675  {
676  return (diagPtr_ && !lowerPtr_ && !upperPtr_);
677  }
678 
679  bool symmetric() const noexcept
680  {
681  return (diagPtr_ && (!lowerPtr_ && upperPtr_));
682  }
683 
684  bool asymmetric() const noexcept
685  {
686  return (diagPtr_ && lowerPtr_ && upperPtr_);
687  }
688 
689 
690  // Operations
691 
692  void sumDiag();
693  void negSumDiag();
694 
695  void sumMagOffDiag(scalarField& sumOff) const;
696 
697  //- Matrix multiplication with updated interfaces.
698  void Amul
699  (
701  const tmp<solveScalarField>&,
702  const FieldField<Field, scalar>&,
704  const direction cmpt
705  ) const;
706 
707  //- Matrix transpose multiplication with updated interfaces.
708  void Tmul
709  (
711  const tmp<solveScalarField>&,
712  const FieldField<Field, scalar>&,
714  const direction cmpt
715  ) const;
716 
717 
718  //- Sum the coefficients on each row of the matrix
719  void sumA
720  (
724  ) const;
725 
727  void residual
728  (
729  solveScalarField& rA,
730  const solveScalarField& psi,
731  const scalarField& source,
732  const FieldField<Field, scalar>& interfaceBouCoeffs,
733  const lduInterfaceFieldPtrsList& interfaces,
734  const direction cmpt
735  ) const;
736 
738  (
739  const solveScalarField& psi,
740  const scalarField& source,
741  const FieldField<Field, scalar>& interfaceBouCoeffs,
742  const lduInterfaceFieldPtrsList& interfaces,
743  const direction cmpt
744  ) const;
745 
746 
747  //- Initialise the update of interfaced interfaces
748  //- for matrix operations
750  (
751  const bool add,
752  const FieldField<Field, scalar>& interfaceCoeffs,
753  const lduInterfaceFieldPtrsList& interfaces,
754  const solveScalarField& psiif,
755  solveScalarField& result,
756  const direction cmpt
757  ) const;
758 
759  //- Update interfaced interfaces for matrix operations
761  (
762  const bool add,
763  const FieldField<Field, scalar>& interfaceCoeffs,
764  const lduInterfaceFieldPtrsList& interfaces,
765  const solveScalarField& psiif,
767  const direction cmpt,
768  const label startRequest // starting request (for non-blocking)
769  ) const;
770 
771  //- Set the residual field using an IOField on the object registry
772  //- if it exists
773  void setResidualField
774  (
775  const scalarField& residual,
776  const word& fieldName,
777  const bool initial
778  ) const;
779 
780  template<class Type>
781  tmp<Field<Type>> H(const Field<Type>&) const;
782 
783  template<class Type>
784  tmp<Field<Type>> H(const tmp<Field<Type>>&) const;
785 
787 
788  template<class Type>
789  tmp<Field<Type>> faceH(const Field<Type>&) const;
790 
791  template<class Type>
792  tmp<Field<Type>> faceH(const tmp<Field<Type>>&) const;
793 
794 
795  // Info
796 
797  //- Return info proxy,
798  //- used to print matrix information to a stream
799  InfoProxy<lduMatrix> info() const noexcept { return *this; }
800 
801 
802  // Member Operators
803 
804  void operator=(const lduMatrix&);
805 
806  void negate();
807 
808  void operator+=(const lduMatrix&);
809  void operator-=(const lduMatrix&);
810 
811  void operator*=(const scalarField&);
812  void operator*=(scalar);
813 
814 
815  // Ostream Operators
816 
817  friend Ostream& operator<<(Ostream&, const lduMatrix&);
818  friend Ostream& operator<<(Ostream&, const InfoProxy<lduMatrix>&);
819 };
820 
821 
822 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
823 
824 } // End namespace Foam
825 
826 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
827 
828 #ifdef NoRepository
829  #include "lduMatrixTemplates.C"
830 #endif
831 
832 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
833 
834 #endif
835 
836 // ************************************************************************* //
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:636
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:396
const FieldField< Field, scalar > & interfaceBouCoeffs() const noexcept
Definition: lduMatrix.H:311
tmp< Field< Type > > faceH(const Field< Type > &) const
const lduInterfaceFieldPtrsList & interfaces() const noexcept
Definition: lduMatrix.H:516
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
Definition: lduMatrix.H:786
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:511
uint8_t direction
Definition: direction.H:46
const FieldField< Field, scalar > & interfaceIntCoeffs_
Definition: lduMatrix.H:138
const lduMatrix & matrix() const noexcept
Definition: lduMatrix.H:501
preconditioner(const solver &sol)
Construct for given solver.
Definition: lduMatrix.H:606
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:558
lduMatrix(const lduMesh &)
Construct given an LDU addressed mesh.
Definition: lduMatrix.C:54
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:327
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:137
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:656
pTraits< solveScalar >::cmptType cmptType
Component type.
Definition: Field.H:123
const lduMatrix & matrix() const noexcept
Definition: lduMatrix.H:306
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:920
List< lduScheduleEntry > lduSchedule
A List of lduSchedule entries.
Definition: lduSchedule.H:46
Triggers for starting/stopping code profiling.
scalarField & upper()
Definition: lduMatrix.C:208
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:118
int log_
Verbosity level for solver output statements.
Definition: lduMatrix.H:149
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:113
label minIter_
Minimum number of iterations in the solver.
Definition: lduMatrix.H:154
lduMatrix::normTypes normType_
The normalisation type.
Definition: lduMatrix.H:164
virtual ~preconditioner()=default
Destructor.
const lduMatrix & matrix_
Definition: lduMatrix.H:393
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:742
Abstract base-class for lduMatrix smoothers.
Definition: lduMatrix.H:386
Forward declarations of the specialisations of Field<T> for scalar, vector and tensor.
scalar tolerance_
Final convergence tolerance.
Definition: lduMatrix.H:169
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:159
static autoPtr< preconditioner > New(const solver &sol, const dictionary &solverControls)
Return a new preconditioner.
Abstract base-class for lduMatrix solvers.
Definition: lduMatrix.H:129
const word & fieldName() const noexcept
Definition: lduMatrix.H:496
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.
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:321
const FieldField< Field, scalar > & interfaceBouCoeffs() const noexcept
Definition: lduMatrix.H:506
static const scalar defaultTolerance
Default (absolute) tolerance (1e-6)
Definition: lduMatrix.H:123
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:394
normTypes
Enumerated matrix normalisation types.
Definition: lduMatrix.H:103
void operator=(const lduMatrix &)
const lduMatrix & matrix_
Definition: lduMatrix.H:136
bool hasUpper() const noexcept
Definition: lduMatrix.H:771
bool asymmetric() const noexcept
Definition: lduMatrix.H:791
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:139
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
Definition: lduMatrix.H:781
virtual ~smoother()=default
Destructor.
virtual const word & type() const =0
Runtime type information.
dictionary controlDict_
Dictionary of solution controls.
Definition: lduMatrix.H:144
"default" norm (== L1_scaled)
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:79
"none" norm (returns 1)
Abstract base-class for lduMatrix preconditioners.
Definition: lduMatrix.H:549
A helper class for outputting values to Ostream.
Definition: ensightCells.H:43
const lduMesh & mesh() const noexcept
Return the LDU mesh from which the addressing is obtained.
Definition: lduMatrix.H:718
const lduAddressing & lduAddr() const
Return the LDU addressing.
Definition: lduMatrix.H:734
scalarField & lower()
Definition: lduMatrix.C:179
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 is obtained.
Definition: lduMatrix.H:726
virtual ~solver()=default
Destructor.
const FieldField< Field, scalar > & interfaceIntCoeffs() const noexcept
Definition: lduMatrix.H:316
const FieldField< Field, scalar > & interfaceIntCoeffs_
Definition: lduMatrix.H:395
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:766
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
scalarField & diag()
Definition: lduMatrix.C:197
profilingTrigger profiling_
Profiling instrumentation.
Definition: lduMatrix.H:179
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:686
bool hasLower() const noexcept
Definition: lduMatrix.H:776
const word & fieldName() const noexcept
Definition: lduMatrix.H:301
virtual void setFinished(const solverPerformance &perf) const
Signal end of solver.
Definition: lduMatrix.H:668
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.
~lduMatrix()
Destructor.
Definition: lduMatrix.C:160
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:174
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))