PatchInjection.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-2023 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 "distribution.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_(this->readDuration(dict, owner)),
43  massFlowRate_(this->readMassFlowRate(dict, owner, duration_)),
44  parcelsPerSecond_(this->readParcelsPerSecond(dict, owner)),
45  U0_(this->coeffDict().lookup("U0")),
46  sizeDistribution_
47  (
49  (
50  this->coeffDict().subDict("sizeDistribution"),
51  owner.rndGen(),
52  this->sizeSampleQ()
53  )
54  )
55 {}
56 
57 
58 template<class CloudType>
60 (
62 )
63 :
66  duration_(im.duration_),
67  massFlowRate_(im.massFlowRate_),
68  parcelsPerSecond_(im.parcelsPerSecond_),
69  U0_(im.U0_),
70  sizeDistribution_(im.sizeDistribution_().clone().ptr())
71 {}
72 
73 
74 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
75 
76 template<class CloudType>
78 {}
79 
80 
81 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
82 
83 template<class CloudType>
85 {
86  patchInjectionBase::topoChange(this->owner().mesh());
87 }
88 
89 
90 template<class CloudType>
92 {
93  return this->SOI_ + duration_;
94 }
95 
96 
97 template<class CloudType>
99 (
100  const scalar time0,
101  const scalar time1
102 )
103 {
104  if (time0 >= 0 && time0 < duration_)
105  {
106  scalar nParcels = parcelsPerSecond_.integral(time0, time1);
107 
108  Random& rnd = this->owner().rndGen();
109 
110  label nParcelsToInject = floor(nParcels);
111 
112  // Inject an additional parcel with a probability based on the
113  // remainder after the floor function
114  if
115  (
116  nParcelsToInject > 0
117  && (nParcels - scalar(nParcelsToInject) > rnd.globalScalar01())
118  )
119  {
120  ++nParcelsToInject;
121  }
122 
123  return nParcelsToInject;
124  }
125  else
126  {
127  return 0;
128  }
129 }
130 
131 
132 template<class CloudType>
134 (
135  const scalar time0,
136  const scalar time1
137 )
138 {
139  if (time0 >= 0 && time0 < duration_)
140  {
141  return massFlowRate_.integral(time0, time1);
142  }
143  else
144  {
145  return 0;
146  }
147 }
148 
149 
150 template<class CloudType>
152 (
153  const label,
154  const label,
155  const scalar,
156  barycentric& coordinates,
157  label& celli,
158  label& tetFacei,
159  label& tetPti,
160  label& facei
161 )
162 {
164  (
165  this->owner().mesh(),
166  this->owner().rndGen(),
167  coordinates,
168  celli,
169  tetFacei,
170  tetPti,
171  facei
172  );
173 }
174 
175 
176 template<class CloudType>
178 (
179  const label,
180  const label,
181  const scalar,
182  typename CloudType::parcelType& parcel
183 )
184 {
185  // set particle velocity
186  parcel.U() = U0_;
187 
188  // set particle diameter
189  parcel.d() = sizeDistribution_->sample();
190 }
191 
192 
193 template<class CloudType>
195 {
196  return false;
197 }
198 
199 
200 // ************************************************************************* //
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:79
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:221
Templated injection model class.
Patch injection.
virtual void topoChange()
Set injector locations when mesh is updated.
virtual ~PatchInjection()
Destructor.
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
PatchInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
virtual void setPositionAndCell(const fvMesh &mesh, Random &rnd, barycentric &coordinates, label &celli, label &tetFacei, label &tetPti, label &facei)
Inherit setPositionAndCell from patchInjectionBase.
virtual label nParcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
virtual scalar timeEnd() const
Return the end-of-injection time.
virtual scalar massToInject(const scalar time0, const scalar time1)
Parcel mass to introduce relative to SOI.
Random number generator.
Definition: Random.H:58
scalar globalScalar01()
Advance the state and return a scalar sample from a uniform.
Definition: Random.C:94
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Accumulating histogram of values. Specified bin resolution automatic generation of bins.
Definition: distribution.H:61
Base class for patch-based injection models.
virtual void setPositionAndCell(const fvMesh &mesh, Random &rnd, barycentric &coordinates, label &celli, label &tetFacei, label &tetPti, label &facei)
Set the injection position and owner cell, tetFace and tetPt.
virtual void topoChange(const polyMesh &mesh)
Update patch geometry and derived info for injection locations.
A class for handling words, derived from string.
Definition: word.H:62
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
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
T clone(const T &t)
Definition: List.H:55
dictionary dict
Random rndGen(label(0))