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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2013-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::twoPhaseSystem
28 
29 Description
30 
31 SourceFiles
32  twoPhaseSystem.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef twoPhaseSystem_H
37 #define twoPhaseSystem_H
38 
39 #include "IOdictionary.H"
40 #include "phaseModel.H"
41 #include "phasePair.H"
42 #include "orderedPhasePair.H"
43 #include "volFields.H"
44 #include "surfaceFields.H"
45 #include "dragModel.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 class virtualMassModel;
53 class heatTransferModel;
54 class liftModel;
55 class wallLubricationModel;
56 class turbulentDispersionModel;
57 
58 class blendingMethod;
59 template<class modelType> class BlendedInterfacialModel;
60 
61 /*---------------------------------------------------------------------------*\
62  Class twoPhaseSystem Declaration
63 \*---------------------------------------------------------------------------*/
64 
65 class twoPhaseSystem
66 :
67  public IOdictionary
68 {
69  // Private data
70 
71  //- Reference to the mesh
72  const fvMesh& mesh_;
73 
74  //- Phase model 1
75  phaseModel phase1_;
76 
77  //- Phase model 2
78  phaseModel phase2_;
79 
80  //- Total volumetric flux
81  surfaceScalarField phi_;
82 
83  //- Dilatation term
84  volScalarField dgdt_;
85 
86  //- Optional dispersion diffusivity
87  tmp<surfaceScalarField> pPrimeByA_;
88 
89  //- Unordered phase pair
90  autoPtr<phasePair> pair_;
91 
92  //- Phase pair for phase 1 dispersed in phase 2
93  autoPtr<orderedPhasePair> pair1In2_;
94 
95  //- Phase pair for phase 2 dispersed in phase 1
96  autoPtr<orderedPhasePair> pair2In1_;
97 
98  //- Blending methods
99  HashTable<autoPtr<blendingMethod>> blendingMethods_;
100 
101  //- Drag model
102  autoPtr<BlendedInterfacialModel<dragModel>> drag_;
103 
104  //- Virtual mass model
105  autoPtr<BlendedInterfacialModel<virtualMassModel>> virtualMass_;
106 
107  //- Heat transfer model
108  autoPtr<BlendedInterfacialModel<heatTransferModel>> heatTransfer_;
109 
110  //- Lift model
111  autoPtr<BlendedInterfacialModel<liftModel>> lift_;
112 
113  //- Wall lubrication model
114  autoPtr<BlendedInterfacialModel<wallLubricationModel>>
115  wallLubrication_;
116 
117  //- Wall lubrication model
118  autoPtr<BlendedInterfacialModel<turbulentDispersionModel>>
119  turbulentDispersion_;
120 
121 
122  // Private member functions
123 
124  //- Return the mixture flux
125  tmp<surfaceScalarField> calcPhi() const;
126 
127 
128 public:
129 
130  // Constructors
131 
132  //- Construct from fvMesh
133  twoPhaseSystem(const fvMesh&, const dimensionedVector& g);
134 
135 
136  //- Destructor
137  virtual ~twoPhaseSystem();
138 
139 
140  // Member Functions
141 
142  //- Return the mixture density
143  tmp<volScalarField> rho() const;
144 
145  //- Return the mixture velocity
146  tmp<volVectorField> U() const;
147 
148  //- Return the drag coefficient
149  tmp<volScalarField> Kd() const;
150 
151  //- Return the face drag coefficient
152  tmp<surfaceScalarField> Kdf() const;
153 
154  //- Return the virtual mass coefficient
155  tmp<volScalarField> Vm() const;
156 
157  //- Return the face virtual mass coefficient
158  tmp<surfaceScalarField> Vmf() const;
159 
160  //- Return the heat transfer coefficient
161  tmp<volScalarField> Kh() const;
162 
163  //- Return the combined force (lift + wall-lubrication)
164  tmp<volVectorField> F() const;
165 
166  //- Return the combined face-force (lift + wall-lubrication)
167  tmp<surfaceScalarField> Ff() const;
168 
169  //- Return the turbulent diffusivity
170  // Multiplies the phase-fraction gradient
171  tmp<volScalarField> D() const;
172 
173  //- Solve for the two-phase-fractions
174  void solve();
175 
176  //- Correct two-phase properties other than turbulence
177  void correct();
178 
179  //- Correct two-phase turbulence
180  void correctTurbulence();
181 
182  //- Read base phaseProperties dictionary
183  bool read();
184 
185  // Access
186 
187  //- Access a sub model between a phase pair
188  template<class modelType>
189  const modelType& lookupSubModel(const phasePair& key) const;
190 
191  //- Access a sub model between two phases
192  template<class modelType>
193  const modelType& lookupSubModel
194  (
195  const phaseModel& dispersed,
196  const phaseModel& continuous
197  ) const;
198 
199  //- Return the surface tension coefficient
200  const dimensionedScalar& sigma() const;
201 
202  //- Return the mesh
203  inline const fvMesh& mesh() const;
204 
205  //- Return phase model 1
206  inline const phaseModel& phase1() const;
207 
208  //- Return non-const access to phase model 1
209  inline phaseModel& phase1();
210 
211  //- Return phase model 2
212  inline const phaseModel& phase2() const;
213 
214  //- Return non-const access to phase model 2
215  inline phaseModel& phase2();
216 
217  //- Return the phase not given as an argument
218  inline const phaseModel& otherPhase(const phaseModel& phase) const;
219 
220  //- Return the mixture flux
221  inline const surfaceScalarField& phi() const;
222 
223  //- Return non-const access to the mixture flux
224  inline surfaceScalarField& phi();
225 
226  //- Return the dilatation term
227  inline const volScalarField& dgdt() const;
228 
229  //- Return non-const access to the dilatation parameter
230  inline volScalarField& dgdt();
231 
232  //- Return non-const access to the dispersion diffusivity
233  inline tmp<surfaceScalarField>& pPrimeByA();
234 };
235 
236 
237 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
238 
239 } // End namespace Foam
240 
241 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
242 
243 #include "twoPhaseSystemI.H"
244 
245 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
246 
247 #ifdef NoRepository
248  #include "twoPhaseSystemTemplates.C"
249 #endif
250 
251 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
252 
253 #endif
254 
255 // ************************************************************************* //
volFields.H
Foam::twoPhaseSystem::F
tmp< volVectorField > F() const
Return the combined force (lift + wall-lubrication)
Definition: twoPhaseSystem.C:334
Foam::twoPhaseSystem::pPrimeByA
tmp< surfaceScalarField > & pPrimeByA()
Return non-const access to the dispersion diffusivity.
Definition: twoPhaseSystemI.H:100
Foam::twoPhaseSystem::mesh
const fvMesh & mesh() const
Return the mesh.
Definition: twoPhaseSystemI.H:30
Foam::twoPhaseSystem::phase1
const phaseModel & phase1() const
Return phase model 1.
Definition: twoPhaseSystemI.H:30
Foam::twoPhaseSystem::correct
void correct()
Correct two-phase properties other than turbulence.
Definition: twoPhaseSystem.C:541
Foam::twoPhaseSystem::D
tmp< volScalarField > D() const
Return the turbulent diffusivity.
Definition: twoPhaseSystem.C:346
Foam::glTF::key
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
Definition: foamGltfBase.H:108
Foam::twoPhaseSystem::sigma
tmp< volScalarField > sigma() const
Return the surface tension coefficient.
Definition: twoPhaseSystem.C:83
Foam::twoPhaseSystem::U
tmp< volVectorField > U() const
Return the mixture velocity.
Definition: twoPhaseSystem.C:290
surfaceFields.H
Foam::surfaceFields.
Foam::twoPhaseSystem::phase2
const phaseModel & phase2() const
Return phase model 2.
Definition: twoPhaseSystemI.H:42
Foam::twoPhaseSystem::Ff
tmp< surfaceScalarField > Ff() const
Return the combined face-force (lift + wall-lubrication)
Definition: twoPhaseSystem.C:340
Foam::twoPhaseSystem::lookupSubModel
const modelType & lookupSubModel(const phasePair &key) const
Access a sub model between a phase pair.
Definition: twoPhaseSystemTemplates.C:41
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
Foam::twoPhaseSystem::dgdt
const volScalarField & dgdt() const
Return the dilatation term.
Definition: twoPhaseSystemI.H:88
Foam::twoPhaseSystem::Vm
tmp< volScalarField > Vm() const
Return the virtual mass coefficient.
Definition: twoPhaseSystem.C:113
Foam::twoPhaseSystem::otherPhase
const phaseModel & otherPhase(const phaseModel &phase) const
Return the phase not given as an argument.
Definition: twoPhaseSystemI.H:55
Foam::twoPhaseSystem::Kdf
tmp< surfaceScalarField > Kdf() const
Return the face drag coefficient.
Definition: twoPhaseSystem.C:103
Foam::twoPhaseSystem::Vmf
tmp< surfaceScalarField > Vmf() const
Return the face virtual mass coefficient.
Definition: twoPhaseSystem.C:322
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::twoPhaseSystem::solve
virtual void solve()
Solve for the phase fractions.
Definition: twoPhaseSystem.C:122
Foam::twoPhaseSystem::Kh
tmp< volScalarField > Kh() const
Return the heat transfer coefficient.
Definition: twoPhaseSystem.C:328
Foam::twoPhaseSystem::twoPhaseSystem
twoPhaseSystem(const fvMesh &)
Construct from fvMesh.
Definition: twoPhaseSystem.C:59
g
const uniformDimensionedVectorField & g
Definition: createFluidFields.H:26
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::twoPhaseSystem::phase2_
phaseModel & phase2_
Phase model 2.
Definition: twoPhaseSystem.H:84
twoPhaseSystemTemplates.C
IOdictionary.H
twoPhaseSystemI.H
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
Foam::twoPhaseSystem::read
bool read()
Read base phaseProperties dictionary.
Definition: twoPhaseSystem.C:555
Foam::twoPhaseSystem::phi
const surfaceScalarField & phi() const
Return the mixture flux.
Definition: twoPhaseSystemI.H:76
Foam::twoPhaseSystem::phase1_
phaseModel & phase1_
Phase model 1.
Definition: twoPhaseSystem.H:81
Foam::twoPhaseSystem::correctTurbulence
void correctTurbulence()
Correct two-phase turbulence.
Definition: twoPhaseSystem.C:548
Foam::twoPhaseSystem::~twoPhaseSystem
virtual ~twoPhaseSystem()
Destructor.
Definition: twoPhaseSystem.C:76
Foam::twoPhaseSystem::Kd
tmp< volScalarField > Kd() const
Return the drag coefficient.
Definition: twoPhaseSystem.C:93
Foam::twoPhaseSystem::rho
tmp< volScalarField > rho() const
Return the mixture density.
Definition: twoPhaseSystem.C:284