TomiyamaLift.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2014-2018 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "TomiyamaLift.H"
29 #include "phasePair.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace liftModels
37 {
38  defineTypeNameAndDebug(TomiyamaLift, 0);
39  addToRunTimeSelectionTable(liftModel, TomiyamaLift, dictionary);
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
47 (
48  const dictionary& dict,
49  const phasePair& pair
50 )
51 :
52  liftModel(dict, pair)
53 {}
54 
55 
56 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
57 
59 {}
60 
61 
62 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
63 
65 {
66  volScalarField EoH(pair_.EoH2());
67 
69  (
70  0.0010422*pow3(EoH) - 0.0159*sqr(EoH) - 0.0204*EoH + 0.474
71  );
72 
73  return
74  neg(EoH - scalar(4))*min(0.288*tanh(0.121*pair_.Re()), f)
75  + pos0(EoH - scalar(4))*neg(EoH - scalar(10.7))*f
76  + pos0(EoH - scalar(10.7))*(-0.288);
77 }
78 
79 
80 // ************************************************************************* //
Foam::phasePair
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:53
Foam::liftModels::TomiyamaLift::~TomiyamaLift
virtual ~TomiyamaLift()
Destructor.
Definition: TomiyamaLift.C:58
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::pos0
dimensionedScalar pos0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:188
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
TomiyamaLift.H
Foam::tanh
dimensionedScalar tanh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:272
Foam::liftModels::defineTypeNameAndDebug
defineTypeNameAndDebug(constantLiftCoefficient, 0)
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::liftModels::TomiyamaLift::TomiyamaLift
TomiyamaLift(const dictionary &dict, const phasePair &pair)
Construct from a dictionary and a phase pair.
Definition: TomiyamaLift.C:47
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::liftModel
Definition: liftModel.H:55
Foam::liftModels::TomiyamaLift::Cl
virtual tmp< volScalarField > Cl() const
Lift coefficient.
Definition: TomiyamaLift.C:64
f
labelList f(nPoints)
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::liftModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(liftModel, constantLiftCoefficient, dictionary)
Foam::neg
dimensionedScalar neg(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:199