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