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