Loading...
Searching...
No Matches
TomiyamaLiftForce.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) 2012-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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 "TomiyamaLiftForce.H"
29
30// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
31
32template<class CloudType>
34(
35 const typename CloudType::parcelType& p,
36 const typename CloudType::parcelType::trackingData& td,
37 const vector& curlUc,
38 const scalar Re,
39 const scalar muc
40) const
41{
42 const vector& g = this->owner().g().value();
43
44 scalar Eo = p.Eo(td, sigma_);
45 scalar dH = p.d()*cbrt(1.0 + 0.163*pow(Eo, 0.757));
46 scalar Eod = p.Eo(g, p.rho(), td.rhoc(), p.U(), dH, sigma_);
47 scalar f = 0.00105*pow3(Eod) - 0.0159*sqr(Eod) - 0.0204*Eod + 0.474;
48
49 if (Eod <= 4)
50 {
51 return min(0.288*tanh(0.121*Re), f);
52 }
53 else if ((Eod > 4) && (Eod <= 10))
54 {
55 return f;
56 }
57 else
58 {
59 return -0.27;
60 }
61}
62
63
64// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
65
66template<class CloudType>
68(
69 CloudType& owner,
70 const fvMesh& mesh,
71 const dictionary& dict,
72 const word& forceType
73)
75 LiftForce<CloudType>(owner, mesh, dict, forceType),
76 sigma_(this->coeffs().getScalar("sigma"))
77{}
78
79
80template<class CloudType>
82(
83 const TomiyamaLiftForce& lf
84)
85:
88{}
89
90
91// * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
92
93template<class CloudType>
95{}
96
97
98// ************************************************************************* //
const uniformDimensionedVectorField & g
Base class for particle lift force models.
Definition LiftForce.H:56
LiftForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict, const word &forceType)
Construct from mesh.
Definition LiftForce.C:46
const fvMesh & mesh() const noexcept
Return the mesh database.
const CloudType & owner() const noexcept
Return const access to the cloud owner.
const dictionary & coeffs() const noexcept
Return the force coefficients dictionary.
Tomiyama particle lift force model applicable to deformable bubbles.
scalar sigma_
Surface tension.
TomiyamaLiftForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict, const word &forceType=typeName)
Construct from mesh.
virtual ~TomiyamaLiftForce()
Destructor.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
dynamicFvMesh & mesh
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
DSMCCloud< dsmcParcel > CloudType
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar tanh(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
dimensionedScalar cbrt(const dimensionedScalar &ds)
labelList f(nPoints)
dictionary dict