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-2017 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 // Forward declaration of friend functions and operators
50 
51 class streamLineParticle;
52 
53 Ostream& operator<<(Ostream&, const streamLineParticle&);
54 
55 /*---------------------------------------------------------------------------*\
56  Class streamLineParticle Declaration
57 \*---------------------------------------------------------------------------*/
58 
60 :
61  public particle
62 {
63 public:
64 
65  class trackingData
66  :
67  public particle::TrackingData<Cloud<streamLineParticle>>
68  {
69  public:
70 
71  // Public data
72 
74 
76 
77  const label UIndex_;
78 
79  const bool trackForward_;
80 
81  const label nSubCycle_;
82 
83  const scalar trackLength_;
84 
86 
88 
90 
91 
92  // Constructors
93 
94  //- Construct from components
96  (
98  const PtrList<interpolation<scalar>>& vsInterp,
99  const PtrList<interpolation<vector>>& vvInterp,
100  const label UIndex,
101  const bool trackForward,
102  const label nSubCycle,
103  const scalar trackLength,
104  DynamicList<List<point>>& allPositions,
105  List<DynamicList<scalarList>>& allScalars,
106  List<DynamicList<vectorList>>& allVectors
107  )
108  :
110  vsInterp_(vsInterp),
111  vvInterp_(vvInterp),
112  UIndex_(UIndex),
113  trackForward_(trackForward),
114  nSubCycle_(nSubCycle),
115  trackLength_(trackLength),
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  //- Interpolate all quantities; return interpolated velocity.
143  vector interpolateFields
144  (
145  const trackingData&,
146  const point&,
147  const label celli,
148  const label facei
149  );
150 
151 
152 public:
153 
154  // Constructors
155 
156  //- Construct from components
158  (
159  const polyMesh& c,
160  const vector& position,
161  const label celli,
162  const label lifeTime
163  );
164 
165  //- Construct from Istream
167  (
168  const polyMesh& c,
169  Istream& is,
170  bool readFields = true
171  );
172 
173  //- Construct copy
175 
176  //- Construct and return a clone
177  autoPtr<particle> clone() const
178  {
179  return autoPtr<particle>(new streamLineParticle(*this));
180  }
181 
182  //- Factory class to read-construct particles used for parallel transfer
183  class iNew
184  {
185  const polyMesh& mesh_;
186 
187  public:
189  iNew(const polyMesh& mesh)
190  :
191  mesh_(mesh)
192  {}
194  autoPtr<streamLineParticle> operator()(Istream& is) const
195  {
197  (
198  new streamLineParticle(mesh_, is, true)
199  );
200  }
201  };
202 
203 
204  // Member Functions
205 
206  // Tracking
207 
208  //- Track all particles to their end point
209  bool move(trackingData&, const scalar);
210 
211  //- Overridable function to handle the particle hitting a patch
212  // Executed before other patch-hitting functions
213  bool hitPatch
214  (
215  const polyPatch&,
216  trackingData& td,
217  const label patchi,
218  const scalar trackFraction,
219  const tetIndices& tetIs
220  );
221 
222  //- Overridable function to handle the particle hitting a wedge
223  void hitWedgePatch
224  (
225  const wedgePolyPatch&,
226  trackingData& td
227  );
228 
229  //- Overridable function to handle the particle hitting a
230  // symmetry plane
232  (
233  const symmetryPlanePolyPatch&,
234  trackingData& td
235  );
236 
237  //- Overridable function to handle the particle hitting a
238  // symmetry patch
239  void hitSymmetryPatch
240  (
241  const symmetryPolyPatch&,
242  trackingData& td
243  );
244 
245  //- Overridable function to handle the particle hitting a cyclic
246  void hitCyclicPatch
247  (
248  const cyclicPolyPatch&,
249  trackingData& td
250  );
251 
252  //- Overridable function to handle the particle hitting a
253  //- processorPatch
254  void hitProcessorPatch
255  (
256  const processorPolyPatch&,
257  trackingData& td
258  );
259 
260  //- Overridable function to handle the particle hitting a wallPatch
261  void hitWallPatch
262  (
263  const wallPolyPatch&,
264  trackingData& td,
265  const tetIndices&
266  );
267 
268  //- Overridable function to handle the particle hitting a polyPatch
269  void hitPatch
270  (
271  const polyPatch&,
272  trackingData& td
273  );
274 
275 
276  // I-O
277 
278  //- Read
280 
281  //- Write
282  static void writeFields(const Cloud<streamLineParticle>&);
283 
284 
285  // Ostream Operator
286 
287  friend Ostream& operator<<(Ostream&, const streamLineParticle&);
288 };
289 
290 
291 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
292 
293 } // End namespace Foam
294 
295 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
296 
297 #endif
298 
299 // ************************************************************************* //
Symmetry patch for non-planar or multi-plane patches.
static void writeFields(const Cloud< streamLineParticle > &)
Write.
bool move(trackingData &, const scalar)
Track all particles to their end point.
const polyMesh & mesh() const
Return the mesh database.
Definition: particleI.H:45
List< DynamicList< vectorList > > & allVectors_
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 parallel transfer.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:60
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.
List< DynamicList< scalarList > > & allScalars_
Base particle class.
Definition: particle.H:81
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.
autoPtr< particle > clone() const
Construct and return a clone.
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
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)
Construct from components.
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.
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:63
const dimensionedScalar c
Speed of light in a vacuum.
Ostream & operator<<(Ostream &, const ensightPart &)
friend Ostream & operator<<(Ostream &, const streamLineParticle &)
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:52
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
void hitWedgePatch(const wedgePolyPatch &, trackingData &td)
Overridable function to handle the particle hitting a wedge.
vector position() const
Return current particle position.
Definition: particleI.H:204
void hitWallPatch(const wallPolyPatch &, trackingData &td, const tetIndices &)
Overridable function to handle the particle hitting a wallPatch.
Namespace for OpenFOAM.