SyamlalOBrien.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) 2011-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 "SyamlalOBrien.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace dragModels
34 {
35  defineTypeNameAndDebug(SyamlalOBrien, 0);
36 
38  (
39  dragModel,
40  SyamlalOBrien,
41  dictionary
42  );
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const dictionary& interfaceDict,
52  const phaseModel& phase1,
53  const phaseModel& phase2
54 )
55 :
56  dragModel(interfaceDict, phase1, phase2)
57 {}
58 
59 
60 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
61 
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67 
69 (
70  const volScalarField& Ur
71 ) const
72 {
73  volScalarField alpha2(max(phase2_, scalar(1.0e-6)));
74  volScalarField A(pow(alpha2, 4.14));
76  (
77  neg(alpha2 - 0.85)*(0.8*pow(alpha2, 1.28))
78  + pos0(alpha2 - 0.85)*(pow(alpha2, 2.65))
79  );
80 
81  volScalarField Re(max(Ur*phase1_.d()/phase2_.nu(), scalar(1.0e-3)));
82 
84  (
85  0.5*
86  (
87  A - 0.06*Re + sqrt(sqr(0.06*Re) + 0.12*Re*(2.0*B - A) + sqr(A))
88  )
89  );
90 
91  volScalarField Cds(sqr(0.63 + 4.8*sqrt(Vr/Re)));
92 
93  return 0.75*Cds*phase2_.rho()*Ur/(phase1_.d()*sqr(Vr));
94 }
95 
96 
97 // ************************************************************************* //
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const phaseModel & phase1_
Definition: dragModel.H:57
dimensionedScalar sqrt(const dimensionedScalar &ds)
tmp< volScalarField > d() const
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
virtual tmp< volScalarField > K() const
Return the drag coefficient K.
const phaseModel & phase2_
Definition: dragModel.H:58
alpha2
Definition: alphaEqn.H:115
phaseModel & phase1
SyamlalOBrien(const dictionary &interfaceDict, const phaseModel &phase1, const phaseModel &phase2)
Construct from components.
dimensionedScalar pos0(const dimensionedScalar &ds)
const dimensionedScalar & nu() const
Return the laminar viscosity.
Definition: phaseModel.H:150
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionedScalar & rho() const
Definition: phaseModel.H:165
virtual ~SyamlalOBrien()
Destructor.
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:98
phaseModel & phase2
A class for managing temporary objects.
Definition: PtrList.H:53
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
Namespace for OpenFOAM.