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-2025 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  volumeRatio_
47  (
48  this->coeffDict().found("concentration")
49  || this->coeffDict().found("volumeRatio")
50  ? Function1<scalar>::New
51  (
52  this->coeffDict().found("volumeRatio")
53  ? "volumeRatio"
54  : "concentration",
55  this->owner().db().time().userUnits(),
56  dimless,
57  this->coeffDict()
58  )
59  : autoPtr<Function1<scalar>>()
60  ),
61  massRatio_
62  (
63  this->coeffDict().found("massRatio")
64  ? Function1<scalar>::New
65  (
66  "massRatio",
67  this->owner().db().time().userUnits(),
68  dimless,
69  this->coeffDict()
70  )
71  : autoPtr<Function1<scalar>>()
72  ),
73  parcelConcentration_
74  (
75  this->coeffDict().template lookup<scalar>("parcelConcentration")
76  ),
77  sizeDistribution_
78  (
80  (
81  dimLength,
82  this->coeffDict().subDict("sizeDistribution"),
83  this->sizeSampleQ(),
84  owner.rndGen().generator()
85  )
86  )
87 {
88  if (volumeRatio_.valid() && massRatio_.valid())
89  {
91  << "keywords volumeRatio (or concentration) and "
92  << "massRatio both defined in dictionary "
93  << this->coeffDict().name() << exit(FatalIOError);
94  }
95 
96  if (!volumeRatio_.valid() && !massRatio_.valid())
97  {
99  << "keyword volumeRatio or massRatio is "
100  << "undefined in dictionary "
101  << this->coeffDict().name() << exit(FatalIOError);
102  }
103 }
104 
105 
106 template<class CloudType>
108 (
110 )
111 :
113  patchInjectionBase(im),
114  phiName_(im.phiName_),
115  rhoName_(im.rhoName_),
116  duration_(im.duration_),
117  volumeRatio_(im.volumeRatio_, false),
118  massRatio_(im.massRatio_, false),
119  parcelConcentration_(im.parcelConcentration_),
120  sizeDistribution_(im.sizeDistribution_, false)
121 {}
122 
123 
124 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
125 
126 template<class CloudType>
128 {}
129 
130 
131 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
132 
133 template<class CloudType>
135 {
136  patchInjectionBase::topoChange(this->owner().mesh());
137 }
138 
139 
140 template<class CloudType>
142 {
143  return this->SOI_ + duration_;
144 }
145 
146 
147 template<class CloudType>
149 {
150  const polyMesh& mesh = this->owner().mesh();
151 
152  const surfaceScalarField& phi =
154  const scalarField& phip = phi.boundaryField()[patchId_];
155 
156  scalar flowRateIn;
157  if (phi.dimensions() == dimVolumetricFlux)
158  {
159  flowRateIn = max(scalar(0), -sum(phip));
160  }
161  else
162  {
163  const volScalarField& rho =
164  mesh.lookupObject<volScalarField>(rhoName_);
165  const scalarField& rhop = rho.boundaryField()[patchId_];
166 
167  flowRateIn = max(scalar(0), -sum(phip/rhop));
168  }
169 
170  reduce(flowRateIn, sumOp<scalar>());
171 
172  return flowRateIn;
173 }
174 
175 
176 template<class CloudType>
178 {
179  const polyMesh& mesh = this->owner().mesh();
180 
181  const surfaceScalarField& phi =
183  const scalarField& phip = phi.boundaryField()[patchId_];
184 
185  scalar flowRateIn;
186  if (phi.dimensions() == dimVolumetricFlux)
187  {
188  const volScalarField& rho =
189  mesh.lookupObject<volScalarField>(rhoName_);
190  const scalarField& rhop = rho.boundaryField()[patchId_];
191 
192  flowRateIn = max(scalar(0), -sum(phip*rhop));
193  }
194  else
195  {
196  flowRateIn = max(scalar(0), -sum(phip));
197  }
198 
199  reduce(flowRateIn, sumOp<scalar>());
200 
201  return flowRateIn;
202 }
203 
204 
205 template<class CloudType>
207 (
208  const scalar t0,
209  const scalar t1
210 )
211 {
212  if (t1 >= 0 && t0 < duration_)
213  {
214  const scalar tMid = (t0 + t1)/2;
215 
216  const scalar volumeFlowRateToInject =
217  volumeRatio_.valid()
218  ? volumeRatio_->value(tMid)
219  *volumetricFlowRate()
220  : massRatio_->value(tMid)
221  *massFlowRate()
222  /this->owner().constProps().rho0();
223 
224  return
225  parcelConcentration_
226  *volumeFlowRateToInject
227  *(min(t1, duration_) - max(t0, 0));
228  }
229  else
230  {
231  return 0;
232  }
233 }
234 
235 
236 template<class CloudType>
238 (
239  const scalar t0,
240  const scalar t1
241 )
242 {
243  if (t1 >= 0 && t0 < duration_)
244  {
245  const scalar tMid = (t0 + t1)/2;
246 
247  const scalar massFlowRateToInject =
248  volumeRatio_.valid()
249  ? volumeRatio_->value(tMid)
250  *volumetricFlowRate()
251  *this->owner().constProps().rho0()
252  : massRatio_->value(tMid)
253  *massFlowRate();
254 
255  return
256  massFlowRateToInject
257  *(min(t1, duration_) - max(t0, 0));
258  }
259  else
260  {
261  return 0;
262  }
263 }
264 
265 
266 template<class CloudType>
268 (
269  const label,
270  const label,
271  const scalar,
273  label& celli,
274  label& tetFacei,
275  label& tetPti,
276  label& facei
277 )
278 {
280  (
281  this->owner().mesh(),
282  this->owner().rndGen(),
283  coordinates,
284  celli,
285  tetFacei,
286  tetPti,
287  facei
288  );
289 }
290 
291 
292 template<class CloudType>
294 (
295  const label,
296  const label,
297  const scalar,
298  typename CloudType::parcelType::trackingData& td,
299  typename CloudType::parcelType& parcel
300 )
301 {
302  // Set particle velocity to carrier velocity
303  parcel.U() = this->owner().U()[parcel.cell()];
304 
305  // Set particle diameter
306  parcel.d() = sizeDistribution_->sample();
307 }
308 
309 
310 template<class CloudType>
312 {
313  return false;
314 }
315 
316 
317 // ************************************************************************* //
bool found
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 in which the quantity of particulate material introduced is specified as a fraction o...
virtual scalar nParcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
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 massFlowRate() const
Return the total mass flow rate across the patch [kg/s].
scalar volumetricFlowRate() const
Return the total volumetric flow rate across the patch [m^3/s].
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.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
const fileName & name() const
Return the dictionary name.
Definition: dictionary.H:111
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Base class for statistical distributions.
Definition: distribution.H:76
const polyMesh & mesh() const
Return reference to polyMesh.
Definition: fvMesh.H:443
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
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
Definition: subModelBase.C:128
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
barycentric coordinates(const polyMesh &mesh, const point &position, const label celli, const label facei, const label faceTrii, const scalar stepFraction)
Return the coordinates given the position and tet topology.
Definition: tracking.C:1259
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
const dimensionSet dimless
const dimensionSet dimLength
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh, PrimitiveField > &df)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
dictionary dict
randomGenerator rndGen(653213)
Foam::surfaceFields.