1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2009-2009 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 option) any later version.
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
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM; if not, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
29 template<class ParcelType>
30 inline const Foam::word& Foam::DsmcCloud<ParcelType>::cloudName() const
36 template<class ParcelType>
37 inline const Foam::fvMesh& Foam::DsmcCloud<ParcelType>::mesh() const
43 template<class ParcelType>
44 inline const Foam::IOdictionary&
45 Foam::DsmcCloud<ParcelType>::particleProperties() const
47 return particleProperties_;
51 template<class ParcelType>
52 inline const Foam::List<Foam::word>&
53 Foam::DsmcCloud<ParcelType>::typeIdList() const
59 template<class ParcelType>
60 inline Foam::scalar Foam::DsmcCloud<ParcelType>::nParticle() const
66 template<class ParcelType>
67 inline const Foam::List<Foam::DynamicList<ParcelType*> >&
68 Foam::DsmcCloud<ParcelType>::cellOccupancy() const
70 return cellOccupancy_;
74 template<class ParcelType>
75 inline Foam::volScalarField& Foam::DsmcCloud<ParcelType>::sigmaTcRMax()
81 template<class ParcelType>
82 inline Foam::scalarField&
83 Foam::DsmcCloud<ParcelType>::collisionSelectionRemainder()
85 return collisionSelectionRemainder_;
89 template<class ParcelType>
90 inline const Foam::List<typename ParcelType::constantProperties>&
91 Foam::DsmcCloud<ParcelType>::constProps() const
97 template<class ParcelType>
98 inline const typename ParcelType::constantProperties&
99 Foam::DsmcCloud<ParcelType>::constProps
104 if (typeId < 0 || typeId >= constProps_.size())
106 FatalErrorIn("Foam::DsmcCloud<ParcelType>::constProps(label typeId)")
107 << "constantProperties for requested typeId index "
108 << typeId << " do not exist" << nl
109 << abort(FatalError);
112 return constProps_[typeId];
116 template<class ParcelType>
117 inline Foam::Random& Foam::DsmcCloud<ParcelType>::rndGen()
123 template<class ParcelType>
124 inline void Foam::DsmcCloud<ParcelType>::storeDeltaT()
126 cachedDeltaT_ = mesh().time().deltaT().value();
130 template<class ParcelType>
131 inline Foam::scalar Foam::DsmcCloud<ParcelType>::cachedDeltaT() const
133 return cachedDeltaT_;
137 template<class ParcelType>
138 inline const Foam::volScalarField& Foam::DsmcCloud<ParcelType>::q() const
144 template<class ParcelType>
145 inline Foam::volScalarField& Foam::DsmcCloud<ParcelType>::q()
151 template<class ParcelType>
152 inline const Foam::volVectorField& Foam::DsmcCloud<ParcelType>::fD() const
158 template<class ParcelType>
159 inline Foam::volVectorField& Foam::DsmcCloud<ParcelType>::fD()
165 template<class ParcelType>
166 inline const Foam::volScalarField& Foam::DsmcCloud<ParcelType>::T() const
172 template<class ParcelType>
173 inline const Foam::volVectorField& Foam::DsmcCloud<ParcelType>::U() const
179 template<class ParcelType>
180 inline const Foam::BinaryCollisionModel<Foam::DsmcCloud<ParcelType> >&
181 Foam::DsmcCloud<ParcelType>::binaryCollision() const
183 return binaryCollisionModel_;
187 template<class ParcelType>
188 inline Foam::BinaryCollisionModel<Foam::DsmcCloud<ParcelType> >&
189 Foam::DsmcCloud<ParcelType>::binaryCollision()
191 return binaryCollisionModel_();
195 template<class ParcelType>
196 inline const Foam::WallInteractionModel<Foam::DsmcCloud<ParcelType> >&
197 Foam::DsmcCloud<ParcelType>::wallInteraction() const
199 return wallInteractionModel_;
203 template<class ParcelType>
204 inline Foam::WallInteractionModel<Foam::DsmcCloud<ParcelType> >&
205 Foam::DsmcCloud<ParcelType>::wallInteraction()
207 return wallInteractionModel_();
211 template<class ParcelType>
212 inline const Foam::InflowBoundaryModel<Foam::DsmcCloud<ParcelType> >&
213 Foam::DsmcCloud<ParcelType>::inflowBoundary() const
215 return inflowBoundaryModel_;
219 template<class ParcelType>
220 inline Foam::InflowBoundaryModel<Foam::DsmcCloud<ParcelType> >&
221 Foam::DsmcCloud<ParcelType>::inflowBoundary()
223 return inflowBoundaryModel_();
227 template<class ParcelType>
228 inline Foam::scalar Foam::DsmcCloud<ParcelType>::massInSystem() const
230 scalar sysMass = 0.0;
232 forAllConstIter(typename DsmcCloud<ParcelType>, *this, iter)
234 const ParcelType& p = iter();
236 const typename ParcelType::constantProperties& cP = constProps
241 sysMass += cP.mass();
244 return nParticle_*sysMass;
248 template<class ParcelType>
249 inline Foam::vector Foam::DsmcCloud<ParcelType>::linearMomentumOfSystem() const
251 vector linearMomentum(vector::zero);
253 forAllConstIter(typename DsmcCloud<ParcelType>, *this, iter)
255 const ParcelType& p = iter();
257 const typename ParcelType::constantProperties& cP = constProps
262 linearMomentum += cP.mass()*p.U();
265 return nParticle_*linearMomentum;
269 template<class ParcelType>
271 Foam::DsmcCloud<ParcelType>::linearKineticEnergyOfSystem() const
273 scalar linearKineticEnergy = 0.0;
275 forAllConstIter(typename DsmcCloud<ParcelType>, *this, iter)
277 const ParcelType& p = iter();
279 const typename ParcelType::constantProperties& cP = constProps
284 linearKineticEnergy += 0.5*cP.mass()*(p.U() & p.U());
287 return nParticle_*linearKineticEnergy;
291 template<class ParcelType>
293 Foam::DsmcCloud<ParcelType>::internalEnergyOfSystem() const
295 scalar internalEnergy = 0.0;
297 forAllConstIter(typename DsmcCloud<ParcelType>, *this, iter)
299 const ParcelType& p = iter();
301 internalEnergy += p.Ei();
304 return nParticle_*internalEnergy;
308 template<class ParcelType>
309 inline Foam::scalar Foam::DsmcCloud<ParcelType>::maxwellianAverageSpeed
315 return 2.0*sqrt(2.0*kb*temperature/(mathematicalConstant::pi*mass));
319 template<class ParcelType>
320 inline Foam::scalarField Foam::DsmcCloud<ParcelType>::maxwellianAverageSpeed
322 scalarField temperature,
326 return 2.0*sqrt(2.0*kb*temperature/(mathematicalConstant::pi*mass));
330 template<class ParcelType>
331 inline Foam::scalar Foam::DsmcCloud<ParcelType>::maxwellianRMSSpeed
337 return sqrt(3.0*kb*temperature/mass);
341 template<class ParcelType>
342 inline Foam::scalarField Foam::DsmcCloud<ParcelType>::maxwellianRMSSpeed
344 scalarField temperature,
348 return sqrt(3.0*kb*temperature/mass);
352 template<class ParcelType>
354 Foam::DsmcCloud<ParcelType>::maxwellianMostProbableSpeed
360 return sqrt(2.0*kb*temperature/mass);
364 template<class ParcelType>
365 inline Foam::scalarField
366 Foam::DsmcCloud<ParcelType>::maxwellianMostProbableSpeed
368 scalarField temperature,
372 return sqrt(2.0*kb*temperature/mass);
376 template<class ParcelType>
377 inline const Foam::tmp<Foam::volScalarField>
378 Foam::DsmcCloud<ParcelType>::rhoN() const
380 tmp<volScalarField> trhoN
386 this->name() + "rhoN",
387 this->db().time().timeName(),
394 dimensionedScalar("zero", dimensionSet(0, -3, 0, 0, 0), VSMALL)
398 scalarField& rhoN = trhoN().internalField();
399 forAllConstIter(typename DsmcCloud<ParcelType>, *this, iter)
401 const ParcelType& p = iter();
402 const label cellI = p.cell();
407 rhoN *= nParticle_/mesh().cellVolumes();
413 template<class ParcelType>
414 inline const Foam::tmp<Foam::volScalarField>
415 Foam::DsmcCloud<ParcelType>::rhoM() const
417 tmp<volScalarField> trhoM
423 this->name() + "rhoM",
424 this->db().time().timeName(),
431 dimensionedScalar("zero", dimensionSet(1, -3, 0, 0, 0), VSMALL)
435 scalarField& rhoM = trhoM().internalField();
436 forAllConstIter(typename DsmcCloud<ParcelType>, *this, iter)
438 const ParcelType& p = iter();
439 const label cellI = p.cell();
441 rhoM[cellI] += constProps(p.typeId()).mass();
444 rhoM *= nParticle_/mesh().cellVolumes();
450 template<class ParcelType>
451 inline const Foam::tmp<Foam::volScalarField>
452 Foam::DsmcCloud<ParcelType>::dsmcRhoN() const
454 tmp<volScalarField> tdsmcRhoN
460 this->name() + "dsmcRhoN",
461 this->db().time().timeName(),
468 dimensionedScalar("zero", dimensionSet(0, -3, 0, 0, 0), 0.0)
472 scalarField& dsmcRhoN = tdsmcRhoN().internalField();
473 forAllConstIter(typename DsmcCloud<ParcelType>, *this, iter)
475 const ParcelType& p = iter();
476 const label cellI = p.cell();
485 template<class ParcelType>
486 inline const Foam::tmp<Foam::volVectorField>
487 Foam::DsmcCloud<ParcelType>::momentum() const
489 tmp<volVectorField> tmomentum
495 this->name() + "momentum",
496 this->db().time().timeName(),
506 dimensionSet(1, -2, -1, 0, 0),
512 vectorField& momentum = tmomentum().internalField();
513 forAllConstIter(typename DsmcCloud<ParcelType>, *this, iter)
515 const ParcelType& p = iter();
516 const label cellI = p.cell();
518 momentum[cellI] += constProps(p.typeId()).mass()*p.U();
521 momentum *= nParticle_/mesh().cellVolumes();
527 template<class ParcelType>
528 inline const Foam::tmp<Foam::volScalarField>
529 Foam::DsmcCloud<ParcelType>::linearKE() const
531 tmp<volScalarField> tlinearKE
537 this->name() + "linearKE",
538 this->db().time().timeName(),
545 dimensionedScalar("zero", dimensionSet(1, -1, -2, 0, 0), 0.0)
549 scalarField& linearKE = tlinearKE().internalField();
550 forAllConstIter(typename DsmcCloud<ParcelType>, *this, iter)
552 const ParcelType& p = iter();
553 const label cellI = p.cell();
555 linearKE[cellI] += 0.5*constProps(p.typeId()).mass()*(p.U() & p.U());
558 linearKE *= nParticle_/mesh().cellVolumes();
564 template<class ParcelType>
565 inline const Foam::tmp<Foam::volScalarField>
566 Foam::DsmcCloud<ParcelType>::internalE() const
568 tmp<volScalarField> tinternalE
574 this->name() + "internalE",
575 this->db().time().timeName(),
582 dimensionedScalar("zero", dimensionSet(1, -1, -2, 0, 0), 0.0)
586 scalarField& internalE = tinternalE().internalField();
587 forAllConstIter(typename DsmcCloud<ParcelType>, *this, iter)
589 const ParcelType& p = iter();
590 const label cellI = p.cell();
592 internalE[cellI] += p.Ei();
595 internalE *= nParticle_/mesh().cellVolumes();
601 template<class ParcelType>
602 inline const Foam::tmp<Foam::volScalarField>
603 Foam::DsmcCloud<ParcelType>::iDof() const
605 tmp<volScalarField> tiDof
611 this->name() + "iDof",
612 this->db().time().timeName(),
619 dimensionedScalar("zero", dimensionSet(0, -3, 0, 0, 0), VSMALL)
623 scalarField& iDof = tiDof().internalField();
625 forAllConstIter(typename DsmcCloud<ParcelType>, *this, iter)
627 const ParcelType& p = iter();
628 const label cellI = p.cell();
630 iDof[cellI] += constProps(p.typeId()).internalDegreesOfFreedom();
633 iDof *= nParticle_/mesh().cellVolumes();
639 template<class ParcelType>
640 inline void Foam::DsmcCloud<ParcelType>::clear()
642 return IDLList<ParcelType>::clear();
646 // ************************************************************************* //