GidaspowErgunWenYuDrag.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) 2025 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 "GidaspowErgunWenYuDrag.H"
27 #include "SchillerNaumannDrag.H"
30 #include "coupledToFluid.H"
31 #include "sphericalCoupled.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace Lagrangian
38 {
41  (
45  );
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
53 Foam::Lagrangian::GidaspowErgunWenYuDrag::calcD
54 (
55  const LagrangianModelRef& model,
56  const LagrangianSubMesh& subMesh
57 ) const
58 {
59  const clouds::spherical& sCloud = cloud<clouds::spherical>();
60  const clouds::sphericalCoupled& scCloud = cloud<clouds::sphericalCoupled>();
61 
62  const LagrangianSubScalarField Re = scCloud.Re(model, subMesh);
63  const LagrangianSubScalarSubField d(sCloud.d(model, subMesh));
64 
65  const LagrangianSubScalarField alpha(min(sCloud.alpha(subMesh), alphaMax_));
66  const LagrangianSubScalarField alphac(1 - alpha);
67 
68  const LagrangianSubScalarField CdRe
69  (
70  // Use Wen-Yu at low particulate fractions (< 20%) ...
71  pos0(alphac - 0.8)
73  *pow(alphac, -2.65)
74 
75  // ... and Ergun at high particulate fractions (> 20%)
76  + neg(alphac - 0.8)*(4.0/3.0)*(150*alpha/alphac + 1.75*Re)
77  );
78 
80  <
81  clouds::coupledToIncompressibleFluid,
82  clouds::coupledToFluid
83  >();
84 
85  tmp<LagrangianSubScalarField> tmucByRhoOrMuc =
86  isCloud<clouds::coupledToIncompressibleFluid>()
87  ? (
88  cloud<clouds::coupledToIncompressibleFluid>().nuc(model, subMesh)
89  /cloud<clouds::coupledToIncompressibleFluid>().rhoByRhoc
90  )
91  : tmp<LagrangianSubScalarField>
92  (
93  cloud<clouds::coupledToFluid>().muc(model, subMesh)
94  );
95 
96  return
98  (
99  "D:" + Foam::name(subMesh.group()),
100  CdRe*(constant::mathematical::pi/8)*d*tmucByRhoOrMuc
101  );
102 }
103 
104 
105 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
106 
108 (
109  const word& name,
110  const LagrangianMesh& mesh,
111  const dictionary& modelDict,
112  const dictionary& stateDict
113 )
114 :
115  drag(name, mesh, modelDict, stateDict),
116  alphaMax_(modelDict.lookup<scalar>("alphaMax", unitFraction))
117 {}
118 
119 
120 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
static tmp< DimensionedField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Mesh &mesh, const dimensionSet &, const PrimitiveField< Type > &)
Return a temporary field constructed from name, mesh,.
Class containing Lagrangian geometry and topology.
Base class for Lagrangian models.
GidaspowErgunWenYuDrag(const word &name, const LagrangianMesh &mesh, const dictionary &modelDict, const dictionary &stateDict)
Construct from components.
static tmp< LagrangianSubScalarField > CdRe(const LagrangianSubScalarField &Re)
Return the drag coefficient times Reynold's number, as a function.
void assertCloud() const
Generate an error if the cloud is not one of the given types.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(LagrangianModel, constantCoefficientVirtualMass, dictionary)
defineTypeNameAndDebug(constantCoefficientVirtualMass, 0)
Namespace for OpenFOAM.
dimensionedScalar pos0(const dimensionedScalar &ds)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
dimensionedScalar neg(const dimensionedScalar &ds)
LagrangianSubSubField< scalar > LagrangianSubScalarSubField
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
LagrangianSubField< scalar > LagrangianSubScalarField
const unitConversion unitFraction
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97