volumeInjection.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "volumeInjection.H"
27 #include "CompactListList.H"
28 #include "LagrangianFields.H"
29 #include "tetIndices.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace Lagrangian
37 {
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
46 void Foam::Lagrangian::volumeInjection::readCoeffs(const dictionary& modelDict)
47 {
48  zone_.read(modelDict);
49 
50  haveNumber_ = modelDict.found("number");
51  const bool haveNumberDensity = modelDict.found("numberDensity");
52 
53  if (haveNumber_ == haveNumberDensity)
54  {
55  FatalIOErrorInFunction(modelDict)
56  << "keywords number and numberDensity are both "
57  << (haveNumber_ ? "" : "un") << "defined in "
58  << "dictionary " << modelDict.name()
59  << exit(FatalIOError);
60  }
61 
62  numberOrNumberDensity_ =
63  haveNumber_
64  ? modelDict.lookup<scalar>("number", dimless)
65  : modelDict.lookup<scalar>("numberDensity", dimless/dimVolume);
66 
67  time_ =
68  modelDict.lookupOrDefault<scalar>
69  (
70  "time",
71  mesh().time().userUnits(),
72  mesh().time().beginTime().value()
73  );
74 }
75 
76 
77 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
78 
80 (
81  const word& name,
82  const LagrangianMesh& mesh,
83  const dictionary& modelDict,
84  const dictionary& stateDict
85 )
86 :
88  zone_(mesh.mesh()),
89  haveNumber_(false),
90  numberOrNumberDensity_(NaN),
91  time_(NaN),
92  globalRndGen_("globalRndGen", stateDict, name, true),
93  localRndGen_("localRndGen", stateDict, name, false),
94  timeIndex_(-1)
95 {
96  readCoeffs(modelDict);
97 }
98 
99 
100 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
101 
103 {
104  // In theory the set can change as a result of motion. We can't handle this
105  // entirely correctly as we can't construct the set at an intermediate
106  // time. We assume that the most likely scenario is injection at the
107  // beginning of the time-step. So, when we are inside the injection
108  // time-step we delay update of the set until after injection.
109 
110  const scalar t1 = mesh().time().value();
111  const scalar t0 = t1 - mesh().time().deltaT().value();
112 
113  if (t0 <= time_ && time_ < t1) return;
114 
115  if (mesh().mesh().moving())
116  {
117  zone_.movePoints();
118  }
119 }
120 
121 
123 (
125  const LagrangianSubMesh&
126 ) const
127 {
128  const scalar t1 = mesh.time().value();
129  const scalar t0 = t1 - mesh.time().deltaT().value();
130 
131  if (!(t0 <= time_ && time_ < t1)) return mesh.subNone();
132 
133  const scalar fraction = (time_ - t0)/(t1 - t0);
134 
135  // Restart the generators if necessary and set the time index up to date
136  localRndGen_.start(timeIndex_ == db().time().timeIndex());
137  globalRndGen_.start(timeIndex_ == db().time().timeIndex());
138  timeIndex_ = db().time().timeIndex();
139 
140  // Reference the mesh faces and the set cells
141  const faceList& faces = mesh.mesh().faces();
142  const labelList& setCellCells = zone_.zone();
143  const UIndirectList<cell> setCells(mesh.mesh().cells(), setCellCells);
144 
145  // Create point and cell-centre fields at the current fraction
146  const pointField points
147  (
148  (1 - fraction)*mesh.mesh().oldPoints()
149  + fraction*mesh.mesh().points()
150  );
151  const pointField cellCentres
152  (
153  (1 - fraction)*mesh.mesh().oldCellCentres()
154  + fraction*mesh.mesh().cellCentres()
155  );
156 
157  // Count the number of tetrahedra in each set cell
158  labelList setCellNTets(setCells.size(), label(0));
159  forAll(setCells, setCelli)
160  {
161  forAll(setCells[setCelli], cellFacei)
162  {
163  setCellNTets[setCelli] +=
164  faces[setCells[setCelli][cellFacei]].nTriangles();
165  }
166  }
167 
168  // Construct indexing for the cell-tetrahedra and cumulative sums of the
169  // volumes in the cells and cell-tetrahedra
170  CompactListList<labelPair> setCellTetFaceAndFaceTri
171  (
172  setCellNTets,
173  labelPair(-1, -1)
174  );
175  scalarList setCellSumVolume(setCells.size(), scalar(0));
176  CompactListList<scalar> setCellTetSumVolume(setCellNTets, scalar(0));
177  forAll(setCells, setCelli)
178  {
179  const cell& c = setCells[setCelli];
180 
181  label cellTeti = 0;
182 
183  forAll(c, cFacei)
184  {
185  const face& f = faces[c[cFacei]];
186 
187  for (label fTrii = 1; fTrii <= f.nTriangles(); ++ fTrii)
188  {
189  setCellTetFaceAndFaceTri[setCelli][cellTeti] =
190  labelPair(c[cFacei], fTrii);
191 
192  setCellSumVolume[setCelli] +=
193  tetIndices(setCellCells[setCelli], c[cFacei], fTrii)
194  .tet(mesh.mesh(), points, cellCentres)
195  .mag();
196 
197  setCellTetSumVolume[setCelli][cellTeti] =
198  setCellSumVolume[setCelli];
199 
200  cellTeti ++;
201  }
202  }
203  }
204  for (label setCelli = 1; setCelli < setCells.size(); ++ setCelli)
205  {
206  setCellSumVolume[setCelli] += setCellSumVolume[setCelli - 1];
207  }
208 
209  // Construct a cumulative sum of the volumes across the processes
210  scalarList procSumVolume(Pstream::nProcs(), -vGreat);
211  procSumVolume[Pstream::myProcNo()] =
212  setCells.size() ? setCellSumVolume.last() : scalar(0);
214  Pstream::listCombineScatter(procSumVolume);
215  for (label proci = 1; proci < Pstream::nProcs(); proci ++)
216  {
217  procSumVolume[proci] += procSumVolume[proci - 1];
218  }
219 
220  // Determine the number of particles to inject
221  const label numberInt =
222  haveNumber_
223  ? round(numberOrNumberDensity_)
224  : round(numberOrNumberDensity_*procSumVolume.last());
225 
226  // Binary search for a given volume within a supplied cumulative volume sum.
227  // Once the interval has been found, subtract the lower value from that
228  // volume so that a subsequent binary search can be done on a sub-set of
229  // volumes. Note that the supplied cumulative volume sum does not include
230  // the leading zero.
231  auto findVolume = []
232  (
233  const scalarUList& volumes,
234  scalar& volume
235  )
236  {
237  label i0 = -1, i1 = volumes.size() - 1;
238 
239  while (i0 + 1 != i1)
240  {
241  const label i = (i0 + i1)/2;
242 
243  (volume < volumes[i] ? i1 : i0) = i;
244  }
245 
246  if (i0 != -1)
247  {
248  volume -= volumes[i0];
249  }
250 
251  return i1;
252  };
253 
254  // Initialise storage for the injection geometry and topology. This is
255  // dynamic as we don't know how much will end up on each processor yet.
256  DynamicList<barycentric> injectCoordinates(numberInt/Pstream::nProcs());
257  DynamicList<label> injectCells(numberInt/Pstream::nProcs());
258  DynamicList<label> injectFaces(numberInt/Pstream::nProcs());
259  DynamicList<label> injectFaceTris(numberInt/Pstream::nProcs());
260 
261  // Create a (global) list of volumes at which to inject. Each volume is
262  // searched for in the cumulative lists to identify a tetrahedron into
263  // which to inject.
264  scalarField volume(globalRndGen_.scalar01(numberInt)*procSumVolume.last());
265  forAll(volume, volumei)
266  {
267  const label proci = findVolume(procSumVolume, volume[volumei]);
268 
269  if (Pstream::myProcNo() == proci)
270  {
271  const label setCelli =
272  findVolume(setCellSumVolume, volume[volumei]);
273  const label celli = setCellCells[setCelli];
274  const label cellTeti =
275  findVolume(setCellTetSumVolume[setCelli], volume[volumei]);
276 
277  injectCoordinates.append(barycentric01(localRndGen_));
278  injectCells.append(celli);
279  injectFaces.append
280  (
281  setCellTetFaceAndFaceTri[setCelli][cellTeti].first()
282  );
283  injectFaceTris.append
284  (
285  setCellTetFaceAndFaceTri[setCelli][cellTeti].second()
286  );
287  }
288  }
289 
290  // Create the particles at the identified locations
291  LagrangianSubMesh injectionMesh =
292  mesh.inject
293  (
294  *this,
295  barycentricField(injectCoordinates),
296  labelField(injectCells),
297  labelField(injectFaces),
298  labelField(injectFaceTris),
300  scalarField(injectCoordinates.size(), fraction)
301  );
302 
303  // See the note in volumeInjection::correct
304  if (mesh.mesh().moving())
305  {
306  zone_.movePoints();
307  }
308 
309  return injectionMesh;
310 }
311 
312 
314 {
315  zone_.topoChange(map);
316 }
317 
318 
320 {
321  zone_.mapMesh(map);
322 }
323 
324 
326 (
327  const polyDistributionMap& map
328 )
329 {
330  zone_.distribute(map);
331 }
332 
333 
335 {
336  if (LagrangianInjection::read(modelDict))
337  {
338  readCoeffs(modelDict);
339  return true;
340  }
341  else
342  {
343  return false;
344  }
345 }
346 
347 
349 {
351 
352  writeEntry(os, "globalRndGen", globalRndGen_);
353 }
354 
355 
357 {
359 
360  writeEntry(os, "localRndGen", localRndGen_);
361 }
362 
363 
364 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
A packed storage unstructured matrix of objects of type <T> using an offset table for access.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
Base class for Lagrangian injections. Minimal wrapper over LagrangianSource. Implements some utility ...
Class containing Lagrangian geometry and topology.
const Time & time() const
Return time.
static const word fractionName
Name of the tracked fraction field.
Base class for Lagrangian models.
virtual bool read(const dictionary &modelDict)
Read dictionary.
const LagrangianMesh & mesh() const
The mesh.
Mesh that relates to a sub-section of a Lagrangian mesh. This is used to construct fields that relate...
virtual void correct()
Correct the LagrangianModel.
virtual void writeState(Ostream &os) const
Write state.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
LagrangianSubMesh modify(LagrangianMesh &mesh, const LagrangianSubMesh &) const
Create new elements in the Lagrangian mesh.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
volumeInjection(const word &name, const LagrangianMesh &mesh, const dictionary &modelDict, const dictionary &stateDict)
Construct from components.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual void writeProcessorState(Ostream &os) const
Write state.
virtual bool read(const dictionary &modelDict)
Read dictionary.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:28
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:46
const unitConversion & userUnits() const
Return the user-time unit conversion.
Definition: Time.C:858
virtual dimensionedScalar beginTime() const
Return begin time (initial start time)
Definition: Time.C:780
A List with indirect addressing.
Definition: UIndirectList.H:60
label size() const
Return the number of elements in the list.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: UList.H:74
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
T & last()
Return the last element of the list.
Definition: UListI.H:128
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:60
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
const Type & value() const
Return const reference to value.
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:1244
const polyMesh & mesh() const
Return reference to polyMesh.
Definition: fvMesh.H:443
bool read(const dictionary &dict)
Read coefficients dictionary.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1344
virtual const pointField & oldPoints() const
Return old points for mesh motion.
Definition: polyMesh.C:1369
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1331
virtual const pointField & oldCellCentres() const
Return old cell centres for mesh motion.
Definition: polyMesh.C:1387
bool moving() const
Is mesh moving.
Definition: polyMesh.H:473
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const vectorField & cellCentres() const
const cellList & cells() const
bool writeState(const bool write) const
Write state.
Definition: stateModel.C:104
virtual void writeProcessorState(Ostream &os) const
Write processor state.
Definition: stateModel.C:132
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:82
tetPointRef tet(const polyMesh &mesh, const pointField &meshPoints, const pointField &cellCentres) const
Return the geometry corresponding to this tet and the given.
Definition: tetIndicesI.H:109
scalar mag() const
Return volume.
Definition: tetrahedronI.H:170
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
const pointField & points
addToRunTimeSelectionTable(LagrangianModel, constantCoefficientVirtualMass, dictionary)
defineTypeNameAndDebug(constantCoefficientVirtualMass, 0)
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
Field< barycentric > barycentricField
Barycentric field.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const dimensionSet dimless
barycentric barycentric01(randomGenerator &rndGen)
Generate a random barycentric coordinate within the unit tetrahedron.
Definition: barycentric.C:31
labelList second(const UList< labelPair > &p)
Definition: patchToPatch.C:49
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const dimensionSet dimVolume
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:49
IOerror FatalIOError
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
labelList f(nPoints)