LarsenBorgnakkeVariableHardSphere.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-2015 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
29#include "constants.H"
30
31using namespace Foam::constant::mathematical;
32
33// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34
35template<class CloudType>
37(
38 scalar ChiA,
39 scalar ChiB
40)
41{
42 CloudType& cloud(this->owner());
43
44 Random& rndGen = cloud.rndGen();
45
46 scalar ChiAMinusOne = ChiA - 1;
47 scalar ChiBMinusOne = ChiB - 1;
48
49 if (ChiAMinusOne < SMALL && ChiBMinusOne < SMALL)
50 {
51 return rndGen.sample01<scalar>();
52 }
53
54 scalar energyRatio;
55 scalar P;
56
57 do
58 {
59 P = 0;
60
61 energyRatio = rndGen.sample01<scalar>();
62
63 if (ChiAMinusOne < SMALL)
64 {
65 P = 1.0 - pow(energyRatio, ChiB);
66 }
67 else if (ChiBMinusOne < SMALL)
68 {
69 P = 1.0 - pow(energyRatio, ChiA);
70 }
71 else
72 {
73 P =
74 pow
75 (
76 (ChiAMinusOne + ChiBMinusOne)*energyRatio/ChiAMinusOne,
77 ChiAMinusOne
78 )
79 *pow
80 (
81 (ChiAMinusOne + ChiBMinusOne)*(1 - energyRatio)
82 /ChiBMinusOne,
83 ChiBMinusOne
84 );
85 }
86 } while (P < rndGen.sample01<scalar>());
87
88 return energyRatio;
89}
90
91
92// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
93
94template<class CloudType>
97(
98 const dictionary& dict,
100)
101:
103 Tref_(this->coeffDict().getScalar("Tref")),
104 relaxationCollisionNumber_
105 (
106 this->coeffDict().getScalar("relaxationCollisionNumber")
107 )
108{}
109
110
111// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
112
113template<class CloudType>
116{}
117
118
119// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
120
121template<class CloudType>
123{
124 return true;
125}
126
127
128template<class CloudType>
130(
131 const typename CloudType::parcelType& pP,
132 const typename CloudType::parcelType& pQ
133) const
134{
135 const CloudType& cloud(this->owner());
136
137 label typeIdP = pP.typeId();
138 label typeIdQ = pQ.typeId();
139
140 scalar dPQ =
141 0.5
142 *(
143 cloud.constProps(typeIdP).d()
144 + cloud.constProps(typeIdQ).d()
145 );
146
147 scalar omegaPQ =
148 0.5
149 *(
150 cloud.constProps(typeIdP).omega()
151 + cloud.constProps(typeIdQ).omega()
152 );
153
154 scalar cR = mag(pP.U() - pQ.U());
155
156 if (cR < VSMALL)
157 {
158 return 0;
159 }
160
161 scalar mP = cloud.constProps(typeIdP).mass();
162 scalar mQ = cloud.constProps(typeIdQ).mass();
163 scalar mR = mP*mQ/(mP + mQ);
164
165 // calculating cross section = pi*dPQ^2, where dPQ is from Bird, eq. 4.79
166 scalar sigmaTPQ =
167 pi*dPQ*dPQ
168 *pow(2.0*physicoChemical::k.value()*Tref_/(mR*cR*cR), omegaPQ - 0.5)
169 /exp(Foam::lgamma(2.5 - omegaPQ));
170
171 return sigmaTPQ*cR;
172}
173
174
175template<class CloudType>
177(
178 typename CloudType::parcelType& pP,
179 typename CloudType::parcelType& pQ
180)
181{
182 CloudType& cloud(this->owner());
183
184 label typeIdP = pP.typeId();
185 label typeIdQ = pQ.typeId();
186 vector& UP = pP.U();
187 vector& UQ = pQ.U();
188 scalar& EiP = pP.Ei();
189 scalar& EiQ = pQ.Ei();
190
191 Random& rndGen = cloud.rndGen();
192
193 scalar inverseCollisionNumber = 1/relaxationCollisionNumber_;
194
195 // Larsen Borgnakke internal energy redistribution part. Using the serial
196 // application of the LB method, as per the INELRS subroutine in Bird's
197 // DSMC0R.FOR
198
199 scalar preCollisionEiP = EiP;
200 scalar preCollisionEiQ = EiQ;
201
202 direction iDofP = cloud.constProps(typeIdP).internalDegreesOfFreedom();
203 direction iDofQ = cloud.constProps(typeIdQ).internalDegreesOfFreedom();
204
205 scalar omegaPQ =
206 0.5
207 *(
208 cloud.constProps(typeIdP).omega()
209 + cloud.constProps(typeIdQ).omega()
210 );
211
212 scalar mP = cloud.constProps(typeIdP).mass();
213 scalar mQ = cloud.constProps(typeIdQ).mass();
214 scalar mR = mP*mQ/(mP + mQ);
215 vector Ucm = (mP*UP + mQ*UQ)/(mP + mQ);
216 scalar cRsqr = magSqr(UP - UQ);
217 scalar availableEnergy = 0.5*mR*cRsqr;
218 scalar ChiB = 2.5 - omegaPQ;
219
220 if (iDofP > 0)
221 {
222 if (inverseCollisionNumber > rndGen.sample01<scalar>())
223 {
224 availableEnergy += preCollisionEiP;
225
226 if (iDofP == 2)
227 {
228 scalar energyRatio =
229 1.0 - pow(rndGen.sample01<scalar>(), (1.0/ChiB));
230 EiP = energyRatio*availableEnergy;
231 }
232 else
233 {
234 scalar ChiA = 0.5*iDofP;
235 EiP = energyRatio(ChiA, ChiB)*availableEnergy;
236 }
237
238 availableEnergy -= EiP;
239 }
240 }
241
242 if (iDofQ > 0)
243 {
244 if (inverseCollisionNumber > rndGen.sample01<scalar>())
245 {
246 availableEnergy += preCollisionEiQ;
247
248 if (iDofQ == 2)
249 {
250 const scalar energyRatio =
251 1.0 - pow(rndGen.sample01<scalar>(), (1.0/ChiB));
252
253 EiQ = energyRatio*availableEnergy;
254 }
255 else
256 {
257 const scalar ChiA = 0.5*iDofQ;
258 EiQ = energyRatio(ChiA, ChiB)*availableEnergy;
259 }
260
261 availableEnergy -= EiQ;
262 }
263 }
264
265 // Rescale the translational energy
266 scalar cR = sqrt(2.0*availableEnergy/mR);
267
268 // Variable Hard Sphere collision part
269 scalar cosTheta = 2.0*rndGen.sample01<scalar>() - 1.0;
270 scalar sinTheta = sqrt(1.0 - cosTheta*cosTheta);
271 scalar phi = twoPi*rndGen.sample01<scalar>();
272
273 vector postCollisionRelU =
274 cR
275 *vector
276 (
277 cosTheta,
278 sinTheta*cos(phi),
279 sinTheta*sin(phi)
280 );
281
282 UP = Ucm + postCollisionRelU*mQ/(mP + mQ);
283 UQ = Ucm - postCollisionRelU*mP/(mP + mQ);
284}
285
286
287// ************************************************************************* //
surfaceScalarField & phi
Templated DSMC particle collision class.
virtual void collide()=0
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
Variable Hard Sphere BinaryCollision Model with Larsen Borgnakke internal energy redistribution....
virtual scalar sigmaTcR(const typename CloudType::parcelType &pP, const typename CloudType::parcelType &pQ) const
Return the collision cross section * relative velocity product.
virtual bool active() const
Flag to indicate whether model activates collision model.
Random number generator.
Definition: Random.H:60
Type sample01()
Return a sample whose components lie in the range [0,1].
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:60
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mathematical constants.
constexpr scalar pi(M_PI)
constexpr scalar twoPi(2 *M_PI)
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar lgamma(const dimensionedScalar &ds)
dimensionedScalar sin(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
uint8_t direction
Definition: direction.H:56
Vector< scalar > vector
Definition: vector.H:61
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensionedScalar cos(const dimensionedScalar &ds)
dictionary dict
Random rndGen
Definition: createFields.H:23