MovingPhaseModel.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2015 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. Provides access to the turbulent quantities.
30 
31  Possible future extensions include separating the turbulent fuctionality
32  into another layer. It should also be possible to replace this layer with a
33  stationary phase model, in order to model packed beds or simple porous
34  media. This would probably require extra functionality, such as returning
35  the inputs into the general pressure equation (A, HbyA, etc ...).
36 
37  Note that this class does not return the turbulence model, it just provides
38  indirect access to the turbulent data. This is so a layer without
39  turbulence modelling (such as a stationary model) could be substituted.
40 
41 SourceFiles
42  MovingPhaseModel.C
43 
44 \*---------------------------------------------------------------------------*/
45 
46 #ifndef MovingPhaseModel_H
47 #define MovingPhaseModel_H
48 
49 #include "phaseModel.H"
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 namespace Foam
55 {
56 
57 /*---------------------------------------------------------------------------*\
58  Class phaseModel Declaration
59 \*---------------------------------------------------------------------------*/
60 
61 template<class BasePhaseModel>
62 class MovingPhaseModel
63 :
64  public BasePhaseModel
65 {
66  // Private 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)
82 
83  //- Dilatation rate
85 
86  //- Turbulence model
88 
89  //- Continuity error
91 
92  //- Phase diffusivity divided by the momentum coefficient.
93  // Used for implicit treatment of the phase pressure and dispersion
95 
96 
97  // Private static member functions
98 
99  //- Calculate and return the flux field
101 
102 
103 public:
104 
105  // Constructors
106 
108  (
109  const phaseSystem& fluid,
110  const word& phaseName,
111  const label index
112  );
113 
114 
115  //- Destructor
116  virtual ~MovingPhaseModel();
117 
118 
119  // Member Functions
120 
121  //- Correct the phase properties other than the thermo and turbulence
122  virtual void correct();
123 
124  //- Correct the kinematics
125  virtual void correctKinematics();
126 
127  //- Correct the turbulence
128  virtual void correctTurbulence();
129 
130  //- Correct the energy transport e.g. alphat
131  virtual void correctEnergyTransport();
132 
133  //- Return the momentum equation
134  virtual tmp<fvVectorMatrix> UEqn();
135 
136 
137  // Implicit phase pressure and dispersion support
138 
139  //- Return the phase diffusivity divided by the momentum coefficient
140  virtual const surfaceScalarField& DbyA() const;
141 
142  //- Set the phase diffusivity divided by the momentum coefficient
143  virtual void DbyA(const tmp<surfaceScalarField>& DbyA);
144 
145 
146  // Momentum
147 
148  //- Constant access the velocity
149  virtual tmp<volVectorField> U() const;
150 
151  //- Access the velocity
152  virtual volVectorField& U();
153 
154  //- Return the substantive acceleration
155  virtual tmp<volVectorField> DUDt() const;
156 
157  //- Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
158  virtual const tmp<volScalarField>& divU() const;
159 
160  //- Set the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
161  virtual void divU(const tmp<volScalarField>& divU);
162 
163  //- Constant access the continuity error
164  virtual tmp<volScalarField> continuityError() const;
165 
166  //- Constant access the volumetric flux
167  virtual tmp<surfaceScalarField> phi() const;
168 
169  //- Access the volumetric flux
170  virtual surfaceScalarField& phi();
171 
172  //- Constant access the volumetric flux of the phase
173  virtual tmp<surfaceScalarField> alphaPhi() const;
174 
175  //- Access the volumetric flux of the phase
176  virtual surfaceScalarField& alphaPhi();
177 
178  //- Constant access the mass flux of the phase
179  virtual tmp<surfaceScalarField> alphaRhoPhi() const;
180 
181  //- Access the mass flux of the phase
182  virtual surfaceScalarField& alphaRhoPhi();
183 
184 
185  // Turbulence
186 
187  //- Return the turbulence model
188  virtual const phaseCompressibleTurbulenceModel& turbulence() const;
189 };
190 
191 
192 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
193 
194 } // End namespace Foam
195 
196 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
197 
198 #ifdef NoRepository
199 # include "MovingPhaseModel.C"
200 #endif
201 
202 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
203 
204 #endif
205 
206 // ************************************************************************* //
Foam::MovingPhaseModel::U_
volVectorField U_
Velocity field.
Definition: MovingPhaseModel.H:68
Foam::MovingPhaseModel::divU_
tmp< volScalarField > divU_
Dilatation rate.
Definition: MovingPhaseModel.H:83
Foam::MovingPhaseModel::alphaRhoPhi_
surfaceScalarField alphaRhoPhi_
Mass flux.
Definition: MovingPhaseModel.H:77
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::MovingPhaseModel::phi_
surfaceScalarField phi_
Flux.
Definition: MovingPhaseModel.H:71
Foam::MovingPhaseModel::alphaPhi
virtual tmp< surfaceScalarField > alphaPhi() const
Constant access the volumetric flux of the phase.
Foam::MovingPhaseModel::U
virtual tmp< volVectorField > U() const
Constant access the velocity.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
fluid
multiphaseSystem & fluid
Definition: createFields.H:10
Foam::MovingPhaseModel::continuityError_
volScalarField continuityError_
Continuity error.
Definition: MovingPhaseModel.H:89
Foam::MovingPhaseModel::alphaPhi_
surfaceScalarField alphaPhi_
Volumetric flux.
Definition: MovingPhaseModel.H:74
Foam::MovingPhaseModel::DbyA_
tmp< surfaceScalarField > DbyA_
Phase diffusivity divided by the momentum coefficient.
Definition: MovingPhaseModel.H:93
Foam::MovingPhaseModel
Class which represents a moving fluid phase. Holds the velocity, fluxes and turbulence model....
Definition: MovingPhaseModel.H:61
Foam::MovingPhaseModel::~MovingPhaseModel
virtual ~MovingPhaseModel()
Destructor.
Foam::MovingPhaseModel::phi
virtual tmp< surfaceScalarField > phi() const
Constant access the volumetric flux.
Foam::MovingPhaseModel::DbyA
virtual const surfaceScalarField & DbyA() const
Return the phase diffusivity divided by the momentum coefficient.
Foam::MovingPhaseModel::correctKinematics
virtual void correctKinematics()
Correct the kinematics.
Foam::MovingPhaseModel::DUDt
virtual tmp< volVectorField > DUDt() const
Return the substantive acceleration.
Foam::MovingPhaseModel::MovingPhaseModel
MovingPhaseModel(const phaseSystem &fluid, const word &phaseName, const label index)
Foam::label
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
Foam::MovingPhaseModel::correctTurbulence
virtual void correctTurbulence()
Correct the turbulence.
Foam::MovingPhaseModel::turbulence_
autoPtr< phaseCompressibleTurbulenceModel > turbulence_
Turbulence model.
Definition: MovingPhaseModel.H:86
Foam::MovingPhaseModel::UEqn
virtual tmp< fvVectorMatrix > UEqn()
Return the momentum equation.
Foam::MovingPhaseModel::turbulence
virtual const phaseCompressibleTurbulenceModel & turbulence() const
Return the turbulence model.
phaseCompressibleTurbulenceModel.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::MovingPhaseModel::alphaRhoPhi
virtual tmp< surfaceScalarField > alphaRhoPhi() const
Constant access the mass flux of the phase.
Foam::MovingPhaseModel::divU
virtual const tmp< volScalarField > & divU() const
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::ThermalDiffusivity
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Definition: ThermalDiffusivity.H:48
MovingPhaseModel.C
Foam::MovingPhaseModel::correct
virtual void correct()
Correct the phase properties other than the thermo and turbulence.
Foam::MovingPhaseModel::correctEnergyTransport
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat.
Foam::MovingPhaseModel::continuityError
virtual tmp< volScalarField > continuityError() const
Constant access the continuity error.
Foam::phaseSystem
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:67
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::MovingPhaseModel::DUDt_
volVectorField DUDt_
Lagrangian acceleration field (needed for virtual-mass)
Definition: MovingPhaseModel.H:80