DILUPreconditioner.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::DILUPreconditioner
26 
27 Description
28  Simplified diagonal-based incomplete LU preconditioner for asymmetric
29  matrices. The reciprocal of the preconditioned diagonal is calculated
30  and stored.
31 
32 SourceFiles
33  DILUPreconditioner.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef DILUPreconditioner_H
38 #define DILUPreconditioner_H
39 
40 #include "lduMatrix.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 /*---------------------------------------------------------------------------*\
48  Class DILUPreconditioner Declaration
49 \*---------------------------------------------------------------------------*/
50 
52 :
54 {
55  // Private Data
56 
57  //- The reciprocal preconditioned diagonal
58  scalarField rD_;
59 
60 
61 public:
62 
63  //- Runtime type information
64  TypeName("DILU");
65 
66 
67  // Constructors
68 
69  //- Construct from matrix components and preconditioner solver controls
71  (
72  const lduMatrix::solver&,
73  const dictionary& solverControlsUnused
74  );
75 
76 
77  //- Destructor
78  virtual ~DILUPreconditioner()
79  {}
80 
81 
82  // Member Functions
83 
84  //- Calculate the reciprocal of the preconditioned diagonal
85  static void calcReciprocalD(scalarField& rD, const lduMatrix& matrix);
86 
87  //- Return wA the preconditioned form of residual rA
88  virtual void precondition
89  (
90  scalarField& wA,
91  const scalarField& rA,
92  const direction cmpt=0
93  ) const;
94 
95  //- Return wT the transpose-matrix preconditioned form of residual rT.
96  virtual void preconditionT
97  (
98  scalarField& wT,
99  const scalarField& rT,
100  const direction cmpt=0
101  ) const;
102 };
103 
104 
105 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
106 
107 } // End namespace Foam
108 
109 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
110 
111 #endif
112 
113 // ************************************************************************* //
virtual void precondition(scalarField &wA, const scalarField &rA, const direction cmpt=0) const
Return wA the preconditioned form of residual rA.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual ~DILUPreconditioner()
Destructor.
uint8_t direction
Definition: direction.H:45
DILUPreconditioner(const lduMatrix::solver &, const dictionary &solverControlsUnused)
Construct from matrix components and preconditioner solver controls.
Abstract base-class for lduMatrix solvers.
Definition: lduMatrix.H:93
virtual void preconditionT(scalarField &wT, const scalarField &rT, const direction cmpt=0) const
Return wT the transpose-matrix preconditioned form of residual rT.
lduMatrix is a general matrix class in which the coefficients are stored as three arrays...
Definition: lduMatrix.H:79
Abstract base-class for lduMatrix preconditioners.
Definition: lduMatrix.H:410
static void calcReciprocalD(scalarField &rD, const lduMatrix &matrix)
Calculate the reciprocal of the preconditioned diagonal.
TypeName("DILU")
Runtime type information.
Namespace for OpenFOAM.
Simplified diagonal-based incomplete LU preconditioner for asymmetric matrices. The reciprocal of the...