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 -------------------------------------------------------------------------------
11 License
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 
34 using namespace Foam::constant;
35 
36 
37 template<class CloudType>
38 const Foam::Enum
39 <
41 >
43 ({
44  { injectionMethod::imPoint, "point" },
45  { injectionMethod::imDisc, "disc" },
46 });
47 
48 template<class CloudType>
49 const Foam::Enum
50 <
52 >
54 ({
55  { flowType::ftConstantVelocity, "constantVelocity" },
56  { flowType::ftPressureDrivenVelocity, "pressureDrivenVelocity" },
57  { flowType::ftFlowRateAndDischarge, "flowRateAndDischarge" },
58 });
59 
60 
61 // * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * * //
62 
63 template<class CloudType>
65 {
66  const auto& mesh = this->owner().mesh();
67 
68  // Position
69  positionVsTime_.reset
70  (
71  Function1<vector>::New("position", this->coeffDict(), &mesh)
72  );
73 
74  positionVsTime_->userTimeToTime(this->owner().time());
75 
76  if (positionVsTime_->constant())
77  {
78  position_ = positionVsTime_->value(0);
79  }
80 
81  // Direction
82  directionVsTime_.reset
83  (
84  Function1<vector>::New("direction", this->coeffDict(), &mesh)
85  );
86 
87  directionVsTime_->userTimeToTime(this->owner().time());
88 
89  if (directionVsTime_->constant())
90  {
91  direction_ = directionVsTime_->value(0);
92  direction_.normalise();
93 
94  Random& rndGen = this->owner().rndGen();
95 
96  // Determine direction vectors tangential to direction
97  vector tangent = Zero;
98  scalar magTangent = 0.0;
99 
100  while(magTangent < SMALL)
101  {
102  vector v = rndGen.globalSample01<vector>();
103 
104  tangent = v - (v & direction_)*direction_;
105  magTangent = mag(tangent);
106  }
107 
108  tanVec1_ = tangent/magTangent;
109  tanVec2_ = direction_^tanVec1_;
110  }
111 }
112 
113 
114 template<class CloudType>
116 {
117  switch (flowType_)
118  {
119  case flowType::ftConstantVelocity:
120  {
121  this->coeffDict().readEntry("UMag", UMag_);
122  break;
123  }
124  case flowType::ftPressureDrivenVelocity:
125  {
126  Pinj_.reset
127  (
128  Function1<scalar>::New
129  (
130  "Pinj",
131  this->coeffDict(),
132  &this->owner().mesh()
133  )
134  );
135  Pinj_->userTimeToTime(this->owner().time());
136  break;
137  }
138  case flowType::ftFlowRateAndDischarge:
139  {
140  Cd_.reset
141  (
142  Function1<scalar>::New
143  (
144  "Cd",
145  this->coeffDict(),
146  &this->owner().mesh()
147  )
148  );
149  Cd_->userTimeToTime(this->owner().time());
150  break;
151  }
152  default:
153  {
155  << "Unhandled flow type "
156  << flowTypeNames[flowType_]
157  << exit(FatalError);
158  }
159  }
160 }
161 
162 
163 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
164 
165 template<class CloudType>
167 (
168  const dictionary& dict,
169  CloudType& owner,
170  const word& modelName
171 )
172 :
173  InjectionModel<CloudType>(dict, owner, modelName, typeName),
174  injectionMethod_
175  (
176  injectionMethodNames.get("injectionMethod", this->coeffDict())
177  ),
178  flowType_(flowTypeNames.get("flowType", this->coeffDict())),
179  outerDiameter_(this->coeffDict().getScalar("outerDiameter")),
180  innerDiameter_(this->coeffDict().getScalar("innerDiameter")),
181  duration_(this->coeffDict().getScalar("duration")),
182  positionVsTime_(nullptr),
183  position_(Zero),
184  injectorCell_(-1),
185  tetFacei_(-1),
186  tetPti_(-1),
187  directionVsTime_(nullptr),
188  direction_(Zero),
189  omegaPtr_
190  (
191  Function1<scalar>::NewIfPresent
192  (
193  "omega",
194  this->coeffDict(),
195  &owner.mesh()
196  )
197  ),
198  parcelsPerSecond_(this->coeffDict().getScalar("parcelsPerSecond")),
199  flowRateProfile_
200  (
201  Function1<scalar>::New
202  (
203  "flowRateProfile",
204  this->coeffDict(),
205  &owner.mesh()
206  )
207  ),
208  thetaInner_
209  (
210  Function1<scalar>::New
211  (
212  "thetaInner",
213  this->coeffDict(),
214  &owner.mesh()
215  )
216  ),
217  thetaOuter_
218  (
219  Function1<scalar>::New
220  (
221  "thetaOuter",
222  this->coeffDict(),
223  &owner.mesh()
224  )
225  ),
226  sizeDistribution_
227  (
229  (
230  this->coeffDict().subDict("sizeDistribution"),
231  owner.rndGen()
232  )
233  ),
234  tanVec1_(Zero),
235  tanVec2_(Zero),
236  normal_(Zero),
237 
238  UMag_(0.0),
239  Cd_(nullptr),
240  Pinj_(nullptr)
241 {
242  if (innerDiameter_ >= outerDiameter_)
243  {
245  << "Inner diameter must be less than the outer diameter:" << nl
246  << " innerDiameter: " << innerDiameter_ << nl
247  << " outerDiameter: " << outerDiameter_
248  << exit(FatalError);
249  }
250 
251  // Convert from user time to reduce the number of time conversion calls
252  const Time& time = owner.db().time();
253  duration_ = time.userTimeToTime(duration_);
254  flowRateProfile_->userTimeToTime(time);
255  thetaInner_->userTimeToTime(time);
256  thetaOuter_->userTimeToTime(time);
257 
258  if (omegaPtr_)
259  {
260  omegaPtr_->userTimeToTime(time);
261  }
262 
263  setInjectionGeometry();
264 
265  setFlowType();
266 
267 
268  // Set total volume to inject
269  this->volumeTotal_ = flowRateProfile_->integrate(0.0, duration_);
271  updateMesh();
272 }
273 
274 
275 template<class CloudType>
277 (
279 )
280 :
282  injectionMethod_(im.injectionMethod_),
283  flowType_(im.flowType_),
284  outerDiameter_(im.outerDiameter_),
285  innerDiameter_(im.innerDiameter_),
286  duration_(im.duration_),
287  positionVsTime_(im.positionVsTime_.clone()),
288  position_(im.position_),
289  injectorCell_(im.injectorCell_),
290  tetFacei_(im.tetFacei_),
291  tetPti_(im.tetPti_),
292  directionVsTime_(im.directionVsTime_.clone()),
293  direction_(im.direction_),
294  omegaPtr_(im.omegaPtr_.clone()),
295  parcelsPerSecond_(im.parcelsPerSecond_),
296  flowRateProfile_(im.flowRateProfile_.clone()),
297  thetaInner_(im.thetaInner_.clone()),
298  thetaOuter_(im.thetaOuter_.clone()),
299  sizeDistribution_(im.sizeDistribution_.clone()),
300  tanVec1_(im.tanVec1_),
301  tanVec2_(im.tanVec2_),
302  normal_(im.normal_),
303  UMag_(im.UMag_),
304  Cd_(im.Cd_.clone()),
305  Pinj_(im.Pinj_.clone())
306 {}
307 
308 
309 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
310 
311 template<class CloudType>
313 {
314  // Set/cache the injector cell info for static methods
315  if (positionVsTime_->constant())
316  {
317  position_ = positionVsTime_->value(0);
318 
319  this->findCellAtPosition
320  (
321  injectorCell_,
322  tetFacei_,
323  tetPti_,
324  position_
325  );
326  }
327 }
328 
329 
330 template<class CloudType>
332 {
333  return this->SOI_ + duration_;
334 }
335 
336 
337 template<class CloudType>
339 (
340  const scalar time0,
341  const scalar time1
342 )
343 {
344  if ((time0 >= 0.0) && (time0 < duration_))
345  {
346  return floor((time1 - time0)*parcelsPerSecond_);
347  }
349  return 0;
350 }
351 
352 
353 template<class CloudType>
355 (
356  const scalar time0,
357  const scalar time1
358 )
359 {
360  if ((time0 >= 0.0) && (time0 < duration_))
361  {
362  return flowRateProfile_->integrate(time0, time1);
363  }
365  return 0.0;
366 }
367 
368 
369 template<class CloudType>
371 (
372  const label,
373  const label,
374  const scalar time,
375  vector& position,
376  label& cellOwner,
377  label& tetFacei,
378  label& tetPti
379 )
380 {
381  Random& rndGen = this->owner().rndGen();
382  const scalar t = time - this->SOI_;
383 
384  if (!directionVsTime_->constant())
385  {
386  direction_ = directionVsTime_->value(t);
387  direction_.normalise();
388 
389  // Determine direction vectors tangential to direction
390  vector tangent = Zero;
391  scalar magTangent = 0.0;
392 
393  while(magTangent < SMALL)
394  {
396 
397  tangent = v - (v & direction_)*direction_;
398  magTangent = mag(tangent);
399  }
400 
401  tanVec1_ = tangent/magTangent;
402  tanVec2_ = direction_^tanVec1_;
403  }
404 
405  scalar beta = mathematical::twoPi*rndGen.globalSample01<scalar>();
406  normal_ = tanVec1_*cos(beta) + tanVec2_*sin(beta);
407 
408  switch (injectionMethod_)
409  {
410  case injectionMethod::imPoint:
411  {
412  if (positionVsTime_->constant())
413  {
414  position = position_;
415  cellOwner = injectorCell_;
416  tetFacei = tetFacei_;
417  tetPti = tetPti_;
418  }
419  else
420  {
421  position = positionVsTime_->value(t);
422 
423  this->findCellAtPosition
424  (
425  cellOwner,
426  tetFacei,
427  tetPti,
428  position
429  );
430  }
431  break;
432  }
433  case injectionMethod::imDisc:
434  {
435  scalar frac = rndGen.globalSample01<scalar>();
436  scalar dr = outerDiameter_ - innerDiameter_;
437  scalar r = 0.5*(innerDiameter_ + frac*dr);
438 
439  position = positionVsTime_->value(t) + r*normal_;
440 
441  this->findCellAtPosition
442  (
443  cellOwner,
444  tetFacei,
445  tetPti,
446  position
447  );
448  break;
449  }
450  default:
451  {
453  << "Unhandled injection method "
454  << injectionMethodNames[injectionMethod_]
455  << exit(FatalError);
456  }
457  }
458 }
459 
460 
461 template<class CloudType>
463 (
464  const label parcelI,
465  const label,
466  const scalar time,
467  typename CloudType::parcelType& parcel
468 )
469 {
470  Random& rndGen = this->owner().rndGen();
471 
472  // Set particle velocity
473  scalar t = time - this->SOI_;
474  scalar ti = thetaInner_->value(t);
475  scalar to = thetaOuter_->value(t);
476  scalar coneAngle = degToRad(rndGen.sample01<scalar>()*(to - ti) + ti);
477 
478  scalar alpha = sin(coneAngle);
479  scalar dcorr = cos(coneAngle);
480 
481  vector normal = alpha*normal_;
482  vector dirVec = dcorr*direction_;
483  dirVec += normal;
484  dirVec.normalise();
485 
486  switch (flowType_)
487  {
488  case flowType::ftConstantVelocity:
489  {
490  parcel.U() = UMag_*dirVec;
491  break;
492  }
493  case flowType::ftPressureDrivenVelocity:
494  {
495  scalar pAmbient = this->owner().pAmbient();
496  scalar rho = parcel.rho();
497  scalar UMag = ::sqrt(2.0*(Pinj_->value(t) - pAmbient)/rho);
498  parcel.U() = UMag*dirVec;
499  break;
500  }
501  case flowType::ftFlowRateAndDischarge:
502  {
503  scalar Ao = 0.25*mathematical::pi*outerDiameter_*outerDiameter_;
504  scalar Ai = 0.25*mathematical::pi*innerDiameter_*innerDiameter_;
505  scalar massFlowRate =
506  this->massTotal()
507  *flowRateProfile_->value(t)
508  /this->volumeTotal();
509 
510  scalar Umag = massFlowRate/(parcel.rho()*Cd_->value(t)*(Ao - Ai));
511  parcel.U() = Umag*dirVec;
512  break;
513  }
514  default:
515  {
517  << "Unhandled injection method "
518  << flowTypeNames[flowType_]
519  << exit(FatalError);
520  }
521  }
522 
523  if (omegaPtr_)
524  {
525  const scalar omega = omegaPtr_->value(t);
526 
527  const vector p0(parcel.position() - positionVsTime_->value(t));
528  const vector r(p0 - (p0 & direction_)*direction_);
529  const scalar rMag = mag(r);
530 
531  const vector d = normalised(normal_ ^ dirVec);
532 
533  parcel.U() += omega*rMag*d;
534  }
536  // Set particle diameter
537  parcel.d() = sizeDistribution_->sample();
538 }
539 
540 
541 template<class CloudType>
543 {
544  return false;
545 }
546 
547 
548 template<class CloudType>
550 {
551  return true;
552 }
553 
554 
555 // ************************************************************************* //
Different types of constants.
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
dictionary dict
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
Templated injection model class.
Unit conversion functions.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
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.
Random rndGen
Definition: createFields.H:23
dimensionedScalar sqrt(const dimensionedScalar &ds)
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.
static const Enum< flowType > flowTypeNames
Vector< Cmpt > & normalise(const scalar tol=ROOTVSMALL)
Inplace normalise the vector by its magnitude.
Definition: VectorI.H:114
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:674
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
const CloudType & owner() const
Return const access to the owner cloud.
injectionMethod
Injection method enumeration.
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
Type sample01()
Return a sample whose components lie in the range [0,1].
dynamicFvMesh & mesh
dimensionedScalar cos(const dimensionedScalar &ds)
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
constexpr scalar twoPi(2 *M_PI)
A class for handling words, derived from Foam::string.
Definition: word.H:63
Random & rndGen()
Return reference to the random object.
Definition: DSMCCloudI.H:117
const Time & time() const noexcept
Return time registry.
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:450
constexpr scalar pi(M_PI)
flowType
Flow type enumeration.
Vector< scalar > vector
Definition: vector.H:57
Random number generator.
Definition: Random.H:55
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:290
dimensionedScalar sin(const dimensionedScalar &ds)
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
scalar timeEnd() const
Return the end-of-injection time.
virtual void updateMesh()
Set injector locations when mesh is updated.
Type globalSample01()
Return a sample whose components lie in the range [0,1].
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
A library of runtime-selectable doubly-truncated probability distribution models. Returns random samp...
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
ConeNozzleInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
scalar volumeTotal_
Total volume of particles introduced by this injector [m^3] Note: scaled to ensure massTotal is achie...
static const Enum< injectionMethod > injectionMethodNames
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
const volScalarField & p0
Definition: EEqn.H:36
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127