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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2017-2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::wallBoundedParticle
29 
30 Description
31  Particle class that tracks on triangles of boundary faces. Use
32  trackToEdge similar to trackToFace on particle.
33 
34 SourceFiles
35  wallBoundedParticle.C
36  wallBoundedParticleTemplates.C
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef Foam_wallBoundedParticle_H
41 #define Foam_wallBoundedParticle_H
42 
43 #include "particle.H"
44 #include "autoPtr.H"
45 #include "InfoProxy.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 // Forward Declarations
53 class wallBoundedParticle;
54 
55 Ostream& operator<<(Ostream&, const wallBoundedParticle&);
56 Ostream& operator<<(Ostream&, const InfoProxy<wallBoundedParticle>&);
57 
58 
59 /*---------------------------------------------------------------------------*\
60  Class wallBoundedParticle Declaration
61 \*---------------------------------------------------------------------------*/
62 
64 :
65  public particle
66 {
67 public:
68 
69  //- Class used to pass tracking data to the trackToFace function
70  class trackingData
71  :
73  {
74  public:
75 
76  const bitSet& isWallPatch_;
77 
78  // Constructors
79 
80  template<class TrackCloudType>
82  (
83  const TrackCloudType& cloud,
84  const bitSet& isWallPatch
85  )
86  :
88  isWallPatch_(isWallPatch)
89  {}
90  };
91 
92 
93 protected:
94 
95  // Protected Data
96 
97  //- Particle position is updated locally as opposed to via track
98  // functions of the base Foam::particle class
100 
101  //- Particle is on mesh edge:
102  // const face& f = mesh.faces()[tetFace()]
103  // const edge e(f[meshEdgeStart_], f.nextLabel(meshEdgeStart_));
104  // Note that this real edge
105  // is also one of the edges of the face-triangle (from
106  // tetFace()+tetPt()).
107  label meshEdgeStart_;
108 
109  //- Particle is on diagonal edge:
110  // const face& f = mesh.faces()[tetFace()]
111  // label faceBasePti = mesh.tetBasePtIs()[facei];
112  // label diagPti = (faceBasePti+diagEdge_)%f.size();
113  // const edge e(f[faceBasePti], f[diagPti]);
114  label diagEdge_;
115 
116 
117  // Protected Member Functions
118 
119  //- Construct current edge
120  edge currentEdge() const;
121 
122  //- Replacement for particle::crossEdgeConnectedFace that avoids bombing
123  // out on invalid tet decomposition (tetBasePtIs = -1)
125  (
126  const label& celli,
127  label& tetFacei,
128  label& tetPti,
129  const edge& e
130  );
131 
132  //- Cross mesh edge into different face on same cell
133  void crossEdgeConnectedFace(const edge& meshEdge);
134 
135  //- Cross diagonal edge into different triangle on same face,cell
136  void crossDiagonalEdge();
137 
138  //- Track through single triangle
139  scalar trackFaceTri(const vector& n, const vector& endPosition, label&);
140 
141  //- Is current triangle in the track direction
142  bool isTriAlongTrack(const vector& n, const point& endPosition) const;
143 
144 
145 public:
146 
147  // Static Data Members
148 
149  //- Size in bytes of the fields
150  static const std::size_t sizeofFields_;
151 
152 
153  // Constructors
154 
155  //- Construct from components
157  (
158  const polyMesh& c,
159  const point& position,
160  const label celli,
161  const label tetFacei,
162  const label tetPti,
163  const label meshEdgeStart,
164  const label diagEdge
165  );
166 
167  //- Construct from Istream
169  (
170  const polyMesh& c,
171  Istream& is,
172  bool readFields = true,
173  bool newFormat = true
174  );
175 
176  //- Construct copy
178 
179  //- Construct and return a clone
180  autoPtr<particle> clone() const
181  {
182  return autoPtr<particle>(new wallBoundedParticle(*this));
183  }
184 
185  //- Factory class to read-construct particles used for
186  // parallel transfer
187  class iNew
188  {
189  const polyMesh& mesh_;
190 
191  public:
192 
193  iNew(const polyMesh& mesh)
194  :
195  mesh_(mesh)
196  {}
197 
198  autoPtr<wallBoundedParticle> operator()
199  (
200  Istream& is
201  ) const
202  {
203  return autoPtr<wallBoundedParticle>
204  (
205  new wallBoundedParticle(mesh_, is, true)
206  );
207  }
208  };
210 
211  // Member Functions
212 
213  // Access
214 
215  //- The mesh edge label or -1
216  label meshEdgeStart() const noexcept { return meshEdgeStart_; }
217 
218  //- The diagonal edge label or -1
219  label diagEdge() const noexcept { return diagEdge_; }
220 
221 
222  // Track
223 
224  //- Equivalent of trackToFace
225  template<class TrackCloudType>
226  scalar trackToEdge
227  (
228  TrackCloudType& cloud,
229  trackingData& td,
230  const vector& endPosition
231  );
232 
233 
234  // Patch interactions
235 
236  //- Do all patch interaction
237  template<class TrackCloudType>
238  void patchInteraction
239  (
240  TrackCloudType& cloud,
241  trackingData& td,
242  const scalar trackFraction
243  );
244 
245  //- Overridable function to handle the particle hitting a
246  //- processorPatch
247  template<class TrackCloudType>
248  void hitProcessorPatch
249  (
250  TrackCloudType& cloud,
251  trackingData& td
252  );
253 
254  //- Overridable function to handle the particle hitting a wallPatch
255  template<class TrackCloudType>
256  void hitWallPatch
257  (
258  TrackCloudType& cloud,
259  trackingData& td
260  );
261 
262 
263  // Info
264 
265  //- Return info proxy,
266  //- used to print particle information to a stream
268  {
269  return *this;
270  }
271 
272 
273  // I-O
274 
275  //- Read
276  template<class CloudType>
277  static void readFields(CloudType&);
278 
279  //- Write
280  template<class CloudType>
281  static void writeFields(const CloudType&);
282 
283 
284  // Ostream Operator
285 
286  friend Ostream& operator<<
287  (
288  Ostream&,
289  const wallBoundedParticle&
290  );
291 
292  friend Ostream& operator<<
293  (
294  Ostream&,
296  );
297 };
298 
299 
300 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
301 
302 } // End namespace Foam
303 
304 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
305 
306 #ifdef NoRepository
308 #endif
309 
310 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
311 
312 #endif
314 // ************************************************************************* //
label meshEdgeStart_
Particle is on mesh edge:
Class used to pass tracking data to the trackToFace function.
trackingData(const TrackCloudType &cloud, const bitSet &isWallPatch)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
label diagEdge() const noexcept
The diagonal edge label or -1.
static const std::size_t sizeofFields_
Size in bytes of the fields.
void patchInteraction(TrackCloudType &cloud, trackingData &td, const scalar trackFraction)
Do all patch interaction.
static void readFields(CloudType &)
Read.
Base particle class.
Definition: particle.H:69
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition: edge.H:59
void hitProcessorPatch(TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a processorPatch.
bool isTriAlongTrack(const vector &n, const point &endPosition) const
Is current triangle in the track direction.
scalar trackToEdge(TrackCloudType &cloud, trackingData &td, const vector &endPosition)
Equivalent of trackToFace.
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:53
point localPosition_
Particle position is updated locally as opposed to via track.
label diagEdge_
Particle is on diagonal edge:
Particle class that tracks on triangles of boundary faces. Use trackToEdge similar to trackToFace on ...
static void writeFields(const CloudType &)
Write.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
const direction noexcept
Definition: Scalar.H:258
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
wallBoundedParticle(const polyMesh &c, const point &position, const label celli, const label tetFacei, const label tetPti, const label meshEdgeStart, const label diagEdge)
Construct from components.
void crossDiagonalEdge()
Cross diagonal edge into different triangle on same face,cell.
void hitWallPatch(TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a wallPatch.
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:77
const polyMesh & mesh() const noexcept
Return the mesh database.
Definition: particleI.H:110
scalar trackFaceTri(const vector &n, const vector &endPosition, label &)
Track through single triangle.
A helper class for outputting values to Ostream.
Definition: ensightCells.H:43
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
const dimensionedScalar c
Speed of light in a vacuum.
label meshEdgeStart() const noexcept
The mesh edge label or -1.
label n
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
volScalarField & p
InfoProxy< wallBoundedParticle > info() const noexcept
Return info proxy, used to print particle information to a stream.
edge currentEdge() const
Construct current edge.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
vector position() const
Return current particle position.
Definition: particleI.H:283
void crossEdgeConnectedFace(const label &celli, label &tetFacei, label &tetPti, const edge &e)
Replacement for particle::crossEdgeConnectedFace that avoids bombing.
autoPtr< particle > clone() const
Construct and return a clone.
Namespace for OpenFOAM.