ConeInjection.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
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 "TimeFunction1.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  tanVec1_[i] = normalised(perpendicular(axis));
111  tanVec2_[i] = normalised(axis^tanVec1_[i]);
112  }
113 
114  // Set total volume to inject
115  this->volumeTotal_ = flowRateProfile_.integrate(0.0, duration_);
116 
117  updateMesh();
118 }
119 
120 
121 template<class CloudType>
123 (
124  const ConeInjection<CloudType>& im
125 )
126 :
128  positionAxis_(im.positionAxis_),
129  injectorCells_(im.injectorCells_),
130  injectorTetFaces_(im.injectorTetFaces_),
131  injectorTetPts_(im.injectorTetPts_),
132  duration_(im.duration_),
133  parcelsPerInjector_(im.parcelsPerInjector_),
134  flowRateProfile_(im.flowRateProfile_),
135  Umag_(im.Umag_),
136  thetaInner_(im.thetaInner_),
137  thetaOuter_(im.thetaOuter_),
138  sizeDistribution_(im.sizeDistribution_().clone().ptr()),
139  nInjected_(im.nInjected_),
140  tanVec1_(im.tanVec1_),
141  tanVec2_(im.tanVec2_)
142 {}
143 
144 
145 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
146 
147 template<class CloudType>
149 {}
150 
151 
152 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
153 
154 template<class CloudType>
156 {
157  // Set/cache the injector cells
158  forAll(positionAxis_, i)
159  {
160  this->findCellAtPosition
161  (
162  injectorCells_[i],
163  injectorTetFaces_[i],
164  injectorTetPts_[i],
165  positionAxis_[i].first()
166  );
167  }
168 }
169 
170 
171 template<class CloudType>
173 {
174  return this->SOI_ + duration_;
175 }
176 
177 
178 template<class CloudType>
180 (
181  const scalar time0,
182  const scalar time1
183 )
184 {
185  if ((time0 >= 0.0) && (time0 < duration_))
186  {
187  const scalar targetVolume = flowRateProfile_.integrate(0, time1);
188 
189  const label targetParcels =
190  parcelsPerInjector_*targetVolume/this->volumeTotal_;
191 
192  const label nToInject = targetParcels - nInjected_;
193 
194  nInjected_ += nToInject;
195 
196  return positionAxis_.size()*nToInject;
197  }
198  else
199  {
200  return 0;
201  }
202 }
203 
204 
205 template<class CloudType>
207 (
208  const scalar time0,
209  const scalar time1
210 )
211 {
212  if ((time0 >= 0.0) && (time0 < duration_))
213  {
214  return flowRateProfile_.integrate(time0, time1);
215  }
216  else
217  {
218  return 0.0;
219  }
220 }
221 
222 
223 template<class CloudType>
225 (
226  const label parcelI,
227  const label,
228  const scalar,
229  vector& position,
230  label& cellOwner,
231  label& tetFacei,
232  label& tetPti
233 )
234 {
235  const label i = parcelI % positionAxis_.size();
236 
237  position = positionAxis_[i].first();
238  cellOwner = injectorCells_[i];
239  tetFacei = injectorTetFaces_[i];
240  tetPti = injectorTetPts_[i];
241 }
242 
243 
244 template<class CloudType>
246 (
247  const label parcelI,
248  const label,
249  const scalar time,
250  typename CloudType::parcelType& parcel
251 )
252 {
253  Random& rnd = this->owner().rndGen();
254 
255  // set particle velocity
256  const label i = parcelI % positionAxis_.size();
257 
258  scalar t = time - this->SOI_;
259  scalar ti = thetaInner_.value(t);
260  scalar to = thetaOuter_.value(t);
261  scalar coneAngle = degToRad(rnd.scalarAB(ti, to));
262 
263  scalar alpha = sin(coneAngle);
264  scalar dcorr = cos(coneAngle);
265  scalar beta = twoPi*rnd.sample01<scalar>();
266 
267  vector normal = alpha*(tanVec1_[i]*cos(beta) + tanVec2_[i]*sin(beta));
268  vector dirVec = dcorr*positionAxis_[i].second();
269  dirVec += normal;
270  dirVec /= mag(dirVec);
271 
272  parcel.U() = Umag_.value(t)*dirVec;
273 
274  // set particle diameter
275  parcel.d() = sizeDistribution_().sample();
276 }
277 
278 
279 template<class CloudType>
281 {
282  return false;
283 }
284 
285 
286 template<class CloudType>
288 {
289  return true;
290 }
291 
292 
293 // ************************************************************************* //
ConeInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
Definition: ConeInjection.C:37
dictionary dict
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:413
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Templated injection model class.
scalar scalarAB(const scalar a, const scalar b)
Advance the state and return a scalar sample from a uniform.
Definition: RandomI.H:63
Unit conversion functions.
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
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.
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
Type sample01()
Advance the state and return a sample of a given type from a.
Definition: RandomI.H:70
stressControl lookup("compactNormalStress") >> compactNormalStress
dimensionedScalar cos(const dimensionedScalar &ds)
mathematical constants.
A class for handling words, derived from string.
Definition: word.H:59
scalar timeEnd() const
Return the end-of-injection time.
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
Random number generator.
Definition: Random.H:57
const scalar twoPi(2 *pi)
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:218
dimensionedScalar sin(const dimensionedScalar &ds)
virtual ~ConeInjection()
Destructor.
Vector< Cmpt > perpendicular(const Vector< Cmpt > &v)
Definition: VectorI.H:166
static autoPtr< distributionModel > New(const dictionary &dict, Random &rndGen)
Selector.
virtual void updateMesh()
Set injector locations when mesh is updated.
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
dimensioned< scalar > mag(const dimensioned< Type > &)
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:69
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
Multi-point cone injection model.
Definition: ConeInjection.H:62