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-2022 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:
67 
68  // Private typedefs
69 
70  typedef HashPtrTable
71  <
75  > KdTable;
76 
77  typedef HashPtrTable
78  <
80  phaseInterfaceKey,
82  > KdfTable;
83 
84  typedef HashPtrTable
85  <
87  phaseInterfaceKey,
89  > VmTable;
90 
91  typedef HashPtrTable
92  <
94  phaseInterfaceKey,
96  > VmfTable;
97 
98  typedef HashTable
99  <
101  phaseInterfaceKey,
103  > dragModelTable;
104 
105  typedef HashTable
106  <
108  phaseInterfaceKey,
111 
112  typedef HashTable
113  <
115  phaseInterfaceKey,
117  > liftModelTable;
118 
119  typedef HashTable
120  <
122  phaseInterfaceKey,
125 
126  typedef HashTable
127  <
129  phaseInterfaceKey,
132 
133 
134  // Private Data
135 
136  //- Drag coefficients
137  KdTable Kds_;
138 
139  //- Face drag coefficients
140  KdfTable Kdfs_;
141 
142  //- Virtual mass coefficients
143  VmTable Vms_;
144 
145  //- Face virtual mass coefficients
146  VmfTable Vmfs_;
147 
148 
149  // Sub Models
150 
151  //- Drag models
152  dragModelTable dragModels_;
153 
154  //- Virtual mass models
155  virtualMassModelTable virtualMassModels_;
156 
157  //- Lift models
158  liftModelTable liftModels_;
159 
160  //- Wall lubrication models
161  wallLubricationModelTable wallLubricationModels_;
162 
163  //- Turbulent dispersion models
164  turbulentDispersionModelTable turbulentDispersionModels_;
165 
166 
167 protected:
168 
169  // Protected Member Functions
170 
171  //- Add momentum transfer terms which result from bulk mass transfers
172  void addDmdtUfs
173  (
174  const phaseSystem::dmdtfTable& dmdtfs,
176  );
177 
178 
179 public:
180 
181  // Constructors
182 
183  //- Construct from fvMesh
185 
186 
187  //- Destructor
189 
190 
191  // Member Functions
192 
193  //- Return the momentum transfer matrices for the cell-based algorithm.
194  // This includes implicit and explicit forces that add into the cell
195  // UEqn in the normal way.
197 
198  //- As momentumTransfer, but for the face-based algorithm
200 
201  //- Return implicit force coefficients on the faces, for the face-based
202  // algorithm.
203  virtual PtrList<surfaceScalarField> AFfs() const;
204 
205  //- Return the explicit force fluxes for the cell-based algorithm, that
206  // do not depend on phase mass/volume fluxes, and can therefore be
207  // evaluated outside the corrector loop. This includes things like
208  // lift, turbulent dispersion, and wall lubrication.
210  (
212  );
213 
214  //- As phiFs, but for the face-based algorithm
216  (
218  );
219 
220  //- Return the explicit drag force fluxes for the cell-based algorithm.
221  // These depend on phase mass/volume fluxes, and must therefore be
222  // evaluated inside the corrector loop.
224  (
225  const PtrList<volScalarField>& rAUs
226  ) const;
227 
228  //- As phiKdPhis, but for the face-based algorithm
230  (
231  const PtrList<surfaceScalarField>& rAUfs
232  ) const;
233 
234  //- Return the explicit part of the drag force for the cell-based
235  // algorithm. This is the cell-equivalent of phiKdPhis. These depend on
236  // phase velocities, and must therefore be evaluated inside the
237  // corrector loop.
239  (
240  const PtrList<volScalarField>& rAUs
241  ) const;
242 
243  //- Solve the drag system for the velocities and fluxes
244  virtual void partialElimination
245  (
246  const PtrList<volScalarField>& rAUs,
247  const PtrList<volVectorField>& KdUByAs,
249  const PtrList<surfaceScalarField>& phiKdPhis
250  );
251 
252  //- As partialElimination, but for the face-based algorithm. Only solves
253  // for the fluxes.
254  virtual void partialEliminationf
255  (
256  const PtrList<surfaceScalarField>& rAUfs,
257  const PtrList<surfaceScalarField>& alphafs,
258  const PtrList<surfaceScalarField>& phiKdPhifs
259  );
260 
261  //- Return the flux corrections for the cell-based algorithm. These
262  // depend on phase mass/volume fluxes, and must therefore be evaluated
263  // inside the corrector loop.
265  (
266  const PtrList<volScalarField>& rAUs,
267  const bool includeVirtualMass = false
268  ) const;
269 
270  //- Returns true if the phase pressure is treated implicitly
271  // in the phase fraction equation
272  virtual bool implicitPhasePressure(const phaseModel& phase) const;
273 
274  //- Returns true if the phase pressure is treated implicitly
275  // in the phase fraction equation for any phase
276  virtual bool implicitPhasePressure() const;
277 
278  //- Return the phase diffusivity
279  // divided by the momentum central coefficient
281  (
282  const PtrList<volScalarField>& rAUs,
283  const PtrList<surfaceScalarField>& rAUfs
284  ) const;
285 
286  //- Read base phaseProperties dictionary
287  virtual bool read();
288 };
289 
290 
291 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
292 
293 } // End namespace Foam
294 
295 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
296 
297 #ifdef NoRepository
299 #endif
300 
301 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
302 
303 #endif
304 
305 // ************************************************************************* //
Class which models interfacial momentum transfer between a number of phases. Drag, virtual mass, lift, wall lubrication and turbulent dispersion are all modelled. The explicit contribution from the drag is omitted from the transfer matrices, as this forms part of the solution of the pressure equation.
virtual PtrList< surfaceScalarField > phiFs(const PtrList< volScalarField > &rAUs)
Return the explicit force fluxes for the cell-based algorithm, that.
PtrList< surfaceScalarField > alphafs(phases.size())
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransferf()
As momentumTransfer, but for the face-based algorithm.
virtual bool read()
Read base phaseProperties dictionary.
Word-pair based class used for keying interface models in hash tables.
MomentumTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
A HashTable specialisation for hashing pointers.
Definition: HashPtrTable.H:50
virtual bool implicitPhasePressure() const
Returns true if the phase pressure is treated implicitly.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
Info<< "Reading field p_rgh\"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);if(fluid.found("pMin")){ IOWarningInFunction(fluid)<< "Pressure limits, pMin and pMax, are now read from "<< pimple.dict().name()<< endl;}pressureReference pressureReference(p, p_rgh, pimple.dict(), fluid.incompressible());if(fluid.incompressible()){ p=p_rgh+fluid.rho() *gh;}if(p_rgh.needReference() &&fluid.incompressible()){ p+=dimensionedScalar("p", p.dimensions(), pressureReference.refValue() - getRefCellValue(p, pressureReference.refCell()));}p_rgh=p - fluid.rho() *gh;mesh.schemes().setFluxRequired(p_rgh.name());PtrList< volScalarField > rAUs
Definition: createFields.H:67
virtual void partialEliminationf(const PtrList< surfaceScalarField > &rAUfs, const PtrList< surfaceScalarField > &alphafs, const PtrList< surfaceScalarField > &phiKdPhifs)
As partialElimination, but for the face-based algorithm. Only solves.
virtual PtrList< volVectorField > KdUByAs(const PtrList< volScalarField > &rAUs) const
Return the explicit part of the drag force for the cell-based.
An STL-conforming hash table.
Definition: HashTable.H:61
virtual PtrList< surfaceScalarField > phiKdPhis(const PtrList< volScalarField > &rAUs) const
Return the explicit drag force fluxes for the cell-based algorithm.
PtrList< surfaceScalarField > rAUfs
Definition: createFields.H:68
virtual PtrList< surfaceScalarField > ddtCorrByAs(const PtrList< volScalarField > &rAUs, const bool includeVirtualMass=false) const
Return the flux corrections for the cell-based algorithm. These.
virtual ~MomentumTransferPhaseSystem()
Destructor.
virtual PtrList< surfaceScalarField > phiFfs(const PtrList< surfaceScalarField > &rAUfs)
As phiFs, but for the face-based algorithm.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
virtual void partialElimination(const PtrList< volScalarField > &rAUs, const PtrList< volVectorField > &KdUByAs, const PtrList< surfaceScalarField > &alphafs, const PtrList< surfaceScalarField > &phiKdPhis)
Solve the drag system for the velocities and fluxes.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:52
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
virtual PtrList< surfaceScalarField > AFfs() const
Return implicit force coefficients on the faces, for the face-based.
virtual PtrList< surfaceScalarField > DByAfs(const PtrList< volScalarField > &rAUs, const PtrList< surfaceScalarField > &rAUfs) const
Return the phase diffusivity.
void addDmdtUfs(const phaseSystem::dmdtfTable &dmdtfs, phaseSystem::momentumTransferTable &eqns)
Add momentum transfer terms which result from bulk mass transfers.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
Namespace for OpenFOAM.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransfer()
Return the momentum transfer matrices for the cell-based algorithm.
virtual PtrList< surfaceScalarField > phiKdPhifs(const PtrList< surfaceScalarField > &rAUfs) const
As phiKdPhis, but for the face-based algorithm.