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  kPtr_(nullptr),
81  turbulence_(this->coeffs().getBool("turbulence")),
82  ownK_(false),
83  useSpherical_(this->coeffs().getOrDefault("spherical", true))
84 {}
85 
86 
87 template<class CloudType>
89 (
90  const BrownianMotionForce& bmf
91 )
92 :
94  rndGen_(bmf.rndGen_),
95  lambda_(bmf.lambda_),
96  kPtr_(nullptr),
97  turbulence_(bmf.turbulence_),
98  ownK_(false),
99  useSpherical_(bmf.useSpherical_)
100 {}
101 
102 
103 // * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
104 
105 template<class CloudType>
107 {
108  if (ownK_)
109  {
110  deleteDemandDrivenData(kPtr_);
111  ownK_ = false;
112  }
113 }
114 
115 
116 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
117 
118 template<class CloudType>
120 {
121  if (turbulence_)
122  {
123  if (store)
124  {
125  tmp<volScalarField> tk = kModel();
126  if (tk.movable())
127  {
128  // Take ownership
129  kPtr_ = tk.ptr();
130  ownK_ = true;
131  }
132  else
133  {
134  kPtr_ = &tk.cref();
135  ownK_ = false;
136  }
137  }
138  else
139  {
140  if (ownK_)
141  {
142  deleteDemandDrivenData(kPtr_);
143  ownK_ = false;
144  }
145  }
146  }
147 }
148 
149 
150 template<class CloudType>
152 (
153  const typename CloudType::parcelType& p,
154  const typename CloudType::parcelType::trackingData& td,
155  const scalar dt,
156  const scalar mass,
157  const scalar Re,
158  const scalar muc
159 ) const
160 {
161  forceSuSp value(Zero);
162 
163  const scalar dp = p.d();
164  const scalar Tc = td.Tc();
165 
166  const scalar alpha = 2.0*lambda_/dp;
167  const scalar cc = 1.0 + alpha*(1.257 + 0.4*exp(-1.1/alpha));
168 
169  // Boltzmann constant
170  const scalar kb = physicoChemical::k.value();
171 
172  scalar f = 0;
173  if (turbulence_)
174  {
175  const label celli = p.cell();
176  const volScalarField& k = *kPtr_;
177  const scalar kc = k[celli];
178  const scalar Dp = kb*Tc*cc/(3*mathematical::pi*muc*dp);
179  f = sqrt(2.0*sqr(kc)*sqr(Tc)/(Dp*dt));
180  }
181  else
182  {
183  const scalar s0 =
184  216*muc*kb*Tc/(sqr(mathematical::pi)*pow5(dp)*sqr(p.rho())*cc);
185  f = mass*sqrt(mathematical::pi*s0/dt);
186  }
187 
188  Random& rnd = this->owner().rndGen();
189 
190  if (useSpherical_)
191  {
192  // To generate a spherical distribution:
193  const scalar theta = rnd.sample01<scalar>()*twoPi;
194  const scalar u = 2*rnd.sample01<scalar>() - 1;
195 
196  const scalar a = sqrt(1 - sqr(u));
197  const vector dir(a*cos(theta), a*sin(theta), u);
198 
199  value.Su() = f*mag(rnd.GaussNormal<scalar>())*dir;
200  }
201  else
202  {
203  // Generate a cubic distribution (3 independent directions)
204  value.Su() = f*rnd.GaussNormal<vector>();
205 
206  // OLD CODE for cubic distribution
207  // const scalar sqrt2 = sqrt(2.0);
208  // for (direction dir = 0; dir < vector::nComponents; dir++)
209  // {
210  // const scalar x = rnd.sample01<scalar>();
211  // const scalar eta = sqrt2*Math::erfInv(2*x - 1.0);
212  // value.Su()[dir] = f*eta;
213  // }
214  }
215 
216  return value;
217 }
218 
219 
220 // ************************************************************************* //
Different types of constants.
const Type & value() const noexcept
Return const reference to value.
dictionary dict
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Random rndGen
Definition: createFields.H:23
compressible::turbulenceModel & turb
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
dimensionedScalar sqrt(const dimensionedScalar &ds)
Abstract base class for particle forces.
Definition: ParticleForce.H:54
label k
Boltzmann constant.
dimensionedScalar pow5(const dimensionedScalar &ds)
Fundamental dimensioned constants.
constexpr const char *const group
Group name for atomic constants.
Helper container for force Su and Sp terms.
Definition: forceSuSp.H:60
const fvMesh & mesh() const noexcept
Return the mesh database.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
dynamicFvMesh & mesh
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar exp(const dimensionedScalar &ds)
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.
static const word propertiesName
Default name of the turbulence properties dictionary.
constexpr scalar twoPi(2 *M_PI)
constexpr scalar pi(M_PI)
Vector< scalar > vector
Definition: vector.H:57
errorManip< error > abort(error &err)
Definition: errorManip.H:139
BrownianMotionForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict)
Construct from mesh.
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:290
dimensionedScalar sin(const dimensionedScalar &ds)
labelList f(nPoints)
virtual void cacheFields(const bool store)
Cache fields.
scalarField Re(const UList< complex > &cmplx)
Extract real component.
Definition: complexField.C:207
const dimensionedScalar k
Boltzmann constant.
Template functions to aid in the implementation of demand driven data.
U
Definition: pEqn.H:72
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual ~BrownianMotionForce()
Destructor.
scalar kb
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
void deleteDemandDrivenData(DataPtr &dataPtr)
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
Calculates particle Brownian motion force.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127