VariableHardSphere.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 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 "VariableHardSphere.H"
27 #include "constants.H"
28 
29 using namespace Foam::constant::mathematical;
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class CloudType>
35 (
36  const dictionary& dict,
38 )
39 :
41  Tref_(this->coeffDict().template lookup<scalar>("Tref"))
42 {}
43 
44 
45 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
46 
47 template<class CloudType>
49 {}
50 
51 
52 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
53 
54 template<class CloudType>
56 (
57  const typename CloudType::parcelType& pP,
58  const typename CloudType::parcelType& pQ
59 ) const
60 {
61  const CloudType& cloud(this->owner());
62 
63  label typeIdP = pP.typeId();
64  label typeIdQ = pQ.typeId();
65 
66  scalar dPQ =
67  0.5
68  *(
69  cloud.constProps(typeIdP).d()
70  + cloud.constProps(typeIdQ).d()
71  );
72 
73  scalar omegaPQ =
74  0.5
75  *(
76  cloud.constProps(typeIdP).omega()
77  + cloud.constProps(typeIdQ).omega()
78  );
79 
80  scalar cR = mag(pP.U() - pQ.U());
81 
82  if (cR < vSmall)
83  {
84  return 0;
85  }
86 
87  scalar mP = cloud.constProps(typeIdP).mass();
88 
89  scalar mQ = cloud.constProps(typeIdQ).mass();
90 
91  scalar mR = mP*mQ/(mP + mQ);
92 
93  // calculating cross section = pi*dPQ^2, where dPQ is from Bird, eq. 4.79
94  scalar sigmaTPQ =
95  pi*dPQ*dPQ
96  *pow(2.0*physicoChemical::k.value()*Tref_/(mR*cR*cR), omegaPQ - 0.5)
97  /exp(Foam::lgamma(2.5 - omegaPQ));
98 
99  return sigmaTPQ*cR;
100 }
101 
102 
103 template<class CloudType>
105 (
106  typename CloudType::parcelType& pP,
107  typename CloudType::parcelType& pQ
108 )
109 {
110  CloudType& cloud(this->owner());
111 
112  label typeIdP = pP.typeId();
113  label typeIdQ = pQ.typeId();
114  vector& UP = pP.U();
115  vector& UQ = pQ.U();
116 
117  randomGenerator& rndGen(cloud.rndGen());
118 
119  scalar mP = cloud.constProps(typeIdP).mass();
120 
121  scalar mQ = cloud.constProps(typeIdQ).mass();
122 
123  vector Ucm = (mP*UP + mQ*UQ)/(mP + mQ);
124 
125  scalar cR = mag(UP - UQ);
126 
127  scalar cosTheta = 2.0*rndGen.scalar01() - 1.0;
128 
129  scalar sinTheta = sqrt(1.0 - cosTheta*cosTheta);
130 
131  scalar phi = twoPi*rndGen.scalar01();
132 
133  vector postCollisionRelU =
134  cR
135  *vector
136  (
137  cosTheta,
138  sinTheta*cos(phi),
139  sinTheta*sin(phi)
140  );
141 
142  UP = Ucm + postCollisionRelU*mQ/(mP + mQ);
143 
144  UQ = Ucm - postCollisionRelU*mP/(mP + mQ);
145 }
146 
147 
148 // ************************************************************************* //
label k
Templated DSMC particle collision class.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:80
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:225
VariableHardSphere(const dictionary &dict, CloudType &cloud)
Construct from dictionary.
virtual scalar sigmaTcR(const typename CloudType::parcelType &pP, const typename CloudType::parcelType &pQ) const
Return the collision cross section * relative velocity product.
virtual ~VariableHardSphere()
Destructor.
virtual void collide(typename CloudType::parcelType &pP, typename CloudType::parcelType &pQ)
Apply collision.
A cloud is a collection of lagrangian particles.
Definition: cloud.H:55
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Random number generator.
scalar scalar01()
Return a scalar uniformly distributed between zero and one.
mathematical constants.
const scalar twoPi(2 *pi)
dimensionedScalar exp(const dimensionedScalar &ds)
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 lgamma(const dimensionedScalar &ds)
dimensionedScalar sin(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
dimensionedScalar cos(const dimensionedScalar &ds)
dictionary dict
randomGenerator rndGen(653213)