BrownianMotionForce.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 "BrownianMotionForce.H"
27 #include "mathematicalConstants.H"
28 #include "fundamentalConstants.H"
29 #include "demandDrivenData.H"
30 #include "momentumTransportModel.H"
31 
32 using namespace Foam::constant;
33 
34 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
35 
36 template<class CloudType>
37 Foam::scalar Foam::BrownianMotionForce<CloudType>::erfInv(const scalar y) const
38 {
39  const scalar a = 0.147;
40  scalar k = 2.0/(mathematical::pi*a) + 0.5*log(1.0 - y*y);
41  scalar h = log(1.0 - y*y)/a;
42  scalar x = sqrt(-k + sqrt(k*k - h));
43 
44  if (y < 0.0)
45  {
46  return -x;
47  }
48  else
49  {
50  return x;
51  }
52 }
53 
54 
55 template<class CloudType>
58 {
59  const objectRegistry& obr = this->owner().mesh();
60 
61  if (obr.foundType<momentumTransportModel>(this->owner().U().group()))
62  {
63  const momentumTransportModel& model =
64  obr.lookupType<momentumTransportModel>(this->owner().U().group());
65 
66  return model.k();
67  }
68  else
69  {
71  << "Turbulence model not found in mesh database" << nl
72  << "Database objects include: " << obr.sortedToc()
73  << abort(FatalError);
74 
75  return tmp<volScalarField>(nullptr);
76  }
77 }
78 
79 
80 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
81 
82 template<class CloudType>
84 (
85  CloudType& owner,
86  const fvMesh& mesh,
87  const dictionary& dict
88 )
89 :
90  ParticleForce<CloudType>(owner, mesh, dict, typeName, true),
91  rndGen_(owner.rndGen()),
92  lambda_(this->coeffs().template lookup<scalar>("lambda")),
93  turbulence_(readBool(this->coeffs().lookup("turbulence"))),
94  kPtr_(nullptr),
95  ownK_(false)
96 {}
97 
98 
99 template<class CloudType>
101 (
102  const BrownianMotionForce& bmf
103 )
104 :
105  ParticleForce<CloudType>(bmf),
106  rndGen_(bmf.rndGen_),
107  lambda_(bmf.lambda_),
108  turbulence_(bmf.turbulence_),
109  kPtr_(nullptr),
110  ownK_(false)
111 {}
112 
113 
114 // * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
115 
116 template<class CloudType>
118 {}
119 
120 
121 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
122 
123 template<class CloudType>
125 {
126  if (turbulence_)
127  {
128  if (store)
129  {
130  tmp<volScalarField> tk = kModel();
131  if (tk.isTmp())
132  {
133  kPtr_ = tk.ptr();
134  ownK_ = true;
135  }
136  else
137  {
138  kPtr_ = &tk();
139  ownK_ = false;
140  }
141  }
142  else
143  {
144  if (ownK_ && kPtr_)
145  {
146  deleteDemandDrivenData(kPtr_);
147  ownK_ = false;
148  }
149  }
150  }
151 }
152 
153 
154 template<class CloudType>
156 (
157  const typename CloudType::parcelType& p,
158  const typename CloudType::parcelType::trackingData& td,
159  const scalar dt,
160  const scalar mass,
161  const scalar Re,
162  const scalar muc
163 ) const
164 {
165  forceSuSp value(Zero, 0.0);
166 
167  const scalar dp = p.d();
168  const scalar Tc = td.Tc();
169 
170  const scalar alpha = 2.0*lambda_/dp;
171  const scalar cc = 1.0 + alpha*(1.257 + 0.4*exp(-1.1/alpha));
172 
173  // Boltzmann constant
174  const scalar kb = physicoChemical::k.value();
175 
176  scalar f = 0;
177  if (turbulence_)
178  {
179  const label celli = p.cell();
180  const volScalarField& k = *kPtr_;
181  const scalar kc = k[celli];
182  const scalar Dp = kb*Tc*cc/(3*mathematical::pi*muc*dp);
183  f = sqrt(2.0*sqr(kc)*sqr(Tc)/(Dp*dt));
184  }
185  else
186  {
187  const scalar s0 =
188  216*muc*kb*Tc/(sqr(mathematical::pi)*pow5(dp)*sqr(p.rho())*cc);
189  f = mass*sqrt(mathematical::pi*s0/dt);
190  }
191 
192 
193  // To generate a cubic distribution (3 independent directions) :
194  // const scalar sqrt2 = sqrt(2.0);
195  // for (direction dir = 0; dir < vector::nComponents; dir++)
196  // {
197  // const scalar x = rndGen_.sample01<scalar>();
198  // const scalar eta = sqrt2*erfInv(2*x - 1.0);
199  // value.Su()[dir] = f*eta;
200  // }
201 
202 
203  // To generate a spherical distribution:
204 
205  Random& rnd = this->owner().rndGen();
206 
207  const scalar theta = rnd.scalar01()*twoPi;
208  const scalar u = 2*rnd.scalar01() - 1;
209 
210  const scalar a = sqrt(1 - sqr(u));
211  const vector dir(a*cos(theta), a*sin(theta), u);
212 
213  value.Su() = f*mag(rnd.scalarNormal())*dir;
214 
215  return value;
216 }
217 
218 
219 // ************************************************************************* //
scalar y
label k
Calculates particle Brownian motion force.
BrownianMotionForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict)
Construct from mesh.
virtual void cacheFields(const bool store)
Cache fields.
virtual ~BrownianMotionForce()
Destructor.
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.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:79
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:221
Generic GeometricField class.
Abstract base class for particle forces.
Definition: ParticleForce.H:54
const fvMesh & mesh() const
Return the mesh database.
Random number generator.
Definition: Random.H:58
scalar scalar01()
Advance the state and return a scalar sample from a uniform.
Definition: RandomI.H:56
scalar scalarNormal()
Advance the state and return a scalar sample from a normal.
Definition: Random.C:31
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Helper container for force Su and Sp terms.
Definition: forceSuSp.H:64
const vector & Su() const
Return const access to the explicit contribution [kg m/s^2].
Definition: forceSuSpI.H:56
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
A class for managing temporary objects.
Definition: tmp.H:55
bool isTmp() const
Return true if this is really a temporary object.
Definition: tmpI.H:153
T * ptr() const
Return tmp pointer for reuse.
Definition: tmpI.H:205
Template functions to aid in the implementation of demand driven data.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
Fundamental dimensioned constants.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
compressibleMomentumTransportModel momentumTransportModel
const scalar twoPi(2 *pi)
const dimensionedScalar h
Planck constant.
Collection of constants.
static const zero Zero
Definition: zero.H:97
dimensionedScalar pow5(const dimensionedScalar &ds)
bool readBool(Istream &)
Definition: boolIO.C:60
dimensionedScalar exp(const dimensionedScalar &ds)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
dimensionedSymmTensor sqr(const dimensionedVector &dv)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
dimensionedScalar sin(const dimensionedScalar &ds)
void deleteDemandDrivenData(DataPtr &dataPtr)
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
error FatalError
static const char nl
Definition: Ostream.H:260
dimensionedScalar cos(const dimensionedScalar &ds)
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
labelList f(nPoints)
dictionary dict
volScalarField & p
Random rndGen(label(0))