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 "ConeInjectionMP.H"
28 #include "DataEntry.H"
29 #include "mathematicalConstants.H"
31 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
33 template<class CloudType>
34 Foam::label Foam::ConeInjectionMP<CloudType>::parcelsToInject
40 if ((time0 >= 0.0) && (time0 < duration_))
42 const scalar targetVolume = volumeFlowRate_().integrate(0, time1);
44 const label targetParcels =
45 parcelsPerInjector_*targetVolume/this->volumeTotal_;
47 const label nToInject = targetParcels - nInjected_;
49 nInjected_ += nToInject;
51 return positions_.size()*nToInject;
60 template<class CloudType>
61 Foam::scalar Foam::ConeInjectionMP<CloudType>::volumeToInject
67 if ((time0 >= 0.0) && (time0 < duration_))
69 return volumeFlowRate_().integrate(time0, time1);
78 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
80 template<class CloudType>
81 Foam::ConeInjectionMP<CloudType>::ConeInjectionMP
83 const dictionary& dict,
87 InjectionModel<CloudType>(dict, owner, typeName),
88 positionsFile_(this->coeffDict().lookup("positionsFile")),
94 owner.db().time().constant(),
100 injectorCells_(positions_.size()),
101 axesFile_(this->coeffDict().lookup("axesFile")),
107 owner.db().time().constant(),
113 duration_(readScalar(this->coeffDict().lookup("duration"))),
116 readScalar(this->coeffDict().lookup("parcelsPerInjector"))
120 DataEntry<scalar>::New
128 DataEntry<scalar>::New
136 DataEntry<scalar>::New
144 DataEntry<scalar>::New
154 this->coeffDict().subDict("parcelPDF"),
158 nInjected_(this->parcelsAddedTotal()),
159 tanVec1_(positions_.size()),
160 tanVec2_(positions_.size())
162 // Normalise direction vector and determine direction vectors
163 // tangential to direction
166 axes_[i] /= mag(axes_[i]);
168 vector tangent = vector::zero;
169 scalar magTangent = 0.0;
171 while (magTangent < SMALL)
173 vector v = this->owner().rndGen().vector01();
175 tangent = v - (v & axes_[i])*axes_[i];
176 magTangent = mag(tangent);
179 tanVec1_[i] = tangent/magTangent;
180 tanVec2_[i] = axes_[i]^tanVec1_[i];
183 // Set total volume to inject
184 this->volumeTotal_ = volumeFlowRate_().integrate(0.0, duration_);
186 // Set/cache the injector cells
187 forAll(positions_, i)
189 this->findCellAtPosition
198 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
200 template<class CloudType>
201 Foam::ConeInjectionMP<CloudType>::~ConeInjectionMP()
205 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
207 template<class CloudType>
208 bool Foam::ConeInjectionMP<CloudType>::active() const
214 template<class CloudType>
215 Foam::scalar Foam::ConeInjectionMP<CloudType>::timeEnd() const
217 return this->SOI_ + duration_;
221 template<class CloudType>
222 void Foam::ConeInjectionMP<CloudType>::setPositionAndCell
231 const label i = parcelI%positions_.size();
233 position = positions_[i];
234 cellOwner = injectorCells_[i];
238 template<class CloudType>
239 void Foam::ConeInjectionMP<CloudType>::setProperties
244 typename CloudType::parcelType& parcel
247 // set particle velocity
248 const label i = parcelI%positions_.size();
250 const scalar deg2Rad = mathematicalConstant::pi/180.0;
252 scalar t = time - this->SOI_;
253 scalar ti = thetaInner_().value(t);
254 scalar to = thetaOuter_().value(t);
255 scalar coneAngle = this->owner().rndGen().scalar01()*(to - ti) + ti;
257 coneAngle *= deg2Rad;
258 scalar alpha = sin(coneAngle);
259 scalar dcorr = cos(coneAngle);
260 scalar beta = mathematicalConstant::twoPi*this->owner().rndGen().scalar01();
262 vector normal = alpha*(tanVec1_[i]*cos(beta) + tanVec2_[i]*sin(beta));
263 vector dirVec = dcorr*axes_[i];
266 dirVec /= mag(dirVec);
268 parcel.U() = Umag_().value(t)*dirVec;
270 // set particle diameter
271 parcel.d() = parcelPDF_().sample();
275 template<class CloudType>
276 bool Foam::ConeInjectionMP<CloudType>::fullyDescribed() const
282 template<class CloudType>
283 bool Foam::ConeInjectionMP<CloudType>::validInjection(const label)
289 // ************************************************************************* //