LagrangianSubMesh.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 "LagrangianSubMesh.H"
27 #include "LagrangianSubFields.H"
28 #include "LagrangianMesh.H"
29 #include "tracking.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
36 }
37 
38 
39 // * * * * * * * * * * * * * * Private Constructors * * * * * * * * * * * * //
40 
41 Foam::LagrangianSubMesh::LagrangianSubMesh
42 (
43  const LagrangianMesh& mesh,
44  const LagrangianGroup group,
45  const label size,
46  const label start,
47  const label index
48 )
49 :
50  GeoMesh<polyMesh>(mesh.mesh()),
51  mesh_(mesh),
52  group_(group),
53  size_(size),
54  start_(start),
55  index_(index)
56 {}
57 
58 
59 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
60 
61 Foam::LagrangianSubMesh::LagrangianSubMesh
62 (
63  const LagrangianMesh& mesh,
64  const LagrangianGroup group,
65  const label size,
66  const label start
67 )
68 :
70  mesh_(mesh),
71  group_(group),
72  size_(size),
73  start_(start),
74  index_(mesh.subMeshIndex())
75 {}
76 
77 
78 Foam::LagrangianSubMesh::LagrangianSubMesh
79 (
80  const LagrangianMesh& mesh,
81  const labelList& groupOffsets,
83 )
84 :
86  mesh_(mesh),
87  group_(group),
88  size_
89  (
90  groupOffsets[static_cast<label>(group) + 1]
91  - groupOffsets[static_cast<label>(group)]
92  ),
93  start_(groupOffsets[static_cast<label>(group)]),
94  index_(mesh.subMeshIndex())
95 {}
96 
97 
98 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
99 
101 {}
102 
103 
104 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
105 
106 template<>
108 (
109  const LagrangianSubScalarSubField& fraction
110 )
111 {
112  const LagrangianSubMesh& subMesh = fraction.mesh();
113 
114  tmp<vectorField> tnf(new vectorField(subMesh.size()));
115  vectorField& nf = tnf.ref();
116 
117  forAll(subMesh, subi)
118  {
119  nf[subi] =
121  (
122  subMesh.mesh().mesh(),
123  subMesh.mesh().coordinates()[subi + subMesh.start()],
124  subMesh.mesh().celli()[subi + subMesh.start()],
125  subMesh.mesh().facei()[subi + subMesh.start()],
126  subMesh.mesh().faceTrii()[subi + subMesh.start()],
127  fraction[subi]
128  ).first();
129  }
130 
131  return tnf;
132 }
133 
134 
135 template<>
137 (
138  const LagrangianSubScalarSubField& fraction
139 )
140 {
141  const LagrangianSubMesh& subMesh = fraction.mesh();
142 
143  return
145  (
146  "nf:" + Foam::name(subMesh.group()),
147  subMesh,
148  dimless,
149  nf<vectorField>(fraction)
150  );
151 }
152 
153 
154 template<>
156 (
158 ) const
159 {
160  tmp<vectorField> tnf(new vectorField(size()));
161  vectorField& nf = tnf.ref();
162 
163  forAll(*this, subi)
164  {
165  nf[subi] =
167  (
168  mesh().mesh(),
169  mesh().coordinates()[subi + start()],
170  mesh().celli()[subi + start()],
171  mesh().facei()[subi + start()],
172  mesh().faceTrii()[subi + start()],
173  fraction[subi + start()]
174  ).first();
175  }
176 
177  return tnf;
178 }
179 
180 
181 template<>
183 (
185 ) const
186 {
187  return
189  (
190  "nf:" + Foam::name(group()),
191  *this,
192  dimless,
193  nf<vectorField>(fraction)
194  );
195 }
196 
197 
198 template<>
200 (
201  const LagrangianSubScalarSubField& fraction
202 )
203 {
204  const LagrangianSubMesh& subMesh = fraction.mesh();
205 
206  tmp<vectorField> tUf(new vectorField(subMesh.size()));
207  vectorField& Uf = tUf.ref();
208 
209  forAll(subMesh, subi)
210  {
211  Uf[subi] =
213  (
214  subMesh.mesh().mesh(),
215  subMesh.mesh().coordinates()[subi + subMesh.start()],
216  subMesh.mesh().celli()[subi + subMesh.start()],
217  subMesh.mesh().facei()[subi + subMesh.start()],
218  subMesh.mesh().faceTrii()[subi + subMesh.start()],
219  fraction[subi]
220  ).second()/subMesh.mesh().time().deltaTValue();
221  }
222 
223  return tUf;
224 }
225 
226 
227 template<>
229 (
230  const LagrangianSubScalarSubField& fraction
231 )
232 {
233  const LagrangianSubMesh& subMesh = fraction.mesh();
234 
235  return
237  (
238  "Uf:" + Foam::name(subMesh.group()),
239  subMesh,
240  dimVelocity,
241  Uf<vectorField>(fraction)
242  );
243 }
244 
245 
246 template<>
248 (
250 ) const
251 {
252  tmp<vectorField> tUf(new vectorField(size()));
253  vectorField& Uf = tUf.ref();
254 
255  forAll(*this, subi)
256  {
257  Uf[subi] =
259  (
260  mesh().mesh(),
261  mesh().coordinates()[subi + start()],
262  mesh().celli()[subi + start()],
263  mesh().facei()[subi + start()],
264  mesh().faceTrii()[subi + start()],
265  fraction[subi + start()]
266  ).second()/mesh().time().deltaTValue();
267  }
268 
269  return tUf;
270 }
271 
272 
273 template<>
275 (
277 ) const
278 {
279  return
281  (
282  "Uf:" + Foam::name(group()),
283  *this,
284  dimVelocity,
285  Uf<vectorField>(fraction)
286  );
287 }
288 
289 
290 // * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * * //
291 
293 {
294  if (&mesh_ != &subMesh.mesh_)
295  {
297  << "Cannot combine sub-meshes which relate to different meshes"
298  << exit(FatalError);
299  }
300 
301  if (group_ != subMesh.group_)
302  {
304  << "Cannot combine sub-meshes with different groups "
305  << exit(FatalError);
306  }
307 
308  if (size_ + start_ != subMesh.start_)
309  {
311  << "Cannot combine sub-meshes that are not contiguous"
312  << exit(FatalError);
313  }
314 
315  size_ += subMesh.size_;
316 }
317 
318 
319 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Mesh & mesh() const
Return mesh.
static tmp< DimensionedField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Mesh &mesh, const dimensionSet &, const PrimitiveField< Type > &)
Return a temporary field constructed from name, mesh,.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:47
Class containing Lagrangian geometry and topology.
const Time & time() const
Return time.
const labelIODynamicField & faceTrii() const
Access the face-tet indices.
const polyMesh & mesh() const
Access the mesh.
const labelIODynamicField & celli() const
Access the cell indices.
const barycentricIODynamicField & coordinates() const
Access the coordinates.
const labelIODynamicField & facei() const
Access the cell-face indices.
Mesh that relates to a sub-section of a Lagrangian mesh. This is used to construct fields that relate...
LagrangianGroup group() const
Return the group.
static tmp< FieldType > nf(const LagrangianSubScalarSubField &fraction)
Return the face normals at the Lagrangian locations.
label size() const
Return size.
void operator+=(const LagrangianSubMesh &)
Add a sub-mesh to this one. Must relate to adjacent elements.
static tmp< FieldType > Uf(const LagrangianSubScalarSubField &fraction)
Return the face velocities at the Lagrangian locations.
const LagrangianMesh & mesh() const
Return the mesh.
label start() const
Return start.
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:34
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
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
const char *const group
Group name for atomic constants.
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
Pair< vector > faceNormalAndDisplacement(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label facei, const label faceTrii, const scalar stepFraction)
Return the normal of the corresponding point on the associated face and.
Definition: tracking.C:1302
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 dimless
defineTypeNameAndDebug(combustionModel, 0)
Field< vector > vectorField
Specialisation of Field<T> for vector.
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimVelocity
LagrangianGroup
Lagrangian group enumeration.
error FatalError