PatchFlowRateInjection.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-2021 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 "PatchFlowRateInjection.H"
27 #include "TimeFunction1.H"
28 #include "distributionModel.H"
29 #include "mathematicalConstants.H"
30 #include "surfaceFields.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
34 template<class CloudType>
36 (
37  const dictionary& dict,
38  CloudType& owner,
39  const word& modelName
40 )
41 :
42  InjectionModel<CloudType>(dict, owner, modelName,typeName),
43  patchInjectionBase(owner.mesh(), this->coeffDict().lookup("patchName")),
44  phiName_(this->coeffDict().template lookupOrDefault<word>("phi", "phi")),
45  rhoName_(this->coeffDict().template lookupOrDefault<word>("rho", "rho")),
46  duration_(this->coeffDict().template lookup<scalar>("duration")),
47  concentration_
48  (
50  (
51  owner.db().time(),
52  "concentration",
53  this->coeffDict()
54  )
55  ),
56  parcelConcentration_
57  (
58  this->coeffDict().template lookup<scalar>("parcelConcentration")
59  ),
60  sizeDistribution_
61  (
63  (
64  this->coeffDict().subDict("sizeDistribution"),
65  owner.rndGen()
66  )
67  )
68 {
69  duration_ = owner.db().time().userTimeToTime(duration_);
70 
71  patchInjectionBase::updateMesh(owner.mesh());
72 
73  // Re-initialise total mass/volume to inject to zero
74  // - will be reset during each injection
75  this->volumeTotal_ = 0.0;
76  this->massTotal_ = 0.0;
77 }
78 
79 
80 template<class CloudType>
82 (
84 )
85 :
88  phiName_(im.phiName_),
89  rhoName_(im.rhoName_),
90  duration_(im.duration_),
91  concentration_(im.concentration_),
92  parcelConcentration_(im.parcelConcentration_),
93  sizeDistribution_(im.sizeDistribution_().clone().ptr())
94 {}
95 
96 
97 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
98 
99 template<class CloudType>
101 {}
102 
103 
104 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
105 
106 template<class CloudType>
108 {
109  patchInjectionBase::updateMesh(this->owner().mesh());
110 }
111 
112 
113 template<class CloudType>
115 {
116  return this->SOI_ + duration_;
117 }
118 
119 
120 template<class CloudType>
122 {
123  const polyMesh& mesh = this->owner().mesh();
124 
125  const surfaceScalarField& phi =
126  mesh.lookupObject<surfaceScalarField>(phiName_);
127 
128  const scalarField& phip = phi.boundaryField()[patchId_];
129 
130  scalar flowRateIn = 0.0;
131  if (phi.dimensions() == dimFlux)
132  {
133  flowRateIn = max(0.0, -sum(phip));
134  }
135  else
136  {
137  const volScalarField& rho =
138  mesh.lookupObject<volScalarField>(rhoName_);
139  const scalarField& rhop = rho.boundaryField()[patchId_];
140 
141  flowRateIn = max(0.0, -sum(phip/rhop));
142  }
143 
144  reduce(flowRateIn, sumOp<scalar>());
145 
146  return flowRateIn;
147 }
148 
149 
150 template<class CloudType>
152 (
153  const scalar time0,
154  const scalar time1
155 )
156 {
157  if ((time0 >= 0.0) && (time0 < duration_))
158  {
159  scalar dt = time1 - time0;
160 
161  scalar c = concentration_.value(0.5*(time0 + time1));
162 
163  scalar nParcels = parcelConcentration_*c*flowRate()*dt;
164 
165  Random& rnd = this->owner().rndGen();
166 
167  label nParcelsToInject = floor(nParcels);
168 
169  // Inject an additional parcel with a probability based on the
170  // remainder after the floor function
171  if
172  (
173  nParcelsToInject > 0
174  && nParcels - scalar(nParcelsToInject) > rnd.globalScalar01()
175  )
176  {
177  ++nParcelsToInject;
178  }
179 
180  return nParcelsToInject;
181  }
182  else
183  {
184  return 0;
185  }
186 }
187 
188 
189 template<class CloudType>
191 (
192  const scalar time0,
193  const scalar time1
194 )
195 {
196  scalar volume = 0.0;
197 
198  if ((time0 >= 0.0) && (time0 < duration_))
199  {
200  scalar c = concentration_.value(0.5*(time0 + time1));
201 
202  volume = c*(time1 - time0)*flowRate();
203  }
204 
205  this->volumeTotal_ = volume;
206  this->massTotal_ = volume*this->owner().constProps().rho0();
207 
208  return volume;
209 }
210 
211 
212 template<class CloudType>
214 (
215  const label,
216  const label,
217  const scalar,
218  vector& position,
219  label& cellOwner,
220  label& tetFacei,
221  label& tetPti
222 )
223 {
224  patchInjectionBase::setPositionAndCell
225  (
226  this->owner().mesh(),
227  this->owner().rndGen(),
228  position,
229  cellOwner,
230  tetFacei,
231  tetPti
232  );
233 }
234 
235 
236 template<class CloudType>
238 (
239  const label,
240  const label,
241  const scalar,
242  typename CloudType::parcelType& parcel
243 )
244 {
245  // Set particle velocity to carrier velocity
246  parcel.U() = this->owner().U()[parcel.cell()];
247 
248  // Set particle diameter
249  parcel.d() = sizeDistribution_->sample();
250 }
251 
252 
253 template<class CloudType>
255 {
256  return false;
257 }
258 
259 
260 template<class CloudType>
262 {
263  return true;
264 }
265 
266 
267 // ************************************************************************* //
Foam::surfaceFields.
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
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Templated injection model class.
virtual void updateMesh()
Set injector locations when mesh is updated.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
const dimensionedScalar c
Speed of light in a vacuum.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Random rndGen(label(0))
const dimensionSet & dimensions() const
Return dimensions.
dynamicFvMesh & mesh
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
A class for handling words, derived from string.
Definition: word.H:59
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
const dimensionSet dimFlux
Base class for patch-based injection models.
const Type & value() const
Return const reference to value.
phi
Definition: correctPhi.H:3
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
scalar globalScalar01()
Advance the state and return a scalar sample from a uniform.
Definition: Random.C:94
Random number generator.
Definition: Random.H:57
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:221
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, by using patch flow rate to determine concentration and velocity.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
scalar timeEnd() const
Return the end-of-injection time.
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
PatchFlowRateInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
virtual scalar flowRate() const
Return the total volumetric flow rate across the patch [m^3/s].
virtual ~PatchFlowRateInjection()
Destructor.