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-2021 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 // Trait for converting the ThermoModel's thermo type to the thermo type needed
57 // for the thermophysical transport model type; i.e., from rho-type thermo to
58 // fluid-type thermo.
59 
60 template<class ThermoModel>
62 
63 template<>
65 {
66  typedef fluidThermo type;
67 };
68 
69 template<>
71 {
72  typedef fluidReactionThermo type;
73 };
74 
75 /*---------------------------------------------------------------------------*\
76  Class MovingPhaseModel Declaration
77 \*---------------------------------------------------------------------------*/
78 
79 template<class BasePhaseModel>
80 class MovingPhaseModel
81 :
82  public BasePhaseModel
83 {
84 protected:
85 
86  // Protected typedefs
87 
88  //- Thermo type for the thermophysical transport model
89  typedef
91  <
92  typename BasePhaseModel::thermoModel
95 
96 
97  // Protected data
98 
99  //- Velocity field
100  volVectorField U_;
101 
102  //- Flux
103  surfaceScalarField phi_;
104 
105  //- Volumetric flux
106  surfaceScalarField alphaPhi_;
107 
108  //- Mass flux
109  surfaceScalarField alphaRhoPhi_;
110 
111  //- Face velocity field
113 
114  //- Lagrangian acceleration field (needed for virtual-mass)
115  mutable tmp<volVectorField> DUDt_;
116 
117  //- Lagrangian acceleration field on the faces (needed for virtual-mass)
118  mutable tmp<surfaceScalarField> DUDtf_;
119 
120  //- Dilatation rate
121  tmp<volScalarField> divU_;
122 
123  //- Turbulence model
125 
126  //- Thermophysical transport model
127  autoPtr
128  <
130  <
132  transportThermoModel
133  >
134  > thermophysicalTransport_;
135 
136  //- Continuity error
137  volScalarField continuityError_;
138 
139  //- Kinetic Energy
140  mutable tmp<volScalarField> K_;
141 
142 
143 private:
144 
145  // Private static member functions
146 
147  //- Calculate and return the flux field
149 
150 
151 public:
152 
153  // Constructors
154 
156  (
157  const phaseSystem& fluid,
158  const word& phaseName,
159  const bool referencePhase,
160  const label index
161  );
162 
163 
164  //- Destructor
165  virtual ~MovingPhaseModel();
166 
167 
168  // Member Functions
169 
170  //- Correct the phase properties other than the thermo and turbulence
171  virtual void correct();
172 
173  //- Correct the continuity error
174  virtual void correctContinuityError(const volScalarField& source);
175 
176  //- Correct the kinematics
177  virtual void correctKinematics();
178 
179  //- Correct the turbulence
180  virtual void correctTurbulence();
181 
182  //- Correct the energy transport e.g. alphat
183  virtual void correctEnergyTransport();
184 
185  //- Correct the face velocity for moving meshes
186  virtual void correctUf();
187 
188 
189  // Momentum
190 
191  //- Return whether the phase is stationary
192  virtual bool stationary() const;
193 
194  //- Return the momentum equation
195  virtual tmp<fvVectorMatrix> UEqn();
196 
197  //- Return the momentum equation for the face-based algorithm
198  virtual tmp<fvVectorMatrix> UfEqn();
199 
200  //- Return the velocity
201  virtual tmp<volVectorField> U() const;
202 
203  //- Access the velocity
204  virtual volVectorField& URef();
205 
206  //- Return the volumetric flux
207  virtual tmp<surfaceScalarField> phi() const;
208 
209  //- Access the volumetric flux
210  virtual surfaceScalarField& phiRef();
211 
212  //- Return the face velocity
213  // Required for moving mesh cases
214  virtual tmp<surfaceVectorField> Uf() const;
215 
216  //- Access the face velocity
217  // Required for moving mesh cases
218  virtual surfaceVectorField& UfRef();
219 
220  //- Return the volumetric flux of the phase
221  virtual tmp<surfaceScalarField> alphaPhi() const;
222 
223  //- Access the volumetric flux of the phase
224  virtual surfaceScalarField& alphaPhiRef();
225 
226  //- Return the mass flux of the phase
227  virtual tmp<surfaceScalarField> alphaRhoPhi() const;
228 
229  //- Access the mass flux of the phase
230  virtual surfaceScalarField& alphaRhoPhiRef();
231 
232  //- Return the substantive acceleration
233  virtual tmp<volVectorField> DUDt() const;
234 
235  //- Return the substantive acceleration on the faces
236  virtual tmp<surfaceScalarField> DUDtf() const;
237 
238  //- Return the continuity error
239  virtual tmp<volScalarField> continuityError() const;
240 
241  //- Return the phase kinetic energy
242  virtual tmp<volScalarField> K() const;
243 
244 
245  // Compressibility (variable density)
246 
247  //- Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
248  virtual tmp<volScalarField> divU() const;
249 
250  //- Set the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
251  virtual void divU(tmp<volScalarField> divU);
252 
253 
254  // Momentum transport
255 
256  //- Return the turbulent kinetic energy
257  virtual tmp<volScalarField> k() const;
258 
259  //- Return the phase-pressure'
260  // (derivative of phase-pressure w.r.t. phase-fraction)
261  virtual tmp<volScalarField> pPrime() const;
262 
263 
264  // Thermophysical transport
265 
266  //- Return the effective thermal conductivity on a patch
267  virtual tmp<scalarField> kappaEff(const label patchi) const;
268 
269  //- Return the source term for the energy equation
270  virtual tmp<fvScalarMatrix> divq(volScalarField& he) const;
271 
272  //- Return the source term for the given specie mass-fraction
273  // equation
274  virtual tmp<fvScalarMatrix> divj(volScalarField& Yi) const;
275 };
276 
277 
278 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
279 
280 } // End namespace Foam
281 
282 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
283 
284 #ifdef NoRepository
285  #include "MovingPhaseModel.C"
286 #endif
287 
288 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
289 
290 #endif
291 
292 // ************************************************************************* //
U
Definition: pEqn.H:72
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-10/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:68
Base-class for multi-component fluid thermodynamic properties.
label k
Boltzmann constant.
Templated base class for multiphase thermophysical transport models.
CGAL::Exact_predicates_exact_constructions_kernel K
Class which represents a moving fluid phase. Holds the velocity, fluxes and turbulence model and can ...
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:68
A class for handling words, derived from string.
Definition: word.H:59
Base-class for multi-component fluid thermodynamic properties based on density.
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:53
MovingPhaseModelTransportThermoModel< typename BasePhaseModel::thermoModel >::type transportThermoModel
Thermo type for the thermophysical transport model.
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))
phaseSystem & fluid
Definition: createFields.H:11
phi
Definition: correctPhi.H:3
autoPtr< surfaceVectorField > Uf
thermo he()
tmp< volScalarField > divU
Definition: alphaSuSp.H:9
phaseCompressibleMomentumTransportModel momentumTransportModel
label patchi
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
void correctUf(autoPtr< surfaceVectorField > &Uf, const volVectorField &U, const surfaceScalarField &phi, const MRFType &MRF)
Base-class for fluid thermodynamic properties based on density.
Definition: rhoThermo.H:52
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
Namespace for OpenFOAM.
fvVectorMatrix & UEqn
Definition: UEqn.H:13