LduMatrix.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::LduMatrix
26 
27 Description
28  LduMatrix is a general matrix class in which the coefficients are
29  stored as three arrays, one for the upper triangle, one for the
30  lower triangle and a third for the diagonal.
31 
32  Addressing arrays must be supplied for the upper and lower triangles.
33 
34  Note:
35  It might be better if this class were organised as a hierarchy starting
36  from an empty matrix, then deriving diagonal, symmetric and asymmetric
37  matrices.
38 
39 SourceFiles
40  LduMatrixATmul.C
41  LduMatrix.C
42  LduMatrixOperations.C
43  LduMatrixSolver.C
44  LduMatrixPreconditioner.C
45  LduMatrixTests.C
46  LduMatrixUpdateMatrixInterfaces.C
47 
48 \*---------------------------------------------------------------------------*/
49 
50 #ifndef LduMatrix_H
51 #define LduMatrix_H
52 
53 #include "lduMesh.H"
54 #include "Field.H"
55 #include "FieldField.H"
57 #include "SolverPerformance.H"
58 #include "typeInfo.H"
59 #include "autoPtr.H"
60 #include "runTimeSelectionTables.H"
61 
62 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
63 
64 namespace Foam
65 {
66 
67 // Forward declaration of friend functions and operators
68 
69 template<class Type, class DType, class LUType>
70 class LduMatrix;
71 
72 template<class Type, class DType, class LUType>
73 Ostream& operator<<
74 (
75  Ostream&,
77 );
78 
79 
80 /*---------------------------------------------------------------------------*\
81  Class LduMatrix Declaration
82 \*---------------------------------------------------------------------------*/
83 
84 template<class Type, class DType, class LUType>
85 class LduMatrix
86 {
87  // private data
88 
89  //- LDU mesh reference
90  const lduMesh& lduMesh_;
91 
92  //- Diagonal coefficients
93  Field<DType> *diagPtr_;
94 
95  //- Off-diagonal coefficients
96  Field<LUType> *upperPtr_, *lowerPtr_;
97 
98  //- Source
99  Field<Type> *sourcePtr_;
100 
101  //- Field interfaces (processor patches etc.)
103 
104  //- Off-diagonal coefficients for interfaces
105  FieldField<Field, LUType> interfacesUpper_, interfacesLower_;
106 
107 
108 public:
110  friend class SolverPerformance<Type>;
111 
112  //- Abstract base-class for LduMatrix solvers
113  class solver
114  {
115  protected:
116 
117  // Protected data
121 
122  //- Dictionary of controls
124 
125  //- Default maximum number of iterations in the solver
126  static const label defaultMaxIter_ = 1000;
127 
128  //- Maximum number of iterations in the solver
129  label maxIter_;
130 
131  //- Minimum number of iterations in the solver
132  label minIter_;
133 
134  //- Final convergence tolerance
135  Type tolerance_;
136 
137  //- Convergence tolerance relative to the initial
138  Type relTol_;
139 
140 
141  // Protected Member Functions
142 
143  //- Read a control parameter from controlDict
144  template<class T>
145  inline void readControl
146  (
147  const dictionary& controlDict,
148  T& control,
149  const word& controlName
150  );
151 
152 
153  //- Read the control parameters from the controlDict_
154  virtual void readControls();
155 
156 
157  public:
158 
159  //- Runtime type information
160  virtual const word& type() const = 0;
161 
162 
163  // Declare run-time constructor selection tables
164 
166  (
167  autoPtr,
168  solver,
169  symMatrix,
170  (
171  const word& fieldName,
172  const LduMatrix<Type, DType, LUType>& matrix,
173  const dictionary& solverDict
174  ),
175  (
176  fieldName,
177  matrix,
178  solverDict
179  )
180  );
181 
183  (
184  autoPtr,
185  solver,
186  asymMatrix,
187  (
188  const word& fieldName,
189  const LduMatrix<Type, DType, LUType>& matrix,
190  const dictionary& solverDict
191  ),
192  (
193  fieldName,
194  matrix,
195  solverDict
196  )
197  );
198 
199 
200  // Constructors
201 
202  solver
203  (
204  const word& fieldName,
205  const LduMatrix<Type, DType, LUType>& matrix,
206  const dictionary& solverDict
207  );
208 
209 
210  // Selectors
211 
212  //- Return a new solver
213  static autoPtr<solver> New
214  (
215  const word& fieldName,
216  const LduMatrix<Type, DType, LUType>& matrix,
217  const dictionary& solverDict
218  );
219 
220 
221  // Destructor
223  virtual ~solver()
224  {}
225 
226 
227  // Member Functions
228 
229  // Access
231  const word& fieldName() const
232  {
233  return fieldName_;
234  }
237  {
238  return matrix_;
239  }
240 
241 
242  //- Read and reset the solver parameters from the given dictionary
243  virtual void read(const dictionary& solverDict);
244 
246  (
248  ) const = 0;
249 
250  //- Return the matrix norm used to normalise the residual for the
251  // stopping criterion
252  Type normFactor
253  (
254  const Field<Type>& psi,
255  const Field<Type>& Apsi,
256  Field<Type>& tmpField
257  ) const;
258  };
259 
260 
261  //- Abstract base-class for LduMatrix smoothers
262  class smoother
263  {
264  protected:
265 
266  // Protected data
270 
271 
272  public:
273 
274  //- Runtime type information
275  virtual const word& type() const = 0;
276 
277 
278  // Declare run-time constructor selection tables
279 
281  (
282  autoPtr,
283  smoother,
284  symMatrix,
285  (
286  const word& fieldName,
287  const LduMatrix<Type, DType, LUType>& matrix
288  ),
289  (
290  fieldName,
291  matrix
292  )
293  );
294 
296  (
297  autoPtr,
298  smoother,
299  asymMatrix,
300  (
301  const word& fieldName,
302  const LduMatrix<Type, DType, LUType>& matrix
303  ),
304  (
305  fieldName,
306  matrix
307  )
308  );
309 
310 
311  // Constructors
312 
313  smoother
314  (
315  const word& fieldName,
316  const LduMatrix<Type, DType, LUType>& matrix
317  );
318 
319 
320  // Selectors
321 
322  //- Return a new smoother
323  static autoPtr<smoother> New
324  (
325  const word& fieldName,
326  const LduMatrix<Type, DType, LUType>& matrix,
327  const dictionary& smootherDict
328  );
329 
330 
331  // Destructor
333  virtual ~smoother()
334  {}
335 
336 
337  // Member Functions
338 
339  // Access
341  const word& fieldName() const
342  {
343  return fieldName_;
344  }
347  {
348  return matrix_;
349  }
350 
351 
352  //- Smooth the solution for a given number of sweeps
353  virtual void smooth
354  (
355  Field<Type>& psi,
356  const label nSweeps
357  ) const = 0;
358  };
359 
360 
361  //- Abstract base-class for LduMatrix preconditioners
362  class preconditioner
363  {
364  protected:
365 
366  // Protected data
367 
368  //- Reference to the base-solver this preconditioner is used with
369  const solver& solver_;
370 
371 
372  public:
373 
374  //- Runtime type information
375  virtual const word& type() const = 0;
376 
377 
378  // Declare run-time constructor selection tables
379 
381  (
382  autoPtr,
384  symMatrix,
385  (
386  const solver& sol,
387  const dictionary& preconditionerDict
388  ),
389  (sol, preconditionerDict)
390  );
391 
393  (
394  autoPtr,
396  asymMatrix,
397  (
398  const solver& sol,
399  const dictionary& preconditionerDict
400  ),
401  (sol, preconditionerDict)
402  );
403 
404 
405  // Constructors
406 
408  (
409  const solver& sol
410  )
411  :
412  solver_(sol)
413  {}
414 
415 
416  // Selectors
417 
418  //- Return a new preconditioner
420  (
421  const solver& sol,
422  const dictionary& preconditionerDict
423  );
424 
425 
426  // Destructor
428  virtual ~preconditioner()
429  {}
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  //- Copy constructor
474 
475  //- Copy constructor 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();
531  {
532  return interfacesUpper_;
533  }
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;
547  {
548  return interfacesUpper_;
549  }
552  {
553  return interfacesLower_;
554  }
555 
557  bool hasDiag() const
558  {
559  return (diagPtr_);
560  }
562  bool hasUpper() const
563  {
564  return (upperPtr_);
565  }
567  bool hasLower() const
568  {
569  return (lowerPtr_);
570  }
572  bool hasSource() const
573  {
574  return (sourcePtr_);
575  }
577  bool diagonal() const
578  {
579  return
580  (
581  diagPtr_
582  && Pstream::parRun()
583  ?
584  !lowerPtr_ && !upperPtr_
585  :
586  !(lowerPtr_ && lowerPtr_->size())
587  && !(upperPtr_ && upperPtr_->size())
588  );
589  }
591  bool symmetric() const
592  {
593  return
594  (
595  diagPtr_
596  && Pstream::parRun()
597  ?
598  !lowerPtr_ && upperPtr_
599  :
600  !(lowerPtr_ && lowerPtr_->size())
601  && (upperPtr_ && upperPtr_->size())
602  );
603  }
605  bool asymmetric() const
606  {
607  return
608  (
609  diagPtr_
610  && Pstream::parRun()
611  ?
612  lowerPtr_ && upperPtr_
613  :
614  (lowerPtr_ && lowerPtr_->size())
615  && (upperPtr_ && upperPtr_->size())
616  );
617  }
618 
619 
620  // operations
621 
622  void sumDiag();
623  void negSumDiag();
624 
625  void sumMagOffDiag(Field<LUType>& sumOff) const;
626 
627  //- Matrix multiplication
628  void Amul(Field<Type>&, const tmp<Field<Type>>&) const;
629 
630  //- Matrix transpose multiplication
631  void Tmul(Field<Type>&, const tmp<Field<Type>>&) const;
632 
633 
634  //- Sum the coefficients on each row of the matrix
635  void sumA(Field<Type>&) const;
636 
637 
638  void residual(Field<Type>& rA, const Field<Type>& psi) const;
639 
640  tmp<Field<Type>> residual(const Field<Type>& psi) const;
641 
642 
643  //- Initialise the update of interfaced interfaces
644  // for matrix operations
646  (
647  const FieldField<Field, LUType>& interfaceCoeffs,
648  const Field<Type>& psiif,
649  Field<Type>& result
650  ) const;
651 
652  //- Update interfaced interfaces for matrix operations
654  (
655  const FieldField<Field, LUType>& interfaceCoeffs,
656  const Field<Type>& psiif,
657  Field<Type>& result
658  ) const;
659 
660 
661  tmp<Field<Type>> H(const Field<Type>&) const;
662  tmp<Field<Type>> H(const tmp<Field<Type>>&) const;
663 
664  tmp<Field<Type>> faceH(const Field<Type>&) const;
665  tmp<Field<Type>> faceH(const tmp<Field<Type>>&) const;
666 
667 
668  // Member Operators
669 
671 
672  void negate();
673 
676 
677  void operator*=(const scalarField&);
678  void operator*=(scalar);
679 
680 
681  // Ostream operator
682 
683  friend Ostream& operator<< <Type, DType, LUType>
684  (
685  Ostream&,
687  );
688 };
689 
690 
691 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
692 
693 } // End namespace Foam
694 
695 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
697 #define makeLduMatrix(Type, DType, LUType) \
698  \
699 typedef Foam::LduMatrix<Type, DType, LUType> \
700  ldu##Type##DType##LUType##Matrix; \
701  \
702 defineNamedTemplateTypeNameAndDebug(ldu##Type##DType##LUType##Matrix, 0); \
703  \
704  \
705 typedef LduMatrix<Type, DType, LUType>::smoother \
706  ldu##Type##DType##LUType##Smoother; \
707  \
708 defineTemplateRunTimeSelectionTable \
709 ( \
710  ldu##Type##DType##LUType##Smoother, \
711  symMatrix \
712 ); \
713  \
714 defineTemplateRunTimeSelectionTable \
715 ( \
716  ldu##Type##DType##LUType##Smoother, \
717  asymMatrix \
718 ); \
719  \
720  \
721 typedef LduMatrix<Type, DType, LUType>::preconditioner \
722  ldu##Type##DType##LUType##Preconditioner; \
723  \
724 defineTemplateRunTimeSelectionTable \
725 ( \
726  ldu##Type##DType##LUType##Preconditioner, \
727  symMatrix \
728 ); \
729  \
730 defineTemplateRunTimeSelectionTable \
731 ( \
732  ldu##Type##DType##LUType##Preconditioner, \
733  asymMatrix \
734 ); \
735  \
736  \
737 typedef LduMatrix<Type, DType, LUType>::solver \
738  ldu##Type##DType##LUType##Solver; \
739  \
740 defineTemplateRunTimeSelectionTable \
741 ( \
742  ldu##Type##DType##LUType##Solver, \
743  symMatrix \
744 ); \
745  \
746 defineTemplateRunTimeSelectionTable \
747 ( \
748  ldu##Type##DType##LUType##Solver, \
749  asymMatrix \
750 );
751 
753 #define makeLduPreconditioner(Precon, Type, DType, LUType) \
754  \
755 typedef Precon<Type, DType, LUType> \
756  Precon##Type##DType##LUType##Preconditioner; \
757 defineNamedTemplateTypeNameAndDebug \
758 ( \
759  Precon##Type##DType##LUType##Preconditioner, \
760  0 \
761 );
763 #define makeLduSymPreconditioner(Precon, Type, DType, LUType) \
764  \
765 LduMatrix<Type, DType, LUType>::preconditioner:: \
766 addsymMatrixConstructorToTable<Precon##Type##DType##LUType##Preconditioner> \
767 add##Precon##Type##DType##LUType##PreconditionerSymMatrixConstructorToTable_;
769 #define makeLduAsymPreconditioner(Precon, Type, DType, LUType) \
770  \
771 LduMatrix<Type, DType, LUType>::preconditioner:: \
772 addasymMatrixConstructorToTable<Precon##Type##DType##LUType##Preconditioner> \
773 add##Precon##Type##DType##LUType##PreconditionerAsymMatrixConstructorToTable_;
774 
776 #define makeLduSmoother(Smoother, Type, DType, LUType) \
777  \
778 typedef Smoother<Type, DType, LUType> \
779  Smoother##Type##DType##LUType##Smoother; \
780  \
781 defineNamedTemplateTypeNameAndDebug \
782 ( \
783  Smoother##Type##DType##LUType##Smoother, \
784  0 \
785 );
787 #define makeLduSymSmoother(Smoother, Type, DType, LUType) \
788  \
789 LduMatrix<Type, DType, LUType>::smoother:: \
790  addsymMatrixConstructorToTable<Smoother##Type##DType##LUType##Smoother> \
791  add##Smoother##Type##DType##LUType##SymMatrixConstructorToTable_;
793 #define makeLduAsymSmoother(Smoother, Type, DType, LUType) \
794  \
795 LduMatrix<Type, DType, LUType>::smoother:: \
796  addasymMatrixConstructorToTable<Smoother##Type##DType##LUType##Smoother> \
797  add##Smoother##Type##DType##LUType##AsymMatrixConstructorToTable_;
798 
800 #define makeLduSolver(Solver, Type, DType, LUType) \
801  \
802 typedef Solver<Type, DType, LUType> \
803  Solver##Type##DType##LUType##Solver; \
804  \
805 defineNamedTemplateTypeNameAndDebug \
806 ( \
807  Solver##Type##DType##LUType##Solver, \
808  0 \
809 );
811 #define makeLduSymSolver(Solver, Type, DType, LUType) \
812  \
813 LduMatrix<Type, DType, LUType>::solver:: \
814  addsymMatrixConstructorToTable<Solver##Type##DType##LUType##Solver> \
815  add##Solver##Type##DType##LUType##SymMatrixConstructorToTable_;
817 #define makeLduAsymSolver(Solver, Type, DType, LUType) \
818  \
819 LduMatrix<Type, DType, LUType>::solver:: \
820  addasymMatrixConstructorToTable<Solver##Type##DType##LUType##Solver> \
821  add##Solver##Type##DType##LUType##AsymMatrixConstructorToTable_;
822 
823 
824 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
825 
826 #ifdef NoRepository
827  #include "LduMatrixI.H"
828  #include "LduMatrix.C"
829 #endif
830 
831 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
832 
833 #endif
834 
835 // ************************************************************************* //
bool hasSource() const
Definition: LduMatrix.H:571
bool symmetric() const
Definition: LduMatrix.H:590
virtual const word & type() const =0
Runtime type information.
runTime controlDict().lookup("adjustTimeStep") >> adjustTimeStep
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
Field< Type > & source()
Definition: LduMatrix.C:248
bool hasLower() const
Definition: LduMatrix.H:566
FieldField< Field, LUType > & interfacesLower()
Definition: LduMatrix.H:534
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
bool hasUpper() const
Definition: LduMatrix.H:561
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
LduMatrix(const lduMesh &)
Construct given an LDU addressed mesh.
Definition: LduMatrix.C:32
void sumMagOffDiag(Field< LUType > &sumOff) const
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:59
Field< LUType > & lower()
Definition: LduMatrix.C:225
label minIter_
Minimum number of iterations in the solver.
Definition: LduMatrix.H:131
const LduMatrix< Type, DType, LUType > & matrix() const
Definition: LduMatrix.H:235
void sumA(Field< Type > &) const
Sum the coefficients on each row of the matrix.
Generic field type.
Definition: FieldField.H:51
void smooth(volScalarField &field, const scalar coeff)
Definition: fvcSmooth.C:35
const word & fieldName() const
Definition: LduMatrix.H:230
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.
dictionary controlDict_
Dictionary of controls.
Definition: LduMatrix.H:122
ClassName("LduMatrix")
void Amul(Field< Type > &, const tmp< Field< Type >> &) const
Matrix multiplication.
void initMatrixInterfaces(const FieldField< Field, LUType > &interfaceCoeffs, const Field< Type > &psiif, Field< Type > &result) const
Initialise the update of interfaced interfaces.
Abstract base-class for LduMatrix smoothers.
Definition: LduMatrix.H:261
A class for handling words, derived from string.
Definition: word.H:59
void readControl(const dictionary &controlDict, T &control, const word &controlName)
Read a control parameter from controlDict.
Definition: LduMatrixI.H:31
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics...
void updateMatrixInterfaces(const FieldField< Field, LUType > &interfaceCoeffs, const Field< Type > &psiif, Field< Type > &result) const
Update interfaced interfaces for matrix operations.
Abstract base-class for LduMatrix preconditioners.
Definition: LduMatrix.H:361
tmp< Field< Type > > H(const Field< Type > &) const
const lduSchedule & patchSchedule() const
Return the patch evaluation schedule.
Definition: LduMatrix.H:504
const volScalarField & psi
void operator*=(const scalarField &)
bool hasDiag() const
Definition: LduMatrix.H:556
Field< LUType > & upper()
Definition: LduMatrix.C:202
void operator-=(const LduMatrix< Type, DType, LUType > &)
Abstract base-class for LduMatrix solvers.
Definition: LduMatrix.H:112
void residual(Field< Type > &rA, const Field< Type > &psi) const
solver(const word &fieldName, const LduMatrix< Type, DType, LUType > &matrix, const dictionary &solverDict)
void Tmul(Field< Type > &, const tmp< Field< Type >> &) const
Matrix transpose multiplication.
const lduAddressing & lduAddr() const
Return the LDU addressing.
Definition: LduMatrix.H:498
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
bool asymmetric() const
Definition: LduMatrix.H:604
bool diagonal() const
Definition: LduMatrix.H:576
const LduInterfaceFieldPtrsList< Type > & interfaces() const
Return interfaces.
Definition: LduMatrix.H:510
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
Field< DType > & diag()
Definition: LduMatrix.C:190
Type relTol_
Convergence tolerance relative to the initial.
Definition: LduMatrix.H:137
FieldField< Field, LUType > & interfacesUpper()
Definition: LduMatrix.H:529
LduMatrix is a general matrix class in which the coefficients are stored as three arrays...
Definition: LduMatrix.H:69
virtual void read(const dictionary &solverDict)
Read and reset the solver parameters from the given dictionary.
label maxIter_
Maximum number of iterations in the solver.
Definition: LduMatrix.H:128
void operator+=(const LduMatrix< Type, DType, LUType > &)
tmp< Field< Type > > faceH(const Field< Type > &) const
static autoPtr< solver > New(const word &fieldName, const LduMatrix< Type, DType, LUType > &matrix, const dictionary &solverDict)
Return a new solver.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
const lduMesh & mesh() const
Return the LDU mesh from which the addressing is obtained.
Definition: LduMatrix.H:492
The class contains the addressing required by the lduMatrix: upper, lower and losort.
virtual SolverPerformance< Type > solve(Field< Type > &psi) const =0
const LduMatrix< Type, DType, LUType > & matrix_
Definition: LduMatrix.H:119
Macros to ease declaration of run-time selection tables.
virtual void readControls()
Read the control parameters from the controlDict_.
A class for managing temporary objects.
Definition: PtrList.H:53
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:353
declareRunTimeSelectionTable(autoPtr, solver, symMatrix,(const word &fieldName, const LduMatrix< Type, DType, LUType > &matrix, const dictionary &solverDict),(fieldName, matrix, solverDict))
void operator=(const LduMatrix< Type, DType, LUType > &)
List of coupled interface fields to be used in coupling.
Type tolerance_
Final convergence tolerance.
Definition: LduMatrix.H:134
Namespace for OpenFOAM.
virtual const lduSchedule & patchSchedule() const =0