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-2024 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 #include "standardNormal.H"
32 
33 using namespace Foam::constant;
34 
35 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
36 
37 template<class CloudType>
40 {
41  const objectRegistry& obr = this->owner().mesh();
42 
43  if (obr.foundType<momentumTransportModel>(this->owner().U().group()))
44  {
45  const momentumTransportModel& model =
46  obr.lookupType<momentumTransportModel>(this->owner().U().group());
47 
48  return model.k();
49  }
50  else
51  {
53  << "Turbulence model not found in mesh database" << nl
54  << "Database objects include: " << obr.sortedToc()
55  << abort(FatalError);
56 
57  return tmp<volScalarField>(nullptr);
58  }
59 }
60 
61 
62 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
63 
64 template<class CloudType>
66 (
67  CloudType& owner,
68  const fvMesh& mesh,
69  const dictionary& dict
70 )
71 :
72  ParticleForce<CloudType>(owner, mesh, dict, typeName, true),
73  lambda_(this->coeffs().template lookup<scalar>("lambda")),
74  turbulence_(readBool(this->coeffs().lookup("turbulence"))),
75  kPtr_(nullptr),
76  ownK_(false)
77 {}
78 
79 
80 template<class CloudType>
82 (
83  const BrownianMotionForce& bmf
84 )
85 :
87  lambda_(bmf.lambda_),
88  turbulence_(bmf.turbulence_),
89  kPtr_(nullptr),
90  ownK_(false)
91 {}
92 
93 
94 // * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
95 
96 template<class CloudType>
98 {}
99 
100 
101 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
102 
103 template<class CloudType>
105 {
106  if (turbulence_)
107  {
108  if (store)
109  {
110  tmp<volScalarField> tk = kModel();
111  if (tk.isTmp())
112  {
113  kPtr_ = tk.ptr();
114  ownK_ = true;
115  }
116  else
117  {
118  kPtr_ = &tk();
119  ownK_ = false;
120  }
121  }
122  else
123  {
124  if (ownK_ && kPtr_)
125  {
126  deleteDemandDrivenData(kPtr_);
127  ownK_ = false;
128  }
129  }
130  }
131 }
132 
133 
134 template<class CloudType>
136 (
137  const typename CloudType::parcelType& p,
138  const typename CloudType::parcelType::trackingData& td,
139  const scalar dt,
140  const scalar mass,
141  const scalar Re,
142  const scalar muc
143 ) const
144 {
145  forceSuSp value(Zero, 0.0);
146 
147  const scalar dp = p.d();
148  const scalar Tc = td.Tc();
149 
150  const scalar alpha = 2.0*lambda_/dp;
151  const scalar cc = 1.0 + alpha*(1.257 + 0.4*exp(-1.1/alpha));
152 
153  // Boltzmann constant
154  const scalar kb = physicoChemical::k.value();
155 
156  scalar f = 0;
157  if (turbulence_)
158  {
159  const label celli = p.cell();
160  const volScalarField& k = *kPtr_;
161  const scalar kc = k[celli];
162  const scalar Dp = kb*Tc*cc/(3*mathematical::pi*muc*dp);
163  f = sqrt(2.0*sqr(kc)*sqr(Tc)/(Dp*dt));
164  }
165  else
166  {
167  const scalar s0 =
168  216*muc*kb*Tc/(sqr(mathematical::pi)*pow5(dp)*sqr(p.rho())*cc);
169  f = mass*sqrt(mathematical::pi*s0/dt);
170  }
171 
172  randomGenerator& rndGen = this->owner().rndGen();
173  distributions::standardNormal& stdNormal = this->owner().stdNormal();
174 
175  // To generate a cubic distribution (i.e., 3 independent directions):
176  // value.Su() = f*stdNormal.sample<vector>();
177 
178  // To generate a spherical distribution:
179  const scalar theta = rndGen.scalar01()*twoPi;
180  const scalar u = 2*rndGen.scalar01() - 1;
181  const scalar a = sqrt(1 - sqr(u));
182  const vector dir(a*cos(theta), a*sin(theta), u);
183  value.Su() = f*mag(stdNormal.sample())*dir;
184 
185  return value;
186 }
187 
188 
189 // ************************************************************************* //
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:80
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:225
Generic GeometricField class.
Abstract base class for particle forces.
Definition: ParticleForce.H:54
const fvMesh & mesh() const
Return the mesh database.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Standard normal distribution. Not selectable.
virtual scalar sample() const
Sample the distribution.
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:99
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Random number generator.
scalar scalar01()
Return a scalar uniformly distributed between zero and one.
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:334
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)
Collection of constants.
static const zero Zero
Definition: zero.H:97
dimensionedScalar pow5(const dimensionedScalar &ds)
bool readBool(Istream &)
Definition: boolIO.C:66
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)
void deleteDemandDrivenData(DataType *&dataPtr)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
dimensionedScalar sin(const dimensionedScalar &ds)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
error FatalError
static const char nl
Definition: Ostream.H:266
dimensionedScalar cos(const dimensionedScalar &ds)
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
labelList f(nPoints)
dictionary dict
randomGenerator rndGen(653213)
volScalarField & p