fixed writing out entries in advective bc
[OpenFOAM-1.6-ext.git] / src / lagrangian / intermediate / clouds / Templates / KinematicCloud / KinematicCloudI.H
blob9a2522e801f203fb92e77f8ac80d008eed497a9b
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 "fvmSup.H"
29 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
31 template<class ParcelType>
32 inline Foam::label Foam::KinematicCloud<ParcelType>::parcelTypeId() const
34     return parcelTypeId_;
38 template<class ParcelType>
39 inline const Foam::fvMesh& Foam::KinematicCloud<ParcelType>::mesh() const
41     return mesh_;
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
57     return constProps_;
61 template<class ParcelType>
62 inline const Foam::Switch Foam::KinematicCloud<ParcelType>::active() const
64     return active_;
68 template<class ParcelType>
69 inline const Foam::Switch Foam::KinematicCloud<ParcelType>::coupled() const
71     return coupled_;
75 template <class ParcelType>
76 inline const Foam::Switch
77 Foam::KinematicCloud<ParcelType>::cellValueSourceCorrection() const
79     return cellValueSourceCorrection_;
83 template<class ParcelType>
84 inline const Foam::volScalarField&
85 Foam::KinematicCloud<ParcelType>::rho() const
87     return rho_;
91 template<class ParcelType>
92 inline const Foam::volVectorField& Foam::KinematicCloud<ParcelType>::U() const
94     return U_;
98 template<class ParcelType>
99 inline const Foam::volScalarField& Foam::KinematicCloud<ParcelType>::mu() const
101     return mu_;
105 template<class ParcelType>
106 inline const Foam::dimensionedVector&
107 Foam::KinematicCloud<ParcelType>::g() const
109     return g_;
113 template<class ParcelType>
114 inline const Foam::particleForces&
115 Foam::KinematicCloud<ParcelType>::forces() const
117     return forces_;
121 template<class ParcelType>
122 inline const Foam::dictionary&
123 Foam::KinematicCloud<ParcelType>::interpolationSchemes() const
125     return interpolationSchemes_;
129 template<class ParcelType>
130 inline const Foam::DispersionModel<Foam::KinematicCloud<ParcelType> >&
131 Foam::KinematicCloud<ParcelType>::dispersion() const
133     return dispersionModel_;
137 template<class ParcelType>
138 inline Foam::DispersionModel<Foam::KinematicCloud<ParcelType> >&
139 Foam::KinematicCloud<ParcelType>::dispersion()
141     return dispersionModel_();
145 template<class ParcelType>
146 inline const Foam::DragModel<Foam::KinematicCloud<ParcelType> >&
147 Foam::KinematicCloud<ParcelType>::drag() const
149     return dragModel_;
153 template<class ParcelType>
154 inline const Foam::InjectionModel<Foam::KinematicCloud<ParcelType> >&
155 Foam::KinematicCloud<ParcelType>::injection() const
157     return injectionModel_;
161 template<class ParcelType>
162 inline const Foam::PatchInteractionModel<Foam::KinematicCloud<ParcelType> >&
163 Foam::KinematicCloud<ParcelType>::patchInteraction() const
165     return patchInteractionModel_;
169 template<class ParcelType>
170 inline Foam::InjectionModel<Foam::KinematicCloud<ParcelType> >&
171 Foam::KinematicCloud<ParcelType>::injection()
173     return injectionModel_();
177 template<class ParcelType>
178 inline Foam::PostProcessingModel<Foam::KinematicCloud<ParcelType> >&
179 Foam::KinematicCloud<ParcelType>::postProcessing()
181     return postProcessingModel_();
185 template<class ParcelType>
186 inline const Foam::vectorIntegrationScheme&
187 Foam::KinematicCloud<ParcelType>::UIntegrator() const
189     return UIntegrator_;
193 template<class ParcelType>
194 inline Foam::scalar Foam::KinematicCloud<ParcelType>::massInSystem() const
196     scalar sysMass = 0.0;
197     forAllConstIter(typename KinematicCloud<ParcelType>, *this, iter)
198     {
199          const ParcelType& p = iter();
200          sysMass += p.mass()*p.nParticle();
201     }
203     return sysMass;
207 template<class ParcelType>
208 inline Foam::Random& Foam::KinematicCloud<ParcelType>::rndGen()
210     return rndGen_;
214 template<class ParcelType>
215 inline Foam::DimensionedField<Foam::vector, Foam::volMesh>&
216 Foam::KinematicCloud<ParcelType>::UTrans()
218     return UTrans_;
222 template<class ParcelType>
223 inline Foam::tmp<Foam::DimensionedField<Foam::vector, Foam::volMesh> >
224 Foam::KinematicCloud<ParcelType>::SU() const
226     tmp<DimensionedField<vector, volMesh> > tSU
227     (
228         new DimensionedField<vector, volMesh>
229         (
230             IOobject
231             (
232                 this->name() + "SU",
233                 this->db().time().timeName(),
234                 this->mesh(),
235                 IOobject::NO_READ,
236                 IOobject::AUTO_WRITE
237             ),
238             this->mesh(),
239             dimensionedVector
240             (
241                  "zero",
242                  dimDensity*dimVelocity/dimTime,
243                  vector::zero
244             )
245         )
246     );
248     vectorField& SU = tSU().field();
249     SU = UTrans_/(mesh_.V()*this->db().time().deltaT());
251     return tSU;
255 template<class ParcelType>
256 inline const Foam::tmp<Foam::volScalarField>
257 Foam::KinematicCloud<ParcelType>::theta() const
259     tmp<volScalarField> ttheta
260     (
261         new volScalarField
262         (
263             IOobject
264             (
265                 this->name() + "Theta",
266                 this->db().time().timeName(),
267                 this->db(),
268                 IOobject::NO_READ,
269                 IOobject::NO_WRITE,
270                 false
271             ),
272             mesh_,
273             dimensionedScalar("zero", dimless, 0.0)
274         )
275     );
277     scalarField& theta = ttheta().internalField();
278     forAllConstIter(typename KinematicCloud<ParcelType>, *this, iter)
279     {
280         const ParcelType& p = iter();
281         const label cellI = p.cell();
283         theta[cellI] += p.nParticle()*p.volume();
284     }
286     theta /= mesh().V();
288     return ttheta;
292 template<class ParcelType>
293 inline const Foam::tmp<Foam::volScalarField>
294 Foam::KinematicCloud<ParcelType>::alpha() const
296     tmp<volScalarField> talpha
297     (
298         new volScalarField
299         (
300             IOobject
301             (
302                 this->name() + "Alpha",
303                 this->db().time().timeName(),
304                 this->db(),
305                 IOobject::NO_READ,
306                 IOobject::NO_WRITE,
307                 false
308             ),
309             mesh_,
310             dimensionedScalar("zero", dimless, 0.0)
311         )
312     );
314     scalarField& alpha = talpha().internalField();
315     forAllConstIter(typename KinematicCloud<ParcelType>, *this, iter)
316     {
317         const ParcelType& p = iter();
318         const label cellI = p.cell();
320         alpha[cellI] += p.nParticle()*p.mass();
321     }
323     alpha /= (mesh().V()*rho_);
325     return talpha;
329 template<class ParcelType>
330 inline const Foam::tmp<Foam::volScalarField>
331 Foam::KinematicCloud<ParcelType>::rhoEff() const
333     tmp<volScalarField> trhoEff
334     (
335         new volScalarField
336         (
337             IOobject
338             (
339                 this->name() + "RhoEff",
340                 this->db().time().timeName(),
341                 this->db(),
342                 IOobject::NO_READ,
343                 IOobject::NO_WRITE,
344                 false
345             ),
346             mesh_,
347             dimensionedScalar("zero", dimDensity, 0.0)
348         )
349     );
351     scalarField& rhoEff = trhoEff().internalField();
352     forAllConstIter(typename KinematicCloud<ParcelType>, *this, iter)
353     {
354         const ParcelType& p = iter();
355         const label cellI = p.cell();
357         rhoEff[cellI] += p.nParticle()*p.mass();
358     }
360     rhoEff /= mesh().V();
362     return trhoEff;
366 // ************************************************************************* //