Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / lagrangian / dsmc / clouds / Templates / DsmcCloud / DsmcCloudI.H
blob57b48986a134fc3972eccb398a2250457bca0f2e
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2009-2010 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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
19     for more details.
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "constants.H"
28 using namespace Foam::constant;
29 using namespace Foam::constant::mathematical;
31 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
33 template<class ParcelType>
34 inline const Foam::word& Foam::DsmcCloud<ParcelType>::cloudName() const
36     return cloudName_;
40 template<class ParcelType>
41 inline const Foam::fvMesh& Foam::DsmcCloud<ParcelType>::mesh() const
43     return mesh_;
47 template<class ParcelType>
48 inline const Foam::IOdictionary&
49 Foam::DsmcCloud<ParcelType>::particleProperties() const
51     return particleProperties_;
55 template<class ParcelType>
56 inline const Foam::List<Foam::word>&
57 Foam::DsmcCloud<ParcelType>::typeIdList() const
59     return typeIdList_;
63 template<class ParcelType>
64 inline Foam::scalar Foam::DsmcCloud<ParcelType>::nParticle() const
66     return nParticle_;
70 template<class ParcelType>
71 inline const Foam::List<Foam::DynamicList<ParcelType*> >&
72 Foam::DsmcCloud<ParcelType>::cellOccupancy() const
74     return cellOccupancy_;
78 template<class ParcelType>
79 inline Foam::volScalarField& Foam::DsmcCloud<ParcelType>::sigmaTcRMax()
81     return sigmaTcRMax_;
85 template<class ParcelType>
86 inline Foam::scalarField&
87 Foam::DsmcCloud<ParcelType>::collisionSelectionRemainder()
89     return collisionSelectionRemainder_;
93 template<class ParcelType>
94 inline const Foam::List<typename ParcelType::constantProperties>&
95 Foam::DsmcCloud<ParcelType>::constProps() const
97     return constProps_;
101 template<class ParcelType>
102 inline const typename ParcelType::constantProperties&
103 Foam::DsmcCloud<ParcelType>::constProps
105     label typeId
106 ) const
108     if (typeId < 0 || typeId >= constProps_.size())
109     {
110         FatalErrorIn("Foam::DsmcCloud<ParcelType>::constProps(label typeId)")
111             << "constantProperties for requested typeId index "
112             << typeId << " do not exist" << nl
113             << abort(FatalError);
114     }
116     return constProps_[typeId];
120 template<class ParcelType>
121 inline Foam::Random& Foam::DsmcCloud<ParcelType>::rndGen()
123     return rndGen_;
127 template<class ParcelType>
128 inline Foam::volScalarField::GeometricBoundaryField&
129 Foam::DsmcCloud<ParcelType>::qBF()
131     return q_.boundaryField();
135 template<class ParcelType>
136 inline Foam::volVectorField::GeometricBoundaryField&
137 Foam::DsmcCloud<ParcelType>::fDBF()
139     return fD_.boundaryField();
143 template<class ParcelType>
144 inline Foam::volScalarField::GeometricBoundaryField&
145 Foam::DsmcCloud<ParcelType>::rhoNBF()
147     return rhoN_.boundaryField();
151 template<class ParcelType>
152 inline Foam::volScalarField::GeometricBoundaryField&
153 Foam::DsmcCloud<ParcelType>::rhoMBF()
155     return rhoM_.boundaryField();
159 template<class ParcelType>
160 inline Foam::volScalarField::GeometricBoundaryField&
161 Foam::DsmcCloud<ParcelType>::linearKEBF()
163     return linearKE_.boundaryField();
167 template<class ParcelType>
168 inline Foam::volScalarField::GeometricBoundaryField&
169 Foam::DsmcCloud<ParcelType>::internalEBF()
171     return internalE_.boundaryField();
175 template<class ParcelType>
176 inline Foam::volScalarField::GeometricBoundaryField&
177 Foam::DsmcCloud<ParcelType>::iDofBF()
179     return iDof_.boundaryField();
183 template<class ParcelType>
184 inline Foam::volVectorField::GeometricBoundaryField&
185 Foam::DsmcCloud<ParcelType>::momentumBF()
187     return momentum_.boundaryField();
191 template<class ParcelType>
192 inline const Foam::volScalarField&
193 Foam::DsmcCloud<ParcelType>::boundaryT() const
195     return boundaryT_;
199 template<class ParcelType>
200 inline const Foam::volVectorField&
201 Foam::DsmcCloud<ParcelType>::boundaryU() const
203     return boundaryU_;
207 template<class ParcelType>
208 inline const Foam::BinaryCollisionModel<Foam::DsmcCloud<ParcelType> >&
209 Foam::DsmcCloud<ParcelType>::binaryCollision() const
211     return binaryCollisionModel_;
215 template<class ParcelType>
216 inline Foam::BinaryCollisionModel<Foam::DsmcCloud<ParcelType> >&
217 Foam::DsmcCloud<ParcelType>::binaryCollision()
219     return binaryCollisionModel_();
223 template<class ParcelType>
224 inline const Foam::WallInteractionModel<Foam::DsmcCloud<ParcelType> >&
225 Foam::DsmcCloud<ParcelType>::wallInteraction() const
227     return wallInteractionModel_;
231 template<class ParcelType>
232 inline Foam::WallInteractionModel<Foam::DsmcCloud<ParcelType> >&
233 Foam::DsmcCloud<ParcelType>::wallInteraction()
235     return wallInteractionModel_();
239 template<class ParcelType>
240 inline const Foam::InflowBoundaryModel<Foam::DsmcCloud<ParcelType> >&
241 Foam::DsmcCloud<ParcelType>::inflowBoundary() const
243     return inflowBoundaryModel_;
247 template<class ParcelType>
248 inline Foam::InflowBoundaryModel<Foam::DsmcCloud<ParcelType> >&
249 Foam::DsmcCloud<ParcelType>::inflowBoundary()
251     return inflowBoundaryModel_();
255 template<class ParcelType>
256 inline Foam::scalar Foam::DsmcCloud<ParcelType>::massInSystem() const
258     scalar sysMass = 0.0;
260     forAllConstIter(typename DsmcCloud<ParcelType>, *this, iter)
261     {
262         const ParcelType& p = iter();
264         const typename ParcelType::constantProperties& cP = constProps
265         (
266             p.typeId()
267         );
269         sysMass += cP.mass();
270     }
272     return nParticle_*sysMass;
276 template<class ParcelType>
277 inline Foam::vector Foam::DsmcCloud<ParcelType>::linearMomentumOfSystem() const
279     vector linearMomentum(vector::zero);
281     forAllConstIter(typename DsmcCloud<ParcelType>, *this, iter)
282     {
283         const ParcelType& p = iter();
285         const typename ParcelType::constantProperties& cP = constProps
286         (
287             p.typeId()
288         );
290         linearMomentum += cP.mass()*p.U();
291     }
293     return nParticle_*linearMomentum;
297 template<class ParcelType>
298 inline Foam::scalar
299 Foam::DsmcCloud<ParcelType>::linearKineticEnergyOfSystem() const
301     scalar linearKineticEnergy = 0.0;
303     forAllConstIter(typename DsmcCloud<ParcelType>, *this, iter)
304     {
305         const ParcelType& p = iter();
307         const typename ParcelType::constantProperties& cP = constProps
308         (
309             p.typeId()
310         );
312         linearKineticEnergy += 0.5*cP.mass()*(p.U() & p.U());
313     }
315     return nParticle_*linearKineticEnergy;
319 template<class ParcelType>
320 inline Foam::scalar
321 Foam::DsmcCloud<ParcelType>::internalEnergyOfSystem() const
323     scalar internalEnergy = 0.0;
325     forAllConstIter(typename DsmcCloud<ParcelType>, *this, iter)
326     {
327         const ParcelType& p = iter();
329         internalEnergy += p.Ei();
330     }
332     return nParticle_*internalEnergy;
336 template<class ParcelType>
337 inline Foam::scalar Foam::DsmcCloud<ParcelType>::maxwellianAverageSpeed
339     scalar temperature,
340     scalar mass
341 ) const
343     return
344         2.0*sqrt(2.0*physicoChemical::k.value()*temperature/(pi*mass));
348 template<class ParcelType>
349 inline Foam::scalarField Foam::DsmcCloud<ParcelType>::maxwellianAverageSpeed
351     scalarField temperature,
352     scalar mass
353 ) const
355     tmp<scalarField> tfld =
356         2.0*sqrt(2.0*physicoChemical::k.value()*temperature/(pi*mass));
357     return tfld();
361 template<class ParcelType>
362 inline Foam::scalar Foam::DsmcCloud<ParcelType>::maxwellianRMSSpeed
364     scalar temperature,
365     scalar mass
366 ) const
368     return sqrt(3.0*physicoChemical::k.value()*temperature/mass);
372 template<class ParcelType>
373 inline Foam::scalarField Foam::DsmcCloud<ParcelType>::maxwellianRMSSpeed
375     scalarField temperature,
376     scalar mass
377 ) const
379     tmp<scalarField> tfld =
380         sqrt(3.0*physicoChemical::k.value()*temperature/mass);
381     return tfld();
385 template<class ParcelType>
386 inline Foam::scalar
387 Foam::DsmcCloud<ParcelType>::maxwellianMostProbableSpeed
389     scalar temperature,
390     scalar mass
391 ) const
393     return sqrt(2.0*physicoChemical::k.value()*temperature/mass);
397 template<class ParcelType>
398 inline Foam::scalarField
399 Foam::DsmcCloud<ParcelType>::maxwellianMostProbableSpeed
401     scalarField temperature,
402     scalar mass
403 ) const
405     tmp<scalarField> tfld =
406         sqrt(2.0*physicoChemical::k.value()*temperature/mass);
407     return tfld();
411 template<class ParcelType>
412 inline const Foam::volScalarField& Foam::DsmcCloud<ParcelType>::q() const
414     return q_;
418 template<class ParcelType>
419 inline const Foam::volVectorField& Foam::DsmcCloud<ParcelType>::fD() const
421     return fD_;
425 template<class ParcelType>
426 inline const Foam::volScalarField&
427 Foam::DsmcCloud<ParcelType>::rhoN() const
429     return rhoN_;
433 template<class ParcelType>
434 inline const Foam::volScalarField& Foam::DsmcCloud<ParcelType>::rhoM() const
436     return rhoM_;
440 template<class ParcelType>
441 inline const Foam::volScalarField&
442 Foam::DsmcCloud<ParcelType>::dsmcRhoN() const
444     return dsmcRhoN_;
448 template<class ParcelType>
449 inline const Foam::volScalarField&
450 Foam::DsmcCloud<ParcelType>::linearKE() const
452     return linearKE_;
456 template<class ParcelType>
457 inline const Foam::volScalarField&
458 Foam::DsmcCloud<ParcelType>::internalE() const
460     return internalE_;
464 template<class ParcelType>
465 inline const Foam::volScalarField&
466 Foam::DsmcCloud<ParcelType>::iDof() const
468     return iDof_;
472 template<class ParcelType>
473 inline const Foam::volVectorField& Foam::DsmcCloud<ParcelType>::momentum() const
475     return momentum_;
479 template<class ParcelType>
480 inline void Foam::DsmcCloud<ParcelType>::clear()
482     return IDLList<ParcelType>::clear();
486 // ************************************************************************* //