Loading...
Searching...
No Matches
ConeNozzleInjection.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2015-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "ConeNozzleInjection.H"
30#include "Function1.H"
31#include "unitConversion.H"
32#include "distributionModel.H"
33#include "axisAngleRotation.H"
34
35using namespace Foam::constant;
36
37
38template<class CloudType>
39const Foam::Enum
40<
42>
44({
45 { injectionMethod::imPoint, "point" },
46 { injectionMethod::imDisc, "disc" },
47 { injectionMethod::imDiscSegments, "discSegments" },
48});
49
50
51template<class CloudType>
52const Foam::Enum
53<
55>
57({
58 { flowType::ftConstantVelocity, "constantVelocity" },
59 { flowType::ftPressureDrivenVelocity, "pressureDrivenVelocity" },
60 { flowType::ftFlowRateAndDischarge, "flowRateAndDischarge" },
61});
62
63
64// * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * * //
65
66template<class CloudType>
67void Foam::ConeNozzleInjection<CloudType>::setInjectionGeometry()
68{
69 const auto& mesh = this->owner().mesh();
70
71 // Position
72 positionVsTime_.reset
73 (
74 Function1<vector>::New("position", this->coeffDict(), &mesh)
75 );
76
77 positionVsTime_->userTimeToTime(this->owner().time());
78
79 if (positionVsTime_->constant())
80 {
81 position_ = positionVsTime_->value(0);
82 }
83
84 // Direction
85 directionVsTime_.reset
86 (
87 Function1<vector>::New("direction", this->coeffDict(), &mesh)
88 );
89
90 directionVsTime_->userTimeToTime(this->owner().time());
91
92 if (directionVsTime_->constant())
93 {
94 direction_ = directionVsTime_->value(0);
95 direction_.normalise();
96
97 Random& rndGen = this->owner().rndGen();
98
99 // Determine direction vectors tangential to direction
100 vector tangent = Zero;
101 scalar magTangent = 0.0;
102
103 while(magTangent < SMALL)
104 {
105 vector v = rndGen.globalSample01<vector>();
106
107 tangent = v - (v & direction_)*direction_;
108 magTangent = mag(tangent);
109 }
110
111 tanVec1_ = tangent/magTangent;
112 tanVec2_ = direction_^tanVec1_;
113 }
114}
115
116
117template<class CloudType>
118void Foam::ConeNozzleInjection<CloudType>::setFlowType()
119{
120 switch (flowType_)
121 {
122 case flowType::ftConstantVelocity:
123 {
124 this->coeffDict().readCompat("Umag", {{"UMag", -2506}}, UMag_);
125 break;
126 }
127 case flowType::ftPressureDrivenVelocity:
128 {
129 Pinj_.reset
130 (
131 Function1<scalar>::New
132 (
133 "Pinj",
134 this->coeffDict(),
135 &this->owner().mesh()
136 )
137 );
138 Pinj_->userTimeToTime(this->owner().time());
139 break;
140 }
141 case flowType::ftFlowRateAndDischarge:
142 {
143 Cd_.reset
144 (
145 Function1<scalar>::New
146 (
147 "Cd",
148 this->coeffDict(),
149 &this->owner().mesh()
150 )
151 );
152 Cd_->userTimeToTime(this->owner().time());
153 break;
154 }
155 default:
156 {
158 << "Unhandled flow type "
159 << flowTypeNames[flowType_]
160 << exit(FatalError);
161 }
163}
164
165
166// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
167
168template<class CloudType>
170(
171 const dictionary& dict,
172 CloudType& owner,
173 const word& modelName
174)
175:
176 InjectionModel<CloudType>(dict, owner, modelName, typeName),
177 injectionMethod_
178 (
179 injectionMethodNames.get("injectionMethod", this->coeffDict())
180 ),
181 flowType_(flowTypeNames.get("flowType", this->coeffDict())),
182 outerDiameter_(this->coeffDict().getScalar("outerDiameter")),
183 innerDiameter_(this->coeffDict().getScalar("innerDiameter")),
184 duration_(this->coeffDict().getScalar("duration")),
185 positionVsTime_(nullptr),
186 position_(Zero),
187 injectorCell_(-1),
188 tetFacei_(-1),
189 tetPti_(-1),
190 directionVsTime_(nullptr),
191 direction_(Zero),
192 omegaPtr_
193 (
194 Function1<scalar>::NewIfPresent
195 (
196 "omega",
197 this->coeffDict(),
198 &owner.mesh()
199 )
200 ),
201 parcelsPerSecond_(this->coeffDict().getScalar("parcelsPerSecond")),
202 flowRateProfile_
203 (
204 Function1<scalar>::New
205 (
206 "flowRateProfile",
207 this->coeffDict(),
208 &owner.mesh()
209 )
210 ),
211 thetaInner_
212 (
213 Function1<scalar>::New
214 (
215 "thetaInner",
216 this->coeffDict(),
217 &owner.mesh()
218 )
219 ),
220 thetaOuter_
221 (
222 Function1<scalar>::New
223 (
224 "thetaOuter",
225 this->coeffDict(),
226 &owner.mesh()
227 )
228 ),
229 sizeDistribution_
230 (
232 (
233 this->coeffDict().subDict("sizeDistribution"),
234 owner.rndGen()
235 )
236 ),
237 t0_(this->template getModelProperty<scalar>("t0")),
238 nInjectors_
239 (
240 this->coeffDict().template getOrDefault<scalar>("nInjectors", 1)
241 ),
242 Uinjector_(Zero),
243 initialInjectorDir_
244 (
245 this->coeffDict().template getOrDefault<vector>
246 (
247 "initialInjectorDir",
248 Zero
249 )
250 ),
251 tanVec1_(Zero),
252 tanVec2_(Zero),
253 normal_(Zero),
254
255 UMag_(0.0),
256 Cd_(nullptr),
257 Pinj_(nullptr)
258{
259 if (innerDiameter_ >= outerDiameter_)
260 {
262 << "Inner diameter must be less than the outer diameter:" << nl
263 << " innerDiameter: " << innerDiameter_ << nl
264 << " outerDiameter: " << outerDiameter_
265 << exit(FatalError);
266 }
267
268 if (nInjectors_ < SMALL)
269 {
271 << "Number of injectors in angular-segmented disc "
272 << "must be positive" << nl
273 << " nInjectors: " << nInjectors_ << nl
274 << exit(FatalIOError);
275 }
276
277 // Convert from user time to reduce the number of time conversion calls
278 const Time& time = owner.db().time();
279 duration_ = time.userTimeToTime(duration_);
280 flowRateProfile_->userTimeToTime(time);
281 thetaInner_->userTimeToTime(time);
282 thetaOuter_->userTimeToTime(time);
283
284 if (omegaPtr_)
285 {
286 omegaPtr_->userTimeToTime(time);
287 }
288
289 setInjectionGeometry();
290
291 setFlowType();
292
293
294 // Set total volume to inject
295 this->volumeTotal_ = flowRateProfile_->integrate(0.0, duration_);
296
297 updateMesh();
298}
299
300
301template<class CloudType>
303(
305)
306:
308 injectionMethod_(im.injectionMethod_),
309 flowType_(im.flowType_),
310 outerDiameter_(im.outerDiameter_),
311 innerDiameter_(im.innerDiameter_),
312 duration_(im.duration_),
313 positionVsTime_(im.positionVsTime_.clone()),
314 position_(im.position_),
315 injectorCell_(im.injectorCell_),
316 tetFacei_(im.tetFacei_),
317 tetPti_(im.tetPti_),
318 directionVsTime_(im.directionVsTime_.clone()),
319 direction_(im.direction_),
320 omegaPtr_(im.omegaPtr_.clone()),
321 parcelsPerSecond_(im.parcelsPerSecond_),
322 flowRateProfile_(im.flowRateProfile_.clone()),
323 thetaInner_(im.thetaInner_.clone()),
324 thetaOuter_(im.thetaOuter_.clone()),
325 sizeDistribution_(im.sizeDistribution_.clone()),
326 t0_(im.t0_),
327 nInjectors_(im.nInjectors_),
328 Uinjector_(im.Uinjector_),
329 initialInjectorDir_(im.initialInjectorDir_),
330 tanVec1_(im.tanVec1_),
331 tanVec2_(im.tanVec2_),
332 normal_(im.normal_),
333 UMag_(im.UMag_),
334 Cd_(im.Cd_.clone()),
335 Pinj_(im.Pinj_.clone())
336{}
337
338
339// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
340
341template<class CloudType>
343{
344 // Set/cache the injector cell info for static methods
345 if (positionVsTime_->constant())
346 {
347 position_ = positionVsTime_->value(0);
348
349 this->findCellAtPosition
350 (
351 injectorCell_,
352 tetFacei_,
353 tetPti_,
354 position_
355 );
356 }
357}
358
359
360template<class CloudType>
362{
363 return this->SOI_ + duration_;
364}
365
366
367template<class CloudType>
369(
370 const scalar time0,
371 const scalar time1
372)
373{
374 if ((time0 >= 0.0) && (time0 < duration_))
375 {
376 return floor((time1 - time0)*parcelsPerSecond_);
378
379 return 0;
380}
381
382
383template<class CloudType>
385(
386 const scalar time0,
387 const scalar time1
388)
389{
390 if ((time0 >= 0.0) && (time0 < duration_))
391 {
392 return flowRateProfile_->integrate(time0, time1);
394
395 return 0.0;
396}
397
398
399template<class CloudType>
401(
402 const label,
403 const label,
404 const scalar time,
405 vector& position,
406 label& cellOwner,
407 label& tetFacei,
408 label& tetPti
409)
410{
411 Random& rndGen = this->owner().rndGen();
412 const scalar t = time - this->SOI_;
413
414 if (!directionVsTime_->constant())
415 {
416 direction_ = directionVsTime_->value(t);
417 direction_.normalise();
418
419 // Determine direction vectors tangential to direction
420 vector tangent = Zero;
421 scalar magTangent = 0.0;
422
423 while(magTangent < SMALL)
424 {
425 vector v = rndGen.globalSample01<vector>();
426
427 tangent = v - (v & direction_)*direction_;
428 magTangent = mag(tangent);
429 }
430
431 tanVec1_ = tangent/magTangent;
432 tanVec2_ = direction_^tanVec1_;
433 }
434
436 normal_ = tanVec1_*cos(beta) + tanVec2_*sin(beta);
437
438 switch (injectionMethod_)
439 {
440 case injectionMethod::imPoint:
441 {
442 if (positionVsTime_->constant())
443 {
444 position = position_;
445 cellOwner = injectorCell_;
446 tetFacei = tetFacei_;
447 tetPti = tetPti_;
448 }
449 else
450 {
451 position = positionVsTime_->value(t);
452
453 // Estimate the moving injector velocity
454 const vector position0(positionVsTime_->value(t0_));
455 const scalar dt = t - t0_;
456 if (dt > 0)
457 {
458 Uinjector_ = (position - position0)/dt;
459 }
460
461 this->findCellAtPosition
462 (
463 cellOwner,
464 tetFacei,
465 tetPti,
466 position
467 );
468 }
469 break;
470 }
471 case injectionMethod::imDisc:
472 {
473 scalar frac = rndGen.globalSample01<scalar>();
474 scalar dr = outerDiameter_ - innerDiameter_;
475 scalar r = 0.5*(innerDiameter_ + frac*dr);
476
477 position = positionVsTime_->value(t) + r*normal_;
478
479 // Estimate the moving injector velocity
480 const vector position0(positionVsTime_->value(t0_) + r*normal_);
481 const scalar dt = t - t0_;
482 if (dt > 0)
483 {
484 Uinjector_ = (position - position0)/dt;
485 }
486
487 this->findCellAtPosition
488 (
489 cellOwner,
490 tetFacei,
491 tetPti,
492 position
493 );
494 break;
495 }
496 case injectionMethod::imDiscSegments:
497 {
498 // Calculate the uniform angular increment in radians
499 const scalar angleIncrement = mathematical::twoPi/nInjectors_;
500
501 // Randomly set the index of injector angles
502 const label injectorIndex =
503 rndGen.globalPosition<label>(0, nInjectors_ - 1);
504
505 // Calculate the angle for the current injection
506 const scalar angle = injectorIndex*angleIncrement;
507
508 // Calculate the rotation tensor around an axis by an angle
509 const tensor R
510 (
512 (
513 direction_,
514 angle,
515 false // radians
516 )
517 );
518
519 // Compute a random radius between innerDiameter_ and outerDiameter_
520 const scalar fraction = rndGen.globalSample01<scalar>();
521 const scalar dr = outerDiameter_ - innerDiameter_;
522 const scalar radius = 0.5*(innerDiameter_ + fraction*dr);
523
524 // Rotate the initial direction to get the normal direction
525 // for this injector
526 normal_ = R & initialInjectorDir_;
527
528 // Compute the particle's injection position
529 const vector radialOffset(radius*normal_);
530 position = positionVsTime_->value(t) + radialOffset;
531
532 // Estimate the injector velocity if time has advanced
533 const scalar dt = t - t0_;
534 if (dt > 0)
535 {
536 const vector position0
537 (
538 positionVsTime_->value(t0_) + radialOffset
539 );
540 Uinjector_ = (position - position0)/dt;
541 }
542
543 this->findCellAtPosition
544 (
545 cellOwner,
546 tetFacei,
547 tetPti,
548 position
549 );
550 break;
551 }
552 default:
553 {
555 << "Unhandled injection method "
556 << injectionMethodNames[injectionMethod_]
557 << exit(FatalError);
558 }
559 }
561 // Update the previous time step
562 t0_ = t;
563}
564
565
566template<class CloudType>
568(
569 const label parcelI,
570 const label,
571 const scalar time,
572 typename CloudType::parcelType& parcel
573)
574{
575 Random& rndGen = this->owner().rndGen();
576
577 // Set particle velocity
578 scalar t = time - this->SOI_;
579 scalar ti = thetaInner_->value(t);
580 scalar to = thetaOuter_->value(t);
581 scalar coneAngle = degToRad(rndGen.sample01<scalar>()*(to - ti) + ti);
582
583 scalar alpha = sin(coneAngle);
584 scalar dcorr = cos(coneAngle);
585
586 vector normal = alpha*normal_;
587 vector dirVec = dcorr*direction_;
588 dirVec += normal;
589 dirVec.normalise();
590
591 switch (flowType_)
592 {
593 case flowType::ftConstantVelocity:
594 {
595 parcel.U() = UMag_*dirVec + Uinjector_;
596 break;
597 }
598 case flowType::ftPressureDrivenVelocity:
599 {
600 scalar pAmbient = this->owner().pAmbient();
601 scalar rho = parcel.rho();
602 scalar UMag = ::sqrt(2.0*(Pinj_->value(t) - pAmbient)/rho);
603 parcel.U() = UMag*dirVec + Uinjector_;
604 break;
605 }
606 case flowType::ftFlowRateAndDischarge:
607 {
608 scalar Ao = 0.25*mathematical::pi*outerDiameter_*outerDiameter_;
609 scalar Ai = 0.25*mathematical::pi*innerDiameter_*innerDiameter_;
610 scalar massFlowRate =
611 this->massTotal()
612 *flowRateProfile_->value(t)
613 /this->volumeTotal();
614
615 scalar Umag = massFlowRate/(parcel.rho()*Cd_->value(t)*(Ao - Ai));
616 parcel.U() = Umag*dirVec + Uinjector_;
617 break;
618 }
619 default:
620 {
622 << "Unhandled injection method "
623 << flowTypeNames[flowType_]
624 << exit(FatalError);
625 }
626 }
627
628
629 if (omegaPtr_)
630 {
631 const scalar omega = omegaPtr_->value(t);
632
633 const vector p0(parcel.position() - positionVsTime_->value(t));
634 const vector r(p0 - (p0 & direction_)*direction_);
635 const scalar rMag = mag(r);
636
637 const vector d = normalised(normal_ ^ dirVec);
638
639 parcel.U() += omega*rMag*d;
640 }
642 // Set particle diameter
643 parcel.d() = sizeDistribution_->sample();
644}
645
646
647template<class CloudType>
649{
650 return false;
651}
652
653
654template<class CloudType>
656{
657 return true;
658}
659
660
661template<class CloudType>
663{
665
666 if (this->writeTime())
667 {
668 this->setModelProperty("t0", t0_);
669 }
670}
671
672
673// ************************************************************************* //
#define R(A, B, C, D, E, F, K, M)
const CloudType & owner() const
Return const access to the owner cloud.
virtual bool writeTime() const
Flag to indicate when to write a property.
ConeNozzleInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
virtual autoPtr< InjectionModel< CloudType > > clone() const
Construct and return a clone.
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
flowType
Flow type enumeration.
virtual void setPositionAndCell(const label parcelI, const label nParcels, const scalar time, vector &position, label &cellOwner, label &tetFacei, label &tetPti)
Set the injection position and owner cell.
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
static const Enum< injectionMethod > injectionMethodNames
injectionMethod
Injection method enumeration.
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
static const Enum< flowType > flowTypeNames
virtual void info()
Write injection info.
virtual void updateMesh()
Set injector locations when mesh is updated.
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
scalar timeEnd() const
Return the end-of-injection time.
Random & rndGen()
Return reference to the random object.
Definition DSMCCloudI.H:117
const vector & U() const
Return const access to velocity.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition Function1.H:92
Templated injection model class.
scalar massTotal() const
Return mass of particles to introduce.
static autoPtr< InjectionModel< CloudType > > New(const dictionary &dict, CloudType &owner)
Selector with lookup from dictionary.
virtual void info()
Write injection info.
scalar volumeTotal() const
Return the total volume to be injected across the event.
scalar volumeTotal_
Total volume of particles introduced by this injector [m^3] Note: scaled to ensure massTotal is achie...
InjectionModel(CloudType &owner)
Construct null from owner.
virtual bool findCellAtPosition(label &celli, label &tetFacei, label &tetPti, vector &position, bool errorOnNotFound=true)
Find the cell that contains the supplied position.
scalar SOI_
Start of injection [s].
Random number generator.
Definition Random.H:56
Type globalSample01()
Return a sample whose components lie in the range [0,1].
Type globalPosition(const Type &start, const Type &end)
Return a sample on the interval [start,end].
virtual scalar userTimeToTime(const scalar theta) const
Convert the user-time (e.g. CA deg) to real-time (s).
Definition TimeState.C:42
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
Vector< Cmpt > & normalise(const scalar tol=ROOTVSMALL)
Inplace normalise the vector by its magnitude.
Definition VectorI.H:114
static tensor rotation(const vector &axis, const scalar angle, bool degrees=false)
The rotation tensor for given axis/angle.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A library of runtime-selectable doubly-truncated probability distribution models. Returns random samp...
vector position() const
Return current particle position.
Definition particleI.H:283
bool getModelProperty(const word &entryName, Type &value) const
Retrieve generic property from the sub-model.
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
const word & modelName() const
Return const access to the name of the sub-model.
void setModelProperty(const word &entryName, const Type &value)
Add generic property to the sub-model.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
const volScalarField & p0
Definition EEqn.H:36
dynamicFvMesh & mesh
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
constexpr scalar pi(M_PI)
constexpr scalar twoPi(2 *M_PI)
Different types of constants.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
DSMCCloud< dsmcParcel > CloudType
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
dimensionedScalar sin(const dimensionedScalar &ds)
Tensor< scalar > tensor
Definition symmTensor.H:57
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensionedScalar cos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & alpha
dictionary dict
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
Unit conversion functions.
Random rndGen