fixed writing out entries in advective bc
[OpenFOAM-1.6-ext.git] / src / lagrangian / intermediate / parcels / Templates / KinematicParcel / KinematicParcelI.H
blob1afaf39c879a208629a4205d30c2ea08ed14d330
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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 #include "mathematicalConstants.H"
29 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
31 template <class ParcelType>
32 inline Foam::KinematicParcel<ParcelType>::constantProperties::constantProperties
34     const dictionary& parentDict
37     dict_(parentDict.subDict("constantProperties")),
38     rhoMin_(dimensionedScalar(dict_.lookup("rhoMin")).value()),
39     rho0_(dimensionedScalar(dict_.lookup("rho0")).value()),
40     minParticleMass_
41     (
42         dimensionedScalar(dict_.lookup("minParticleMass")).value()
43     )
47 template <class ParcelType>
48 inline Foam::KinematicParcel<ParcelType>::trackData::trackData
50     KinematicCloud<ParcelType>& cloud,
51     const constantProperties& constProps,
52     const interpolation<scalar>& rhoInterp,
53     const interpolation<vector>& UInterp,
54     const interpolation<scalar>& muInterp,
55     const vector& g
58     Particle<ParcelType>::trackData(cloud),
59     cloud_(cloud),
60     constProps_(constProps),
61     rhoInterp_(rhoInterp),
62     UInterp_(UInterp),
63     muInterp_(muInterp),
64     g_(g)
68 template <class ParcelType>
69 inline Foam::KinematicParcel<ParcelType>::KinematicParcel
71     KinematicCloud<ParcelType>& owner,
72     const vector& position,
73     const label cellI
76     Particle<ParcelType>(owner, position, cellI),
77     typeId_(owner.parcelTypeId()),
78     nParticle_(0),
79     d_(0.0),
80     U_(vector::zero),
81     rho_(0.0),
82     tTurb_(0.0),
83     UTurb_(vector::zero),
84     rhoc_(0.0),
85     Uc_(vector::zero),
86     muc_(0.0)
90 template <class ParcelType>
91 inline Foam::KinematicParcel<ParcelType>::KinematicParcel
93     KinematicCloud<ParcelType>& owner,
94     const vector& position,
95     const label cellI,
96     const label typeId,
97     const scalar nParticle0,
98     const scalar d0,
99     const vector& U0,
100     const constantProperties& constProps
103     Particle<ParcelType>(owner, position, cellI),
104     typeId_(typeId),
105     nParticle_(nParticle0),
106     d_(d0),
107     U_(U0),
108     rho_(constProps.rho0()),
109     tTurb_(0.0),
110     UTurb_(vector::zero),
111     rhoc_(0.0),
112     Uc_(vector::zero),
113     muc_(0.0)
117 // * * * * * * * * * constantProperties Member Functions * * * * * * * * * * //
119 template <class ParcelType>
120 inline const Foam::dictionary&
121 Foam::KinematicParcel<ParcelType>::constantProperties::dict() const
123     return dict_;
127 template <class ParcelType>
128 inline Foam::scalar
129 Foam::KinematicParcel<ParcelType>::constantProperties::rhoMin() const
131     return rhoMin_;
135 template <class ParcelType>
136 inline Foam::scalar
137 Foam::KinematicParcel<ParcelType>::constantProperties::rho0() const
139     return rho0_;
143 template <class ParcelType>
144 inline Foam::scalar
145 Foam::KinematicParcel<ParcelType>::constantProperties::minParticleMass() const
147     return minParticleMass_;
151 // * * * * * * * * * * * trackData Member Functions  * * * * * * * * * * * * //
153 template <class ParcelType>
154 inline Foam::KinematicCloud<ParcelType>&
155 Foam::KinematicParcel<ParcelType>::trackData::cloud()
157     return cloud_;
161 template <class ParcelType>
162 inline const typename Foam::KinematicParcel<ParcelType>::constantProperties&
163 Foam::KinematicParcel<ParcelType>::trackData::constProps() const
165     return constProps_;
169 template<class ParcelType>
170 inline const Foam::interpolation<Foam::scalar>&
171 Foam::KinematicParcel<ParcelType>::trackData::rhoInterp() const
173     return rhoInterp_;
177 template <class ParcelType>
178 inline const Foam::interpolation<Foam::vector>&
179 Foam::KinematicParcel<ParcelType>::trackData::UInterp() const
181     return UInterp_;
185 template<class ParcelType>
186 inline const Foam::interpolation<Foam::scalar>&
187 Foam::KinematicParcel<ParcelType>::trackData::muInterp() const
189     return muInterp_;
193 template<class ParcelType>
194 inline const Foam::vector&
195 Foam::KinematicParcel<ParcelType>::trackData::g() const
197     return g_;
201 // * * * * * * * * * * KinematicParcel Member Functions  * * * * * * * * * * //
203 template <class ParcelType>
204 inline Foam::label Foam::KinematicParcel<ParcelType>::typeId() const
206     return typeId_;
210 template <class ParcelType>
211 inline Foam::scalar Foam::KinematicParcel<ParcelType>::nParticle() const
213     return nParticle_;
217 template <class ParcelType>
218 inline Foam::scalar Foam::KinematicParcel<ParcelType>::d() const
220     return d_;
224 template <class ParcelType>
225 inline const Foam::vector& Foam::KinematicParcel<ParcelType>::U() const
227     return U_;
231 template <class ParcelType>
232 inline Foam::scalar Foam::KinematicParcel<ParcelType>::rho() const
234     return rho_;
238 template <class ParcelType>
239 inline Foam::scalar Foam::KinematicParcel<ParcelType>::tTurb() const
241     return tTurb_;
245 template <class ParcelType>
246 inline const Foam::vector& Foam::KinematicParcel<ParcelType>::UTurb() const
248     return UTurb_;
252 template <class ParcelType>
253 inline Foam::label Foam::KinematicParcel<ParcelType>::typeId()
255     return typeId_;
259 template <class ParcelType>
260 inline Foam::scalar& Foam::KinematicParcel<ParcelType>::nParticle()
262     return nParticle_;
266 template <class ParcelType>
267 inline Foam::scalar& Foam::KinematicParcel<ParcelType>::d()
269     return d_;
273 template <class ParcelType>
274 inline Foam::vector& Foam::KinematicParcel<ParcelType>::U()
276     return U_;
280 template <class ParcelType>
281 inline Foam::scalar& Foam::KinematicParcel<ParcelType>::rho()
283     return rho_;
287 template <class ParcelType>
288 inline Foam::scalar& Foam::KinematicParcel<ParcelType>::tTurb()
290     return tTurb_;
294 template <class ParcelType>
295 inline Foam::vector& Foam::KinematicParcel<ParcelType>::UTurb()
297     return UTurb_;
301 template <class ParcelType>
302 inline Foam::scalar Foam::KinematicParcel<ParcelType>::wallImpactDistance
304     const vector&
305 ) const
307     return 0.5*d_;
311 template <class ParcelType>
312 inline Foam::label Foam::KinematicParcel<ParcelType>::faceInterpolation() const
314     // Use volume-based interpolation if dealing with external faces
315     if (this->cloud().internalFace(this->face()))
316     {
317         return this->face();
318     }
319     else
320     {
321         return -1;
322     }
326 template <class ParcelType>
327 inline Foam::scalar Foam::KinematicParcel<ParcelType>::massCell
329     const label cellI
330 ) const
332     return rhoc_*this->cloud().pMesh().cellVolumes()[cellI];
336 template <class ParcelType>
337 inline Foam::scalar Foam::KinematicParcel<ParcelType>::mass() const
339     return rho_*volume();
343 template <class ParcelType>
344 inline Foam::scalar Foam::KinematicParcel<ParcelType>::volume() const
346     return volume(d_);
350 template <class ParcelType>
351 inline Foam::scalar
352 Foam::KinematicParcel<ParcelType>::volume(const scalar d) const
354     return mathematicalConstant::pi/6.0*pow3(d);
358 template <class ParcelType>
359 inline Foam::scalar Foam::KinematicParcel<ParcelType>::areaP() const
361     return areaP(d_);
365 template <class ParcelType>
366 inline Foam::scalar
367 Foam::KinematicParcel<ParcelType>::areaP(const scalar d) const
369     return 0.25*areaS(d);
373 template <class ParcelType>
374 inline Foam::scalar Foam::KinematicParcel<ParcelType>::areaS() const
376     return areaS(d_);
380 template <class ParcelType>
381 inline Foam::scalar
382 Foam::KinematicParcel<ParcelType>::areaS(const scalar d) const
384     return mathematicalConstant::pi*d*d;
388 template <class ParcelType>
389 inline Foam::scalar
390 Foam::KinematicParcel<ParcelType>::Re
392     const vector& U,
393     const scalar d,
394     const scalar rhoc,
395     const scalar muc
396 ) const
398     return rhoc*mag(U - Uc_)*d/muc;
402 // ************************************************************************* //