PatchInjection.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 "PatchInjection.H"
27 #include "distribution.H"
28 
29 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
30 
31 template<class CloudType>
33 (
34  typename CloudType::parcelType::trackingData& td
35 )
36 {
38 
39  if (U0Name_ == word::null)
40  {
41  U0InterpPtr_.clear();
42  }
43  else if (U0Name_ == this->owner().U().name())
44  {
45  U0InterpPtr_ = tmpNrc<interpolation<vector>>(td.UInterp());
46  }
47  else
48  {
49  U0InterpPtr_ =
51  (
53  (
54  this->owner().solution().interpolationSchemes(),
55  this->owner().mesh().template lookupObject<volVectorField>
56  (
57  U0Name_
58  )
59  ).ptr()
60  );
61  }
62 }
63 
64 
65 template<class CloudType>
67 (
68  const label parcelsAdded,
69  const scalar massAdded,
70  typename CloudType::parcelType::trackingData& td
71 )
72 {
73  InjectionModel<CloudType>::postInject(parcelsAdded, massAdded, td);
74 
75  U0InterpPtr_.clear();
76 }
77 
78 
79 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
80 
81 template<class CloudType>
83 (
84  const dictionary& dict,
85  CloudType& owner,
86  const word& modelName
87 )
88 :
89  InjectionModel<CloudType>(dict, owner, modelName, typeName),
90  patchInjectionBase(owner.mesh(), this->coeffDict().lookup("patchName")),
91  duration_(this->readDuration(dict, owner)),
92  massFlowRate_(this->readMassFlowRate(dict, owner, duration_)),
93  parcelsPerSecond_(this->readParcelsPerSecond(dict, owner)),
94  U0_(vector::uniform(NaN)),
95  U0Name_(word::null),
96  U0InterpPtr_(nullptr),
97  sizeDistribution_
98  (
100  (
101  dimLength,
102  this->coeffDict().subDict("sizeDistribution"),
103  this->sizeSampleQ(),
104  owner.rndGen().generator()
105  )
106  )
107 {
108  ITstream& is = this->coeffDict().lookup("U0");
109 
110  token t(is);
111  is.putBack(t);
112 
113  if (t.isWord())
114  {
115  U0Name_ = word(is);
116  }
117  else
118  {
119  U0_ = vector(is);
120  }
121 }
122 
123 
124 template<class CloudType>
126 (
127  const PatchInjection<CloudType>& im
128 )
129 :
131  patchInjectionBase(im),
132  duration_(im.duration_),
133  massFlowRate_(im.massFlowRate_, false),
134  parcelsPerSecond_(im.parcelsPerSecond_, false),
135  U0_(im.U0_),
136  U0Name_(word::null),
137  U0InterpPtr_(nullptr),
138  sizeDistribution_(im.sizeDistribution_, false)
139 {}
140 
141 
142 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
143 
144 template<class CloudType>
146 {}
147 
148 
149 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
150 
151 template<class CloudType>
153 {
154  patchInjectionBase::topoChange(this->owner().mesh());
155 }
156 
157 
158 template<class CloudType>
160 {
161  return this->SOI_ + duration_;
162 }
163 
164 
165 template<class CloudType>
167 (
168  const scalar t0,
169  const scalar t1
170 )
171 {
172  if (t1 >= 0 && t0 < duration_)
173  {
174  return parcelsPerSecond_->integral(max(t0, 0), min(t1, duration_));
175  }
176  else
177  {
178  return 0;
179  }
180 }
181 
182 
183 template<class CloudType>
185 (
186  const scalar t0,
187  const scalar t1
188 )
189 {
190  if (t1 >= 0 && t0 < duration_)
191  {
192  return massFlowRate_->integral(max(t0, 0), min(t1, duration_));
193  }
194  else
195  {
196  return 0;
197  }
198 }
199 
200 
201 template<class CloudType>
203 (
204  const label,
205  const label,
206  const scalar,
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::trackingData& td,
234  typename CloudType::parcelType& parcel
235 )
236 {
237  // set particle velocity
238  parcel.U() =
239  U0InterpPtr_.valid()
240  ? U0InterpPtr_().interpolate
241  (
242  parcel.coordinates(),
243  parcel.currentTetIndices(td.mesh)
244  )
245  : U0_;
246 
247  // set particle diameter
248  parcel.d() = sizeDistribution_->sample();
249 }
250 
251 
252 template<class CloudType>
254 {
255  return false;
256 }
257 
258 
259 // ************************************************************************* //
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:80
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:225
Input token stream.
Definition: ITstream.H:53
Templated injection model class.
void putBack(const token &)
Put back token.
Definition: Istream.C:30
Patch injection.
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 ~PatchInjection()
Destructor.
virtual void preInject(typename CloudType::parcelType::trackingData &td)
Pre injection hook.
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType::trackingData &td, typename CloudType::parcelType &parcel)
virtual void postInject(const label parcelsAdded, const scalar massAdded, typename CloudType::parcelType::trackingData &td)
Post injection hook.
PatchInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
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 keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:740
Base class for statistical distributions.
Definition: distribution.H:76
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.
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:51
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
Definition: subModelBase.C:128
A class for managing temporary objects without reference counting.
Definition: tmpNrc.H:53
A token holds items read from Istream.
Definition: token.H:73
bool isWord() const
Definition: tokenI.H:298
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)
U
Definition: pEqn.H:72
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
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 dimLength
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
dictionary dict
randomGenerator rndGen(653213)