VariableHardSphere.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 -------------------------------------------------------------------------------
10 License
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 
28 #include "VariableHardSphere.H"
29 #include "constants.H"
30 
31 using namespace Foam::constant::mathematical;
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
35 template<class CloudType>
37 (
38  const dictionary& dict,
40 )
41 :
43  Tref_(this->coeffDict().getScalar("Tref"))
44 {}
45 
46 
47 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
48 
49 template<class CloudType>
51 {}
52 
53 
54 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
55 
56 template<class CloudType>
58 {
59  return true;
60 }
61 
62 
63 template<class CloudType>
65 (
66  const typename CloudType::parcelType& pP,
67  const typename CloudType::parcelType& pQ
68 ) const
69 {
70  const CloudType& cloud(this->owner());
71 
72  label typeIdP = pP.typeId();
73  label typeIdQ = pQ.typeId();
74 
75  scalar dPQ =
76  0.5
77  *(
78  cloud.constProps(typeIdP).d()
79  + cloud.constProps(typeIdQ).d()
80  );
81 
82  scalar omegaPQ =
83  0.5
84  *(
85  cloud.constProps(typeIdP).omega()
86  + cloud.constProps(typeIdQ).omega()
87  );
88 
89  scalar cR = mag(pP.U() - pQ.U());
90 
91  if (cR < VSMALL)
92  {
93  return 0;
94  }
95 
96  scalar mP = cloud.constProps(typeIdP).mass();
97 
98  scalar mQ = cloud.constProps(typeIdQ).mass();
99 
100  scalar mR = mP*mQ/(mP + mQ);
101 
102  // calculating cross section = pi*dPQ^2, where dPQ is from Bird, eq. 4.79
103  scalar sigmaTPQ =
104  pi*dPQ*dPQ
105  *pow(2.0*physicoChemical::k.value()*Tref_/(mR*cR*cR), omegaPQ - 0.5)
106  /exp(Foam::lgamma(2.5 - omegaPQ));
107 
108  return sigmaTPQ*cR;
109 }
110 
111 
112 template<class CloudType>
114 (
115  typename CloudType::parcelType& pP,
116  typename CloudType::parcelType& pQ
117 )
118 {
119  CloudType& cloud(this->owner());
120 
121  label typeIdP = pP.typeId();
122  label typeIdQ = pQ.typeId();
123  vector& UP = pP.U();
124  vector& UQ = pQ.U();
125 
126  Random& rndGen = cloud.rndGen();
127 
128  scalar mP = cloud.constProps(typeIdP).mass();
129 
130  scalar mQ = cloud.constProps(typeIdQ).mass();
131 
132  vector Ucm = (mP*UP + mQ*UQ)/(mP + mQ);
133 
134  scalar cR = mag(UP - UQ);
135 
136  scalar cosTheta = 2.0*rndGen.sample01<scalar>() - 1.0;
137 
138  scalar sinTheta = sqrt(1.0 - cosTheta*cosTheta);
139 
140  scalar phi = twoPi*rndGen.sample01<scalar>();
141 
142  vector postCollisionRelU =
143  cR
144  *vector
145  (
146  cosTheta,
147  sinTheta*cos(phi),
148  sinTheta*sin(phi)
149  );
150 
151  UP = Ucm + postCollisionRelU*mQ/(mP + mQ);
152 
153  UQ = Ucm - postCollisionRelU*mP/(mP + mQ);
154 }
155 
156 
157 // ************************************************************************* //
Foam::Random
Random number generator.
Definition: Random.H:59
Foam::constant::physicoChemical::k
const dimensionedScalar k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::BinaryCollisionModel
Templated DSMC particle collision class.
Definition: DSMCCloud.H:58
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
VariableHardSphere.H
Foam::lgamma
dimensionedScalar lgamma(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:278
Foam::VariableHardSphere::sigmaTcR
virtual scalar sigmaTcR(const typename CloudType::parcelType &pP, const typename CloudType::parcelType &pQ) const
Return the collision cross section * relative velocity product.
Definition: VariableHardSphere.C:65
Foam::VariableHardSphere::VariableHardSphere
VariableHardSphere(const dictionary &dict, CloudType &cloud)
Construct from dictionary.
Definition: VariableHardSphere.C:37
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:71
Foam::constant::mathematical::twoPi
constexpr scalar twoPi(2 *M_PI)
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::VariableHardSphere::collide
virtual void collide(typename CloudType::parcelType &pP, typename CloudType::parcelType &pQ)
Apply collision.
Definition: VariableHardSphere.C:114
constants.H
Foam::cloud
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:57
Foam::constant::mathematical
Mathematical constants.
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::Vector< scalar >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
rndGen
Random rndGen
Definition: createFields.H:23
Foam::DSMCCloud::parcelType
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
Foam::VariableHardSphere::~VariableHardSphere
virtual ~VariableHardSphere()
Destructor.
Definition: VariableHardSphere.C:50
Foam::VariableHardSphere::active
virtual bool active() const
Flag to indicate whether model activates collision model.
Definition: VariableHardSphere.C:57
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265