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 "TimeFunction1.H"
31 #include "distributionModel.H"
32 #include "mathematicalConstants.H"
33 #include "surfaceFields.H"
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
37 template<class CloudType>
39 (
40  const dictionary& dict,
41  CloudType& owner,
42  const word& modelName
43 )
44 :
45  InjectionModel<CloudType>(dict, owner, modelName,typeName),
46  patchInjectionBase(owner.mesh(), this->coeffDict().getWord("patch")),
47  phiName_(this->coeffDict().template getOrDefault<word>("phi", "phi")),
48  rhoName_(this->coeffDict().template getOrDefault<word>("rho", "rho")),
49  duration_(this->coeffDict().getScalar("duration")),
50  concentration_
51  (
53  (
54  owner.db().time(),
55  "concentration",
56  this->coeffDict()
57  )
58  ),
59  parcelConcentration_
60  (
61  this->coeffDict().getScalar("parcelConcentration")
62  ),
63  sizeDistribution_
64  (
66  (
67  this->coeffDict().subDict("sizeDistribution"),
68  owner.rndGen()
69  )
70  )
71 {
72  duration_ = owner.db().time().userTimeToTime(duration_);
73 
74  patchInjectionBase::updateMesh(owner.mesh());
75 
76  // Re-initialise total mass/volume to inject to zero
77  // - will be reset during each injection
78  this->volumeTotal_ = 0.0;
79  this->massTotal_ = 0.0;
80 }
81 
82 
83 template<class CloudType>
85 (
87 )
88 :
91  phiName_(im.phiName_),
92  rhoName_(im.rhoName_),
93  duration_(im.duration_),
94  concentration_(im.concentration_),
95  parcelConcentration_(im.parcelConcentration_),
96  sizeDistribution_(im.sizeDistribution_.clone())
97 {}
98 
99 
100 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
101 
102 template<class CloudType>
104 {}
105 
106 
107 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
108 
109 template<class CloudType>
111 {
112  patchInjectionBase::updateMesh(this->owner().mesh());
113 }
114 
115 
116 template<class CloudType>
118 {
119  return this->SOI_ + duration_;
120 }
121 
122 
123 template<class CloudType>
125 {
126  const polyMesh& mesh = this->owner().mesh();
127 
128  const surfaceScalarField& phi =
129  mesh.lookupObject<surfaceScalarField>(phiName_);
130 
131  const scalarField& phip = phi.boundaryField()[patchId_];
132 
133  scalar flowRateIn = 0.0;
134  if (phi.dimensions() == dimVelocity*dimArea)
135  {
136  flowRateIn = max(0.0, -sum(phip));
137  }
138  else
139  {
140  const volScalarField& rho =
141  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>());
148 
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  }
188 
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();
211 
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 {
228  patchInjectionBase::setPositionAndCell
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()];
251 
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::PatchFlowRateInjection
Patch injection, by using patch flow rate to determine concentration and velocity.
Definition: PatchFlowRateInjection.H:70
Foam::Random
Random number generator.
Definition: Random.H:59
Foam::PatchFlowRateInjection::flowRate
virtual scalar flowRate() const
Return the total volumetric flow rate across the patch [m3/s].
Definition: PatchFlowRateInjection.C:124
mathematicalConstants.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::DSMCCloud::constProps
const List< typename ParcelType::constantProperties > & constProps() const
Return all of the constant properties.
Definition: DSMCCloudI.H:98
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::PatchFlowRateInjection::setProperties
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
Definition: PatchFlowRateInjection.C:242
Foam::InjectionModel
Templated injection model class.
Definition: InjectionModel.H:73
Foam::PatchFlowRateInjection::~PatchFlowRateInjection
virtual ~PatchFlowRateInjection()
Destructor.
Definition: PatchFlowRateInjection.C:103
Foam::TimeFunction1< scalar >
surfaceFields.H
Foam::surfaceFields.
Foam::fieldTypes::volume
const wordList volume
Standard volume field types (scalar, vector, tensor, etc)
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
rho
rho
Definition: readInitialConditions.H:88
Foam::PatchFlowRateInjection::updateMesh
virtual void updateMesh()
Set injector locations when mesh is updated.
Definition: PatchFlowRateInjection.C:110
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::sumOp
Definition: ops.H:213
Foam::DSMCCloud::rndGen
Random & rndGen()
Return reference to the random object.
Definition: DSMCCloudI.H:124
Foam::PatchFlowRateInjection::timeEnd
scalar timeEnd() const
Return the end-of-injection time.
Definition: PatchFlowRateInjection.C:117
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:60
Foam::Field< scalar >
Foam::PatchFlowRateInjection::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: PatchFlowRateInjection.C:218
Foam::DSMCCloud::mesh
const fvMesh & mesh() const
Return reference to the mesh.
Definition: DSMCCloudI.H:44
Foam::PatchFlowRateInjection::validInjection
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
Definition: PatchFlowRateInjection.C:265
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:71
Foam::PatchFlowRateInjection::volumeToInject
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
Definition: PatchFlowRateInjection.C:195
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
reduce
reduce(hasMovingMesh, orOp< bool >())
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
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:1500
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
Foam::PatchFlowRateInjection::PatchFlowRateInjection
PatchFlowRateInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
Definition: PatchFlowRateInjection.C:39
Foam::Vector< scalar >
Foam::patchInjectionBase
Definition: patchInjectionBase.H:64
Foam::PatchFlowRateInjection::parcelsToInject
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
Definition: PatchFlowRateInjection.C:155
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
PatchFlowRateInjection.H
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
rndGen
Random rndGen
Definition: createFields.H:23
Foam::DSMCCloud::parcelType
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
Foam::PatchFlowRateInjection::fullyDescribed
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
Definition: PatchFlowRateInjection.C:258
TimeFunction1.H
Foam::GeometricField< scalar, fvsPatchField, surfaceMesh >
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