ReactingMultiphaseLookupTableInjection.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2019-2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class CloudType>
36(
37 const dictionary& dict,
38 CloudType& owner,
39 const word& modelName
40)
41:
42 InjectionModel<CloudType>(dict, owner, modelName, typeName),
43 inputFileName_(this->coeffDict().lookup("inputFile")),
44 duration_(this->coeffDict().getScalar("duration")),
45 parcelsPerSecond_(this->coeffDict().getScalar("parcelsPerSecond")),
46 randomise_(this->coeffDict().getBool("randomise")),
47 injectors_
48 (
50 (
51 inputFileName_,
52 owner.db().time().constant(),
53 owner.db(),
54 IOobject::MUST_READ,
55 IOobject::NO_WRITE
56 )
57 ),
58 injectorCells_(injectors_.size()),
59 injectorTetFaces_(injectors_.size()),
60 injectorTetPts_(injectors_.size())
61{
62 updateMesh();
63
64 duration_ = owner.db().time().userTimeToTime(duration_);
65
66 // Determine volume of particles to inject
67 this->volumeTotal_ = 0.0;
68 forAll(injectors_, i)
69 {
70 this->volumeTotal_ += injectors_[i].mDot()/injectors_[i].rho();
71 }
72 this->volumeTotal_ *= duration_;
73}
74
75
76template<class CloudType>
79(
81)
82:
84 inputFileName_(im.inputFileName_),
85 duration_(im.duration_),
86 parcelsPerSecond_(im.parcelsPerSecond_),
87 randomise_(im.randomise_),
88 injectors_(im.injectors_),
89 injectorCells_(im.injectorCells_),
90 injectorTetFaces_(im.injectorTetFaces_),
91 injectorTetPts_(im.injectorTetPts_)
92{}
93
94
95// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
96
97template<class CloudType>
99{
100 // Set/cache the injector cells
101 bitSet reject(injectors_.size());
102
103 // Set/cache the injector cells
104 forAll(injectors_, i)
105 {
106 if
107 (
108 !this->findCellAtPosition
109 (
110 injectorCells_[i],
111 injectorTetFaces_[i],
112 injectorTetPts_[i],
113 injectors_[i].x(),
114 !this->ignoreOutOfBounds_
115
116 )
117 )
118 {
119 reject.set(i);
120 }
121 }
122
123 label nRejected = reject.count();
124
125 if (nRejected)
126 {
127 reject.flip();
128 inplaceSubset(reject, injectorCells_);
129 inplaceSubset(reject, injectorTetFaces_);
130 inplaceSubset(reject, injectorTetPts_);
131 inplaceSubset(reject, injectors_);
132
133 Info<< " " << nRejected
134 << " positions rejected, out of bounds" << endl;
135 }
136}
137
138
139template<class CloudType>
140Foam::scalar
142{
143 return this->SOI_ + duration_;
144}
145
146
147template<class CloudType>
148Foam::label
150(
151 const scalar time0,
152 const scalar time1
153)
154{
155 if ((time0 >= 0.0) && (time0 < duration_))
156 {
157 return floor(injectorCells_.size()*(time1 - time0)*parcelsPerSecond_);
158 }
159
160 return 0;
161}
162
163
164template<class CloudType>
165Foam::scalar
167(
168 const scalar time0,
169 const scalar time1
170)
171{
172 scalar volume = 0.0;
173 if ((time0 >= 0.0) && (time0 < duration_))
174 {
175 forAll(injectors_, i)
176 {
177 volume += injectors_[i].mDot()/injectors_[i].rho()*(time1 - time0);
178 }
179 }
180
181 return volume;
182}
183
184
185template<class CloudType>
187(
188 const label parcelI,
189 const label nParcels,
190 const scalar time,
191 vector& position,
192 label& cellOwner,
193 label& tetFacei,
194 label& tetPti
195)
196{
197 label injectorI = 0;
198 if (randomise_)
199 {
200 Random& rnd = this->owner().rndGen();
201 injectorI = rnd.position<label>(0, injectorCells_.size() - 1);
202 }
203 else
204 {
205 injectorI = parcelI*injectorCells_.size()/nParcels;
206 }
207
208 position = injectors_[injectorI].x();
209 cellOwner = injectorCells_[injectorI];
210 tetFacei = injectorTetFaces_[injectorI];
211 tetPti = injectorTetPts_[injectorI];
212}
213
214
215template<class CloudType>
217(
218 const label parcelI,
219 const label nParcels,
220 const scalar,
221 typename CloudType::parcelType& parcel
222)
223{
224 label injectorI = parcelI*injectorCells_.size()/nParcels;
225
226 // set particle velocity
227 parcel.U() = injectors_[injectorI].U();
228
229 // set particle diameter
230 parcel.d() = injectors_[injectorI].d();
231
232 // set particle density
233 parcel.rho() = injectors_[injectorI].rho();
234
235 // set particle temperature
236 parcel.T() = injectors_[injectorI].T();
237
238 // set particle specific heat capacity
239 parcel.Cp() = injectors_[injectorI].Cp();
240
241 // set particle component mass fractions
242 parcel.Y() = injectors_[injectorI].Y();
243
244 // set particle gaseous component mass fractions
245 parcel.YGas() = injectors_[injectorI].YGas();
246
247 // set particle liquid component mass fractions
248 parcel.YLiquid() = injectors_[injectorI].YLiquid();
249
250 // set particle solid component mass fractions
251 parcel.YSolid() = injectors_[injectorI].YSolid();
252}
253
254
255template<class CloudType>
256bool
258{
259 return true;
260}
261
262
263template<class CloudType>
265(
266 const label
267)
268{
269 return true;
270}
271
272
273// ************************************************************************* //
const CloudType & owner() const
Return const access to the owner cloud.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:500
Templated injection model class.
Random number generator.
Definition: Random.H:60
Type position(const Type &start, const Type &end)
Return a sample on the interval [start,end].
Particle injection sources read from look-up table. Each row corresponds to an injection site.
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
virtual void setPositionAndCell(const label parcelI, const label nParcels, const scalar time, vector &position, label &cellOwner, label &tetFacei, label &tetPti)
Set the injection position and owner cell, tetFace and tetPt.
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
virtual void updateMesh()
Set injector locations when mesh is updated.
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
scalar timeEnd() const
Return the end-of-injection time.
virtual scalar userTimeToTime(const scalar theta) const
Convert the user-time (e.g. CA deg) to real-time (s).
Definition: TimeState.C:49
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:66
unsigned int count(const bool on=true) const
Count number of bits set.
Definition: bitSetI.H:500
void flip()
Invert all bits in the addressable region.
Definition: bitSetI.H:634
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:590
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const Time & time() const noexcept
Return time registry.
constant condensation/saturation model.
Lookup type of boundary radiation properties.
Definition: lookup.H:66
A class for handling words, derived from Foam::string.
Definition: word.H:68
void inplaceSubset(const BoolListType &select, ListType &input, const bool invert=false)
Inplace extract elements of the input list when select is true.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333