1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright held by original author
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 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
39 template<class ParcelType>
40 void Foam::KinematicCloud<ParcelType>::preEvolve()
42 this->dispersion().cacheFields(true);
43 forces_.cacheFields(true);
47 template<class ParcelType>
48 void Foam::KinematicCloud<ParcelType>::evolveCloud()
50 autoPtr<interpolation<scalar> > rhoInterpolator =
51 interpolation<scalar>::New
53 interpolationSchemes_,
57 autoPtr<interpolation<vector> > UInterpolator =
58 interpolation<vector>::New
60 interpolationSchemes_,
64 autoPtr<interpolation<scalar> > muInterpolator =
65 interpolation<scalar>::New
67 interpolationSchemes_,
71 typename ParcelType::trackData td
81 this->injection().inject(td);
88 Cloud<ParcelType>::move(td);
92 template<class ParcelType>
93 void Foam::KinematicCloud<ParcelType>::postEvolve()
97 this->writePositions();
100 this->dispersion().cacheFields(false);
101 forces_.cacheFields(false);
103 this->postProcessing().post();
107 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
109 template<class ParcelType>
110 Foam::KinematicCloud<ParcelType>::KinematicCloud
112 const word& cloudName,
113 const volScalarField& rho,
114 const volVectorField& U,
115 const volScalarField& mu,
116 const dimensionedVector& g,
120 Cloud<ParcelType>(rho.mesh(), cloudName, false),
127 cloudName + "Properties",
128 rho.mesh().time().constant(),
134 constProps_(particleProperties_),
135 active_(particleProperties_.lookup("active")),
136 parcelTypeId_(readLabel(particleProperties_.lookup("parcelTypeId"))),
137 coupled_(particleProperties_.lookup("coupled")),
138 cellValueSourceCorrection_
140 particleProperties_.lookup("cellValueSourceCorrection")
147 forces_(mesh_, particleProperties_, g_.value()),
148 interpolationSchemes_(particleProperties_.subDict("interpolationSchemes")),
151 DispersionModel<KinematicCloud<ParcelType> >::New
159 DragModel<KinematicCloud<ParcelType> >::New
167 InjectionModel<KinematicCloud<ParcelType> >::New
173 patchInteractionModel_
175 PatchInteractionModel<KinematicCloud<ParcelType> >::New
183 PostProcessingModel<KinematicCloud<ParcelType> >::New
185 this->particleProperties_,
191 vectorIntegrationScheme::New
194 particleProperties_.subDict("integrationSchemes")
201 this->name() + "UTrans",
202 this->db().time().timeName(),
209 dimensionedVector("zero", dimMass*dimVelocity, vector::zero)
214 ParcelType::readFields(*this);
219 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
221 template<class ParcelType>
222 Foam::KinematicCloud<ParcelType>::~KinematicCloud()
226 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
228 template<class ParcelType>
229 void Foam::KinematicCloud<ParcelType>::checkParcelProperties
232 const scalar lagrangianDt,
233 const bool fullyDescribed
238 parcel.rho() = constProps_.rho0();
241 scalar carrierDt = this->db().time().deltaTValue();
242 parcel.stepFraction() = (carrierDt - lagrangianDt)/carrierDt;
246 template<class ParcelType>
247 void Foam::KinematicCloud<ParcelType>::resetSourceTerms()
249 UTrans_.field() = vector::zero;
253 template<class ParcelType>
254 void Foam::KinematicCloud<ParcelType>::evolve()
270 template<class ParcelType>
271 void Foam::KinematicCloud<ParcelType>::info() const
273 Info<< "Cloud: " << this->name() << nl
274 << " Total number of parcels added = "
275 << this->injection().parcelsAddedTotal() << nl
276 << " Total mass introduced = "
277 << this->injection().massInjected() << nl
278 << " Current number of parcels = "
279 << returnReduce(this->size(), sumOp<label>()) << nl
280 << " Current mass in system = "
281 << returnReduce(massInSystem(), sumOp<scalar>()) << nl;
285 // ************************************************************************* //