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  << "Turbulence model not found in mesh database" << nl
76  << "Database objects include: " << obr.sortedToc()
77  << abort(FatalError);
78 
79  return tmp<volScalarField>(NULL);
80  }
81 }
82 
83 
84 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
85 
86 template<class CloudType>
88 (
89  CloudType& owner,
90  const fvMesh& mesh,
91  const dictionary& dict
92 )
93 :
94  ParticleForce<CloudType>(owner, mesh, dict, typeName, true),
95  rndGen_(owner.rndGen()),
96  lambda_(readScalar(this->coeffs().lookup("lambda"))),
97  turbulence_(readBool(this->coeffs().lookup("turbulence"))),
98  kPtr_(NULL),
99  ownK_(false)
100 {}
101 
102 
103 template<class CloudType>
105 (
106  const BrownianMotionForce& bmf
107 )
108 :
110  rndGen_(bmf.rndGen_),
111  lambda_(bmf.lambda_),
112  turbulence_(bmf.turbulence_),
113  kPtr_(NULL),
114  ownK_(false)
115 {}
116 
117 
118 // * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
119 
120 template<class CloudType>
122 {}
123 
124 
125 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
126 
127 template<class CloudType>
129 {
130  if (turbulence_)
131  {
132  if (store)
133  {
134  tmp<volScalarField> tk = kModel();
135  if (tk.isTmp())
136  {
137  kPtr_ = tk.ptr();
138  ownK_ = true;
139  }
140  else
141  {
142  kPtr_ = tk.operator->();
143  ownK_ = false;
144  }
145  }
146  else
147  {
148  if (ownK_ && kPtr_)
149  {
150  deleteDemandDrivenData(kPtr_);
151  ownK_ = false;
152  }
153  }
154  }
155 }
156 
157 
158 template<class CloudType>
160 (
161  const typename CloudType::parcelType& p,
162  const scalar dt,
163  const scalar mass,
164  const scalar Re,
165  const scalar muc
166 ) const
167 {
168  forceSuSp value(vector::zero, 0.0);
169 
170  const scalar dp = p.d();
171  const scalar Tc = p.Tc();
172 
173  const scalar eta = rndGen_.sample01<scalar>();
174  const scalar alpha = 2.0*lambda_/dp;
175  const scalar cc = 1.0 + alpha*(1.257 + 0.4*exp(-1.1/alpha));
176 
177  const scalar sigma = physicoChemical::sigma.value();
178 
179  scalar f = 0.0;
180  if (turbulence_)
181  {
182  const label cellI = p.cell();
183  const volScalarField& k = *kPtr_;
184  const scalar kc = k[cellI];
185  const scalar Dp = sigma*Tc*cc/(3*mathematical::pi*muc*dp);
186  f = eta/mass*sqrt(2.0*sqr(kc)*sqr(Tc)/(Dp*dt));
187  }
188  else
189  {
190  const scalar rhoRatio = p.rho()/p.rhoc();
191  const scalar s0 =
192  216*muc*sigma*Tc/(sqr(mathematical::pi)*pow5(dp)*(rhoRatio)*cc);
193  f = eta*sqrt(mathematical::pi*s0/dt);
194  }
195 
196  const scalar sqrt2 = sqrt(2.0);
197  for (label i = 0; i < 3; i++)
198  {
199  const scalar x = rndGen_.sample01<scalar>();
200  const scalar eta = sqrt2*erfInv(2*x - 1.0);
201  value.Su()[i] = mass*f*eta;
202  }
203 
204  return value;
205 }
206 
207 
208 // ************************************************************************* //
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
mathematicalConstants.H
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::BrownianMotionForce::BrownianMotionForce
BrownianMotionForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict)
Construct from mesh.
Definition: BrownianMotionForce.C:88
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::constant
Collection of constants.
Definition: atomicConstants.C:37
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:216
Foam::BrownianMotionForce::~BrownianMotionForce
virtual ~BrownianMotionForce()
Destructor.
Definition: BrownianMotionForce.C:121
Foam::Re
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:97
Foam::BrownianMotionForce::cacheFields
virtual void cacheFields(const bool store)
Cache fields.
Definition: BrownianMotionForce.C:128
Foam::constant::atomic::group
const char *const group
Group name for atomic constants.
Definition: atomicConstants.C:40
Foam::BrownianMotionForce::erfInv
scalar erfInv(const scalar y) const
Inverse erf for Gaussian distribution.
Definition: BrownianMotionForce.C:36
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:252
Foam::DSMCCloud::rndGen
Random & rndGen()
Return refernce to the random object.
Definition: DSMCCloudI.H:121
U
U
Definition: pEqn.H:46
Foam::tmp::isTmp
bool isTmp() const
Return true if this is really a temporary object.
Definition: tmpI.H:125
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:50
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:40
Foam::label
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
Foam::constant::universal::h
const dimensionedScalar h
Planck constant.
Definition: createFields.H:6
Foam::DSMCCloud::mesh
const fvMesh & mesh() const
Return refernce to the mesh.
Definition: DSMCCloudI.H:41
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::forceSuSp
Helper container for force Su and Sp terms.
Definition: forceSuSp.H:58
Foam::ParticleForce
Abstract base class for particle forces.
Definition: ParticleForce.H:53
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:60
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::pow5
dimensionedScalar pow5(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:109
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:253
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::turbulenceModel::k
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::tmp::ptr
T * ptr() const
Return tmp pointer for reuse.
Definition: tmpI.H:146
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::HashTable::sortedToc
List< Key > sortedToc() const
Return the table of contents as a sorted list.
Definition: HashTable.C:216
Foam::BrownianMotionForce
Calculates particle Brownian motion force.
Definition: BrownianMotionForce.H:51
Foam::BrownianMotionForce::calcCoupled
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.
Definition: BrownianMotionForce.C:160
Foam::objectRegistry::foundObject
bool foundObject(const word &name) const
Is the named Type found?
Definition: objectRegistryTemplates.C:142
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::BrownianMotionForce::lambda_
const scalar lambda_
Molecular free path length [m].
Definition: BrownianMotionForce.H:61
Foam::readScalar
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
f
labelList f(nPoints)
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::constant::physicoChemical::sigma
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Foam::constant::mathematical::pi
const scalar pi(M_PI)
Foam::DSMCCloud::parcelType
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:217
Foam::readBool
bool readBool(Istream &)
Definition: boolIO.C:60
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Definition: objectRegistryTemplates.C:165
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::forceSuSp::Su
const vector & Su() const
Return const access to the explicit contribution [kg.m/s2].
Definition: forceSuSpI.H:56
Foam::BrownianMotionForce::turbulence_
bool turbulence_
Turbulence flag.
Definition: BrownianMotionForce.H:64
turbulenceModel.H
Foam::BrownianMotionForce::rndGen_
cachedRandom & rndGen_
Reference to the cloud random number generator.
Definition: BrownianMotionForce.H:58
Foam::BrownianMotionForce::kModel
tmp< volScalarField > kModel() const
Return the k field from the turbulence model.
Definition: BrownianMotionForce.C:56
BrownianMotionForce.H
y
scalar y
Definition: LISASMDCalcMethod1.H:14
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress