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 | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::wallBoundedStreamLineParticle
26 
27 Description
28  Particle class that samples fields as it passes through. Used in streamline
29  calculation.
30 
31 SourceFiles
32  wallBoundedStreamLineParticle.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef wallBoundedStreamLineParticle_H
37 #define wallBoundedStreamLineParticle_H
38 
39 #include "wallBoundedParticle.H"
40 #include "autoPtr.H"
41 #include "interpolation.H"
42 #include "vectorList.H"
43 #include "InfoProxy.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 class wallBoundedStreamLineParticleCloud;
51 
52 
53 // Forward declaration of friend functions and operators
54 
55 class wallBoundedStreamLineParticle;
56 
57 Ostream& operator<<(Ostream&, const wallBoundedStreamLineParticle&);
58 
59 
60 /*---------------------------------------------------------------------------*\
61  Class wallBoundedStreamLineParticle Declaration
62 \*---------------------------------------------------------------------------*/
63 
65 :
66  public wallBoundedParticle
67 {
68 
69 public:
70 
71  //- Class used to pass tracking data to the trackToEdge function
72  class trackingData
73  :
75  <
76  Cloud<wallBoundedStreamLineParticle>
77  >
78  {
79 
80  public:
81 
82 
85  const label UIndex_;
86  const bool trackForward_;
87  const scalar trackLength_;
88 
92 
93 
94  // Constructors
95 
97  (
99  const PtrList<interpolation<scalar>>& vsInterp,
100  const PtrList<interpolation<vector>>& vvInterp,
101  const label UIndex,
102  const bool trackForward,
103  const scalar trackLength,
104  const PackedBoolList& isWallPatch,
105 
106  DynamicList<List<point>>& allPositions,
107  List<DynamicList<scalarList>>& allScalars,
108  List<DynamicList<vectorList>>& allVectors
109  )
110  :
112  <
114  >
115  (
116  cloud,
117  isWallPatch
118  ),
119  vsInterp_(vsInterp),
120  vvInterp_(vvInterp),
121  UIndex_(UIndex),
122  trackForward_(trackForward),
123  trackLength_(trackLength),
124 
125  allPositions_(allPositions),
126  allScalars_(allScalars),
127  allVectors_(allVectors)
128  {}
129  };
130 
131 
132 private:
133 
134  // Private data
135 
136  //- Lifetime of particle. Particle dies when reaches 0.
137  label lifeTime_;
138 
139  //- Sampled positions
140  DynamicList<point> sampledPositions_;
141 
142  //- Sampled scalars
143  List<DynamicList<scalar>> sampledScalars_;
144 
145  //- Sampled vectors
146  List<DynamicList<vector>> sampledVectors_;
147 
148 
149  // Private Member Functions
150 
151  vector interpolateFields
152  (
153  const trackingData& td,
154  const point& position,
155  const label celli,
156  const label facei
157  );
158 
159  vector sample(trackingData& td);
160 
161 
162 public:
163 
164  // Constructors
165 
166  //- Construct from components
168  (
169  const polyMesh& c,
170  const vector& position,
171  const label celli,
172  const label tetFacei,
173  const label tetPtI,
174  const label meshEdgeStart,
175  const label diagEdge,
176  const label lifeTime
177  );
178 
179  //- Construct from Istream
181  (
182  const polyMesh& c,
183  Istream& is,
184  bool readFields = true
185  );
186 
187  //- Construct copy
189 
190  //- Construct and return a clone
191  autoPtr<particle> clone() const
192  {
194  }
195 
196  //- Factory class to read-construct particles used for
197  // parallel transfer
198  class iNew
199  {
200  const polyMesh& mesh_;
201 
202  public:
204  iNew(const polyMesh& mesh)
205  :
206  mesh_(mesh)
207  {}
208 
210  (
211  Istream& is
212  ) const
213  {
215  (
216  new wallBoundedStreamLineParticle(mesh_, is, true)
217  );
218  }
219  };
220 
221 
222  // Member Functions
223 
224  // Tracking
225 
226  //- Track all particles to their end point
227  bool move(trackingData&, const scalar trackTime);
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 #endif
259 
260 // ************************************************************************* //
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
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.
const PtrList< interpolation< vector > > & vvInterp_
Factory class to read-construct particles used for.
autoPtr< particle > clone() const
Construct and return a clone.
wallBoundedStreamLineParticle(const polyMesh &c, const vector &position, const label celli, const label tetFacei, const label tetPtI, const label meshEdgeStart, const label diagEdge, const label lifeTime)
Construct from components.
bool move(trackingData &, const scalar trackTime)
Track all particles to their end point.
const vector & position() const
Return current particle position.
Definition: particleI.H:586
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
const polyMesh & mesh_
Reference to the polyMesh database.
Definition: particle.H:137
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:53
const PtrList< interpolation< scalar > > & vsInterp_
static void readFields(Cloud< wallBoundedStreamLineParticle > &)
Read.
label meshEdgeStart() const
-1 or label of mesh edge
A bit-packed bool list.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:62
const dimensionedScalar c
Speed of light in a vacuum.
Ostream & operator<<(Ostream &, const ensightPart &)
trackingData(Cloud< wallBoundedStreamLineParticle > &cloud, const PtrList< interpolation< scalar >> &vsInterp, const PtrList< interpolation< vector >> &vvInterp, const label UIndex, const bool trackForward, const scalar trackLength, const PackedBoolList &isWallPatch, DynamicList< List< point >> &allPositions, List< DynamicList< scalarList >> &allScalars, List< DynamicList< vectorList >> &allVectors)
const polyMesh & mesh() const
Return the mesh database.
Definition: particleI.H:580
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Class used to pass tracking data to the trackToFace function.
label diagEdge() const
-1 or diagonal edge
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.
Namespace for OpenFOAM.