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-2022 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 "autoPtr.H"
44 #include "interpolation.H"
45 #include "vectorList.H"
46 #include "InfoProxy.H"
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 namespace Foam
51 {
52 
53 // Forward Declarations
54 class wallBoundedStreamLineParticle;
55 class wallBoundedStreamLineParticleCloud;
56 
57 Ostream& operator<<(Ostream&, const wallBoundedStreamLineParticle&);
58 
59 
60 /*---------------------------------------------------------------------------*\
61  Class wallBoundedStreamLineParticle Declaration
62 \*---------------------------------------------------------------------------*/
63 
65 :
66  public wallBoundedParticle
67 {
68 public:
69 
70  //- Class used to pass tracking data to the trackToEdge function
71  class trackingData
72  :
74  {
75  public:
76 
77  // Public Data
78 
82  const scalar trackLength_;
83 
87 
88 
89  // Constructors
90 
91  template<class TrackCloudType>
93  (
94  TrackCloudType& cloud,
95  const PtrList<interpolation<scalar>>& vsInterp,
96  const PtrList<interpolation<vector>>& vvInterp,
97  const interpolation<vector>& UInterp,
98  const scalar trackLength,
99  const bitSet& isWallPatch,
100 
101  DynamicList<List<point>>& allPositions,
102  List<DynamicList<scalarList>>& allScalars,
103  List<DynamicList<vectorList>>& allVectors
104  )
105  :
106  wallBoundedParticle::trackingData(cloud, isWallPatch),
107  vsInterp_(vsInterp),
108  vvInterp_(vvInterp),
109  UInterp_(UInterp),
110  trackLength_(trackLength),
111 
112  allPositions_(allPositions),
113  allScalars_(allScalars),
114  allVectors_(allVectors)
115  {}
116 
117  virtual ~trackingData() = default;
118  };
119 
120 
121 protected:
122 
123  // Protected Data
124 
125  //- Track with +U or -U
126  bool trackForward_;
127 
128  //- Lifetime of particle. Particle dies when reaches 0.
129  label lifeTime_;
131  //- Sampled positions
133 
134  //- Sampled scalars
136 
137  //- Sampled vectors
139 
141  // Protected Member Functions
142 
144  (
145  const trackingData& td,
146  const point& position,
147  const label celli,
148  const label facei
149  );
150 
152 
153 
154 public:
155 
156  // Constructors
157 
158  //- Construct from components
160  (
161  const polyMesh& c,
162  const point& position,
163  const label celli,
164  const label tetFacei,
165  const label tetPti,
166  const label meshEdgeStart,
167  const label diagEdge,
168  const bool trackForward,
169  const label lifeTime
170  );
171 
172  //- Construct from Istream
174  (
175  const polyMesh& c,
176  Istream& is,
177  bool readFields = true,
178  bool newFormat = true
179  );
180 
181  //- Construct copy
183 
184  //- Construct and return a clone
185  autoPtr<particle> clone() const
186  {
188  }
189 
190  //- Factory class to read-construct particles used for
191  // parallel transfer
192  class iNew
193  {
194  const polyMesh& mesh_;
195 
196  public:
197 
198  iNew(const polyMesh& mesh)
199  :
200  mesh_(mesh)
201  {}
202 
204  (
205  Istream& is
206  ) const
207  {
209  (
210  new wallBoundedStreamLineParticle(mesh_, is, true)
211  );
212  }
213  };
214 
215 
216  // Member Functions
217 
218  // Tracking
219 
220  //- Track all particles to their end point
221  template<class TrackCloudType>
222  bool move
223  (
224  TrackCloudType& cloud,
225  trackingData& td,
226  const scalar trackTime
227  );
228 
229 
230  // I-O
231 
232  //- Read
234 
235  //- Write
236  static void writeFields
237  (
239  );
240 
241 
242  // Ostream Operator
243 
244  friend Ostream& operator<<
245  (
246  Ostream&,
248  );
249 };
250 
251 
252 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
253 
254 } // End namespace Foam
255 
256 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
257 
258 #ifdef NoRepository
260 #endif
261 
262 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
263 
264 #endif
265 
266 // ************************************************************************* //
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
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 ...
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
Construct and 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:74
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.