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-2023 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  momentumTransport model and can generate the momentum equation. The
30  interface is quite restrictive as it also has to support an equivalent
31  stationary model, which does not store motion fields or a momentumTransport
32  model.
33 
34  Possible future extensions include separating the turbulent functionality
35  into another layer.
36 
37 See also
38  StationaryPhaseModel
39 
40 SourceFiles
41  MovingPhaseModel.C
42 
43 \*---------------------------------------------------------------------------*/
44 
45 #ifndef MovingPhaseModel_H
46 #define MovingPhaseModel_H
47 
48 #include "phaseModel.H"
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 namespace Foam
55 {
56 
57 // Trait for converting the ThermoModel's thermo type to the thermo type needed
58 // for the thermophysical transport model type; i.e., from rho-type thermo to
59 // fluid-type thermo.
60 
61 template<class ThermoModel>
63 
64 template<>
66 {
67  typedef fluidThermo type;
68 };
69 
70 template<>
72 {
74 };
75 
76 /*---------------------------------------------------------------------------*\
77  Class MovingPhaseModel Declaration
78 \*---------------------------------------------------------------------------*/
79 
80 template<class BasePhaseModel>
81 class MovingPhaseModel
82 :
83  public BasePhaseModel
84 {
85 protected:
86 
87  // Protected typedefs
88 
89  //- Thermo type for the thermophysical transport model
90  typedef
92  <
93  typename BasePhaseModel::thermoModel
96 
97 
98  // Protected data
99 
100  //- Velocity field
102 
103  //- Flux
105 
106  //- Volumetric flux
108 
109  //- Mass flux
111 
112  //- Face velocity field
114 
115  //- Lagrangian acceleration field (needed for virtual-mass)
116  mutable tmp<volVectorField> DUDt_;
117 
118  //- Lagrangian acceleration field on the faces (needed for virtual-mass)
120 
121  //- Dilatation rate
123 
124  //- Turbulence model
126 
127  //- Thermophysical transport model
128  autoPtr
129  <
131  <
134  >
136 
137  //- Continuity error
139 
140  //- Kinetic Energy
141  mutable tmp<volScalarField> K_;
142 
143 
144 private:
145 
146  // Private static member functions
147 
148  //- Calculate and return the flux field
150 
151 
152 public:
153 
154  // Constructors
155 
157  (
158  const phaseSystem& fluid,
159  const word& phaseName,
160  const bool referencePhase,
161  const label index
162  );
163 
164 
165  //- Destructor
166  virtual ~MovingPhaseModel();
167 
168 
169  // Member Functions
170 
171  //- Correct the phase properties other than the thermo
172  // and momentumTransport
173  virtual void correct();
174 
175  //- Correct the continuity error
176  virtual void correctContinuityError(const volScalarField& source);
177 
178  //- Correct the kinematics
179  virtual void correctKinematics();
180 
181  //- Predict the momentumTransport
182  virtual void predictMomentumTransport();
183 
184  //- Predict the energy transport e.g. alphat
185  virtual void predictThermophysicalTransport();
186 
187  //- Correct the momentumTransport
188  virtual void correctMomentumTransport();
189 
190  //- Correct the energy transport e.g. alphat
191  virtual void correctThermophysicalTransport();
192 
193  //- Correct the face velocity for moving meshes
194  virtual void correctUf();
195 
196 
197  // Momentum
198 
199  //- Return whether the phase is stationary
200  virtual bool stationary() const;
201 
202  //- Return the momentum equation
203  virtual tmp<fvVectorMatrix> UEqn();
204 
205  //- Return the momentum equation for the face-based algorithm
206  virtual tmp<fvVectorMatrix> UfEqn();
207 
208  //- Return the velocity
209  virtual tmp<volVectorField> U() const;
210 
211  //- Access the velocity
212  virtual volVectorField& URef();
213 
214  //- Access the velocity
215  virtual const volVectorField& URef() const;
216 
217  //- Return the volumetric flux
218  virtual tmp<surfaceScalarField> phi() const;
219 
220  //- Access the volumetric flux
221  virtual surfaceScalarField& phiRef();
222 
223  //- Access the volumetric flux
224  virtual const surfaceScalarField& phiRef() const;
225 
226  //- Return the face velocity
227  // Required for moving mesh cases
228  virtual const autoPtr<surfaceVectorField>& Uf() const;
229 
230  //- Access the face velocity
231  // Required for moving mesh cases
232  virtual surfaceVectorField& UfRef();
233 
234  //- Access the face velocity
235  // Required for moving mesh cases
236  virtual const surfaceVectorField& UfRef() const;
237 
238  //- Return the volumetric flux of the phase
239  virtual tmp<surfaceScalarField> alphaPhi() const;
240 
241  //- Access the volumetric flux of the phase
242  virtual surfaceScalarField& alphaPhiRef();
243 
244  //- Access the volumetric flux of the phase
245  virtual const surfaceScalarField& alphaPhiRef() const;
246 
247  //- Return the mass flux of the phase
248  virtual tmp<surfaceScalarField> alphaRhoPhi() const;
249 
250  //- Access the mass flux of the phase
252 
253  //- Access the mass flux of the phase
254  virtual const surfaceScalarField& alphaRhoPhiRef() const;
255 
256  //- Return the substantive acceleration
257  virtual tmp<volVectorField> DUDt() const;
258 
259  //- Return the substantive acceleration on the faces
260  virtual tmp<surfaceScalarField> DUDtf() const;
261 
262  //- Return the continuity error
263  virtual tmp<volScalarField> continuityError() const;
264 
265  //- Return the phase kinetic energy
266  virtual tmp<volScalarField> K() const;
267 
268 
269  // Compressibility (variable density)
270 
271  //- Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
272  virtual const autoPtr<volScalarField>& divU() const;
273 
274  //- Set the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
275  virtual void divU(tmp<volScalarField> divU);
276 
277 
278  // Momentum transport
279 
280  //- Return the turbulent kinetic energy
281  virtual tmp<volScalarField> k() const;
282 
283  //- Return the phase-pressure'
284  // (derivative of phase-pressure w.r.t. phase-fraction)
285  virtual tmp<volScalarField> pPrime() const;
286 
287 
288  // Thermophysical transport
289 
290  //- Return the effective thermal conductivity on a patch
291  virtual tmp<scalarField> kappaEff(const label patchi) const;
292 
293  //- Return the source term for the energy equation
294  virtual tmp<fvScalarMatrix> divq(volScalarField& he) const;
295 
296  //- Return the source term for the given specie mass-fraction
297  // equation
298  virtual tmp<fvScalarMatrix> divj(volScalarField& Yi) const;
299 };
300 
301 
302 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
303 
304 } // End namespace Foam
305 
306 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
307 
308 #ifdef NoRepository
309  #include "MovingPhaseModel.C"
310 #endif
311 
312 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
313 
314 #endif
315 
316 // ************************************************************************* //
Generic GeometricField class.
Class which represents a moving fluid phase. Holds the velocity, fluxes and momentumTransport model a...
virtual void correctUf()
Correct the face velocity for moving meshes.
virtual void correctKinematics()
Correct the kinematics.
autoPtr< volScalarField > divU_
Dilatation rate.
virtual tmp< fvVectorMatrix > UEqn()
Return the momentum equation.
virtual tmp< surfaceScalarField > DUDtf() const
Return the substantive acceleration on the faces.
MovingPhaseModelTransportThermoModel< typename BasePhaseModel::thermoModel >::type transportThermoModel
Thermo type for the thermophysical transport model.
volVectorField U_
Velocity field.
virtual void correct()
Correct the phase properties other than the thermo.
virtual tmp< volScalarField > k() const
Return the turbulent kinetic energy.
virtual void predictThermophysicalTransport()
Predict the energy transport e.g. alphat.
virtual ~MovingPhaseModel()
Destructor.
virtual tmp< volVectorField > U() const
Return the velocity.
virtual void correctMomentumTransport()
Correct the momentumTransport.
tmp< surfaceScalarField > DUDtf_
Lagrangian acceleration field on the faces (needed for virtual-mass)
autoPtr< PhaseThermophysicalTransportModel< phaseCompressible::momentumTransportModel, transportThermoModel > > thermophysicalTransport_
Thermophysical transport model.
volScalarField continuityError_
Continuity error.
virtual tmp< fvScalarMatrix > divj(volScalarField &Yi) const
Return the source term for the given specie mass-fraction.
virtual surfaceScalarField & alphaRhoPhiRef()
Access the mass flux of the phase.
surfaceScalarField alphaRhoPhi_
Mass flux.
virtual const autoPtr< volScalarField > & divU() const
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
virtual void correctThermophysicalTransport()
Correct the energy transport e.g. alphat.
virtual surfaceScalarField & phiRef()
Access the volumetric flux.
virtual tmp< surfaceScalarField > alphaPhi() const
Return the volumetric flux of the phase.
virtual tmp< volScalarField > K() const
Return the phase kinetic energy.
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure'.
surfaceScalarField phi_
Flux.
MovingPhaseModel(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
surfaceScalarField alphaPhi_
Volumetric flux.
virtual const autoPtr< surfaceVectorField > & Uf() const
Return the face velocity.
virtual surfaceVectorField & UfRef()
Access the face velocity.
tmp< volVectorField > DUDt_
Lagrangian acceleration field (needed for virtual-mass)
virtual surfaceScalarField & alphaPhiRef()
Access the volumetric flux of the phase.
virtual tmp< surfaceScalarField > alphaRhoPhi() const
Return the mass flux of the phase.
virtual tmp< surfaceScalarField > phi() const
Return the volumetric flux.
virtual tmp< fvVectorMatrix > UfEqn()
Return the momentum equation for the face-based algorithm.
virtual volVectorField & URef()
Access the velocity.
virtual tmp< volVectorField > DUDt() const
Return the substantive acceleration.
virtual bool stationary() const
Return whether the phase is stationary.
tmp< volScalarField > K_
Kinetic Energy.
virtual void predictMomentumTransport()
Predict the momentumTransport.
virtual tmp< volScalarField > continuityError() const
Return the continuity error.
virtual tmp< scalarField > kappaEff(const label patchi) const
Return the effective thermal conductivity on a patch.
autoPtr< phaseCompressible::momentumTransportModel > momentumTransport_
Turbulence model.
virtual void correctContinuityError(const volScalarField &source)
Correct the continuity error.
autoPtr< surfaceVectorField > Uf_
Face velocity field.
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
Templated base class for multiphase thermophysical transport models.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Base-class for multi-component fluid thermodynamic properties.
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:57
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:73
Base-class for multi-component fluid thermodynamic properties based on density.
Base-class for fluid thermodynamic properties based on density.
Definition: rhoThermo.H:55
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
label patchi
phaseCompressibleMomentumTransportModel momentumTransportModel
Namespace for OpenFOAM.
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
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
thermo he()