LagrangianSubMesh.H
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 Class
25  Foam::LagrangianSubMesh
26 
27 Description
28  Mesh that relates to a sub-section of a Lagrangian mesh. This is used to
29  construct fields that relate to a contiguous sub-set of the Lagrangian
30  elements. This class only stores references and the indices defining the
31  range of the sub-set, so it is very lightweight and can be constructed and
32  thrown away largely without consideration of expense.
33 
34 SourceFiles
35  LagrangianSubMesh.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef LagrangianSubMesh_H
40 #define LagrangianSubMesh_H
41 
42 #include "GeoMesh.H"
43 #include "LagrangianState.H"
44 #include "LagrangianFieldsFwd.H"
45 #include "LagrangianSubFieldsFwd.H"
46 #include "polyMesh.H"
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 namespace Foam
51 {
52 
53 class LagrangianMesh;
54 
55 /*---------------------------------------------------------------------------*\
56  Class LagrangianSubMesh Declaration
57 \*---------------------------------------------------------------------------*/
58 
60 :
61  public GeoMesh<polyMesh>
62 {
63  // Private Data
64 
65  //- Reference to the Lagrangian mesh
66  const LagrangianMesh& mesh_;
67 
68  //- Group of the elements in this subset
69  const LagrangianGroup group_;
70 
71  //- Size of the subset
72  label size_;
73 
74  //- Start index of the subset
75  label start_;
76 
77  //- Index with which to compare instances
78  const label index_;
79 
80 
81  // Private Constructors
82 
83  //- Construct from components
84  explicit LagrangianSubMesh
85  (
86  const LagrangianMesh& mesh,
87  const LagrangianGroup group,
88  const label size,
89  const label start,
90  const label index
91  );
92 
93 
94 public:
95 
96  // Friend classes
97 
98  //- Allow the Lagrangian mesh to construct with a specified index
99  friend class LagrangianMesh;
100 
101 
102  //- Runtime type information
103  ClassName("LagrangianSubMesh");
104 
105 
106  // Public Type Definitions
107 
108  //- Mesh type
109  typedef LagrangianSubMesh Mesh;
110 
111 
112  // Constructors
113 
114  //- Construct from components, except for the index which is
115  // automatically generated from the complete mesh
116  explicit LagrangianSubMesh
117  (
118  const LagrangianMesh& mesh,
119  const LagrangianGroup group,
120  const label size,
121  const label start
122  );
123 
124  //- Construct from Lagrangian mesh, group offsets and group
125  explicit LagrangianSubMesh
126  (
127  const LagrangianMesh& mesh,
128  const labelList& groupOffsets,
129  const LagrangianGroup group
130  );
131 
132 
133  //- Destructor
135 
136 
137  // Member functions
138 
139  // Access
140 
141  //- Return the mesh
142  inline const LagrangianMesh& mesh() const;
143 
144  //- Return the group
145  inline LagrangianGroup group() const;
146 
147  //- Return size
148  inline label size() const;
149 
150  //- Return size
151  inline label globalSize() const;
152 
153  //- Return size
154  static inline label size(const LagrangianSubMesh& subMesh);
155 
156  //- Return whether or not the mesh is empty
157  inline bool empty() const;
158 
159  //- Return start
160  inline label start() const;
161 
162  //- Return end
163  inline label end() const;
164 
165  //- Return the index
166  inline label index() const;
167 
168 
169  // Sub-setting
170 
171  //- Return a sub-list corresponding to this sub-mesh
172  template<class Type>
173  SubList<Type> sub(const List<Type>& list) const;
174 
175  //- Return a sub-field corresponding to this sub-mesh
176  template<class Type>
177  SubField<Type> sub(const Field<Type>& field) const;
178 
179  //- Return a sub-dimensioned-field corresponding to this sub-mesh
180  template<class Type, template<class> class PrimitiveField>
182  (
184  ) const;
185 
186 
187  // Geometry
188 
189  //- Return the face normals at the Lagrangian locations
190  template<class FieldType>
191  static tmp<FieldType> nf
192  (
193  const LagrangianSubScalarSubField& fraction
194  );
195 
196  //- Return the face normals at the Lagrangian locations
197  template<class FieldType>
199  (
201  ) const;
202 
203  //- Return the face velocities at the Lagrangian locations
204  template<class FieldType>
205  static tmp<FieldType> Uf
206  (
207  const LagrangianSubScalarSubField& fraction
208  );
209 
210  //- Return the face velocities at the Lagrangian locations
211  template<class FieldType>
213  (
215  ) const;
216 
217 
218  // Member Operators
219 
220  //- Add a sub-mesh to this one. Must relate to adjacent elements.
221  void operator+=(const LagrangianSubMesh&);
222 };
223 
224 
225 // Template Specialisations
226 
227 template<>
229 (
230  const LagrangianSubScalarSubField& fraction
231 );
232 
233 template<>
235 (
236  const LagrangianSubScalarSubField& fraction
237 );
238 
239 template<>
241 (
243 ) const;
244 
245 template<>
247 (
249 ) const;
250 
251 template<>
253 (
254  const LagrangianSubScalarSubField& fraction
255 );
256 
257 template<>
259 (
260  const LagrangianSubScalarSubField& fraction
261 );
262 
263 template<>
265 (
267 ) const;
268 
269 template<>
271 (
273 ) const;
274 
275 
276 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
277 
278 } // End namespace Foam
279 
280 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
281 
282 #include "LagrangianSubMeshI.H"
283 
284 #ifdef NoRepository
286 #endif
287 
288 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
289 
290 #endif
291 
292 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Pre-declare SubField and related Field type.
Definition: Field.H:83
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:47
Class containing Lagrangian geometry and topology.
Mesh that relates to a sub-section of a Lagrangian mesh. This is used to construct fields that relate...
SubList< Type > sub(const List< Type > &list) const
Return a sub-list corresponding to this sub-mesh.
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.
LagrangianSubMesh Mesh
Mesh type.
label end() const
Return end.
bool empty() const
Return whether or not the mesh is empty.
void operator+=(const LagrangianSubMesh &)
Add a sub-mesh to this one. Must relate to adjacent elements.
ClassName("LagrangianSubMesh")
Runtime type information.
label index() const
Return the index.
static tmp< FieldType > Uf(const LagrangianSubScalarSubField &fraction)
Return the face velocities at the Lagrangian locations.
const LagrangianMesh & mesh() const
Return the mesh.
label globalSize() const
Return size.
label start() const
Return start.
Pre-declare related SubField type.
Definition: SubField.H:63
A class for managing temporary objects.
Definition: tmp.H:55
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
LagrangianGroup
Lagrangian group enumeration.