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 wallBoundedParticle_H
41 #define 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 declaration of friend functions and operators
53 
54 class wallBoundedParticle;
55 
56 Ostream& operator<<(Ostream&, const wallBoundedParticle&);
57 Ostream& operator<<(Ostream&, const InfoProxy<wallBoundedParticle>&);
58 
59 
60 /*---------------------------------------------------------------------------*\
61  Class wallBoundedParticle Declaration
62 \*---------------------------------------------------------------------------*/
63 
65 :
66  public particle
67 {
68 public:
69 
70  //- Class used to pass tracking data to the trackToFace function
71  class trackingData
72  :
74  {
75 
76  public:
77 
78  const bitSet& isWallPatch_;
79 
80  // Constructors
81 
82  template <class TrackCloudType>
83  inline trackingData
84  (
85  const TrackCloudType& cloud,
86  const bitSet& isWallPatch
87  )
88  :
90  isWallPatch_(isWallPatch)
91  {}
92  };
93 
94 
95 protected:
96 
97  // Protected data
98 
99  //- Particle position is updated locally as opposed to via track
100  // functions of the base Foam::particle class
102 
103  //- Particle is on mesh edge:
104  // const face& f = mesh.faces()[tetFace()]
105  // const edge e(f[meshEdgeStart_], f.nextLabel(meshEdgeStart_));
106  // Note that this real edge
107  // is also one of the edges of the face-triangle (from
108  // tetFace()+tetPt()).
109  label meshEdgeStart_;
110 
111  //- Particle is on diagonal edge:
112  // const face& f = mesh.faces()[tetFace()]
113  // label faceBasePti = mesh.tetBasePtIs()[facei];
114  // label diagPti = (faceBasePti+diagEdge_)%f.size();
115  // const edge e(f[faceBasePti], f[diagPti]);
116  label diagEdge_;
117 
118 
119  // Protected Member Functions
120 
121  //- Construct current edge
122  edge currentEdge() const;
123 
124  //- Replacement for particle::crossEdgeConnectedFace that avoids bombing
125  // out on invalid tet decomposition (tetBasePtIs = -1)
127  (
128  const label& celli,
129  label& tetFacei,
130  label& tetPti,
131  const edge& e
132  );
133 
134  //- Cross mesh edge into different face on same cell
135  void crossEdgeConnectedFace(const edge& meshEdge);
136 
137  //- Cross diagonal edge into different triangle on same face,cell
138  void crossDiagonalEdge();
139 
140  //- Track through single triangle
141  scalar trackFaceTri(const vector& n, const vector& endPosition, label&);
142 
143  //- Is current triangle in the track direction
144  bool isTriAlongTrack(const vector& n, const point& endPosition) const;
145 
146 
147 public:
148 
149  // Static Data Members
150 
151  //- Size in bytes of the fields
152  static const std::size_t sizeofFields_;
153 
154 
155  // Constructors
156 
157  //- Construct from components
159  (
160  const polyMesh& c,
161  const point& position,
162  const label celli,
163  const label tetFacei,
164  const label tetPti,
165  const label meshEdgeStart,
166  const label diagEdge
167  );
168 
169  //- Construct from Istream
171  (
172  const polyMesh& c,
173  Istream& is,
174  bool readFields = true,
175  bool newFormat = true
176  );
177 
178  //- Construct copy
180 
181  //- Construct and return a clone
182  autoPtr<particle> clone() const
183  {
184  return autoPtr<particle>(new wallBoundedParticle(*this));
185  }
186 
187  //- Factory class to read-construct particles used for
188  // parallel transfer
189  class iNew
190  {
191  const polyMesh& mesh_;
192 
193  public:
194 
195  iNew(const polyMesh& mesh)
196  :
197  mesh_(mesh)
198  {}
199 
201  (
202  Istream& is
203  ) const
204  {
206  (
207  new wallBoundedParticle(mesh_, is, true)
208  );
209  }
210  };
211 
212 
213  // Member Functions
214 
215  // Access
216 
217  //- -1 or label of mesh edge
218  inline label meshEdgeStart() const
219  {
220  return meshEdgeStart_;
221  }
222 
223  //- -1 or diagonal edge
224  inline label diagEdge() const
225  {
226  return diagEdge_;
227  }
228 
229 
230  // Track
231 
232  //- Equivalent of trackToFace
233  template<class TrackCloudType>
234  scalar trackToEdge
235  (
236  TrackCloudType& cloud,
237  trackingData& td,
238  const vector& endPosition
239  );
240 
241 
242  // Patch interactions
243 
244  //- Do all patch interaction
245  template<class TrackCloudType>
246  void patchInteraction
247  (
248  TrackCloudType& cloud,
249  trackingData& td,
250  const scalar trackFraction
251  );
252 
253  //- Overridable function to handle the particle hitting a
254  //- processorPatch
255  template<class TrackCloudType>
256  void hitProcessorPatch
257  (
258  TrackCloudType& cloud,
259  trackingData& td
260  );
261 
262  //- Overridable function to handle the particle hitting a wallPatch
263  template<class TrackCloudType>
264  void hitWallPatch
265  (
266  TrackCloudType& cloud,
267  trackingData& td
268  );
269 
270 
271  // Info
272 
273  //- Return info proxy.
274  // Used to print particle information to a stream
275  inline InfoProxy<wallBoundedParticle> info() const
276  {
277  return *this;
278  }
279 
280 
281  // I-O
282 
283  //- Read
284  template<class CloudType>
285  static void readFields(CloudType&);
286 
287  //- Write
288  template<class CloudType>
289  static void writeFields(const CloudType&);
290 
291 
292  // Ostream Operator
293 
294  friend Ostream& operator<<
295  (
296  Ostream&,
297  const wallBoundedParticle&
298  );
299 
300  friend Ostream& operator<<
301  (
302  Ostream&,
304  );
305 };
306 
307 
308 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
309 
310 } // End namespace Foam
311 
312 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
313 
314 #ifdef NoRepository
316 #endif
317 
318 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
319 
320 #endif
321 
322 // ************************************************************************* //
Foam::wallBoundedParticle::trackingData::isWallPatch_
const bitSet & isWallPatch_
Definition: wallBoundedParticle.H:77
wallBoundedParticleTemplates.C
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::wallBoundedParticle::sizeofFields_
static const std::size_t sizeofFields_
Size in bytes of the fields.
Definition: wallBoundedParticle.H:151
Foam::InfoProxy
A helper class for outputting values to Ostream.
Definition: InfoProxy.H:47
Foam::wallBoundedParticle::hitProcessorPatch
void hitProcessorPatch(TrackCloudType &cloud, trackingData &td)
Definition: wallBoundedParticleTemplates.C:382
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:63
Foam::wallBoundedParticle::patchInteraction
void patchInteraction(TrackCloudType &cloud, trackingData &td, const scalar trackFraction)
Do all patch interaction.
Definition: wallBoundedParticleTemplates.C:35
Foam::wallBoundedParticle::wallBoundedParticle
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.
Definition: wallBoundedParticle.C:412
Foam::wallBoundedParticle::readFields
static void readFields(CloudType &)
Read.
InfoProxy.H
Foam::wallBoundedParticle::hitWallPatch
void hitWallPatch(TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a wallPatch.
Definition: wallBoundedParticleTemplates.C:411
Foam::wallBoundedParticle::writeFields
static void writeFields(const CloudType &)
Write.
Foam::wallBoundedParticle::iNew
Factory class to read-construct particles used for.
Definition: wallBoundedParticle.H:188
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::wallBoundedParticle::meshEdgeStart_
label meshEdgeStart_
Particle is on mesh edge:
Definition: wallBoundedParticle.H:108
Foam::wallBoundedParticle::crossDiagonalEdge
void crossDiagonalEdge()
Cross diagonal edge into different triangle on same face,cell.
Definition: wallBoundedParticle.C:239
Foam::wallBoundedParticle::info
InfoProxy< wallBoundedParticle > info() const
Return info proxy.
Definition: wallBoundedParticle.H:274
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::wallBoundedParticle::trackToEdge
scalar trackToEdge(TrackCloudType &cloud, trackingData &td, const vector &endPosition)
Equivalent of trackToFace.
Foam::operator<<
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::wallBoundedParticle::meshEdgeStart
label meshEdgeStart() const
-1 or label of mesh edge
Definition: wallBoundedParticle.H:217
Foam::wallBoundedParticle::crossEdgeConnectedFace
void crossEdgeConnectedFace(const label &celli, label &tetFacei, label &tetPti, const edge &e)
Replacement for particle::crossEdgeConnectedFace that avoids bombing.
Definition: wallBoundedParticle.C:83
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::wallBoundedParticle
Particle class that tracks on triangles of boundary faces. Use trackToEdge similar to trackToFace on ...
Definition: wallBoundedParticle.H:63
Foam::particle::position
vector position() const
Return current particle position.
Definition: particleI.H:314
Foam::particle::mesh
const polyMesh & mesh() const
Return the mesh database.
Definition: particleI.H:137
Foam::CloudType
DSMCCloud< dsmcParcel > CloudType
Definition: makeDSMCParcelBinaryCollisionModels.C:38
Foam::wallBoundedParticle::iNew::iNew
iNew(const polyMesh &mesh)
Definition: wallBoundedParticle.H:194
Foam::wallBoundedParticle::trackingData::trackingData
trackingData(const TrackCloudType &cloud, const bitSet &isWallPatch)
Definition: wallBoundedParticle.H:83
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:71
Foam::wallBoundedParticle::isTriAlongTrack
bool isTriAlongTrack(const vector &n, const point &endPosition) const
Is current triangle in the track direction.
Definition: wallBoundedParticle.C:365
Foam::wallBoundedParticle::clone
autoPtr< particle > clone() const
Construct and return a clone.
Definition: wallBoundedParticle.H:181
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::wallBoundedParticle::trackingData
Class used to pass tracking data to the trackToFace function.
Definition: wallBoundedParticle.H:70
Foam::wallBoundedParticle::localPosition_
point localPosition_
Particle position is updated locally as opposed to via track.
Definition: wallBoundedParticle.H:100
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::cloud
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:57
Foam::Vector< scalar >
Foam::particle
Base particle class.
Definition: particle.H:76
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
particle.H
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:223
Foam::wallBoundedParticle::currentEdge
edge currentEdge() const
Construct current edge.
Definition: wallBoundedParticle.C:41
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::wallBoundedParticle::trackFaceTri
scalar trackFaceTri(const vector &n, const vector &endPosition, label &)
Track through single triangle.
Definition: wallBoundedParticle.C:286
Foam::particle::trackingData
Definition: particle.H:95
Foam::wallBoundedParticle::diagEdge_
label diagEdge_
Particle is on diagonal edge:
Definition: wallBoundedParticle.H:115
autoPtr.H