MomentumTransferPhaseSystem.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2015-2023 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 momentum 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 class blendingMethod;
51 class blendedDragModel;
52 class blendedVirtualMassModel;
53 class blendedLiftModel;
54 class blendedWallLubricationModel;
55 class blendedTurbulentDispersionModel;
56 
57 /*---------------------------------------------------------------------------*\
58  Class MomentumTransferPhaseSystem Declaration
59 \*---------------------------------------------------------------------------*/
60 
61 template<class BasePhaseSystem>
63 :
64  public BasePhaseSystem
65 {
66  // Private typedefs
67 
68  typedef HashPtrTable
69  <
73  > KdTable;
74 
75  typedef HashPtrTable
76  <
80  > KdfTable;
81 
82  typedef HashPtrTable
83  <
87  > VmTable;
88 
89  typedef HashTable
90  <
95 
96  typedef HashTable
97  <
102 
103  typedef HashTable
104  <
108  > liftModelTable;
109 
110  typedef HashTable
111  <
116 
117  typedef HashTable
118  <
123 
124 
125  // Private Data
126 
127  //- Drag coefficients
128  KdTable Kds_;
129 
130  //- Face drag coefficients
131  KdfTable Kdfs_;
132 
133  //- Virtual mass coefficients
134  VmTable Vms_;
135 
136 
137  // Sub Models
138 
139  //- Drag models
140  dragModelTable dragModels_;
141 
142  //- Virtual mass models
143  virtualMassModelTable virtualMassModels_;
144 
145  //- Lift models
146  liftModelTable liftModels_;
147 
148  //- Wall lubrication models
149  wallLubricationModelTable wallLubricationModels_;
150 
151  //- Turbulent dispersion models
152  turbulentDispersionModelTable turbulentDispersionModels_;
153 
154 
155 protected:
156 
157  // Protected Member Functions
158 
159  //- Add momentum transfer terms which result from bulk mass transfers
160  void addDmdtUfs
161  (
162  const phaseSystem::dmdtfTable& dmdtfs,
164  );
165 
166  void addTmpField
167  (
168  tmp<surfaceScalarField>& result,
169  const tmp<surfaceScalarField>& field
170  ) const;
171 
172 
173 public:
174 
175  // Constructors
176 
177  //- Construct from fvMesh
179 
180 
181  //- Destructor
183 
184 
185  // Member Functions
186 
187  //- Return the momentum transfer matrices for the cell-based algorithm.
188  // This includes implicit and explicit forces that add into the cell
189  // UEqn in the normal way.
191 
192  //- As momentumTransfer, but for the face-based algorithm
194 
195  //- Return implicit force coefficients on the faces, for the face-based
196  // algorithm.
197  virtual PtrList<surfaceScalarField> KdVmfs() const;
198 
199  //- Return the explicit force fluxes for the cell-based algorithm, that
200  // do not depend on phase mass/volume fluxes, and can therefore be
201  // evaluated outside the corrector loop. This includes things like
202  // lift, turbulent dispersion, and wall lubrication.
203  virtual PtrList<surfaceScalarField> Fs() const;
204 
205  //- As Fs, but for the face-based algorithm
206  virtual PtrList<surfaceScalarField> Ffs() const;
207 
208  //- Return the explicit drag force fluxes for the cell-based algorithm.
209  // These depend on phase mass/volume fluxes, and must therefore be
210  // evaluated inside the corrector loop.
211  virtual PtrList<surfaceScalarField> KdPhis() const;
212 
213  //- As KdPhis, but for the face-based algorithm
214  virtual PtrList<surfaceScalarField> KdPhifs() const;
215 
216  //- Return the implicit part of the drag force
217  virtual PtrList<volScalarField> Kds() const;
218 
219  //- Return the explicit part of the drag force for the cell-based
220  // algorithm. This is the cell-equivalent of KdPhis. These depend on
221  // phase velocities, and must therefore be evaluated inside the
222  // corrector loop.
223  virtual PtrList<volVectorField> KdUs() const;
224 
225  //- Returns true if the phase pressure is treated implicitly
226  // in the phase fraction equation
227  virtual bool implicitPhasePressure(const phaseModel& phase) const;
228 
229  //- Returns true if the phase pressure is treated implicitly
230  // in the phase fraction equation for any phase
231  virtual bool implicitPhasePressure() const;
232 
233  //- Return the phase diffusivity
234  // divided by the momentum central coefficient
236  (
237  const PtrList<volScalarField>& rAUs,
238  const PtrList<surfaceScalarField>& rAUfs
239  ) const;
240 
241  //- Return the flux corrections for the cell-based algorithm. These
242  // depend on phase mass/volume fluxes, and must therefore be evaluated
243  // inside the corrector loop.
244  virtual PtrList<surfaceScalarField> ddtCorrs() const;
245 
246  //- Set the cell and faces drag correction fields
247  virtual void dragCorrs
248  (
250  PtrList<surfaceScalarField>& dragCorrf
251  ) const;
252 
253  //- Solve the drag system for the velocities and fluxes
254  virtual void partialElimination
255  (
256  const PtrList<volScalarField>& rAUs,
258  const PtrList<surfaceScalarField>& alphafs,
259  const PtrList<surfaceScalarField>& rAUfs,
261  );
262 
263  //- As partialElimination, but for the face-based algorithm. Only solves
264  // for the fluxes.
265  virtual void partialEliminationf
266  (
267  const PtrList<surfaceScalarField>& rAUfs,
268  const PtrList<surfaceScalarField>& alphafs,
270  );
271 
272  //- Read base phaseProperties dictionary
273  virtual bool read();
274 };
275 
276 
277 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
278 
279 } // End namespace Foam
280 
281 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
282 
283 #ifdef NoRepository
285 #endif
286 
287 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
288 
289 #endif
290 
291 // ************************************************************************* //
A HashTable specialisation for hashing pointers.
Definition: HashPtrTable.H:67
An STL-conforming hash table.
Definition: HashTable.H:127
Class which models interfacial momentum transfer between a number of phases. Drag,...
virtual void partialElimination(const PtrList< volScalarField > &rAUs, const PtrList< volVectorField > &KdUs, const PtrList< surfaceScalarField > &alphafs, const PtrList< surfaceScalarField > &rAUfs, const PtrList< surfaceScalarField > &KdPhis)
Solve the drag system for the velocities and fluxes.
virtual void dragCorrs(PtrList< volVectorField > &dragCorrs, PtrList< surfaceScalarField > &dragCorrf) const
Set the cell and faces drag correction fields.
virtual PtrList< surfaceScalarField > ddtCorrs() const
Return the flux corrections for the cell-based algorithm. These.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransfer()
Return the momentum transfer matrices for the cell-based algorithm.
virtual tmp< surfaceScalarField > alphaDByAf(const PtrList< volScalarField > &rAUs, const PtrList< surfaceScalarField > &rAUfs) const
Return the phase diffusivity.
virtual PtrList< volVectorField > KdUs() const
Return the explicit part of the drag force for the cell-based.
virtual void partialEliminationf(const PtrList< surfaceScalarField > &rAUfs, const PtrList< surfaceScalarField > &alphafs, const PtrList< surfaceScalarField > &KdPhifs)
As partialElimination, but for the face-based algorithm. Only solves.
virtual PtrList< volScalarField > Kds() const
Return the implicit part of the drag force.
void addDmdtUfs(const phaseSystem::dmdtfTable &dmdtfs, phaseSystem::momentumTransferTable &eqns)
Add momentum transfer terms which result from bulk mass transfers.
virtual PtrList< surfaceScalarField > Ffs() const
As Fs, but for the face-based algorithm.
virtual bool implicitPhasePressure() const
Returns true if the phase pressure is treated implicitly.
virtual PtrList< surfaceScalarField > Fs() const
Return the explicit force fluxes for the cell-based algorithm, that.
virtual PtrList< surfaceScalarField > KdVmfs() const
Return implicit force coefficients on the faces, for the face-based.
MomentumTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransferf()
As momentumTransfer, but for the face-based algorithm.
virtual PtrList< surfaceScalarField > KdPhifs() const
As KdPhis, but for the face-based algorithm.
virtual PtrList< surfaceScalarField > KdPhis() const
Return the explicit drag force fluxes for the cell-based algorithm.
virtual bool read()
Read base phaseProperties dictionary.
void addTmpField(tmp< surfaceScalarField > &result, const tmp< surfaceScalarField > &field) const
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Word-pair based class used for keying interface models in hash tables.
A class for managing temporary objects.
Definition: tmp.H:55
Namespace for OpenFOAM.
SurfaceField< scalar > surfaceScalarField
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61