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-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 "PatchInjection.H"
27 #include "TimeDataEntry.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_, false)
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  cachedRandom& rnd = this->owner().rndGen();
125  scalar rndPos = rnd.globalPosition(scalar(0), scalar(1));
126  label nParcelsToInject = floor(nParcels);
127 
128  // Inject an additional parcel with a probability based on the
129  // remainder after the floor function
130  if
131  (
132  nParcelsToInject > 0
133  && (nParcels - scalar(nParcelsToInject) > rndPos)
134  )
135  {
136  ++nParcelsToInject;
137  }
138 
139  return nParcelsToInject;
140  }
141  else
142  {
143  return 0;
144  }
145 }
146 
147 
148 template<class CloudType>
150 (
151  const scalar time0,
152  const scalar time1
153 )
154 {
155  if ((time0 >= 0.0) && (time0 < duration_))
156  {
157  return flowRateProfile_.integrate(time0, time1);
158  }
159  else
160  {
161  return 0.0;
162  }
163 }
164 
165 
166 template<class CloudType>
168 (
169  const label,
170  const label,
171  const scalar,
172  vector& position,
173  label& cellOwner,
174  label& tetFaceI,
175  label& tetPtI
176 )
177 {
178  patchInjectionBase::setPositionAndCell
179  (
180  this->owner().mesh(),
181  this->owner().rndGen(),
182  position,
183  cellOwner,
184  tetFaceI,
185  tetPtI
186  );
187 }
188 
189 
190 template<class CloudType>
192 (
193  const label,
194  const label,
195  const scalar,
196  typename CloudType::parcelType& parcel
197 )
198 {
199  // set particle velocity
200  parcel.U() = U0_;
201 
202  // set particle diameter
203  parcel.d() = sizeDistribution_->sample();
204 }
205 
206 
207 template<class CloudType>
209 {
210  return false;
211 }
212 
213 
214 template<class CloudType>
216 {
217  return true;
218 }
219 
220 
221 // ************************************************************************* //
Foam::PatchInjection::fullyDescribed
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
Definition: PatchInjection.C:208
Foam::PatchInjection::U0_
const vector U0_
Initial parcel velocity [m/s].
Definition: PatchInjection.H:80
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::PatchInjection::parcelsToInject
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
Definition: PatchInjection.C:116
Foam::PatchInjection::duration_
scalar duration_
Injection duration [s].
Definition: PatchInjection.H:74
Foam::PatchInjection::validInjection
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
Definition: PatchInjection.C:215
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::InjectionModel
Templated injection model class.
Definition: InjectionModel.H:67
Foam::cachedRandom::globalPosition
Type globalPosition(const Type &start, const Type &end)
Return a sample between start and end.
Definition: cachedRandomTemplates.C:111
PatchInjection.H
Foam::TimeDataEntry< scalar >
Foam::DSMCCloud::rndGen
Random & rndGen()
Return refernce to the random object.
Definition: DSMCCloudI.H:121
Foam::PatchInjection::PatchInjection
PatchInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
Definition: PatchInjection.C:34
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::DSMCCloud::mesh
const fvMesh & mesh() const
Return refernce to the mesh.
Definition: DSMCCloudI.H:41
Foam::PatchInjection::parcelsPerSecond_
const label parcelsPerSecond_
Number of parcels to introduce per second [].
Definition: PatchInjection.H:77
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
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::PatchInjection::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: PatchInjection.C:168
Foam::PatchInjection
Patch injection.
Definition: PatchInjection.H:66
readScalar
#define readScalar
Definition: doubleScalar.C:38
Foam::PatchInjection::timeEnd
virtual scalar timeEnd() const
Return the end-of-injection time.
Definition: PatchInjection.C:108
Foam::PatchInjection::flowRateProfile_
const TimeDataEntry< scalar > flowRateProfile_
Flow rate profile relative to SOI [].
Definition: PatchInjection.H:83
Foam::PatchInjection::sizeDistribution_
const autoPtr< distributionModels::distributionModel > sizeDistribution_
Parcel size distribution model.
Definition: PatchInjection.H:86
Foam::Vector< scalar >
Foam::patchInjectionBase
Definition: patchInjectionBase.H:61
Foam::PatchInjection::~PatchInjection
virtual ~PatchInjection()
Destructor.
Definition: PatchInjection.C:94
Foam::PatchInjection::setProperties
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
Definition: PatchInjection.C:192
Foam::PatchInjection::volumeToInject
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
Definition: PatchInjection.C:150
Foam::DSMCCloud::parcelType
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:217
rndGen
cachedRandom rndGen(label(0), -1)
Foam::PatchInjection::updateMesh
virtual void updateMesh()
Set injector locations when mesh is updated.
Definition: PatchInjection.C:101
distributionModel.H
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress