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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2015-2021 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 "PatchInjection.H"
30 #include "distributionModel.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 template<class CloudType>
37 ({
38  { vtFixedValue, "fixedValue" },
39  { vtPatchValue, "patchValue" },
40  { vtZeroGradient, "zeroGradient" },
41 });
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
45 template<class CloudType>
47 (
48  const dictionary& dict,
49  CloudType& owner,
50  const word& modelName
51 )
52 :
53  InjectionModel<CloudType>(dict, owner, modelName, typeName),
54  patchInjectionBase(owner.mesh(), this->coeffDict().getWord("patch")),
55  duration_(this->coeffDict().getScalar("duration")),
56  parcelsPerSecond_
57  (
58  this->coeffDict().getScalar("parcelsPerSecond")
59  ),
60  velocityType_
61  (
62  velocityTypeNames_.getOrDefault
63  (
64  "velocityType",
65  this->coeffDict(),
66  vtFixedValue
67  )
68  ),
69  U0_
70  (
71  velocityType_ == vtFixedValue
72  ? this->coeffDict().template get<vector>("U0")
73  : Zero
74  ),
75  flowRateProfile_
76  (
78  (
79  "flowRateProfile",
80  this->coeffDict(),
81  &owner.mesh()
82  )
83  ),
84  sizeDistribution_
85  (
87  (
88  this->coeffDict().subDict("sizeDistribution"),
89  owner.rndGen()
90  )
91  ),
92  currentParceli_(-1),
93  currentFacei_(-1)
94 {
95  // Convert from user time to reduce the number of time conversion calls
96  const Time& time = owner.db().time();
97  duration_ = time.userTimeToTime(duration_);
98  flowRateProfile_->userTimeToTime(time);
99 
100  patchInjectionBase::updateMesh(owner.mesh());
101 
102  // Set total volume/mass to inject
103  this->volumeTotal_ = flowRateProfile_->integrate(0.0, duration_);
104 }
105 
106 
107 template<class CloudType>
109 (
110  const PatchInjection<CloudType>& im
111 )
112 :
114  patchInjectionBase(im),
115  duration_(im.duration_),
116  parcelsPerSecond_(im.parcelsPerSecond_),
117  velocityType_(im.velocityType_),
118  U0_(im.U0_),
119  flowRateProfile_(im.flowRateProfile_.clone()),
120  sizeDistribution_(im.sizeDistribution_.clone()),
121  currentParceli_(im.currentParceli_),
122  currentFacei_(im.currentFacei_)
123 {}
124 
125 
126 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
127 
128 template<class CloudType>
130 {
131  patchInjectionBase::updateMesh(this->owner().mesh());
132 }
133 
134 
135 template<class CloudType>
137 {
138  return this->SOI_ + duration_;
139 }
140 
141 
142 template<class CloudType>
144 (
145  const scalar time0,
146  const scalar time1
147 )
148 {
149  if ((time0 >= 0.0) && (time0 < duration_))
150  {
151  scalar nParcels = (time1 - time0)*parcelsPerSecond_;
152  Random& rnd = this->owner().rndGen();
153  scalar rndPos = rnd.globalPosition(scalar(0), scalar(1));
154  label nParcelsToInject = floor(nParcels);
155 
156  // Inject an additional parcel with a probability based on the
157  // remainder after the floor function
158  if
159  (
160  nParcelsToInject > 0
161  && (nParcels - scalar(nParcelsToInject) > rndPos)
162  )
163  {
164  ++nParcelsToInject;
165  }
166 
167  return nParcelsToInject;
168  }
169 
170  return 0;
171 }
172 
173 
174 template<class CloudType>
176 (
177  const scalar time0,
178  const scalar time1
179 )
180 {
181  if ((time0 >= 0.0) && (time0 < duration_))
182  {
183  return flowRateProfile_->integrate(time0, time1);
184  }
185 
186  return 0.0;
187 }
188 
189 
190 template<class CloudType>
192 (
193  const label parcelI,
194  const label nParcels,
195  const scalar time,
196  vector& position,
197  label& cellOwner,
198  label& tetFacei,
199  label& tetPti
200 )
201 {
202  currentParceli_ = parcelI;
203 
204  currentFacei_ = patchInjectionBase::setPositionAndCell
205  (
206  this->owner().mesh(),
207  this->owner().rndGen(),
208  position,
209  cellOwner,
210  tetFacei,
211  tetPti
212  );
213 }
214 
215 
216 template<class CloudType>
218 (
219  const label parcelI,
220  const label nParcels,
221  const scalar time,
222  typename CloudType::parcelType& parcel
223 )
224 {
225  // Set particle velocity
226  switch (velocityType_)
227  {
228  case vtFixedValue:
229  {
230  parcel.U() = U0_;
231  break;
232  }
233  case vtPatchValue:
234  {
235  if (parcelI != currentParceli_)
236  {
238  << "Synchronisation problem: "
239  << "attempting to set injected parcel " << parcelI
240  << " properties using cached parcel " << currentParceli_
241  << " properties" << endl;
242  }
243 
244  const label patchFacei = currentFacei_;
245 
246  if (patchFacei < 0)
247  {
249  << "Unable to set parcel velocity using patch value "
250  << "due to missing face index: patchFacei=" << patchFacei
251  << abort(FatalError);
252  }
253 
254  const volVectorField& U = this->owner().U();
255  const label patchi = this->patchId_;
256  parcel.U() = U.boundaryField()[patchi][patchFacei];
257  break;
258  }
259  case vtZeroGradient:
260  {
261  const label celli = parcel.cell();
262 
263  if (celli < 0)
264  {
266  << "Unable to set parcel velocity using zeroGradient "
267  << "due to missing cell index"
268  << abort(FatalError);
269  }
270 
271  const volVectorField& U = this->owner().U();
272  parcel.U() = U[celli];
273  break;
274  }
275  default:
276  {
278  << "Unhandled velocityType "
279  << velocityTypeNames_[velocityType_]
280  << ". Available options are:"
281  << velocityTypeNames_.sortedToc()
282  << abort(FatalError);
283  }
284  }
285 
286  // Set particle diameter
287  parcel.d() = sizeDistribution_->sample();
288 }
289 
290 
291 template<class CloudType>
293 {
294  return false;
295 }
296 
297 
298 template<class CloudType>
300 {
301  return true;
302 }
303 
304 
305 // ************************************************************************* //
Foam::PatchInjection::fullyDescribed
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
Definition: PatchInjection.C:292
Foam::Random
Random number generator.
Definition: Random.H:59
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:57
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::PatchInjection::parcelsToInject
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
Definition: PatchInjection.C:144
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::PatchInjection::validInjection
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
Definition: PatchInjection.C:299
Foam::TimeState::userTimeToTime
virtual scalar userTimeToTime(const scalar theta) const
Convert the user-time (e.g. CA deg) to real-time (s).
Definition: TimeState.C:49
Foam::InjectionModel
Templated injection model class.
Definition: InjectionModel.H:73
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
PatchInjection.H
Foam::Function1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: propellerInfo.H:291
Foam::DSMCCloud::rndGen
Random & rndGen()
Return reference to the random object.
Definition: DSMCCloudI.H:124
Foam::PatchInjection::PatchInjection
PatchInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
Definition: PatchInjection.C:47
Foam::DSMCCloud::mesh
const fvMesh & mesh() const
Return reference to the mesh.
Definition: DSMCCloudI.H:44
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:192
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:71
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::PatchInjection::velocityTypeNames_
static const Enum< velocityType > velocityTypeNames_
Velocity type names.
Definition: PatchInjection.H:87
Foam::dictionary::getWord
word getWord(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Same as get< word >(const word&, keyType::option)
Definition: dictionary.H:1514
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::PatchInjection
Patch injection.
Definition: PatchInjection.H:69
Foam::PatchInjection::timeEnd
virtual scalar timeEnd() const
Return the end-of-injection time.
Definition: PatchInjection.C:136
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
U
U
Definition: pEqn.H:72
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::Vector< scalar >
Foam::patchInjectionBase
Definition: patchInjectionBase.H:64
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:218
rndGen
Random rndGen
Definition: createFields.H:23
Foam::PatchInjection::volumeToInject
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
Definition: PatchInjection.C:176
Foam::DSMCCloud::parcelType
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
Foam::GeometricField< vector, fvPatchField, volMesh >
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::PatchInjection::updateMesh
virtual void updateMesh()
Set injector locations when mesh is updated.
Definition: PatchInjection.C:129
Foam::Random::globalPosition
Type globalPosition(const Type &start, const Type &end)
Return a sample on the interval [start,end].
Definition: RandomTemplates.C:126
distributionModel.H