1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-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 \*---------------------------------------------------------------------------*/
29 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
31 template<class ParcelType>
32 inline Foam::label Foam::KinematicCloud<ParcelType>::parcelTypeId() const
38 template<class ParcelType>
39 inline const Foam::fvMesh& Foam::KinematicCloud<ParcelType>::mesh() const
45 template<class ParcelType>
46 inline const Foam::IOdictionary&
47 Foam::KinematicCloud<ParcelType>::particleProperties() const
49 return particleProperties_;
53 template<class ParcelType>
54 inline const typename ParcelType::constantProperties&
55 Foam::KinematicCloud<ParcelType>::constProps() const
61 template<class ParcelType>
62 inline const Foam::Switch Foam::KinematicCloud<ParcelType>::coupled() const
68 template <class ParcelType>
69 inline const Foam::Switch
70 Foam::KinematicCloud<ParcelType>::cellValueSourceCorrection() const
72 return cellValueSourceCorrection_;
76 template<class ParcelType>
77 inline const Foam::volScalarField&
78 Foam::KinematicCloud<ParcelType>::rho() const
84 template<class ParcelType>
85 inline const Foam::volVectorField& Foam::KinematicCloud<ParcelType>::U() const
91 template<class ParcelType>
92 inline const Foam::volScalarField& Foam::KinematicCloud<ParcelType>::mu() const
98 template<class ParcelType>
99 inline const Foam::dimensionedVector&
100 Foam::KinematicCloud<ParcelType>::g() const
106 template<class ParcelType>
107 inline const Foam::particleForces&
108 Foam::KinematicCloud<ParcelType>::forces() const
114 template<class ParcelType>
115 inline const Foam::dictionary&
116 Foam::KinematicCloud<ParcelType>::interpolationSchemes() const
118 return interpolationSchemes_;
122 template<class ParcelType>
123 inline const Foam::DispersionModel<Foam::KinematicCloud<ParcelType> >&
124 Foam::KinematicCloud<ParcelType>::dispersion() const
126 return dispersionModel_;
130 template<class ParcelType>
131 inline Foam::DispersionModel<Foam::KinematicCloud<ParcelType> >&
132 Foam::KinematicCloud<ParcelType>::dispersion()
134 return dispersionModel_();
138 template<class ParcelType>
139 inline const Foam::DragModel<Foam::KinematicCloud<ParcelType> >&
140 Foam::KinematicCloud<ParcelType>::drag() const
146 template<class ParcelType>
147 inline const Foam::InjectionModel<Foam::KinematicCloud<ParcelType> >&
148 Foam::KinematicCloud<ParcelType>::injection() const
150 return injectionModel_;
154 template<class ParcelType>
155 inline const Foam::PatchInteractionModel<Foam::KinematicCloud<ParcelType> >&
156 Foam::KinematicCloud<ParcelType>::patchInteraction() const
158 return patchInteractionModel_;
162 template<class ParcelType>
163 inline Foam::InjectionModel<Foam::KinematicCloud<ParcelType> >&
164 Foam::KinematicCloud<ParcelType>::injection()
166 return injectionModel_();
170 template<class ParcelType>
171 inline Foam::PostProcessingModel<Foam::KinematicCloud<ParcelType> >&
172 Foam::KinematicCloud<ParcelType>::postProcessing()
174 return postProcessingModel_();
178 template<class ParcelType>
179 inline const Foam::vectorIntegrationScheme&
180 Foam::KinematicCloud<ParcelType>::UIntegrator() const
186 template<class ParcelType>
187 inline Foam::scalar Foam::KinematicCloud<ParcelType>::massInSystem() const
189 scalar sysMass = 0.0;
190 forAllConstIter(typename KinematicCloud<ParcelType>, *this, iter)
192 const ParcelType& p = iter();
193 sysMass += p.mass()*p.nParticle();
200 template<class ParcelType>
201 inline Foam::Random& Foam::KinematicCloud<ParcelType>::rndGen()
207 template<class ParcelType>
208 inline Foam::DimensionedField<Foam::vector, Foam::volMesh>&
209 Foam::KinematicCloud<ParcelType>::UTrans()
215 template<class ParcelType>
216 inline Foam::tmp<Foam::DimensionedField<Foam::vector, Foam::volMesh> >
217 Foam::KinematicCloud<ParcelType>::SU() const
219 tmp<DimensionedField<vector, volMesh> > tSU
221 new DimensionedField<vector, volMesh>
226 this->db().time().timeName(),
235 dimDensity*dimVelocity/dimTime,
241 vectorField& SU = tSU().field();
242 SU = UTrans_/(mesh_.V()*this->db().time().deltaT());
248 template<class ParcelType>
249 inline const Foam::tmp<Foam::volScalarField>
250 Foam::KinematicCloud<ParcelType>::theta() const
252 tmp<volScalarField> ttheta
258 this->name() + "Theta",
259 this->db().time().timeName(),
266 dimensionedScalar("zero", dimless, 0.0)
270 scalarField& theta = ttheta().internalField();
271 forAllConstIter(typename KinematicCloud<ParcelType>, *this, iter)
273 const ParcelType& p = iter();
274 const label cellI = p.cell();
276 theta[cellI] += p.nParticle()*p.volume();
285 template<class ParcelType>
286 inline const Foam::tmp<Foam::volScalarField>
287 Foam::KinematicCloud<ParcelType>::alpha() const
289 tmp<volScalarField> talpha
295 this->name() + "Alpha",
296 this->db().time().timeName(),
303 dimensionedScalar("zero", dimless, 0.0)
307 scalarField& alpha = talpha().internalField();
308 forAllConstIter(typename KinematicCloud<ParcelType>, *this, iter)
310 const ParcelType& p = iter();
311 const label cellI = p.cell();
313 alpha[cellI] += p.nParticle()*p.mass();
316 alpha /= (mesh().V()*rho_);
322 template<class ParcelType>
323 inline const Foam::tmp<Foam::volScalarField>
324 Foam::KinematicCloud<ParcelType>::rhoEff() const
326 tmp<volScalarField> trhoEff
332 this->name() + "RhoEff",
333 this->db().time().timeName(),
340 dimensionedScalar("zero", dimDensity, 0.0)
344 scalarField& rhoEff = trhoEff().internalField();
345 forAllConstIter(typename KinematicCloud<ParcelType>, *this, iter)
347 const ParcelType& p = iter();
348 const label cellI = p.cell();
350 rhoEff[cellI] += p.nParticle()*p.mass();
353 rhoEff /= mesh().V();
359 // ************************************************************************* //