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  parcelsPerSecond_(this->coeffDict().getScalar("parcelsPerSecond")),
190  flowRateProfile_
191  (
193  (
194  "flowRateProfile",
195  this->coeffDict(),
196  &owner.mesh()
197  )
198  ),
199  thetaInner_
200  (
202  (
203  "thetaInner",
204  this->coeffDict(),
205  &owner.mesh()
206  )
207  ),
208  thetaOuter_
209  (
211  (
212  "thetaOuter",
213  this->coeffDict(),
214  &owner.mesh()
215  )
216  ),
217  sizeDistribution_
218  (
220  (
221  this->coeffDict().subDict("sizeDistribution"),
222  owner.rndGen()
223  )
224  ),
225  tanVec1_(Zero),
226  tanVec2_(Zero),
227  normal_(Zero),
228 
229  UMag_(0.0),
230  Cd_(nullptr),
231  Pinj_(nullptr)
232 {
233  if (innerDiameter_ >= outerDiameter_)
234  {
236  << "Inner diameter must be less than the outer diameter:" << nl
237  << " innerDiameter: " << innerDiameter_ << nl
238  << " outerDiameter: " << outerDiameter_
239  << exit(FatalError);
240  }
241 
242  // Convert from user time to reduce the number of time conversion calls
243  const Time& time = owner.db().time();
244  duration_ = time.userTimeToTime(duration_);
245  flowRateProfile_->userTimeToTime(time);
246  thetaInner_->userTimeToTime(time);
247  thetaOuter_->userTimeToTime(time);
248 
249  setInjectionGeometry();
250 
251  setFlowType();
252 
253 
254  // Set total volume to inject
255  this->volumeTotal_ = flowRateProfile_->integrate(0.0, duration_);
256 
257  updateMesh();
258 }
259 
260 
261 template<class CloudType>
263 (
265 )
266 :
268  injectionMethod_(im.injectionMethod_),
269  flowType_(im.flowType_),
270  outerDiameter_(im.outerDiameter_),
271  innerDiameter_(im.innerDiameter_),
272  duration_(im.duration_),
273  positionVsTime_(im.positionVsTime_.clone()),
274  position_(im.position_),
275  injectorCell_(im.injectorCell_),
276  tetFacei_(im.tetFacei_),
277  tetPti_(im.tetPti_),
278  directionVsTime_(im.directionVsTime_.clone()),
279  direction_(im.direction_),
280  parcelsPerSecond_(im.parcelsPerSecond_),
281  flowRateProfile_(im.flowRateProfile_.clone()),
282  thetaInner_(im.thetaInner_.clone()),
283  thetaOuter_(im.thetaOuter_.clone()),
284  sizeDistribution_(im.sizeDistribution_.clone()),
285  tanVec1_(im.tanVec1_),
286  tanVec2_(im.tanVec2_),
287  normal_(im.normal_),
288  UMag_(im.UMag_),
289  Cd_(im.Cd_.clone()),
290  Pinj_(im.Pinj_.clone())
291 {}
292 
293 
294 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
295 
296 template<class CloudType>
298 {
299  // Set/cache the injector cell info for static methods
300  if (positionVsTime_->constant())
301  {
302  position_ = positionVsTime_->value(0);
303 
304  this->findCellAtPosition
305  (
306  injectorCell_,
307  tetFacei_,
308  tetPti_,
309  position_
310  );
311  }
312 }
313 
314 
315 template<class CloudType>
317 {
318  return this->SOI_ + duration_;
319 }
320 
321 
322 template<class CloudType>
324 (
325  const scalar time0,
326  const scalar time1
327 )
328 {
329  if ((time0 >= 0.0) && (time0 < duration_))
330  {
331  return floor((time1 - time0)*parcelsPerSecond_);
332  }
333 
334  return 0;
335 }
336 
337 
338 template<class CloudType>
340 (
341  const scalar time0,
342  const scalar time1
343 )
344 {
345  if ((time0 >= 0.0) && (time0 < duration_))
346  {
347  return flowRateProfile_->integrate(time0, time1);
348  }
349 
350  return 0.0;
351 }
352 
353 
354 template<class CloudType>
356 (
357  const label,
358  const label,
359  const scalar time,
360  vector& position,
361  label& cellOwner,
362  label& tetFacei,
363  label& tetPti
364 )
365 {
366  Random& rndGen = this->owner().rndGen();
367  const scalar t = time - this->SOI_;
368 
369  if (!directionVsTime_->constant())
370  {
371  direction_ = directionVsTime_->value(t);
372  direction_.normalise();
373 
374  // Determine direction vectors tangential to direction
375  vector tangent = Zero;
376  scalar magTangent = 0.0;
377 
378  while(magTangent < SMALL)
379  {
380  vector v = rndGen.globalSample01<vector>();
381 
382  tangent = v - (v & direction_)*direction_;
383  magTangent = mag(tangent);
384  }
385 
386  tanVec1_ = tangent/magTangent;
387  tanVec2_ = direction_^tanVec1_;
388  }
389 
390  scalar beta = mathematical::twoPi*rndGen.globalSample01<scalar>();
391  normal_ = tanVec1_*cos(beta) + tanVec2_*sin(beta);
392 
393  switch (injectionMethod_)
394  {
395  case injectionMethod::imPoint:
396  {
397  if (positionVsTime_->constant())
398  {
399  position = position_;
400  cellOwner = injectorCell_;
401  tetFacei = tetFacei_;
402  tetPti = tetPti_;
403  }
404  else
405  {
406  position = positionVsTime_->value(t);
407 
408  this->findCellAtPosition
409  (
410  cellOwner,
411  tetFacei,
412  tetPti,
413  position
414  );
415  }
416  break;
417  }
418  case injectionMethod::imDisc:
419  {
420  scalar frac = rndGen.globalSample01<scalar>();
421  scalar dr = outerDiameter_ - innerDiameter_;
422  scalar r = 0.5*(innerDiameter_ + frac*dr);
423 
424  position = positionVsTime_->value(t) + r*normal_;
425 
426  this->findCellAtPosition
427  (
428  cellOwner,
429  tetFacei,
430  tetPti,
431  position
432  );
433  break;
434  }
435  default:
436  {
438  << "Unhandled injection method "
439  << injectionMethodNames[injectionMethod_]
440  << exit(FatalError);
441  }
442  }
443 }
444 
445 
446 template<class CloudType>
448 (
449  const label parcelI,
450  const label,
451  const scalar time,
452  typename CloudType::parcelType& parcel
453 )
454 {
455  Random& rndGen = this->owner().rndGen();
456 
457  // Set particle velocity
458  scalar t = time - this->SOI_;
459  scalar ti = thetaInner_->value(t);
460  scalar to = thetaOuter_->value(t);
461  scalar coneAngle = degToRad(rndGen.sample01<scalar>()*(to - ti) + ti);
462 
463  scalar alpha = sin(coneAngle);
464  scalar dcorr = cos(coneAngle);
465 
466  vector normal = alpha*normal_;
467  vector dirVec = dcorr*direction_;
468  dirVec += normal;
469  dirVec.normalise();
470 
471  switch (flowType_)
472  {
473  case flowType::ftConstantVelocity:
474  {
475  parcel.U() = UMag_*dirVec;
476  break;
477  }
478  case flowType::ftPressureDrivenVelocity:
479  {
480  scalar pAmbient = this->owner().pAmbient();
481  scalar rho = parcel.rho();
482  scalar UMag = ::sqrt(2.0*(Pinj_->value(t) - pAmbient)/rho);
483  parcel.U() = UMag*dirVec;
484  break;
485  }
486  case flowType::ftFlowRateAndDischarge:
487  {
488  scalar Ao = 0.25*mathematical::pi*outerDiameter_*outerDiameter_;
489  scalar Ai = 0.25*mathematical::pi*innerDiameter_*innerDiameter_;
490  scalar massFlowRate =
491  this->massTotal()
492  *flowRateProfile_->value(t)
493  /this->volumeTotal();
494 
495  scalar Umag = massFlowRate/(parcel.rho()*Cd_->value(t)*(Ao - Ai));
496  parcel.U() = Umag*dirVec;
497  break;
498  }
499  default:
500  {
502  << "Unhandled injection method "
503  << flowTypeNames[flowType_]
504  << exit(FatalError);
505  }
506  }
507 
508  // Set particle diameter
509  parcel.d() = sizeDistribution_->sample();
510 }
511 
512 
513 template<class CloudType>
515 {
516  return false;
517 }
518 
519 
520 template<class CloudType>
522 {
523  return true;
524 }
525 
526 
527 // ************************************************************************* //
Foam::ConeNozzleInjection::injectionMethod
injectionMethod
Injection method enumeration.
Definition: ConeNozzleInjection.H:96
Foam::Random
Random number generator.
Definition: Random.H:59
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:57
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::ConeNozzleInjection::parcelsToInject
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
Definition: ConeNozzleInjection.C:324
Foam::constant
Different types of constants.
Definition: atomicConstants.C:38
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
Function1.H
Foam::TimeState::userTimeToTime
virtual scalar userTimeToTime(const scalar theta) const
Convert the user-time (e.g. CA deg) to real-time (s).
Definition: TimeState.C:49
ConeNozzleInjection.H
Foam::ConeNozzleInjection::validInjection
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
Definition: ConeNozzleInjection.C:521
unitConversion.H
Unit conversion functions.
Foam::InjectionModel
Templated injection model class.
Definition: InjectionModel.H:73
Foam::Vector::normalise
Vector< Cmpt > & normalise()
Normalise the vector by its magnitude.
Definition: VectorI.H:123
rho
rho
Definition: readInitialConditions.H:88
Foam::Function1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: propellerInfo.H:291
Foam::DSMCCloud::rndGen
Random & rndGen()
Return reference to the random object.
Definition: DSMCCloudI.H:124
Foam::ConeNozzleInjection::volumeToInject
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
Definition: ConeNozzleInjection.C:340
Foam::ConeNozzleInjection::injectionMethodNames
static const Enum< injectionMethod > injectionMethodNames
Definition: ConeNozzleInjection.H:102
Foam::DSMCCloud::mesh
const fvMesh & mesh() const
Return reference to the mesh.
Definition: DSMCCloudI.H:44
Foam::ConeNozzleInjection::flowTypeNames
static const Enum< flowType > flowTypeNames
Definition: ConeNozzleInjection.H:112
Foam::ConeNozzleInjection::ConeNozzleInjection
ConeNozzleInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
Definition: ConeNozzleInjection.C:167
beta
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
Foam::ConeNozzleInjection
Cone injection.
Definition: ConeNozzleInjection.H:89
Foam::ConeNozzleInjection::updateMesh
virtual void updateMesh()
Set injector locations when mesh is updated.
Definition: ConeNozzleInjection.C:297
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:71
Foam::constant::mathematical::twoPi
constexpr scalar twoPi(2 *M_PI)
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::degToRad
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
Definition: unitConversion.H:48
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::distributionModel::New
static autoPtr< distributionModel > New(const dictionary &dict, Random &rndGen)
Selector.
Definition: distributionModelNew.C:34
Foam::ConeNozzleInjection::timeEnd
scalar timeEnd() const
Return the end-of-injection time.
Definition: ConeNozzleInjection.C:316
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::Vector< scalar >
Foam::ConeNozzleInjection::flowType
flowType
Flow type enumeration.
Definition: ConeNozzleInjection.H:105
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::ConeNozzleInjection::setProperties
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
Definition: ConeNozzleInjection.C:448
Foam::primitiveMesh::reset
void reset(const label nPoints, const label nInternalFaces, const label nFaces, const label nCells)
Reset this primitiveMesh given the primitive array sizes.
Definition: primitiveMesh.C:207
Foam::ConeNozzleInjection::fullyDescribed
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
Definition: ConeNozzleInjection.C:514
rndGen
Random rndGen
Definition: createFields.H:23
Foam::DSMCCloud::parcelType
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
distributionModel.H
Foam::ConeNozzleInjection::setPositionAndCell
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.
Definition: ConeNozzleInjection.C:356
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265