FieldActivatedInjection.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 
27 #include "volFields.H"
28 #include "mathematicalConstants.H"
29 
30 using namespace Foam::constant::mathematical;
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  factor_(this->coeffDict().template lookup<scalar>("factor")),
44  referenceField_
45  (
46  owner.db().objectRegistry::template lookupObject<volScalarField>
47  (
48  this->coeffDict().lookup("referenceField")
49  )
50  ),
51  thresholdField_
52  (
53  owner.db().objectRegistry::template lookupObject<volScalarField>
54  (
55  this->coeffDict().lookup("thresholdField")
56  )
57  ),
58  positionsFile_(this->coeffDict().lookup("positionsFile")),
59  positions_
60  (
61  IOobject
62  (
63  positionsFile_,
64  owner.db().time().constant(),
65  owner.mesh(),
66  IOobject::MUST_READ,
67  IOobject::NO_WRITE
68  )
69  ),
70  injectorCoordinates_(positions_.size()),
71  injectorCells_(positions_.size()),
72  injectorTetFaces_(positions_.size()),
73  injectorTetPts_(positions_.size()),
74  massTotal_(this->readMassTotal(dict, owner)),
75  nParcelsPerInjector_
76  (
77  this->coeffDict().template lookup<label>("parcelsPerInjector")
78  ),
79  nParcelsInjected_(positions_.size(), 0),
80  U0_(this->coeffDict().lookup("U0")),
81  diameters_(positions_.size()),
82  sizeDistribution_
83  (
85  (
86  dimLength,
87  this->coeffDict().subDict("sizeDistribution"),
88  this->sizeSampleQ(),
89  owner.rndGen().generator()
90  )
91  )
92 {
93  // Construct parcel diameters - one per injector cell
94  forAll(diameters_, i)
95  {
96  diameters_[i] = sizeDistribution_->sample();
97  }
98 
99  topoChange();
100 }
101 
102 
103 template<class CloudType>
105 (
107 )
108 :
110  factor_(im.factor_),
111  referenceField_(im.referenceField_),
112  thresholdField_(im.thresholdField_),
113  positionsFile_(im.positionsFile_),
114  positions_(im.positions_),
115  injectorCoordinates_(im.injectorCoordinates_),
116  injectorCells_(im.injectorCells_),
117  injectorTetFaces_(im.injectorTetFaces_),
118  injectorTetPts_(im.injectorTetPts_),
119  massTotal_(im.massTotal_),
120  nParcelsPerInjector_(im.nParcelsPerInjector_),
121  nParcelsInjected_(im.nParcelsInjected_),
122  U0_(im.U0_),
123  diameters_(im.diameters_),
124  sizeDistribution_(im.sizeDistribution_, false)
125 {}
126 
127 
128 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
129 
130 template<class CloudType>
132 {}
133 
134 
135 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
136 
137 template<class CloudType>
139 {
140  // Set/cache the injector cells
141  forAll(positions_, i)
142  {
143  this->findCellAtPosition
144  (
145  positions_[i],
146  injectorCoordinates_[i],
147  injectorCells_[i],
148  injectorTetFaces_[i],
149  injectorTetPts_[i]
150  );
151  }
152 }
153 
154 
155 template<class CloudType>
157 {
158  return great;
159 }
160 
161 
162 template<class CloudType>
164 (
165  const scalar time0,
166  const scalar time1
167 )
168 {
169  if (sum(nParcelsInjected_) < nParcelsPerInjector_*positions_.size())
170  {
171  return positions_.size();
172  }
173  else
174  {
175  return 0;
176  }
177 }
178 
179 
180 template<class CloudType>
182 (
183  const scalar time0,
184  const scalar time1
185 )
186 {
187  if (sum(nParcelsInjected_) < nParcelsPerInjector_*positions_.size())
188  {
189  return massTotal_/nParcelsPerInjector_;
190  }
191  else
192  {
193  return 0;
194  }
195 }
196 
197 
198 template<class CloudType>
200 (
201  const label parceli,
202  const label,
203  const scalar,
204  barycentric& coordinates,
205  label& celli,
206  label& tetFacei,
207  label& tetPti,
208  label& facei
209 )
210 {
211  const label injectorCelli = injectorCells_[parceli];
212 
213  if
214  (
215  nParcelsInjected_[parceli] < nParcelsPerInjector_
216  && factor_*referenceField_[injectorCelli] > thresholdField_[injectorCelli]
217  )
218  {
219  coordinates = injectorCoordinates_[parceli];
220  celli = injectorCells_[parceli];
221  tetFacei = injectorTetFaces_[parceli];
222  tetPti = injectorTetPts_[parceli];
223 
224  nParcelsInjected_[parceli]++;
225  }
226 }
227 
228 
229 template<class CloudType>
231 (
232  const label parceli,
233  const label,
234  const scalar,
235  typename CloudType::parcelType::trackingData& td,
236  typename CloudType::parcelType& parcel
237 )
238 {
239  // set particle velocity
240  parcel.U() = U0_;
241 
242  // set particle diameter
243  parcel.d() = diameters_[parceli];
244 }
245 
246 
247 template<class CloudType>
249 {
250  return false;
251 }
252 
253 
254 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:80
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:225
Injection at specified positions, with the conditions:
virtual void topoChange()
Set injector locations when mesh is updated.
virtual ~FieldActivatedInjection()
Destructor.
virtual void setProperties(const label parceli, const label nParcels, const scalar time, typename CloudType::parcelType::trackingData &td, typename CloudType::parcelType &parcel)
Set the parcel properties.
virtual label nParcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
FieldActivatedInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
virtual void setPositionAndCell(const label parceli, const label nParcels, const scalar time, barycentric &coordinates, label &celli, label &tetFacei, label &tetPti, label &facei)
Set the injection position and owner cell, tetFace and tetPt.
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
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.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Templated injection model class.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Base class for statistical distributions.
Definition: distribution.H:76
Registry of regIOobjects.
A class for handling words, derived from string.
Definition: word.H:62
mathematical constants.
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)
const dimensionSet dimLength
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
dictionary dict
randomGenerator rndGen(653213)