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

DsmcCloudI_.H

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 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00027 
00028 template<class ParcelType>
00029 inline const Foam::word& Foam::DsmcCloud<ParcelType>::cloudName() const
00030 {
00031     return cloudName_;
00032 }
00033 
00034 
00035 template<class ParcelType>
00036 inline const Foam::fvMesh& Foam::DsmcCloud<ParcelType>::mesh() const
00037 {
00038     return mesh_;
00039 }
00040 
00041 
00042 template<class ParcelType>
00043 inline const Foam::IOdictionary&
00044 Foam::DsmcCloud<ParcelType>::particleProperties() const
00045 {
00046     return particleProperties_;
00047 }
00048 
00049 
00050 template<class ParcelType>
00051 inline const Foam::List<Foam::word>&
00052 Foam::DsmcCloud<ParcelType>::typeIdList() const
00053 {
00054     return typeIdList_;
00055 }
00056 
00057 
00058 template<class ParcelType>
00059 inline Foam::scalar Foam::DsmcCloud<ParcelType>::nParticle() const
00060 {
00061     return nParticle_;
00062 }
00063 
00064 
00065 template<class ParcelType>
00066 inline const Foam::List<Foam::DynamicList<ParcelType*> >&
00067 Foam::DsmcCloud<ParcelType>::cellOccupancy() const
00068 {
00069     return cellOccupancy_;
00070 }
00071 
00072 
00073 template<class ParcelType>
00074 inline Foam::volScalarField& Foam::DsmcCloud<ParcelType>::sigmaTcRMax()
00075 {
00076     return sigmaTcRMax_;
00077 }
00078 
00079 
00080 template<class ParcelType>
00081 inline Foam::scalarField&
00082 Foam::DsmcCloud<ParcelType>::collisionSelectionRemainder()
00083 {
00084     return collisionSelectionRemainder_;
00085 }
00086 
00087 
00088 template<class ParcelType>
00089 inline const Foam::List<typename ParcelType::constantProperties>&
00090 Foam::DsmcCloud<ParcelType>::constProps() const
00091 {
00092     return constProps_;
00093 }
00094 
00095 
00096 template<class ParcelType>
00097 inline const typename ParcelType::constantProperties&
00098 Foam::DsmcCloud<ParcelType>::constProps
00099 (
00100     label typeId
00101 ) const
00102 {
00103     if (typeId < 0 || typeId >= constProps_.size())
00104     {
00105         FatalErrorIn("Foam::DsmcCloud<ParcelType>::constProps(label typeId)")
00106             << "constantProperties for requested typeId index "
00107             << typeId << " do not exist" << nl
00108             << abort(FatalError);
00109     }
00110 
00111     return constProps_[typeId];
00112 }
00113 
00114 
00115 template<class ParcelType>
00116 inline Foam::Random& Foam::DsmcCloud<ParcelType>::rndGen()
00117 {
00118     return rndGen_;
00119 }
00120 
00121 
00122 template<class ParcelType>
00123 inline Foam::volScalarField::GeometricBoundaryField&
00124 Foam::DsmcCloud<ParcelType>::qBF()
00125 {
00126     return q_.boundaryField();
00127 }
00128 
00129 
00130 template<class ParcelType>
00131 inline Foam::volVectorField::GeometricBoundaryField&
00132 Foam::DsmcCloud<ParcelType>::fDBF()
00133 {
00134     return fD_.boundaryField();
00135 }
00136 
00137 
00138 template<class ParcelType>
00139 inline Foam::volScalarField::GeometricBoundaryField&
00140 Foam::DsmcCloud<ParcelType>::rhoNBF()
00141 {
00142     return rhoN_.boundaryField();
00143 }
00144 
00145 
00146 template<class ParcelType>
00147 inline Foam::volScalarField::GeometricBoundaryField&
00148 Foam::DsmcCloud<ParcelType>::rhoMBF()
00149 {
00150     return rhoM_.boundaryField();
00151 }
00152 
00153 
00154 template<class ParcelType>
00155 inline Foam::volScalarField::GeometricBoundaryField&
00156 Foam::DsmcCloud<ParcelType>::linearKEBF()
00157 {
00158     return linearKE_.boundaryField();
00159 }
00160 
00161 
00162 template<class ParcelType>
00163 inline Foam::volScalarField::GeometricBoundaryField&
00164 Foam::DsmcCloud<ParcelType>::internalEBF()
00165 {
00166     return internalE_.boundaryField();
00167 }
00168 
00169 
00170 template<class ParcelType>
00171 inline Foam::volScalarField::GeometricBoundaryField&
00172 Foam::DsmcCloud<ParcelType>::iDofBF()
00173 {
00174     return iDof_.boundaryField();
00175 }
00176 
00177 
00178 template<class ParcelType>
00179 inline Foam::volVectorField::GeometricBoundaryField&
00180 Foam::DsmcCloud<ParcelType>::momentumBF()
00181 {
00182     return momentum_.boundaryField();
00183 }
00184 
00185 
00186 template<class ParcelType>
00187 inline const Foam::volScalarField&
00188 Foam::DsmcCloud<ParcelType>::boundaryT() const
00189 {
00190     return boundaryT_;
00191 }
00192 
00193 
00194 template<class ParcelType>
00195 inline const Foam::volVectorField&
00196 Foam::DsmcCloud<ParcelType>::boundaryU() const
00197 {
00198     return boundaryU_;
00199 }
00200 
00201 
00202 template<class ParcelType>
00203 inline const Foam::BinaryCollisionModel<Foam::DsmcCloud<ParcelType> >&
00204 Foam::DsmcCloud<ParcelType>::binaryCollision() const
00205 {
00206     return binaryCollisionModel_;
00207 }
00208 
00209 
00210 template<class ParcelType>
00211 inline Foam::BinaryCollisionModel<Foam::DsmcCloud<ParcelType> >&
00212 Foam::DsmcCloud<ParcelType>::binaryCollision()
00213 {
00214     return binaryCollisionModel_();
00215 }
00216 
00217 
00218 template<class ParcelType>
00219 inline const Foam::WallInteractionModel<Foam::DsmcCloud<ParcelType> >&
00220 Foam::DsmcCloud<ParcelType>::wallInteraction() const
00221 {
00222     return wallInteractionModel_;
00223 }
00224 
00225 
00226 template<class ParcelType>
00227 inline Foam::WallInteractionModel<Foam::DsmcCloud<ParcelType> >&
00228 Foam::DsmcCloud<ParcelType>::wallInteraction()
00229 {
00230     return wallInteractionModel_();
00231 }
00232 
00233 
00234 template<class ParcelType>
00235 inline const Foam::InflowBoundaryModel<Foam::DsmcCloud<ParcelType> >&
00236 Foam::DsmcCloud<ParcelType>::inflowBoundary() const
00237 {
00238     return inflowBoundaryModel_;
00239 }
00240 
00241 
00242 template<class ParcelType>
00243 inline Foam::InflowBoundaryModel<Foam::DsmcCloud<ParcelType> >&
00244 Foam::DsmcCloud<ParcelType>::inflowBoundary()
00245 {
00246     return inflowBoundaryModel_();
00247 }
00248 
00249 
00250 template<class ParcelType>
00251 inline Foam::scalar Foam::DsmcCloud<ParcelType>::massInSystem() const
00252 {
00253     scalar sysMass = 0.0;
00254 
00255     forAllConstIter(typename DsmcCloud<ParcelType>, *this, iter)
00256     {
00257         const ParcelType& p = iter();
00258 
00259         const typename ParcelType::constantProperties& cP = constProps
00260         (
00261             p.typeId()
00262         );
00263 
00264         sysMass += cP.mass();
00265     }
00266 
00267     return nParticle_*sysMass;
00268 }
00269 
00270 
00271 template<class ParcelType>
00272 inline Foam::vector Foam::DsmcCloud<ParcelType>::linearMomentumOfSystem() const
00273 {
00274     vector linearMomentum(vector::zero);
00275 
00276     forAllConstIter(typename DsmcCloud<ParcelType>, *this, iter)
00277     {
00278         const ParcelType& p = iter();
00279 
00280         const typename ParcelType::constantProperties& cP = constProps
00281         (
00282             p.typeId()
00283         );
00284 
00285         linearMomentum += cP.mass()*p.U();
00286     }
00287 
00288     return nParticle_*linearMomentum;
00289 }
00290 
00291 
00292 template<class ParcelType>
00293 inline Foam::scalar
00294 Foam::DsmcCloud<ParcelType>::linearKineticEnergyOfSystem() const
00295 {
00296     scalar linearKineticEnergy = 0.0;
00297 
00298     forAllConstIter(typename DsmcCloud<ParcelType>, *this, iter)
00299     {
00300         const ParcelType& p = iter();
00301 
00302         const typename ParcelType::constantProperties& cP = constProps
00303         (
00304             p.typeId()
00305         );
00306 
00307         linearKineticEnergy += 0.5*cP.mass()*(p.U() & p.U());
00308     }
00309 
00310     return nParticle_*linearKineticEnergy;
00311 }
00312 
00313 
00314 template<class ParcelType>
00315 inline Foam::scalar
00316 Foam::DsmcCloud<ParcelType>::internalEnergyOfSystem() const
00317 {
00318     scalar internalEnergy = 0.0;
00319 
00320     forAllConstIter(typename DsmcCloud<ParcelType>, *this, iter)
00321     {
00322         const ParcelType& p = iter();
00323 
00324         internalEnergy += p.Ei();
00325     }
00326 
00327     return nParticle_*internalEnergy;
00328 }
00329 
00330 
00331 template<class ParcelType>
00332 inline Foam::scalar Foam::DsmcCloud<ParcelType>::maxwellianAverageSpeed
00333 (
00334     scalar temperature,
00335     scalar mass
00336 ) const
00337 {
00338     return
00339         2.0*sqrt(2.0*kb*temperature/(mathematicalConstant::pi*mass));
00340 }
00341 
00342 
00343 template<class ParcelType>
00344 inline Foam::scalarField Foam::DsmcCloud<ParcelType>::maxwellianAverageSpeed
00345 (
00346     scalarField temperature,
00347     scalar mass
00348 ) const
00349 {
00350     return
00351         2.0*sqrt(2.0*kb*temperature/(mathematicalConstant::pi*mass));
00352 }
00353 
00354 
00355 template<class ParcelType>
00356 inline Foam::scalar Foam::DsmcCloud<ParcelType>::maxwellianRMSSpeed
00357 (
00358     scalar temperature,
00359     scalar mass
00360 ) const
00361 {
00362     return sqrt(3.0*kb*temperature/mass);
00363 }
00364 
00365 
00366 template<class ParcelType>
00367 inline Foam::scalarField Foam::DsmcCloud<ParcelType>::maxwellianRMSSpeed
00368 (
00369     scalarField temperature,
00370     scalar mass
00371 ) const
00372 {
00373     return sqrt(3.0*kb*temperature/mass);
00374 }
00375 
00376 
00377 template<class ParcelType>
00378 inline Foam::scalar
00379 Foam::DsmcCloud<ParcelType>::maxwellianMostProbableSpeed
00380 (
00381     scalar temperature,
00382     scalar mass
00383 ) const
00384 {
00385     return sqrt(2.0*kb*temperature/mass);
00386 }
00387 
00388 
00389 template<class ParcelType>
00390 inline Foam::scalarField
00391 Foam::DsmcCloud<ParcelType>::maxwellianMostProbableSpeed
00392 (
00393     scalarField temperature,
00394     scalar mass
00395 ) const
00396 {
00397     return sqrt(2.0*kb*temperature/mass);
00398 }
00399 
00400 
00401 template<class ParcelType>
00402 inline const Foam::volScalarField& Foam::DsmcCloud<ParcelType>::q() const
00403 {
00404     return q_;
00405 }
00406 
00407 
00408 template<class ParcelType>
00409 inline const Foam::volVectorField& Foam::DsmcCloud<ParcelType>::fD() const
00410 {
00411     return fD_;
00412 }
00413 
00414 
00415 template<class ParcelType>
00416 inline const Foam::volScalarField&
00417 Foam::DsmcCloud<ParcelType>::rhoN() const
00418 {
00419     return rhoN_;
00420 }
00421 
00422 
00423 template<class ParcelType>
00424 inline const Foam::volScalarField& Foam::DsmcCloud<ParcelType>::rhoM() const
00425 {
00426     return rhoM_;
00427 }
00428 
00429 
00430 template<class ParcelType>
00431 inline const Foam::volScalarField&
00432 Foam::DsmcCloud<ParcelType>::dsmcRhoN() const
00433 {
00434     return dsmcRhoN_;
00435 }
00436 
00437 
00438 template<class ParcelType>
00439 inline const Foam::volScalarField&
00440 Foam::DsmcCloud<ParcelType>::linearKE() const
00441 {
00442     return linearKE_;
00443 }
00444 
00445 
00446 template<class ParcelType>
00447 inline const Foam::volScalarField&
00448 Foam::DsmcCloud<ParcelType>::internalE() const
00449 {
00450     return internalE_;
00451 }
00452 
00453 
00454 template<class ParcelType>
00455 inline const Foam::volScalarField&
00456 Foam::DsmcCloud<ParcelType>::iDof() const
00457 {
00458     return iDof_;
00459 }
00460 
00461 
00462 template<class ParcelType>
00463 inline const Foam::volVectorField& Foam::DsmcCloud<ParcelType>::momentum() const
00464 {
00465     return momentum_;
00466 }
00467 
00468 
00469 template<class ParcelType>
00470 inline void Foam::DsmcCloud<ParcelType>::clear()
00471 {
00472     return IDLList<ParcelType>::clear();
00473 }
00474 
00475 
00476 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines