streamLineParticle.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 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 class streamLineParticleCloud;
53 
54 Ostream& operator<<(Ostream&, const streamLineParticle&);
55 
56 /*---------------------------------------------------------------------------*\
57  Class streamLineParticle Declaration
58 \*---------------------------------------------------------------------------*/
59 
61 :
62  public particle
63 {
64 public:
65 
66  class trackingData
67  :
69  {
70  public:
71 
72  // Public data
73 
75 
77 
78  const label UIndex_;
79 
80  bool trackForward_;
81 
82  const label nSubCycle_;
83 
84  const scalar trackLength_;
85 
87 
89 
91 
92 
93  // Constructors
94 
95  //- Construct from components
97  (
99  const PtrList<interpolation<scalar>>& vsInterp,
100  const PtrList<interpolation<vector>>& vvInterp,
101  const label UIndex,
102  const bool trackForward,
103  const label nSubCycle,
104  const scalar trackLength,
105  DynamicList<List<point>>& allPositions,
106  List<DynamicList<scalarList>>& allScalars,
107  List<DynamicList<vectorList>>& allVectors
108  )
109  :
110  particle::trackingData(cloud),
111  vsInterp_(vsInterp),
112  vvInterp_(vvInterp),
113  UIndex_(UIndex),
114  trackForward_(trackForward),
115  nSubCycle_(nSubCycle),
116  trackLength_(trackLength),
117  allPositions_(allPositions),
118  allScalars_(allScalars),
119  allVectors_(allVectors)
120  {}
121  };
122 
123 
124 private:
125 
126  // Private Data
127 
128  //- Lifetime of particle. Particle dies when reaches 0.
129  label lifeTime_;
130 
131  //- Sampled positions
132  DynamicList<point> sampledPositions_;
133 
134  //- Sampled scalars
135  List<DynamicList<scalar>> sampledScalars_;
136 
137  //- Sampled vectors
138  List<DynamicList<vector>> sampledVectors_;
139 
140 
141  // Private Member Functions
142 
143  //- Interpolate all quantities; return interpolated velocity.
144  vector interpolateFields
145  (
146  const trackingData&,
147  const point&,
148  const label celli,
149  const label facei
150  );
151 
152 
153 public:
154 
155  // Constructors
156 
157  //- Construct from components
159  (
160  const polyMesh& c,
161  const vector& position,
162  const label celli,
163  const label lifeTime
164  );
165 
166  //- Construct from Istream
168  (
169  const polyMesh& c,
170  Istream& is,
171  bool readFields = true
172  );
173 
174  //- Construct copy
176 
177  //- Construct and return a clone
178  autoPtr<particle> clone() const
179  {
180  return autoPtr<particle>(new streamLineParticle(*this));
181  }
182 
183  //- Factory class to read-construct particles used for parallel transfer
184  class iNew
185  {
186  const polyMesh& mesh_;
187 
188  public:
190  iNew(const polyMesh& mesh)
191  :
192  mesh_(mesh)
193  {}
195  autoPtr<streamLineParticle> operator()(Istream& is) const
196  {
198  (
199  new streamLineParticle(mesh_, is, true)
200  );
201  }
202  };
203 
204 
205  // Member Functions
206 
207  // Tracking
208 
209  //- Track all particles to their end point
210  bool move(streamLineParticleCloud&, trackingData&, const scalar);
211 
212  //- Overridable function to handle the particle hitting a patch
213  // Executed before other patch-hitting functions
215 
216  //- Overridable function to handle the particle hitting a wedge
218 
219  //- Overridable function to handle the particle hitting a
220  // symmetry plane
222 
223  //- Overridable function to handle the particle hitting a
224  // symmetry patch
226 
227  //- Overridable function to handle the particle hitting a cyclic
229 
230  //- Overridable function to handle the particle hitting a
231  // cyclicAMIPatch
232  void hitCyclicAMIPatch
233  (
234  const vector&,
235  const scalar,
237  trackingData&
238  );
239 
240  //- Overridable function to handle the particle hitting a
241  // cyclicACMIPatch
242  void hitCyclicACMIPatch
243  (
244  const vector&,
245  const scalar,
247  trackingData&
248  );
249 
250  //- Overridable function to handle the particle hitting a
251  // cyclicRepeatAMIPatch
253  (
254  const vector&,
255  const scalar,
257  trackingData&
258  );
259 
260  //- Overridable function to handle the particle hitting a
261  //- processorPatch
263 
264  //- Overridable function to handle the particle hitting a wallPatch
266 
267 
268  // I-O
269 
270  //- Read
272 
273  //- Write
274  static void writeFields(const Cloud<streamLineParticle>&);
275 
276 
277  // Ostream Operator
278 
279  friend Ostream& operator<<(Ostream&, const streamLineParticle&);
280 };
281 
282 
283 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
284 
285 } // End namespace Foam
286 
287 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
288 
289 #endif
290 
291 // ************************************************************************* //
void hitProcessorPatch(streamLineParticleCloud &, trackingData &)
static void writeFields(const Cloud< streamLineParticle > &)
Write.
const polyMesh & mesh() const
Return the mesh database.
Definition: particleI.H:125
List< DynamicList< vectorList > > & allVectors_
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
void hitWedgePatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a wedge.
Factory class to read-construct particles used for parallel transfer.
trackingData(streamLineParticleCloud &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.
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
const dimensionedScalar & c
Speed of light in a vacuum.
bool hitPatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a patch.
List< DynamicList< scalarList > > & allScalars_
Base particle class.
Definition: particle.H:84
autoPtr< particle > clone() const
Construct and return a clone.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
void hitSymmetryPlanePatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a.
DynamicList< vectorList > & allPositions_
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.
void hitCyclicAMIPatch(const vector &, const scalar, streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a.
static void readFields(Cloud< streamLineParticle > &)
Read.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
const PtrList< interpolation< scalar > > & vsInterp_
void hitCyclicRepeatAMIPatch(const vector &, const scalar, streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a.
streamLineParticle(const polyMesh &c, const vector &position, const label celli, const label lifeTime)
Construct from components.
A Cloud of streamLine particles.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
Ostream & operator<<(Ostream &, const ensightPart &)
bool move(streamLineParticleCloud &, trackingData &, const scalar)
Track all particles to their end point.
friend Ostream & operator<<(Ostream &, const streamLineParticle &)
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
void hitCyclicACMIPatch(const vector &, const scalar, streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a.
volScalarField & p
void hitWallPatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a wallPatch.
vector position() const
Return current particle position.
Definition: particleI.H:278
void hitSymmetryPatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a.
void hitCyclicPatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a cyclic.
Namespace for OpenFOAM.