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-2024 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 "distribution.H"
28 #include "mathematicalConstants.H"
29 #include "surfaceFields.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class CloudType>
35 (
36  const dictionary& dict,
37  CloudType& owner,
38  const word& modelName
39 )
40 :
41  InjectionModel<CloudType>(dict, owner, modelName,typeName),
42  patchInjectionBase(owner.mesh(), this->coeffDict().lookup("patchName")),
43  phiName_(this->coeffDict().template lookupOrDefault<word>("phi", "phi")),
44  rhoName_(this->coeffDict().template lookupOrDefault<word>("rho", "rho")),
45  duration_(this->readDuration(dict, owner)),
46  concentration_
47  (
48  Function1<scalar>::New
49  (
50  "concentration",
51  this->owner().db().time().userUnits(),
52  dimless,
53  this->coeffDict()
54  )
55  ),
56  parcelConcentration_
57  (
58  this->coeffDict().template lookup<scalar>("parcelConcentration")
59  ),
60  sizeDistribution_
61  (
63  (
64  dimLength,
65  this->coeffDict().subDict("sizeDistribution"),
66  this->sizeSampleQ(),
67  owner.rndGen().generator()
68  )
69  )
70 {}
71 
72 
73 template<class CloudType>
75 (
77 )
78 :
81  phiName_(im.phiName_),
82  rhoName_(im.rhoName_),
83  duration_(im.duration_),
84  concentration_(im.concentration_, false),
85  parcelConcentration_(im.parcelConcentration_),
86  sizeDistribution_(im.sizeDistribution_, false)
87 {}
88 
89 
90 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
91 
92 template<class CloudType>
94 {}
95 
96 
97 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
98 
99 template<class CloudType>
101 {
102  patchInjectionBase::topoChange(this->owner().mesh());
103 }
104 
105 
106 template<class CloudType>
108 {
109  return this->SOI_ + duration_;
110 }
111 
112 
113 template<class CloudType>
115 {
116  const polyMesh& mesh = this->owner().mesh();
117 
118  const surfaceScalarField& phi =
119  mesh.lookupObject<surfaceScalarField>(phiName_);
120 
121  const scalarField& phip = phi.boundaryField()[patchId_];
122 
123  scalar flowRateIn = 0.0;
124  if (phi.dimensions() == dimVolumetricFlux)
125  {
126  flowRateIn = max(0.0, -sum(phip));
127  }
128  else
129  {
130  const volScalarField& rho =
131  mesh.lookupObject<volScalarField>(rhoName_);
132  const scalarField& rhop = rho.boundaryField()[patchId_];
133 
134  flowRateIn = max(0.0, -sum(phip/rhop));
135  }
136 
137  reduce(flowRateIn, sumOp<scalar>());
138 
139  return flowRateIn;
140 }
141 
142 
143 template<class CloudType>
145 (
146  const scalar time0,
147  const scalar time1
148 )
149 {
150  if (time0 >= 0 && time0 < duration_)
151  {
152  scalar dt = time1 - time0;
153 
154  scalar c = concentration_->value(0.5*(time0 + time1));
155 
156  scalar nParcels = parcelConcentration_*c*flowRate()*dt;
157 
158  randomGenerator& rndGen = this->owner().rndGen();
159 
160  label nParcelsToInject = floor(nParcels);
161 
162  // Inject an additional parcel with a probability based on the
163  // remainder after the floor function
164  if (nParcels - scalar(nParcelsToInject) > this->globalScalar01(rndGen))
165  {
166  ++nParcelsToInject;
167  }
168 
169  return nParcelsToInject;
170  }
171  else
172  {
173  return 0;
174  }
175 }
176 
177 
178 template<class CloudType>
180 (
181  const scalar time0,
182  const scalar time1
183 )
184 {
185  scalar volume = 0;
186 
187  if (time0 >= 0 && time0 < duration_)
188  {
189  scalar c = concentration_->value(0.5*(time0 + time1));
190 
191  volume = c*(time1 - time0)*flowRate();
192  }
193 
194  return volume*this->owner().constProps().rho0();
195 }
196 
197 
198 template<class CloudType>
200 (
201  const label,
202  const label,
203  const scalar,
204  barycentric& coordinates,
205  label& celli,
206  label& tetFacei,
207  label& tetPti,
208  label& facei
209 )
210 {
212  (
213  this->owner().mesh(),
214  this->owner().rndGen(),
215  coordinates,
216  celli,
217  tetFacei,
218  tetPti,
219  facei
220  );
221 }
222 
223 
224 template<class CloudType>
226 (
227  const label,
228  const label,
229  const scalar,
230  typename CloudType::parcelType::trackingData& td,
231  typename CloudType::parcelType& parcel
232 )
233 {
234  // Set particle velocity to carrier velocity
235  parcel.U() = this->owner().U()[parcel.cell()];
236 
237  // Set particle diameter
238  parcel.d() = sizeDistribution_->sample();
239 }
240 
241 
242 template<class CloudType>
244 {
245  return false;
246 }
247 
248 
249 // ************************************************************************* //
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:80
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:225
const dimensionSet & dimensions() const
Return dimensions.
Run-time selectable general function of one variable.
Definition: Function1.H:125
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Templated injection model class.
Patch injection, by using patch flow rate to determine concentration and velocity.
virtual void topoChange()
Set injector locations when mesh is updated.
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType::trackingData &td, typename CloudType::parcelType &parcel)
virtual ~PatchFlowRateInjection()
Destructor.
PatchFlowRateInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
scalar flowRate() const
Return the total volumetric flow rate across the patch [m^3/s].
virtual label nParcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
virtual void setPositionAndCell(const fvMesh &mesh, randomGenerator &rndGen, barycentric &coordinates, label &celli, label &tetFacei, label &tetPti, label &facei)
Inherit setPositionAndCell from patchInjectionBase.
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.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const Type & value() const
Return const reference to value.
Base class for statistical distributions.
Definition: distribution.H:76
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
Base class for patch-based injection models.
virtual void topoChange(const polyMesh &mesh)
Update patch geometry and derived info for injection locations.
virtual void setPositionAndCell(const fvMesh &mesh, randomGenerator &rndGen, barycentric &coordinates, label &celli, label &tetFacei, label &tetPti, label &facei)
Set the injection position and owner cell, tetFace and tetPt.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
Random number generator.
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)
const dimensionedScalar c
Speed of light in a vacuum.
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
const dimensionSet dimVolumetricFlux
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
const dimensionSet dimless
const dimensionSet dimLength
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dictionary dict
randomGenerator rndGen(653213)
Foam::surfaceFields.