trackedParticle.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2017 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::trackedParticle
26 
27 Description
28  Particle class that marks cells it passes through. Used to mark cells
29  visited by feature edges.
30 
31 SourceFiles
32  trackedParticle.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef trackedParticle_H
37 #define trackedParticle_H
38 
39 #include "particle.H"
40 #include "autoPtr.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 class trackedParticleCloud;
48 
49 
50 // Forward declaration of friend functions and operators
51 
52 class trackedParticle;
53 
54 Ostream& operator<<(Ostream&, const trackedParticle&);
55 
56 
57 /*---------------------------------------------------------------------------*\
58  Class trackedParticle Declaration
59 \*---------------------------------------------------------------------------*/
60 
61 class trackedParticle
62 :
63  public particle
64 {
65  // Private data
66 
67  //- Start point to track from
68  point start_;
69 
70  //- End point to track to
71  point end_;
72 
73  //- Level of this particle
74  label level_;
75 
76  //- Passive label (used to store feature edge mesh)
77  label i_;
78 
79  //- Passive label (used to store feature edge point)
80  label j_;
81 
82  //- Passive label (used to store feature edge label)
83  label k_;
84 
85 
86 public:
87 
88  friend class Cloud<trackedParticle>;
89 
90  //- Class used to pass tracking data to the trackToFace function
91  class trackingData
92  :
93  public particle::TrackingData<Cloud<trackedParticle>>
94  {
95  public:
96 
98 
100 
101 
102  // Constructors
103 
105  (
107  labelList& maxLevel,
108  List<PackedBoolList>& featureEdgeVisited
109  )
110  :
112  maxLevel_(maxLevel),
113  featureEdgeVisited_(featureEdgeVisited)
114  {}
115 
116  };
117 
118 
119 
120  // Constructors
121 
122  //- Construct from components
124  (
125  const polyMesh& mesh,
126  const barycentric& coordinates,
127  const label celli,
128  const label tetFacei,
129  const label tetPtI,
130  const point& end,
131  const label level,
132  const label i,
133  const label j,
134  const label k
135  );
136 
137  //- Construct from a position and a cell, searching for the rest of the
138  // required topology
140  (
141  const polyMesh& mesh,
142  const vector& position,
143  const label celli,
144  const point& end,
145  const label level,
146  const label i,
147  const label j,
148  const label k
149  );
150 
151  //- Construct from Istream
153  (
154  const polyMesh& mesh,
155  Istream& is,
156  bool readFields = true
157  );
158 
159  //- Construct and return a clone
160  autoPtr<particle> clone() const
161  {
162  return autoPtr<particle>(new trackedParticle(*this));
163  }
164 
165  //- Factory class to read-construct particles used for
166  // parallel transfer
167  class iNew
168  {
169  const polyMesh& mesh_;
170 
171  public:
173  iNew(const polyMesh& mesh)
174  :
175  mesh_(mesh)
176  {}
178  autoPtr<trackedParticle> operator()(Istream& is) const
179  {
181  (
182  new trackedParticle(mesh_, is, true)
183  );
184  }
185  };
186 
187 
188  // Member Functions
189 
190  //- Point to track from
191  point& start()
192  {
193  return start_;
194  }
195 
196  //- Point to track to
197  point& end()
198  {
199  return end_;
200  }
201 
202  //- Transported label
203  label i() const
204  {
205  return i_;
206  }
207 
208  //- Transported label
209  label& i()
210  {
211  return i_;
212  }
213 
214  //- Transported label
215  label j() const
216  {
217  return j_;
218  }
219 
220  //- Transported label
221  label& j()
222  {
223  return j_;
224  }
225 
226  //- Transported label
227  label k() const
228  {
229  return k_;
230  }
231 
232  //- Transported label
233  label& k()
234  {
235  return k_;
236  }
237 
238 
239 
240  // Tracking
241 
242  //- Track all particles to their end point
243  bool move(trackingData&, const scalar);
244 
245 
246  //- Overridable function to handle the particle hitting a patch
247  // Executed before other patch-hitting functions
248  bool hitPatch
249  (
250  const polyPatch&,
251  trackingData& td,
252  const label patchi,
253  const scalar trackFraction,
254  const tetIndices& tetIs
255  );
256 
257  //- Overridable function to handle the particle hitting a wedge
258  void hitWedgePatch
259  (
260  const wedgePolyPatch&,
261  trackingData& td
262  );
263 
264  //- Overridable function to handle the particle hitting a
265  // symmetry plane
267  (
268  const symmetryPlanePolyPatch&,
269  trackingData& td
270  );
271 
272  //- Overridable function to handle the particle hitting a
273  // symmetry patch
274  void hitSymmetryPatch
275  (
276  const symmetryPolyPatch&,
277  trackingData& td
278  );
279 
280  //- Overridable function to handle the particle hitting a cyclic
281  void hitCyclicPatch
282  (
283  const cyclicPolyPatch&,
284  trackingData& td
285  );
286 
287  //- Overridable function to handle the particle hitting a cyclicAMI
288  void hitCyclicAMIPatch
289  (
290  const cyclicAMIPolyPatch&,
291  trackingData& td,
292  const vector&
293  );
294 
295  //- Overridable function to handle the particle hitting a
296  //- processorPatch
297  void hitProcessorPatch
298  (
299  const processorPolyPatch&,
300  trackingData& td
301  );
302 
303  //- Overridable function to handle the particle hitting a wallPatch
304  void hitWallPatch
305  (
306  const wallPolyPatch&,
307  trackingData& td,
308  const tetIndices&
309  );
310 
311  //- Overridable function to handle the particle hitting a polyPatch
312  void hitPatch
313  (
314  const polyPatch&,
315  trackingData& td
316  );
317 
318  //- Convert processor patch addressing to the global equivalents
319  // and set the celli to the face-neighbour
321 
322 
323  // Ostream Operator
324 
325  friend Ostream& operator<<(Ostream&, const trackedParticle&);
326 };
327 
328 
329 template<>
330 inline bool contiguous<trackedParticle>()
331 {
332  return true;
333 }
334 
335 //template<>
336 //void Cloud<trackedParticle>::readFields();
337 //
338 //template<>
339 //void Cloud<trackedParticle>::writeFields() const;
340 
341 
342 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
343 
344 } // End namespace Foam
345 
346 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
347 
348 #endif
349 
350 // ************************************************************************* //
Symmetry patch for non-planar or multi-plane patches.
const polyMesh & mesh() const
Return the mesh database.
Definition: particleI.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
void hitCyclicAMIPatch(const cyclicAMIPolyPatch &, trackingData &td, const vector &)
Overridable function to handle the particle hitting a cyclicAMI.
List< PackedBoolList > & featureEdgeVisited_
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
autoPtr< particle > clone() const
Construct and return a clone.
void hitCyclicPatch(const cyclicPolyPatch &, trackingData &td)
Overridable function to handle the particle hitting a cyclic.
static void readFields(CloudType &c)
Read the fields associated with the owner cloud.
Base particle class.
Definition: particle.H:81
Neighbour processor patch.
trackedParticle(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPtI, const point &end, const label level, const label i, const label j, const label k)
Construct from components.
void hitProcessorPatch(const processorPolyPatch &, trackingData &td)
Particle class that marks cells it passes through. Used to mark cells visited by feature edges...
void hitWallPatch(const wallPolyPatch &, trackingData &td, const tetIndices &)
Overridable function to handle the particle hitting a wallPatch.
Cyclic plane patch.
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
Foam::wallPolyPatch.
Definition: wallPolyPatch.H:48
point & end()
Point to track to.
Wedge front and back plane patch.
point & start()
Point to track from.
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:81
label k() const
Transported label.
Cyclic patch for Arbitrary Mesh Interface (AMI)
Factory class to read-construct particles used for.
Base cloud calls templated on particle type.
Definition: Cloud.H:52
bool contiguous< trackedParticle >()
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
void hitSymmetryPlanePatch(const symmetryPlanePolyPatch &, trackingData &td)
Overridable function to handle the particle hitting a.
void correctAfterParallelTransfer(const label, trackingData &)
Convert processor patch addressing to the global equivalents.
Class used to pass tracking data to the trackToFace function.
label patchi
void hitSymmetryPatch(const symmetryPolyPatch &, trackingData &td)
Overridable function to handle the particle hitting a.
label i() const
Transported label.
Ostream & operator<<(Ostream &, const ensightPart &)
trackingData(Cloud< trackedParticle > &cloud, labelList &maxLevel, List< PackedBoolList > &featureEdgeVisited)
const barycentric & coordinates() const
Return current particle coordinates.
Definition: particleI.H:51
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
label j() const
Transported label.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
void hitWedgePatch(const wedgePolyPatch &, trackingData &td)
Overridable function to handle the particle hitting a wedge.
friend Ostream & operator<<(Ostream &, const trackedParticle &)
vector position() const
Return current particle position.
Definition: particleI.H:204
bool move(trackingData &, const scalar)
Track all particles to their end point.
bool hitPatch(const polyPatch &, trackingData &td, const label patchi, const scalar trackFraction, const tetIndices &tetIs)
Overridable function to handle the particle hitting a patch.
Namespace for OpenFOAM.