ConeInjection.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-2012 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 "ConeInjection.H"
27 #include "TimeDataEntry.H"
28 #include "mathematicalConstants.H"
29 #include "unitConversion.H"
30 
31 using namespace Foam::constant::mathematical;
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
35 template<class CloudType>
37 (
38  const dictionary& dict,
39  CloudType& owner,
40  const word& modelName
41 )
42 :
43  InjectionModel<CloudType>(dict, owner, modelName, typeName),
44  positionAxis_(this->coeffDict().lookup("positionAxis")),
45  injectorCells_(positionAxis_.size()),
46  injectorTetFaces_(positionAxis_.size()),
47  injectorTetPts_(positionAxis_.size()),
48  duration_(readScalar(this->coeffDict().lookup("duration"))),
49  parcelsPerInjector_
50  (
51  readScalar(this->coeffDict().lookup("parcelsPerInjector"))
52  ),
53  flowRateProfile_
54  (
56  (
57  owner.db().time(),
58  "flowRateProfile",
59  this->coeffDict()
60  )
61  ),
62  Umag_
63  (
65  (
66  owner.db().time(),
67  "Umag",
68  this->coeffDict()
69  )
70  ),
71  thetaInner_
72  (
74  (
75  owner.db().time(),
76  "thetaInner",
77  this->coeffDict()
78  )
79  ),
80  thetaOuter_
81  (
83  (
84  owner.db().time(),
85  "thetaOuter",
86  this->coeffDict()
87  )
88  ),
89  sizeDistribution_
90  (
92  (
93  this->coeffDict().subDict("sizeDistribution"), owner.rndGen()
94  )
95  ),
96  nInjected_(this->parcelsAddedTotal()),
97  tanVec1_(positionAxis_.size()),
98  tanVec2_(positionAxis_.size())
99 {
100  duration_ = owner.db().time().userTimeToTime(duration_);
101 
102  // Normalise direction vector and determine direction vectors
103  // tangential to injector axis direction
104  forAll(positionAxis_, i)
105  {
106  vector& axis = positionAxis_[i].second();
107 
108  axis /= mag(axis);
109 
110  vector tangent = vector::zero;
111  scalar magTangent = 0.0;
112 
113  cachedRandom& rnd = this->owner().rndGen();
114  while (magTangent < SMALL)
115  {
116  vector v = rnd.sample01<vector>();
117 
118  tangent = v - (v & axis)*axis;
119  magTangent = mag(tangent);
120  }
121 
122  tanVec1_[i] = tangent/magTangent;
123  tanVec2_[i] = axis^tanVec1_[i];
124  }
125 
126  // Set total volume to inject
127  this->volumeTotal_ = flowRateProfile_.integrate(0.0, duration_);
128 
129  updateMesh();
130 }
131 
132 
133 template<class CloudType>
135 (
136  const ConeInjection<CloudType>& im
137 )
138 :
140  positionAxis_(im.positionAxis_),
141  injectorCells_(im.injectorCells_),
142  injectorTetFaces_(im.injectorTetFaces_),
143  injectorTetPts_(im.injectorTetPts_),
144  duration_(im.duration_),
145  parcelsPerInjector_(im.parcelsPerInjector_),
146  flowRateProfile_(im.flowRateProfile_),
147  Umag_(im.Umag_),
148  thetaInner_(im.thetaInner_),
149  thetaOuter_(im.thetaOuter_),
150  sizeDistribution_(im.sizeDistribution_, false),
151  nInjected_(im.nInjected_),
152  tanVec1_(im.tanVec1_),
153  tanVec2_(im.tanVec2_)
154 {}
155 
156 
157 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
158 
159 template<class CloudType>
161 {}
162 
163 
164 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
165 
166 template<class CloudType>
168 {
169  // Set/cache the injector cells
170  forAll(positionAxis_, i)
171  {
172  this->findCellAtPosition
173  (
174  injectorCells_[i],
175  injectorTetFaces_[i],
176  injectorTetPts_[i],
177  positionAxis_[i].first()
178  );
179  }
180 }
181 
182 
183 template<class CloudType>
185 {
186  return this->SOI_ + duration_;
187 }
188 
189 
190 template<class CloudType>
192 (
193  const scalar time0,
194  const scalar time1
195 )
196 {
197  if ((time0 >= 0.0) && (time0 < duration_))
198  {
199  const scalar targetVolume = flowRateProfile_.integrate(0, time1);
200 
201  const label targetParcels =
202  parcelsPerInjector_*targetVolume/this->volumeTotal_;
203 
204  const label nToInject = targetParcels - nInjected_;
205 
206  return positionAxis_.size()*nToInject;
207  }
208  else
209  {
210  return 0;
211  }
212 }
213 
214 
215 template<class CloudType>
217 (
218  const scalar time0,
219  const scalar time1
220 )
221 {
222  if ((time0 >= 0.0) && (time0 < duration_))
223  {
224  return flowRateProfile_.integrate(time0, time1);
225  }
226  else
227  {
228  return 0.0;
229  }
230 }
231 
232 
233 template<class CloudType>
235 (
236  const label parcelI,
237  const label,
238  const scalar,
239  vector& position,
240  label& cellOwner,
241  label& tetFaceI,
242  label& tetPtI
243 )
244 {
245  const label i = parcelI % positionAxis_.size();
246 
247  position = positionAxis_[i].first();
248  cellOwner = injectorCells_[i];
249  tetFaceI = injectorTetFaces_[i];
250  tetPtI = injectorTetPts_[i];
251 }
252 
253 
254 template<class CloudType>
256 (
257  const label parcelI,
258  const label,
259  const scalar time,
260  typename CloudType::parcelType& parcel
261 )
262 {
263  cachedRandom& rnd = this->owner().rndGen();
264 
265  // Set particle velocity
266  const label i = parcelI % positionAxis_.size();
267 
268  scalar t = time - this->SOI_;
269  scalar ti = thetaInner_.value(t);
270  scalar to = thetaOuter_.value(t);
271  scalar coneAngle = degToRad(rnd.position<scalar>(ti, to));
272 
273  scalar alpha = sin(coneAngle);
274  scalar dcorr = cos(coneAngle);
275  scalar beta = twoPi*rnd.sample01<scalar>();
276 
277  vector normal = alpha*(tanVec1_[i]*cos(beta) + tanVec2_[i]*sin(beta));
278  vector dirVec = dcorr*positionAxis_[i].second();
279  dirVec += normal;
280  dirVec /= mag(dirVec);
281 
282  parcel.U() = Umag_.value(t)*dirVec;
283 
284  // Set particle diameter
285  parcel.d() = sizeDistribution_().sample();
286 
287  // Increment number of particles injected
288  nInjected_++;
289 }
290 
291 
292 template<class CloudType>
294 {
295  return false;
296 }
297 
298 
299 template<class CloudType>
301 {
302  return true;
303 }
304 
305 
306 // ************************************************************************* //
beta
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
mathematicalConstants.H
Foam::ConeInjection::~ConeInjection
virtual ~ConeInjection()
Destructor.
Definition: ConeInjection.C:160
Foam::ConeInjection::parcelsPerInjector_
const label parcelsPerInjector_
Number of parcels to introduce per injector.
Definition: ConeInjection.H:84
Foam::cachedRandom::sample01
Type sample01()
Return a sample whose components lie in the range 0-1.
Definition: cachedRandomTemplates.C:32
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::ConeInjection::injectorCells_
labelList injectorCells_
List of cell labels corresponding to injector positions.
Definition: ConeInjection.H:72
Foam::ConeInjection::thetaInner_
const TimeDataEntry< scalar > thetaInner_
Inner half-cone angle relative to SOI [deg].
Definition: ConeInjection.H:93
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::ConeInjection::injectorTetFaces_
labelList injectorTetFaces_
List of tetFace labels corresponding to injector positions.
Definition: ConeInjection.H:75
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::ConeInjection::flowRateProfile_
const TimeDataEntry< scalar > flowRateProfile_
Flow rate profile relative to SOI [].
Definition: ConeInjection.H:87
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:255
Foam::ConeInjection::validInjection
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
Definition: ConeInjection.C:300
unitConversion.H
Unit conversion functions.
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::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::TimeDataEntry< scalar >
Foam::DSMCCloud::rndGen
Random & rndGen()
Return refernce to the random object.
Definition: DSMCCloudI.H:121
Foam::ConeInjection::tanVec2_
vectorList tanVec2_
Second tangential vector.
Definition: ConeInjection.H:111
Foam::ConeInjection::thetaOuter_
const TimeDataEntry< scalar > thetaOuter_
Outer half-cone angle relative to SOI [deg].
Definition: ConeInjection.H:96
ConeInjection.H
Foam::ConeInjection::volumeToInject
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
Definition: ConeInjection.C:217
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::cachedRandom
Random number generator.
Definition: cachedRandom.H:63
Foam::ConeInjection::nInjected_
label nInjected_
Number of parcels per injector already injected.
Definition: ConeInjection.H:102
Foam::cachedRandom::position
Type position(const Type &start, const Type &end)
Return a sample between start and end.
Definition: cachedRandomTemplates.C:58
Foam::ConeInjection::fullyDescribed
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
Definition: ConeInjection.C:293
Foam::ConeInjection::ConeInjection
ConeInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
Definition: ConeInjection.C:37
Foam::ConeInjection::duration_
scalar duration_
Injection duration [s].
Definition: ConeInjection.H:81
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
TimeDataEntry.H
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::ConeInjection::Umag_
const TimeDataEntry< scalar > Umag_
Parcel velocity magnitude relative to SOI [m/s].
Definition: ConeInjection.H:90
Foam::ConeInjection::sizeDistribution_
const autoPtr< distributionModels::distributionModel > sizeDistribution_
Parcel size distribution model.
Definition: ConeInjection.H:99
Foam::constant::mathematical
mathematical constants.
Foam::ConeInjection::parcelsToInject
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
Definition: ConeInjection.C:192
Foam::ConeInjection
Multi-point cone injection model.
Definition: ConeInjection.H:62
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::ConeInjection::setProperties
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
Definition: ConeInjection.C:256
Foam::ConeInjection::injectorTetPts_
labelList injectorTetPts_
List of tetPt labels corresponding to injector positions.
Definition: ConeInjection.H:78
Foam::ConeInjection::updateMesh
virtual void updateMesh()
Set injector locations when mesh is updated.
Definition: ConeInjection.C:167
Foam::DSMCCloud::parcelType
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:217
Foam::ConeInjection::tanVec1_
vectorList tanVec1_
First tangential vector.
Definition: ConeInjection.H:108
Foam::ConeInjection::timeEnd
scalar timeEnd() const
Return the end-of-injection time.
Definition: ConeInjection.C:184
normal
A normal distribution model.
Foam::ConeInjection::positionAxis_
List< Tuple2< vector, vector > > positionAxis_
List of position and axis for each injector.
Definition: ConeInjection.H:69
Foam::ConeInjection::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, tetFace and tetPt.
Definition: ConeInjection.C:235
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress
Foam::degToRad
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
Definition: unitConversion.H:45
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:256