manualInjection.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 "manualInjection.H"
27 #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::manualInjection::readCoeffs(const dictionary& modelDict)
45 {
46  fileName_ = modelDict.lookup<fileName>("file");
47 
48  units_.readIfPresent("units", modelDict);
49 
50  time_ =
51  modelDict.lookupOrDefault<scalar>
52  (
53  "time",
54  mesh().time().userUnits(),
55  mesh().time().beginTime().value()
56  );
57 }
58 
59 
60 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
61 
63 (
64  const word& name,
65  const LagrangianMesh& mesh,
66  const dictionary& modelDict,
67  const dictionary& stateDict
68 )
69 :
71  fileName_(fileName::null),
72  units_(dimLength),
73  time_(NaN)
74 {
75  readCoeffs(modelDict);
76 }
77 
78 
79 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
80 
82 (
84  const LagrangianSubMesh&
85 ) const
86 {
87  const scalar t1 = mesh.time().value();
88  const scalar t0 = t1 - mesh.time().deltaT().value();
89 
90  if (!(t0 <= time_ && time_ < t1)) return mesh.subNone();
91 
92  // Read the positions file
93  GlobalIOField<point> positions
94  (
95  IOobject
96  (
97  fileName_,
98  mesh.time().constant(),
99  mesh,
102  )
103  );
104 
105  units_.makeStandard(positions);
106 
107  const label number = positions.size();
108 
109  scalarField fraction(number, (time_ - t0)/(t1 - t0));
110 
111  // Locate within the mesh
113  labelField celli(number, -1), facei(number), faceTrii(number);
114  const List<LagrangianMesh::location> locations =
115  mesh.locate
116  (
117  positions,
118  coordinates,
119  celli,
120  facei,
121  faceTrii,
122  fraction
123  );
124 
125  // Check for any failures
126  checkLocation(locations, positions);
127 
128  // Remove particles not on this process
129  filter(coordinates, celli, facei, faceTrii, fraction);
130 
131  // Inject particles
132  return
133  mesh.inject
134  (
135  *this,
136  coordinates,
137  celli,
138  facei,
139  faceTrii,
141  fraction
142  );
143 }
144 
145 
147 {
148  if (LagrangianInjection::read(modelDict))
149  {
150  readCoeffs(modelDict);
151  return true;
152  }
153  else
154  {
155  return false;
156  }
157 }
158 
159 
160 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
A primitive field of type <Type> with automated input and output.
Definition: GlobalIOField.H:50
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
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...
LagrangianSubMesh modify(LagrangianMesh &mesh, const LagrangianSubMesh &) const
Create new elements in the Lagrangian mesh.
virtual bool read(const dictionary &modelDict)
Read dictionary.
manualInjection(const word &name, const LagrangianMesh &mesh, const dictionary &modelDict, const dictionary &stateDict)
Construct from components.
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
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
static const word & constant()
Return constant name.
Definition: TimePaths.H:122
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 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 class for handling file names.
Definition: fileName.H:82
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
bool readIfPresent(const word &keyword, const dictionary &)
Update if found in the dictionary.
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)
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.
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 dimLength
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.