ThermoLookupTableInjection.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#include "scalarIOList.H"
31
32// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33
34template<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>
78(
80)
81:
83 inputFileName_(im.inputFileName_),
84 duration_(im.duration_),
85 parcelsPerSecond_(im.parcelsPerSecond_),
86 randomise_(im.randomise_),
87 injectors_(im.injectors_),
88 injectorCells_(im.injectorCells_),
89 injectorTetFaces_(im.injectorTetFaces_),
90 injectorTetPts_(im.injectorTetPts_)
91{}
92
93
94// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95
96template<class CloudType>
98{
99 // Set/cache the injector cells
100 bitSet reject(injectors_.size());
101
102 // Set/cache the injector cells
103 forAll(injectors_, i)
104 {
105 if
106 (
107 !this->findCellAtPosition
108 (
109 injectorCells_[i],
110 injectorTetFaces_[i],
111 injectorTetPts_[i],
112 injectors_[i].x(),
113 !this->ignoreOutOfBounds_
114
115 )
116 )
117 {
118 reject.set(i);
119 }
120 }
121
122 const label nRejected = reject.count();
123
124 if (nRejected)
125 {
126 reject.flip();
127 inplaceSubset(reject, injectorCells_);
128 inplaceSubset(reject, injectorTetFaces_);
129 inplaceSubset(reject, injectorTetPts_);
130 inplaceSubset(reject, injectors_);
131
132 Info<< " " << nRejected
133 << " positions rejected, out of bounds" << endl;
134 }
135}
136
137
138template<class CloudType>
140{
141 return this->SOI_ + duration_;
142}
143
144
145template<class CloudType>
147(
148 const scalar time0,
149 const scalar time1
150)
151{
152 if ((time0 >= 0.0) && (time0 < duration_))
153 {
154 return floor(injectorCells_.size()*(time1 - time0)*parcelsPerSecond_);
155 }
156
157 return 0;
158}
159
160
161template<class CloudType>
163(
164 const scalar time0,
165 const scalar time1
166)
167{
168 scalar volume = 0.0;
169 if ((time0 >= 0.0) && (time0 < duration_))
170 {
171 forAll(injectors_, i)
172 {
173 volume += injectors_[i].mDot()/injectors_[i].rho()*(time1 - time0);
174 }
175 }
176
177 return volume;
178}
179
180
181template<class CloudType>
183(
184 const label parcelI,
185 const label nParcels,
186 const scalar time,
187 vector& position,
188 label& cellOwner,
189 label& tetFacei,
190 label& tetPti
191)
192{
193 label injectorI = 0;
194 if (randomise_)
195 {
196 Random& rnd = this->owner().rndGen();
197 injectorI = rnd.position<label>(0, injectorCells_.size() - 1);
198 }
199 else
200 {
201 injectorI = parcelI*injectorCells_.size()/nParcels;
202 }
203
204 position = injectors_[injectorI].x();
205 cellOwner = injectorCells_[injectorI];
206 tetFacei = injectorTetFaces_[injectorI];
207 tetPti = injectorTetPts_[injectorI];
208}
209
210
211template<class CloudType>
213(
214 const label parcelI,
215 const label nParcels,
216 const scalar,
217 typename CloudType::parcelType* pPtr
218)
219{
220 label injectorI = parcelI*injectorCells_.size()/nParcels;
221
222 // set particle velocity
223 parcel.U() = injectors_[injectorI].U();
224
225 // set particle diameter
226 parcel.d() = injectors_[injectorI].d();
227
228 // set particle density
229 parcel.rho() = injectors_[injectorI].rho();
230
231 // set particle temperature
232 parcel.T() = injectors_[injectorI].T();
233
234 // set particle specific heat capacity
235 parcel.cp() = injectors_[injectorI].cp();
236}
237
238
239template<class CloudType>
241{
242 return true;
243}
244
245
246template<class CloudType>
248{
249 return true;
250}
251
252
253// ************************************************************************* //
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 bool validInjection(const label parcelI)
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