initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / lagrangian / intermediate / parcels / Templates / KinematicParcel / KinematicParcelI.H
blobe0d927153780af56e42ac11af2db03636bd583b0
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-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 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
29 template <class ParcelType>
30 inline Foam::KinematicParcel<ParcelType>::constantProperties::constantProperties
32     const dictionary& parentDict
35     dict_(parentDict.subDict("constantProperties")),
36     rhoMin_(dimensionedScalar(dict_.lookup("rhoMin")).value()),
37     rho0_(dimensionedScalar(dict_.lookup("rho0")).value()),
38     minParticleMass_
39     (
40         dimensionedScalar(dict_.lookup("minParticleMass")).value()
41     )
45 template <class ParcelType>
46 inline Foam::KinematicParcel<ParcelType>::trackData::trackData
48     KinematicCloud<ParcelType>& cloud,
49     const constantProperties& constProps,
50     const interpolation<scalar>& rhoInterp,
51     const interpolation<vector>& UInterp,
52     const interpolation<scalar>& muInterp,
53     const vector& g
56     Particle<ParcelType>::trackData(cloud),
57     cloud_(cloud),
58     constProps_(constProps),
59     rhoInterp_(rhoInterp),
60     UInterp_(UInterp),
61     muInterp_(muInterp),
62     g_(g)
66 template <class ParcelType>
67 inline Foam::KinematicParcel<ParcelType>::KinematicParcel
69     KinematicCloud<ParcelType>& owner,
70     const vector& position,
71     const label cellI
74     Particle<ParcelType>(owner, position, cellI),
75     typeId_(owner.parcelTypeId()),
76     nParticle_(0),
77     d_(0.0),
78     U_(vector::zero),
79     rho_(0.0),
80     tTurb_(0.0),
81     UTurb_(vector::zero),
82     rhoc_(0.0),
83     Uc_(vector::zero),
84     muc_(0.0)
88 template <class ParcelType>
89 inline Foam::KinematicParcel<ParcelType>::KinematicParcel
91     KinematicCloud<ParcelType>& owner,
92     const vector& position,
93     const label cellI,
94     const label typeId,
95     const scalar nParticle0,
96     const scalar d0,
97     const vector& U0,
98     const constantProperties& constProps
101     Particle<ParcelType>(owner, position, cellI),
102     typeId_(typeId),
103     nParticle_(nParticle0),
104     d_(d0),
105     U_(U0),
106     rho_(constProps.rho0()),
107     tTurb_(0.0),
108     UTurb_(vector::zero),
109     rhoc_(0.0),
110     Uc_(vector::zero),
111     muc_(0.0)
115 // * * * * * * * * * constantProperties Member Functions * * * * * * * * * * //
117 template <class ParcelType>
118 inline const Foam::dictionary&
119 Foam::KinematicParcel<ParcelType>::constantProperties::dict() const
121     return dict_;
125 template <class ParcelType>
126 inline Foam::scalar
127 Foam::KinematicParcel<ParcelType>::constantProperties::rhoMin() const
129     return rhoMin_;
133 template <class ParcelType>
134 inline Foam::scalar
135 Foam::KinematicParcel<ParcelType>::constantProperties::rho0() const
137     return rho0_;
141 template <class ParcelType>
142 inline Foam::scalar
143 Foam::KinematicParcel<ParcelType>::constantProperties::minParticleMass() const
145     return minParticleMass_;
149 // * * * * * * * * * * * trackData Member Functions  * * * * * * * * * * * * //
151 template <class ParcelType>
152 inline Foam::KinematicCloud<ParcelType>&
153 Foam::KinematicParcel<ParcelType>::trackData::cloud()
155     return cloud_;
159 template <class ParcelType>
160 inline const typename Foam::KinematicParcel<ParcelType>::constantProperties&
161 Foam::KinematicParcel<ParcelType>::trackData::constProps() const
163     return constProps_;
167 template<class ParcelType>
168 inline const Foam::interpolation<Foam::scalar>&
169 Foam::KinematicParcel<ParcelType>::trackData::rhoInterp() const
171     return rhoInterp_;
175 template <class ParcelType>
176 inline const Foam::interpolation<Foam::vector>&
177 Foam::KinematicParcel<ParcelType>::trackData::UInterp() const
179     return UInterp_;
183 template<class ParcelType>
184 inline const Foam::interpolation<Foam::scalar>&
185 Foam::KinematicParcel<ParcelType>::trackData::muInterp() const
187     return muInterp_;
191 template<class ParcelType>
192 inline const Foam::vector&
193 Foam::KinematicParcel<ParcelType>::trackData::g() const
195     return g_;
199 // * * * * * * * * * * KinematicParcel Member Functions  * * * * * * * * * * //
201 template <class ParcelType>
202 inline Foam::label Foam::KinematicParcel<ParcelType>::typeId() const
204     return typeId_;
208 template <class ParcelType>
209 inline Foam::scalar Foam::KinematicParcel<ParcelType>::nParticle() const
211     return nParticle_;
215 template <class ParcelType>
216 inline Foam::scalar Foam::KinematicParcel<ParcelType>::d() const
218     return d_;
222 template <class ParcelType>
223 inline const Foam::vector& Foam::KinematicParcel<ParcelType>::U() const
225     return U_;
229 template <class ParcelType>
230 inline Foam::scalar Foam::KinematicParcel<ParcelType>::rho() const
232     return rho_;
236 template <class ParcelType>
237 inline Foam::scalar Foam::KinematicParcel<ParcelType>::tTurb() const
239     return tTurb_;
243 template <class ParcelType>
244 inline const Foam::vector& Foam::KinematicParcel<ParcelType>::UTurb() const
246     return UTurb_;
250 template <class ParcelType>
251 inline Foam::label Foam::KinematicParcel<ParcelType>::typeId()
253     return typeId_;
257 template <class ParcelType>
258 inline Foam::scalar& Foam::KinematicParcel<ParcelType>::nParticle()
260     return nParticle_;
264 template <class ParcelType>
265 inline Foam::scalar& Foam::KinematicParcel<ParcelType>::d()
267     return d_;
271 template <class ParcelType>
272 inline Foam::vector& Foam::KinematicParcel<ParcelType>::U()
274     return U_;
278 template <class ParcelType>
279 inline Foam::scalar& Foam::KinematicParcel<ParcelType>::rho()
281     return rho_;
285 template <class ParcelType>
286 inline Foam::scalar& Foam::KinematicParcel<ParcelType>::tTurb()
288     return tTurb_;
292 template <class ParcelType>
293 inline Foam::vector& Foam::KinematicParcel<ParcelType>::UTurb()
295     return UTurb_;
299 template <class ParcelType>
300 inline Foam::scalar Foam::KinematicParcel<ParcelType>::wallImpactDistance
302     const vector&
303 ) const
305     return 0.5*d_;
309 template <class ParcelType>
310 inline Foam::label Foam::KinematicParcel<ParcelType>::faceInterpolation() const
312     // Use volume-based interpolation if dealing with external faces
313     if (this->cloud().internalFace(this->face()))
314     {
315         return this->face();
316     }
317     else
318     {
319         return -1;
320     }
324 template <class ParcelType>
325 inline Foam::scalar Foam::KinematicParcel<ParcelType>::massCell
327     const label cellI
328 ) const
330     return rhoc_*this->cloud().pMesh().cellVolumes()[cellI];
334 template <class ParcelType>
335 inline Foam::scalar Foam::KinematicParcel<ParcelType>::mass() const
337     return rho_*volume();
341 template <class ParcelType>
342 inline Foam::scalar Foam::KinematicParcel<ParcelType>::volume() const
344     return volume(d_);
348 template <class ParcelType>
349 inline Foam::scalar
350 Foam::KinematicParcel<ParcelType>::volume(const scalar d) const
352     return mathematicalConstant::pi/6.0*pow3(d);
356 template <class ParcelType>
357 inline Foam::scalar Foam::KinematicParcel<ParcelType>::areaP() const
359     return areaP(d_);
363 template <class ParcelType>
364 inline Foam::scalar
365 Foam::KinematicParcel<ParcelType>::areaP(const scalar d) const
367     return 0.25*areaS(d);
371 template <class ParcelType>
372 inline Foam::scalar Foam::KinematicParcel<ParcelType>::areaS() const
374     return areaS(d_);
378 template <class ParcelType>
379 inline Foam::scalar
380 Foam::KinematicParcel<ParcelType>::areaS(const scalar d) const
382     return mathematicalConstant::pi*d*d;
386 // ************************************************************************* //