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-2020 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 functionality
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"
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 namespace Foam
54 {
55 
56 /*---------------------------------------------------------------------------*\
57  Class MovingPhaseModel Declaration
58 \*---------------------------------------------------------------------------*/
59 
60 template<class BasePhaseModel>
61 class MovingPhaseModel
62 :
63  public BasePhaseModel
64 {
65 protected:
66 
67  // Protected data
68 
69  //- Velocity field
71 
72  //- Flux
74 
75  //- Volumetric flux
77 
78  //- Mass flux
80 
81  //- Lagrangian acceleration field (needed for virtual-mass)
82  mutable tmp<volVectorField> DUDt_;
83 
84  //- Lagrangian acceleration field on the faces (needed for virtual-mass)
86 
87  //- Dilatation rate
89 
90  //- Turbulence model
92 
93  //- Thermophysical transport model
94  autoPtr
95  <
97  <
99  typename BasePhaseModel::thermoModel
100  >
102 
103  //- Continuity error
105 
106  //- Kinetic Energy
107  mutable tmp<volScalarField> K_;
108 
109 
110 private:
111 
112  // Private static member functions
113 
114  //- Calculate and return the flux field
116 
117 
118 public:
119 
120  // Constructors
121 
123  (
124  const phaseSystem& fluid,
125  const word& phaseName,
126  const bool referencePhase,
127  const label index
128  );
129 
130 
131  //- Destructor
132  virtual ~MovingPhaseModel();
133 
134 
135  // Member Functions
136 
137  //- Correct the phase properties other than the thermo and turbulence
138  virtual void correct();
139 
140  //- Correct the continuity error
141  virtual void correctContinuityError(const volScalarField& source);
142 
143  //- Correct the kinematics
144  virtual void correctKinematics();
145 
146  //- Correct the turbulence
147  virtual void correctTurbulence();
148 
149  //- Correct the energy transport e.g. alphat
150  virtual void correctEnergyTransport();
151 
152 
153  // Momentum
154 
155  //- Return whether the phase is stationary
156  virtual bool stationary() const;
157 
158  //- Return the momentum equation
159  virtual tmp<fvVectorMatrix> UEqn();
160 
161  //- Return the momentum equation for the face-based algorithm
162  virtual tmp<fvVectorMatrix> UfEqn();
163 
164  //- Return the velocity
165  virtual tmp<volVectorField> U() const;
166 
167  //- Access the velocity
168  virtual volVectorField& URef();
169 
170  //- Return the volumetric flux
171  virtual tmp<surfaceScalarField> phi() const;
172 
173  //- Access the volumetric flux
174  virtual surfaceScalarField& phiRef();
175 
176  //- Return the volumetric flux of the phase
177  virtual tmp<surfaceScalarField> alphaPhi() const;
178 
179  //- Access the volumetric flux of the phase
180  virtual surfaceScalarField& alphaPhiRef();
181 
182  //- Return the mass flux of the phase
183  virtual tmp<surfaceScalarField> alphaRhoPhi() const;
184 
185  //- Access the mass flux of the phase
187 
188  //- Return the substantive acceleration
189  virtual tmp<volVectorField> DUDt() const;
190 
191  //- Return the substantive acceleration on the faces
192  virtual tmp<surfaceScalarField> DUDtf() const;
193 
194  //- Return the continuity error
195  virtual tmp<volScalarField> continuityError() const;
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  // Momentum transport
211 
212  //- Return the effective thermal conductivity on a patch
213  virtual tmp<scalarField> kappaEff(const label patchi) const;
214 
215  //- Return the turbulent kinetic energy
216  virtual tmp<volScalarField> k() const;
217 
218  //- Return the phase-pressure'
219  // (derivative of phase-pressure w.r.t. phase-fraction)
220  virtual tmp<volScalarField> pPrime() const;
221 
222 
223  // Thermophysical transport
224 
225  //- Return the source term for the energy equation
226  virtual tmp<fvScalarMatrix> divq(volScalarField& he) const;
227 
228  //- Return the source term for the given specie mass-fraction
229  // equation
230  virtual tmp<fvScalarMatrix> divj(volScalarField& Yi) const;
231 };
232 
233 
234 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
235 
236 } // End namespace Foam
237 
238 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239 
240 #ifdef NoRepository
241  #include "MovingPhaseModel.C"
242 #endif
243 
244 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
245 
246 #endif
247 
248 // ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:50
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
PhaseCompressibleMomentumTransportModel< phaseModel > phaseCompressibleMomentumTransportModel
Typedef for phaseCompressibleMomentumTransportModel.
virtual surfaceScalarField & alphaRhoPhiRef()
Access the mass flux of the phase.
volVectorField U_
Velocity field.
autoPtr< phaseCompressibleMomentumTransportModel > turbulence_
Turbulence model.
virtual bool stationary() const
Return whether the phase is stationary.
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 tmp< fvScalarMatrix > divj(volScalarField &Yi) const
Return the source term for the given specie mass-fraction.
virtual void correctTurbulence()
Correct the turbulence.
Templated base class for multiphase thermophysical transport models.
volScalarField continuityError_
Continuity error.
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 ...
MovingPhaseModel(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
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.
surfaceScalarField phi_
Flux.
tmp< surfaceScalarField > DUDtf_
Lagrangian acceleration field on the faces (needed for virtual-mass)
phaseSystem & fluid
Definition: createFields.H:11
tmp< volScalarField > K_
Kinetic Energy.
virtual volVectorField & URef()
Access the velocity.
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure&#39;.
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< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
autoPtr< PhaseThermophysicalTransportModel< phaseCompressibleMomentumTransportModel, typename BasePhaseModel::thermoModel > > thermophysicalTransport_
Thermophysical transport model.
virtual tmp< volScalarField > K() const
Return the phase kinetic energy.
label patchi
virtual void correctContinuityError(const volScalarField &source)
Correct the continuity error.
virtual surfaceScalarField & alphaPhiRef()
Access the volumetric flux of the phase.
virtual tmp< fvVectorMatrix > UEqn()
Return the momentum equation.
virtual tmp< scalarField > kappaEff(const label patchi) const
Return the effective thermal conductivity on a patch.
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.
Namespace for OpenFOAM.