TomiyamaDragForce.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) 2024 OpenCFD Ltd.
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 "TomiyamaDragForce.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 template<class CloudType>
35 {
36  { contaminationType::PURE, "pure" },
37  { contaminationType::SLIGHT, "slight" },
38  { contaminationType::FULL, "full" },
39 };
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 template<class CloudType>
45 Foam::scalar Foam::TomiyamaDragForce<CloudType>::CdRe(const scalar Re) const
46 {
47  const scalar f = 1 + 0.15*pow(Re, 0.687);
48 
49  switch (contaminationType_)
50  {
51  case contaminationType::PURE:
52  {
53  // Eq. 31 pure system
54  return min(16*f, 48);
55  }
56  case contaminationType::SLIGHT:
57  {
58  // Eq. 32 slightly contaminated system
59  return min(24*f, 72);
60  }
61  case contaminationType::FULL:
62  {
63  // Eq. 33 fully contaminated system
64  return 24*f;
65  }
66  default:
67  {}
68  }
69 
70  return Zero;
71 }
72 
73 
74 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
75 
76 template<class CloudType>
78 (
79  CloudType& owner,
80  const fvMesh& mesh,
81  const dictionary& dict
82 )
83 :
84  ParticleForce<CloudType>(owner, mesh, dict, typeName, true),
85  sigma_(this->coeffs().getScalar("sigma")),
86  contaminationType_
87  (
88  contaminationTypeNames.get("contamination", this->coeffs())
89  )
90 {}
91 
92 
93 template<class CloudType>
95 (
97 )
98 :
100  sigma_(df.sigma_),
101  contaminationType_(df.contaminationType_)
102 {}
103 
104 
105 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106 
107 template<class CloudType>
109 (
110  const typename CloudType::parcelType& p,
111  const typename CloudType::parcelType::trackingData& td,
112  const scalar dt,
113  const scalar mass,
114  const scalar Re,
115  const scalar muc
116 ) const
117 {
118  const scalar Eo = p.Eo(td, sigma_);
119  const scalar CdRe = max(this->CdRe(Re), Re*8/3*Eo/(Eo + 4));
120 
121  return forceSuSp(Zero, mass*0.75*muc*CdRe/(p.rho()*sqr(p.d())));
122 }
123 
124 
125 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Abstract base class for particle forces.
Definition: ParticleForce.H:54
TomiyamaDragForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict)
Construct from mesh.
Helper container for force Su and Sp terms.
Definition: forceSuSp.H:60
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.
dynamicFvMesh & mesh
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:290
labelList f(nPoints)
scalarField Re(const UList< complex > &cmplx)
Extract real component.
Definition: complexField.C:207
Particle-drag model wherein drag forces (per unit carrier-fluid velocity) are dynamically computed us...
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
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