ErgunWenYuDragForce.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2013-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "ErgunWenYuDragForce.H"
29 #include "volFields.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 template<class CloudType>
34 Foam::scalar Foam::ErgunWenYuDragForce<CloudType>::CdRe(const scalar alphacRe)
35 const
36 {
37  // (ZZB:Eq. 14, GLSLR:Table 3)
38  if (alphacRe < 1000.0)
39  {
40  return 24.0*(1.0 + 0.15*pow(alphacRe, 0.687));
41  }
42 
43  return 0.44*alphacRe;
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
49 template<class CloudType>
51 (
52  CloudType& owner,
53  const fvMesh& mesh,
54  const dictionary& dict
55 )
56 :
57  ParticleForce<CloudType>(owner, mesh, dict, typeName, true),
58  alphac_
59  (
60  this->mesh().template lookupObject<volScalarField>
61  (
62  this->coeffs().getWord("alphac")
63  )
64  )
65 {}
66 
67 
68 template<class CloudType>
70 (
72 )
73 :
75  alphac_
76  (
77  this->mesh().template lookupObject<volScalarField>
78  (
79  this->coeffs().getWord("alphac")
80  )
81  )
82 {}
83 
84 
85 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
86 
87 template<class CloudType>
89 (
90  const typename CloudType::parcelType& p,
91  const typename CloudType::parcelType::trackingData& td,
92  const scalar dt,
93  const scalar mass,
94  const scalar Re,
95  const scalar muc
96 ) const
97 {
98  const scalar alphac = alphac_[p.cell()];
99 
100  if (alphac < 0.8)
101  {
102  return forceSuSp
103  (
104  Zero,
105  (mass/p.rho())
106  *(150.0*(1.0 - alphac)/alphac + 1.75*Re)*muc/(alphac*sqr(p.d()))
107  );
108  }
109  else
110  {
111  return forceSuSp
112  (
113  Zero,
114  (mass/p.rho())
115  *0.75*CdRe(alphac*Re)*muc*pow(alphac, -2.65)/(alphac*sqr(p.d()))
116  );
117  }
118 }
119 
120 
121 // ************************************************************************* //
virtual forceSuSp calcCoupled(const typename CloudType::parcelType &p, const typename CloudType::parcelType::trackingData &td, const scalar dt, const scalar mass, const scalar Re, const scalar muc) const
Calculate the coupled force.
dictionary dict
ErgunWenYuDragForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict)
Construct from mesh.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Particle-drag model wherein drag forces (per unit carrier-fluid velocity) are dynamically computed ba...
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Abstract base class for particle forces.
Definition: ParticleForce.H:54
Helper container for force Su and Sp terms.
Definition: forceSuSp.H:60
dynamicFvMesh & mesh
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:290
scalarField Re(const UList< complex > &cmplx)
Extract real component.
Definition: complexField.C:207
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
volScalarField & p
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127