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-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 "SyamlalOBrien.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace dragModels
34 {
35  defineTypeNameAndDebug(SyamlalOBrien, 0);
36  addToRunTimeSelectionTable(dragModel, SyamlalOBrien, 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 {}
52 
53 
54 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
55 
57 {}
58 
59 
60 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
61 
63 {
65  (
66  max(1 - interface_.dispersed(), interface_.continuous().residualAlpha())
67  );
68 
69  const volScalarField A(pow(alpha2, 4.14));
70  const volScalarField B
71  (
72  neg(alpha2 - 0.85)*(0.8*pow(alpha2, 1.28))
73  + pos0(alpha2 - 0.85)*(pow(alpha2, 2.65))
74  );
75  const volScalarField Re(interface_.Re());
76  const volScalarField Vr
77  (
78  0.5
79  *(
80  A - 0.06*Re + sqrt(sqr(0.06*Re) + 0.12*Re*(2*B - A) + sqr(A))
81  )
82  );
83  volScalarField CdsRe(sqr(0.63*sqrt(Re) + 4.8*sqrt(Vr)));
84 
85  return
86  CdsRe
87  *max(interface_.continuous(), interface_.continuous().residualAlpha())
88  /sqr(Vr);
89 }
90 
91 
92 // ************************************************************************* //
dictionary dict
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
alpha2
Definition: alphaEqn.H:115
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
SyamlalOBrien(const dictionary &dict, const phaseInterface &interface, const bool registerObject)
Construct from a dictionary and an interface.
dimensionedScalar pos0(const dimensionedScalar &ds)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
virtual ~SyamlalOBrien()
Destructor.
virtual tmp< volScalarField > CdRe() const
Drag coefficient.
A class for managing temporary objects.
Definition: PtrList.H:53
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
Namespace for OpenFOAM.