LagrangianMeshI.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 \*---------------------------------------------------------------------------*/
25 
26 #include "LagrangianMesh.H"
27 #include "DimensionedField.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 template<class Field>
32 inline void Foam::LagrangianMesh::checkFieldSize(const Field& field) const
33 {
34  if (field.size() != size())
35  {
37  << "size of " << field.objectPath() << " = "
38  << field.size() << " is not the same as the size of the "
39  << "mesh = " << size() << abort(FatalError);
40  }
41 }
42 
43 
44 template<class Type>
45 inline void Foam::LagrangianMesh::checkPtr
46 (
47  const autoPtr<Type>& ptr,
48  const word& name
49 ) const
50 {
51  if (!ptr.valid())
52  {
54  << name << " requested for a Lagrangian mesh in which "
55  << "they have not been constructed"
56  << exit(FatalError);
57  }
58 }
59 
60 
61 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
62 
64 {
65  return mesh_;
66 }
67 
68 
71 {
72  return coordinates_;
73 }
74 
75 
77 {
78  return celli_;
79 }
80 
81 
83 {
84  return facei_;
85 }
86 
87 
89 {
90  return faceTrii_;
91 }
92 
93 
94 inline const Foam::LagrangianBoundaryMesh&
96 {
97  return boundary_;
98 }
99 
100 
102 {
103  return statesPtr_.valid();
104 }
105 
106 
109 {
110  checkPtr(statesPtr_, "States");
111  return statesPtr_();
112 }
113 
114 
116 {
117  checkPtr(statesPtr_, "States");
118  return statesPtr_();
119 }
120 
121 
123 {
124  return statesPtr_.valid() ? statesPtr_()[i] : LagrangianState::none;
125 }
126 
127 
129 (
130  const LagrangianSubMesh& subMesh,
131  const label subi
132 ) const
133 {
134  return state(subMesh.start() + subi);
135 }
136 
137 
139 {
140  if (statesPtr_.valid())
141  {
142  subMeshIndex_ ++;
143  }
144  else
145  {
146  subMeshIndex_ --;
147  }
148 
149  if (subMeshIndex_ == labelMax || subMeshIndex_ == labelMin)
150  {
152  << "Sub mesh index counter has overflowed."
153  << exit(FatalError);
154  }
155 
156  return subMeshIndex_;
157 }
158 
159 
161 {
162  return coordinates().size();
163 }
164 
165 
167 {
168  return returnReduce(size(), sumOp<label>());
169 }
170 
171 
173 {
174  return mesh.coordinates().size();
175 }
176 
177 
179 {
180  return
181  static_cast<label>(LagrangianState::onPatchZero)
182  + mesh().boundaryMesh().size() + 1;
183 }
184 
185 
187 {
188  return
189  static_cast<label>(LagrangianGroup::onPatchZero)
190  + mesh().boundaryMesh().size() + 1;
191 }
192 
193 
195 (
196  const LagrangianState state
197 ) const
198 {
199  return
200  min
201  (
202  static_cast<label>(stateToGroup(state)),
203  nGroups() - 1
204  );
205 }
206 
207 
209 (
210  const LagrangianGroup group
211 ) const
212 {
213  checkPtr(offsetsPtr_, "Offsets");
214  return LagrangianSubMesh(*this, offsetsPtr_(), group);
215 }
216 
217 
219 {
220  return LagrangianSubMesh(*this, LagrangianGroup::none, 0, size());
221 }
222 
223 
225 {
226  return subAll_;
227 }
228 
229 
231 {
232  checkPtr(offsetsPtr_, "Offsets");
233  return
235  (
236  *this,
238  size()
239  - offsetsPtr_()[static_cast<label>(LagrangianGroup::inInternalMesh)],
240  offsetsPtr_()[static_cast<label>(LagrangianGroup::inInternalMesh)]
241  );
242 }
243 
244 
245 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
246 
247 inline const Foam::vector&
249 {
250  return linear_[i];
251 }
252 
253 
254 inline const Foam::vector& Foam::LagrangianMesh::linearDisplacement::operator()
255 (
256  const label i,
257  const scalar f
258 ) const
259 {
260  return linear_[i];
261 }
262 
263 
266 {
267  return Pair<vector>(linear_[i], quadratic_[i]);
268 }
269 
270 
271 inline Foam::vector Foam::LagrangianMesh::parabolicDisplacement::operator()
272 (
273  const label i,
274  const scalar f
275 ) const
276 {
277  return linear_[i] + (1 - f)*quadratic_[i];
278 }
279 
280 
281 // ************************************************************************* //
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
Boundary part of a Lagrangian mesh. Just a list of Lagrangian patches with some added convenience fun...
const vector & operator()(const label i) const
The displacement for a given index.
Pair< vector > operator()(const label i) const
The displacements for a given index.
Class containing Lagrangian geometry and topology.
LagrangianState state(const label i) const
Return the state for an element of the mesh, or a none state.
const labelIODynamicField & faceTrii() const
Access the face-tet indices.
const LagrangianSubMesh & subAll() const
Return a sub-mesh for all elements.
const polyMesh & mesh() const
Access the mesh.
label size() const
Get the number of elements.
bool changing() const
Return whether or not the mesh is changing.
const labelIODynamicField & celli() const
Access the cell indices.
const LagrangianBoundaryMesh & boundary() const
Return reference to boundary mesh.
label nGroups() const
Return the number of groups.
LagrangianSubMesh subNone() const
Return a sub-mesh for no elements.
LagrangianSubMesh subIncomplete() const
Return a sub-mesh for all incomplete elements.
label nStates() const
Return the number of states.
label subMeshIndex() const
Get a sub-mesh index.
label stateToGroupi(const LagrangianState state) const
Convert a state to a group index.
label globalSize() const
Get the global number of elements.
LagrangianSubMesh sub(const LagrangianGroup group) const
Return a sub-mesh for the given group.
const List< LagrangianState > & states() const
Return the states.
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...
label start() const
Return start.
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 ordered pair of two objects of type <Type> with first() and second() elements.
Definition: Pair.H:66
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
static direction size()
Return the number of elements in the VectorSpace = Ncmpts.
Definition: VectorSpaceI.H:83
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
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
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
errorManip< error > abort(error &err)
Definition: errorManip.H:131
LagrangianState
Lagrangian state enumeration.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
LagrangianGroup stateToGroup(const LagrangianState &state)
Convert from a state enumeration to the corresponding group enumerations.
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
LagrangianGroup
Lagrangian group enumeration.
error FatalError
static const label labelMax
Definition: label.H:62
static const label labelMin
Definition: label.H:61
labelList f(nPoints)