wallBoundedStreamLineParticle.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 -------------------------------------------------------------------------------
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::wallBoundedStreamLineParticle
29 
30 Description
31  Particle class that samples fields as it passes through. Used in streamline
32  calculation.
33 
34 SourceFiles
35  wallBoundedStreamLineParticle.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef Foam_wallBoundedStreamLineParticle_H
40 #define Foam_wallBoundedStreamLineParticle_H
41 
42 #include "wallBoundedParticle.H"
43 #include "interpolation.H"
44 #include "vectorList.H"
45 #include "InfoProxy.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 // Forward Declarations
53 class wallBoundedStreamLineParticle;
54 class wallBoundedStreamLineParticleCloud;
55 
56 Ostream& operator<<(Ostream&, const wallBoundedStreamLineParticle&);
57 
58 
59 /*---------------------------------------------------------------------------*\
60  Class wallBoundedStreamLineParticle Declaration
61 \*---------------------------------------------------------------------------*/
62 
64 :
65  public wallBoundedParticle
66 {
67 public:
68 
69  //- Class used to pass tracking data to the trackToEdge function
70  class trackingData
71  :
73  {
74  public:
75 
76  // Public Data
77 
81  const scalar trackLength_;
82 
86 
87 
88  // Constructors
89 
90  template<class TrackCloudType>
92  (
93  TrackCloudType& cloud,
94  const PtrList<interpolation<scalar>>& vsInterp,
95  const PtrList<interpolation<vector>>& vvInterp,
96  const interpolation<vector>& UInterp,
97  const scalar trackLength,
98  const bitSet& isWallPatch,
99 
100  DynamicList<List<point>>& allPositions,
101  List<DynamicList<scalarList>>& allScalars,
102  List<DynamicList<vectorList>>& allVectors
103  )
104  :
105  wallBoundedParticle::trackingData(cloud, isWallPatch),
106  vsInterp_(vsInterp),
107  vvInterp_(vvInterp),
108  UInterp_(UInterp),
109  trackLength_(trackLength),
110 
111  allPositions_(allPositions),
112  allScalars_(allScalars),
113  allVectors_(allVectors)
114  {}
115 
116  virtual ~trackingData() = default;
117  };
118 
119 
120 protected:
121 
122  // Protected Data
123 
124  //- Track with +U or -U
125  bool trackForward_;
126 
127  //- Lifetime of particle. Particle dies when reaches 0.
128  label lifeTime_;
130  //- Sampled positions
132 
133  //- Sampled scalars
135 
136  //- Sampled vectors
138 
140  // Protected Member Functions
141 
143  (
144  const trackingData& td,
145  const point& position,
146  const label celli,
147  const label facei
148  );
149 
151 
152 
153 public:
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  const bool trackForward,
168  const label lifeTime
169  );
170 
171  //- Construct from Istream
173  (
174  const polyMesh& c,
175  Istream& is,
176  bool readFields = true,
177  bool newFormat = true
178  );
179 
180  //- Construct copy
182 
183  //- Return a clone
184  autoPtr<particle> clone() const
185  {
186  return particle::Clone(*this);
187  }
188 
189  //- Factory class to read-construct particles (for parallel transfer)
190  class iNew
191  {
192  const polyMesh& mesh_;
193 
194  public:
195 
196  iNew(const polyMesh& mesh)
197  :
198  mesh_(mesh)
199  {}
200 
202  (
203  Istream& is
204  ) const
205  {
207  (
208  new wallBoundedStreamLineParticle(mesh_, is, true)
209  );
210  }
211  };
212 
214  // Member Functions
215 
216  // Tracking
217 
218  //- Track all particles to their end point
219  template<class TrackCloudType>
220  bool move
221  (
222  TrackCloudType& cloud,
223  trackingData& td,
224  const scalar trackTime
225  );
226 
227 
228  // I-O
229 
230  //- Read
232 
233  //- Write
234  static void writeFields
235  (
237  );
238 
239 
240  // Ostream Operator
241 
242  friend Ostream& operator<<
243  (
244  Ostream&,
246  );
247 };
248 
249 
250 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
251 
252 } // End namespace Foam
253 
254 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
255 
256 #ifdef NoRepository
258 #endif
259 
260 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
261 
262 #endif
263 
264 // ************************************************************************* //
wallBoundedStreamLineParticle(const polyMesh &c, const point &position, const label celli, const label tetFacei, const label tetPti, const label meshEdgeStart, const label diagEdge, const bool trackForward, const label lifeTime)
Construct from components.
Class used to pass tracking data to the trackToFace function.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
static void writeFields(const Cloud< wallBoundedStreamLineParticle > &)
Write.
label diagEdge() const noexcept
The diagonal edge label or -1.
const PtrList< interpolation< vector > > & vvInterp_
label lifeTime_
Lifetime of particle. Particle dies when reaches 0.
List< DynamicList< scalar > > sampledScalars_
Sampled scalars.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:51
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:53
bool move(TrackCloudType &cloud, trackingData &td, const scalar trackTime)
Track all particles to their end point.
Particle class that tracks on triangles of boundary faces. Use trackToEdge similar to trackToFace on ...
static autoPtr< particle > Clone(const Derived &p)
Clone a particle.
Definition: particle.H:552
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const PtrList< interpolation< scalar > > & vsInterp_
DynamicList< point > sampledPositions_
Sampled positions.
autoPtr< particle > clone() const
Return a clone.
static void readFields(Cloud< wallBoundedStreamLineParticle > &)
Read.
trackingData(TrackCloudType &cloud, const PtrList< interpolation< scalar >> &vsInterp, const PtrList< interpolation< vector >> &vvInterp, const interpolation< vector > &UInterp, const scalar trackLength, const bitSet &isWallPatch, DynamicList< List< point >> &allPositions, List< DynamicList< scalarList >> &allScalars, List< DynamicList< vectorList >> &allVectors)
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
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
const dimensionedScalar c
Speed of light in a vacuum.
label meshEdgeStart() const noexcept
The mesh edge label or -1.
Abstract base class for volume field interpolation.
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:75
Class used to pass tracking data to the trackToEdge function.
volScalarField & p
Particle class that samples fields as it passes through. Used in streamline calculation.
List< DynamicList< vector > > sampledVectors_
Sampled vectors.
vector interpolateFields(const trackingData &td, const point &position, const label celli, const label facei)
vector position() const
Return current particle position.
Definition: particleI.H:283
Namespace for OpenFOAM.