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-2020 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 "TimeFunction1.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  { injectionMethod::imMovingPoint, "movingPoint" },
47 });
48 
49 template<class CloudType>
50 const Foam::Enum
51 <
53 >
55 ({
56  { flowType::ftConstantVelocity, "constantVelocity" },
57  { flowType::ftPressureDrivenVelocity, "pressureDrivenVelocity" },
58  { flowType::ftFlowRateAndDischarge, "flowRateAndDischarge" },
59 });
60 
61 
62 // * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * * //
63 
64 template<class CloudType>
66 {
67  switch (injectionMethod_)
68  {
69  case injectionMethod::imPoint:
70  case injectionMethod::imDisc:
71  {
72  this->coeffDict().readEntry("position", position_);
73  break;
74  }
75  case injectionMethod::imMovingPoint:
76  {
77  positionVsTime_.reset(this->coeffDict());
78  break;
79  }
80  default:
81  {
83  << "Unhandled injection method "
84  << injectionMethodNames[injectionMethod_]
85  << exit(FatalError);
86  }
87  }
88 }
89 
90 
91 template<class CloudType>
93 {
94  switch (flowType_)
95  {
96  case flowType::ftConstantVelocity:
97  {
98  this->coeffDict().readEntry("UMag", UMag_);
99  break;
100  }
101  case flowType::ftPressureDrivenVelocity:
102  {
103  Pinj_.reset(this->coeffDict());
104  break;
105  }
106  case flowType::ftFlowRateAndDischarge:
107  {
108  Cd_.reset(this->coeffDict());
109  break;
110  }
111  default:
112  {
114  << "Unhandled flow type "
115  << flowTypeNames[flowType_]
116  << exit(FatalError);
117  }
118  }
119 }
120 
121 
122 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
123 
124 template<class CloudType>
126 (
127  const dictionary& dict,
128  CloudType& owner,
129  const word& modelName
130 )
131 :
132  InjectionModel<CloudType>(dict, owner, modelName, typeName),
133  injectionMethod_
134  (
135  injectionMethodNames.get("injectionMethod", this->coeffDict())
136  ),
137  flowType_(flowTypeNames.get("flowType", this->coeffDict())),
138  outerDiameter_(this->coeffDict().getScalar("outerDiameter")),
139  innerDiameter_(this->coeffDict().getScalar("innerDiameter")),
140  duration_(this->coeffDict().getScalar("duration")),
141  positionVsTime_(owner.db().time(), "position"),
142  position_(Zero),
143  injectorCell_(-1),
144  tetFacei_(-1),
145  tetPti_(-1),
146  direction_(this->coeffDict().lookup("direction")),
147  parcelsPerSecond_(this->coeffDict().getScalar("parcelsPerSecond")),
148  flowRateProfile_
149  (
151  (
152  owner.db().time(),
153  "flowRateProfile",
154  this->coeffDict()
155  )
156  ),
157  thetaInner_
158  (
160  (
161  owner.db().time(),
162  "thetaInner",
163  this->coeffDict()
164  )
165  ),
166  thetaOuter_
167  (
169  (
170  owner.db().time(),
171  "thetaOuter",
172  this->coeffDict()
173  )
174  ),
175  sizeDistribution_
176  (
178  (
179  this->coeffDict().subDict("sizeDistribution"),
180  owner.rndGen()
181  )
182  ),
183  tanVec1_(Zero),
184  tanVec2_(Zero),
185  normal_(Zero),
186 
187  UMag_(0.0),
188  Cd_(owner.db().time(), "Cd"),
189  Pinj_(owner.db().time(), "Pinj")
190 {
191  if (innerDiameter_ >= outerDiameter_)
192  {
194  << "Inner diameter must be less than the outer diameter:" << nl
195  << " innerDiameter: " << innerDiameter_ << nl
196  << " outerDiameter: " << outerDiameter_
197  << exit(FatalError);
198  }
199 
200  duration_ = owner.db().time().userTimeToTime(duration_);
201 
202  setInjectionMethod();
203 
204  setFlowType();
205 
206  Random& rndGen = this->owner().rndGen();
207 
208  // Normalise direction vector
209  direction_.normalise();
210 
211  // Determine direction vectors tangential to direction
212  vector tangent = Zero;
213  scalar magTangent = 0.0;
214 
215  while(magTangent < SMALL)
216  {
217  vector v = rndGen.globalSample01<vector>();
218 
219  tangent = v - (v & direction_)*direction_;
220  magTangent = mag(tangent);
221  }
222 
223  tanVec1_ = tangent/magTangent;
224  tanVec2_ = direction_^tanVec1_;
225 
226  // Set total volume to inject
227  this->volumeTotal_ = flowRateProfile_.integrate(0.0, duration_);
228 
229  updateMesh();
230 }
231 
232 
233 template<class CloudType>
235 (
237 )
238 :
240  injectionMethod_(im.injectionMethod_),
241  flowType_(im.flowType_),
242  outerDiameter_(im.outerDiameter_),
243  innerDiameter_(im.innerDiameter_),
244  duration_(im.duration_),
245  positionVsTime_(im.positionVsTime_),
246  position_(im.position_),
247  injectorCell_(im.injectorCell_),
248  tetFacei_(im.tetFacei_),
249  tetPti_(im.tetPti_),
250  direction_(im.direction_),
251  parcelsPerSecond_(im.parcelsPerSecond_),
252  flowRateProfile_(im.flowRateProfile_),
253  thetaInner_(im.thetaInner_),
254  thetaOuter_(im.thetaOuter_),
255  sizeDistribution_(im.sizeDistribution_.clone()),
256  tanVec1_(im.tanVec1_),
257  tanVec2_(im.tanVec2_),
258  normal_(im.normal_),
259  UMag_(im.UMag_),
260  Cd_(im.Cd_),
261  Pinj_(im.Pinj_)
262 {}
263 
264 
265 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
266 
267 template<class CloudType>
269 {}
270 
271 
272 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
273 
274 template<class CloudType>
276 {
277  // Set/cache the injector cell info for static methods
278 
279  switch (injectionMethod_)
280  {
281  case injectionMethod::imPoint:
282  {
283  this->findCellAtPosition
284  (
285  injectorCell_,
286  tetFacei_,
287  tetPti_,
288  position_
289  );
290  break;
291  }
292  default:
293  {
294  // Do nothing for the other methods
295  }
296  }
297 }
298 
299 
300 template<class CloudType>
302 {
303  return this->SOI_ + duration_;
304 }
305 
306 
307 template<class CloudType>
309 (
310  const scalar time0,
311  const scalar time1
312 )
313 {
314  if ((time0 >= 0.0) && (time0 < duration_))
315  {
316  return floor((time1 - time0)*parcelsPerSecond_);
317  }
318 
319  return 0;
320 }
321 
322 
323 template<class CloudType>
325 (
326  const scalar time0,
327  const scalar time1
328 )
329 {
330  if ((time0 >= 0.0) && (time0 < duration_))
331  {
332  return flowRateProfile_.integrate(time0, time1);
333  }
334 
335  return 0.0;
336 }
337 
338 
339 template<class CloudType>
341 (
342  const label,
343  const label,
344  const scalar time,
345  vector& position,
346  label& cellOwner,
347  label& tetFacei,
348  label& tetPti
349 )
350 {
351  Random& rndGen = this->owner().rndGen();
352 
353  scalar beta = mathematical::twoPi*rndGen.globalSample01<scalar>();
354  normal_ = tanVec1_*cos(beta) + tanVec2_*sin(beta);
355 
356  switch (injectionMethod_)
357  {
358  case injectionMethod::imPoint:
359  {
360  position = position_;
361  cellOwner = injectorCell_;
362  tetFacei = tetFacei_;
363  tetPti = tetPti_;
364 
365  break;
366  }
367  case injectionMethod::imMovingPoint:
368  {
369  position = positionVsTime_.value(time - this->SOI_);
370 
371  this->findCellAtPosition
372  (
373  cellOwner,
374  tetFacei,
375  tetPti,
376  position
377  );
378 
379  break;
380  }
381  case injectionMethod::imDisc:
382  {
383  scalar frac = rndGen.globalSample01<scalar>();
384  scalar dr = outerDiameter_ - innerDiameter_;
385  scalar r = 0.5*(innerDiameter_ + frac*dr);
386 
387  position = position_ + r*normal_;
388 
389  this->findCellAtPosition
390  (
391  cellOwner,
392  tetFacei,
393  tetPti,
394  position
395  );
396  break;
397  }
398  default:
399  {
401  << "Unhandled injection method "
402  << injectionMethodNames[injectionMethod_]
403  << exit(FatalError);
404  }
405  }
406 }
407 
408 
409 template<class CloudType>
411 (
412  const label parcelI,
413  const label,
414  const scalar time,
415  typename CloudType::parcelType& parcel
416 )
417 {
418  Random& rndGen = this->owner().rndGen();
419 
420  // Set particle velocity
421  scalar t = time - this->SOI_;
422  scalar ti = thetaInner_.value(t);
423  scalar to = thetaOuter_.value(t);
424  scalar coneAngle = degToRad(rndGen.sample01<scalar>()*(to - ti) + ti);
425 
426  scalar alpha = sin(coneAngle);
427  scalar dcorr = cos(coneAngle);
428 
429  vector normal = alpha*normal_;
430  vector dirVec = dcorr*direction_;
431  dirVec += normal;
432  dirVec.normalise();
433 
434  switch (flowType_)
435  {
436  case flowType::ftConstantVelocity:
437  {
438  parcel.U() = UMag_*dirVec;
439  break;
440  }
441  case flowType::ftPressureDrivenVelocity:
442  {
443  scalar pAmbient = this->owner().pAmbient();
444  scalar rho = parcel.rho();
445  scalar UMag = ::sqrt(2.0*(Pinj_.value(t) - pAmbient)/rho);
446  parcel.U() = UMag*dirVec;
447  break;
448  }
449  case flowType::ftFlowRateAndDischarge:
450  {
451  scalar Ao = 0.25*mathematical::pi*outerDiameter_*outerDiameter_;
452  scalar Ai = 0.25*mathematical::pi*innerDiameter_*innerDiameter_;
453  scalar massFlowRate =
454  this->massTotal()
455  *flowRateProfile_.value(t)
456  /this->volumeTotal();
457 
458  scalar Umag = massFlowRate/(parcel.rho()*Cd_.value(t)*(Ao - Ai));
459  parcel.U() = Umag*dirVec;
460  break;
461  }
462  default:
463  {
465  << "Unhandled injection method "
466  << flowTypeNames[flowType_]
467  << exit(FatalError);
468  }
469  }
470 
471  // Set particle diameter
472  parcel.d() = sizeDistribution_->sample();
473 }
474 
475 
476 template<class CloudType>
478 {
479  return false;
480 }
481 
482 
483 template<class CloudType>
485 {
486  return true;
487 }
488 
489 
490 // ************************************************************************* //
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::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
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:309
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
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:484
unitConversion.H
Unit conversion functions.
Foam::InjectionModel
Templated injection model class.
Definition: InjectionModel.H:73
Foam::TimeFunction1< scalar >
Foam::Vector::normalise
Vector< Cmpt > & normalise()
Normalise the vector by its magnitude.
Definition: VectorI.H:123
rho
rho
Definition: readInitialConditions.H:88
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:325
Foam::ConeNozzleInjection::injectionMethodNames
static const Enum< injectionMethod > injectionMethodNames
Definition: ConeNozzleInjection.H:103
Foam::ConeNozzleInjection::flowTypeNames
static const Enum< flowType > flowTypeNames
Definition: ConeNozzleInjection.H:113
Foam::ConeNozzleInjection::ConeNozzleInjection
ConeNozzleInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
Definition: ConeNozzleInjection.C:126
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:275
Foam::radiation::lookup
Lookup type of boundary radiation properties.
Definition: lookup.H:63
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:121
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:301
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::Vector< scalar >
Foam::ConeNozzleInjection::flowType
flowType
Flow type enumeration.
Definition: ConeNozzleInjection.H:106
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:411
Foam::ConeNozzleInjection::fullyDescribed
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
Definition: ConeNozzleInjection.C:477
rndGen
Random rndGen
Definition: createFields.H:23
Foam::DSMCCloud::parcelType
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
TimeFunction1.H
Foam::ConeNozzleInjection::~ConeNozzleInjection
virtual ~ConeNozzleInjection()
Destructor.
Definition: ConeNozzleInjection.C:268
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:341
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265