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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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
19  for more details.
20 
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/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "ConeNozzleInjection.H"
27 #include "TimeDataEntry.H"
28 #include "mathematicalConstants.H"
29 #include "distributionModel.H"
30 
31 using namespace Foam::constant;
32 
33 // * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * * //
34 
35 template<class CloudType>
37 {
38  word injectionMethodType = this->coeffDict().lookup("injectionMethod");
39  if (injectionMethodType == "disc")
40  {
41  injectionMethod_ = imDisc;
42  }
43  else if (injectionMethodType == "point")
44  {
45  injectionMethod_ = imPoint;
46 
47  // Set/cache the injector cell
48  this->findCellAtPosition
49  (
50  injectorCell_,
51  tetFaceI_,
52  tetPtI_,
53  position_,
54  false
55  );
56  }
57  else
58  {
60  << "injectionMethod must be either 'point' or 'disc'"
61  << exit(FatalError);
62  }
63 }
64 
65 
66 template<class CloudType>
68 {
69  word flowType = this->coeffDict().lookup("flowType");
70  if (flowType == "constantVelocity")
71  {
72  this->coeffDict().lookup("UMag") >> UMag_;
73  flowType_ = ftConstantVelocity;
74  }
75  else if (flowType == "pressureDrivenVelocity")
76  {
77  Pinj_.reset(this->coeffDict());
78  flowType_ = ftPressureDrivenVelocity;
79  }
80  else if (flowType == "flowRateAndDischarge")
81  {
82  Cd_.reset(this->coeffDict());
83  flowType_ = ftFlowRateAndDischarge;
84  }
85  else
86  {
88  << "flowType must be either 'constantVelocity', "
89  <<"'pressureDrivenVelocity' or 'flowRateAndDischarge'"
90  << exit(FatalError);
91  }
92 }
93 
94 
95 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
96 
97 template<class CloudType>
99 (
100  const dictionary& dict,
101  CloudType& owner,
102  const word& modelName
103 )
104 :
105  InjectionModel<CloudType>(dict, owner, modelName, typeName),
106  injectionMethod_(imPoint),
107  flowType_(ftConstantVelocity),
108  outerDiameter_(readScalar(this->coeffDict().lookup("outerDiameter"))),
109  innerDiameter_(readScalar(this->coeffDict().lookup("innerDiameter"))),
110  duration_(readScalar(this->coeffDict().lookup("duration"))),
111  position_(this->coeffDict().lookup("position")),
112  injectorCell_(-1),
113  tetFaceI_(-1),
114  tetPtI_(-1),
115  direction_(this->coeffDict().lookup("direction")),
116  parcelsPerSecond_
117  (
118  readScalar(this->coeffDict().lookup("parcelsPerSecond"))
119  ),
120  flowRateProfile_
121  (
123  (
124  owner.db().time(),
125  "flowRateProfile",
126  this->coeffDict()
127  )
128  ),
129  thetaInner_
130  (
132  (
133  owner.db().time(),
134  "thetaInner",
135  this->coeffDict()
136  )
137  ),
138  thetaOuter_
139  (
141  (
142  owner.db().time(),
143  "thetaOuter",
144  this->coeffDict()
145  )
146  ),
147  sizeDistribution_
148  (
150  (
151  this->coeffDict().subDict("sizeDistribution"),
152  owner.rndGen()
153  )
154  ),
155  tanVec1_(vector::zero),
156  tanVec2_(vector::zero),
157  normal_(vector::zero),
158 
159  UMag_(0.0),
160  Cd_(owner.db().time(), "Cd"),
161  Pinj_(owner.db().time(), "Pinj")
162 {
163  if (innerDiameter_ >= outerDiameter_)
164  {
166  << "Inner diameter must be less than the outer diameter:" << nl
167  << " innerDiameter: " << innerDiameter_ << nl
168  << " outerDiameter: " << outerDiameter_
169  << exit(FatalError);
170  }
171 
172  duration_ = owner.db().time().userTimeToTime(duration_);
173 
174  setInjectionMethod();
175 
176  setFlowType();
177 
178  cachedRandom& rndGen = this->owner().rndGen();
179 
180  // Normalise direction vector
181  direction_ /= mag(direction_);
182 
183  // Determine direction vectors tangential to direction
184  vector tangent = vector::zero;
185  scalar magTangent = 0.0;
186 
187  while(magTangent < SMALL)
188  {
189  vector v = rndGen.sample01<vector>();
190 
191  tangent = v - (v & direction_)*direction_;
192  magTangent = mag(tangent);
193  }
194 
195  tanVec1_ = tangent/magTangent;
196  tanVec2_ = direction_^tanVec1_;
197 
198  // Set total volume to inject
199  this->volumeTotal_ = flowRateProfile_.integrate(0.0, duration_);
200 
201  updateMesh();
202 }
203 
204 
205 template<class CloudType>
207 (
209 )
210 :
212  injectionMethod_(im.injectionMethod_),
213  flowType_(im.flowType_),
214  outerDiameter_(im.outerDiameter_),
215  innerDiameter_(im.innerDiameter_),
216  duration_(im.duration_),
217  position_(im.position_),
218  injectorCell_(im.injectorCell_),
219  tetFaceI_(im.tetFaceI_),
220  tetPtI_(im.tetPtI_),
221  direction_(im.direction_),
222  parcelsPerSecond_(im.parcelsPerSecond_),
223  flowRateProfile_(im.flowRateProfile_),
224  thetaInner_(im.thetaInner_),
225  thetaOuter_(im.thetaOuter_),
226  sizeDistribution_(im.sizeDistribution_, false),
227  tanVec1_(im.tanVec1_),
228  tanVec2_(im.tanVec1_),
229  normal_(im.normal_),
230  UMag_(im.UMag_),
231  Cd_(im.Cd_),
232  Pinj_(im.Pinj_)
233 {}
234 
235 
236 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
237 
238 template<class CloudType>
240 {}
241 
242 
243 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
244 
245 template<class CloudType>
247 {
248  // Set/cache the injector cells
249  switch (injectionMethod_)
250  {
251  case imPoint:
252  {
253  this->findCellAtPosition
254  (
255  injectorCell_,
256  tetFaceI_,
257  tetPtI_,
258  position_
259  );
260  }
261  default:
262  {
263  // Do nothing
264  }
265  }
266 }
267 
268 
269 template<class CloudType>
271 {
272  return this->SOI_ + duration_;
273 }
274 
275 
276 template<class CloudType>
278 (
279  const scalar time0,
280  const scalar time1
281 )
282 {
283  if ((time0 >= 0.0) && (time0 < duration_))
284  {
285  return floor((time1 - time0)*parcelsPerSecond_);
286  }
287  else
288  {
289  return 0;
290  }
291 }
292 
293 
294 template<class CloudType>
296 (
297  const scalar time0,
298  const scalar time1
299 )
300 {
301  if ((time0 >= 0.0) && (time0 < duration_))
302  {
303  return flowRateProfile_.integrate(time0, time1);
304  }
305  else
306  {
307  return 0.0;
308  }
309 }
310 
311 
312 template<class CloudType>
314 (
315  const label,
316  const label,
317  const scalar,
318  vector& position,
319  label& cellOwner,
320  label& tetFaceI,
321  label& tetPtI
322 )
323 {
324  cachedRandom& rndGen = this->owner().rndGen();
325 
326  scalar beta = mathematical::twoPi*rndGen.sample01<scalar>();
327  normal_ = tanVec1_*cos(beta) + tanVec2_*sin(beta);
328 
329  switch (injectionMethod_)
330  {
331  case imPoint:
332  {
333  position = position_;
334  cellOwner = injectorCell_;
335  tetFaceI = tetFaceI_;
336  tetPtI = tetPtI_;
337 
338  break;
339  }
340  case imDisc:
341  {
342  scalar frac = rndGen.globalSample01<scalar>();
343  scalar dr = outerDiameter_ - innerDiameter_;
344  scalar r = 0.5*(innerDiameter_ + frac*dr);
345 
346  position = position_ + r*normal_;
347 
348  this->findCellAtPosition
349  (
350  cellOwner,
351  tetFaceI,
352  tetPtI,
353  position
354  );
355  break;
356  }
357  default:
358  {
360  << "Unknown injectionMethod type" << nl
361  << exit(FatalError);
362  }
363  }
364 }
365 
366 
367 template<class CloudType>
369 (
370  const label parcelI,
371  const label,
372  const scalar time,
373  typename CloudType::parcelType& parcel
374 )
375 {
376  cachedRandom& rndGen = this->owner().rndGen();
377 
378  // Set particle velocity
379  const scalar deg2Rad = mathematical::pi/180.0;
380 
381  scalar t = time - this->SOI_;
382  scalar ti = thetaInner_.value(t);
383  scalar to = thetaOuter_.value(t);
384  scalar coneAngle = rndGen.sample01<scalar>()*(to - ti) + ti;
385 
386  coneAngle *= deg2Rad;
387  scalar alpha = sin(coneAngle);
388  scalar dcorr = cos(coneAngle);
389 
390  vector normal = alpha*normal_;
391  vector dirVec = dcorr*direction_;
392  dirVec += normal;
393  dirVec /= mag(dirVec);
394 
395  switch (flowType_)
396  {
397  case ftConstantVelocity:
398  {
399  parcel.U() = UMag_*dirVec;
400  break;
401  }
402  case ftPressureDrivenVelocity:
403  {
404  scalar pAmbient = this->owner().pAmbient();
405  scalar rho = parcel.rho();
406  scalar UMag = ::sqrt(2.0*(Pinj_.value(t) - pAmbient)/rho);
407  parcel.U() = UMag*dirVec;
408  break;
409  }
410  case ftFlowRateAndDischarge:
411  {
412  scalar Ao = 0.25*mathematical::pi*outerDiameter_*outerDiameter_;
413  scalar Ai = 0.25*mathematical::pi*innerDiameter_*innerDiameter_;
414  scalar massFlowRate =
415  this->massTotal()
416  *flowRateProfile_.value(t)
417  /this->volumeTotal();
418 
419  scalar Umag = massFlowRate/(parcel.rho()*Cd_.value(t)*(Ao - Ai));
420  parcel.U() = Umag*dirVec;
421  break;
422  }
423  default:
424  {
425  }
426  }
427 
428  // Set particle diameter
429  parcel.d() = sizeDistribution_->sample();
430 }
431 
432 
433 template<class CloudType>
435 {
436  return false;
437 }
438 
439 
440 template<class CloudType>
442 {
443  return true;
444 }
445 
446 
447 // ************************************************************************* //
beta
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
mathematicalConstants.H
Foam::ConeNozzleInjection::flowRateProfile_
const TimeDataEntry< scalar > flowRateProfile_
Flow rate profile relative to SOI [].
Definition: ConeNozzleInjection.H:133
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::ConeNozzleInjection::UMag_
scalar UMag_
Constant velocity [m/s].
Definition: ConeNozzleInjection.H:160
Foam::ConeNozzleInjection::position_
vector position_
Injector position [m].
Definition: ConeNozzleInjection.H:115
Foam::ConeNozzleInjection::parcelsToInject
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
Definition: ConeNozzleInjection.C:278
Foam::ConeNozzleInjection::direction_
vector direction_
Injector direction [].
Definition: ConeNozzleInjection.H:127
Foam::constant
Collection of constants.
Definition: atomicConstants.C:37
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:255
ConeNozzleInjection.H
Foam::ConeNozzleInjection::normal_
vector normal_
Injection vector orthogonal to direction.
Definition: ConeNozzleInjection.H:154
Foam::ConeNozzleInjection::validInjection
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
Definition: ConeNozzleInjection.C:441
Foam::InjectionModel
Templated injection model class.
Definition: InjectionModel.H:67
Foam::distributionModels::distributionModel::New
static autoPtr< distributionModel > New(const dictionary &dict, cachedRandom &rndGen)
Selector.
Definition: distributionModelNew.C:32
Foam::constant::mathematical::twoPi
const scalar twoPi(2 *pi)
Foam::ConeNozzleInjection::tetFaceI_
label tetFaceI_
Index of tet face for injector cell.
Definition: ConeNozzleInjection.H:121
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::ConeNozzleInjection::injectorCell_
label injectorCell_
Cell containing injector position [].
Definition: ConeNozzleInjection.H:118
Foam::ConeNozzleInjection::tanVec1_
vector tanVec1_
First tangential vector.
Definition: ConeNozzleInjection.H:148
Foam::TimeDataEntry< scalar >
Foam::DSMCCloud::rndGen
Random & rndGen()
Return refernce to the random object.
Definition: DSMCCloudI.H:121
Foam::ConeNozzleInjection::setInjectionMethod
void setInjectionMethod()
Set the injection type.
Definition: ConeNozzleInjection.C:36
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::ConeNozzleInjection::volumeToInject
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
Definition: ConeNozzleInjection.C:296
Foam::cachedRandom
Random number generator.
Definition: cachedRandom.H:63
Foam::ConeNozzleInjection::setFlowType
void setFlowType()
Set the injection flow type.
Definition: ConeNozzleInjection.C:67
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::ConeNozzleInjection::parcelsPerSecond_
const label parcelsPerSecond_
Number of parcels to introduce per second [].
Definition: ConeNozzleInjection.H:130
Foam::ConeNozzleInjection::ConeNozzleInjection
ConeNozzleInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
Definition: ConeNozzleInjection.C:99
Foam::ConeNozzleInjection::Cd_
TimeDataEntry< scalar > Cd_
Discharge coefficient, relative to SOI [m/s].
Definition: ConeNozzleInjection.H:163
Foam::ConeNozzleInjection::flowType_
flowType flowType_
Flow type.
Definition: ConeNozzleInjection.H:103
Foam::ConeNozzleInjection
Cone injection.
Definition: ConeNozzleInjection.H:73
Foam::ConeNozzleInjection::updateMesh
virtual void updateMesh()
Set injector locations when mesh is updated.
Definition: ConeNozzleInjection.C:246
Foam::ConeNozzleInjection::thetaInner_
const TimeDataEntry< scalar > thetaInner_
Inner half-cone angle relative to SOI [deg].
Definition: ConeNozzleInjection.H:136
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
TimeDataEntry.H
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::ConeNozzleInjection::Pinj_
TimeDataEntry< scalar > Pinj_
Injection pressure [Pa].
Definition: ConeNozzleInjection.H:166
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
rho
rho
Definition: pEqn.H:3
Foam::ConeNozzleInjection::injectionMethod_
injectionMethod injectionMethod_
Point/disc injection method.
Definition: ConeNozzleInjection.H:100
Foam::ConeNozzleInjection::innerDiameter_
const scalar innerDiameter_
Inner nozzle diameter [m].
Definition: ConeNozzleInjection.H:109
Foam::ConeNozzleInjection::timeEnd
scalar timeEnd() const
Return the end-of-injection time.
Definition: ConeNozzleInjection.C:270
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::readScalar
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
Foam::Vector< scalar >
Foam::ConeNozzleInjection::tetPtI_
label tetPtI_
Index of tet point for injector cell.
Definition: ConeNozzleInjection.H:124
Foam::ConeNozzleInjection::flowType
flowType
Flow type enumeration.
Definition: ConeNozzleInjection.H:87
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
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:314
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:369
Foam::constant::mathematical::pi
const scalar pi(M_PI)
Foam::ConeNozzleInjection::fullyDescribed
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
Definition: ConeNozzleInjection.C:434
Foam::ConeNozzleInjection::thetaOuter_
const TimeDataEntry< scalar > thetaOuter_
Outer half-cone angle relative to SOI [deg].
Definition: ConeNozzleInjection.H:139
Foam::DSMCCloud::parcelType
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:217
Foam::ConeNozzleInjection::duration_
scalar duration_
Injection duration [s].
Definition: ConeNozzleInjection.H:112
rndGen
cachedRandom rndGen(label(0), -1)
Foam::ConeNozzleInjection::sizeDistribution_
const autoPtr< distributionModels::distributionModel > sizeDistribution_
Parcel size PDF model.
Definition: ConeNozzleInjection.H:142
Foam::ConeNozzleInjection::~ConeNozzleInjection
virtual ~ConeNozzleInjection()
Destructor.
Definition: ConeNozzleInjection.C:239
Foam::ConeNozzleInjection::outerDiameter_
const scalar outerDiameter_
Outer nozzle diameter [m].
Definition: ConeNozzleInjection.H:106
normal
A normal distribution model.
distributionModel.H
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:256