Beetstra.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) 2016-2018 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 "Beetstra.H"
27 #include "phasePair.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace dragModels
35 {
36  defineTypeNameAndDebug(Beetstra, 0);
37  addToRunTimeSelectionTable(dragModel, Beetstra, dictionary);
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
45 (
46  const dictionary& dict,
47  const phasePair& pair,
48  const bool registerObject
49 )
50 :
51  dragModel(dict, pair, registerObject),
52  residualRe_("residualRe", dimless, dict.lookup("residualRe"))
53 {}
54 
55 
56 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
57 
59 {}
60 
61 
62 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
63 
65 {
67  (
69  );
70 
72  (
73  max(scalar(1) - pair_.dispersed(), pair_.continuous().residualAlpha())
74  );
75 
77 
78  volScalarField ResLim
79  (
80  "ReLim",
81  max(Res, residualRe_)
82  );
83 
85  (
86  "F0",
87  10*alpha1/sqr(alpha2) + sqr(alpha2)*(1 + 1.5*sqrt(alpha1))
88  );
89 
91  (
92  "F1",
93  0.413*Res/(24*sqr(alpha2))*(1.0/alpha2
94  + 3*alpha1*alpha2 + 8.4*pow(ResLim, -0.343))
95  /(1 + pow(10.0, 3*alpha1)*pow(ResLim, -(1 + 4*alpha1)/2.0))
96  );
97 
98  return 24*alpha2*(F0 + F1);
99 }
100 
101 
102 // ************************************************************************* //
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 sqrt(const dimensionedScalar &ds)
const phasePair & pair_
Phase pair.
Definition: dragModel.H:62
Macros for easy insertion into run-time selection tables.
virtual ~Beetstra()
Destructor.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Beetstra(const dictionary &dict, const phasePair &pair, const bool registerObject)
Construct from a dictionary and a phase pair.
stressControl lookup("compactNormalStress") >> compactNormalStress
virtual const phaseModel & continuous() const
Continuous phase.
const volScalarField & alpha1
tmp< volScalarField > Re() const
Reynolds number.
alpha2
Definition: alphaEqn.H:115
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
virtual const phaseModel & dispersed() const
Dispersed phase.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
virtual tmp< volScalarField > CdRe() const
Drag coefficient.
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.