Loading...
Searching...
No Matches
CellZoneInjection.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-2025 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 "CellZoneInjection.H"
32#include "globalIndex.H"
33#include "Pstream.H"
34
35// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36
37template<class CloudType>
38void Foam::CellZoneInjection<CloudType>::setPositions
39(
40 const labelList& cellZoneCells
41)
42{
43 const fvMesh& mesh = this->owner().mesh();
44 const scalarField& V = mesh.V();
45 const label nCells = cellZoneCells.size();
46 Random& rnd = this->owner().rndGen();
47
48 DynamicList<vector> positions(nCells); // initial size only
49 DynamicList<label> injectorCells(nCells); // initial size only
50 DynamicList<label> injectorTetFaces(nCells); // initial size only
51 DynamicList<label> injectorTetPts(nCells); // initial size only
52
53 scalar newParticlesTotal = 0.0;
54 label addParticlesTotal = 0;
55
56 forAll(cellZoneCells, i)
57 {
58 const label celli = cellZoneCells[i];
59
60 // Calc number of particles to add
61 const scalar newParticles = V[celli]*numberDensity_;
62 newParticlesTotal += newParticles;
63 label addParticles = floor(newParticles);
64 addParticlesTotal += addParticles;
65
66 const scalar diff = newParticlesTotal - addParticlesTotal;
67 if (diff > 1)
68 {
69 label corr = floor(diff);
70 addParticles += corr;
71 addParticlesTotal += corr;
72 }
73
74 // Construct cell tet indices
75 const List<tetIndices> cellTetIs =
76 polyMeshTetDecomposition::cellTetIndices(mesh, celli);
77
78 // Construct cell tet volume fractions
79 scalarList cTetVFrac(cellTetIs.size(), Zero);
80 for (label tetI = 1; tetI < cellTetIs.size() - 1; tetI++)
81 {
82 cTetVFrac[tetI] =
83 cTetVFrac[tetI-1] + cellTetIs[tetI].tet(mesh).mag()/V[celli];
84 }
85 cTetVFrac.last() = 1.0;
86
87 // Set new particle position and cellId
88 for (label pI = 0; pI < addParticles; pI++)
89 {
90 const scalar volFrac = rnd.sample01<scalar>();
91 label tetI = 0;
92 forAll(cTetVFrac, vfI)
93 {
94 if (cTetVFrac[vfI] > volFrac)
95 {
96 tetI = vfI;
97 break;
98 }
99 }
100 positions.append(cellTetIs[tetI].tet(mesh).randomPoint(rnd));
101
102 injectorCells.append(celli);
103 injectorTetFaces.append(cellTetIs[tetI].face());
104 injectorTetPts.append(cellTetIs[tetI].tetPt());
105 }
106 }
107
108 // Parallel operation manipulations
109 const label myProci = UPstream::myProcNo();
110 globalIndex globalPositions(positions.size());
111
112 List<vector> allPositions(globalPositions.totalSize(), point::max);
113 List<label> allInjectorCells(globalPositions.totalSize(), -1);
114 List<label> allInjectorTetFaces(globalPositions.totalSize(), -1);
115 List<label> allInjectorTetPts(globalPositions.totalSize(), -1);
116
117 // Gather all positions on to all processors
118 SubList<vector>
119 (
120 allPositions,
121 globalPositions.range(myProci)
122 ) = positions;
123
124 Pstream::listReduce(allPositions, minOp<point>());
125
126 // Gather local cell tet and tet-point Ids, but leave non-local ids set -1
127 SubList<label>
128 (
129 allInjectorCells,
130 globalPositions.range(myProci)
131 ) = injectorCells;
132 SubList<label>
133 (
134 allInjectorTetFaces,
135 globalPositions.range(myProci)
136 ) = injectorTetFaces;
137 SubList<label>
138 (
139 allInjectorTetPts,
140 globalPositions.range(myProci)
141 ) = injectorTetPts;
142
143 // Transfer data
144 positions_.transfer(allPositions);
145 injectorCells_.transfer(allInjectorCells);
146 injectorTetFaces_.transfer(allInjectorTetFaces);
147 injectorTetPts_.transfer(allInjectorTetPts);
148
149
150 if (debug)
151 {
152 OFstream points("points.obj");
153 forAll(positions_, i)
154 {
155 meshTools::writeOBJ(points, positions_[i]);
156 }
158}
159
160
161// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
162
163template<class CloudType>
165(
166 const dictionary& dict,
167 CloudType& owner,
168 const word& modelName
169)
170:
171 InjectionModel<CloudType>(dict, owner, modelName, typeName),
172 cellZoneName_(this->coeffDict().lookup("cellZone")),
173 numberDensity_(this->coeffDict().getScalar("numberDensity")),
174 positions_(),
175 injectorCells_(),
176 injectorTetFaces_(),
177 injectorTetPts_(),
178 diameters_(),
179 U0_(this->coeffDict().lookup("U0")),
180 sizeDistribution_
181 (
182 distributionModel::New
183 (
184 this->coeffDict().subDict("sizeDistribution"), owner.rndGen()
185 )
187{
188 updateMesh();
189}
190
191
192template<class CloudType>
194(
196)
197:
199 cellZoneName_(im.cellZoneName_),
200 numberDensity_(im.numberDensity_),
201 positions_(im.positions_),
202 injectorCells_(im.injectorCells_),
203 injectorTetFaces_(im.injectorTetFaces_),
204 injectorTetPts_(im.injectorTetPts_),
205 diameters_(im.diameters_),
206 U0_(im.U0_),
207 sizeDistribution_(im.sizeDistribution_.clone())
208{}
209
210
211// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
212
213template<class CloudType>
215{}
216
217
218// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
219
220template<class CloudType>
222{
223 // Set/cache the injector cells
224 const fvMesh& mesh = this->owner().mesh();
225 const label zoneI = mesh.cellZones().findZoneID(cellZoneName_);
226
227 if (zoneI < 0)
228 {
230 << "Unknown cell zone name: " << cellZoneName_
231 << ". Valid cell zones are: " << mesh.cellZones().names()
232 << nl << exit(FatalError);
233 }
234
235 const labelList& cellZoneCells = mesh.cellZones()[zoneI];
236 const label nCells = cellZoneCells.size();
237 const scalar nCellsTotal = returnReduce(nCells, sumOp<label>());
238 const scalar VCells = sum(scalarField(mesh.V(), cellZoneCells));
239 const scalar VCellsTotal = returnReduce(VCells, sumOp<scalar>());
240 Info<< " cell zone size = " << nCellsTotal << endl;
241 Info<< " cell zone volume = " << VCellsTotal << endl;
242
243 if ((nCellsTotal == 0) || (VCellsTotal*numberDensity_ < 1))
244 {
246 << "Number of particles to be added to cellZone " << cellZoneName_
247 << " is zero" << endl;
248 }
249 else
250 {
251 setPositions(cellZoneCells);
252
253 Info<< " number density = " << numberDensity_ << nl
254 << " number of particles = " << positions_.size() << endl;
255
256 // Construct parcel diameters
257 diameters_.setSize(positions_.size());
258 forAll(diameters_, i)
259 {
260 diameters_[i] = sizeDistribution_->sample();
261 }
262 }
264 // Determine volume of particles to inject
265 this->volumeTotal_ = sum(pow3(diameters_))*constant::mathematical::pi/6.0;
266}
267
268
269template<class CloudType>
272 // Not used
273 return this->SOI_;
274}
275
276
277template<class CloudType>
279(
280 const scalar time0,
281 const scalar time1
282)
283{
284 if ((0.0 >= time0) && (0.0 < time1))
285 {
286 return positions_.size();
288
289 return 0;
290}
291
292
293template<class CloudType>
295(
296 const scalar time0,
297 const scalar time1
298)
299{
300 // All parcels introduced at SOI
301 if ((0.0 >= time0) && (0.0 < time1))
302 {
303 return this->volumeTotal_;
305
306 return 0.0;
307}
308
309
310template<class CloudType>
312(
313 const label parcelI,
314 const label nParcels,
315 const scalar time,
316 vector& position,
317 label& cellOwner,
318 label& tetFacei,
319 label& tetPti
320)
321{
322 position = positions_[parcelI];
323 cellOwner = injectorCells_[parcelI];
324 tetFacei = injectorTetFaces_[parcelI];
325 tetPti = injectorTetPts_[parcelI];
326}
327
328
329template<class CloudType>
331(
332 const label parcelI,
333 const label,
334 const scalar,
335 typename CloudType::parcelType& parcel
336)
337{
338 // set particle velocity
339 parcel.U() = U0_;
341 // set particle diameter
342 parcel.d() = diameters_[parcelI];
343}
344
345
346template<class CloudType>
348{
349 return false;
350}
351
352
353template<class CloudType>
355{
356 return true;
357}
358
359
360// ************************************************************************* //
Injection positions specified by a particle number density within a cell set.
virtual autoPtr< InjectionModel< CloudType > > clone() const
Construct and return a clone.
scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
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 ~CellZoneInjection()
Destructor.
CellZoneInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
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.
const CloudType & owner() const
Return const access to the owner cloud.
const vector & U() const
Return const access to velocity.
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.
scalar SOI_
Start of injection [s].
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
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...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
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
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const pointField & points
#define WarningInFunction
Report a warning using Foam::Warning.
const expr V(m.psi().mesh().V())
constexpr scalar pi(M_PI)
DSMCCloud< dsmcParcel > CloudType
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedScalar pow3(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition triad.C:373
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Random rndGen