twoPhaseSystem.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) 2013-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::twoPhaseSystem
26 
27 Description
28 
29 SourceFiles
30  twoPhaseSystem.C
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #ifndef twoPhaseSystem_H
35 #define twoPhaseSystem_H
36 
37 #include "IOdictionary.H"
38 #include "phaseModel.H"
39 #include "phasePair.H"
40 #include "orderedPhasePair.H"
41 #include "volFields.H"
42 #include "surfaceFields.H"
43 #include "dragModel.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 class virtualMassModel;
51 class heatTransferModel;
52 class liftModel;
53 class wallLubricationModel;
54 class turbulentDispersionModel;
55 
56 class blendingMethod;
57 template <class modelType> class BlendedInterfacialModel;
58 
59 /*---------------------------------------------------------------------------*\
60  Class twoPhaseSystem Declaration
61 \*---------------------------------------------------------------------------*/
62 
63 class twoPhaseSystem
64 :
65  public IOdictionary
66 {
67  // Private data
68 
69  //- Reference to the mesh
70  const fvMesh& mesh_;
71 
72  //- Phase model 1
74 
75  //- Phase model 2
77 
78  //- Total volumetric flux
80 
81  //- Dilatation term
83 
84  //- Optional dispersion diffusivity
86 
87  //- Unordered phase pair
89 
90  //- Phase pair for phase 1 dispersed in phase 2
92 
93  //- Phase pair for phase 2 dispersed in phase 1
95 
96  //- Blending methods
98 
99  //- Drag model
101 
102  //- Virtual mass model
104 
105  //- Heat transfer model
107 
108  //- Lift model
110 
111  //- Wall lubrication model
114 
115  //- Wall lubrication model
118 
119 
120  // Private member functions
121 
122  //- Return the mixture flux
124 
125 
126 public:
127 
128  // Constructors
129 
130  //- Construct from fvMesh
131  twoPhaseSystem(const fvMesh&, const dimensionedVector& g);
132 
133 
134  //- Destructor
135  virtual ~twoPhaseSystem();
136 
137 
138  // Member Functions
139 
140  //- Return the mixture density
141  tmp<volScalarField> rho() const;
142 
143  //- Return the mixture velocity
144  tmp<volVectorField> U() const;
145 
146  //- Return the drag coefficient
147  tmp<volScalarField> Kd() const;
148 
149  //- Return the face drag coefficient
151 
152  //- Return the virtual mass coefficient
153  tmp<volScalarField> Vm() const;
154 
155  //- Return the face virtual mass coefficient
157 
158  //- Return the heat transfer coefficient
159  tmp<volScalarField> Kh() const;
160 
161  //- Return the combined force (lift + wall-lubrication)
162  tmp<volVectorField> F() const;
163 
164  //- Return the combined face-force (lift + wall-lubrication)
165  tmp<surfaceScalarField> Ff() const;
166 
167  //- Return the turbulent diffusivity
168  // Multiplies the phase-fraction gradient
169  tmp<volScalarField> D() const;
170 
171  //- Solve for the two-phase-fractions
172  void solve();
173 
174  //- Correct two-phase properties other than turbulence
175  void correct();
176 
177  //- Correct two-phase turbulence
178  void correctTurbulence();
179 
180  //- Read base phaseProperties dictionary
181  bool read();
182 
183  // Access
184 
185  //- Return the drag model for the given phase
186  const dragModel& drag(const phaseModel& phase) const;
187 
188  //- Return the virtual mass model for the given phase
189  const virtualMassModel& virtualMass(const phaseModel& phase) const;
190 
191  //- Return the surface tension coefficient
192  const dimensionedScalar& sigma() const;
193 
194  //- Return the mesh
195  inline const fvMesh& mesh() const;
196 
197  //- Return phase model 1
198  inline const phaseModel& phase1() const;
199 
200  //- Return non-const access to phase model 1
201  inline phaseModel& phase1();
202 
203  //- Return phase model 2
204  inline const phaseModel& phase2() const;
205 
206  //- Return non-const access to phase model 2
207  inline phaseModel& phase2();
208 
209  //- Return the phase not given as an argument
210  inline const phaseModel& otherPhase(const phaseModel& phase) const;
211 
212  //- Return the mixture flux
213  inline const surfaceScalarField& phi() const;
214 
215  //- Return non-const access to the the mixture flux
216  inline surfaceScalarField& phi();
217 
218  //- Return the dilatation term
219  inline const volScalarField& dgdt() const;
220 
221  //- Return non-const access to the dilatation parameter
222  inline volScalarField& dgdt();
223 
224  //- Return non-const access to the dispersion diffusivity
226 };
227 
228 
229 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
230 
231 } // End namespace Foam
232 
233 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
234 
235 #include "twoPhaseSystemI.H"
236 
237 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
238 
239 #endif
240 
241 // ************************************************************************* //
Foam::virtualMassModel
Definition: virtualMassModel.H:52
volFields.H
Foam::phaseModel
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
Foam::twoPhaseSystem::pPrimeByA
tmp< surfaceScalarField > & pPrimeByA()
Return non-const access to the dispersion diffusivity.
Definition: twoPhaseSystemI.H:98
Foam::twoPhaseSystem::mesh
const fvMesh & mesh() const
Return the mesh.
Definition: twoPhaseSystemI.H:28
Foam::twoPhaseSystem::pair1In2_
autoPtr< orderedPhasePair > pair1In2_
Phase pair for phase 1 dispersed in phase 2.
Definition: twoPhaseSystem.H:90
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::phase
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:52
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::twoPhaseSystem::rho
tmp< volScalarField > rho() const
Return the mixture density.
Foam::twoPhaseSystem::Vm
tmp< volScalarField > Vm() const
Return the virtual mass coefficient.
Foam::twoPhaseSystem::virtualMass_
autoPtr< BlendedInterfacialModel< virtualMassModel > > virtualMass_
Virtual mass model.
Definition: twoPhaseSystem.H:102
Foam::twoPhaseSystem::correct
void correct()
Correct two-phase properties other than turbulence.
Foam::twoPhaseSystem::phi_
surfaceScalarField phi_
Total volumetric flux.
Definition: twoPhaseSystem.H:78
Foam::twoPhaseSystem::heatTransfer_
autoPtr< BlendedInterfacialModel< heatTransferModel > > heatTransfer_
Heat transfer model.
Definition: twoPhaseSystem.H:105
Foam::twoPhaseSystem::Kh
tmp< volScalarField > Kh() const
Return the heat transfer coefficient.
Foam::twoPhaseSystem::sigma
tmp< volScalarField > sigma() const
Return the surface tension coefficient.
g
const dimensionedVector & g
Definition: setRegionFluidFields.H:33
Foam::twoPhaseSystem::drag_
autoPtr< BlendedInterfacialModel< dragModel > > drag_
Drag model.
Definition: twoPhaseSystem.H:99
surfaceFields.H
Foam::surfaceFields.
Foam::twoPhaseSystem::Kd
tmp< volScalarField > Kd() const
Return the drag coefficient.
Foam::twoPhaseSystem::phase1_
phaseModel phase1_
Phase model 1.
Definition: twoPhaseSystem.H:72
Foam::twoPhaseSystem::phase2
const phaseModel & phase2() const
Constant access phase model 2.
Definition: twoPhaseSystemI.H:40
Foam::twoPhaseSystem::pair2In1_
autoPtr< orderedPhasePair > pair2In1_
Phase pair for phase 2 dispersed in phase 1.
Definition: twoPhaseSystem.H:93
Foam::twoPhaseSystem::pair_
autoPtr< phasePair > pair_
Unordered phase pair.
Definition: twoPhaseSystem.H:87
Foam::twoPhaseSystem::calcPhi
tmp< surfaceScalarField > calcPhi() const
Return the mixture flux.
Foam::twoPhaseSystem::dgdt_
volScalarField dgdt_
Dilatation term.
Definition: twoPhaseSystem.H:81
Foam::twoPhaseSystem::otherPhase
const phaseModel & otherPhase(const phaseModel &phase) const
Constant access the phase not given as an argument.
Definition: twoPhaseSystemI.H:53
Foam::twoPhaseSystem::pPrimeByA_
tmp< surfaceScalarField > pPrimeByA_
Optional dispersion diffusivity.
Definition: twoPhaseSystem.H:84
Foam::twoPhaseSystem::phase1
const phaseModel & phase1() const
Constant access phase model 1.
Definition: twoPhaseSystemI.H:28
Foam::twoPhaseSystem::wallLubrication_
autoPtr< BlendedInterfacialModel< wallLubricationModel > > wallLubrication_
Wall lubrication model.
Definition: twoPhaseSystem.H:112
Foam::twoPhaseSystem::phi
const surfaceScalarField & phi() const
Return the mixture flux.
Definition: twoPhaseSystemI.H:74
twoPhaseSystemI.H
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::twoPhaseSystem::twoPhaseSystem
twoPhaseSystem(const fvMesh &)
Construct from fvMesh.
Foam::twoPhaseSystem::lift_
autoPtr< BlendedInterfacialModel< liftModel > > lift_
Lift model.
Definition: twoPhaseSystem.H:108
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::twoPhaseSystem::phase2_
phaseModel phase2_
Phase model 2.
Definition: twoPhaseSystem.H:75
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::twoPhaseSystem::phase2_
phaseModel & phase2_
Phase model 2.
Definition: twoPhaseSystem.H:95
Foam::twoPhaseSystem::Vmf
tmp< surfaceScalarField > Vmf() const
Return the face virtual mass coefficient.
Foam::twoPhaseSystem::virtualMass
const virtualMassModel & virtualMass(const phaseModel &phase) const
Return the virtual mass model for the given phase.
Foam::HashTable
An STL-conforming hash table.
Definition: HashTable.H:61
IOdictionary.H
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::string::hash
Hashing function class, shared by all the derived classes.
Definition: string.H:90
Foam::twoPhaseSystem::turbulentDispersion_
autoPtr< BlendedInterfacialModel< turbulentDispersionModel > > turbulentDispersion_
Wall lubrication model.
Definition: twoPhaseSystem.H:116
Foam::twoPhaseSystem::Kdf
tmp< surfaceScalarField > Kdf() const
Return the face drag coefficient.
Foam::twoPhaseSystem::D
tmp< volScalarField > D() const
Return the turbulent diffusivity.
Foam::twoPhaseSystem::read
bool read()
Read base phaseProperties dictionary.
Foam::twoPhaseSystem::dgdt
const volScalarField & dgdt() const
Return the dilatation term.
Definition: twoPhaseSystemI.H:86
Foam::twoPhaseSystem::phase1_
phaseModel & phase1_
Phase model 1.
Definition: twoPhaseSystem.H:92
Foam::twoPhaseSystem::correctTurbulence
void correctTurbulence()
Correct two-phase turbulence.
Foam::dragModel
Definition: dragModel.H:50
Foam::twoPhaseSystem::F
tmp< volVectorField > F() const
Return the combined force (lift + wall-lubrication)
Foam::twoPhaseSystem::solve
virtual void solve()
Solve for the phase fractions.
Foam::twoPhaseSystem::U
tmp< volVectorField > U() const
Return the mixture velocity.
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::twoPhaseSystem::mesh_
const fvMesh & mesh_
Reference to the mesh.
Definition: twoPhaseSystem.H:69
Foam::twoPhaseSystem::drag
const dragModel & drag(const phaseModel &phase) const
Return the drag model for the given phase.
Foam::twoPhaseSystem::~twoPhaseSystem
virtual ~twoPhaseSystem()
Destructor.
Foam::twoPhaseSystem::blendingMethods_
HashTable< autoPtr< blendingMethod >, word, word::hash > blendingMethods_
Blending methods.
Definition: twoPhaseSystem.H:96
Foam::twoPhaseSystem::Ff
tmp< surfaceScalarField > Ff() const
Return the combined face-force (lift + wall-lubrication)