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