lduMatrix.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2016-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  It might be better if this class were organised as a hierarchy starting
38  from an empty matrix, then deriving diagonal, symmetric and asymmetric
39  matrices.
40 
41 SourceFiles
42  lduMatrixATmul.C
43  lduMatrix.C
44  lduMatrixTemplates.C
45  lduMatrixOperations.C
46  lduMatrixSolver.C
47  lduMatrixPreconditioner.C
48  lduMatrixTests.C
49  lduMatrixUpdateMatrixInterfaces.C
50 
51 \*---------------------------------------------------------------------------*/
52 
53 #ifndef lduMatrix_H
54 #define lduMatrix_H
55 
56 #include "lduMesh.H"
57 #include "primitiveFieldsFwd.H"
58 #include "FieldField.H"
60 #include "typeInfo.H"
61 #include "autoPtr.H"
62 #include "runTimeSelectionTables.H"
63 #include "solverPerformance.H"
64 #include "InfoProxy.H"
65 #include "profilingTrigger.H"
66 
67 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68 
69 namespace Foam
70 {
71 
72 // Forward declaration of friend functions and operators
73 
74 class lduMatrix;
75 
76 Ostream& operator<<(Ostream&, const lduMatrix&);
77 Ostream& operator<<(Ostream&, const InfoProxy<lduMatrix>&);
78 
79 
80 /*---------------------------------------------------------------------------*\
81  Class lduMatrix Declaration
82 \*---------------------------------------------------------------------------*/
83 
84 class lduMatrix
85 {
86  // private data
87 
88  //- LDU mesh reference
89  const lduMesh& lduMesh_;
90 
91  //- Coefficients (not including interfaces)
92  scalarField *lowerPtr_, *diagPtr_, *upperPtr_;
93 
94 
95 public:
96 
97  //- Abstract base-class for lduMatrix solvers
98  class solver
99  {
100  protected:
101 
102  // Protected data
103 
109 
110  //- Dictionary of controls
112 
113  //- Default maximum number of iterations in the solver
114  static const label defaultMaxIter_;
115 
116  //- Maximum number of iterations in the solver
117  label maxIter_;
118 
119  //- Minimum number of iterations in the solver
120  label minIter_;
121 
122  //- Final convergence tolerance
123  scalar tolerance_;
124 
125  //- Convergence tolerance relative to the initial
126  scalar relTol_;
127 
129 
130 
131  // Protected Member Functions
132 
133  //- Read the control parameters from the controlDict_
134  virtual void readControls();
135 
136 
137  public:
138 
139  //- Runtime type information
140  virtual const word& type() const = 0;
141 
142 
143  // Declare run-time constructor selection tables
144 
146  (
147  autoPtr,
148  solver,
149  symMatrix,
150  (
151  const word& fieldName,
152  const lduMatrix& matrix,
156  const dictionary& solverControls
157  ),
158  (
159  fieldName,
160  matrix,
163  interfaces,
164  solverControls
165  )
166  );
167 
169  (
170  autoPtr,
171  solver,
172  asymMatrix,
173  (
174  const word& fieldName,
175  const lduMatrix& matrix,
179  const dictionary& solverControls
180  ),
181  (
182  fieldName,
183  matrix,
186  interfaces,
187  solverControls
188  )
189  );
190 
191 
192  // Constructors
193 
194  solver
195  (
196  const word& fieldName,
197  const lduMatrix& matrix,
201  const dictionary& solverControls
202  );
203 
204  // Selectors
205 
206  //- Return a new solver
207  static autoPtr<solver> New
208  (
209  const word& fieldName,
210  const lduMatrix& matrix,
214  const dictionary& solverControls
215  );
216 
217 
218 
219  //- Destructor
220  virtual ~solver() = default;
221 
222 
223  // Member functions
224 
225  // Access
226 
227  const word& fieldName() const
228  {
229  return fieldName_;
230  }
231 
232  const lduMatrix& matrix() const
233  {
234  return matrix_;
235  }
236 
238  {
239  return interfaceBouCoeffs_;
240  }
241 
243  {
244  return interfaceIntCoeffs_;
245  }
246 
248  {
249  return interfaces_;
250  }
251 
252 
253  //- Read and reset the solver parameters from the given stream
254  virtual void read(const dictionary&);
255 
256  //- Solve with given field and rhs
257  virtual solverPerformance solve
258  (
259  scalarField& psi,
260  const scalarField& source,
261  const direction cmpt=0
262  ) const = 0;
263 
264  //- Solve with given field and rhs (in solveScalar precision).
265  // Default is to call solve routine
267  (
269  const solveScalarField& source,
270  const direction cmpt=0
271  ) const;
272 
273  //- Return the matrix norm used to normalise the residual for the
274  //- stopping criterion
276  (
277  const solveScalarField& psi,
278  const solveScalarField& source,
279  const solveScalarField& Apsi,
280  solveScalarField& tmpField
281  ) const;
282  };
283 
284 
285  //- Abstract base-class for lduMatrix smoothers
286  class smoother
287  {
288  protected:
289 
290  // Protected data
291 
297 
298 
299  public:
300 
301  //- Find the smoother name (directly or from a sub-dictionary)
302  static word getName(const dictionary&);
303 
304  //- Runtime type information
305  virtual const word& type() const = 0;
306 
307 
308  // Declare run-time constructor selection tables
309 
311  (
312  autoPtr,
313  smoother,
314  symMatrix,
315  (
316  const word& fieldName,
317  const lduMatrix& matrix,
321  ),
322  (
323  fieldName,
324  matrix,
327  interfaces
328  )
329  );
330 
332  (
333  autoPtr,
334  smoother,
335  asymMatrix,
336  (
337  const word& fieldName,
338  const lduMatrix& matrix,
342  ),
343  (
344  fieldName,
345  matrix,
348  interfaces
349  )
350  );
351 
352 
353  // Constructors
354 
355  smoother
356  (
357  const word& fieldName,
358  const lduMatrix& matrix,
362  );
363 
364 
365  // Selectors
366 
367  //- Return a new smoother
368  static autoPtr<smoother> New
369  (
370  const word& fieldName,
371  const lduMatrix& matrix,
375  const dictionary& solverControls
376  );
377 
378 
379  //- Destructor
380  virtual ~smoother() = default;
381 
382 
383  // Member functions
384 
385  // Access
386 
387  const word& fieldName() const
388  {
389  return fieldName_;
390  }
391 
392  const lduMatrix& matrix() const
393  {
394  return matrix_;
395  }
396 
398  {
399  return interfaceBouCoeffs_;
400  }
401 
403  {
404  return interfaceIntCoeffs_;
405  }
406 
408  {
409  return interfaces_;
410  }
411 
412 
413  //- Smooth the solution for a given number of sweeps
414  virtual void smooth
415  (
417  const scalarField& source,
418  const direction cmpt,
419  const label nSweeps
420  ) const = 0;
421 
422  //- Smooth the solution for a given number of sweeps
423  virtual void scalarSmooth
424  (
426  const solveScalarField& source,
427  const direction cmpt,
428  const label nSweeps
429  ) const = 0;
430  };
431 
432 
433  //- Abstract base-class for lduMatrix preconditioners
434  class preconditioner
435  {
436  protected:
437 
438  // Protected data
439 
440  //- Reference to the base-solver this preconditioner is used with
441  const solver& solver_;
442 
443 
444  public:
445 
446  //- Find the preconditioner name (directly or from a sub-dictionary)
447  static word getName(const dictionary&);
448 
449  //- Runtime type information
450  virtual const word& type() const = 0;
451 
452 
453  // Declare run-time constructor selection tables
454 
456  (
457  autoPtr,
459  symMatrix,
460  (
461  const solver& sol,
462  const dictionary& solverControls
463  ),
464  (sol, solverControls)
465  );
466 
468  (
469  autoPtr,
471  asymMatrix,
472  (
473  const solver& sol,
474  const dictionary& solverControls
475  ),
476  (sol, solverControls)
477  );
478 
479 
480  // Constructors
481 
483  (
484  const solver& sol
485  )
486  :
487  solver_(sol)
488  {}
489 
490 
491  // Selectors
492 
493  //- Return a new preconditioner
495  (
496  const solver& sol,
497  const dictionary& solverControls
498  );
499 
500 
501  //- Destructor
502  virtual ~preconditioner() = default;
503 
504 
505  // Member functions
506 
507  //- Read and reset the preconditioner parameters
508  //- from the given stream
509  virtual void read(const dictionary&)
510  {}
511 
512  //- Return wA the preconditioned form of residual rA
513  virtual void precondition
514  (
515  solveScalarField& wA,
516  const solveScalarField& rA,
517  const direction cmpt=0
518  ) const = 0;
519 
520  //- Return wT the transpose-matrix preconditioned form of
521  //- residual rT.
522  // This is only required for preconditioning asymmetric matrices.
523  virtual void preconditionT
524  (
525  solveScalarField& wT,
526  const solveScalarField& rT,
527  const direction cmpt=0
528  ) const
529  {
531  }
532  };
533 
534 
535  // Static data
536 
537  // Declare name of the class and its debug switch
538  ClassName("lduMatrix");
539 
540 
541  // Constructors
542 
543  //- Construct given an LDU addressed mesh.
544  // The coefficients are initially empty for subsequent setting.
545  lduMatrix(const lduMesh&);
546 
547  //- Construct as copy
548  lduMatrix(const lduMatrix&);
549 
550  //- Construct as copy or re-use as specified.
551  lduMatrix(lduMatrix&, bool reuse);
552 
553  //- Construct given an LDU addressed mesh and an Istream
554  //- from which the coefficients are read
555  lduMatrix(const lduMesh&, Istream&);
556 
557 
558  //- Destructor
559  ~lduMatrix();
560 
561 
562  // Member functions
563 
564  // Access to addressing
565 
566  //- Return the LDU mesh from which the addressing is obtained
567  const lduMesh& mesh() const
568  {
569  return lduMesh_;
570  }
571 
572  //- Return the LDU addressing
573  const lduAddressing& lduAddr() const
574  {
575  return lduMesh_.lduAddr();
576  }
577 
578  //- Return the patch evaluation schedule
579  const lduSchedule& patchSchedule() const
580  {
581  return lduAddr().patchSchedule();
582  }
583 
584 
585  // Access to coefficients
586 
587  scalarField& lower();
588  scalarField& diag();
589  scalarField& upper();
590 
591  // Size with externally provided sizes (for constructing with 'fake'
592  // mesh in GAMG)
593 
594  scalarField& lower(const label size);
595  scalarField& diag(const label nCoeffs);
596  scalarField& upper(const label nCoeffs);
597 
598 
599  const scalarField& lower() const;
600  const scalarField& diag() const;
601  const scalarField& upper() const;
602 
603  bool hasDiag() const
604  {
605  return (diagPtr_);
606  }
607 
608  bool hasUpper() const
609  {
610  return (upperPtr_);
611  }
612 
613  bool hasLower() const
614  {
615  return (lowerPtr_);
616  }
617 
618  bool diagonal() const
619  {
620  return (diagPtr_ && !lowerPtr_ && !upperPtr_);
621  }
622 
623  bool symmetric() const
624  {
625  return (diagPtr_ && (!lowerPtr_ && upperPtr_));
626  }
627 
628  bool asymmetric() const
629  {
630  return (diagPtr_ && lowerPtr_ && upperPtr_);
631  }
632 
633 
634  // operations
635 
636  void sumDiag();
637  void negSumDiag();
638 
639  void sumMagOffDiag(scalarField& sumOff) const;
640 
641  //- Matrix multiplication with updated interfaces.
642  void Amul
643  (
645  const tmp<solveScalarField>&,
648  const direction cmpt
649  ) const;
650 
651  //- Matrix transpose multiplication with updated interfaces.
652  void Tmul
653  (
655  const tmp<solveScalarField>&,
658  const direction cmpt
659  ) const;
660 
661 
662  //- Sum the coefficients on each row of the matrix
663  void sumA
664  (
668  ) const;
669 
670 
671  void residual
672  (
673  solveScalarField& rA,
674  const solveScalarField& psi,
675  const scalarField& source,
676  const FieldField<Field, scalar>& interfaceBouCoeffs,
677  const lduInterfaceFieldPtrsList& interfaces,
678  const direction cmpt
679  ) const;
680 
682  (
683  const solveScalarField& psi,
684  const scalarField& source,
685  const FieldField<Field, scalar>& interfaceBouCoeffs,
686  const lduInterfaceFieldPtrsList& interfaces,
687  const direction cmpt
688  ) const;
689 
690 
691  //- Initialise the update of interfaced interfaces
692  //- for matrix operations
694  (
695  const bool add,
696  const FieldField<Field, scalar>& interfaceCoeffs,
697  const lduInterfaceFieldPtrsList& interfaces,
698  const solveScalarField& psiif,
699  solveScalarField& result,
700  const direction cmpt
701  ) const;
702 
703  //- Update interfaced interfaces for matrix operations
705  (
706  const bool add,
707  const FieldField<Field, scalar>& interfaceCoeffs,
708  const lduInterfaceFieldPtrsList& interfaces,
709  const solveScalarField& psiif,
710  solveScalarField& result,
711  const direction cmpt,
712  const label startRequest // starting request (for non-blocking)
713  ) const;
714 
715  //- Set the residual field using an IOField on the object registry
716  //- if it exists
717  void setResidualField
718  (
719  const scalarField& residual,
720  const word& fieldName,
721  const bool initial
722  ) const;
723 
724  template<class Type>
725  tmp<Field<Type>> H(const Field<Type>&) const;
726 
727  template<class Type>
728  tmp<Field<Type>> H(const tmp<Field<Type>>&) const;
729 
730  tmp<scalarField> H1() const;
731 
732  template<class Type>
733  tmp<Field<Type>> faceH(const Field<Type>&) const;
734 
735  template<class Type>
736  tmp<Field<Type>> faceH(const tmp<Field<Type>>&) const;
737 
738 
739  // Info
740 
741  //- Return info proxy.
742  // Used to print matrix information to a stream
743  InfoProxy<lduMatrix> info() const
744  {
745  return *this;
746  }
747 
748 
749  // Member operators
750 
751  void operator=(const lduMatrix&);
752 
753  void negate();
754 
755  void operator+=(const lduMatrix&);
756  void operator-=(const lduMatrix&);
757 
758  void operator*=(const scalarField&);
759  void operator*=(scalar);
760 
761 
762  // Ostream operator
763 
764  friend Ostream& operator<<(Ostream&, const lduMatrix&);
766 };
767 
768 
769 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
770 
771 } // End namespace Foam
772 
773 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
774 
775 #ifdef NoRepository
776  #include "lduMatrixTemplates.C"
777 #endif
778 
779 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
780 
781 #endif
782 
783 // ************************************************************************* //
Foam::lduMatrix::solver::solve
virtual solverPerformance solve(scalarField &psi, const scalarField &source, const direction cmpt=0) const =0
Solve with given field and rhs.
Foam::lduAddressing
The class contains the addressing required by the lduMatrix: upper, lower and losort.
Definition: lduAddressing.H:114
Foam::lduMatrix::preconditioner::getName
static word getName(const dictionary &)
Find the preconditioner name (directly or from a sub-dictionary)
Definition: lduMatrixPreconditioner.C:43
Foam::lduMatrix::operator-=
void operator-=(const lduMatrix &)
Definition: lduMatrixOperations.C:223
Foam::lduMatrix::lduAddr
const lduAddressing & lduAddr() const
Return the LDU addressing.
Definition: lduMatrix.H:572
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
primitiveFieldsFwd.H
Forward declarations of the specialisations of Field<T> for scalar, vector and tensor.
Foam::lduMatrix::solver::interfaces
const lduInterfaceFieldPtrsList & interfaces() const
Definition: lduMatrix.H:246
Foam::lduMatrix::info
InfoProxy< lduMatrix > info() const
Return info proxy.
Definition: lduMatrix.H:742
Foam::lduMatrix::hasLower
bool hasLower() const
Definition: lduMatrix.H:612
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::InfoProxy
A helper class for outputting values to Ostream.
Definition: InfoProxy.H:47
Foam::FieldField
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:53
typeInfo.H
FieldField.H
Foam::lduMatrix::smoother::scalarSmooth
virtual void scalarSmooth(solveScalarField &psi, const solveScalarField &source, const direction cmpt, const label nSweeps) const =0
Smooth the solution for a given number of sweeps.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::lduMatrix::smoother::declareRunTimeSelectionTable
declareRunTimeSelectionTable(autoPtr, smoother, symMatrix,(const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces),(fieldName, matrix, interfaceBouCoeffs, interfaceIntCoeffs, interfaces))
Foam::lduMatrix::operator+=
void operator+=(const lduMatrix &)
Definition: lduMatrixOperations.C:144
Foam::lduMatrix::~lduMatrix
~lduMatrix()
Destructor.
Definition: lduMatrix.C:155
Foam::lduMatrix::preconditioner::read
virtual void read(const dictionary &)
Definition: lduMatrix.H:508
Foam::lduMatrix::ClassName
ClassName("lduMatrix")
InfoProxy.H
Foam::lduMatrix::upper
scalarField & upper()
Definition: lduMatrix.C:203
Foam::lduMatrix::solver::profiling_
profilingTrigger profiling_
Definition: lduMatrix.H:127
Foam::lduMatrix
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition: lduMatrix.H:83
Foam::lduMatrix::sumMagOffDiag
void sumMagOffDiag(scalarField &sumOff) const
Definition: lduMatrixOperations.C:71
Foam::lduMatrix::hasDiag
bool hasDiag() const
Definition: lduMatrix.H:602
Foam::lduMatrix::smoother::interfaceBouCoeffs
const FieldField< Field, scalar > & interfaceBouCoeffs() const
Definition: lduMatrix.H:396
Foam::lduMatrix::symmetric
bool symmetric() const
Definition: lduMatrix.H:622
Foam::lduMatrix::solver::scalarSolve
virtual solverPerformance scalarSolve(solveScalarField &psi, const solveScalarField &source, const direction cmpt=0) const
Solve with given field and rhs (in solveScalar precision).
Definition: lduMatrixSolver.C:180
Foam::lduMatrix::preconditioner::New
static autoPtr< preconditioner > New(const solver &sol, const dictionary &solverControls)
Return a new preconditioner.
Definition: lduMatrixPreconditioner.C:68
Foam::lduMatrix::solver::matrix
const lduMatrix & matrix() const
Definition: lduMatrix.H:231
Foam::lduMatrix::solver
Abstract base-class for lduMatrix solvers.
Definition: lduMatrix.H:97
lduMatrixTemplates.C
lduMatrix member H operations.
Foam::lduMatrix::smoother::fieldName
const word & fieldName() const
Definition: lduMatrix.H:386
Foam::lduMatrix::solver::normFactor
solveScalarField::cmptType normFactor(const solveScalarField &psi, const solveScalarField &source, const solveScalarField &Apsi, solveScalarField &tmpField) const
Definition: lduMatrixSolver.C:197
Foam::lduMatrix::updateMatrixInterfaces
void updateMatrixInterfaces(const bool add, const FieldField< Field, scalar > &interfaceCoeffs, const lduInterfaceFieldPtrsList &interfaces, const solveScalarField &psiif, solveScalarField &result, const direction cmpt, const label startRequest) const
Update interfaced interfaces for matrix operations.
Definition: lduMatrixUpdateMatrixInterfaces.C:103
Foam::lduMatrix::smoother::matrix_
const lduMatrix & matrix_
Definition: lduMatrix.H:292
lduInterfaceFieldPtrsList.H
Foam::lduMatrix::smoother::type
virtual const word & type() const =0
Runtime type information.
Foam::lduMatrix::Amul
void Amul(solveScalarField &, const tmp< solveScalarField > &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &, const direction cmpt) const
Matrix multiplication with updated interfaces.
Definition: lduMatrixATmul.C:38
Foam::lduMatrix::lduMatrix
lduMatrix(const lduMesh &)
Construct given an LDU addressed mesh.
Definition: lduMatrix.C:49
Foam::lduMatrix::smoother::interfaceIntCoeffs
const FieldField< Field, scalar > & interfaceIntCoeffs() const
Definition: lduMatrix.H:401
Foam::lduMatrix::solver::interfaceBouCoeffs_
const FieldField< Field, scalar > & interfaceBouCoeffs_
Definition: lduMatrix.H:105
Foam::lduMatrix::negSumDiag
void negSumDiag()
Definition: lduMatrixOperations.C:53
Foam::lduMatrix::solver::minIter_
label minIter_
Minimum number of iterations in the solver.
Definition: lduMatrix.H:119
Foam::operator<<
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:445
Foam::Field< solveScalar >::cmptType
pTraits< solveScalar >::cmptType cmptType
Component type.
Definition: Field.H:86
Foam::lduMatrix::smoother::interfaceBouCoeffs_
const FieldField< Field, scalar > & interfaceBouCoeffs_
Definition: lduMatrix.H:293
Foam::lduMatrix::setResidualField
void setResidualField(const scalarField &residual, const word &fieldName, const bool initial) const
Definition: lduMatrix.C:322
Foam::lduMatrix::solver::read
virtual void read(const dictionary &)
Read and reset the solver parameters from the given stream.
Definition: lduMatrixSolver.C:172
Foam::Field< scalar >
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::lduMatrix::smoother::matrix
const lduMatrix & matrix() const
Definition: lduMatrix.H:391
Foam::lduMatrix::solver::maxIter_
label maxIter_
Maximum number of iterations in the solver.
Definition: lduMatrix.H:116
Foam::UPtrList< const lduInterfaceField >
solverPerformance.H
Foam::lduAddressing::patchSchedule
virtual const lduSchedule & patchSchedule() const =0
Foam::lduMatrix::diagonal
bool diagonal() const
Definition: lduMatrix.H:617
Foam::lduMatrix::preconditioner::precondition
virtual void precondition(solveScalarField &wA, const solveScalarField &rA, const direction cmpt=0) const =0
Return wA the preconditioned form of residual rA.
Foam::lduMatrix::solver::fieldName_
word fieldName_
Definition: lduMatrix.H:103
Foam::lduMatrix::H1
tmp< scalarField > H1() const
Definition: lduMatrixATmul.C:306
Foam::lduMatrix::preconditioner::~preconditioner
virtual ~preconditioner()=default
Destructor.
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::readControls
virtual void readControls()
Read the control parameters from the controlDict_.
Definition: lduMatrixSolver.C:163
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::smoother
Abstract base-class for lduMatrix smoothers.
Definition: lduMatrix.H:285
Foam::lduMatrix::smoother::getName
static word getName(const dictionary &)
Find the smoother name (directly or from a sub-dictionary)
Definition: lduMatrixSmoother.C:43
Foam::lduMatrix::solver::New
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.
Definition: lduMatrixSolver.C:45
Foam::lduMatrix::negate
void negate()
Definition: lduMatrixOperations.C:125
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::lduMatrix::H
tmp< Field< Type > > H(const Field< Type > &) const
Foam::lduMatrix::solver::relTol_
scalar relTol_
Convergence tolerance relative to the initial.
Definition: lduMatrix.H:125
Foam::lduMatrix::smoother::New
static autoPtr< smoother > New(const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces, const dictionary &solverControls)
Return a new smoother.
Definition: lduMatrixSmoother.C:67
Foam::lduMatrix::diag
scalarField & diag()
Definition: lduMatrix.C:192
Foam::lduMatrix::solver::matrix_
const lduMatrix & matrix_
Definition: lduMatrix.H:104
Foam::lduMatrix::solver::interfaceIntCoeffs_
const FieldField< Field, scalar > & interfaceIntCoeffs_
Definition: lduMatrix.H:106
Foam::lduMatrix::initMatrixInterfaces
void initMatrixInterfaces(const bool add, const FieldField< Field, scalar > &interfaceCoeffs, const lduInterfaceFieldPtrsList &interfaces, const solveScalarField &psiif, solveScalarField &result, const direction cmpt) const
Definition: lduMatrixUpdateMatrixInterfaces.C:34
Foam::lduMatrix::sumA
void sumA(solveScalarField &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &) const
Sum the coefficients on each row of the matrix.
Definition: lduMatrixATmul.C:168
Foam::lduMatrix::patchSchedule
const lduSchedule & patchSchedule() const
Return the patch evaluation schedule.
Definition: lduMatrix.H:578
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::lduMatrix::solver::controlDict_
dictionary controlDict_
Dictionary of controls.
Definition: lduMatrix.H:110
Foam::lduMatrix::solver::tolerance_
scalar tolerance_
Final convergence tolerance.
Definition: lduMatrix.H:122
Foam::lduMatrix::solver::interfaces_
lduInterfaceFieldPtrsList interfaces_
Definition: lduMatrix.H:107
Foam::profilingTrigger
Triggers for starting/stopping code profiling.
Definition: profilingTrigger.H:54
Foam::lduMatrix::operator=
void operator=(const lduMatrix &)
Definition: lduMatrixOperations.C:91
Foam::lduMatrix::sumDiag
void sumDiag()
Definition: lduMatrixOperations.C:36
Foam::lduMatrix::smoother::interfaces_
const lduInterfaceFieldPtrsList & interfaces_
Definition: lduMatrix.H:295
runTimeSelectionTables.H
Macros to ease declaration of run-time selection tables.
Foam::lduMatrix::smoother::interfaces
const lduInterfaceFieldPtrsList & interfaces() const
Definition: lduMatrix.H:406
Foam::lduMatrix::smoother::smoother
smoother(const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces)
Definition: lduMatrixSmoother.C:161
Foam::List< lduScheduleEntry >
Foam::lduMatrix::solver::declareRunTimeSelectionTable
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))
Foam::lduMatrix::solver::interfaceBouCoeffs
const FieldField< Field, scalar > & interfaceBouCoeffs() const
Definition: lduMatrix.H:236
Foam::lduMatrix::operator*=
void operator*=(const scalarField &)
Definition: lduMatrixOperations.C:302
Foam::lduMatrix::faceH
tmp< Field< Type > > faceH(const Field< Type > &) const
Foam::lduMatrix::operator<<
friend Ostream & operator<<(Ostream &, const lduMatrix &)
Foam::lduMatrix::mesh
const lduMesh & mesh() const
Return the LDU mesh from which the addressing is obtained.
Definition: lduMatrix.H:566
Foam::direction
uint8_t direction
Definition: direction.H:52
Foam::lduMatrix::preconditioner::preconditioner
preconditioner(const solver &sol)
Definition: lduMatrix.H:482
Foam::lduMatrix::solver::fieldName
const word & fieldName() const
Definition: lduMatrix.H:226
Foam::lduMatrix::smoother::fieldName_
word fieldName_
Definition: lduMatrix.H:291
Foam::lduMatrix::asymmetric
bool asymmetric() const
Definition: lduMatrix.H:627
Foam::lduMatrix::hasUpper
bool hasUpper() const
Definition: lduMatrix.H:607
Foam::lduMatrix::smoother::smooth
virtual void smooth(solveScalarField &psi, const scalarField &source, const direction cmpt, const label nSweeps) const =0
Smooth the solution for a given number of sweeps.
Foam::lduMatrix::preconditioner::type
virtual const word & type() const =0
Runtime type information.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::lduMatrix::preconditioner
Abstract base-class for lduMatrix preconditioners.
Definition: lduMatrix.H:433
Foam::lduMatrix::preconditioner::solver_
const solver & solver_
Reference to the base-solver this preconditioner is used with.
Definition: lduMatrix.H:440
psi
const volScalarField & psi
Definition: createFieldRefs.H:1
lduMesh.H
Foam::lduMatrix::Tmul
void Tmul(solveScalarField &, const tmp< solveScalarField > &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &, const direction cmpt) const
Matrix transpose multiplication with updated interfaces.
Definition: lduMatrixATmul.C:104
Foam::lduMatrix::residual
void residual(solveScalarField &rA, const solveScalarField &psi, const scalarField &source, const FieldField< Field, scalar > &interfaceBouCoeffs, const lduInterfaceFieldPtrsList &interfaces, const direction cmpt) const
Definition: lduMatrixATmul.C:217
Foam::lduMatrix::solver::defaultMaxIter_
static const label defaultMaxIter_
Default maximum number of iterations in the solver.
Definition: lduMatrix.H:113
Foam::lduMatrix::lower
scalarField & lower()
Definition: lduMatrix.C:174
Foam::lduMatrix::solver::interfaceIntCoeffs
const FieldField< Field, scalar > & interfaceIntCoeffs() const
Definition: lduMatrix.H:241
Foam::SolverPerformance
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
Definition: SolverPerformance.H:51
profilingTrigger.H
Foam::lduMatrix::solver::type
virtual const word & type() const =0
Runtime type information.
Foam::lduMatrix::preconditioner::preconditionT
virtual void preconditionT(solveScalarField &wT, const solveScalarField &rT, const direction cmpt=0) const
Definition: lduMatrix.H:523
Foam::lduMesh::lduAddr
virtual const lduAddressing & lduAddr() const =0
Return ldu addressing.
Foam::lduMatrix::smoother::~smoother
virtual ~smoother()=default
Destructor.
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::solver
solver(const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces, const dictionary &solverControls)
Definition: lduMatrixSolver.C:140
Foam::lduMatrix::preconditioner::declareRunTimeSelectionTable
declareRunTimeSelectionTable(autoPtr, preconditioner, symMatrix,(const solver &sol, const dictionary &solverControls),(sol, solverControls))
Foam::lduMatrix::smoother::interfaceIntCoeffs_
const FieldField< Field, scalar > & interfaceIntCoeffs_
Definition: lduMatrix.H:294
Foam::lduMatrix::solver::~solver
virtual ~solver()=default
Destructor.
autoPtr.H