MomentumTransferPhaseSystem.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-2016 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::MomentumTransferPhaseSystem
26 
27 Description
28  Class which models interfacial momenum transfer between a number of phases.
29  Drag, virtual mass, lift, wall lubrication and turbulent dispersion are all
30  modelled. The explicit contribution from the drag is omitted from the
31  transfer matrices, as this forms part of the solution of the pressure
32  equation.
33 
34 SourceFiles
35  MomentumTransferPhaseSystem.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef MomentumTransferPhaseSystem_H
40 #define MomentumTransferPhaseSystem_H
41 
42 #include "phaseSystem.H"
43 #include "HashPtrTable.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 template<class modelType>
51 class BlendedInterfacialModel;
52 
53 class blendingMethod;
54 class dragModel;
55 class virtualMassModel;
56 class liftModel;
57 class wallLubricationModel;
58 class turbulentDispersionModel;
59 
60 /*---------------------------------------------------------------------------*\
61  Class MomentumTransferPhaseSystem Declaration
62 \*---------------------------------------------------------------------------*/
63 
64 template<class BasePhaseSystem>
66 :
67  public BasePhaseSystem
68 {
69 protected:
70 
71  // Protected typedefs
72 
73  typedef HashTable
74  <
79 
80  typedef HashTable
81  <
83  phasePairKey,
86 
87  typedef HashTable
88  <
90  phasePairKey,
93 
94  typedef HashTable
95  <
97  phasePairKey,
100 
101  typedef HashTable
102  <
104  phasePairKey,
107 
108 
109 private:
110 
111  // Private data
112 
113  //- Drag coefficients
115 
116  //- Virtual mass coefficients
118 
119  // Sub Models
120 
121  //- Drag models
122  dragModelTable dragModels_;
123 
124  //- Virtual mass models
125  virtualMassModelTable virtualMassModels_;
126 
127  //- Lift models
128  liftModelTable liftModels_;
129 
130  //- Wall lubrication models
131  wallLubricationModelTable wallLubricationModels_;
132 
133  //- Turbulent dispersion models
134  turbulentDispersionModelTable turbulentDispersionModels_;
135 
136  //- Construct element phasei of Fs if not set and return
137  // Used by Fs()
138  volVectorField& setF
139  (
141  ) const;
142 
143  //- Construct element phasei of phiDs if not set and return
144  // Used by phiDs()
145  surfaceScalarField& setPhiD
146  (
148  ) const;
149 
150 
151 public:
152 
153  // Constructors
154 
155  //- Construct from fvMesh
157 
158 
159  //- Destructor
161 
162 
163  // Member Functions
164 
165  //- Constant access to drag coefficients
166  virtual const phaseSystem::KdTable& Kds() const
167  {
168  return Kds_;
169  }
170 
171  //- Return the drag coefficient
172  virtual tmp<volScalarField> Kd(const phasePairKey& key) const;
173 
174  //- Return the face drag coefficient
175  virtual tmp<surfaceScalarField> Kdf(const phasePairKey& key) const;
176 
177  //- Return the drag coefficient for phase
178  virtual tmp<volScalarField> Kd(const phaseModel& phase) const;
179 
180  //- Return the virtual mass coefficient
181  virtual tmp<volScalarField> Vm(const phasePairKey& key) const;
182 
183  //- Return the face virtual mass coefficient
184  virtual tmp<surfaceScalarField> Vmf(const phasePairKey& key) const;
185 
186  //- Return the combined force (lift + wall-lubrication)
187  virtual tmp<volVectorField> F(const phasePairKey& key) const;
188 
189  //- Return the combined force (lift + wall-lubrication)
190  virtual autoPtr<PtrList<volVectorField>> Fs() const;
191 
192  //- Return the turbulent dispersion force on faces for phase pair
194  (
196  ) const;
197 
198  //- Return the combined face-force (lift + wall-lubrication)
199  virtual tmp<surfaceScalarField> Ff(const phasePairKey& key) const;
200 
201  //- Return the turbulent diffusivity
202  // Multiplies the phase-fraction gradient
203  virtual tmp<volScalarField> D(const phasePairKey& key) const;
204 
205  //- Return the momentum transfer matrices
207  momentumTransfer() const;
208 
209  //- Read base phaseProperties dictionary
210  virtual bool read();
211 };
212 
213 
214 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
215 
216 } // End namespace Foam
217 
218 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219 
220 #ifdef NoRepository
222 #endif
223 
224 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
225 
226 #endif
227 
228 // ************************************************************************* //
Class which models interfacial momenum transfer between a number of phases. Drag, virtual mass...
virtual tmp< volScalarField > Kd(const phasePairKey &key) const
Return the drag coefficient.
virtual tmp< surfaceScalarField > Kdf(const phasePairKey &key) const
Return the face drag coefficient.
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
virtual const phaseSystem::KdTable & Kds() const
Constant access to drag coefficients.
label phasei
Definition: pEqn.H:27
HashTable< autoPtr< BlendedInterfacialModel< dragModel > >, phasePairKey, phasePairKey::hash > dragModelTable
virtual tmp< surfaceScalarField > Ff(const phasePairKey &key) const
Return the combined face-force (lift + wall-lubrication)
virtual bool read()
Read base phaseProperties dictionary.
virtual autoPtr< PtrList< surfaceScalarField > > phiDs(const PtrList< volScalarField > &rAUs) const
Return the turbulent dispersion force on faces for phase pair.
MomentumTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
HashTable< autoPtr< BlendedInterfacialModel< liftModel > >, phasePairKey, phasePairKey::hash > liftModelTable
PtrList< volScalarField > rAUs(fluid.phases().size())
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransfer() const
Return the momentum transfer matrices.
HashTable< autoPtr< BlendedInterfacialModel< virtualMassModel > >, phasePairKey, phasePairKey::hash > virtualMassModelTable
virtual tmp< volScalarField > Vm(const phasePairKey &key) const
Return the virtual mass coefficient.
An STL-conforming hash table.
Definition: HashTable.H:61
virtual tmp< volScalarField > D(const phasePairKey &key) const
Return the turbulent diffusivity.
virtual ~MomentumTransferPhaseSystem()
Destructor.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:62
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:52
virtual tmp< volVectorField > F(const phasePairKey &key) const
Return the combined force (lift + wall-lubrication)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
virtual autoPtr< PtrList< volVectorField > > Fs() const
Return the combined force (lift + wall-lubrication)
A class for managing temporary objects.
Definition: PtrList.H:54
virtual tmp< surfaceScalarField > Vmf(const phasePairKey &key) const
Return the face virtual mass coefficient.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
HashTable< autoPtr< BlendedInterfacialModel< turbulentDispersionModel > >, phasePairKey, phasePairKey::hash > turbulentDispersionModelTable
Namespace for OpenFOAM.
HashTable< autoPtr< BlendedInterfacialModel< wallLubricationModel > >, phasePairKey, phasePairKey::hash > wallLubricationModelTable