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