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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::streamLineParticle
28 
29 Description
30  Particle class that samples fields as it passes through. Used in streamline
31  calculation.
32 
33 SourceFiles
34  streamLineParticle.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef streamLineParticle_H
39 #define streamLineParticle_H
40 
41 #include "particle.H"
42 #include "autoPtr.H"
43 #include "interpolation.H"
44 #include "vectorList.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 // Forward declaration of friend functions and operators
52 
53 class streamLineParticle;
54 class streamLineParticleCloud;
55 
56 Ostream& operator<<(Ostream&, const streamLineParticle&);
57 
58 /*---------------------------------------------------------------------------*\
59  Class streamLineParticle Declaration
60 \*---------------------------------------------------------------------------*/
61 
63 :
64  public particle
65 {
66 public:
67 
68  class trackingData
69  :
71  {
72  public:
73 
74  // Public data
75 
77 
79 
80  const label UIndex_;
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 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  nSubCycle_(nSubCycle),
114  trackLength_(trackLength),
115  allPositions_(allPositions),
116  allScalars_(allScalars),
117  allVectors_(allVectors)
118  {}
119  };
120 
121 
122 private:
123 
124  // Private data
125 
126  //- Whether particle transports with +U or -U
127  bool trackForward_;
128 
129  //- Lifetime of particle. Particle dies when reaches 0.
130  label lifeTime_;
131 
132  //- Sampled positions
133  DynamicList<point> sampledPositions_;
134 
135  //- Sampled scalars
136  List<DynamicList<scalar>> sampledScalars_;
137 
138  //- Sampled vectors
139  List<DynamicList<vector>> sampledVectors_;
140 
141 
142  // Private Member Functions
143 
144  //- Interpolate all quantities; return interpolated velocity.
145  vector interpolateFields
146  (
147  const trackingData&,
148  const point&,
149  const label celli,
150  const label facei
151  );
152 
153 
154 public:
155 
156  // Constructors
157 
158  //- Construct from components
160  (
161  const polyMesh& c,
162  const vector& position,
163  const label celli,
164  const bool trackForward,
165  const label lifeTime
166  );
167 
168  //- Construct from Istream
170  (
171  const polyMesh& c,
172  Istream& is,
173  bool readFields = true,
174  bool newFormat = true
175  );
176 
177  //- Construct copy
179 
180  //- Construct and return a clone
181  autoPtr<particle> clone() const
182  {
183  return autoPtr<particle>(new streamLineParticle(*this));
184  }
185 
186  //- Factory class to read-construct particles used for parallel transfer
187  class iNew
188  {
189  const polyMesh& mesh_;
190 
191  public:
192 
193  iNew(const polyMesh& mesh)
194  :
195  mesh_(mesh)
196  {}
197 
199  {
201  (
202  new streamLineParticle(mesh_, is, true)
203  );
204  }
205  };
206 
207 
208  // Member Functions
209 
210  // Tracking
211 
212  //- Track all particles to their end point
213  bool move(streamLineParticleCloud&, trackingData&, const scalar);
214 
215  //- Overridable function to handle the particle hitting a patch
216  // Executed before other patch-hitting functions
217  bool hitPatch(streamLineParticleCloud&, trackingData&);
218 
219  //- Overridable function to handle the particle hitting a wedge
220  void hitWedgePatch(streamLineParticleCloud&, trackingData&);
221 
222  //- Overridable function to handle the particle hitting a
223  // symmetry plane
224  void hitSymmetryPlanePatch(streamLineParticleCloud&, trackingData&);
225 
226  //- Overridable function to handle the particle hitting a
227  // symmetry patch
228  void hitSymmetryPatch(streamLineParticleCloud&, trackingData&);
229 
230  //- Overridable function to handle the particle hitting a cyclic
231  void hitCyclicPatch(streamLineParticleCloud&, trackingData&);
232 
233  //- Overridable function to handle the particle hitting a
234  // cyclicAMIPatch
235  void hitCyclicAMIPatch
236  (
238  trackingData&,
239  const vector& direction
240  );
241 
242  //- Overridable function to handle the particle hitting a
243  // cyclicACMIPatch
244  void hitCyclicACMIPatch
245  (
247  trackingData&,
248  const vector& direction
249  );
250 
251  //- Overridable function to handle the particle hitting a
252  //- processorPatch
253  void hitProcessorPatch(streamLineParticleCloud&, trackingData&);
254 
255  //- Overridable function to handle the particle hitting a wallPatch
256  void hitWallPatch(streamLineParticleCloud&, trackingData&);
257 
258 
259  // I-O
260 
261  //- Read
263 
264  //- Write
265  static void writeFields(const Cloud<streamLineParticle>&);
266 
267 
268  // Ostream Operator
269 
270  friend Ostream& operator<<(Ostream&, const streamLineParticle&);
271 };
272 
273 
274 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
275 
276 } // End namespace Foam
277 
278 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
279 
280 #endif
281 
282 // ************************************************************************* //
Foam::streamLineParticle::streamLineParticle
streamLineParticle(const polyMesh &c, const vector &position, const label celli, const bool trackForward, const label lifeTime)
Construct from components.
Definition: streamLineParticle.C:86
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::streamLineParticle::hitSymmetryPatch
void hitSymmetryPatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a.
Definition: streamLineParticle.C:307
Foam::streamLineParticle::hitProcessorPatch
void hitProcessorPatch(streamLineParticleCloud &, trackingData &)
Definition: streamLineParticle.C:353
Foam::streamLineParticle::trackingData
Definition: streamLineParticle.H:67
vectorList.H
Foam::streamLineParticleCloud
A Cloud of streamLine particles.
Definition: streamLineParticleCloud.H:52
Foam::streamLineParticle
Particle class that samples fields as it passes through. Used in streamline calculation.
Definition: streamLineParticle.H:61
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:55
Foam::streamLineParticle::writeFields
static void writeFields(const Cloud< streamLineParticle > &)
Write.
Definition: streamLineParticle.C:404
interpolation.H
Foam::streamLineParticle::trackingData::UIndex_
const label UIndex_
Definition: streamLineParticle.H:79
Foam::streamLineParticle::hitCyclicACMIPatch
void hitCyclicACMIPatch(streamLineParticleCloud &, trackingData &, const vector &direction)
Overridable function to handle the particle hitting a.
Definition: streamLineParticle.C:341
Foam::streamLineParticle::hitWedgePatch
void hitWedgePatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a wedge.
Definition: streamLineParticle.C:285
Foam::streamLineParticle::trackingData::allScalars_
List< DynamicList< scalarList > > & allScalars_
Definition: streamLineParticle.H:87
Foam::streamLineParticle::trackingData::nSubCycle_
const label nSubCycle_
Definition: streamLineParticle.H:81
Foam::streamLineParticle::trackingData::allPositions_
DynamicList< vectorList > & allPositions_
Definition: streamLineParticle.H:85
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::operator<<
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
Foam::streamLineParticle::operator<<
friend Ostream & operator<<(Ostream &, const streamLineParticle &)
Foam::streamLineParticle::trackingData::trackLength_
const scalar trackLength_
Definition: streamLineParticle.H:83
Foam::streamLineParticle::trackingData::trackingData
trackingData(streamLineParticleCloud &cloud, const PtrList< interpolation< scalar >> &vsInterp, const PtrList< interpolation< vector >> &vvInterp, const label UIndex, const label nSubCycle, const scalar trackLength, DynamicList< List< point >> &allPositions, List< DynamicList< scalarList >> &allScalars, List< DynamicList< vectorList >> &allVectors)
Construct from components.
Definition: streamLineParticle.H:96
Foam::streamLineParticle::move
bool move(streamLineParticleCloud &, trackingData &, const scalar)
Track all particles to their end point.
Definition: streamLineParticle.C:152
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::streamLineParticle::clone
autoPtr< particle > clone() const
Construct and return a clone.
Definition: streamLineParticle.H:180
Foam::streamLineParticle::iNew::iNew
iNew(const polyMesh &mesh)
Definition: streamLineParticle.H:192
Foam::particle::position
vector position() const
Return current particle position.
Definition: particleI.H:314
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
Foam::particle::mesh
const polyMesh & mesh() const
Return the mesh database.
Definition: particleI.H:137
Foam::streamLineParticle::trackingData::allVectors_
List< DynamicList< vectorList > > & allVectors_
Definition: streamLineParticle.H:89
Foam::streamLineParticle::hitWallPatch
void hitWallPatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a wallPatch.
Definition: streamLineParticle.C:364
Foam::streamLineParticle::trackingData::vsInterp_
const PtrList< interpolation< scalar > > & vsInterp_
Definition: streamLineParticle.H:75
Foam::interpolation< scalar >
Foam::streamLineParticle::hitSymmetryPlanePatch
void hitSymmetryPlanePatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a.
Definition: streamLineParticle.C:296
Foam::streamLineParticle::iNew
Factory class to read-construct particles used for parallel transfer.
Definition: streamLineParticle.H:186
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::streamLineParticle::hitCyclicAMIPatch
void hitCyclicAMIPatch(streamLineParticleCloud &, trackingData &, const vector &direction)
Overridable function to handle the particle hitting a.
Definition: streamLineParticle.C:329
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::cloud
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:57
Foam::streamLineParticle::hitPatch
bool hitPatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a patch.
Definition: streamLineParticle.C:277
Foam::Vector< scalar >
Foam::streamLineParticle::iNew::operator()
autoPtr< streamLineParticle > operator()(Istream &is) const
Definition: streamLineParticle.H:197
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::particle
Base particle class.
Definition: particle.H:76
Foam::streamLineParticle::hitCyclicPatch
void hitCyclicPatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a cyclic.
Definition: streamLineParticle.C:318
Foam::Cloud< streamLineParticle >
Foam::direction
uint8_t direction
Definition: direction.H:52
Foam::streamLineParticle::trackingData::vvInterp_
const PtrList< interpolation< vector > > & vvInterp_
Definition: streamLineParticle.H:77
particle.H
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::particle::trackingData
Definition: particle.H:95
autoPtr.H
Foam::streamLineParticle::readFields
static void readFields(Cloud< streamLineParticle > &)
Read.
Definition: streamLineParticle.C:374