PatchFlowRateInjection.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-2020 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 "PatchFlowRateInjection.H"
30 #include "distributionModel.H"
32 #include "surfaceFields.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
36 template<class CloudType>
38 (
39  const dictionary& dict,
40  CloudType& owner,
41  const word& modelName
42 )
43 :
44  InjectionModel<CloudType>(dict, owner, modelName,typeName),
45  patchInjectionBase(owner.mesh(), this->coeffDict().getWord("patch")),
46  phiName_(this->coeffDict().template getOrDefault<word>("phi", "phi")),
47  rhoName_(this->coeffDict().template getOrDefault<word>("rho", "rho")),
48  duration_(this->coeffDict().getScalar("duration")),
49  concentration_
50  (
51  Function1<scalar>::New
52  (
53  "concentration",
54  this->coeffDict(),
55  &owner.mesh()
56  )
57  ),
58  parcelConcentration_
59  (
60  this->coeffDict().getScalar("parcelConcentration")
61  ),
62  sizeDistribution_
63  (
65  (
66  this->coeffDict().subDict("sizeDistribution"),
67  owner.rndGen()
68  )
69  )
70 {
71  // Convert from user time to reduce the number of time conversion calls
72  const Time& time = owner.db().time();
73  duration_ = time.userTimeToTime(duration_);
74  concentration_->userTimeToTime(time);
75 
77 
78  // Re-initialise total mass/volume to inject to zero
79  // - will be reset during each injection
80  this->volumeTotal_ = 0.0;
81  this->massTotal_ = 0.0;
82 }
83 
84 
85 template<class CloudType>
87 (
89 )
90 :
93  phiName_(im.phiName_),
94  rhoName_(im.rhoName_),
95  duration_(im.duration_),
96  concentration_(im.concentration_.clone()),
97  parcelConcentration_(im.parcelConcentration_),
98  sizeDistribution_(im.sizeDistribution_.clone())
99 {}
100 
101 
102 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
103 
104 template<class CloudType>
106 {}
107 
108 
109 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
110 
111 template<class CloudType>
113 {
114  patchInjectionBase::updateMesh(this->owner().mesh());
115 }
116 
117 
118 template<class CloudType>
120 {
121  return this->SOI_ + duration_;
122 }
123 
124 
125 template<class CloudType>
127 {
128  const polyMesh& mesh = this->owner().mesh();
129 
130  const auto& phi = mesh.lookupObject<surfaceScalarField>(phiName_);
131 
132  const scalarField& phip = phi.boundaryField()[patchId_];
133 
134  scalar flowRateIn = 0.0;
135  if (phi.dimensions() == dimVolume/dimTime)
136  {
137  flowRateIn = max(0.0, -sum(phip));
138  }
139  else
140  {
141  const auto& rho = mesh.lookupObject<volScalarField>(rhoName_);
142  const scalarField& rhop = rho.boundaryField()[patchId_];
143 
144  flowRateIn = max(0.0, -sum(phip/rhop));
145  }
146 
147  reduce(flowRateIn, sumOp<scalar>());
149  return flowRateIn;
150 }
151 
152 
153 template<class CloudType>
155 (
156  const scalar time0,
157  const scalar time1
158 )
159 {
160  if ((time0 >= 0.0) && (time0 < duration_))
161  {
162  scalar dt = time1 - time0;
163 
164  scalar c = concentration_->value(0.5*(time0 + time1));
165 
166  scalar nParcels = parcelConcentration_*c*flowRate()*dt;
167 
168  Random& rnd = this->owner().rndGen();
169 
170  label nParcelsToInject = floor(nParcels);
171 
172  // Inject an additional parcel with a probability based on the
173  // remainder after the floor function
174  if
175  (
176  nParcelsToInject > 0
177  && (
178  nParcels - scalar(nParcelsToInject)
179  > rnd.globalPosition(scalar(0), scalar(1))
180  )
181  )
182  {
183  ++nParcelsToInject;
184  }
185 
186  return nParcelsToInject;
187  }
189  return 0;
190 }
191 
192 
193 template<class CloudType>
195 (
196  const scalar time0,
197  const scalar time1
198 )
199 {
200  scalar volume = 0.0;
201 
202  if ((time0 >= 0.0) && (time0 < duration_))
203  {
204  scalar c = concentration_->value(0.5*(time0 + time1));
205 
206  volume = c*(time1 - time0)*flowRate();
207  }
208 
209  this->volumeTotal_ = volume;
210  this->massTotal_ = volume*this->owner().constProps().rho0();
212  return volume;
213 }
214 
215 
216 template<class CloudType>
218 (
219  const label,
220  const label,
221  const scalar,
222  vector& position,
223  label& cellOwner,
224  label& tetFacei,
225  label& tetPti
226 )
227 {
229  (
230  this->owner().mesh(),
231  this->owner().rndGen(),
232  position,
233  cellOwner,
234  tetFacei,
235  tetPti
236  );
237 }
238 
239 
240 template<class CloudType>
242 (
243  const label,
244  const label,
245  const scalar,
246  typename CloudType::parcelType& parcel
247 )
248 {
249  // Set particle velocity to carrier velocity
250  parcel.U() = this->owner().U()[parcel.cell()];
252  // Set particle diameter
253  parcel.d() = sizeDistribution_->sample();
254 }
255 
256 
257 template<class CloudType>
259 {
260  return false;
261 }
262 
263 
264 template<class CloudType>
266 {
267  return true;
268 }
269 
270 
271 // ************************************************************************* //
Foam::surfaceFields.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
const Type & value() const noexcept
Return const reference to value.
dictionary dict
scalar massTotal_
Total mass to inject [kg].
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1)
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Templated injection model class.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
Random rndGen
Definition: createFields.H:23
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
virtual void updateMesh(const polyMesh &mesh)
Update patch geometry and derived info for injection locations.
virtual void updateMesh()
Set injector locations when mesh is updated.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
const CloudType & owner() const
Return const access to the owner cloud.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
dynamicFvMesh & mesh
virtual scalar userTimeToTime(const scalar theta) const
Convert the user-time (e.g. CA deg) to real-time (s).
Definition: TimeState.C:42
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
const Time & time() const noexcept
Return time registry.
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:450
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
Random number generator.
Definition: Random.H:55
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:290
const wordList volume
Standard volume field types (scalar, vector, tensor, etc)
Type globalPosition(const Type &start, const Type &end)
Return a sample on the interval [start,end].
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.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
Patch injection, by using patch flow rate to determine concentration and velocity.
const fvMesh & mesh() const
Return reference to the mesh.
Definition: DSMCCloudI.H:37
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.
label setPositionAndCell(const fvMesh &mesh, const scalar fraction01, Random &rnd, vector &position, label &cellOwner, label &tetFacei, label &tetPti)
Set the injection position and owner cell, tetFace and tetPt.
PatchFlowRateInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
const dimensionedScalar c
Speed of light in a vacuum.
A library of runtime-selectable doubly-truncated probability distribution models. Returns random samp...
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
scalar volumeTotal_
Total volume of particles introduced by this injector [m^3] Note: scaled to ensure massTotal is achie...
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:67
virtual scalar flowRate() const
Return the total volumetric flow rate across the patch [m3/s].
virtual ~PatchFlowRateInjection()
Destructor.