tracking.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) 2024 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 Namespace
25  Foam::tracking
26 
27 Description
28  Functions for tracking locations through a mesh
29 
30 SourceFiles
31  tracking.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef tracking_H
36 #define tracking_H
37 
38 #include "polyMesh.H"
39 #include "barycentric.H"
40 #include "barycentricTensor.H"
41 #include "vector.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 
48 class wedgePolyPatch;
49 class cyclicPolyPatch;
50 class processorPolyPatch;
51 
52 namespace tracking
53 {
54 
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 
57 // Geometry
58 
59  //- Return the position given the coordinates and tet topology
60  inline point position
61  (
62  const polyMesh& mesh,
63  const barycentric& coordinates,
64  const label celli,
65  const label facei,
66  const label faceTrii,
67  const scalar stepFraction
68  );
69 
70  //- Return the coordinates given the position and tet topology
72  (
73  const polyMesh& mesh,
74  const point& position,
75  const label celli,
76  const label facei,
77  const label faceTrii,
78  const scalar stepFraction
79  );
80 
81  //- Return the normal of the corresponding point on the associated face and
82  // the displacement of that point over the time-step
83  Pair<vector> faceNormalAndDisplacement
84  (
85  const polyMesh& mesh,
86  const barycentric& coordinates,
87  const label celli,
88  const label facei,
89  const label faceTrii,
90  const scalar stepFraction
91  );
92 
93 
94 // Tracking
95 
96  //- Track along the displacement for a given fraction of the overall
97  // time-step. End when the track is complete or when a face is hit. Return
98  // whether or not a face was hit (true) or the track completed (false) and
99  // the proportion of the displacement still to be completed. Displacement
100  // can be a single vector, d, in which case the path is linear and equal
101  // to x0 + lambda*d, where x0 is the starting position and lambda is the
102  // local coordinate along the path. Or, displacement can be a pair of
103  // vectors, d1 and d2, in which case the path is parabolic and equal to
104  // x0 + lambda*d1 + lambda*lambda*d2.
105  template<class Displacement>
107  (
108  const polyMesh& mesh,
109  const Displacement& displacement,
110  const scalar fraction,
112  label& celli,
113  label& facei,
114  label& faceTrii,
115  scalar& stepFraction,
116  scalar& stepFractionBehind,
117  label& nTracksBehind,
118  const string& debugPrefix = NullObjectRef<string>()
119  );
120 
121  //- As toFace, except that if the track ends on an internal face then this
122  // face will be crossed
123  template<class Displacement>
125  (
126  const polyMesh& mesh,
127  const Displacement& displacement,
128  const scalar fraction,
130  label& celli,
131  label& facei,
132  label& faceTrii,
133  scalar& stepFraction,
134  scalar& stepFractionBehind,
135  label& nTracksBehind,
136  const string& debugPrefix = NullObjectRef<string>()
137  );
138 
139  //- As toFace, except that the track continues across multiple cells until
140  // it ends or until a boundary face is hit
141  template<class Displacement>
143  (
144  const polyMesh& mesh,
145  const Displacement& displacement,
146  const scalar fraction,
148  label& celli,
149  label& facei,
150  label& faceTrii,
151  scalar& stepFraction,
152  scalar& stepFractionBehind,
153  label& nTracksBehind,
154  const string& debugPrefix = NullObjectRef<string>()
155  );
156 
157 
158 // Initialisation
159 
160  //- Initialise the location at the given position. Returns whether or not a
161  // boundary was hit. The cell index should be initialised to -1, or the
162  // index of the cell.
163  bool locate
164  (
165  const polyMesh& mesh,
166  const point& position,
168  label& celli,
169  label& facei,
170  label& faceTrii,
171  const scalar stepFraction,
172  const string& debugPrefix = NullObjectRef<string>()
173  );
174 
175 
176 // Topology Changes
177 
178  //- Cross an internal face
179  void crossInternalFace
180  (
181  const polyMesh& mesh,
183  label& celli,
184  label& facei,
185  label& faceTrii
186  );
187 
188  //- Cross a wedge patch
189  void crossWedge
190  (
191  const wedgePolyPatch& inPatch,
193  label& celli,
194  label& facei,
195  label& faceTrii,
196  const scalar stepFraction
197  );
198 
199  //- Cross a cyclic patch
200  void crossCyclic
201  (
202  const cyclicPolyPatch& inPatch,
204  label& celli,
205  label& facei,
206  label& faceTrii
207  );
208 
209  //- Initialise crossing of a processor patch. Breaks the topology in order
210  // to store the destination patch face in advance of communication.
211  void inProcessor
212  (
213  const processorPolyPatch& inPatch,
214  label& celli,
215  label& facei
216  );
217 
218  //- Complete crossing of a processor patch. Restore the topology.
219  void outProcessor
220  (
221  const processorPolyPatch& outPatch,
223  label& celli,
224  label& facei,
225  label& faceTrii
226  );
227 
228 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
229 
230 } // End namespace tracking
231 } // End namespace Foam
232 
233 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
234 
235 #include "trackingI.H"
236 
237 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
238 
239 #endif
240 
241 // ************************************************************************* //
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:66
Cyclic plane patch.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
Neighbour processor patch.
Wedge front and back plane patch.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
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
Tuple2< bool, scalar > toCell(const polyMesh &mesh, const Displacement &displacement, const scalar fraction, barycentric &coordinates, label &celli, label &facei, label &faceTrii, scalar &stepFraction, scalar &stepFractionBehind, label &nTracksBehind, const string &debugPrefix=NullObjectRef< string >())
As toFace, except that if the track ends on an internal face then this.
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
Tuple2< bool, scalar > toBoundary(const polyMesh &mesh, const Displacement &displacement, const scalar fraction, barycentric &coordinates, label &celli, label &facei, label &faceTrii, scalar &stepFraction, scalar &stepFractionBehind, label &nTracksBehind, const string &debugPrefix=NullObjectRef< string >())
As toFace, except that the track continues across multiple cells until.
void crossInternalFace(const polyMesh &mesh, barycentric &coordinates, label &celli, label &facei, label &faceTrii)
Cross an internal face.
Definition: tracking.C:1681
void crossCyclic(const cyclicPolyPatch &inPatch, barycentric &coordinates, label &celli, label &facei, label &faceTrii)
Cross a cyclic patch.
Definition: tracking.C:1736
void outProcessor(const processorPolyPatch &outPatch, barycentric &coordinates, label &celli, label &facei, label &faceTrii)
Complete crossing of a processor patch. Restore the topology.
Definition: tracking.C:1777
void inProcessor(const processorPolyPatch &inPatch, label &celli, label &facei)
Initialise crossing of a processor patch. Breaks the topology in order.
Definition: tracking.C:1762
void crossWedge(const wedgePolyPatch &inPatch, barycentric &coordinates, label &celli, label &facei, label &faceTrii, const scalar stepFraction)
Cross a wedge patch.
Definition: tracking.C:1700
point position(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label facei, const label faceTrii, const scalar stepFraction)
Return the position given the coordinates and tet topology.
Definition: trackingI.H:224
Tuple2< bool, scalar > toFace(const polyMesh &mesh, const Displacement &displacement, const scalar fraction, barycentric &coordinates, label &celli, label &facei, label &faceTrii, scalar &stepFraction, scalar &stepFractionBehind, label &nTracksBehind, const string &debugPrefix=NullObjectRef< string >())
Track along the displacement for a given fraction of the overall.
bool locate(const polyMesh &mesh, const point &position, barycentric &coordinates, label &celli, label &facei, label &faceTrii, const scalar stepFraction, const string &debugPrefix=NullObjectRef< string >())
Initialise the location at the given position. Returns whether or not a.
Definition: tracking.C:1593
Namespace for OpenFOAM.
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition: barycentric.H:45
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
vector point
Point is a vector.
Definition: point.H:41