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-2018 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  Class which solves the volume fraction equations for two phases.
31 
32 SourceFiles
33  twoPhaseSystem.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef twoPhaseSystem_H
38 #define twoPhaseSystem_H
39 
40 #include "phaseSystem.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 class dragModel;
48 class virtualMassModel;
49 
50 /*---------------------------------------------------------------------------*\
51  Class twoPhaseSystem Declaration
52 \*---------------------------------------------------------------------------*/
53 
54 class twoPhaseSystem
55 :
56  public phaseSystem
57 {
58 private:
59 
60  // Private member functions
61 
62  //- Return the drag coefficient for phase pair
63  virtual tmp<volScalarField> Kd(const phasePairKey& key) const = 0;
64 
65  //- Return the face drag coefficient for phase pair
66  virtual tmp<surfaceScalarField> Kdf(const phasePairKey& key) const = 0;
67 
68  //- Return the virtual mass coefficient for phase pair
69  virtual tmp<volScalarField> Vm(const phasePairKey& key) const = 0;
70 
71 
72 protected:
73 
74  // Protected data
75 
76  //- Flag to indicate that returned lists of fields are "complete"; i.e.,
77  // that an absence of force is returned as a constructed list of zeros,
78  // rather than a null pointer
79  static const bool fillFields_ = true;
80 
81  //- Phase model 1
83 
84  //- Phase model 2
86 
87 
88 public:
89 
90  //- Runtime type information
91  TypeName("twoPhaseSystem");
92 
93  // Declare runtime construction
94 
96  (
97  autoPtr,
99  dictionary,
100  (
101  const fvMesh& mesh
102  ),
103  (mesh)
104  );
105 
106 
107  // Constructors
108 
109  //- Construct from fvMesh
110  twoPhaseSystem(const fvMesh&);
111 
112 
113  //- Destructor
114  virtual ~twoPhaseSystem();
115 
116 
117  // Selectors
118 
120  (
121  const fvMesh& mesh
122  );
123 
124 
125  // Member Functions
126 
127  using phaseSystem::sigma;
128  using phaseSystem::dmdts;
129 
130  //- Return phase model 1
131  const phaseModel& phase1() const;
132 
133  //- Access phase model 1
134  phaseModel& phase1();
135 
136  //- Return phase model 2
137  const phaseModel& phase2() const;
138 
139  //- Access phase model 2
140  phaseModel& phase2();
141 
142  //- Return the phase not given as an argument
143  const phaseModel& otherPhase(const phaseModel& phase) const;
144 
145  //- Return the surface tension coefficient
146  tmp<volScalarField> sigma() const;
147 
148  //- Return the drag coefficient
149  tmp<volScalarField> Kd() const;
150 
151  //- Return the face drag coefficient
153 
154  //- Return the virtual mass coefficient
155  tmp<volScalarField> Vm() const;
156 
157  //- Solve for the phase fractions
158  virtual void solve();
159 };
160 
161 
162 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163 
164 } // End namespace Foam
165 
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167 
168 #include "twoPhaseSystemI.H"
169 
170 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171 
172 #endif
173 
174 // ************************************************************************* //
Foam::phaseModel
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:54
Foam::twoPhaseSystem::New
static autoPtr< twoPhaseSystem > New(const fvMesh &mesh)
Definition: twoPhaseSystemNew.C:35
Foam::twoPhaseSystem::declareRunTimeSelectionTable
declareRunTimeSelectionTable(autoPtr, twoPhaseSystem, dictionary,(const fvMesh &mesh),(mesh))
Foam::twoPhaseSystem::mesh
const fvMesh & mesh() const
Return the mesh.
Definition: twoPhaseSystemI.H:30
twoPhaseSystemI.H
Foam::twoPhaseSystem::phase1
const phaseModel & phase1() const
Return phase model 1.
Definition: twoPhaseSystemI.H:30
Foam::twoPhaseSystem
Class which solves the volume fraction equations for two phases.
Definition: twoPhaseSystem.H:53
Foam::phase
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:54
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::phaseSystem::sigma
tmp< volScalarField > sigma(const phasePairKey &key) const
Return the surface tension coefficient for a pair.
Definition: phaseSystem.C:319
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::TypeName
TypeName("twoPhaseSystem")
Runtime type information.
Foam::twoPhaseSystem::phase2
const phaseModel & phase2() const
Return phase model 2.
Definition: twoPhaseSystemI.H:42
Foam::phaseSystem::dmdts
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
Definition: phaseSystem.C:351
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::fillFields_
static const bool fillFields_
Flag to indicate that returned lists of fields are "complete"; i.e.,.
Definition: twoPhaseSystem.H:78
Foam::phasePairKey
An ordered or unorder pair of phase names. Typically specified as follows.
Definition: phasePairKey.H:65
Foam::twoPhaseSystem::solve
virtual void solve()
Solve for the phase fractions.
Definition: twoPhaseSystem.C:122
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::twoPhaseSystem::twoPhaseSystem
twoPhaseSystem(const fvMesh &)
Construct from fvMesh.
Definition: twoPhaseSystem.C:59
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::twoPhaseSystem::phase2_
phaseModel & phase2_
Phase model 2.
Definition: twoPhaseSystem.H:84
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::twoPhaseSystem::phase1_
phaseModel & phase1_
Phase model 1.
Definition: twoPhaseSystem.H:81
Foam::twoPhaseSystem::~twoPhaseSystem
virtual ~twoPhaseSystem()
Destructor.
Definition: twoPhaseSystem.C:76
Foam::phaseSystem
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:66
Foam::twoPhaseSystem::Kd
tmp< volScalarField > Kd() const
Return the drag coefficient.
Definition: twoPhaseSystem.C:93