BrownianMotionForce.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) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "BrownianMotionForce.H"
30 #include "mathematicalConstants.H"
31 #include "fundamentalConstants.H"
32 #include "demandDrivenData.H"
33 #include "turbulenceModel.H"
34 
35 using namespace Foam::constant;
36 
37 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
38 
39 template<class CloudType>
42 {
43  const objectRegistry& obr = this->owner().mesh();
44  const word turbName =
46  (
48  this->owner().U().group()
49  );
50 
51  const turbulenceModel* turb = obr.findObject<turbulenceModel>(turbName);
52 
53  if (turb)
54  {
55  return turb->k();
56  }
57 
59  << "Turbulence model not found in mesh database" << nl
60  << "Database objects include: " << obr.sortedToc()
61  << abort(FatalError);
62 
63  return nullptr;
64 }
65 
66 
67 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
68 
69 template<class CloudType>
71 (
72  CloudType& owner,
73  const fvMesh& mesh,
74  const dictionary& dict
75 )
76 :
77  ParticleForce<CloudType>(owner, mesh, dict, typeName, true),
78  rndGen_(owner.rndGen()),
79  lambda_(this->coeffs().getScalar("lambda")),
80  turbulence_(this->coeffs().getBool("turbulence")),
81  kPtr_(nullptr),
82  ownK_(false)
83 {}
84 
85 
86 template<class CloudType>
88 (
89  const BrownianMotionForce& bmf
90 )
91 :
93  rndGen_(bmf.rndGen_),
94  lambda_(bmf.lambda_),
95  turbulence_(bmf.turbulence_),
96  kPtr_(nullptr),
97  ownK_(false)
98 {}
99 
100 
101 // * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
102 
103 template<class CloudType>
105 {}
106 
107 
108 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
109 
110 template<class CloudType>
112 {
113  if (turbulence_)
114  {
115  if (store)
116  {
117  tmp<volScalarField> tk = kModel();
118  if (tk.isTmp())
119  {
120  kPtr_ = tk.ptr();
121  ownK_ = true;
122  }
123  else
124  {
125  kPtr_ = &tk();
126  ownK_ = false;
127  }
128  }
129  else
130  {
131  if (ownK_ && kPtr_)
132  {
133  deleteDemandDrivenData(kPtr_);
134  ownK_ = false;
135  }
136  }
137  }
138 }
139 
140 
141 template<class CloudType>
143 (
144  const typename CloudType::parcelType& p,
145  const typename CloudType::parcelType::trackingData& td,
146  const scalar dt,
147  const scalar mass,
148  const scalar Re,
149  const scalar muc
150 ) const
151 {
152  forceSuSp value(Zero);
153 
154  const scalar dp = p.d();
155  const scalar Tc = td.Tc();
156 
157  const scalar alpha = 2.0*lambda_/dp;
158  const scalar cc = 1.0 + alpha*(1.257 + 0.4*exp(-1.1/alpha));
159 
160  // Boltzmann constant
161  const scalar kb = physicoChemical::k.value();
162 
163  scalar f = 0;
164  if (turbulence_)
165  {
166  const label celli = p.cell();
167  const volScalarField& k = *kPtr_;
168  const scalar kc = k[celli];
169  const scalar Dp = kb*Tc*cc/(3*mathematical::pi*muc*dp);
170  f = sqrt(2.0*sqr(kc)*sqr(Tc)/(Dp*dt));
171  }
172  else
173  {
174  const scalar s0 =
175  216*muc*kb*Tc/(sqr(mathematical::pi)*pow5(dp)*sqr(p.rho())*cc);
176  f = mass*sqrt(mathematical::pi*s0/dt);
177  }
178 
179 
180  // To generate a cubic distribution (3 independent directions) :
181  // const scalar sqrt2 = sqrt(2.0);
182  // for (direction dir = 0; dir < vector::nComponents; dir++)
183  // {
184  // const scalar x = rndGen_.sample01<scalar>();
185  // const scalar eta = sqrt2*Math::erfInv(2*x - 1.0);
186  // value.Su()[dir] = f*eta;
187  // }
188 
189 
190  // To generate a spherical distribution:
191 
192  Random& rnd = this->owner().rndGen();
193 
194  const scalar theta = rnd.sample01<scalar>()*twoPi;
195  const scalar u = 2*rnd.sample01<scalar>() - 1;
196 
197  const scalar a = sqrt(1 - sqr(u));
198  const vector dir(a*cos(theta), a*sin(theta), u);
199 
200  value.Su() = f*mag(rnd.GaussNormal<scalar>())*dir;
201 
202  return value;
203 }
204 
205 
206 // ************************************************************************* //
Foam::Random
Random number generator.
Definition: Random.H:59
mathematicalConstants.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::Random::GaussNormal
Type GaussNormal()
Definition: RandomTemplates.C:49
Foam::BrownianMotionForce::BrownianMotionForce
BrownianMotionForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict)
Construct from mesh.
Definition: BrownianMotionForce.C:71
Foam::constant::atomic::group
constexpr const char *const group
Group name for atomic constants.
Definition: atomicConstants.H:52
Foam::constant::physicoChemical::k
const dimensionedScalar k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::BrownianMotionForce::calcCoupled
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.
Definition: BrownianMotionForce.C:143
Foam::constant
Different types of constants.
Definition: atomicConstants.C:38
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
Foam::BrownianMotionForce::~BrownianMotionForce
virtual ~BrownianMotionForce()
Destructor.
Definition: BrownianMotionForce.C:104
Foam::Random::sample01
Type sample01()
Return a sample whose components lie in the range [0,1].
Definition: RandomTemplates.C:36
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
Foam::BrownianMotionForce::cacheFields
virtual void cacheFields(const bool store)
Cache fields.
Definition: BrownianMotionForce.C:111
Foam::tmp::isTmp
bool isTmp() const noexcept
Identical to is_pointer()
Definition: tmp.H:295
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::forceSuSp::Su
const vector & Su() const
Return const access to the explicit contribution [kg.m/s2].
Definition: forceSuSpI.H:61
Foam::DSMCCloud::rndGen
Random & rndGen()
Return reference to the random object.
Definition: DSMCCloudI.H:124
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
Foam::compressible::turbulenceModel
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Definition: turbulentFluidThermoModel.H:63
Foam::DSMCCloud::mesh
const fvMesh & mesh() const
Return reference to the mesh.
Definition: DSMCCloudI.H:44
Foam::forceSuSp
Helper container for force Su and Sp terms.
Definition: forceSuSp.H:64
Foam::ParticleForce
Abstract base class for particle forces.
Definition: ParticleForce.H:59
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:71
Foam::constant::mathematical::twoPi
constexpr scalar twoPi(2 *M_PI)
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::pow5
dimensionedScalar pow5(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:111
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam::tmp::ptr
T * ptr() const
Return managed pointer for reuse, or clone() the object reference.
Definition: tmpI.H:255
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::BrownianMotionForce
Calculates particle Brownian motion force.
Definition: BrownianMotionForce.H:63
U
U
Definition: pEqn.H:72
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::nl
constexpr char nl
Definition: Ostream.H:404
f
labelList f(nPoints)
Foam::Vector< scalar >
fundamentalConstants.H
Fundamental dimensioned constants.
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::Re
scalarField Re(const UList< complex > &cf)
Extract real component.
Definition: complexField.C:159
Foam::DSMCCloud::parcelType
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
Foam::ParticleForce::mesh
const fvMesh & mesh() const
Return the mesh database.
Definition: ParticleForceI.H:45
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
kb
scalar kb
Definition: solveBulkSurfactant.H:7
Foam::GeometricField< scalar, fvPatchField, volMesh >
turb
compressible::turbulenceModel & turb
Definition: setRegionFluidFields.H:10
turbulenceModel.H
BrownianMotionForce.H
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265