FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

LarsenBorgnakkeVariableHardSphere.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 2009-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "LarsenBorgnakkeVariableHardSphere.H"
00027 
00028 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00029 
00030 template <class CloudType>
00031 Foam::scalar Foam::LarsenBorgnakkeVariableHardSphere<CloudType>::energyRatio
00032 (
00033     scalar ChiA,
00034     scalar ChiB
00035 )
00036 {
00037     CloudType& cloud(this->owner());
00038 
00039     Random& rndGen(cloud.rndGen());
00040 
00041     scalar ChiAMinusOne = ChiA - 1;
00042 
00043     scalar ChiBMinusOne = ChiB - 1;
00044 
00045     if (ChiAMinusOne < SMALL && ChiBMinusOne < SMALL)
00046     {
00047         return rndGen.scalar01();
00048     }
00049 
00050     scalar energyRatio;
00051 
00052     scalar P;
00053 
00054     do
00055     {
00056         P = 0;
00057 
00058         energyRatio = rndGen.scalar01();
00059 
00060         if (ChiAMinusOne < SMALL)
00061         {
00062             P = 1.0 - pow(energyRatio, ChiB);
00063         }
00064         else if (ChiBMinusOne < SMALL)
00065         {
00066             P = 1.0 - pow(energyRatio, ChiA);
00067         }
00068         else
00069         {
00070             P =
00071                 pow
00072                 (
00073                     (ChiAMinusOne + ChiBMinusOne)*energyRatio/ChiAMinusOne,
00074                     ChiAMinusOne
00075                 )
00076                *pow
00077                 (
00078                     (ChiAMinusOne + ChiBMinusOne)*(1 - energyRatio)
00079                     /ChiBMinusOne,
00080                     ChiBMinusOne
00081                 );
00082         }
00083     } while (P < rndGen.scalar01());
00084 
00085     return energyRatio;
00086 }
00087 
00088 
00089 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00090 
00091 template <class CloudType>
00092 Foam::LarsenBorgnakkeVariableHardSphere<CloudType>::LarsenBorgnakkeVariableHardSphere
00093 (
00094     const dictionary& dict,
00095     CloudType& cloud
00096 )
00097 :
00098     BinaryCollisionModel<CloudType>(dict, cloud, typeName),
00099     Tref_(readScalar(this->coeffDict().lookup("Tref"))),
00100     relaxationCollisionNumber_
00101     (
00102         readScalar(this->coeffDict().lookup("relaxationCollisionNumber"))
00103     )
00104 {}
00105 
00106 
00107 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00108 
00109 template <class CloudType>
00110 Foam::LarsenBorgnakkeVariableHardSphere<CloudType>::
00111 ~LarsenBorgnakkeVariableHardSphere()
00112 {}
00113 
00114 
00115 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00116 
00117 
00118 template <class CloudType>
00119 Foam::scalar Foam::LarsenBorgnakkeVariableHardSphere<CloudType>::sigmaTcR
00120 (
00121     label typeIdP,
00122     label typeIdQ,
00123     const vector& UP,
00124     const vector& UQ
00125 ) const
00126 {
00127     const CloudType& cloud(this->owner());
00128 
00129     scalar dPQ =
00130         0.5
00131        *(
00132             cloud.constProps(typeIdP).d()
00133           + cloud.constProps(typeIdQ).d()
00134         );
00135 
00136     scalar omegaPQ =
00137         0.5
00138        *(
00139             cloud.constProps(typeIdP).omega()
00140           + cloud.constProps(typeIdQ).omega()
00141         );
00142 
00143     scalar cR = mag(UP - UQ);
00144 
00145     if (cR < VSMALL)
00146     {
00147         return 0;
00148     }
00149 
00150     scalar mP = cloud.constProps(typeIdP).mass();
00151 
00152     scalar mQ = cloud.constProps(typeIdQ).mass();
00153 
00154     scalar mR = mP*mQ/(mP + mQ);
00155 
00156     // calculating cross section = pi*dPQ^2, where dPQ is from Bird, eq. 4.79
00157     scalar sigmaTPQ =
00158         mathematicalConstant::pi*dPQ*dPQ
00159        *pow(2.0*CloudType::kb*Tref_/(mR*cR*cR), omegaPQ - 0.5)
00160        /exp(Foam::lgamma(2.5 - omegaPQ));
00161 
00162     return sigmaTPQ*cR;
00163 }
00164 
00165 
00166 template <class CloudType>
00167 void Foam::LarsenBorgnakkeVariableHardSphere<CloudType>::collide
00168 (
00169     label typeIdP,
00170     label typeIdQ,
00171     vector& UP,
00172     vector& UQ,
00173     scalar& EiP,
00174     scalar& EiQ
00175 )
00176 {
00177     CloudType& cloud(this->owner());
00178 
00179     Random& rndGen(cloud.rndGen());
00180 
00181     scalar inverseCollisionNumber = 1/relaxationCollisionNumber_;
00182 
00183     // Larsen Borgnakke internal energy redistribution part.  Using the serial
00184     // application of the LB method, as per the INELRS subroutine in Bird's
00185     // DSMC0R.FOR
00186 
00187     scalar preCollisionEiP = EiP;
00188 
00189     scalar preCollisionEiQ = EiQ;
00190 
00191     scalar iDofP = cloud.constProps(typeIdP).internalDegreesOfFreedom();
00192 
00193     scalar iDofQ = cloud.constProps(typeIdQ).internalDegreesOfFreedom();
00194 
00195     scalar omegaPQ =
00196         0.5
00197        *(
00198             cloud.constProps(typeIdP).omega()
00199           + cloud.constProps(typeIdQ).omega()
00200         );
00201 
00202     scalar mP = cloud.constProps(typeIdP).mass();
00203 
00204     scalar mQ = cloud.constProps(typeIdQ).mass();
00205 
00206     scalar mR = mP*mQ/(mP + mQ);
00207 
00208     vector Ucm = (mP*UP + mQ*UQ)/(mP + mQ);
00209 
00210     scalar cRsqr = magSqr(UP - UQ);
00211 
00212     scalar availableEnergy = 0.5*mR*cRsqr;
00213 
00214     scalar ChiB = 2.5 - omegaPQ;
00215 
00216     if (iDofP > 0)
00217     {
00218         if (inverseCollisionNumber > rndGen.scalar01())
00219         {
00220             availableEnergy += preCollisionEiP;
00221 
00222             scalar ChiA = 0.5*iDofP;
00223 
00224             EiP = energyRatio(ChiA, ChiB)*availableEnergy;
00225 
00226             availableEnergy -= EiP;
00227         }
00228     }
00229 
00230     if (iDofQ > 0)
00231     {
00232         if (inverseCollisionNumber > rndGen.scalar01())
00233         {
00234             availableEnergy += preCollisionEiQ;
00235 
00236             // Change to general LB ratio calculation
00237             scalar energyRatio = 1.0 - pow(rndGen.scalar01(),(1.0/ChiB));
00238 
00239             EiQ = energyRatio*availableEnergy;
00240 
00241             availableEnergy -= EiQ;
00242         }
00243     }
00244 
00245     // Rescale the translational energy
00246     scalar cR = sqrt(2.0*availableEnergy/mR);
00247 
00248     // Variable Hard Sphere collision part
00249 
00250     scalar cosTheta = 2.0*rndGen.scalar01() - 1.0;
00251 
00252     scalar sinTheta = sqrt(1.0 - cosTheta*cosTheta);
00253 
00254     scalar phi = 2.0*mathematicalConstant::pi*rndGen.scalar01();
00255 
00256     vector postCollisionRelU =
00257         cR
00258        *vector
00259         (
00260             cosTheta,
00261             sinTheta*cos(phi),
00262             sinTheta*sin(phi)
00263         );
00264 
00265     UP = Ucm + postCollisionRelU*mQ/(mP + mQ);
00266 
00267     UQ = Ucm - postCollisionRelU*mP/(mP + mQ);
00268 }
00269 
00270 
00271 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines