1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
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
13 the Free Software Foundation, either version 3 of the License, or
14 (at your 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, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "multiHoleInjector.H"
27 #include "addToRunTimeSelectionTable.H"
28 #include "unitConversion.H"
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 defineTypeNameAndDebug(multiHoleInjector, 0);
35 addToRunTimeSelectionTable
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 Foam::multiHoleInjector::multiHoleInjector
48 const Foam::dictionary& dict
51 injectorType(t, dict),
52 propsDict_(dict.subDict(typeName + "Props")),
53 centerPosition_(propsDict_.lookup("position")),
54 xyAngle_(readScalar(propsDict_.lookup("xyAngle"))),
55 zAngle_(readScalar(propsDict_.lookup("zAngle"))),
56 nHoles_(readLabel(propsDict_.lookup("nHoles"))),
57 umbrellaAngle_(readScalar(propsDict_.lookup("umbrellaAngle"))),
58 nozzleTipDiameter_(readScalar(propsDict_.lookup("nozzleTipDiameter"))),
59 angleSpacing_(propsDict_.lookup("angleSpacing")),
60 d_(readScalar(propsDict_.lookup("diameter"))),
61 Cd_(readScalar(propsDict_.lookup("Cd"))),
62 mass_(readScalar(propsDict_.lookup("mass"))),
63 nParcels_(readLabel(propsDict_.lookup("nParcels"))),
64 X_(propsDict_.lookup("X")),
65 massFlowRateProfile_(propsDict_.lookup("massFlowRateProfile")),
66 velocityProfile_(massFlowRateProfile_),
67 injectionPressureProfile_(massFlowRateProfile_),
68 CdProfile_(massFlowRateProfile_),
69 TProfile_(propsDict_.lookup("temperatureProfile")),
70 averageParcelMass_(nHoles_*mass_/nParcels_),
73 pressureIndependentVelocity_(true),
74 tangentialInjectionVector1_(nHoles_),
75 tangentialInjectionVector2_(nHoles_)
77 // check if time entries for soi and eoi match
78 if (mag(massFlowRateProfile_[0][0] - TProfile_[0][0]) > SMALL)
82 "multiHoleInjector::multiHoleInjector"
83 "(const time& t, const dictionary dict)"
84 ) << "Start-times do not match for TemperatureProfile and "
85 << "massFlowRateProfile."
91 mag(massFlowRateProfile_.last()[0] - TProfile_.last()[0])
97 "multiHoleInjector::multiHoleInjector"
98 "(const time& t, const dictionary dict)"
99 ) << "End-times do not match for TemperatureProfile and "
100 << "massFlowRateProfile."
101 << abort(FatalError);
104 // convert CA to real time
105 forAll(massFlowRateProfile_, i)
107 massFlowRateProfile_[i][0] =
108 t.userTimeToTime(massFlowRateProfile_[i][0]);
109 velocityProfile_[i][0] = massFlowRateProfile_[i][0];
110 injectionPressureProfile_[i][0] = massFlowRateProfile_[i][0];
115 TProfile_[i][0] = t.userTimeToTime(TProfile_[i][0]);
118 scalar integratedMFR = integrateTable(massFlowRateProfile_);
120 forAll(massFlowRateProfile_, i)
122 // correct the massFlowRateProfile to match the injected mass
123 massFlowRateProfile_[i][1] *= mass_/integratedMFR;
125 CdProfile_[i][0] = massFlowRateProfile_[i][0];
126 CdProfile_[i][1] = Cd_;
129 setTangentialVectors();
131 // check molar fractions
138 if (mag(Xsum - 1.0) > SMALL)
142 "multiHoleInjector::multiHoleInjector"
143 "(const time& t, const dictionary dict)"
144 ) << "X does not add up to 1.0, correcting molar fractions."
154 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
156 Foam::multiHoleInjector::~multiHoleInjector()
160 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
162 void Foam::multiHoleInjector::setTangentialVectors()
164 scalar alpha = degToRad(xyAngle_);
165 scalar phi = degToRad(zAngle_);
167 vector xp(cos(alpha), sin(alpha), 0.0);
168 vector zp(cos(alpha)*sin(phi), sin(alpha)*sin(phi), cos(phi));
169 if (mag(zp-xp) < 1.0e-15)
171 xp = vector(0.0, 0.0, -1.0);
177 // Info<< "xp = " << xp << endl;
178 // Info<< "yp = " << yp << endl;
179 // Info<< "zp = " << zp << endl;
182 scalar u = degToRad(umbrellaAngle_/2.0);
183 for (label i=0; i<nHoles_; i++)
185 angle += angleSpacing_[i];
186 scalar v = degToRad(angle);
187 direction_[i] = cos(v)*sin(u)*xp + sin(v)*sin(u)*yp + cos(u)*zp;
188 vector dp = direction_[i] - (direction_[i] & zp)*direction_[i];
193 position_[i] = centerPosition_ + 0.5*nozzleTipDiameter_*dp;
196 cachedRandom rndGen(label(0), -1);
198 for (label i=0; i<nHoles_; i++)
200 vector tangent(vector::zero);
204 vector testThis = rndGen.sample01<vector>();
206 tangent = testThis - (testThis & direction_[i])*direction_[i];
210 tangentialInjectionVector1_[i] = tangent/magV;
211 tangentialInjectionVector2_[i] =
212 direction_[i] ^ tangentialInjectionVector1_[i];
217 Foam::label Foam::multiHoleInjector::nParcelsToInject
224 mass_*(fractionOfInjection(time1) - fractionOfInjection(time0));
225 label nParcels = label(mInj/averageParcelMass_ + 0.49);
231 const Foam::vector Foam::multiHoleInjector::position(const label n) const
237 Foam::vector Foam::multiHoleInjector::position
242 const scalar angleOfWedge,
243 const vector& axisOfSymmetry,
244 const vector& axisOfWedge,
245 const vector& axisOfWedgeNormal,
251 scalar is = position_[n] & axisOfSymmetry;
252 scalar magInj = mag(position_[n] - is*axisOfSymmetry);
255 axisOfWedge*cos(0.5*angleOfWedge)
256 + axisOfWedgeNormal*sin(0.5*angleOfWedge);
257 halfWedge /= mag(halfWedge);
259 return (is*axisOfSymmetry + magInj*halfWedge);
263 // otherwise, disc injection
264 scalar iRadius = d_*rndGen.sample01<scalar>();
265 scalar iAngle = constant::mathematical::twoPi*rndGen.sample01<scalar>();
272 tangentialInjectionVector1_[n]*cos(iAngle)
273 + tangentialInjectionVector2_[n]*sin(iAngle)
281 Foam::label Foam::multiHoleInjector::nHoles() const
287 Foam::scalar Foam::multiHoleInjector::d() const
293 const Foam::vector& Foam::multiHoleInjector::direction
299 return direction_[i];
303 Foam::scalar Foam::multiHoleInjector::mass
308 const scalar angleOfWedge
312 mass_*(fractionOfInjection(time1) - fractionOfInjection(time0));
314 // correct mass if calculation is 2D
317 mInj *= 0.5*angleOfWedge/constant::mathematical::pi;
324 Foam::scalar Foam::multiHoleInjector::mass() const
330 const Foam::scalarField& Foam::multiHoleInjector::X() const
336 Foam::List<Foam::multiHoleInjector::pair> Foam::multiHoleInjector::T() const
342 Foam::scalar Foam::multiHoleInjector::T(const scalar time) const
344 return getTableValue(TProfile_, time);
348 Foam::scalar Foam::multiHoleInjector::tsoi() const
350 return massFlowRateProfile_.first()[0];
354 Foam::scalar Foam::multiHoleInjector::teoi() const
356 return massFlowRateProfile_.last()[0];
360 Foam::scalar Foam::multiHoleInjector::massFlowRate
365 return getTableValue(massFlowRateProfile_, time);
369 Foam::scalar Foam::multiHoleInjector::injectionPressure
374 return getTableValue(injectionPressureProfile_, time);
378 Foam::scalar Foam::multiHoleInjector::velocity
383 return getTableValue(velocityProfile_, time);
387 Foam::List<Foam::multiHoleInjector::pair> Foam::multiHoleInjector::CdProfile()
394 Foam::scalar Foam::multiHoleInjector::Cd
403 Foam::scalar Foam::multiHoleInjector::fractionOfInjection
408 return integrateTable(massFlowRateProfile_, time)/mass_;
412 Foam::scalar Foam::multiHoleInjector::injectedMass
417 return mass_*fractionOfInjection(t);
421 void Foam::multiHoleInjector::correctProfiles
423 const liquidMixtureProperties& fuel,
424 const scalar referencePressure
427 scalar A = nHoles_*0.25*constant::mathematical::pi*sqr(d_);
429 forAll(velocityProfile_, i)
431 scalar time = velocityProfile_[i][0];
432 scalar rho = fuel.rho(referencePressure, T(time), X_);
433 scalar v = massFlowRateProfile_[i][1]/(Cd_*rho*A);
434 velocityProfile_[i][1] = v;
435 injectionPressureProfile_[i][1] = referencePressure + 0.5*rho*v*v;
440 Foam::vector Foam::multiHoleInjector::tan1(const label n) const
442 return tangentialInjectionVector1_[n];
446 Foam::vector Foam::multiHoleInjector::tan2(const label n) const
448 return tangentialInjectionVector2_[n];
452 // ************************************************************************* //