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-2022 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"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace dragModels
34 {
35  defineTypeNameAndDebug(Beetstra, 0);
36  addToRunTimeSelectionTable(dragModel, Beetstra, dictionary);
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 (
45  const dictionary& dict,
46  const phaseInterface& interface,
47  const bool registerObject
48 )
49 :
50  dispersedDragModel(dict, interface, registerObject),
51  residualRe_("residualRe", dimless, dict.lookup("residualRe"))
52 {}
53 
54 
55 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
56 
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
62 
64 {
66  (
67  max(interface_.dispersed(), interface_.continuous().residualAlpha())
68  );
69 
71  (
72  max(1 - interface_.dispersed(), interface_.continuous().residualAlpha())
73  );
74 
75  const volScalarField Res(alpha2*interface_.Re());
76 
77  const volScalarField ResLim
78  (
79  "ReLim",
80  max(Res, residualRe_)
81  );
82 
83  const volScalarField F0
84  (
85  "F0",
86  10*alpha1/sqr(alpha2) + sqr(alpha2)*(1 + 1.5*sqrt(alpha1))
87  );
88 
89  const volScalarField F1
90  (
91  "F1",
92  0.413*Res/(24*sqr(alpha2))*(1.0/alpha2
93  + 3*alpha1*alpha2 + 8.4*pow(ResLim, -0.343))
94  /(1 + pow(10.0, 3*alpha1)*pow(ResLim, -(1 + 4*alpha1)/2.0))
95  );
96 
97  return 24*alpha2*(F0 + F1);
98 }
99 
100 
101 // ************************************************************************* //
dictionary dict
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
#define F1(B, C, D)
Definition: SHA1.C:173
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
volScalarField & alpha1(mixture.alpha1())
const dimensionSet dimless
Beetstra(const dictionary &dict, const phaseInterface &interface, const bool registerObject)
Construct from a dictionary and an interface.
Macros for easy insertion into run-time selection tables.
alpha2
Definition: alphaEqn.H:115
virtual ~Beetstra()
Destructor.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
stressControl lookup("compactNormalStress") >> compactNormalStress
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
virtual tmp< volScalarField > CdRe() const
Drag coefficient.
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.