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-2016 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  //- End point to track to
68  point end_;
69 
70  //- Level of this particle
71  label level_;
72 
73  //- Passive label (used to store feature edge mesh)
74  label i_;
75 
76  //- Passive label (used to store feature edge point)
77  label j_;
78 
79  //- Passive label (used to store feature edge label)
80  label k_;
81 
82 
83 public:
84 
85  friend class Cloud<trackedParticle>;
86 
87  //- Class used to pass tracking data to the trackToFace function
88  class trackingData
89  :
90  public particle::TrackingData<Cloud<trackedParticle>>
91  {
92  public:
93 
95 
97 
98 
99  // Constructors
100 
102  (
104  labelList& maxLevel,
105  List<PackedBoolList>& featureEdgeVisited
106  )
107  :
109  maxLevel_(maxLevel),
110  featureEdgeVisited_(featureEdgeVisited)
111  {}
112 
113  };
114 
115 
116 
117  // Constructors
118 
119  //- Construct from components
121  (
122  const polyMesh& mesh,
123  const vector& position,
124  const label celli,
125  const label tetFacei,
126  const label tetPtI,
127  const point& end,
128  const label level,
129  const label i,
130  const label j,
131  const label k
132  );
133 
134  //- Construct from Istream
136  (
137  const polyMesh& mesh,
138  Istream& is,
139  bool readFields = true
140  );
141 
142  //- Construct and return a clone
143  autoPtr<particle> clone() const
144  {
145  return autoPtr<particle>(new trackedParticle(*this));
146  }
147 
148  //- Factory class to read-construct particles used for
149  // parallel transfer
150  class iNew
151  {
152  const polyMesh& mesh_;
153 
154  public:
156  iNew(const polyMesh& mesh)
157  :
158  mesh_(mesh)
159  {}
161  autoPtr<trackedParticle> operator()(Istream& is) const
162  {
164  (
165  new trackedParticle(mesh_, is, true)
166  );
167  }
168  };
169 
170 
171  // Member Functions
172 
173  //- Point to track to
174  point& end()
175  {
176  return end_;
177  }
178 
179  //- Transported label
180  label i() const
181  {
182  return i_;
183  }
184 
185  //- Transported label
186  label& i()
187  {
188  return i_;
189  }
190 
191  //- Transported label
192  label j() const
193  {
194  return j_;
195  }
196 
197  //- Transported label
198  label& j()
199  {
200  return j_;
201  }
202 
203  //- Transported label
204  label k() const
205  {
206  return k_;
207  }
208 
209  //- Transported label
210  label& k()
211  {
212  return k_;
213  }
214 
215 
216 
217  // Tracking
218 
219  //- Track all particles to their end point
220  bool move(trackingData&, const scalar);
221 
222 
223  //- Overridable function to handle the particle hitting a patch
224  // Executed before other patch-hitting functions
225  bool hitPatch
226  (
227  const polyPatch&,
228  trackingData& td,
229  const label patchi,
230  const scalar trackFraction,
231  const tetIndices& tetIs
232  );
233 
234  //- Overridable function to handle the particle hitting a wedge
235  void hitWedgePatch
236  (
237  const wedgePolyPatch&,
238  trackingData& td
239  );
240 
241  //- Overridable function to handle the particle hitting a
242  // symmetry plane
244  (
245  const symmetryPlanePolyPatch&,
246  trackingData& td
247  );
248 
249  //- Overridable function to handle the particle hitting a
250  // symmetry patch
251  void hitSymmetryPatch
252  (
253  const symmetryPolyPatch&,
254  trackingData& td
255  );
256 
257  //- Overridable function to handle the particle hitting a cyclic
258  void hitCyclicPatch
259  (
260  const cyclicPolyPatch&,
261  trackingData& td
262  );
263 
264  //- Overridable function to handle the particle hitting a
265  //- processorPatch
266  void hitProcessorPatch
267  (
268  const processorPolyPatch&,
269  trackingData& td
270  );
271 
272  //- Overridable function to handle the particle hitting a wallPatch
273  void hitWallPatch
274  (
275  const wallPolyPatch&,
276  trackingData& td,
277  const tetIndices&
278  );
279 
280  //- Overridable function to handle the particle hitting a polyPatch
281  void hitPatch
282  (
283  const polyPatch&,
284  trackingData& td
285  );
286 
287  //- Convert processor patch addressing to the global equivalents
288  // and set the celli to the face-neighbour
290 
291 
292  // Ostream Operator
293 
294  friend Ostream& operator<<(Ostream&, const trackedParticle&);
295 };
296 
297 
298 template<>
299 inline bool contiguous<trackedParticle>()
300 {
301  return true;
302 }
303 
304 //template<>
305 //void Cloud<trackedParticle>::readFields();
306 //
307 //template<>
308 //void Cloud<trackedParticle>::writeFields() const;
309 
310 
311 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
312 
313 } // End namespace Foam
314 
315 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
316 
317 #endif
318 
319 // ************************************************************************* //
Symmetry patch for non-planar or multi-plane patches.
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
List< PackedBoolList > & featureEdgeVisited_
trackedParticle(const polyMesh &mesh, const vector &position, 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.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
label i() const
Transported label.
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:78
Neighbour processor patch.
void hitProcessorPatch(const processorPolyPatch &, trackingData &td)
const vector & position() const
Return current particle position.
Definition: particleI.H:586
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.
const polyMesh & mesh_
Reference to the polyMesh database.
Definition: particle.H:137
Wedge front and back plane patch.
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:81
autoPtr< particle > clone() const
Construct and return a clone.
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 j() const
Transported label.
Ostream & operator<<(Ostream &, const ensightPart &)
trackingData(Cloud< trackedParticle > &cloud, labelList &maxLevel, List< PackedBoolList > &featureEdgeVisited)
label k() const
Transported label.
const polyMesh & mesh() const
Return the mesh database.
Definition: particleI.H:580
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
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 &)
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.