PatchInjection.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-2016 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 "PatchInjection.H"
27 #include "TimeFunction1.H"
28 #include "distributionModel.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class CloudType>
34 (
35  const dictionary& dict,
36  CloudType& owner,
37  const word& modelName
38 )
39 :
40  InjectionModel<CloudType>(dict, owner, modelName, typeName),
41  patchInjectionBase(owner.mesh(), this->coeffDict().lookup("patchName")),
42  duration_(readScalar(this->coeffDict().lookup("duration"))),
43  parcelsPerSecond_
44  (
45  readScalar(this->coeffDict().lookup("parcelsPerSecond"))
46  ),
47  U0_(this->coeffDict().lookup("U0")),
48  flowRateProfile_
49  (
51  (
52  owner.db().time(),
53  "flowRateProfile",
54  this->coeffDict()
55  )
56  ),
57  sizeDistribution_
58  (
60  (
61  this->coeffDict().subDict("sizeDistribution"),
62  owner.rndGen()
63  )
64  )
65 {
66  duration_ = owner.db().time().userTimeToTime(duration_);
67 
68  patchInjectionBase::updateMesh(owner.mesh());
69 
70  // Set total volume/mass to inject
71  this->volumeTotal_ = flowRateProfile_.integrate(0.0, duration_);
72 }
73 
74 
75 template<class CloudType>
77 (
79 )
80 :
83  duration_(im.duration_),
84  parcelsPerSecond_(im.parcelsPerSecond_),
85  U0_(im.U0_),
86  flowRateProfile_(im.flowRateProfile_),
87  sizeDistribution_(im.sizeDistribution_().clone().ptr())
88 {}
89 
90 
91 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
92 
93 template<class CloudType>
95 {}
96 
97 
98 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
99 
100 template<class CloudType>
102 {
103  patchInjectionBase::updateMesh(this->owner().mesh());
104 }
105 
106 
107 template<class CloudType>
109 {
110  return this->SOI_ + duration_;
111 }
112 
113 
114 template<class CloudType>
116 (
117  const scalar time0,
118  const scalar time1
119 )
120 {
121  if ((time0 >= 0.0) && (time0 < duration_))
122  {
123  scalar nParcels = (time1 - time0)*parcelsPerSecond_;
124 
125  cachedRandom& rnd = this->owner().rndGen();
126 
127  label nParcelsToInject = floor(nParcels);
128 
129  // Inject an additional parcel with a probability based on the
130  // remainder after the floor function
131  if
132  (
133  nParcelsToInject > 0
134  && (
135  nParcels - scalar(nParcelsToInject)
136  > rnd.globalPosition(scalar(0), scalar(1))
137  )
138  )
139  {
140  ++nParcelsToInject;
141  }
142 
143  return nParcelsToInject;
144  }
145  else
146  {
147  return 0;
148  }
149 }
150 
151 
152 template<class CloudType>
154 (
155  const scalar time0,
156  const scalar time1
157 )
158 {
159  if ((time0 >= 0.0) && (time0 < duration_))
160  {
161  return flowRateProfile_.integrate(time0, time1);
162  }
163  else
164  {
165  return 0.0;
166  }
167 }
168 
169 
170 template<class CloudType>
172 (
173  const label,
174  const label,
175  const scalar,
176  vector& position,
177  label& cellOwner,
178  label& tetFacei,
179  label& tetPtI
180 )
181 {
182  patchInjectionBase::setPositionAndCell
183  (
184  this->owner().mesh(),
185  this->owner().rndGen(),
186  position,
187  cellOwner,
188  tetFacei,
189  tetPtI
190  );
191 }
192 
193 
194 template<class CloudType>
196 (
197  const label,
198  const label,
199  const scalar,
200  typename CloudType::parcelType& parcel
201 )
202 {
203  // set particle velocity
204  parcel.U() = U0_;
205 
206  // set particle diameter
207  parcel.d() = sizeDistribution_->sample();
208 }
209 
210 
211 template<class CloudType>
213 {
214  return false;
215 }
216 
217 
218 template<class CloudType>
220 {
221  return true;
222 }
223 
224 
225 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
cachedRandom rndGen(label(0),-1)
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
dictionary dict
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
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Templated injection model class.
Type globalPosition(const Type &start, const Type &end)
Return a sample between start and end.
Random number generator.
Definition: cachedRandom.H:63
virtual void updateMesh()
Set injector locations when mesh is updated.
virtual scalar timeEnd() const
Return the end-of-injection time.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
virtual ~PatchInjection()
Destructor.
PatchInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:217
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
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.
Patch injection.
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68