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-2020 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 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 private:
70 
71  // Private typedefs
72 
73  typedef HashPtrTable
74  <
78  > KdTable;
79 
80  typedef HashPtrTable
81  <
83  phasePairKey,
85  > KdfTable;
86 
87  typedef HashPtrTable
88  <
90  phasePairKey,
92  > VmTable;
93 
94  typedef HashPtrTable
95  <
97  phasePairKey,
99  > VmfTable;
100 
101  typedef HashTable
102  <
104  phasePairKey,
106  > dragModelTable;
107 
108  typedef HashTable
109  <
111  phasePairKey,
114 
115  typedef HashTable
116  <
118  phasePairKey,
120  > liftModelTable;
121 
122  typedef HashTable
123  <
125  phasePairKey,
128 
129  typedef HashTable
130  <
132  phasePairKey,
135 
136 
137  // Private Data
138 
139  //- Drag coefficients
140  KdTable Kds_;
141 
142  //- Face drag coefficients
143  KdfTable Kdfs_;
144 
145  //- Virtual mass coefficients
146  VmTable Vms_;
147 
148  //- Face virtual mass coefficients
149  VmfTable Vmfs_;
150 
151 
152  // Sub Models
153 
154  //- Drag models
155  dragModelTable dragModels_;
156 
157  //- Virtual mass models
158  virtualMassModelTable virtualMassModels_;
159 
160  //- Lift models
161  liftModelTable liftModels_;
162 
163  //- Wall lubrication models
164  wallLubricationModelTable wallLubricationModels_;
165 
166  //- Turbulent dispersion models
167  turbulentDispersionModelTable turbulentDispersionModels_;
168 
169 
170 protected:
171 
172  // Protected Member Functions
173 
174  //- Add momentum transfer terms which result from bulk mass transfers
175  void addDmdtUfs
176  (
177  const phaseSystem::dmdtfTable& dmdtfs,
179  );
180 
181  //- Return the drag coefficient for the phase pair
182  virtual tmp<volScalarField> Kd(const phasePairKey& key) const;
183 
184  //- Return the face drag coefficient for the phase pair
185  virtual tmp<surfaceScalarField> Kdf(const phasePairKey& key) const;
186 
187  //- Return the virtual mass coefficient for the phase pair
188  virtual tmp<volScalarField> Vm(const phasePairKey& key) const;
189 
190 
191 public:
192 
193  // Constructors
194 
195  //- Construct from fvMesh
197 
198 
199  //- Destructor
201 
202 
203  // Member Functions
204 
205  //- Return the momentum transfer matrices for the cell-based algorithm.
206  // This includes implicit and explicit forces that add into the cell
207  // UEqn in the normal way.
209 
210  //- As momentumTransfer, but for the face-based algorithm
212 
213  //- Return implicit force coefficients on the faces, for the face-based
214  // algorithm.
215  virtual PtrList<surfaceScalarField> AFfs() const;
216 
217  //- Return the explicit force fluxes for the cell-based algorithm, that
218  // do not depend on phase mass/volume fluxes, and can therefore be
219  // evaluated outside the corrector loop. This includes things like
220  // lift, turbulent dispersion, and wall lubrication.
222  (
224  );
225 
226  //- As phiFs, but for the face-based algorithm
228  (
230  );
231 
232  //- Return the explicit drag force fluxes for the cell-based algorithm.
233  // These depend on phase mass/volume fluxes, and must therefore be
234  // evaluated inside the corrector loop.
236  (
237  const PtrList<volScalarField>& rAUs
238  ) const;
239 
240  //- As phiKdPhis, but for the face-based algorithm
242  (
243  const PtrList<surfaceScalarField>& rAUfs
244  ) const;
245 
246  //- Return the explicit part of the drag force for the cell-based
247  // algorithm. This is the cell-equivalent of phiKdPhis. These depend on
248  // phase velocities, and must therefore be evaluated inside the
249  // corrector loop.
251  (
252  const PtrList<volScalarField>& rAUs
253  ) const;
254 
255  //- Solve the drag system for the velocities and fluxes
256  virtual void partialElimination
257  (
258  const PtrList<volScalarField>& rAUs,
259  const PtrList<volVectorField>& KdUByAs,
261  const PtrList<surfaceScalarField>& phiKdPhis
262  );
263 
264  //- As partialElimination, but for the face-based algorithm. Only solves
265  // for the fluxes.
266  virtual void partialEliminationf
267  (
268  const PtrList<surfaceScalarField>& rAUfs,
269  const PtrList<surfaceScalarField>& alphafs,
270  const PtrList<surfaceScalarField>& phiKdPhifs
271  );
272 
273  //- Return the flux corrections for the cell-based algorithm. These
274  // depend on phase mass/volume fluxes, and must therefore be evaluated
275  // inside the corrector loop.
277  (
278  const PtrList<volScalarField>& rAUs,
279  const bool includeVirtualMass = false
280  ) const;
281 
282  //- Returns true if the phase pressure is treated implicitly
283  // in the phase fraction equation
284  virtual bool implicitPhasePressure(const phaseModel& phase) const;
285 
286  //- Returns true if the phase pressure is treated implicitly
287  // in the phase fraction equation for any phase
288  virtual bool implicitPhasePressure() const;
289 
290  //- Return the phase diffusivity
291  // divided by the momentum central coefficient
293  (
294  const PtrList<volScalarField>& rAUs,
295  const PtrList<surfaceScalarField>& rAUfs
296  ) const;
297 
298  //- Read base phaseProperties dictionary
299  virtual bool read();
300 };
301 
302 
303 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
304 
305 } // End namespace Foam
306 
307 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
308 
309 #ifdef NoRepository
311 #endif
312 
313 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
314 
315 #endif
316 
317 // ************************************************************************* //
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.
virtual tmp< surfaceScalarField > Kdf(const phasePairKey &key) const
Return the face drag coefficient for the phase pair.
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.
MomentumTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
A HashTable specialization 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:57
Info<< "Reading field p_rgh\"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);label pRefCell=0;scalar pRefValue=0.0;if(fluid.incompressible()){ p=max(p_rgh+fluid.rho() *gh, pMin);if(p_rgh.needReference()) { setRefCell(p, p_rgh, pimple.dict(), pRefCell, pRefValue);p+=dimensionedScalar("p", p.dimensions(), pRefValue - getRefCellValue(p, pRefCell));p_rgh=p - fluid.rho() *gh;}}mesh.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 tmp< volScalarField > Vm(const phasePairKey &key) const
Return the virtual mass coefficient for the phase pair.
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:78
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
virtual tmp< volScalarField > Kd(const phasePairKey &key) const
Return the drag coefficient for the phase pair.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
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.