streamLineParticle.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::streamLineParticle
26 
27 Description
28  Particle class that samples fields as it passes through. Used in streamline
29  calculation.
30 
31 SourceFiles
32  streamLineParticle.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef streamLineParticle_H
37 #define streamLineParticle_H
38 
39 #include "particle.H"
40 #include "autoPtr.H"
41 #include "interpolation.H"
42 #include "vectorList.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 class streamLineParticleCloud;
50 
51 
52 // Forward declaration of friend functions and operators
53 
54 class streamLineParticle;
55 
56 Ostream& operator<<(Ostream&, const streamLineParticle&);
57 
58 
59 /*---------------------------------------------------------------------------*\
60  Class streamLineParticle Declaration
61 \*---------------------------------------------------------------------------*/
62 
64 :
65  public particle
66 {
67 
68 public:
69 
70  //- Class used to pass tracking data to the trackToFace function
71  class trackingData
72  :
73  public particle::TrackingData<Cloud<streamLineParticle>>
74  {
75 
76  public:
77 
78 
81  const label UIndex_;
82  const bool trackForward_;
84  const scalar trackLength_;
85 
89 
90 
91  // Constructors
92 
94  (
96  const PtrList<interpolation<scalar>>& vsInterp,
97  const PtrList<interpolation<vector>>& vvInterp,
98  const label UIndex,
99  const bool trackForward,
100  const label nSubCycle,
101  const scalar trackLength,
102 
103  DynamicList<List<point>>& allPositions,
104  List<DynamicList<scalarList>>& allScalars,
105  List<DynamicList<vectorList>>& allVectors
106  )
107  :
109  vsInterp_(vsInterp),
110  vvInterp_(vvInterp),
111  UIndex_(UIndex),
112  trackForward_(trackForward),
113  nSubCycle_(nSubCycle),
114  trackLength_(trackLength),
115 
116  allPositions_(allPositions),
117  allScalars_(allScalars),
118  allVectors_(allVectors)
119  {}
120  };
121 
122 
123 private:
124 
125  // Private data
126 
127  //- Lifetime of particle. Particle dies when reaches 0.
128  label lifeTime_;
129 
130  //- Sampled positions
131  DynamicList<point> sampledPositions_;
132 
133  //- Sampled scalars
134  List<DynamicList<scalar>> sampledScalars_;
135 
136  //- Sampled vectors
137  List<DynamicList<vector>> sampledVectors_;
138 
139 
140  // Private Member Functions
141 
142  //- Estimate dt to cross from current face to next one in nSubCycle
143  // steps.
144  scalar calcSubCycleDeltaT
145  (
146  trackingData& td,
147  const scalar dt,
148  const vector& U
149  ) const;
150 
151  void constrainVelocity
152  (
153  trackingData& td,
154  const scalar dt,
155  vector& U
156  );
157 
158  //- Interpolate all quantities; return interpolated velocity.
159  vector interpolateFields
160  (
161  const trackingData&,
162  const point&,
163  const label celli,
164  const label facei
165  );
166 
167 
168 public:
169 
170  // Constructors
171 
172  //- Construct from components
174  (
175  const polyMesh& c,
176  const vector& position,
177  const label celli,
178  const label lifeTime
179  );
180 
181  //- Construct from Istream
183  (
184  const polyMesh& c,
185  Istream& is,
186  bool readFields = true
187  );
188 
189  //- Construct copy
191 
192  //- Construct and return a clone
193  autoPtr<particle> clone() const
194  {
195  return autoPtr<particle>(new streamLineParticle(*this));
196  }
197 
198  //- Factory class to read-construct particles used for
199  // parallel transfer
200  class iNew
201  {
202  const polyMesh& mesh_;
203 
204  public:
206  iNew(const polyMesh& mesh)
207  :
208  mesh_(mesh)
209  {}
211  autoPtr<streamLineParticle> operator()(Istream& is) const
212  {
214  (
215  new streamLineParticle(mesh_, is, true)
216  );
217  }
218  };
219 
220 
221  // Member Functions
222 
223  // Tracking
224 
225  //- Track all particles to their end point
226  bool move(trackingData&, const scalar trackTime);
227 
228 
229  //- Overridable function to handle the particle hitting a patch
230  // Executed before other patch-hitting functions
231  bool hitPatch
232  (
233  const polyPatch&,
234  trackingData& td,
235  const label patchi,
236  const scalar trackFraction,
237  const tetIndices& tetIs
238  );
239 
240  //- Overridable function to handle the particle hitting a wedge
241  void hitWedgePatch
242  (
243  const wedgePolyPatch&,
244  trackingData& td
245  );
246 
247  //- Overridable function to handle the particle hitting a
248  // symmetry plane
250  (
251  const symmetryPlanePolyPatch&,
252  trackingData& td
253  );
254 
255  //- Overridable function to handle the particle hitting a
256  // symmetry patch
257  void hitSymmetryPatch
258  (
259  const symmetryPolyPatch&,
260  trackingData& td
261  );
262 
263  //- Overridable function to handle the particle hitting a cyclic
264  void hitCyclicPatch
265  (
266  const cyclicPolyPatch&,
267  trackingData& td
268  );
269 
270  //- Overridable function to handle the particle hitting a
271  //- processorPatch
272  void hitProcessorPatch
273  (
274  const processorPolyPatch&,
275  trackingData& td
276  );
277 
278  //- Overridable function to handle the particle hitting a wallPatch
279  void hitWallPatch
280  (
281  const wallPolyPatch&,
282  trackingData& td,
283  const tetIndices&
284  );
285 
286  //- Overridable function to handle the particle hitting a polyPatch
287  void hitPatch
288  (
289  const polyPatch&,
290  trackingData& td
291  );
292 
293 
294  // I-O
295 
296  //- Read
298 
299  //- Write
300  static void writeFields(const Cloud<streamLineParticle>&);
301 
302 
303  // Ostream Operator
304 
305  friend Ostream& operator<<(Ostream&, const streamLineParticle&);
306 };
307 
308 
309 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
310 
311 } // End namespace Foam
312 
313 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
314 
315 #endif
316 
317 // ************************************************************************* //
Symmetry patch for non-planar or multi-plane patches.
static void writeFields(const Cloud< streamLineParticle > &)
Write.
List< DynamicList< vectorList > > & allVectors_
U
Definition: pEqn.H:83
void hitSymmetryPlanePatch(const symmetryPlanePolyPatch &, trackingData &td)
Overridable function to handle the particle hitting a.
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
Factory class to read-construct particles used for.
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
const PtrList< interpolation< vector > > & vvInterp_
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
void hitCyclicPatch(const cyclicPolyPatch &, trackingData &td)
Overridable function to handle the particle hitting a cyclic.
Class used to pass tracking data to the trackToFace function.
List< DynamicList< scalarList > > & allScalars_
Base particle class.
Definition: particle.H:78
bool hitPatch(const polyPatch &, trackingData &td, const label patchi, const scalar trackFraction, const tetIndices &tetIs)
Overridable function to handle the particle hitting a patch.
Neighbour processor patch.
const vector & position() const
Return current particle position.
Definition: particleI.H:586
void hitSymmetryPatch(const symmetryPolyPatch &, trackingData &td)
Overridable function to handle the particle hitting a.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
DynamicList< vectorList > & allPositions_
Cyclic plane patch.
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
Particle class that samples fields as it passes through. Used in streamline calculation.
Foam::wallPolyPatch.
Definition: wallPolyPatch.H:48
const polyMesh & mesh_
Reference to the polyMesh database.
Definition: particle.H:137
Wedge front and back plane patch.
trackingData(Cloud< streamLineParticle > &cloud, const PtrList< interpolation< scalar >> &vsInterp, const PtrList< interpolation< vector >> &vvInterp, const label UIndex, const bool trackForward, const label nSubCycle, const scalar trackLength, DynamicList< List< point >> &allPositions, List< DynamicList< scalarList >> &allScalars, List< DynamicList< vectorList >> &allVectors)
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:81
static void readFields(Cloud< streamLineParticle > &)
Read.
autoPtr< particle > clone() const
Construct and return a clone.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
const PtrList< interpolation< scalar > > & vsInterp_
streamLineParticle(const polyMesh &c, const vector &position, const label celli, const label lifeTime)
Construct from components.
label patchi
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 &)
friend Ostream & operator<<(Ostream &, const streamLineParticle &)
const polyMesh & mesh() const
Return the mesh database.
Definition: particleI.H:580
void hitProcessorPatch(const processorPolyPatch &, trackingData &td)
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
volScalarField & p
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
bool move(trackingData &, const scalar trackTime)
Track all particles to their end point.
void hitWedgePatch(const wedgePolyPatch &, trackingData &td)
Overridable function to handle the particle hitting a wedge.
void hitWallPatch(const wallPolyPatch &, trackingData &td, const tetIndices &)
Overridable function to handle the particle hitting a wallPatch.
Namespace for OpenFOAM.