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 | Copyright (C) 2011-2015 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 "demandDrivenData.H"
29 #include "turbulenceModel.H"
30 
31 using namespace Foam::constant;
32 
33 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
34 
35 template<class CloudType>
36 Foam::scalar Foam::BrownianMotionForce<CloudType>::erfInv(const scalar y) const
37 {
38  const scalar a = 0.147;
39  scalar k = 2.0/(mathematical::pi*a) + 0.5*log(1.0 - y*y);
40  scalar h = log(1.0 - y*y)/a;
41  scalar x = sqrt(-k + sqrt(k*k - h));
42 
43  if (y < 0.0)
44  {
45  return -x;
46  }
47  else
48  {
49  return x;
50  }
51 }
52 
53 
54 template<class CloudType>
57 {
58  const objectRegistry& obr = this->owner().mesh();
59  const word turbName =
61  (
63  this->owner().U().group()
64  );
65 
66  if (obr.foundObject<turbulenceModel>(turbName))
67  {
68  const turbulenceModel& model =
69  obr.lookupObject<turbulenceModel>(turbName);
70  return model.k();
71  }
72  else
73  {
75  (
76  "Foam::tmp<Foam::volScalarField>"
77  "Foam::DispersionRASModel<CloudType>::kModel() const"
78  )
79  << "Turbulence model not found in mesh database" << nl
80  << "Database objects include: " << obr.sortedToc()
81  << abort(FatalError);
82 
83  return tmp<volScalarField>(NULL);
84  }
85 }
86 
87 
88 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
89 
90 template<class CloudType>
92 (
93  CloudType& owner,
94  const fvMesh& mesh,
95  const dictionary& dict
96 )
97 :
98  ParticleForce<CloudType>(owner, mesh, dict, typeName, true),
99  rndGen_(owner.rndGen()),
100  lambda_(readScalar(this->coeffs().lookup("lambda"))),
101  turbulence_(readBool(this->coeffs().lookup("turbulence"))),
102  kPtr_(NULL),
103  ownK_(false)
104 {}
105 
106 
107 template<class CloudType>
109 (
110  const BrownianMotionForce& bmf
111 )
112 :
114  rndGen_(bmf.rndGen_),
115  lambda_(bmf.lambda_),
116  turbulence_(bmf.turbulence_),
117  kPtr_(NULL),
118  ownK_(false)
119 {}
120 
121 
122 // * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
123 
124 template<class CloudType>
126 {}
127 
128 
129 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
130 
131 template<class CloudType>
133 {
134  if (turbulence_)
135  {
136  if (store)
137  {
138  tmp<volScalarField> tk = kModel();
139  if (tk.isTmp())
140  {
141  kPtr_ = tk.ptr();
142  ownK_ = true;
143  }
144  else
145  {
146  kPtr_ = tk.operator->();
147  ownK_ = false;
148  }
149  }
150  else
151  {
152  if (ownK_ && kPtr_)
153  {
154  deleteDemandDrivenData(kPtr_);
155  ownK_ = false;
156  }
157  }
158  }
159 }
160 
161 
162 template<class CloudType>
164 (
165  const typename CloudType::parcelType& p,
166  const scalar dt,
167  const scalar mass,
168  const scalar Re,
169  const scalar muc
170 ) const
171 {
172  forceSuSp value(vector::zero, 0.0);
173 
174  const scalar dp = p.d();
175  const scalar Tc = p.Tc();
176 
177  const scalar eta = rndGen_.sample01<scalar>();
178  const scalar alpha = 2.0*lambda_/dp;
179  const scalar cc = 1.0 + alpha*(1.257 + 0.4*exp(-1.1/alpha));
180 
181  const scalar sigma = physicoChemical::sigma.value();
182 
183  scalar f = 0.0;
184  if (turbulence_)
185  {
186  const label cellI = p.cell();
187  const volScalarField& k = *kPtr_;
188  const scalar kc = k[cellI];
189  const scalar Dp = sigma*Tc*cc/(3*mathematical::pi*muc*dp);
190  f = eta/mass*sqrt(2.0*sqr(kc)*sqr(Tc)/(Dp*dt));
191  }
192  else
193  {
194  const scalar rhoRatio = p.rho()/p.rhoc();
195  const scalar s0 =
196  216*muc*sigma*Tc/(sqr(mathematical::pi)*pow5(dp)*(rhoRatio)*cc);
197  f = eta*sqrt(mathematical::pi*s0/dt);
198  }
199 
200  const scalar sqrt2 = sqrt(2.0);
201  for (label i = 0; i < 3; i++)
202  {
203  const scalar x = rndGen_.sample01<scalar>();
204  const scalar eta = sqrt2*erfInv(2*x - 1.0);
205  value.Su()[i] = mass*f*eta;
206  }
207 
208  return value;
209 }
210 
211 
212 // ************************************************************************* //
dimensionedScalar sqrt(const dimensionedScalar &ds)
bool readBool(Istream &)
Definition: boolIO.C:60
T * ptr() const
Return tmp pointer for reuse.
Definition: tmpI.H:148
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Collection of constants.
labelList f(nPoints)
BrownianMotionForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict)
Construct from mesh.
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
const char *const group
Group name for atomic constants.
bool isTmp() const
Return true if this is really a temporary object.
Definition: tmpI.H:127
void deleteDemandDrivenData(DataPtr &dataPtr)
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
dimensionedScalar exp(const dimensionedScalar &ds)
dynamicFvMesh & mesh
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:217
static word groupName(Name name, const word &group)
const vector & Su() const
Return const access to the explicit contribution [kg.m/s2].
Definition: forceSuSpI.H:56
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
dictionary dict
Calculates particle Brownian motion force.
dimensionedScalar log(const dimensionedScalar &ds)
static const char nl
Definition: Ostream.H:260
virtual forceSuSp calcCoupled(const typename CloudType::parcelType &p, const scalar dt, const scalar mass, const scalar Re, const scalar muc) const
Calculate the non-coupled force.
stressControl lookup("compactNormalStress") >> compactNormalStress
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
virtual ~BrownianMotionForce()
Destructor.
Template functions to aid in the implementation of demand driven data.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
virtual void cacheFields(const bool store)
Cache fields.
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
const fvMesh & mesh() const
Return the mesh database.
error FatalError
static const Vector zero
Definition: Vector.H:80
Abstract base class for particle forces.
Definition: ParticleForce.H:53
dimensionedScalar pow5(const dimensionedScalar &ds)
static const word propertiesName
Default name of the turbulence properties dictionary.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
dimensionedSymmTensor sqr(const dimensionedVector &dv)
A class for managing temporary objects.
Definition: PtrList.H:118
U
Definition: pEqn.H:82
const Type & value() const
Return const reference to value.
Helper container for force Su and Sp terms.
Definition: forceSuSp.H:58