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-------------------------------------------------------------------------------
11License
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"
32#include "demandDrivenData.H"
33#include "turbulenceModel.H"
34
35using namespace Foam::constant;
36
37// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
38
39template<class CloudType>
42{
43 const objectRegistry& obr = this->owner().mesh();
44 const word turbName =
45 IOobject::groupName
46 (
47 turbulenceModel::propertiesName,
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
69template<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 turbulence_(this->coeffs().getBool("turbulence")),
81 kPtr_(nullptr),
82 ownK_(false)
83{}
84
85
86template<class CloudType>
88(
89 const BrownianMotionForce& bmf
90)
91:
93 rndGen_(bmf.rndGen_),
94 lambda_(bmf.lambda_),
95 turbulence_(bmf.turbulence_),
96 kPtr_(nullptr),
97 ownK_(false)
98{}
99
100
101// * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
102
103template<class CloudType>
105{}
106
107
108// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
109
110template<class CloudType>
112{
113 if (turbulence_)
114 {
115 if (store)
116 {
117 tmp<volScalarField> tk = kModel();
118 if (tk.isTmp())
119 {
120 kPtr_ = tk.ptr();
121 ownK_ = true;
122 }
123 else
124 {
125 kPtr_ = &tk();
126 ownK_ = false;
127 }
128 }
129 else
130 {
131 if (ownK_ && kPtr_)
132 {
134 ownK_ = false;
135 }
136 }
137 }
138}
139
140
141template<class CloudType>
143(
144 const typename CloudType::parcelType& p,
145 const typename CloudType::parcelType::trackingData& td,
146 const scalar dt,
147 const scalar mass,
148 const scalar Re,
149 const scalar muc
150) const
151{
152 forceSuSp value(Zero);
153
154 const scalar dp = p.d();
155 const scalar Tc = td.Tc();
156
157 const scalar alpha = 2.0*lambda_/dp;
158 const scalar cc = 1.0 + alpha*(1.257 + 0.4*exp(-1.1/alpha));
159
160 // Boltzmann constant
161 const scalar kb = physicoChemical::k.value();
162
163 scalar f = 0;
164 if (turbulence_)
165 {
166 const label celli = p.cell();
167 const volScalarField& k = *kPtr_;
168 const scalar kc = k[celli];
169 const scalar Dp = kb*Tc*cc/(3*mathematical::pi*muc*dp);
170 f = sqrt(2.0*sqr(kc)*sqr(Tc)/(Dp*dt));
171 }
172 else
173 {
174 const scalar s0 =
175 216*muc*kb*Tc/(sqr(mathematical::pi)*pow5(dp)*sqr(p.rho())*cc);
176 f = mass*sqrt(mathematical::pi*s0/dt);
177 }
178
179
180 // To generate a cubic distribution (3 independent directions) :
181 // const scalar sqrt2 = sqrt(2.0);
182 // for (direction dir = 0; dir < vector::nComponents; dir++)
183 // {
184 // const scalar x = rndGen_.sample01<scalar>();
185 // const scalar eta = sqrt2*Math::erfInv(2*x - 1.0);
186 // value.Su()[dir] = f*eta;
187 // }
188
189
190 // To generate a spherical distribution:
191
192 Random& rnd = this->owner().rndGen();
193
194 const scalar theta = rnd.sample01<scalar>()*twoPi;
195 const scalar u = 2*rnd.sample01<scalar>() - 1;
196
197 const scalar a = sqrt(1 - sqr(u));
198 const vector dir(a*cos(theta), a*sin(theta), u);
199
200 value.Su() = f*mag(rnd.GaussNormal<scalar>())*dir;
201
202 return value;
203}
204
205
206// ************************************************************************* //
label k
compressible::turbulenceModel & turb
Calculates particle Brownian motion force.
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:75
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:137
Abstract base class for particle forces.
Definition: ParticleForce.H:60
Random number generator.
Definition: Random.H:60
Type GaussNormal()
Type sample01()
Return a sample whose components lie in the range [0,1].
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const Type & value() const
Return const reference to value.
Class used to pass data into container.
Helper container for force Su and Sp terms.
Definition: forceSuSp.H:67
const vector & Su() const
Return const access to the explicit contribution [kg.m/s2].
Definition: forceSuSpI.H:61
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Registry of regIOobjects.
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
A class for managing temporary objects.
Definition: tmp.H:65
bool isTmp() const noexcept
Identical to is_pointer()
Definition: tmp.H:295
T * ptr() const
Return managed pointer for reuse, or clone() the object reference.
Definition: tmpI.H:255
Abstract base class for turbulence models (RAS, LES and laminar).
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
volScalarField & p
dynamicFvMesh & mesh
Template functions to aid in the implementation of demand driven data.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Fundamental dimensioned constants.
constexpr const char *const group
Group name for atomic constants.
constexpr scalar pi(M_PI)
const dimensionedScalar k
Boltzmann constant.
Different types of constants.
dimensionedScalar pow5(const dimensionedScalar &ds)
dimensionedScalar exp(const dimensionedScalar &ds)
scalarField Re(const UList< complex > &cf)
Extract real component.
Definition: complexField.C:159
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sin(const dimensionedScalar &ds)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
void deleteDemandDrivenData(DataPtr &dataPtr)
dimensionedScalar cos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList f(nPoints)
volScalarField & alpha
dictionary dict
scalar kb
Random rndGen
Definition: createFields.H:23