TomiyamaLift.C
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) 2014-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 \*---------------------------------------------------------------------------*/
25 
26 #include "TomiyamaLift.H"
27 #include "aspectRatioModel.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace liftModels
35 {
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
45 (
46  const dictionary& dict,
47  const phaseInterface& interface
48 )
49 :
50  dispersedLiftModel(dict, interface),
51  aspectRatio_(aspectRatioModel::New(dict.subDict("aspectRatio"), interface))
52 {}
53 
54 
55 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
56 
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
62 
64 {
65  const volScalarField EoH
66  (
67  interface_.Eo(interface_.dispersed().d()/cbrt(aspectRatio_->E()))
68  );
69 
70  const volScalarField f
71  (
72  0.0010422*pow3(EoH) - 0.0159*sqr(EoH) - 0.0204*EoH + 0.474
73  );
74 
75  return
76  neg(EoH - 4)*min(0.288*tanh(0.121*interface_.Re()), f)
77  + pos0(EoH - 4)*neg(EoH - 10.7)*f
78  + pos0(EoH - 10.7)*(-0.288);
79 }
80 
81 
82 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
Model for deviations in the shape of the dispersed phase from spherical. Just a sub-model modifier,...
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Model for the lift force between two phases.
Definition: liftModel.H:53
Lift model of Tomiyama et al.
Definition: TomiyamaLift.H:67
virtual tmp< volScalarField > Cl() const
Lift coefficient.
Definition: TomiyamaLift.C:63
virtual ~TomiyamaLift()
Destructor.
Definition: TomiyamaLift.C:57
TomiyamaLift(const dictionary &dict, const phaseInterface &interface)
Construct from a dictionary and an interface.
Definition: TomiyamaLift.C:45
Class to represent an interface between phases. Derivations can further specify the configuration of ...
A class for managing temporary objects.
Definition: tmp.H:55
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
addToRunTimeSelectionTable(liftModel, constantLiftCoefficient, dictionary)
defineTypeNameAndDebug(constantLiftCoefficient, 0)
Namespace for OpenFOAM.
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar tanh(const dimensionedScalar &ds)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar neg(const dimensionedScalar &ds)
dimensionedScalar cbrt(const dimensionedScalar &ds)
labelList f(nPoints)
dictionary dict