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-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 "PatchFlowRateInjection.H"
27 #include "TimeFunction1.H"
28 #include "distribution.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->readDuration(dict, owner)),
47  concentration_
48  (
49  TimeFunction1<scalar>
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  this->sizeSampleQ()
67  )
68  )
69 {}
70 
71 
72 template<class CloudType>
74 (
76 )
77 :
80  phiName_(im.phiName_),
81  rhoName_(im.rhoName_),
82  duration_(im.duration_),
83  concentration_(im.concentration_),
84  parcelConcentration_(im.parcelConcentration_),
85  sizeDistribution_(im.sizeDistribution_().clone().ptr())
86 {}
87 
88 
89 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
90 
91 template<class CloudType>
93 {}
94 
95 
96 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
97 
98 template<class CloudType>
100 {
101  patchInjectionBase::topoChange(this->owner().mesh());
102 }
103 
104 
105 template<class CloudType>
107 {
108  return this->SOI_ + duration_;
109 }
110 
111 
112 template<class CloudType>
114 {
115  const polyMesh& mesh = this->owner().mesh();
116 
117  const surfaceScalarField& phi =
118  mesh.lookupObject<surfaceScalarField>(phiName_);
119 
120  const scalarField& phip = phi.boundaryField()[patchId_];
121 
122  scalar flowRateIn = 0.0;
123  if (phi.dimensions() == dimFlux)
124  {
125  flowRateIn = max(0.0, -sum(phip));
126  }
127  else
128  {
129  const volScalarField& rho =
130  mesh.lookupObject<volScalarField>(rhoName_);
131  const scalarField& rhop = rho.boundaryField()[patchId_];
132 
133  flowRateIn = max(0.0, -sum(phip/rhop));
134  }
135 
136  reduce(flowRateIn, sumOp<scalar>());
137 
138  return flowRateIn;
139 }
140 
141 
142 template<class CloudType>
144 (
145  const scalar time0,
146  const scalar time1
147 )
148 {
149  if (time0 >= 0 && time0 < duration_)
150  {
151  scalar dt = time1 - time0;
152 
153  scalar c = concentration_.value(0.5*(time0 + time1));
154 
155  scalar nParcels = parcelConcentration_*c*flowRate()*dt;
156 
157  Random& rnd = this->owner().rndGen();
158 
159  label nParcelsToInject = floor(nParcels);
160 
161  // Inject an additional parcel with a probability based on the
162  // remainder after the floor function
163  if
164  (
165  nParcelsToInject > 0
166  && nParcels - scalar(nParcelsToInject) > rnd.globalScalar01()
167  )
168  {
169  ++nParcelsToInject;
170  }
171 
172  return nParcelsToInject;
173  }
174  else
175  {
176  return 0;
177  }
178 }
179 
180 
181 template<class CloudType>
183 (
184  const scalar time0,
185  const scalar time1
186 )
187 {
188  scalar volume = 0;
189 
190  if (time0 >= 0 && time0 < duration_)
191  {
192  scalar c = concentration_.value(0.5*(time0 + time1));
193 
194  volume = c*(time1 - time0)*flowRate();
195  }
196 
197  return volume*this->owner().constProps().rho0();
198 }
199 
200 
201 template<class CloudType>
203 (
204  const label,
205  const label,
206  const scalar,
207  barycentric& coordinates,
208  label& celli,
209  label& tetFacei,
210  label& tetPti,
211  label& facei
212 )
213 {
215  (
216  this->owner().mesh(),
217  this->owner().rndGen(),
218  coordinates,
219  celli,
220  tetFacei,
221  tetPti,
222  facei
223  );
224 }
225 
226 
227 template<class CloudType>
229 (
230  const label,
231  const label,
232  const scalar,
233  typename CloudType::parcelType& parcel
234 )
235 {
236  // Set particle velocity to carrier velocity
237  parcel.U() = this->owner().U()[parcel.cell()];
238 
239  // Set particle diameter
240  parcel.d() = sizeDistribution_->sample();
241 }
242 
243 
244 template<class CloudType>
246 {
247  return false;
248 }
249 
250 
251 // ************************************************************************* //
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:79
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:221
const dimensionSet & dimensions() const
Return dimensions.
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 ~PatchFlowRateInjection()
Destructor.
PatchFlowRateInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
scalar flowRate() const
Return the total volumetric flow rate across the patch [m^3/s].
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
Light wrapper around Function1 to provide a mechanism to update time-based entries.
Definition: TimeFunction1.H:61
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
const Type & value() const
Return const reference to value.
Accumulating histogram of values. Specified bin resolution automatic generation of bins.
Definition: distribution.H:61
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 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.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
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
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
T clone(const T &t)
Definition: List.H:55
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimFlux
dictionary dict
Foam::surfaceFields.
Random rndGen(label(0))