Added measurement of surface values of velocity and temperature to allow slip
[OpenFOAM-1.6.x.git] / src / lagrangian / dsmc / clouds / Templates / DsmcCloud / DsmcCloudI.H
blobafb7f642f20a3d67e4d7845fd56dd7fdf0c3b577
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2009-2009 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 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
19     for more details.
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
32     return cloudName_;
36 template<class ParcelType>
37 inline const Foam::fvMesh& Foam::DsmcCloud<ParcelType>::mesh() const
39     return mesh_;
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
55     return typeIdList_;
59 template<class ParcelType>
60 inline Foam::scalar Foam::DsmcCloud<ParcelType>::nParticle() const
62     return nParticle_;
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()
77     return 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
93     return constProps_;
97 template<class ParcelType>
98 inline const typename ParcelType::constantProperties&
99 Foam::DsmcCloud<ParcelType>::constProps
101     label typeId
102 ) const
104     if (typeId < 0 || typeId >= constProps_.size())
105     {
106         FatalErrorIn("Foam::DsmcCloud<ParcelType>::constProps(label typeId)")
107             << "constantProperties for requested typeId index "
108             << typeId << " do not exist" << nl
109             << abort(FatalError);
110     }
112     return constProps_[typeId];
116 template<class ParcelType>
117 inline Foam::Random& Foam::DsmcCloud<ParcelType>::rndGen()
119     return rndGen_;
123 template<class ParcelType>
124 inline Foam::volScalarField::GeometricBoundaryField&
125 Foam::DsmcCloud<ParcelType>::qBF()
127     return q_.boundaryField();
131 template<class ParcelType>
132 inline Foam::volVectorField::GeometricBoundaryField&
133 Foam::DsmcCloud<ParcelType>::fDBF()
135     return fD_.boundaryField();
139 template<class ParcelType>
140 inline Foam::volScalarField::GeometricBoundaryField&
141 Foam::DsmcCloud<ParcelType>::rhoNBF()
143     return rhoN_.boundaryField();
147 template<class ParcelType>
148 inline Foam::volScalarField::GeometricBoundaryField&
149 Foam::DsmcCloud<ParcelType>::rhoMBF()
151     return rhoM_.boundaryField();
155 template<class ParcelType>
156 inline Foam::volScalarField::GeometricBoundaryField&
157 Foam::DsmcCloud<ParcelType>::linearKEBF()
159     return linearKE_.boundaryField();
163 template<class ParcelType>
164 inline Foam::volScalarField::GeometricBoundaryField&
165 Foam::DsmcCloud<ParcelType>::internalEBF()
167     return internalE_.boundaryField();
171 template<class ParcelType>
172 inline Foam::volScalarField::GeometricBoundaryField&
173 Foam::DsmcCloud<ParcelType>::iDofBF()
175     return iDof_.boundaryField();
179 template<class ParcelType>
180 inline Foam::volVectorField::GeometricBoundaryField&
181 Foam::DsmcCloud<ParcelType>::momentumBF()
183     return momentum_.boundaryField();
187 template<class ParcelType>
188 inline const Foam::volScalarField&
189 Foam::DsmcCloud<ParcelType>::boundaryT() const
191     return boundaryT_;
195 template<class ParcelType>
196 inline const Foam::volVectorField&
197 Foam::DsmcCloud<ParcelType>::boundaryU() const
199     return boundaryU_;
203 template<class ParcelType>
204 inline const Foam::BinaryCollisionModel<Foam::DsmcCloud<ParcelType> >&
205 Foam::DsmcCloud<ParcelType>::binaryCollision() const
207     return binaryCollisionModel_;
211 template<class ParcelType>
212 inline Foam::BinaryCollisionModel<Foam::DsmcCloud<ParcelType> >&
213 Foam::DsmcCloud<ParcelType>::binaryCollision()
215     return binaryCollisionModel_();
219 template<class ParcelType>
220 inline const Foam::WallInteractionModel<Foam::DsmcCloud<ParcelType> >&
221 Foam::DsmcCloud<ParcelType>::wallInteraction() const
223     return wallInteractionModel_;
227 template<class ParcelType>
228 inline Foam::WallInteractionModel<Foam::DsmcCloud<ParcelType> >&
229 Foam::DsmcCloud<ParcelType>::wallInteraction()
231     return wallInteractionModel_();
235 template<class ParcelType>
236 inline const Foam::InflowBoundaryModel<Foam::DsmcCloud<ParcelType> >&
237 Foam::DsmcCloud<ParcelType>::inflowBoundary() const
239     return inflowBoundaryModel_;
243 template<class ParcelType>
244 inline Foam::InflowBoundaryModel<Foam::DsmcCloud<ParcelType> >&
245 Foam::DsmcCloud<ParcelType>::inflowBoundary()
247     return inflowBoundaryModel_();
251 template<class ParcelType>
252 inline Foam::scalar Foam::DsmcCloud<ParcelType>::massInSystem() const
254     scalar sysMass = 0.0;
256     forAllConstIter(typename DsmcCloud<ParcelType>, *this, iter)
257     {
258         const ParcelType& p = iter();
260         const typename ParcelType::constantProperties& cP = constProps
261         (
262             p.typeId()
263         );
265         sysMass += cP.mass();
266     }
268     return nParticle_*sysMass;
272 template<class ParcelType>
273 inline Foam::vector Foam::DsmcCloud<ParcelType>::linearMomentumOfSystem() const
275     vector linearMomentum(vector::zero);
277     forAllConstIter(typename DsmcCloud<ParcelType>, *this, iter)
278     {
279         const ParcelType& p = iter();
281         const typename ParcelType::constantProperties& cP = constProps
282         (
283             p.typeId()
284         );
286         linearMomentum += cP.mass()*p.U();
287     }
289     return nParticle_*linearMomentum;
293 template<class ParcelType>
294 inline Foam::scalar
295 Foam::DsmcCloud<ParcelType>::linearKineticEnergyOfSystem() const
297     scalar linearKineticEnergy = 0.0;
299     forAllConstIter(typename DsmcCloud<ParcelType>, *this, iter)
300     {
301         const ParcelType& p = iter();
303         const typename ParcelType::constantProperties& cP = constProps
304         (
305             p.typeId()
306         );
308         linearKineticEnergy += 0.5*cP.mass()*(p.U() & p.U());
309     }
311     return nParticle_*linearKineticEnergy;
315 template<class ParcelType>
316 inline Foam::scalar
317 Foam::DsmcCloud<ParcelType>::internalEnergyOfSystem() const
319     scalar internalEnergy = 0.0;
321     forAllConstIter(typename DsmcCloud<ParcelType>, *this, iter)
322     {
323         const ParcelType& p = iter();
325         internalEnergy += p.Ei();
326     }
328     return nParticle_*internalEnergy;
332 template<class ParcelType>
333 inline Foam::scalar Foam::DsmcCloud<ParcelType>::maxwellianAverageSpeed
335     scalar temperature,
336     scalar mass
337 ) const
339     return
340         2.0*sqrt(2.0*kb*temperature/(mathematicalConstant::pi*mass));
344 template<class ParcelType>
345 inline Foam::scalarField Foam::DsmcCloud<ParcelType>::maxwellianAverageSpeed
347     scalarField temperature,
348     scalar mass
349 ) const
351     return
352         2.0*sqrt(2.0*kb*temperature/(mathematicalConstant::pi*mass));
356 template<class ParcelType>
357 inline Foam::scalar Foam::DsmcCloud<ParcelType>::maxwellianRMSSpeed
359     scalar temperature,
360     scalar mass
361 ) const
363     return sqrt(3.0*kb*temperature/mass);
367 template<class ParcelType>
368 inline Foam::scalarField Foam::DsmcCloud<ParcelType>::maxwellianRMSSpeed
370     scalarField temperature,
371     scalar mass
372 ) const
374     return sqrt(3.0*kb*temperature/mass);
378 template<class ParcelType>
379 inline Foam::scalar
380 Foam::DsmcCloud<ParcelType>::maxwellianMostProbableSpeed
382     scalar temperature,
383     scalar mass
384 ) const
386     return sqrt(2.0*kb*temperature/mass);
390 template<class ParcelType>
391 inline Foam::scalarField
392 Foam::DsmcCloud<ParcelType>::maxwellianMostProbableSpeed
394     scalarField temperature,
395     scalar mass
396 ) const
398     return sqrt(2.0*kb*temperature/mass);
402 template<class ParcelType>
403 inline const Foam::volScalarField& Foam::DsmcCloud<ParcelType>::q() const
405     return q_;
409 template<class ParcelType>
410 inline const Foam::volVectorField& Foam::DsmcCloud<ParcelType>::fD() const
412     return fD_;
416 template<class ParcelType>
417 inline const Foam::volScalarField&
418 Foam::DsmcCloud<ParcelType>::rhoN() const
420     return rhoN_;
424 template<class ParcelType>
425 inline const Foam::volScalarField& Foam::DsmcCloud<ParcelType>::rhoM() const
427     return rhoM_;
431 template<class ParcelType>
432 inline const Foam::volScalarField&
433 Foam::DsmcCloud<ParcelType>::dsmcRhoN() const
435     return dsmcRhoN_;
439 template<class ParcelType>
440 inline const Foam::volScalarField&
441 Foam::DsmcCloud<ParcelType>::linearKE() const
443     return linearKE_;
447 template<class ParcelType>
448 inline const Foam::volScalarField&
449 Foam::DsmcCloud<ParcelType>::internalE() const
451     return internalE_;
455 template<class ParcelType>
456 inline const Foam::volScalarField&
457 Foam::DsmcCloud<ParcelType>::iDof() const
459     return iDof_;
463 template<class ParcelType>
464 inline const Foam::volVectorField& Foam::DsmcCloud<ParcelType>::momentum() const
466     return momentum_;
470 template<class ParcelType>
471 inline void Foam::DsmcCloud<ParcelType>::clear()
473     return IDLList<ParcelType>::clear();
477 // ************************************************************************* //