BUGFIX: Seg-fault in multiphaseInterFoam. Author: Henrik Rusche. Merge: Hrvoje Jasak
[foam-extend-3.2.git] / src / dbns / numericFlux / fineNumericFlux.H
blob921b827afe96e6e2ec2aba36ce472cac83b295e4
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend is free software: you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation, either version 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 Class
25     fineNumericFlux
27 Description
28     Fine-level numeric flux class for density-based solvers
30 Author
31     Aleksandar Jemcov
32     Rewrite by Hrvoje Jasak
34 SourceFiles
35     fineNumericFlux.H
36     fineNumericFlux.C
38 \*---------------------------------------------------------------------------*/
40 #ifndef fineNumericFlux_H
41 #define fineNumericFlux_H
43 #include "fvMesh.H"
44 #include "volFields.H"
45 #include "surfaceFields.H"
46 #include "basicThermo.H"
48 #include "mgMeshLevel.H"
49 #include "mgFieldLevel.H"
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 namespace Foam
56 /*---------------------------------------------------------------------------*\
57                        Class fineNumericFlux Declaration
58 \*---------------------------------------------------------------------------*/
60 template<class Flux, class Limiter>
61 class fineNumericFlux
63     public numericFluxBase<Flux>
65     // Private data
67         //- Reference to fine mesh level
68         const mgMeshLevel& meshLevel_;
70         //- Reference to fine fields level
71         const mgFieldLevel& fieldLevel_;
74         // Reference to primitive fields
76             //- Static pressure
77             const volScalarField& p_;
79             //- Velocity
80             const volVectorField& U_;
82             //- Static temperature
83             const volScalarField& T_;
85             //- Reference to the thermophysicalModel
86             basicThermo& thermo_;
89         // Fluxes
91             //- Density flux
92             surfaceScalarField& rhoFlux_;
94             //- Velocity flux
95             surfaceVectorField& rhoUFlux_;
97             //- Energy flux
98             surfaceScalarField& rhoEFlux_;
101         // Gradients
103             //- Static pressure gradient
104             volVectorField gradP_;
106             //- Velocity gradient
107             volTensorField gradU_;
109             //- Static temperature gradient
110             volVectorField gradT_;
113     // Private Member Functions
115         //- Disallow default bitwise copy construct
116         fineNumericFlux(const fineNumericFlux&);
118         //- Disallow default bitwise assignment
119         void operator=(const fineNumericFlux&);
122         //- Return internal field of mass flux
123         const scalarField& rhoFluxI() const
124         {
125             return rhoFlux_.internalField();
126         }
128         //- Return access to internal field of mass flux
129         scalarField& rhoFluxI()
130         {
131             return rhoFlux_.internalField();
132         }
134         //- Return internal field of momentum flux
135         const vectorField& rhoUFluxI() const
136         {
137             return rhoUFlux_.internalField();
138         }
140         //- Return access to internal field of momentum flux
141         vectorField& rhoUFluxI()
142         {
143             return rhoUFlux_.internalField();
144         }
146          //- Return access to internal field of energy flux
147         const scalarField& rhoEFluxI() const
148         {
149             return rhoEFlux_.internalField();
150         }
152         //- Return access to internal field of energy flux
153         scalarField& rhoEFluxI()
154         {
155             return rhoEFlux_.internalField();
156         }
159 public:
161     // Constructors
163         //- Construct from components
164         fineNumericFlux
165         (
166             const mgMeshLevel& meshLevel,
167             const mgFieldLevel& fieldLevel,
168             basicThermo& thermo
169         );
173     //- Destructor
174     virtual ~fineNumericFlux()
175     {}
178     // Member Functions
180         //- Return mesh reference
181         const fvMesh& mesh() const
182         {
183             return meshLevel_.mesh();
184         }
187         // Return fluxes
189             //- Return density flux
190             virtual const surfaceScalarField& rhoFlux() const
191             {
192                 return rhoFlux_;
193             }
195             //- Return velocity flux
196             virtual const surfaceVectorField& rhoUFlux() const
197             {
198                 return rhoUFlux_;
199             }
201             //- Return energy flux
202             virtual const surfaceScalarField& rhoEFlux() const
203             {
204                 return rhoEFlux_;
205             }
208        // Return residuals
210             //- Return density equation residual
211             virtual tmp<scalarField> rhoResidual() const
212             {
213                 return fvc::div(rhoFlux_)().internalField();
214             }
216             //- Return momentum equation flux
217             virtual tmp<vectorField> rhoUResidual() const
218             {
219                 return fvc::div(rhoUFlux_)().internalField();
220             }
222             //- Return energy equation flux
223             virtual tmp<scalarField> rhoEResidual() const
224             {
225                 return fvc::div(rhoEFlux_)().internalField();
226             }
229         // Return Gradients
231             //- Return pressure gradient
232             const volVectorField& gradP() const
233             {
234                 return gradP_;
235             }
237             //- Return pressure gradient
238             const volTensorField& gradU() const
239             {
240                 return gradU_;
241             }
243             //- Return Temperature gradient
244             const volVectorField& gradT() const
245             {
246                 return gradT_;
247             }
250         // Update fluxes based on current state
252             //- Compute flux
253             virtual void computeFlux();
257 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
259 } // End namespace Foam
261 #ifdef NoRepository
262 #   include "fineNumericFlux.C"
263 #endif
265 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
267 #endif
269 // ************************************************************************* //