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