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  Copyright (C) 2021 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::MovingPhaseModel
29 
30 Description
31  Class which represents a moving fluid phase. Holds the velocity, fluxes and
32  turbulence model and can generate the momentum equation. The interface is
33  quite restrictive as it also has to support an equivalent stationary model,
34  which does not store motion fields or a turbulence model.
35 
36  Possible future extensions include separating the turbulent fuctionality
37  into another layer.
38 
39 See also
40  StationaryPhaseModel
41 
42 SourceFiles
43  MovingPhaseModel.C
44 
45 \*---------------------------------------------------------------------------*/
46 
47 #ifndef MovingPhaseModel_H
48 #define MovingPhaseModel_H
49 
50 #include "phaseModel.H"
51 #include "phaseCompressibleTurbulenceModel.H"
52 
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55 namespace Foam
56 {
57 
58 /*---------------------------------------------------------------------------*\
59  Class MovingPhaseModel Declaration
60 \*---------------------------------------------------------------------------*/
61 
62 template<class BasePhaseModel>
63 class MovingPhaseModel
64 :
65  public BasePhaseModel
66 {
67 protected:
68 
69  // Protected Data
70 
71  //- Velocity field
72  volVectorField U_;
73 
74  //- Flux
75  surfaceScalarField phi_;
76 
77  //- Volumetric flux
78  surfaceScalarField alphaPhi_;
79 
80  //- Mass flux
82 
83  //- Lagrangian acceleration field (needed for virtual-mass)
85 
86  //- Lagrangian acceleration field on the faces (needed for virtual-mass)
88 
89  //- Dilatation rate
91 
92  //- Turbulence model
94 
95  //- Continuity error due to the flow
97 
98  //- Continuity error due to any sources
100 
101  //- Kinetic Energy
102  mutable tmp<volScalarField> K_;
103 
105 private:
106 
107  // Private Member Functions
108 
109  //- Calculate and return the flux field
111 
112 
113 public:
115  // Constructors
116 
117  //- Construct from phase system and phase name
119  (
120  const phaseSystem& fluid,
121  const word& phaseName,
122  const label index
123  );
124 
125 
126  //- Destructor
127  virtual ~MovingPhaseModel() = default;
128 
129 
130  // Member Functions
131 
132  //- Correct the phase properties other than the thermo and turbulence
133  virtual void correct();
134 
135  //- Correct the kinematics
136  virtual void correctKinematics();
137 
138  //- Correct the turbulence
139  virtual void correctTurbulence();
140 
141  //- Correct the energy transport e.g. alphat
142  virtual void correctEnergyTransport();
143 
144 
145  // Momentum
146 
147  //- Return whether the phase is stationary
148  virtual bool stationary() const;
149 
150  //- Return the momentum equation
151  virtual tmp<fvVectorMatrix> UEqn();
152 
153  //- Return the momentum equation for the face-based algorithm
154  virtual tmp<fvVectorMatrix> UfEqn();
155 
156  //- Return the velocity
157  virtual tmp<volVectorField> U() const;
158 
159  //- Access the velocity
160  virtual volVectorField& URef();
161 
162  //- Return the volumetric flux
163  virtual tmp<surfaceScalarField> phi() const;
164 
165  //- Access the volumetric flux
166  virtual surfaceScalarField& phiRef();
167 
168  //- Return the volumetric flux of the phase
169  virtual tmp<surfaceScalarField> alphaPhi() const;
170 
171  //- Access the volumetric flux of the phase
172  virtual surfaceScalarField& alphaPhiRef();
173 
174  //- Return the mass flux of the phase
175  virtual tmp<surfaceScalarField> alphaRhoPhi() const;
176 
177  //- Access the mass flux of the phase
179 
180  //- Return the substantive acceleration
181  virtual tmp<volVectorField> DUDt() const;
182 
183  //- Return the substantive acceleration on the faces
184  virtual tmp<surfaceScalarField> DUDtf() const;
185 
186  //- Return the continuity error
187  virtual tmp<volScalarField> continuityError() const;
188 
189  //- Return the continuity error due to the flow field
191 
192  //- Return the continuity error due to any sources
194 
195  //- Return the phase kinetic energy
196  virtual tmp<volScalarField> K() const;
197 
198 
199  // Compressibility (variable density)
200 
201  //- Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
202  virtual tmp<volScalarField> divU() const;
203 
204  //- Set the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
205  virtual void divU(tmp<volScalarField> divU);
206 
207 
208  // Turbulence
209 
210  //- Return the turbulent dynamic viscosity
211  virtual tmp<volScalarField> mut() const;
212 
213  //- Return the effective dynamic viscosity
214  virtual tmp<volScalarField> muEff() const;
215 
216  //- Return the turbulent kinematic viscosity
217  virtual tmp<volScalarField> nut() const;
218 
219  //- Return the effective kinematic viscosity
220  virtual tmp<volScalarField> nuEff() const;
221 
222  //- Effective thermal turbulent diffusivity for temperature
223  // of mixture for patch [J/m/s/K]
225 
226  //- Return the effective thermal conductivity
227  virtual tmp<volScalarField> kappaEff() const;
228 
229  //- Return the effective thermal conductivity on a patch
230  virtual tmp<scalarField> kappaEff(const label patchi) const;
231 
232  //- Effective thermal turbulent diffusivity of mixture [kg/m/s]
234 
235  //- Return the effective thermal diffusivity
236  virtual tmp<volScalarField> alphaEff() const;
237 
238  //- Return the effective thermal conductivity on a patch
239  virtual tmp<scalarField> alphaEff(const label patchi) const;
240 
241  //- Return the turbulent kinetic energy
242  virtual tmp<volScalarField> k() const;
243 
244  //- Return the phase-pressure'
245  // (derivative of phase-pressure w.r.t. phase-fraction)
246  virtual tmp<volScalarField> pPrime() const;
247 };
248 
249 
250 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
251 
252 } // End namespace Foam
253 
254 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
255 
256 #ifdef NoRepository
257  #include "MovingPhaseModel.C"
258 #endif
259 
260 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
261 
262 #endif
263 
264 // ************************************************************************* //
virtual tmp< volVectorField > U() const
Access const reference to U.
twoPhaseSystem & fluid
autoPtr< phaseCompressibleTurbulenceModel > turbulence_
Turbulence model.
virtual tmp< fvVectorMatrix > UEqn()
Return the momentum equation.
virtual ~MovingPhaseModel()=default
Destructor.
virtual tmp< volScalarField > mut() const
Return the turbulent dynamic viscosity.
virtual tmp< surfaceScalarField > DUDtf() const
Return the substantive acceleration on the faces.
virtual tmp< fvVectorMatrix > UfEqn()
Return the momentum equation for the face-based algorithm.
virtual tmp< volScalarField > nuEff() const
Return the effective kinematic viscosity.
MovingPhaseModel(const multiphaseInterSystem &fluid, const word &phaseName)
virtual void correctKinematics()
Correct the kinematics.
virtual void correct()
Correct the phase properties other than the thermo and turbulence.
volScalarField continuityErrorFlow_
Continuity error due to the flow.
virtual surfaceScalarField & alphaPhiRef()
Access the volumetric flux of the phase.
tmp< volVectorField > DUDt_
Lagrangian acceleration field (needed for virtual-mass)
kappaEff
Definition: TEqn.H:10
virtual tmp< volScalarField > continuityErrorFlow() const
Return the continuity error due to the flow field.
virtual tmp< volScalarField > K() const
Return the phase kinetic energy.
virtual tmp< volVectorField > DUDt() const
Return the substantive acceleration.
virtual tmp< surfaceScalarField > alphaRhoPhi() const
Return the mass flux of the phase.
virtual volVectorField & URef()
Access the velocity.
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:69
A class for handling words, derived from Foam::string.
Definition: word.H:63
tmp< volScalarField > divU_
Dilatation rate.
virtual surfaceScalarField & phiRef()
Access the volumetric flux.
tmp< surfaceScalarField > DUDtf_
Lagrangian acceleration field on the faces (needed for virtual-mass)
virtual tmp< volScalarField > kappaEff() const
Return the effective thermal conductivity.
tmp< volScalarField > K_
Kinetic Energy.
virtual tmp< volScalarField > divU() const
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
virtual tmp< volScalarField > k() const
Return the turbulent kinetic energy.
virtual tmp< volScalarField > muEff() const
Return the effective dynamic viscosity.
virtual tmp< volScalarField > alphaEff() const
Return the effective thermal diffusivity.
virtual tmp< volScalarField > continuityError() const
Return the continuity error.
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
volScalarField alphaEff("alphaEff", turbulence->nu()/Pr+alphat)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
surfaceScalarField alphaRhoPhi_
Mass flux.
virtual void correctTurbulence()
Correct the turbulence.
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure&#39;.
virtual tmp< volScalarField > continuityErrorSources() const
Return the continuity error due to any sources.
virtual surfaceScalarField & alphaRhoPhiRef()
Access the mass flux of the phase.
virtual bool stationary() const
Return whether the phase is stationary.
volScalarField continuityErrorSources_
Continuity error due to any sources.
Namespace for OpenFOAM.
virtual tmp< volScalarField > nut() const
Return the turbulent kinematic viscosity.
virtual tmp< surfaceScalarField > phi() const
Constant access the volumetric flux.
virtual tmp< surfaceScalarField > alphaPhi() const
Constant access the volumetric flux of the phase.