FieldActivatedInjection.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-2017 OpenFOAM Foundation
9 Copyright (C) 2015-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 "volFields.H"
32
33using namespace Foam::constant::mathematical;
34
35// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36
37template<class CloudType>
39(
40 const dictionary& dict,
41 CloudType& owner,
42 const word& modelName
43)
44:
45 InjectionModel<CloudType>(dict, owner, modelName, typeName),
46 factor_(this->coeffDict().getScalar("factor")),
47 referenceField_
48 (
49 owner.db().objectRegistry::template lookupObject<volScalarField>
50 (
51 this->coeffDict().getWord("referenceField")
52 )
53 ),
54 thresholdField_
55 (
56 owner.db().objectRegistry::template lookupObject<volScalarField>
57 (
58 this->coeffDict().getWord("thresholdField")
59 )
60 ),
61 positionsFile_(this->coeffDict().getWord("positionsFile")),
62 positions_
63 (
65 (
66 positionsFile_,
67 owner.db().time().constant(),
68 owner.mesh(),
69 IOobject::MUST_READ,
70 IOobject::NO_WRITE
71 )
72 ),
73 injectorCells_(positions_.size()),
74 injectorTetFaces_(positions_.size()),
75 injectorTetPts_(positions_.size()),
76 nParcelsPerInjector_
77 (
78 this->coeffDict().getLabel("parcelsPerInjector")
79 ),
80 nParcelsInjected_(),
81 U0_
82 (
83 this->coeffDict().template get<vector>("U0")
84 ),
85 diameters_(),
86 sizeDistribution_
87 (
89 (
90 this->coeffDict().subDict("sizeDistribution"),
91 owner.rndGen()
92 )
93 )
94{
95 updateMesh();
96
97 nParcelsInjected_.setSize(positions_.size(), Zero);
98
99 // Construct parcel diameters - one per injector cell
100 diameters_.setSize(positions_.size());
101 forAll(diameters_, i)
102 {
103 diameters_[i] = sizeDistribution_->sample();
104 }
105
106 // Determine total volume of particles to inject
107 this->volumeTotal_ =
108 nParcelsPerInjector_*sum(pow3(diameters_))*pi/6.0;
109}
110
111
112template<class CloudType>
114(
116)
117:
119 factor_(im.factor_),
120 referenceField_(im.referenceField_),
121 thresholdField_(im.thresholdField_),
122 positionsFile_(im.positionsFile_),
123 positions_(im.positions_),
124 injectorCells_(im.injectorCells_),
125 injectorTetFaces_(im.injectorTetFaces_),
126 injectorTetPts_(im.injectorTetPts_),
127 nParcelsPerInjector_(im.nParcelsPerInjector_),
128 nParcelsInjected_(im.nParcelsInjected_),
129 U0_(im.U0_),
130 diameters_(im.diameters_),
131 sizeDistribution_(im.sizeDistribution_.clone())
132{}
133
134
135// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
136
137template<class CloudType>
139{
140 bitSet reject(positions_.size());
141
142 // Set/cache the injector cells
143 forAll(positions_, i)
144 {
145 if
146 (
147 !this->findCellAtPosition
148 (
149 injectorCells_[i],
150 injectorTetFaces_[i],
151 injectorTetPts_[i],
152 positions_[i],
153 !this->ignoreOutOfBounds_
154
155 )
156 )
157 {
158 reject.set(i);
159 }
160 }
161
162 const label nRejected = reject.count();
163
164 if (nRejected)
165 {
166 reject.flip();
167 inplaceSubset(reject, injectorCells_);
168 inplaceSubset(reject, injectorTetFaces_);
169 inplaceSubset(reject, injectorTetPts_);
170 inplaceSubset(reject, positions_);
171
172 Info<< " " << nRejected
173 << " positions rejected, out of bounds" << endl;
174 }
175}
176
177
178template<class CloudType>
180{
181 return GREAT;
182}
183
184
185template<class CloudType>
187(
188 const scalar time0,
189 const scalar time1
190)
191{
192 if (sum(nParcelsInjected_) < nParcelsPerInjector_*positions_.size())
193 {
194 return positions_.size();
195 }
196
197 return 0;
198}
199
200
201template<class CloudType>
203(
204 const scalar time0,
205 const scalar time1
206)
207{
208 if (sum(nParcelsInjected_) < nParcelsPerInjector_*positions_.size())
209 {
210 return this->volumeTotal_/nParcelsPerInjector_;
211 }
212
213 return 0.0;
214}
215
216
217template<class CloudType>
219(
220 const label parcelI,
221 const label,
222 const scalar,
223 vector& position,
224 label& cellOwner,
225 label& tetFacei,
226 label& tetPti
227)
228{
229 position = positions_[parcelI];
230 cellOwner = injectorCells_[parcelI];
231 tetFacei = injectorTetFaces_[parcelI];
232 tetPti = injectorTetPts_[parcelI];
233}
234
235
236template<class CloudType>
238(
239 const label parcelI,
240 const label,
241 const scalar,
242 typename CloudType::parcelType& parcel
243)
244{
245 // set particle velocity
246 parcel.U() = U0_;
247
248 // set particle diameter
249 parcel.d() = diameters_[parcelI];
250}
251
252
253template<class CloudType>
255{
256 return false;
257}
258
259
260template<class CloudType>
262(
263 const label parcelI
264)
265{
266 const label celli = injectorCells_[parcelI];
267
268 if
269 (
270 nParcelsInjected_[parcelI] < nParcelsPerInjector_
271 && factor_*referenceField_[celli] > thresholdField_[celli]
272 )
273 {
274 nParcelsInjected_[parcelI]++;
275 return true;
276 }
277
278 return false;
279}
280
281
282// ************************************************************************* //
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
Conditional injection at specified positions.
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)
Return flag to identify whether or not injection of parcelI is.
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.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
Templated injection model class.
void setSize(const label n)
Alias for resize()
Definition: List.H:218
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
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
A library of runtime-selectable doubly-truncated probability distribution models. Returns random samp...
Registry of regIOobjects.
constant condensation/saturation model.
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
Mathematical constants.
constexpr scalar pi(M_PI)
void inplaceSubset(const BoolListType &select, ListType &input, const bool invert=false)
Inplace extract elements of the input list when select is true.
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Random rndGen
Definition: createFields.H:23