pointInjection.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 "pointInjection.H"
28 #include "LagrangianFields.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace Lagrangian
35 {
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 void Foam::Lagrangian::pointInjection::readCoeffs(const dictionary& modelDict)
45 {
46  point_.reset
47  (
49  (
50  "point",
51  mesh().time().userUnits(),
52  dimLength,
53  modelDict
54  ).ptr()
55  );
56 
57  numberRate_.reset
58  (
60  (
61  "numberRate",
62  mesh().time().userUnits(),
63  dimRate,
64  modelDict
65  ).ptr()
66  );
67 
68  numberDeferred_ = 0;
69 
70  injectionLocation_ = injectionLocation::unset;
71 }
72 
73 
74 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
75 
77 (
78  const word& name,
79  const LagrangianMesh& mesh,
80  const dictionary& modelDict,
81  const dictionary& stateDict
82 )
83 :
85  point_(nullptr),
86  numberRate_(nullptr),
87  numberDeferred_(stateDict.lookupOrDefault<scalar>("numberDeferred", 0)),
88  rndGen_("rndGen", stateDict, name, true),
89  timeIndex_(-1),
90  injectionLocation_(injectionLocation::unset),
91  coordinates_(barycentric::uniform(NaN)),
92  celli_(-1),
93  facei_(-1),
94  faceTrii_(-1)
95 {
96  readCoeffs(modelDict);
97 }
98 
99 
100 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
101 
103 {
104  // If anything is moving then all coordinates have to be constructed on the
105  // fly, so don't store anything
106  if (mesh().mesh().moving() || !point_->constant())
107  {
108  injectionLocation_ =
109  injectionLocation::multiplePoints;
110  }
111 
112  // If nothing is moving, and coordinates are not already stored, then
113  // construct and store the coordinates for repeated usage
114  if
115  (
116  injectionLocation_ == injectionLocation::unset
117  || injectionLocation_ == injectionLocation::multiplePoints
118  )
119  {
120  const scalar t1 = mesh().time().value();
121  const scalar t0 = t1 - mesh().time().deltaT().value();
122 
123  // Evaluate the position
124  const point p = point_->value(t0);
125 
126  // Locate within the mesh
127  const LagrangianMesh::location l =
128  mesh().locate(p, coordinates_, celli_, facei_, faceTrii_, 0);
129 
130  // Check for failure
131  checkLocation(l, p);
132 
133  // Set the location flag
134  injectionLocation_ =
135  celli_ >= 0
136  ? injectionLocation::fixedPointOnThisProcessor
137  : injectionLocation::fixedPointOnAnotherProcessor;
138  }
139 }
140 
141 
143 (
145  const LagrangianSubMesh&
146 ) const
147 {
148  const scalar t1 = mesh.time().value();
149  const scalar t0 = t1 - mesh.time().deltaT().value();
150 
151  // Restart the generator if necessary and set the time index up to date
152  rndGen_.start(timeIndex_ == db().time().timeIndex());
153  timeIndex_ = db().time().timeIndex();
154 
155  // Calculate the number of particles to inject. Round down to get an
156  // integer number. Store the excess to apply at a later time.
157  const scalar number = numberRate_->integral(t0, t1) + numberDeferred_;
158  const label numberInt = floor(number);
159  numberDeferred_ = number - numberInt;
160 
161  // Inject at random times throughout the time-step
162  scalarField fraction(rndGen_.scalar01(numberInt));
163 
164  // Create particles in the Lagrangian mesh
165  switch (injectionLocation_)
166  {
167  case injectionLocation::unset:
169  << "Injection location has not been set"
170  << exit(FatalError);
172 
173  case injectionLocation::fixedPointOnThisProcessor:
174  {
175  // Use stored coordinates to inject the particles
176  return
177  mesh.inject
178  (
179  *this,
180  barycentricField(numberInt, coordinates_),
181  labelField(numberInt, celli_),
182  labelField(numberInt, facei_),
183  labelField(numberInt, faceTrii_),
185  fraction
186  );
187  }
188 
189  case injectionLocation::fixedPointOnAnotherProcessor:
190  {
191  // Inject nothing. We need this because source conditions evaluated
192  // on the process containing the point may communicate; e.g., to
193  // sum the flow rate across the processes.
194  return
195  mesh.inject
196  (
197  *this,
199  labelField(),
200  labelField(),
201  labelField(),
203  scalarField()
204  );
205  }
206 
207  case injectionLocation::multiplePoints:
208  {
209  // Evaluate the variable positions
210  const pointField positions(point_->value(t0 + fraction*(t1 - t0)));
211 
212  // Locate within the moving mesh
213  barycentricField coordinates(numberInt);
214  labelField celli(numberInt), facei(numberInt), faceTrii(numberInt);
215  const List<LagrangianMesh::location> locations =
216  mesh.locate
217  (
218  positions,
219  coordinates,
220  celli,
221  facei,
222  faceTrii,
223  fraction
224  );
225 
226  // Check for any failures
227  checkLocation(locations, positions);
228 
229  // Remove particles not on this process
230  filter(coordinates, celli, facei, faceTrii, fraction);
231 
232  // Inject particles
233  return
234  mesh.inject
235  (
236  *this,
237  coordinates,
238  celli,
239  facei,
240  faceTrii,
242  fraction
243  );
244  }
245  }
246 
248 }
249 
250 
252 {
253  injectionLocation_ = injectionLocation::unset;
254 }
255 
256 
258 {
259  injectionLocation_ = injectionLocation::unset;
260 }
261 
262 
264 {
265  injectionLocation_ = injectionLocation::unset;
266 }
267 
268 
270 {
271  if (LagrangianInjection::read(modelDict))
272  {
273  readCoeffs(modelDict);
274  return true;
275  }
276  else
277  {
278  return false;
279  }
280 }
281 
282 
284 {
286 
287  writeEntry(os, "numberDeferred", numberDeferred_);
288  writeEntry(os, "rndGen", rndGen_);
289 }
290 
291 
292 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
static autoPtr< Function1< Type > > New(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Select from dictionary.
Definition: Function1New.C:32
Base class for Lagrangian injections. Minimal wrapper over LagrangianSource. Implements some utility ...
Class containing Lagrangian geometry and topology.
location
Enumeration for the locations of searched positions.
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...
pointInjection(const word &name, const LagrangianMesh &mesh, const dictionary &modelDict, const dictionary &stateDict)
Construct from components.
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.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual bool read(const dictionary &modelDict)
Read dictionary.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:46
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.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
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
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
bool writeState(const bool write) const
Write state.
Definition: stateModel.C:104
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 FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
addToRunTimeSelectionTable(LagrangianModel, constantCoefficientVirtualMass, dictionary)
defineTypeNameAndDebug(constantCoefficientVirtualMass, 0)
barycentric coordinates(const polyMesh &mesh, const point &position, const label celli, const label facei, const label faceTrii, const scalar stepFraction)
Return the coordinates given the position and tet topology.
Definition: tracking.C:1259
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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 dimRate
const dimensionSet dimLength
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:49
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
label timeIndex
Definition: getTimeIndex.H:4
volScalarField & p