streamlinesParticle.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-2021 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::streamlinesParticle
26 
27 Description
28  Particle class that samples fields as it passes through. Used in streamlines
29  calculation.
30 
31 SourceFiles
32  streamlinesParticle.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef streamlinesParticle_H
37 #define streamlinesParticle_H
38 
39 #include "particle.H"
40 #include "Cloud.H"
41 #include "autoPtr.H"
42 #include "interpolation.H"
43 #include "vectorList.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 // Forward declaration of friend functions and operators
51 
52 class streamlinesParticle;
53 class streamlinesCloud;
54 
55 Ostream& operator<<(Ostream&, const streamlinesParticle&);
56 
57 /*---------------------------------------------------------------------------*\
58  Class streamlinesParticle Declaration
59 \*---------------------------------------------------------------------------*/
60 
62 :
63  public particle
64 {
65 public:
66 
67  class trackingData
68  :
70  {
71  public:
72 
73  // Public data
74 
76 
78 
79  const label UIndex_;
80 
81  bool trackForward_;
82 
83  const label nSubCycle_;
84 
85  const scalar trackLength_;
86 
88 
90 
92 
94 
95 
96  // Constructors
97 
98  //- Construct from components
100  (
102  const PtrList<interpolation<scalar>>& vsInterp,
103  const PtrList<interpolation<vector>>& vvInterp,
104  const label UIndex,
105  const bool trackForward,
106  const label nSubCycle,
107  const scalar trackLength,
108  DynamicList<List<point>>& allPositions,
109  DynamicList<List<scalar>>& allTimes,
110  List<DynamicList<scalarList>>& allScalars,
111  List<DynamicList<vectorList>>& allVectors
112  )
113  :
114  particle::trackingData(cloud),
115  vsInterp_(vsInterp),
116  vvInterp_(vvInterp),
117  UIndex_(UIndex),
118  trackForward_(trackForward),
119  nSubCycle_(nSubCycle),
120  trackLength_(trackLength),
121  allPositions_(allPositions),
122  allTimes_(allTimes),
123  allScalars_(allScalars),
124  allVectors_(allVectors)
125  {}
126  };
127 
128 
129 private:
130 
131  // Private Data
132 
133  //- Lifetime of particle. Particle dies when reaches 0.
134  label lifeTime_;
135 
136  //- Age of the particle
137  scalar age_;
138 
139  //- Sampled positions
140  DynamicList<point> sampledPositions_;
141 
142  //- Sampled times
143  DynamicList<scalar> sampledTimes_;
144 
145  //- Sampled scalars
146  List<DynamicList<scalar>> sampledScalars_;
147 
148  //- Sampled vectors
149  List<DynamicList<vector>> sampledVectors_;
150 
151 
152  // Private Member Functions
153 
154  //- Interpolate all quantities; return interpolated velocity.
155  vector interpolateFields
156  (
157  const trackingData&,
158  const point&,
159  const label celli,
160  const label facei
161  );
162 
163 
164 public:
165 
166  // Constructors
167 
168  //- Construct from components
170  (
171  const polyMesh& c,
172  const vector& position,
173  const label celli,
174  const label lifeTime
175  );
176 
177  //- Construct from Istream
179  (
180  const polyMesh& c,
181  Istream& is,
182  bool readFields = true
183  );
184 
185  //- Construct copy
187 
188  //- Construct and return a clone
189  autoPtr<particle> clone() const
190  {
191  return autoPtr<particle>(new streamlinesParticle(*this));
192  }
193 
194  //- Factory class to read-construct particles used for parallel transfer
195  class iNew
196  {
197  const polyMesh& mesh_;
198 
199  public:
201  iNew(const polyMesh& mesh)
202  :
203  mesh_(mesh)
204  {}
206  autoPtr<streamlinesParticle> operator()(Istream& is) const
207  {
209  (
210  new streamlinesParticle(mesh_, is, true)
211  );
212  }
213  };
214 
215 
216  // Member Functions
217 
218  // Tracking
219 
220  //- Track all particles to their end point
221  bool move(streamlinesCloud&, trackingData&, const scalar);
222 
223  //- Overridable function to handle the particle hitting a patch
224  // Executed before other patch-hitting functions
226 
227  //- Overridable function to handle the particle hitting a wedge
229 
230  //- Overridable function to handle the particle hitting a
231  // symmetry plane
233 
234  //- Overridable function to handle the particle hitting a
235  // symmetry patch
237 
238  //- Overridable function to handle the particle hitting a cyclic
240 
241  //- Overridable function to handle the particle hitting a
242  // cyclicAMIPatch
243  void hitCyclicAMIPatch
244  (
245  const vector&,
246  const scalar,
248  trackingData&
249  );
250 
251  //- Overridable function to handle the particle hitting a
252  // cyclicACMIPatch
253  void hitCyclicACMIPatch
254  (
255  const vector&,
256  const scalar,
258  trackingData&
259  );
260 
261  //- Overridable function to handle the particle hitting a
262  // cyclicRepeatAMIPatch
264  (
265  const vector&,
266  const scalar,
268  trackingData&
269  );
270 
271  //- Overridable function to handle the particle hitting a
272  //- processorPatch
274 
275  //- Overridable function to handle the particle hitting a wallPatch
277 
278 
279  // I-O
280 
281  //- Read
283 
284  //- Write
285  static void writeFields(const Cloud<streamlinesParticle>&);
286 
287 
288  // Ostream Operator
289 
290  friend Ostream& operator<<(Ostream&, const streamlinesParticle&);
291 };
292 
293 
294 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
295 
296 } // End namespace Foam
297 
298 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
299 
300 #endif
301 
302 // ************************************************************************* //
void hitCyclicACMIPatch(const vector &, const scalar, streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a.
void hitProcessorPatch(streamlinesCloud &, trackingData &)
const polyMesh & mesh() const
Return the mesh database.
Definition: particleI.H:125
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
bool move(streamlinesCloud &, trackingData &, const scalar)
Track all particles to their end point.
void hitCyclicAMIPatch(const vector &, const scalar, streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a.
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
void hitWallPatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a wallPatch.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
trackingData(streamlinesCloud &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, DynamicList< List< scalar >> &allTimes, List< DynamicList< scalarList >> &allScalars, List< DynamicList< vectorList >> &allVectors)
Construct from components.
const dimensionedScalar c
Speed of light in a vacuum.
Base particle class.
Definition: particle.H:83
static void writeFields(const Cloud< streamlinesParticle > &)
Write.
DynamicList< vectorList > & allPositions_
List< DynamicList< scalarList > > & allScalars_
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
A Cloud of streamlines particles.
bool hitPatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a patch.
void hitWedgePatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a wedge.
Particle class that samples fields as it passes through. Used in streamlines calculation.
const PtrList< interpolation< scalar > > & vsInterp_
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void hitCyclicPatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a cyclic.
const PtrList< interpolation< vector > > & vvInterp_
void hitSymmetryPatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a.
Factory class to read-construct particles used for parallel transfer.
List< DynamicList< vectorList > > & allVectors_
autoPtr< particle > clone() const
Construct and return a clone.
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
void hitCyclicRepeatAMIPatch(const vector &, const scalar, streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a.
Ostream & operator<<(Ostream &, const ensightPart &)
static void readFields(Cloud< streamlinesParticle > &)
Read.
friend Ostream & operator<<(Ostream &, const streamlinesParticle &)
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
void hitSymmetryPlanePatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a.
streamlinesParticle(const polyMesh &c, const vector &position, const label celli, const label lifeTime)
Construct from components.
vector position() const
Return current particle position.
Definition: particleI.H:278
Namespace for OpenFOAM.