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"
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  //- Face velocity field
83 
84  //- Dilatation rate
86 
87  //- Turbulence model
89 
90  //- Continuity error
92 
93  //- Kinetic Energy
94  mutable tmp<volScalarField> K_;
95 
96 
97 private:
98 
99  // Private static member functions
100 
101  //- Calculate and return the flux field
103 
104 
105 public:
106 
107  // Constructors
108 
110  (
111  const phaseSystem& fluid,
112  const word& phaseName,
113  const bool referencePhase,
114  const label index
115  );
116 
117 
118  //- Destructor
119  virtual ~MovingPhaseModel();
120 
121 
122  // Member Functions
123 
124  //- Correct the phase properties other than the thermo
125  // and momentumTransport
126  virtual void correct();
127 
128  //- Correct the continuity error
129  virtual void correctContinuityError(const volScalarField& source);
130 
131  //- Correct the kinematics
132  virtual void correctKinematics();
133 
134  //- Predict the momentumTransport
135  virtual void predictMomentumTransport();
136 
137  //- Correct the momentumTransport
138  virtual void correctMomentumTransport();
139 
140  //- Correct the face velocity for moving meshes
141  virtual void correctUf();
142 
143 
144  // Momentum
145 
146  //- Return whether the phase is stationary
147  virtual bool stationary() const;
148 
149  //- Return the momentum equation
150  virtual tmp<fvVectorMatrix> UEqn();
151 
152  //- Return the momentum equation for the face-based algorithm
153  virtual tmp<fvVectorMatrix> UfEqn();
154 
155  //- Return the velocity
156  virtual tmp<volVectorField> U() const;
157 
158  //- Access the velocity
159  virtual volVectorField& URef();
160 
161  //- Access the velocity
162  virtual const volVectorField& URef() const;
163 
164  //- Return the volumetric flux
165  virtual tmp<surfaceScalarField> phi() const;
166 
167  //- Access the volumetric flux
168  virtual surfaceScalarField& phiRef();
169 
170  //- Access the volumetric flux
171  virtual const surfaceScalarField& phiRef() const;
172 
173  //- Return the face velocity
174  // Required for moving mesh cases
175  virtual const autoPtr<surfaceVectorField>& Uf() const;
176 
177  //- Access the face velocity
178  // Required for moving mesh cases
179  virtual surfaceVectorField& UfRef();
180 
181  //- Access the face velocity
182  // Required for moving mesh cases
183  virtual const surfaceVectorField& UfRef() const;
184 
185  //- Return the volumetric flux of the phase
186  virtual tmp<surfaceScalarField> alphaPhi() const;
187 
188  //- Access the volumetric flux of the phase
189  virtual surfaceScalarField& alphaPhiRef();
190 
191  //- Access the volumetric flux of the phase
192  virtual const surfaceScalarField& alphaPhiRef() const;
193 
194  //- Return the mass flux of the phase
195  virtual tmp<surfaceScalarField> alphaRhoPhi() const;
196 
197  //- Access the mass flux of the phase
199 
200  //- Access the mass flux of the phase
201  virtual const surfaceScalarField& alphaRhoPhiRef() const;
202 
203  //- Return the velocity transport matrix
204  virtual tmp<fvVectorMatrix> UgradU() const;
205 
206  //- Return the substantive acceleration matrix
207  virtual tmp<fvVectorMatrix> DUDt() const;
208 
209  //- Return the continuity error
210  virtual tmp<volScalarField> continuityError() const;
211 
212  //- Return the phase kinetic energy
213  virtual tmp<volScalarField> K() const;
214 
215 
216  // Compressibility (variable density)
217 
218  //- Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
219  virtual const autoPtr<volScalarField>& divU() const;
220 
221  //- Set the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
222  virtual void divU(tmp<volScalarField> divU);
223 
224 
225  // Momentum transport
226 
227  //- Return the turbulent kinetic energy
228  virtual tmp<volScalarField> k() const;
229 
230  //- Return the face-phase-pressure'
231  // (derivative of phase-pressure w.r.t. phase-fraction)
232  virtual tmp<surfaceScalarField> pPrimef() const;
233 };
234 
235 
236 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
237 
238 } // End namespace Foam
239 
240 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
241 
242 #ifdef NoRepository
243  #include "MovingPhaseModel.C"
244 #endif
245 
246 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
247 
248 #endif
249 
250 // ************************************************************************* //
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.
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 tmp< fvVectorMatrix > DUDt() const
Return the substantive acceleration matrix.
virtual ~MovingPhaseModel()
Destructor.
virtual tmp< volVectorField > U() const
Return the velocity.
virtual void correctMomentumTransport()
Correct the momentumTransport.
volScalarField continuityError_
Continuity error.
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 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.
surfaceScalarField phi_
Flux.
MovingPhaseModel(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
virtual tmp< surfaceScalarField > pPrimef() const
Return the face-phase-pressure'.
surfaceScalarField alphaPhi_
Volumetric flux.
virtual const autoPtr< surfaceVectorField > & Uf() const
Return the face velocity.
virtual surfaceVectorField & UfRef()
Access the face velocity.
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 bool stationary() const
Return whether the phase is stationary.
tmp< volScalarField > K_
Kinetic Energy.
virtual void predictMomentumTransport()
Predict the momentumTransport.
virtual tmp< fvVectorMatrix > UgradU() const
Return the velocity transport matrix.
virtual tmp< volScalarField > continuityError() const
Return the continuity error.
autoPtr< phaseCompressible::momentumTransportModel > momentumTransport_
Turbulence model.
virtual void correctContinuityError(const volScalarField &source)
Correct the continuity error.
autoPtr< surfaceVectorField > Uf_
Face velocity field.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:73
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
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