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