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 \*---------------------------------------------------------------------------*/
27 #include "KinematicCloud.H"
28 #include "IntegrationScheme.H"
29 #include "interpolation.H"
31 #include "DispersionModel.H"
32 #include "DragModel.H"
33 #include "InjectionModel.H"
34 #include "PatchInteractionModel.H"
35 #include "PostProcessingModel.H"
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 template<class ParcelType>
40 Foam::KinematicCloud<ParcelType>::KinematicCloud
42 const word& cloudName,
43 const volScalarField& rho,
44 const volVectorField& U,
45 const volScalarField& mu,
46 const dimensionedVector& g
49 Cloud<ParcelType>(rho.mesh(), cloudName, false),
56 cloudName + "Properties",
57 rho.mesh().time().constant(),
63 constProps_(particleProperties_),
64 parcelTypeId_(readLabel(particleProperties_.lookup("parcelTypeId"))),
65 coupled_(particleProperties_.lookup("coupled")),
66 cellValueSourceCorrection_
68 particleProperties_.lookup("cellValueSourceCorrection")
75 forces_(mesh_, particleProperties_, g_.value()),
76 interpolationSchemes_(particleProperties_.subDict("interpolationSchemes")),
79 DispersionModel<KinematicCloud<ParcelType> >::New
87 DragModel<KinematicCloud<ParcelType> >::New
95 InjectionModel<KinematicCloud<ParcelType> >::New
101 patchInteractionModel_
103 PatchInteractionModel<KinematicCloud<ParcelType> >::New
111 PostProcessingModel<KinematicCloud<ParcelType> >::New
113 this->particleProperties_,
119 vectorIntegrationScheme::New
122 particleProperties_.subDict("integrationSchemes")
129 this->name() + "UTrans",
130 this->db().time().timeName(),
137 dimensionedVector("zero", dimMass*dimVelocity, vector::zero)
142 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
144 template<class ParcelType>
145 Foam::KinematicCloud<ParcelType>::~KinematicCloud()
149 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
151 template<class ParcelType>
152 void Foam::KinematicCloud<ParcelType>::checkParcelProperties
155 const scalar lagrangianDt,
156 const bool fullyDescribed
161 parcel.rho() = constProps_.rho0();
164 scalar carrierDt = this->db().time().deltaT().value();
165 parcel.stepFraction() = (carrierDt - lagrangianDt)/carrierDt;
169 template<class ParcelType>
170 void Foam::KinematicCloud<ParcelType>::resetSourceTerms()
172 UTrans_.field() = vector::zero;
176 template<class ParcelType>
177 void Foam::KinematicCloud<ParcelType>::preEvolve()
179 this->dispersion().cacheFields(true);
180 forces_.cacheFields(true);
184 template<class ParcelType>
185 void Foam::KinematicCloud<ParcelType>::postEvolve()
189 this->writePositions();
192 this->dispersion().cacheFields(false);
193 forces_.cacheFields(false);
195 this->postProcessing().post();
199 template<class ParcelType>
200 void Foam::KinematicCloud<ParcelType>::evolve()
204 autoPtr<interpolation<scalar> > rhoInterpolator =
205 interpolation<scalar>::New
207 interpolationSchemes_,
211 autoPtr<interpolation<vector> > UInterpolator =
212 interpolation<vector>::New
214 interpolationSchemes_,
218 autoPtr<interpolation<scalar> > muInterpolator =
219 interpolation<scalar>::New
221 interpolationSchemes_,
225 typename ParcelType::trackData td
235 this->injection().inject(td);
242 Cloud<ParcelType>::move(td);
248 template<class ParcelType>
249 void Foam::KinematicCloud<ParcelType>::info() const
251 Info<< "Cloud: " << this->name() << nl
252 << " Total number of parcels added = "
253 << returnReduce(this->injection().parcelsAddedTotal(), sumOp<label>())
255 << " Total mass introduced = "
256 << returnReduce(this->injection().massInjected(), sumOp<scalar>())
258 << " Current number of parcels = "
259 << returnReduce(this->size(), sumOp<label>()) << nl
260 << " Current mass in system = "
261 << returnReduce(massInSystem(), sumOp<scalar>()) << nl;
265 // ************************************************************************* //