Tenneti.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 | Copyright (C) 2016-2017 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 "Tenneti.H"
27 #include "phasePair.H"
28 #include "SchillerNaumann.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace dragModels
36 {
37  defineTypeNameAndDebug(Tenneti, 0);
38  addToRunTimeSelectionTable(dragModel, Tenneti, dictionary);
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
45 Foam::dragModels::Tenneti::Tenneti
46 (
47  const dictionary& dict,
48  const phasePair& pair,
49  const bool registerObject
50 )
51 :
52  dragModel(dict, pair, registerObject),
53  residualRe_("residualRe", dimless, dict.lookup("residualRe"))
54 {}
55 
56 
57 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
58 
59 Foam::dragModels::Tenneti::~Tenneti()
60 {}
61 
62 
63 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
64 
65 Foam::tmp<Foam::volScalarField> Foam::dragModels::Tenneti::CdRe() const
66 {
68  (
70  );
71 
73  (
74  max(scalar(1) - pair_.dispersed(), pair_.continuous().residualAlpha())
75  );
76 
78 
79  volScalarField CdReIsolated
80  (
81  neg(Res - 1000)*24.0*(1.0 + 0.15*pow(Res, 0.687))
82  + pos0(Res - 1000)*0.44*max(Res, residualRe_)
83  );
84 
86  (
87  5.81*alpha1/pow3(alpha2) + 0.48*pow(alpha1, 1.0/3.0)/pow4(alpha2)
88  );
89 
91  (
92  pow3(alpha1)*Res*(0.95 + 0.61*pow3(alpha1)/sqr(alpha2))
93  );
94 
95  // Tenneti et al. correlation includes the mean pressure drag.
96  // This was removed here by multiplying F by alpha2 for consistency with
97  // the formulation used in OpenFOAM
98  return
99  CdReIsolated + 24.0*sqr(alpha2)*(F0 + F1);
100 }
101 
102 
103 // ************************************************************************* //
dictionary dict
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define F1(B, C, D)
Definition: SHA1.C:173
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar neg(const dimensionedScalar &ds)
const phasePair & pair_
Phase pair.
Definition: dragModel.H:62
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
stressControl lookup("compactNormalStress") >> compactNormalStress
virtual const phaseModel & continuous() const
Continuous phase.
tmp< volScalarField > Re() const
Reynolds number.
alpha2
Definition: alphaEqn.H:112
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
virtual const phaseModel & dispersed() const
Dispersed phase.
dimensionedScalar pos0(const dimensionedScalar &ds)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
volScalarField & alpha1
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensionedScalar pow4(const dimensionedScalar &ds)
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.