diskInjection.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 "LagrangianMesh.H"
27 #include "diskInjection.H"
29 #include "LagrangianFields.H"
30 #include "mathematicalConstants.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::diskInjection::readCoeffs(const dictionary& modelDict)
47 {
48  centre_.reset
49  (
51  (
52  "centre",
53  mesh().time().userUnits(),
54  dimLength,
55  modelDict
56  ).ptr()
57  );
58 
59  axis_.reset
60  (
62  (
63  "axis",
64  mesh().time().userUnits(),
65  dimless,
66  modelDict
67  ).ptr()
68  );
69 
70  const bool haveDiemeter = modelDict.found("diameter");
71  const bool haveInnerDiemeter = modelDict.found("innerDiameter");
72  const bool haveOuterDiemeter = modelDict.found("outerDiameter");
73 
74  if (haveDiemeter == (haveInnerDiemeter || haveOuterDiemeter))
75  {
76  FatalIOErrorInFunction(modelDict)
77  << "keywords diameter and innerDiameter/outerDiameter are both "
78  << (haveDiemeter ? "" : "un") << "defined in "
79  << "dictionary " << modelDict.name()
80  << exit(FatalIOError);
81  }
82 
83  if (haveInnerDiemeter != haveOuterDiemeter)
84  {
85  FatalIOErrorInFunction(modelDict)
86  << "keywords innerDiameter and outerDiameter are not both defined "
87  << "in dictionary " << modelDict.name()
88  << exit(FatalIOError);
89  }
90 
91  if (haveDiemeter)
92  {
93  innerDiameter_ = 0;
94  outerDiameter_ = modelDict.lookup<scalar>("diameter", dimLength);
95  }
96  else
97  {
98  innerDiameter_ = modelDict.lookup<scalar>("innerDiameter", dimLength);
99  outerDiameter_ = modelDict.lookup<scalar>("outerDiameter", dimLength);
100  }
101 
102  numberRate_.reset
103  (
105  (
106  "numberRate",
107  mesh().time().userUnits(),
108  dimRate,
109  modelDict
110  ).ptr()
111  );
112 
113  numberDeferred_ = 0;
114 }
115 
116 
117 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
118 
120 (
121  const word& name,
122  const LagrangianMesh& mesh,
123  const dictionary& modelDict,
124  const dictionary& stateDict
125 )
126 :
128  centre_(nullptr),
129  axis_(nullptr),
130  innerDiameter_(NaN),
131  outerDiameter_(NaN),
132  numberRate_(nullptr),
133  numberDeferred_(stateDict.lookupOrDefault<scalar>("numberDeferred", 0)),
134  rndGen_("rndGen", stateDict, name, true),
135  timeIndex_(-1)
136 {
137  readCoeffs(modelDict);
138 }
139 
140 
141 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
142 
144 {
145  if (!rFracPtr_.valid())
146  {
148  << "Radius fraction requested outside of the injection"
149  << exit(FatalError);
150  }
151 
152  return rFracPtr_();
153 }
154 
155 
157 {
158  if (!axisPtr_.valid())
159  {
161  << "Axis requested outside of the injection"
162  << exit(FatalError);
163  }
164 
165  return axisPtr_();
166 }
167 
168 
170 {
171  if (!radialPtr_.valid())
172  {
174  << "Axis requested outside of the injection"
175  << exit(FatalError);
176  }
177 
178  return radialPtr_();
179 }
180 
181 
183 (
185  const LagrangianSubMesh&
186 ) const
187 {
188  const scalar t1 = mesh.time().value();
189  const scalar t0 = t1 - mesh.time().deltaT().value();
190 
191  // Restart the generator if necessary and set the time index up to date
192  rndGen_.start(timeIndex_ == db().time().timeIndex());
193  timeIndex_ = db().time().timeIndex();
194 
195  // Calculate the number of particles to inject. Round down to get an
196  // integer number. Store the excess to apply at a later time.
197  const scalar number = numberRate_->integral(t0, t1) + numberDeferred_;
198  const label numberInt = floor(number);
199  numberDeferred_ = number - numberInt;
200 
201  // Inject at random times throughout the time-step
202  scalarField fraction(rndGen_.scalar01(numberInt));
203 
204  // Evaluate the variable centre and axis, and create radial vectors to
205  // complete the local coordinate system
206  const pointField centre(centre_->value(t0 + fraction*(t1 - t0)));
207  axisPtr_.set(normalised(axis_->value(t0 + fraction*(t1 - t0))).ptr());
208  const vectorField radial1(normalised(perpendicular(axis())));
209  const vectorField radial2(axis() ^ radial1);
210 
211  // Create random radii within and angles around the disk. Store temporarily
212  // so that source conditions can use it.
213  rFracPtr_.set(rndGen_.scalar01(numberInt).ptr());
214  const scalarField r
215  (
216  sqrt
217  (
218  (1 - rFrac())*sqr(innerDiameter_/2)
219  + rFrac()*sqr(outerDiameter_/2)
220  )
221  );
222  const scalarField phi
223  (
224  constant::mathematical::twoPi*rndGen_.scalar01(numberInt)
225  );
226 
227  // Evaluate the radial vector
228  radialPtr_.set((cos(phi)*radial1 + sin(phi)*radial2).ptr());
229 
230  // Evaluate the positions
231  const pointField positions(centre + r*radial());
232 
233  // Locate within the mesh
235  labelField celli(number, -1), facei(number), faceTrii(number);
236  const List<LagrangianMesh::location> locations =
237  mesh.locate
238  (
239  positions,
240  coordinates,
241  celli,
242  facei,
243  faceTrii,
244  fraction
245  );
246 
247  // Check for any failures
248  checkLocation(locations, positions);
249 
250  // Remove particles not on this process
251  filter(coordinates, celli, facei, faceTrii, fraction);
252 
253  // Inject particles
254  LagrangianSubMesh injectionMesh =
255  mesh.inject
256  (
257  *this,
258  coordinates,
259  celli,
260  facei,
261  faceTrii,
263  fraction
264  );
265 
266  // Clean up
267  rFracPtr_.clear();
268  axisPtr_.clear();
269  radialPtr_.clear();
270 
271  return injectionMesh;
272 }
273 
274 
276 {
277  if (LagrangianInjection::read(modelDict))
278  {
279  readCoeffs(modelDict);
280  return true;
281  }
282  else
283  {
284  return false;
285  }
286 }
287 
288 
290 {
292 
293  writeEntry(os, "numberDeferred", numberDeferred_);
294  writeEntry(os, "rndGen", rndGen_);
295 }
296 
297 
298 // ************************************************************************* //
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.
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 writeState(Ostream &os) const
Write state.
LagrangianSubMesh modify(LagrangianMesh &mesh, const LagrangianSubMesh &) const
Create new elements in the Lagrangian mesh.
diskInjection(const word &name, const LagrangianMesh &mesh, const dictionary &modelDict, const dictionary &stateDict)
Construct from components.
const scalarField & rFrac() const
Access the cached radius fractions. Only valid during injection.
const vectorField & radial() const
Access the cached radial vectors. Only valid during injection.
virtual bool read(const dictionary &modelDict)
Read dictionary.
const vectorField & axis() const
Access the cached axes. Only valid during injection.
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
const unitConversion & userUnits() const
Return the user-time unit conversion.
Definition: Time.C:858
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
void clear()
Remove all regIOobject owned by the registry.
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 FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
addToRunTimeSelectionTable(LagrangianModel, constantCoefficientVirtualMass, dictionary)
defineTypeNameAndDebug(constantCoefficientVirtualMass, 0)
const scalar twoPi(2 *pi)
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
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 dimless
dimensionedScalar sin(const dimensionedScalar &ds)
const dimensionSet dimLength
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionSet normalised(const dimensionSet &)
Definition: dimensionSet.C:507
IOerror FatalIOError
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
dimensionSet perpendicular(const dimensionSet &)
Definition: dimensionSet.C:513
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dimensionedScalar cos(const dimensionedScalar &ds)
label timeIndex
Definition: getTimeIndex.H:4