Loading...
Searching...
No Matches
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>
41Foam::BrownianMotionForce<CloudType>::kModel() const
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
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 kPtr_(nullptr),
81 turbulence_(this->coeffs().getBool("turbulence")),
82 ownK_(false),
83 useSpherical_(this->coeffs().getOrDefault("spherical", true))
84{}
85
86
87template<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
105template<class CloudType>
107{
108 if (ownK_)
109 {
111 ownK_ = false;
113}
114
115
116// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
117
118template<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 {
143 ownK_ = false;
145 }
146 }
147}
148
149
150template<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// ************************************************************************* //
label k
compressible::turbulenceModel & turb
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.
Foam::List< Key > sortedToc(const Compare &comp) const
Definition HashTable.C:168
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Abstract base class for particle forces.
const fvMesh & mesh() const noexcept
Return the mesh database.
const CloudType & owner() const noexcept
Return const access to the cloud owner.
ParticleForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict, const word &forceType, const bool readCoeffs)
Construct from mesh.
const dictionary & coeffs() const noexcept
Return the force coefficients dictionary.
Random number generator.
Definition Random.H:56
Type GaussNormal()
Return a sample whose components are normally distributed with zero mean and unity variance N(0,...
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:133
Helper container for force Su and Sp terms.
Definition forceSuSp.H:63
const vector & Su() const
Return const access to the explicit contribution [kg.m/s2].
Definition forceSuSpI.H:54
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
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:75
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
Definition tmpI.H:221
bool movable() const noexcept
True if this is a non-null managed pointer with a unique ref-count.
Definition tmpI.H:214
T * ptr() const
Return managed pointer for reuse, or clone() the object reference.
Definition tmpI.H:256
Abstract base class for turbulence models (RAS, LES and laminar).
static const word propertiesName
Default name of the turbulence properties dictionary.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
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:600
Fundamental dimensioned constants.
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
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)
DSMCCloud< dsmcParcel > CloudType
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar sin(const dimensionedScalar &ds)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition errorManip.H:139
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
scalarField Re(const UList< complex > &cmplx)
Extract real component.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
Vector< scalar > vector
Definition vector.H:57
void deleteDemandDrivenData(DataPtr &dataPtr)
dimensionedScalar cos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
labelList f(nPoints)
volScalarField & alpha
dictionary dict
scalar kb
Random rndGen