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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2015-2019 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 "ConeInjection.H"
30 #include "TimeFunction1.H"
31 #include "mathematicalConstants.H"
32 #include "unitConversion.H"
33 
34 using namespace Foam::constant::mathematical;
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
38 template<class CloudType>
40 (
41  const dictionary& dict,
42  CloudType& owner,
43  const word& modelName
44 )
45 :
46  InjectionModel<CloudType>(dict, owner, modelName, typeName),
47  positionAxis_(this->coeffDict().lookup("positionAxis")),
48  injectorCells_(positionAxis_.size()),
49  injectorTetFaces_(positionAxis_.size()),
50  injectorTetPts_(positionAxis_.size()),
51  duration_(this->coeffDict().getScalar("duration")),
52  parcelsPerInjector_
53  (
54  this->coeffDict().getScalar("parcelsPerInjector")
55  ),
56  flowRateProfile_
57  (
59  (
60  owner.db().time(),
61  "flowRateProfile",
62  this->coeffDict()
63  )
64  ),
65  Umag_
66  (
68  (
69  owner.db().time(),
70  "Umag",
71  this->coeffDict()
72  )
73  ),
74  thetaInner_
75  (
77  (
78  owner.db().time(),
79  "thetaInner",
80  this->coeffDict()
81  )
82  ),
83  thetaOuter_
84  (
86  (
87  owner.db().time(),
88  "thetaOuter",
89  this->coeffDict()
90  )
91  ),
92  sizeDistribution_
93  (
95  (
96  this->coeffDict().subDict("sizeDistribution"), owner.rndGen()
97  )
98  ),
99  nInjected_(this->parcelsAddedTotal()),
100  tanVec1_(),
101  tanVec2_()
102 {
103  updateMesh();
104 
105  tanVec1_.setSize(positionAxis_.size());
106  tanVec2_.setSize(positionAxis_.size());
107 
108  duration_ = owner.db().time().userTimeToTime(duration_);
109 
110  // Normalise direction vector and determine direction vectors
111  // tangential to injector axis direction
112  forAll(positionAxis_, i)
113  {
114  vector& axis = positionAxis_[i].second();
115  axis.normalise();
116 
117  vector tangent = Zero;
118  scalar magTangent = 0.0;
119 
120  Random& rnd = this->owner().rndGen();
121  while (magTangent < SMALL)
122  {
123  vector v = rnd.sample01<vector>();
124 
125  tangent = v - (v & axis)*axis;
126  magTangent = mag(tangent);
127  }
128 
129  tanVec1_[i] = tangent/magTangent;
130  tanVec2_[i] = axis^tanVec1_[i];
131  }
132 
133  // Set total volume to inject
134  this->volumeTotal_ = flowRateProfile_.integrate(0.0, duration_);
135 }
136 
137 
138 template<class CloudType>
140 (
141  const ConeInjection<CloudType>& im
142 )
143 :
145  positionAxis_(im.positionAxis_),
146  injectorCells_(im.injectorCells_),
147  injectorTetFaces_(im.injectorTetFaces_),
148  injectorTetPts_(im.injectorTetPts_),
149  duration_(im.duration_),
150  parcelsPerInjector_(im.parcelsPerInjector_),
151  flowRateProfile_(im.flowRateProfile_),
152  Umag_(im.Umag_),
153  thetaInner_(im.thetaInner_),
154  thetaOuter_(im.thetaOuter_),
155  sizeDistribution_(im.sizeDistribution_.clone()),
156  nInjected_(im.nInjected_),
157  tanVec1_(im.tanVec1_),
158  tanVec2_(im.tanVec2_)
159 {}
160 
161 
162 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
163 
164 template<class CloudType>
166 {
167  bitSet reject(positionAxis_.size());
168 
169  // Set/cache the injector cells
170  forAll(positionAxis_, i)
171  {
172  if
173  (
174  !this->findCellAtPosition
175  (
176  injectorCells_[i],
177  injectorTetFaces_[i],
178  injectorTetPts_[i],
179  positionAxis_[i].first(),
180  !this->ignoreOutOfBounds_
181 
182  )
183  )
184  {
185  reject.set(i);
186  }
187  }
188 
189  const label nRejected = reject.count();
190 
191  if (nRejected)
192  {
193  reject.flip();
194  inplaceSubset(reject, injectorCells_);
195  inplaceSubset(reject, injectorTetFaces_);
196  inplaceSubset(reject, injectorTetPts_);
197  inplaceSubset(reject, positionAxis_);
198 
199  Info<< " " << nRejected
200  << " positions rejected, out of bounds" << endl;
201  }
202 }
203 
204 
205 template<class CloudType>
207 {
208  return this->SOI_ + duration_;
209 }
210 
211 
212 template<class CloudType>
214 (
215  const scalar time0,
216  const scalar time1
217 )
218 {
219  if ((time0 >= 0.0) && (time0 < duration_))
220  {
221  const scalar targetVolume = flowRateProfile_.integrate(0, time1);
222 
223  const scalar volumeFraction = targetVolume/this->volumeTotal_;
224 
225  const label targetParcels =
226  ceil(positionAxis_.size()*parcelsPerInjector_*volumeFraction);
227 
228  return targetParcels - nInjected_;
229  }
230 
231  return 0;
232 }
233 
234 
235 template<class CloudType>
237 (
238  const scalar time0,
239  const scalar time1
240 )
241 {
242  if ((time0 >= 0.0) && (time0 < duration_))
243  {
244  return flowRateProfile_.integrate(time0, time1);
245  }
246 
247  return 0.0;
248 }
249 
250 
251 template<class CloudType>
253 (
254  const label parcelI,
255  const label,
256  const scalar,
257  vector& position,
258  label& cellOwner,
259  label& tetFacei,
260  label& tetPti
261 )
262 {
263  const label i = parcelI % positionAxis_.size();
264 
265  position = positionAxis_[i].first();
266  cellOwner = injectorCells_[i];
267  tetFacei = injectorTetFaces_[i];
268  tetPti = injectorTetPts_[i];
269 }
270 
271 
272 template<class CloudType>
274 (
275  const label parcelI,
276  const label,
277  const scalar time,
278  typename CloudType::parcelType& parcel
279 )
280 {
281  Random& rnd = this->owner().rndGen();
282 
283  // Set particle velocity
284  const label i = parcelI % positionAxis_.size();
285 
286  scalar t = time - this->SOI_;
287  scalar ti = thetaInner_.value(t);
288  scalar to = thetaOuter_.value(t);
289  scalar coneAngle = degToRad(rnd.position<scalar>(ti, to));
290 
291  scalar alpha = sin(coneAngle);
292  scalar dcorr = cos(coneAngle);
293  scalar beta = twoPi*rnd.sample01<scalar>();
294 
295  vector normal = alpha*(tanVec1_[i]*cos(beta) + tanVec2_[i]*sin(beta));
296  vector dirVec = dcorr*positionAxis_[i].second();
297  dirVec += normal;
298  dirVec.normalise();
299 
300  parcel.U() = Umag_.value(t)*dirVec;
301 
302  // Set particle diameter
303  parcel.d() = sizeDistribution_().sample();
304 
305  // Increment number of particles injected
306  nInjected_++;
307 }
308 
309 
310 template<class CloudType>
312 {
313  return false;
314 }
315 
316 
317 template<class CloudType>
319 {
320  return true;
321 }
322 
323 
324 // ************************************************************************* //
Foam::Random
Random number generator.
Definition: Random.H:59
mathematicalConstants.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:63
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
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
Foam::Random::sample01
Type sample01()
Return a sample whose components lie in the range [0,1].
Definition: RandomTemplates.C:36
Foam::ConeInjection::validInjection
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
Definition: ConeInjection.C:318
unitConversion.H
Unit conversion functions.
Foam::InjectionModel
Templated injection model class.
Definition: InjectionModel.H:73
Foam::TimeFunction1< scalar >
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::Vector::normalise
Vector< Cmpt > & normalise()
Normalise the vector by its magnitude.
Definition: VectorI.H:123
Foam::DSMCCloud::rndGen
Random & rndGen()
Return reference to the random object.
Definition: DSMCCloudI.H:124
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:253
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
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:237
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
beta
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
Foam::ConeInjection::fullyDescribed
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
Definition: ConeInjection.C:311
Foam::ConeInjection::ConeInjection
ConeInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
Definition: ConeInjection.C:40
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::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::inplaceSubset
void inplaceSubset(const BoolListType &select, ListType &input, const bool invert=false)
Inplace extract elements of the input list when select is true.
Definition: ListOpsTemplates.C:590
Foam::degToRad
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
Definition: unitConversion.H:48
Foam::distributionModel::New
static autoPtr< distributionModel > New(const dictionary &dict, Random &rndGen)
Selector.
Definition: distributionModelNew.C:34
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:214
Foam::ConeInjection
Multi-point cone injection model.
Definition: ConeInjection.H:67
Foam::Random::position
Type position(const Type &start, const Type &end)
Return a sample on the interval [start,end].
Definition: RandomTemplates.C:62
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:274
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::ConeInjection::updateMesh
virtual void updateMesh()
Set injector locations when mesh is updated.
Definition: ConeInjection.C:165
Foam::DSMCCloud::parcelType
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
TimeFunction1.H
Foam::ConeInjection::timeEnd
scalar timeEnd() const
Return the end-of-injection time.
Definition: ConeInjection.C:206
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265