MovingPhaseModel.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) 2015-2018 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::MovingPhaseModel
28 
29 Description
30  Class which represents a moving fluid phase. Holds the velocity, fluxes and
31  turbulence model and can generate the momentum equation. The interface is
32  quite restrictive as it also has to support an equivalent stationary model,
33  which does not store motion fields or a turbulence model.
34 
35  Possible future extensions include separating the turbulent fuctionality
36  into another layer.
37 
38 See also
39  StationaryPhaseModel
40 
41 SourceFiles
42  MovingPhaseModel.C
43 
44 \*---------------------------------------------------------------------------*/
45 
46 #ifndef MovingPhaseModel_H
47 #define MovingPhaseModel_H
48 
49 #include "phaseModel.H"
50 #include "phaseCompressibleTurbulenceModel.H"
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 namespace Foam
55 {
56 
57 /*---------------------------------------------------------------------------*\
58  Class MovingPhaseModel Declaration
59 \*---------------------------------------------------------------------------*/
60 
61 template<class BasePhaseModel>
62 class MovingPhaseModel
63 :
64  public BasePhaseModel
65 {
66 protected:
67 
68  // Protected data
69 
70  //- Velocity field
71  volVectorField U_;
72 
73  //- Flux
74  surfaceScalarField phi_;
75 
76  //- Volumetric flux
77  surfaceScalarField alphaPhi_;
78 
79  //- Mass flux
81 
82  //- Lagrangian acceleration field (needed for virtual-mass)
83  mutable tmp<volVectorField> DUDt_;
84 
85  //- Lagrangian acceleration field on the faces (needed for virtual-mass)
87 
88  //- Dilatation rate
90 
91  //- Turbulence model
93 
94  //- Continuity error due to the flow
96 
97  //- Continuity error due to any sources
99 
100  //- Kinetic Energy
101  mutable tmp<volScalarField> K_;
102 
103 
104 private:
105 
106  // Private static member functions
107 
108  //- Calculate and return the flux field
110 
111 
112 public:
113 
114  // Constructors
115 
117  (
118  const phaseSystem& fluid,
119  const word& phaseName,
120  const label index
121  );
122 
123 
124  //- Destructor
125  virtual ~MovingPhaseModel();
126 
127 
128  // Member Functions
129 
130  //- Correct the phase properties other than the thermo and turbulence
131  virtual void correct();
132 
133  //- Correct the kinematics
134  virtual void correctKinematics();
135 
136  //- Correct the turbulence
137  virtual void correctTurbulence();
138 
139  //- Correct the energy transport e.g. alphat
140  virtual void correctEnergyTransport();
141 
142 
143  // Momentum
144 
145  //- Return whether the phase is stationary
146  virtual bool stationary() const;
147 
148  //- Return the momentum equation
149  virtual tmp<fvVectorMatrix> UEqn();
150 
151  //- Return the momentum equation for the face-based algorithm
152  virtual tmp<fvVectorMatrix> UfEqn();
153 
154  //- Return the velocity
155  virtual tmp<volVectorField> U() const;
156 
157  //- Access the velocity
158  virtual volVectorField& URef();
159 
160  //- Return the volumetric flux
161  virtual tmp<surfaceScalarField> phi() const;
162 
163  //- Access the volumetric flux
164  virtual surfaceScalarField& phiRef();
165 
166  //- Return the volumetric flux of the phase
167  virtual tmp<surfaceScalarField> alphaPhi() const;
168 
169  //- Access the volumetric flux of the phase
170  virtual surfaceScalarField& alphaPhiRef();
171 
172  //- Return the mass flux of the phase
173  virtual tmp<surfaceScalarField> alphaRhoPhi() const;
174 
175  //- Access the mass flux of the phase
177 
178  //- Return the substantive acceleration
179  virtual tmp<volVectorField> DUDt() const;
180 
181  //- Return the substantive acceleration on the faces
182  virtual tmp<surfaceScalarField> DUDtf() const;
183 
184  //- Return the continuity error
185  virtual tmp<volScalarField> continuityError() const;
186 
187  //- Return the continuity error due to the flow field
189 
190  //- Return the continuity error due to any sources
192 
193  //- Return the phase kinetic energy
194  virtual tmp<volScalarField> K() const;
195 
196 
197  // Compressibility (variable density)
198 
199  //- Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
200  virtual tmp<volScalarField> divU() const;
201 
202  //- Set the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
203  virtual void divU(tmp<volScalarField> divU);
204 
205 
206  // Turbulence
207 
208  //- Return the turbulent dynamic viscosity
209  virtual tmp<volScalarField> mut() const;
210 
211  //- Return the effective dynamic viscosity
212  virtual tmp<volScalarField> muEff() const;
213 
214  //- Return the turbulent kinematic viscosity
215  virtual tmp<volScalarField> nut() const;
216 
217  //- Return the effective kinematic viscosity
218  virtual tmp<volScalarField> nuEff() const;
219 
220  //- Effective thermal turbulent diffusivity for temperature
221  // of mixture for patch [J/m/s/K]
223 
224  //- Return the effective thermal conductivity
225  virtual tmp<volScalarField> kappaEff() const;
226 
227  //- Return the effective thermal conductivity on a patch
228  virtual tmp<scalarField> kappaEff(const label patchi) const;
229 
230  //- Effective thermal turbulent diffusivity of mixture [kg/m/s]
232 
233  //- Return the effective thermal diffusivity
234  virtual tmp<volScalarField> alphaEff() const;
235 
236  //- Return the effective thermal conductivity on a patch
237  virtual tmp<scalarField> alphaEff(const label patchi) const;
238 
239  //- Return the turbulent kinetic energy
240  virtual tmp<volScalarField> k() const;
241 
242  //- Return the phase-pressure'
243  // (derivative of phase-pressure w.r.t. phase-fraction)
244  virtual tmp<volScalarField> pPrime() const;
245 };
246 
247 
248 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
249 
250 } // End namespace Foam
251 
252 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
253 
254 #ifdef NoRepository
255  #include "MovingPhaseModel.C"
256 #endif
257 
258 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
259 
260 #endif
261 
262 // ************************************************************************* //
Foam::MovingPhaseModel::U_
volVectorField U_
Velocity field.
Definition: MovingPhaseModel.H:70
Foam::MovingPhaseModel::divU_
tmp< volScalarField > divU_
Dilatation rate.
Definition: MovingPhaseModel.H:88
Foam::MovingPhaseModel::alphaRhoPhi_
surfaceScalarField alphaRhoPhi_
Mass flux.
Definition: MovingPhaseModel.H:79
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::MovingPhaseModel::kappaEff
virtual tmp< volScalarField > kappaEff() const
Return the effective thermal conductivity.
Definition: MovingPhaseModel.C:506
Foam::MovingPhaseModel::phi_
surfaceScalarField phi_
Flux.
Definition: MovingPhaseModel.H:73
Foam::MovingPhaseModel::DUDt_
tmp< volVectorField > DUDt_
Lagrangian acceleration field (needed for virtual-mass)
Definition: MovingPhaseModel.H:82
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::MovingPhaseModel::continuityErrorSources
virtual tmp< volScalarField > continuityErrorSources() const
Return the continuity error due to any sources.
Definition: MovingPhaseModel.C:434
Foam::MovingPhaseModel::alphaPhi
virtual tmp< surfaceScalarField > alphaPhi() const
Constant access the volumetric flux of the phase.
Definition: MovingPhaseModel.C:97
Foam::MovingPhaseModel::K_
tmp< volScalarField > K_
Kinetic Energy.
Definition: MovingPhaseModel.H:100
Foam::MovingPhaseModel::stationary
virtual bool stationary() const
Return whether the phase is stationary.
Definition: MovingPhaseModel.C:282
kappaEff
kappaEff
Definition: TEqn.H:10
Foam::MovingPhaseModel::continuityError
virtual tmp< volScalarField > continuityError() const
Return the continuity error.
Definition: MovingPhaseModel.C:418
fluid
twoPhaseSystem & fluid
Definition: setRegionFluidFields.H:3
Foam::MovingPhaseModel::mut
virtual tmp< volScalarField > mut() const
Return the turbulent dynamic viscosity.
Definition: MovingPhaseModel.C:474
Foam::MovingPhaseModel::~MovingPhaseModel
virtual ~MovingPhaseModel()=default
Destructor.
Definition: MovingPhaseModel.C:217
Foam::MovingPhaseModel::nuEff
virtual tmp< volScalarField > nuEff() const
Return the effective kinematic viscosity.
Definition: MovingPhaseModel.C:498
Foam::MovingPhaseModel::phiRef
virtual surfaceScalarField & phiRef()
Access the volumetric flux.
Definition: MovingPhaseModel.C:352
Foam::MovingPhaseModel::DUDtf_
tmp< surfaceScalarField > DUDtf_
Lagrangian acceleration field on the faces (needed for virtual-mass)
Definition: MovingPhaseModel.H:85
Foam::MovingPhaseModel::alphaPhiRef
virtual surfaceScalarField & alphaPhiRef()
Access the volumetric flux of the phase.
Definition: MovingPhaseModel.C:368
Foam::MovingPhaseModel::phi
virtual tmp< surfaceScalarField > phi() const
Constant access the volumetric flux.
Definition: MovingPhaseModel.C:81
Foam::MovingPhaseModel::correctTurbulence
virtual void correctTurbulence()
Correct the turbulence.
Definition: MovingPhaseModel.C:264
Foam::MovingPhaseModel::k
virtual tmp< volScalarField > k() const
Return the turbulent kinetic energy.
Definition: MovingPhaseModel.C:538
Foam::MovingPhaseModel::U
virtual tmp< volVectorField > U() const
Access const reference to U.
Definition: MovingPhaseModel.C:113
Foam::MovingPhaseModel::alphaEff
virtual tmp< volScalarField > alphaEff() const
Return the effective thermal diffusivity.
Definition: MovingPhaseModel.C:522
Foam::MovingPhaseModel::turbulence_
autoPtr< phaseCompressibleTurbulenceModel > turbulence_
Turbulence model.
Definition: MovingPhaseModel.H:91
Foam::MovingPhaseModel::continuityErrorSources_
volScalarField continuityErrorSources_
Continuity error due to any sources.
Definition: MovingPhaseModel.H:97
Foam::MovingPhaseModel::UfEqn
virtual tmp< fvVectorMatrix > UfEqn()
Return the momentum equation for the face-based algorithm.
Definition: MovingPhaseModel.C:308
Foam::MovingPhaseModel::MovingPhaseModel
MovingPhaseModel(const phaseSystem &fluid, const word &phaseName)
Definition: MovingPhaseModel.C:48
Foam::MovingPhaseModel::UEqn
virtual tmp< fvVectorMatrix > UEqn()
Return the momentum equation.
Definition: MovingPhaseModel.C:290
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::MovingPhaseModel::DUDtf
virtual tmp< surfaceScalarField > DUDtf() const
Return the substantive acceleration on the faces.
Definition: MovingPhaseModel.C:405
Foam::MovingPhaseModel::correctKinematics
virtual void correctKinematics()
Correct the kinematics.
Definition: MovingPhaseModel.C:239
Foam::MovingPhaseModel::URef
virtual volVectorField & URef()
Access the velocity.
Definition: MovingPhaseModel.C:336
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::MovingPhaseModel::K
virtual tmp< volScalarField > K() const
Return the phase kinetic energy.
Definition: MovingPhaseModel.C:442
Foam::MovingPhaseModel::correct
virtual void correct()
Correct the phase properties other than the thermo and turbulence.
Definition: MovingPhaseModel.C:73
Foam::MovingPhaseModel::alphaRhoPhi
virtual tmp< surfaceScalarField > alphaRhoPhi() const
Return the mass flux of the phase.
Definition: MovingPhaseModel.C:376
Foam::MovingPhaseModel::continuityErrorFlow_
volScalarField continuityErrorFlow_
Continuity error due to the flow.
Definition: MovingPhaseModel.H:94
Foam::MovingPhaseModel::pPrime
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure'.
Definition: MovingPhaseModel.C:546
Foam::MovingPhaseModel::divU
virtual tmp< volScalarField > divU() const
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
Definition: MovingPhaseModel.C:459
Foam::MovingPhaseModel::alphaRhoPhiRef
virtual surfaceScalarField & alphaRhoPhiRef()
Access the mass flux of the phase.
Definition: MovingPhaseModel.C:384
alphaEff
volScalarField alphaEff("alphaEff", turbulence->nu()/Pr+alphat)
Foam::MovingPhaseModel::nut
virtual tmp< volScalarField > nut() const
Return the turbulent kinematic viscosity.
Definition: MovingPhaseModel.C:490
Foam::MovingPhaseModel::continuityErrorFlow
virtual tmp< volScalarField > continuityErrorFlow() const
Return the continuity error due to the flow field.
Definition: MovingPhaseModel.C:426
Foam::MovingPhaseModel::muEff
virtual tmp< volScalarField > muEff() const
Return the effective dynamic viscosity.
Definition: MovingPhaseModel.C:482
Foam::phaseSystem
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:63
Foam::MovingPhaseModel::DUDt
virtual tmp< volVectorField > DUDt() const
Return the substantive acceleration.
Definition: MovingPhaseModel.C:392
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:53
Foam::MovingPhaseModel::correctEnergyTransport
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat.
Definition: MovingPhaseModel.C:273