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-2020 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  It might be better if this class were organised as a hierarchy starting
35  from an empty matrix, then deriving diagonal, symmetric and asymmetric
36  matrices.
37 
38 SourceFiles
39  lduMatrixATmul.C
40  lduMatrix.C
41  lduMatrixTemplates.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 "primitiveFieldsFwd.H"
55 #include "FieldField.H"
57 #include "typeInfo.H"
58 #include "autoPtr.H"
59 #include "runTimeSelectionTables.H"
60 #include "solverPerformance.H"
61 #include "InfoProxy.H"
62 
63 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 
65 namespace Foam
66 {
67 
68 // Forward declaration of friend functions and operators
69 
70 class lduMatrix;
71 
72 Ostream& operator<<(Ostream&, const lduMatrix&);
73 Ostream& operator<<(Ostream&, const InfoProxy<lduMatrix>&);
74 
75 
76 /*---------------------------------------------------------------------------*\
77  Class lduMatrix Declaration
78 \*---------------------------------------------------------------------------*/
79 
80 class lduMatrix
81 {
82  // private data
83 
84  //- LDU mesh reference
85  const lduMesh& lduMesh_;
86 
87  //- Coefficients (not including interfaces)
88  scalarField *lowerPtr_, *diagPtr_, *upperPtr_;
89 
90 
91 public:
92 
93  //- Abstract base-class for lduMatrix solvers
94  class solver
95  {
96  protected:
97 
98  // Protected data
99 
105 
106  //- Dictionary of controls
108 
109  //- Default maximum number of iterations in the solver
110  static const label defaultMaxIter_;
111 
112  //- Maximum number of iterations in the solver
113  label maxIter_;
114 
115  //- Minimum number of iterations in the solver
116  label minIter_;
117 
118  //- Final convergence tolerance
119  scalar tolerance_;
120 
121  //- Convergence tolerance relative to the initial
122  scalar relTol_;
123 
124 
125  // Protected Member Functions
126 
127  //- Read the control parameters from the controlDict_
128  virtual void readControls();
129 
130 
131  public:
132 
133  //- Runtime type information
134  virtual const word& type() const = 0;
135 
136 
137  // Declare run-time constructor selection tables
138 
140  (
141  autoPtr,
142  solver,
143  symMatrix,
144  (
145  const word& fieldName,
146  const lduMatrix& matrix,
147  const FieldField<Field, scalar>& interfaceBouCoeffs,
148  const FieldField<Field, scalar>& interfaceIntCoeffs,
149  const lduInterfaceFieldPtrsList& interfaces,
150  const dictionary& solverControls
151  ),
152  (
153  fieldName,
154  matrix,
155  interfaceBouCoeffs,
156  interfaceIntCoeffs,
157  interfaces,
158  solverControls
159  )
160  );
161 
163  (
164  autoPtr,
165  solver,
166  asymMatrix,
167  (
168  const word& fieldName,
169  const lduMatrix& matrix,
170  const FieldField<Field, scalar>& interfaceBouCoeffs,
171  const FieldField<Field, scalar>& interfaceIntCoeffs,
172  const lduInterfaceFieldPtrsList& interfaces,
173  const dictionary& solverControls
174  ),
175  (
176  fieldName,
177  matrix,
178  interfaceBouCoeffs,
179  interfaceIntCoeffs,
180  interfaces,
181  solverControls
182  )
183  );
184 
185 
186  // Constructors
187 
188  solver
189  (
190  const word& fieldName,
191  const lduMatrix& matrix,
192  const FieldField<Field, scalar>& interfaceBouCoeffs,
193  const FieldField<Field, scalar>& interfaceIntCoeffs,
194  const lduInterfaceFieldPtrsList& interfaces,
195  const dictionary& solverControls
196  );
197 
198  // Selectors
199 
200  //- Return a new solver
201  static autoPtr<solver> New
202  (
203  const word& fieldName,
204  const lduMatrix& matrix,
205  const FieldField<Field, scalar>& interfaceBouCoeffs,
206  const FieldField<Field, scalar>& interfaceIntCoeffs,
207  const lduInterfaceFieldPtrsList& interfaces,
208  const dictionary& solverControls
209  );
210 
211 
212 
213  //- Destructor
214  virtual ~solver()
215  {}
216 
217 
218  // Member Functions
219 
220  // Access
222  const word& fieldName() const
223  {
224  return fieldName_;
225  }
227  const lduMatrix& matrix() const
228  {
229  return matrix_;
230  }
233  {
234  return interfaceBouCoeffs_;
235  }
238  {
239  return interfaceIntCoeffs_;
240  }
243  {
244  return interfaces_;
245  }
246 
247 
248  //- Read and reset the solver parameters from the given stream
249  virtual void read(const dictionary&);
250 
251  virtual solverPerformance solve
252  (
253  scalarField& psi,
254  const scalarField& source,
255  const direction cmpt=0
256  ) const = 0;
257 
258  //- Return the matrix norm used to normalise the residual for the
259  // stopping criterion
260  scalar normFactor
261  (
262  const scalarField& psi,
263  const scalarField& source,
264  const scalarField& Apsi,
265  scalarField& tmpField
266  ) const;
267  };
268 
269 
270  //- Abstract base-class for lduMatrix smoothers
271  class smoother
272  {
273  protected:
274 
275  // Protected data
282 
283 
284  public:
285 
286  //- Find the smoother name (directly or from a sub-dictionary)
287  static word getName(const dictionary&);
288 
289  //- Runtime type information
290  virtual const word& type() const = 0;
291 
292 
293  // Declare run-time constructor selection tables
294 
296  (
297  autoPtr,
298  smoother,
299  symMatrix,
300  (
301  const word& fieldName,
302  const lduMatrix& matrix,
303  const FieldField<Field, scalar>& interfaceBouCoeffs,
304  const FieldField<Field, scalar>& interfaceIntCoeffs,
305  const lduInterfaceFieldPtrsList& interfaces
306  ),
307  (
308  fieldName,
309  matrix,
310  interfaceBouCoeffs,
311  interfaceIntCoeffs,
312  interfaces
313  )
314  );
315 
317  (
318  autoPtr,
319  smoother,
320  asymMatrix,
321  (
322  const word& fieldName,
323  const lduMatrix& matrix,
324  const FieldField<Field, scalar>& interfaceBouCoeffs,
325  const FieldField<Field, scalar>& interfaceIntCoeffs,
326  const lduInterfaceFieldPtrsList& interfaces
327  ),
328  (
329  fieldName,
330  matrix,
331  interfaceBouCoeffs,
332  interfaceIntCoeffs,
333  interfaces
334  )
335  );
336 
337 
338  // Constructors
339 
340  smoother
341  (
342  const word& fieldName,
343  const lduMatrix& matrix,
344  const FieldField<Field, scalar>& interfaceBouCoeffs,
345  const FieldField<Field, scalar>& interfaceIntCoeffs,
346  const lduInterfaceFieldPtrsList& interfaces
347  );
348 
349 
350  // Selectors
351 
352  //- Return a new smoother
353  static autoPtr<smoother> New
354  (
355  const word& fieldName,
356  const lduMatrix& matrix,
357  const FieldField<Field, scalar>& interfaceBouCoeffs,
358  const FieldField<Field, scalar>& interfaceIntCoeffs,
359  const lduInterfaceFieldPtrsList& interfaces,
360  const dictionary& solverControls
361  );
362 
363 
364  //- Destructor
365  virtual ~smoother()
366  {}
367 
368 
369  // Member Functions
370 
371  // Access
373  const word& fieldName() const
374  {
375  return fieldName_;
376  }
378  const lduMatrix& matrix() const
379  {
380  return matrix_;
381  }
384  {
385  return interfaceBouCoeffs_;
386  }
389  {
390  return interfaceIntCoeffs_;
391  }
394  {
395  return interfaces_;
396  }
397 
398 
399  //- Smooth the solution for a given number of sweeps
400  virtual void smooth
401  (
402  scalarField& psi,
403  const scalarField& source,
404  const direction cmpt,
405  const label nSweeps
406  ) const = 0;
407  };
408 
409 
410  //- Abstract base-class for lduMatrix preconditioners
411  class preconditioner
412  {
413  protected:
414 
415  // Protected data
416 
417  //- Reference to the base-solver this preconditioner is used with
418  const solver& solver_;
419 
420 
421  public:
422 
423  //- Find the preconditioner name (directly or from a sub-dictionary)
424  static word getName(const dictionary&);
425 
426  //- Runtime type information
427  virtual const word& type() const = 0;
428 
429 
430  // Declare run-time constructor selection tables
431 
433  (
434  autoPtr,
436  symMatrix,
437  (
438  const solver& sol,
439  const dictionary& solverControls
440  ),
441  (sol, solverControls)
442  );
443 
445  (
446  autoPtr,
448  asymMatrix,
449  (
450  const solver& sol,
451  const dictionary& solverControls
452  ),
453  (sol, solverControls)
454  );
455 
456 
457  // Constructors
458 
460  (
461  const solver& sol
462  )
463  :
464  solver_(sol)
465  {}
466 
467 
468  // Selectors
469 
470  //- Return a new preconditioner
472  (
473  const solver& sol,
474  const dictionary& solverControls
475  );
476 
477 
478  //- Destructor
479  virtual ~preconditioner()
480  {}
481 
482 
483  // Member Functions
484 
485  //- Read and reset the preconditioner parameters
486  // from the given stream
487  virtual void read(const dictionary&)
488  {}
489 
490  //- Return wA the preconditioned form of residual rA
491  virtual void precondition
492  (
493  scalarField& wA,
494  const scalarField& rA,
495  const direction cmpt=0
496  ) const = 0;
497 
498  //- Return wT the transpose-matrix preconditioned form of
499  // residual rT.
500  // This is only required for preconditioning asymmetric matrices.
501  virtual void preconditionT
502  (
503  scalarField& wT,
504  const scalarField& rT,
505  const direction cmpt=0
506  ) const
507  {
509  }
510  };
511 
512 
513  // Static data
514 
515  // Declare name of the class and its debug switch
516  ClassName("lduMatrix");
517 
518 
519  // Constructors
520 
521  //- Construct given an LDU addressed mesh.
522  // The coefficients are initially empty for subsequent setting.
523  lduMatrix(const lduMesh&);
524 
525  //- Copy constructor
526  lduMatrix(const lduMatrix&);
527 
528  //- Copy constructor or re-use as specified.
529  lduMatrix(lduMatrix&, bool reuse);
530 
531  //- Construct given an LDU addressed mesh and an Istream
532  // from which the coefficients are read
533  lduMatrix(const lduMesh&, Istream&);
534 
535 
536  //- Destructor
537  ~lduMatrix();
538 
539 
540  // Member Functions
541 
542  // Access to addressing
543 
544  //- Return the LDU mesh from which the addressing is obtained
545  const lduMesh& mesh() const
546  {
547  return lduMesh_;
548  }
549 
550  //- Return the LDU addressing
551  const lduAddressing& lduAddr() const
552  {
553  return lduMesh_.lduAddr();
554  }
555 
556  //- Return the patch evaluation schedule
557  const lduSchedule& patchSchedule() const
558  {
559  return lduAddr().patchSchedule();
560  }
561 
562 
563  // Access to coefficients
564 
565  scalarField& lower();
566  scalarField& diag();
567  scalarField& upper();
568 
569  // Size with externally provided sizes (for constructing with 'fake'
570  // mesh in GAMG)
571 
572  scalarField& lower(const label size);
573  scalarField& diag(const label nCoeffs);
574  scalarField& upper(const label nCoeffs);
575 
576 
577  const scalarField& lower() const;
578  const scalarField& diag() const;
579  const scalarField& upper() const;
581  bool hasDiag() const
582  {
583  return (diagPtr_);
584  }
586  bool hasUpper() const
587  {
588  return (upperPtr_);
589  }
591  bool hasLower() const
592  {
593  return (lowerPtr_);
594  }
596  bool diagonal() const
597  {
598  return
599  (
600  diagPtr_
601  && Pstream::parRun()
602  ?
603  !lowerPtr_ && !upperPtr_
604  :
605  !(lowerPtr_ && lowerPtr_->size())
606  && !(upperPtr_ && upperPtr_->size())
607  );
608  }
610  bool symmetric() const
611  {
612  return
613  (
614  diagPtr_
615  && Pstream::parRun()
616  ?
617  !lowerPtr_ && upperPtr_
618  :
619  !(lowerPtr_ && lowerPtr_->size())
620  && (upperPtr_ && upperPtr_->size())
621  );
622  }
624  bool asymmetric() const
625  {
626  return
627  (
628  diagPtr_
629  && Pstream::parRun()
630  ?
631  lowerPtr_ && upperPtr_
632  :
633  (lowerPtr_ && lowerPtr_->size())
634  && (upperPtr_ && upperPtr_->size())
635  );
636  }
637 
638 
639  // operations
640 
641  void sumDiag();
642  void negSumDiag();
643 
644  void sumMagOffDiag(scalarField& sumOff) const;
645 
646  //- Matrix multiplication with updated interfaces.
647  void Amul
648  (
649  scalarField&,
650  const tmp<scalarField>&,
653  const direction cmpt
654  ) const;
655 
656  //- Matrix transpose multiplication with updated interfaces.
657  void Tmul
658  (
659  scalarField&,
660  const tmp<scalarField>&,
663  const direction cmpt
664  ) const;
665 
666 
667  //- Sum the coefficients on each row of the matrix
668  void sumA
669  (
670  scalarField&,
673  ) const;
674 
675 
676  void residual
677  (
678  scalarField& rA,
679  const scalarField& psi,
680  const scalarField& source,
681  const FieldField<Field, scalar>& interfaceBouCoeffs,
682  const lduInterfaceFieldPtrsList& interfaces,
683  const direction cmpt
684  ) const;
685 
687  (
688  const scalarField& psi,
689  const scalarField& source,
690  const FieldField<Field, scalar>& interfaceBouCoeffs,
691  const lduInterfaceFieldPtrsList& interfaces,
692  const direction cmpt
693  ) const;
694 
695 
696  //- Initialise the update of interfaced interfaces
697  // for matrix operations
699  (
700  const FieldField<Field, scalar>& interfaceCoeffs,
701  const lduInterfaceFieldPtrsList& interfaces,
702  const scalarField& psiif,
703  scalarField& result,
704  const direction cmpt
705  ) const;
706 
707  //- Update interfaced interfaces for matrix operations
709  (
710  const FieldField<Field, scalar>& interfaceCoeffs,
711  const lduInterfaceFieldPtrsList& interfaces,
712  const scalarField& psiif,
713  scalarField& result,
714  const direction cmpt
715  ) const;
716 
717 
718  template<class Type>
719  tmp<Field<Type>> H(const Field<Type>&) const;
720 
721  template<class Type>
722  tmp<Field<Type>> H(const tmp<Field<Type>>&) const;
723 
724  tmp<scalarField> H1() const;
725 
726  template<class Type>
727  tmp<Field<Type>> faceH(const Field<Type>&) const;
728 
729  template<class Type>
730  tmp<Field<Type>> faceH(const tmp<Field<Type>>&) const;
731 
732 
733  // Info
734 
735  //- Return info proxy.
736  // Used to print matrix information to a stream
737  InfoProxy<lduMatrix> info() const
738  {
739  return *this;
740  }
741 
742 
743  // Member Operators
744 
745  void operator=(const lduMatrix&);
746 
747  void negate();
748 
749  void operator+=(const lduMatrix&);
750  void operator-=(const lduMatrix&);
751 
752  void operator*=(const scalarField&);
753  void operator*=(scalar);
754 
755  void operator/=(const scalarField&);
756  void operator/=(scalar);
757 
758 
759  // Ostream operator
760 
761  friend Ostream& operator<<(Ostream&, const lduMatrix&);
762  friend Ostream& operator<<(Ostream&, const InfoProxy<lduMatrix>&);
763 };
764 
765 
766 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
767 
768 } // End namespace Foam
769 
770 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
771 
772 #ifdef NoRepository
773  #include "lduMatrixTemplates.C"
774 #endif
775 
776 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
777 
778 #endif
779 
780 // ************************************************************************* //
tmp< Field< Type > > faceH(const Field< Type > &) const
virtual void read(const dictionary &)
Read and reset the solver parameters from the given stream.
const FieldField< Field, scalar > & interfaceIntCoeffs_
Definition: lduMatrix.H:102
bool hasUpper() const
Definition: lduMatrix.H:585
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
lduMatrix(const lduMesh &)
Construct given an LDU addressed mesh.
Definition: lduMatrix.C:43
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
const FieldField< Field, scalar > & interfaceBouCoeffs_
Definition: lduMatrix.H:101
uint8_t direction
Definition: direction.H:45
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void initMatrixInterfaces(const FieldField< Field, scalar > &interfaceCoeffs, const lduInterfaceFieldPtrsList &interfaces, const scalarField &psiif, scalarField &result, const direction cmpt) const
Initialise the update of interfaced interfaces.
virtual ~solver()
Destructor.
Definition: lduMatrix.H:213
void Tmul(scalarField &, const tmp< scalarField > &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &, const direction cmpt) const
Matrix transpose multiplication with updated interfaces.
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:59
const lduMatrix & matrix() const
Definition: lduMatrix.H:226
InfoProxy< lduMatrix > info() const
Return info proxy.
Definition: lduMatrix.H:736
bool hasLower() const
Definition: lduMatrix.H:590
const FieldField< Field, scalar > & interfaceBouCoeffs() const
Definition: lduMatrix.H:231
scalarField & upper()
Definition: lduMatrix.C:197
ClassName("lduMatrix")
tmp< scalarField > H1() const
void sumMagOffDiag(scalarField &sumOff) const
Generic field type.
Definition: FieldField.H:51
void smooth(volScalarField &field, const scalar coeff)
Definition: fvcSmooth.C:35
void Amul(scalarField &, const tmp< scalarField > &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &, const direction cmpt) const
Matrix multiplication with updated interfaces.
label minIter_
Minimum number of iterations in the solver.
Definition: lduMatrix.H:115
const FieldField< Field, scalar > & interfaceIntCoeffs() const
Definition: lduMatrix.H:236
void residual(scalarField &rA, const scalarField &psi, const scalarField &source, const FieldField< Field, scalar > &interfaceBouCoeffs, const lduInterfaceFieldPtrsList &interfaces, const direction cmpt) const
void sumA(scalarField &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &) const
Sum the coefficients on each row of the matrix.
const lduSchedule & patchSchedule() const
Return the patch evaluation schedule.
Definition: lduMatrix.H:556
Abstract base-class for lduMatrix smoothers.
Definition: lduMatrix.H:270
Forward declarations of the specialisations of Field<T> for scalar, vector and tensor.
scalar tolerance_
Final convergence tolerance.
Definition: lduMatrix.H:118
bool hasDiag() const
Definition: lduMatrix.H:580
A class for handling words, derived from string.
Definition: word.H:59
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics...
label maxIter_
Maximum number of iterations in the solver.
Definition: lduMatrix.H:112
Abstract base-class for lduMatrix solvers.
Definition: lduMatrix.H:93
lduMatrix member H operations.
const lduInterfaceFieldPtrsList & interfaces() const
Definition: lduMatrix.H:241
const volScalarField & psi
friend Ostream & operator<<(Ostream &, const lduMatrix &)
void updateMatrixInterfaces(const FieldField< Field, scalar > &interfaceCoeffs, const lduInterfaceFieldPtrsList &interfaces, const scalarField &psiif, scalarField &result, const direction cmpt) const
Update interfaced interfaces for matrix operations.
virtual const lduAddressing & lduAddr() const =0
Return ldu addressing.
solver(const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces, const dictionary &solverControls)
virtual const word & type() const =0
Runtime type information.
static autoPtr< solver > 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 solver.
void operator=(const lduMatrix &)
const lduMatrix & matrix_
Definition: lduMatrix.H:100
bool asymmetric() const
Definition: lduMatrix.H:623
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void operator*=(const scalarField &)
void operator/=(const scalarField &)
lduInterfaceFieldPtrsList interfaces_
Definition: lduMatrix.H:103
scalar normFactor(const scalarField &psi, const scalarField &source, const scalarField &Apsi, scalarField &tmpField) const
Return the matrix norm used to normalise the residual for the.
dictionary controlDict_
Dictionary of controls.
Definition: lduMatrix.H:106
lduMatrix is a general matrix class in which the coefficients are stored as three arrays...
Definition: lduMatrix.H:79
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
Abstract base-class for lduMatrix preconditioners.
Definition: lduMatrix.H:410
A helper class for outputting values to Ostream.
Definition: InfoProxy.H:45
const lduAddressing & lduAddr() const
Return the LDU addressing.
Definition: lduMatrix.H:550
scalarField & lower()
Definition: lduMatrix.C:168
bool diagonal() const
Definition: lduMatrix.H:595
Ostream & operator<<(Ostream &, const ensightPart &)
const word & fieldName() const
Definition: lduMatrix.H:221
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
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: PtrList.H:53
scalarField & diag()
Definition: lduMatrix.C:186
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:353
const lduMesh & mesh() const
Return the LDU mesh from which the addressing is obtained.
Definition: lduMatrix.H:544
bool symmetric() const
Definition: lduMatrix.H:609
Namespace for OpenFOAM.
~lduMatrix()
Destructor.
Definition: lduMatrix.C:149
virtual void readControls()
Read the control parameters from the controlDict_.
virtual const lduSchedule & patchSchedule() const =0
static const label defaultMaxIter_
Default maximum number of iterations in the solver.
Definition: lduMatrix.H:109
scalar relTol_
Convergence tolerance relative to the initial.
Definition: lduMatrix.H:121