wallBoundedParticle.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::wallBoundedParticle
26 
27 Description
28  Particle class that tracks on triangles of boundary faces. Use
29  trackToEdge similar to trackToFace on particle.
30 
31 SourceFiles
32  wallBoundedParticle.C
33  wallBoundedParticleTemplates.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef wallBoundedParticle_H
38 #define wallBoundedParticle_H
39 
40 #include "particle.H"
41 #include "autoPtr.H"
42 #include "InfoProxy.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 // Forward declaration of friend functions and operators
50 
51 class wallBoundedParticle;
52 
53 Ostream& operator<<(Ostream&, const wallBoundedParticle&);
54 Ostream& operator<<(Ostream&, const InfoProxy<wallBoundedParticle>&);
55 
56 
57 /*---------------------------------------------------------------------------*\
58  Class wallBoundedParticle Declaration
59 \*---------------------------------------------------------------------------*/
60 
62 :
63  public particle
64 {
65  // Private data
66 
67  //- Size in bytes of the fields
68  static const std::size_t sizeofFields_;
69 
70 
71 public:
72 
73  //- Class used to pass tracking data to the trackToFace function
74  template<class CloudType>
75  class TrackingData
76  :
77  public particle::TrackingData<CloudType>
78  {
79 
80  public:
81 
83 
84  // Constructors
85 
86  inline TrackingData
87  (
89  const PackedBoolList& isWallPatch
90  )
91  :
93  (
94  cloud
95  ),
96  isWallPatch_(isWallPatch)
97  {}
98  };
99 
100 
101 protected:
102 
103  // Protected data
104 
105  //- Particle is on mesh edge:
106  // const face& f = mesh.faces()[tetFace()]
107  // const edge e(f[meshEdgeStart_], f.nextLabel(meshEdgeStart_));
108  // Note that this real edge
109  // is also one of the edges of the face-triangle (from
110  // tetFace()+tetPt()).
112 
113  //- Particle is on diagonal edge:
114  // const face& f = mesh.faces()[tetFace()]
115  // label faceBasePtI = mesh.tetBasePtIs()[facei];
116  // label diagPtI = (faceBasePtI+diagEdge_)%f.size();
117  // const edge e(f[faceBasePtI], f[diagPtI]);
119 
120 
121  // Protected Member Functions
122 
123  //- Construct current edge
124  edge currentEdge() const;
125 
126  //- Check if inside current tet
127  //void checkInside() const;
128 
129  //- Check if on current edge
130  //void checkOnEdge() const;
131 
132  //- Check if point on triangle
133  //void checkOnTriangle(const point&) const;
134 
135  //- Cross mesh edge into different face on same cell
136  void crossEdgeConnectedFace(const edge& meshEdge);
137 
138  //- Cross diagonal edge into different triangle on same face,cell
139  void crossDiagonalEdge();
140 
141  //- Track through single triangle
142  scalar trackFaceTri(const vector& endPosition, label& minEdgeI);
143 
144  //- Is current triangle in the track direction
145  bool isTriAlongTrack(const point& endPosition) const;
146 
147 
148  // Patch interactions
149 
150  //- Do all patch interaction
151  template<class TrackData>
152  void patchInteraction(TrackData& td, const scalar trackFraction);
153 
154  //- Overridable function to handle the particle hitting a patch
155  // Executed before other patch-hitting functions
156  template<class TrackData>
157  bool hitPatch
158  (
159  const polyPatch&,
160  TrackData& td,
161  const label patchi,
162  const scalar trackFraction,
163  const tetIndices& tetIs
164  );
165 
166  //- Overridable function to handle the particle hitting a wedge
167  template<class TrackData>
168  void hitWedgePatch
169  (
170  const wedgePolyPatch&,
171  TrackData& td
172  );
173 
174  //- Overridable function to handle the particle hitting a
175  // symmetry plane
176  template<class TrackData>
178  (
179  const symmetryPlanePolyPatch&,
180  TrackData& td
181  );
182 
183  //- Overridable function to handle the particle hitting a
184  // symmetry patch
185  template<class TrackData>
186  void hitSymmetryPatch
187  (
188  const symmetryPolyPatch&,
189  TrackData& td
190  );
191 
192  //- Overridable function to handle the particle hitting a cyclic
193  template<class TrackData>
194  void hitCyclicPatch
195  (
196  const cyclicPolyPatch&,
197  TrackData& td
198  );
199 
200  //- Overridable function to handle the particle hitting a
201  //- processorPatch
202  template<class TrackData>
203  void hitProcessorPatch
204  (
205  const processorPolyPatch&,
206  TrackData& td
207  );
208 
209  //- Overridable function to handle the particle hitting a wallPatch
210  template<class TrackData>
211  void hitWallPatch
212  (
213  const wallPolyPatch&,
214  TrackData& td,
215  const tetIndices&
216  );
217 
218  //- Overridable function to handle the particle hitting a polyPatch
219  template<class TrackData>
220  void hitPatch
221  (
222  const polyPatch&,
223  TrackData& td
224  );
225 
226 public:
227 
228  // Constructors
229 
230  //- Construct from components
232  (
233  const polyMesh& c,
234  const vector& position,
235  const label celli,
236  const label tetFacei,
237  const label tetPtI,
238  const label meshEdgeStart,
239  const label diagEdge
240  );
241 
242  //- Construct from Istream
244  (
245  const polyMesh& c,
246  Istream& is,
247  bool readFields = true
248  );
249 
250  //- Construct copy
252 
253  //- Construct and return a clone
254  autoPtr<particle> clone() const
255  {
256  return autoPtr<particle>(new wallBoundedParticle(*this));
257  }
258 
259  //- Factory class to read-construct particles used for
260  // parallel transfer
261  class iNew
262  {
263  const polyMesh& mesh_;
264 
265  public:
267  iNew(const polyMesh& mesh)
268  :
269  mesh_(mesh)
270  {}
271 
273  (
274  Istream& is
275  ) const
276  {
278  (
279  new wallBoundedParticle(mesh_, is, true)
280  );
281  }
282  };
283 
284 
285  // Member Functions
286 
287  // Access
288 
289  //- -1 or label of mesh edge
290  inline label meshEdgeStart() const
291  {
292  return meshEdgeStart_;
293  }
294 
295  //- -1 or diagonal edge
296  inline label diagEdge() const
297  {
298  return diagEdge_;
299  }
300 
301 
302  // Track
303 
304  //- Equivalent of trackToFace
305  template<class TrackData>
306  scalar trackToEdge
307  (
308  TrackData& td,
309  const vector& endPosition
310  );
311 
312 
313  // Info
314 
315  //- Return info proxy.
316  // Used to print particle information to a stream
317  inline InfoProxy<wallBoundedParticle> info() const
318  {
319  return *this;
320  }
321 
322 
323  // I-O
324 
325  //- Read
326  template<class CloudType>
327  static void readFields(CloudType&);
328 
329  //- Write
330  template<class CloudType>
331  static void writeFields(const CloudType&);
332 
333 
334  // Ostream Operator
335 
336  friend Ostream& operator<<
337  (
338  Ostream&,
339  const wallBoundedParticle&
340  );
341 
342  friend Ostream& operator<<
343  (
344  Ostream&,
346  );
347 };
348 
349 
350 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
351 
352 } // End namespace Foam
353 
354 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
355 
356 #ifdef NoRepository
358 #endif
359 
360 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
361 
362 #endif
363 
364 // ************************************************************************* //
Symmetry patch for non-planar or multi-plane patches.
label meshEdgeStart_
Particle is on mesh edge:
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 hitSymmetryPatch(const symmetryPolyPatch &, TrackData &td)
Overridable function to handle the particle hitting a.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
Factory class to read-construct particles used for.
bool hitPatch(const polyPatch &, TrackData &td, const label patchi, const scalar trackFraction, const tetIndices &tetIs)
Overridable function to handle the particle hitting a patch.
edge currentEdge() const
Construct current edge.
void hitWedgePatch(const wedgePolyPatch &, TrackData &td)
Overridable function to handle the particle hitting a wedge.
Base particle class.
Definition: particle.H:78
Neighbour processor patch.
void hitProcessorPatch(const processorPolyPatch &, TrackData &td)
const vector & position() const
Return current particle position.
Definition: particleI.H:586
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Cyclic plane patch.
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
Foam::wallPolyPatch.
Definition: wallPolyPatch.H:48
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
label diagEdge_
Particle is on diagonal edge:
Particle class that tracks on triangles of boundary faces. Use trackToEdge similar to trackToFace on ...
bool isTriAlongTrack(const point &endPosition) const
Is current triangle in the track direction.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static void readFields(CloudType &)
Read.
void crossDiagonalEdge()
Cross diagonal edge into different triangle on same face,cell.
void hitWallPatch(const wallPolyPatch &, TrackData &td, const tetIndices &)
Overridable function to handle the particle hitting a wallPatch.
label meshEdgeStart() const
-1 or label of mesh edge
static void writeFields(const CloudType &)
Write.
void hitCyclicPatch(const cyclicPolyPatch &, TrackData &td)
Overridable function to handle the particle hitting a cyclic.
A bit-packed bool list.
void crossEdgeConnectedFace(const edge &meshEdge)
Check if inside current tet.
label patchi
A helper class for outputting values to Ostream.
Definition: InfoProxy.H:45
wallBoundedParticle(const polyMesh &c, const vector &position, const label celli, const label tetFacei, const label tetPtI, const label meshEdgeStart, const label diagEdge)
Construct from components.
scalar trackFaceTri(const vector &endPosition, label &minEdgeI)
Track through single triangle.
void patchInteraction(TrackData &td, const scalar trackFraction)
Do all patch interaction.
scalar trackToEdge(TrackData &td, const vector &endPosition)
Equivalent of trackToFace.
const dimensionedScalar c
Speed of light in a vacuum.
Ostream & operator<<(Ostream &, const ensightPart &)
autoPtr< particle > clone() const
Construct and return a clone.
const polyMesh & mesh() const
Return the mesh database.
Definition: particleI.H:580
TrackingData(CloudType &cloud, const PackedBoolList &isWallPatch)
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
Class used to pass tracking data to the trackToFace function.
label diagEdge() const
-1 or diagonal edge
volScalarField & p
void hitSymmetryPlanePatch(const symmetryPlanePolyPatch &, TrackData &td)
Overridable function to handle the particle hitting a.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
InfoProxy< wallBoundedParticle > info() const
Return info proxy.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
Namespace for OpenFOAM.