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-2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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"
30 #include "mathematicalConstants.H"
32 #include "globalIndex.H"
33 #include "Pstream.H"
34 
35 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36 
37 template<class CloudType>
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  globalIndex globalPositions(positions.size());
110  List<vector> allPositions(globalPositions.size(), point::max);
111  List<label> allInjectorCells(globalPositions.size(), -1);
112  List<label> allInjectorTetFaces(globalPositions.size(), -1);
113  List<label> allInjectorTetPts(globalPositions.size(), -1);
114 
115  // Gather all positions on to all processors
116  SubList<vector>
117  (
118  allPositions,
119  globalPositions.localSize(Pstream::myProcNo()),
120  globalPositions.offset(Pstream::myProcNo())
121  ) = positions;
122 
123  Pstream::listCombineGather(allPositions, minEqOp<point>());
124  Pstream::listCombineScatter(allPositions);
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.localSize(Pstream::myProcNo()),
131  globalPositions.offset(Pstream::myProcNo())
132  ) = injectorCells;
133  SubList<label>
134  (
135  allInjectorTetFaces,
136  globalPositions.localSize(Pstream::myProcNo()),
137  globalPositions.offset(Pstream::myProcNo())
138  ) = injectorTetFaces;
139  SubList<label>
140  (
141  allInjectorTetPts,
142  globalPositions.localSize(Pstream::myProcNo()),
143  globalPositions.offset(Pstream::myProcNo())
144  ) = injectorTetPts;
145 
146  // Transfer data
147  positions_.transfer(allPositions);
148  injectorCells_.transfer(allInjectorCells);
149  injectorTetFaces_.transfer(allInjectorTetFaces);
150  injectorTetPts_.transfer(allInjectorTetPts);
151 
152 
153  if (debug)
154  {
155  OFstream points("points.obj");
156  forAll(positions_, i)
157  {
158  meshTools::writeOBJ(points, positions_[i]);
159  }
160  }
161 }
162 
163 
164 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
165 
166 template<class CloudType>
168 (
169  const dictionary& dict,
170  CloudType& owner,
171  const word& modelName
172 )
173 :
174  InjectionModel<CloudType>(dict, owner, modelName, typeName),
175  cellZoneName_(this->coeffDict().lookup("cellZone")),
176  numberDensity_(this->coeffDict().getScalar("numberDensity")),
177  positions_(),
178  injectorCells_(),
179  injectorTetFaces_(),
180  injectorTetPts_(),
181  diameters_(),
182  U0_(this->coeffDict().lookup("U0")),
183  sizeDistribution_
184  (
186  (
187  this->coeffDict().subDict("sizeDistribution"), owner.rndGen()
188  )
189  )
190 {
191  updateMesh();
192 }
193 
194 
195 template<class CloudType>
197 (
199 )
200 :
202  cellZoneName_(im.cellZoneName_),
203  numberDensity_(im.numberDensity_),
204  positions_(im.positions_),
205  injectorCells_(im.injectorCells_),
206  injectorTetFaces_(im.injectorTetFaces_),
207  injectorTetPts_(im.injectorTetPts_),
208  diameters_(im.diameters_),
209  U0_(im.U0_),
210  sizeDistribution_(im.sizeDistribution_.clone())
211 {}
212 
213 
214 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
215 
216 template<class CloudType>
218 {}
219 
220 
221 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
222 
223 template<class CloudType>
225 {
226  // Set/cache the injector cells
227  const fvMesh& mesh = this->owner().mesh();
228  const label zoneI = mesh.cellZones().findZoneID(cellZoneName_);
229 
230  if (zoneI < 0)
231  {
233  << "Unknown cell zone name: " << cellZoneName_
234  << ". Valid cell zones are: " << mesh.cellZones().names()
235  << nl << exit(FatalError);
236  }
237 
238  const labelList& cellZoneCells = mesh.cellZones()[zoneI];
239  const label nCells = cellZoneCells.size();
240  const scalar nCellsTotal = returnReduce(nCells, sumOp<label>());
241  const scalar VCells = sum(scalarField(mesh.V(), cellZoneCells));
242  const scalar VCellsTotal = returnReduce(VCells, sumOp<scalar>());
243  Info<< " cell zone size = " << nCellsTotal << endl;
244  Info<< " cell zone volume = " << VCellsTotal << endl;
245 
246  if ((nCellsTotal == 0) || (VCellsTotal*numberDensity_ < 1))
247  {
249  << "Number of particles to be added to cellZone " << cellZoneName_
250  << " is zero" << endl;
251  }
252  else
253  {
254  setPositions(cellZoneCells);
255 
256  Info<< " number density = " << numberDensity_ << nl
257  << " number of particles = " << positions_.size() << endl;
258 
259  // Construct parcel diameters
260  diameters_.setSize(positions_.size());
261  forAll(diameters_, i)
262  {
263  diameters_[i] = sizeDistribution_->sample();
264  }
265  }
266 
267  // Determine volume of particles to inject
268  this->volumeTotal_ = sum(pow3(diameters_))*constant::mathematical::pi/6.0;
269 }
270 
271 
272 template<class CloudType>
274 {
275  // Not used
276  return this->SOI_;
277 }
278 
279 
280 template<class CloudType>
282 (
283  const scalar time0,
284  const scalar time1
285 )
286 {
287  if ((0.0 >= time0) && (0.0 < time1))
288  {
289  return positions_.size();
290  }
291 
292  return 0;
293 }
294 
295 
296 template<class CloudType>
298 (
299  const scalar time0,
300  const scalar time1
301 )
302 {
303  // All parcels introduced at SOI
304  if ((0.0 >= time0) && (0.0 < time1))
305  {
306  return this->volumeTotal_;
307  }
308 
309  return 0.0;
310 }
311 
312 
313 template<class CloudType>
315 (
316  const label parcelI,
317  const label nParcels,
318  const scalar time,
319  vector& position,
320  label& cellOwner,
321  label& tetFacei,
322  label& tetPti
323 )
324 {
325  position = positions_[parcelI];
326  cellOwner = injectorCells_[parcelI];
327  tetFacei = injectorTetFaces_[parcelI];
328  tetPti = injectorTetPts_[parcelI];
329 }
330 
331 
332 template<class CloudType>
334 (
335  const label parcelI,
336  const label,
337  const scalar,
338  typename CloudType::parcelType& parcel
339 )
340 {
341  // set particle velocity
342  parcel.U() = U0_;
343 
344  // set particle diameter
345  parcel.d() = diameters_[parcelI];
346 }
347 
348 
349 template<class CloudType>
351 {
352  return false;
353 }
354 
355 
356 template<class CloudType>
358 {
359  return true;
360 }
361 
362 
363 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
CellZoneInjection.H
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::scalarList
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:64
mathematicalConstants.H
Foam::CellZoneInjection::setPositionAndCell
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.
Definition: CellZoneInjection.C:315
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::CellZoneInjection
Injection positions specified by a particle number density within a cell set.
Definition: CellZoneInjection.H:65
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::meshTools::writeOBJ
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:203
globalIndex.H
Foam::InjectionModel
Templated injection model class.
Definition: InjectionModel.H:73
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
polyMeshTetDecomposition.H
Foam::sumOp
Definition: ops.H:213
Foam::DSMCCloud::rndGen
Random & rndGen()
Return reference to the random object.
Definition: DSMCCloudI.H:124
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::diff
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:378
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::DSMCCloud::mesh
const fvMesh & mesh() const
Return reference to the mesh.
Definition: DSMCCloudI.H:44
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::CellZoneInjection::timeEnd
scalar timeEnd() const
Return the end-of-injection time.
Definition: CellZoneInjection.C:273
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::radiation::lookup
Lookup type of boundary radiation properties.
Definition: lookup.H:63
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:71
Foam::CellZoneInjection::parcelsToInject
label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
Definition: CellZoneInjection.C:282
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::CellZoneInjection::~CellZoneInjection
virtual ~CellZoneInjection()
Destructor.
Definition: CellZoneInjection.C:217
Pstream.H
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam::CellZoneInjection::CellZoneInjection
CellZoneInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
Definition: CellZoneInjection.C:168
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::CellZoneInjection::setProperties
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
Definition: CellZoneInjection.C:334
Foam::CellZoneInjection::volumeToInject
scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
Definition: CellZoneInjection.C:298
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::CellZoneInjection::updateMesh
virtual void updateMesh()
Set injector locations when mesh is updated.
Definition: CellZoneInjection.C:224
Foam::Vector< scalar >
Foam::List< label >
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::CellZoneInjection::validInjection
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
Definition: CellZoneInjection.C:357
Foam::DSMCCloud::parcelType
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
Foam::CellZoneInjection::fullyDescribed
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
Definition: CellZoneInjection.C:350
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328