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