192319f8ced148c3004fceb550623053e638705b
[OpenFOAM-1.5.x.git] / src / lagrangian / dieselSpray / injector / swirlInjector / swirlInjector.C
blob192319f8ced148c3004fceb550623053e638705b
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
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 "swirlInjector.H"
28 #include "addToRunTimeSelectionTable.H"
29 #include "Random.H"
30 #include "mathematicalConstants.H"
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 namespace Foam
36 defineTypeNameAndDebug(swirlInjector, 0);
38 addToRunTimeSelectionTable
40     injectorType,
41     swirlInjector,
42     dictionary
45 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
48 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
50 // Construct from components
51 Foam::swirlInjector::swirlInjector
53     const Foam::Time& t,
54     const Foam::dictionary& dict
57     injectorType(t, dict),
58     propsDict_(dict.subDict(typeName + "Props")),
59     position_(propsDict_.lookup("position")),
60     direction_(propsDict_.lookup("direction")),
61     d_(readScalar(propsDict_.lookup("diameter"))),
62     mass_(readScalar(propsDict_.lookup("mass"))),
63     injectionPressure_(readScalar(propsDict_.lookup("injectionPressure"))),
64     T_(readScalar(propsDict_.lookup("temperature"))),
65     nParcels_(readLabel(propsDict_.lookup("nParcels"))),
66     X_(propsDict_.lookup("X")),
67     massFlowRateProfile_(propsDict_.lookup("massFlowRateProfile")),
68     injectionPressureProfile_(propsDict_.lookup("injectionPressureProfile")),
69     velocityProfile_(injectionPressureProfile_),
70     CdProfile_(injectionPressureProfile_),
71     TProfile_(injectionPressureProfile_),
72     averageParcelMass_(mass_/nParcels_),
73     pressureIndependentVelocity_(false)
75     // convert CA to real time
76     forAll(massFlowRateProfile_, i)
77     {
78         massFlowRateProfile_[i][0] = t.userTimeToTime(massFlowRateProfile_[i][0]);
79     }
80     forAll(injectionPressureProfile_, i)
81     {
82         injectionPressureProfile_[i][0] = t.userTimeToTime(injectionPressureProfile_[i][0]);
83     }
85     // check if time entries match
86     if (mag(massFlowRateProfile_[0][0]-injectionPressureProfile_[0][0]) > SMALL)
87     {
88         FatalError << "swirlInjector::swirlInjector(const time& t, const dictionary dict) " << endl
89             << " start-times do not match for injectionPressureProfile and massFlowRateProfile."
90             << abort(FatalError);
91     }
93     // check if time entries match
94     if (mag(massFlowRateProfile_[massFlowRateProfile_.size()-1][0]-injectionPressureProfile_[injectionPressureProfile_.size()-1][0]) > SMALL)
95     {
96         FatalError << "swirlInjector::swirlInjector(const time& t, const dictionary dict) " << endl
97             << " end-times do not match for injectionPressureProfile and massFlowRateProfile."
98             << abort(FatalError);
99     }
101     scalar integratedMFR = integrateTable(massFlowRateProfile_);
102     scalar integratedPressure = integrateTable(injectionPressureProfile_)/(teoi()-tsoi());
104     forAll(massFlowRateProfile_, i)
105     {
106         // correct the massFlowRateProfile to match the injected mass
107         massFlowRateProfile_[i][1] *= mass_/integratedMFR;
109         velocityProfile_[i][0] = massFlowRateProfile_[i][0];
111         TProfile_[i][0] = massFlowRateProfile_[i][0];
112         TProfile_[i][1] = T_;
114         CdProfile_[i][0] = massFlowRateProfile_[i][0];
116     }
118     forAll(injectionPressureProfile_, i)
119     {
120         // correct the pressureProfile to match the injection pressure
121         injectionPressureProfile_[i][1] *= injectionPressure_/integratedPressure;
122     }
124     // Normalize the direction vector
125     direction_ /= mag(direction_);
126     
127     setTangentialVectors();
129     // check molar fractions
130     scalar Xsum = 0.0;
131     forAll(X_, i)
132     {
133         Xsum += X_[i];
134     }
136     if (mag(Xsum - 1.0) > SMALL)
137     {
138         Info << "Warning!!!\n swirlInjector::swirlInjector(const time& t, Istream& is)"
139             << "X does not add up to 1.0, correcting molar fractions."
140             << endl;
141         forAll(X_, i)
142         {
143             X_[i] /= Xsum;
144         }
145     }
149 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
151 Foam::swirlInjector::~swirlInjector()
155 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
157 void Foam::swirlInjector::setTangentialVectors()
159     Random rndGen(label(0));
160     scalar magV = 0.0;
161     vector tangent;
163     while (magV < SMALL)
164     {
165         vector testThis = rndGen.vector01();
167         tangent = testThis - (testThis & direction_)*direction_;
168         magV = mag(tangent);
169     }
171     tangentialInjectionVector1_ = tangent/magV;
172     tangentialInjectionVector2_ = direction_ ^ tangentialInjectionVector1_;
177 Foam::label Foam::swirlInjector::nParcelsToInject
179     const scalar time0,
180     const scalar time1
181 ) const
184     scalar mInj = mass_*(fractionOfInjection(time1)-fractionOfInjection(time0));
185     label nParcels = label(mInj/averageParcelMass_ + 0.49);
186     
187     return nParcels;
190 const Foam::vector Foam::swirlInjector::position(const label n) const
192     return position_;
195 Foam::vector Foam::swirlInjector::position
197     const label n,
198     const scalar time,
199     const bool twoD,
200     const scalar angleOfWedge,
201     const vector& axisOfSymmetry,
202     const vector& axisOfWedge,
203     const vector& axisOfWedgeNormal,
204     Random& rndGen
205 ) const
207     if (twoD)
208     {
209         scalar is = position_ & axisOfSymmetry;
210         scalar magInj = mag(position_ - is*axisOfSymmetry);
212         vector halfWedge =
213             axisOfWedge*cos(0.5*angleOfWedge)
214           + axisOfWedgeNormal*sin(0.5*angleOfWedge);
215         halfWedge /= mag(halfWedge);
217         return (is*axisOfSymmetry + magInj*halfWedge);
218     }
219     else
220     {
221         // otherwise, disc injection
222         scalar iRadius = d_*rndGen.scalar01();
223         scalar iAngle = 2.0*mathematicalConstant::pi*rndGen.scalar01();
225         return
226         ( 
227             position_
228           + iRadius
229           * (
230               tangentialInjectionVector1_*cos(iAngle)
231             + tangentialInjectionVector2_*sin(iAngle)
232           )
233         );
234         
235     }
237     return position_;
240 Foam::label Foam::swirlInjector::nHoles() const
242     return 1;
245 Foam::scalar Foam::swirlInjector::d() const
247     return d_;
250 const Foam::vector& Foam::swirlInjector::direction
252     const label i,
253     const scalar time
254 ) const
256     return direction_;
259 Foam::scalar Foam::swirlInjector::mass
261     const scalar time0,
262     const scalar time1,
263     const bool twoD,
264     const scalar angleOfWedge
265 ) const
267     scalar mInj = mass_*(fractionOfInjection(time1)-fractionOfInjection(time0));
269     // correct mass if calculation is 2D 
270     if (twoD)
271     {
272         mInj *= 0.5*angleOfWedge/mathematicalConstant::pi;
273     }
275     return mInj;
278 Foam::scalar Foam::swirlInjector::mass() const
280     return mass_;
283 Foam::List<Foam::swirlInjector::pair> Foam::swirlInjector::massFlowRateProfile() const
285     return massFlowRateProfile_;
288 Foam::scalar Foam::swirlInjector::massFlowRate(const scalar time) const
290     return getTableValue(massFlowRateProfile_, time);
293 Foam::List<Foam::swirlInjector::pair> Foam::swirlInjector::injectionPressureProfile() const
295     return injectionPressureProfile_;
298 Foam::scalar Foam::swirlInjector::injectionPressure(const scalar time) const
300     return getTableValue(injectionPressureProfile_, time);
303 Foam::List<Foam::swirlInjector::pair> Foam::swirlInjector::velocityProfile() const
305     return velocityProfile_;
308 Foam::scalar Foam::swirlInjector::velocity(const scalar time) const
310     return getTableValue(velocityProfile_, time);
312     
313 Foam::List<Foam::swirlInjector::pair> Foam::swirlInjector::CdProfile() const
315     return CdProfile_;
318 Foam::scalar Foam::swirlInjector::Cd(const scalar time) const
320     return getTableValue(CdProfile_, time);
323 const Foam::scalarField& Foam::swirlInjector::X() const
325     return X_;
328 Foam::List<Foam::swirlInjector::pair> Foam::swirlInjector::T() const
330     return TProfile_;
333 Foam::scalar Foam::swirlInjector::T(const scalar time) const
335     return T_;
338 Foam::scalar Foam::swirlInjector::tsoi() const
340     return massFlowRateProfile_[0][0];
343 Foam::scalar Foam::swirlInjector::teoi() const
345     return massFlowRateProfile_[massFlowRateProfile_.size()-1][0];
348 Foam::scalar Foam::swirlInjector::fractionOfInjection(const scalar time) const
350     return integrateTable(massFlowRateProfile_, time)/mass_;
354 Foam::scalar Foam::swirlInjector::injectedMass
356     const scalar t
357 ) const
359     return mass_*fractionOfInjection(t);
362 void Foam::swirlInjector::correctProfiles
364     const liquidMixture& fuel,
365     const scalar referencePressure
369     scalar A = 0.25*mathematicalConstant::pi*pow(d_, 2.0);
370     scalar pDummy = 1.0e+5;
371     scalar rho = fuel.rho(pDummy, T_, X_);
373     forAll(velocityProfile_, i)
374     {
375         scalar Pinj = getTableValue(injectionPressureProfile_, massFlowRateProfile_[i][0]);
376         scalar mfr = massFlowRateProfile_[i][1]/(rho*A);
377         scalar v = sqrt(2.0*(Pinj - referencePressure)/rho);
378         velocityProfile_[i][1] = v;
379         CdProfile_[i][1] = mfr/v;
380     }
383 Foam::vector Foam::swirlInjector::tan1(const label n) const
385     return tangentialInjectionVector1_;
388 Foam::vector Foam::swirlInjector::tan2(const label n) const
390     return tangentialInjectionVector2_;
393 // ************************************************************************* //