Loading...
Searching...
No Matches
ManualInjection.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-2022 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
29#include "ManualInjection.H"
31#include "bitSet.H"
32
33using namespace Foam::constant::mathematical;
34
35// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36
37template<class CloudType>
39(
40 const dictionary& dict,
42 const word& modelName
43)
44:
46 positionsFile_(this->coeffDict().lookup("positionsFile")),
47 positions_
48 (
50 (
51 positionsFile_,
52 owner.db().time().constant(),
53 owner.mesh(),
54 IOobject::MUST_READ,
55 IOobject::NO_WRITE
56 )
57 ),
58 diameters_(positions_.size()),
59 injectorCells_(positions_.size(), -1),
60 injectorTetFaces_(positions_.size(), -1),
61 injectorTetPts_(positions_.size(), -1),
62 U0_(this->coeffDict().lookup("U0")),
63 sizeDistribution_
64 (
66 (
67 this->coeffDict().subDict("sizeDistribution"),
69 )
70 ),
71 ignoreOutOfBounds_
72 (
73 this->coeffDict().getOrDefault("ignoreOutOfBounds", false)
74 )
75{
76 updateMesh();
77
78 // Construct parcel diameters
79 forAll(diameters_, i)
80 {
81 diameters_[i] = sizeDistribution_->sample();
82 }
84 // Determine volume of particles to inject
85 this->volumeTotal_ = sum(pow3(diameters_))*pi/6.0;
86}
87
88
89template<class CloudType>
91(
93)
94:
96 positionsFile_(im.positionsFile_),
97 positions_(im.positions_),
98 diameters_(im.diameters_),
99 injectorCells_(im.injectorCells_),
100 injectorTetFaces_(im.injectorTetFaces_),
101 injectorTetPts_(im.injectorTetPts_),
102 U0_(im.U0_),
103 sizeDistribution_(im.sizeDistribution_.clone()),
104 ignoreOutOfBounds_(im.ignoreOutOfBounds_)
105{}
106
107
108// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
109
110template<class CloudType>
112{}
113
114
115// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116
117template<class CloudType>
119{
120 label nRejected = 0;
121
122 bitSet keep(positions_.size(), true);
123
124 forAll(positions_, pI)
125 {
126 if
127 (
128 !this->findCellAtPosition
129 (
130 injectorCells_[pI],
131 injectorTetFaces_[pI],
132 injectorTetPts_[pI],
133 positions_[pI],
134 !ignoreOutOfBounds_
135 )
136 )
137 {
138 keep.unset(pI);
139 nRejected++;
140 }
141 }
142
143 if (nRejected > 0)
144 {
145 inplaceSubset(keep, positions_);
146 inplaceSubset(keep, diameters_);
147 inplaceSubset(keep, injectorCells_);
148 inplaceSubset(keep, injectorTetFaces_);
149 inplaceSubset(keep, injectorTetPts_);
150
151 Info<< " " << nRejected
152 << " particles ignored, out of bounds" << endl;
153 }
154}
155
156
157template<class CloudType>
159{
160 // Injection is instantaneous - but allow for a finite interval to
161 // avoid numerical issues when interval is zero
162 return this->SOI_ + SMALL;
163}
164
165
166template<class CloudType>
168(
169 const scalar time0,
170 const scalar time1
171)
172{
173 if ((0.0 >= time0) && (0.0 < time1))
174 {
175 return positions_.size();
177
178 return 0;
179}
180
181
182template<class CloudType>
184(
185 const scalar time0,
186 const scalar time1
187)
188{
189 // All parcels introduced at SOI
190 if ((0.0 >= time0) && (0.0 < time1))
191 {
192 return this->volumeTotal_;
194
195 return 0.0;
196}
197
198
199template<class CloudType>
201(
202 const label parcelI,
203 const label,
204 const scalar,
205 vector& position,
206 label& cellOwner,
207 label& tetFacei,
208 label& tetPti
209)
210{
211 position = positions_[parcelI];
212 cellOwner = injectorCells_[parcelI];
213 tetFacei = injectorTetFaces_[parcelI];
214 tetPti = injectorTetPts_[parcelI];
215}
216
217
218template<class CloudType>
220(
221 const label parcelI,
222 const label,
223 const scalar,
224 typename CloudType::parcelType& parcel
225)
226{
227 // set particle velocity
228 parcel.U() = U0_;
230 // set particle diameter
231 parcel.d() = diameters_[parcelI];
232}
233
234
235template<class CloudType>
237{
238 return false;
239}
240
241
242template<class CloudType>
244{
245 return true;
246}
247
248
249// ************************************************************************* //
const CloudType & owner() const
Return const access to the owner cloud.
const vector & U() const
Return const access to velocity.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
Templated injection model class.
static autoPtr< InjectionModel< CloudType > > New(const dictionary &dict, CloudType &owner)
Selector with lookup from dictionary.
scalar volumeTotal_
Total volume of particles introduced by this injector [m^3] Note: scaled to ensure massTotal is achie...
InjectionModel(CloudType &owner)
Construct null from owner.
virtual bool findCellAtPosition(label &celli, label &tetFacei, label &tetPti, vector &position, bool errorOnNotFound=true)
Find the cell that contains the supplied position.
scalar SOI_
Start of injection [s].
Manual injection.
virtual autoPtr< InjectionModel< CloudType > > clone() const
Construct and return a clone.
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
virtual ~ManualInjection()
Destructor.
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
ManualInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
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.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
bitSet & unset(const bitSet &other)
Unset (subtract) the bits specified in the other bitset, which is a set difference corresponds to the...
Definition bitSetI.H:540
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A library of runtime-selectable doubly-truncated probability distribution models. Returns random samp...
Lookup type of boundary radiation properties.
Definition lookup.H:60
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
const word & modelName() const
Return const access to the name of the sub-model.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
dynamicFvMesh & mesh
Mathematical constants.
constexpr scalar pi(M_PI)
Different types of constants.
void inplaceSubset(const BoolListType &select, ListType &input, const bool invert=false)
Inplace extract elements of the input list when select is true.
DSMCCloud< dsmcParcel > CloudType
dimensionedScalar pow3(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
Vector< scalar > vector
Definition vector.H:57
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Random rndGen