wallBoundedStreamLineParticle.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::wallBoundedStreamLineParticle
26 
27 Description
28  Particle class that samples fields as it passes through. Used in streamline
29  calculation.
30 
31 SourceFiles
32  wallBoundedStreamLineParticle.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef wallBoundedStreamLineParticle_H
37 #define wallBoundedStreamLineParticle_H
38 
39 #include "wallBoundedParticle.H"
40 #include "autoPtr.H"
41 #include "interpolation.H"
42 #include "vectorList.H"
43 #include "InfoProxy.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 class wallBoundedStreamLineParticleCloud;
51 
52 /*---------------------------------------------------------------------------*\
53  Class wallBoundedStreamLineParticle Declaration
54 \*---------------------------------------------------------------------------*/
55 
57 :
58  public wallBoundedParticle
59 {
60 
61 public:
62 
63  //- Class used to pass tracking data to the trackToEdge function
64  class trackingData
65  :
67  <
68  Cloud<wallBoundedStreamLineParticle>
69  >
70  {
71 
72  public:
73 
74 
77  const label UIndex_;
78  const bool trackForward_;
79  const scalar trackLength_;
80 
84 
85 
86  // Constructors
87 
89  (
91  const PtrList<interpolation<scalar> >& vsInterp,
92  const PtrList<interpolation<vector> >& vvInterp,
93  const label UIndex,
94  const bool trackForward,
95  const scalar trackLength,
96  const PackedBoolList& isWallPatch,
97 
98  DynamicList<List<point> >& allPositions,
99  List<DynamicList<scalarList> >& allScalars,
100  List<DynamicList<vectorList> >& allVectors
101  )
102  :
104  <
106  >
107  (
108  cloud,
109  isWallPatch
110  ),
111  vsInterp_(vsInterp),
112  vvInterp_(vvInterp),
113  UIndex_(UIndex),
114  trackForward_(trackForward),
115  trackLength_(trackLength),
116 
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.
130 
131  //- Sampled positions
133 
134  //- Sampled scalars
136 
137  //- Sampled vectors
139 
140 
141  // Private Member Functions
142 
144  (
145  const trackingData& td,
146  const point& position,
147  const label cellI,
148  const label faceI
149  );
150 
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 label tetFaceI,
165  const label tetPtI,
166  const label meshEdgeStart,
167  const label diagEdge,
168  const label lifeTime
169  );
170 
171  //- Construct from Istream
173  (
174  const polyMesh& c,
175  Istream& is,
176  bool readFields = true
177  );
178 
179  //- Construct copy
181 
182  //- Construct and return a clone
183  autoPtr<particle> clone() const
184  {
186  }
187 
188  //- Factory class to read-construct particles used for
189  // parallel transfer
190  class iNew
191  {
192  const polyMesh& mesh_;
193 
194  public:
195 
196  iNew(const polyMesh& mesh)
197  :
198  mesh_(mesh)
199  {}
200 
202  (
203  Istream& is
204  ) const
205  {
207  (
208  new wallBoundedStreamLineParticle(mesh_, is, true)
209  );
210  }
211  };
212 
213 
214  // Member Functions
215 
216  // Tracking
217 
218  //- Track all particles to their end point
219  bool move(trackingData&, const scalar trackTime);
220 
221 
222  // I-O
223 
224  //- Read
226 
227  //- Write
228  static void writeFields
229  (
231  );
232 
233 
234  // Ostream Operator
235 
236  friend Ostream& operator<<
237  (
238  Ostream&,
240  );
241 };
242 
243 
244 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
245 
246 } // End namespace Foam
247 
248 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
249 
250 #endif
251 
252 // ************************************************************************* //
Foam::wallBoundedStreamLineParticle::clone
autoPtr< particle > clone() const
Construct and return a clone.
Definition: wallBoundedStreamLineParticle.H:182
Foam::PackedBoolList
A bit-packed bool list.
Definition: PackedBoolList.H:63
vectorList.H
p
p
Definition: pEqn.H:62
Foam::wallBoundedStreamLineParticle::move
bool move(trackingData &, const scalar trackTime)
Track all particles to their end point.
Definition: wallBoundedStreamLineParticle.C:213
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:56
InfoProxy.H
Foam::wallBoundedStreamLineParticle::sampledPositions_
DynamicList< point > sampledPositions_
Sampled positions.
Definition: wallBoundedStreamLineParticle.H:131
Foam::wallBoundedStreamLineParticle::trackingData::allScalars_
List< DynamicList< scalarList > > & allScalars_
Definition: wallBoundedStreamLineParticle.H:81
Foam::particle::TrackingData< Cloud< wallBoundedStreamLineParticle > >::cloud
Cloud< wallBoundedStreamLineParticle > & cloud()
Return a reference to the cloud.
Definition: particle.H:125
interpolation.H
Foam::wallBoundedStreamLineParticle::iNew::iNew
iNew(const polyMesh &mesh)
Definition: wallBoundedStreamLineParticle.H:195
Foam::wallBoundedStreamLineParticle::trackingData::trackingData
trackingData(Cloud< wallBoundedStreamLineParticle > &cloud, const PtrList< interpolation< scalar > > &vsInterp, const PtrList< interpolation< vector > > &vvInterp, const label UIndex, const bool trackForward, const scalar trackLength, const PackedBoolList &isWallPatch, DynamicList< List< point > > &allPositions, List< DynamicList< scalarList > > &allScalars, List< DynamicList< vectorList > > &allVectors)
Definition: wallBoundedStreamLineParticle.H:88
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::wallBoundedStreamLineParticle::trackingData::UIndex_
const label UIndex_
Definition: wallBoundedStreamLineParticle.H:76
Foam::wallBoundedStreamLineParticle::lifeTime_
label lifeTime_
Lifetime of particle. Particle dies when reaches 0.
Definition: wallBoundedStreamLineParticle.H:128
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::wallBoundedParticle::TrackingData
Class used to pass tracking data to the trackToFace function.
Definition: wallBoundedParticle.H:66
Foam::wallBoundedStreamLineParticle::sample
vector sample(trackingData &td)
Definition: wallBoundedStreamLineParticle.C:106
Foam::wallBoundedStreamLineParticle::trackingData::allVectors_
List< DynamicList< vectorList > > & allVectors_
Definition: wallBoundedStreamLineParticle.H:82
Foam::wallBoundedStreamLineParticle::trackingData::allPositions_
DynamicList< vectorList > & allPositions_
Definition: wallBoundedStreamLineParticle.H:80
Foam::wallBoundedStreamLineParticle::trackingData
Class used to pass tracking data to the trackToEdge function.
Definition: wallBoundedStreamLineParticle.H:63
wallBoundedParticle.H
Foam::wallBoundedParticle::meshEdgeStart
label meshEdgeStart() const
-1 or label of mesh edge
Definition: wallBoundedParticle.H:272
Foam::wallBoundedStreamLineParticle::trackingData::trackForward_
const bool trackForward_
Definition: wallBoundedStreamLineParticle.H:77
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::wallBoundedParticle
Particle class that tracks on triangles of boundary faces. Use trackToEdge similar to trackToFace on ...
Definition: wallBoundedParticle.H:52
Foam::wallBoundedStreamLineParticle
Particle class that samples fields as it passes through. Used in streamline calculation.
Definition: wallBoundedStreamLineParticle.H:55
Foam::wallBoundedStreamLineParticle::trackingData::vvInterp_
const PtrList< interpolation< vector > > & vvInterp_
Definition: wallBoundedStreamLineParticle.H:75
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::wallBoundedStreamLineParticle::writeFields
static void writeFields(const Cloud< wallBoundedStreamLineParticle > &)
Write.
Definition: wallBoundedStreamLineParticle.C:379
Foam::interpolation< scalar >
Foam::wallBoundedStreamLineParticle::interpolateFields
vector interpolateFields(const trackingData &td, const point &position, const label cellI, const label faceI)
Definition: wallBoundedStreamLineParticle.C:32
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::wallBoundedStreamLineParticle::sampledVectors_
List< DynamicList< vector > > sampledVectors_
Sampled vectors.
Definition: wallBoundedStreamLineParticle.H:137
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::wallBoundedStreamLineParticle::readFields
static void readFields(Cloud< wallBoundedStreamLineParticle > &)
Read.
Definition: wallBoundedStreamLineParticle.C:345
Foam::Vector< scalar >
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::wallBoundedStreamLineParticle::iNew
Factory class to read-construct particles used for.
Definition: wallBoundedStreamLineParticle.H:189
Foam::Cloud< wallBoundedStreamLineParticle >
Foam::wallBoundedStreamLineParticle::trackingData::trackLength_
const scalar trackLength_
Definition: wallBoundedStreamLineParticle.H:78
Foam::wallBoundedStreamLineParticle::wallBoundedStreamLineParticle
wallBoundedStreamLineParticle(const polyMesh &c, const vector &position, const label cellI, const label tetFaceI, const label tetPtI, const label meshEdgeStart, const label diagEdge, const label lifeTime)
Construct from components.
Definition: wallBoundedStreamLineParticle.C:134
Foam::wallBoundedStreamLineParticle::trackingData::vsInterp_
const PtrList< interpolation< scalar > > & vsInterp_
Definition: wallBoundedStreamLineParticle.H:74
Foam::particle::position
const vector & position() const
Return current particle position.
Definition: particleI.H:586
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::wallBoundedParticle::diagEdge
label diagEdge() const
-1 or diagonal edge
Definition: wallBoundedParticle.H:278
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::wallBoundedStreamLineParticle::iNew::mesh_
const polyMesh & mesh_
Definition: wallBoundedStreamLineParticle.H:191
Foam::wallBoundedStreamLineParticle::sampledScalars_
List< DynamicList< scalar > > sampledScalars_
Sampled scalars.
Definition: wallBoundedStreamLineParticle.H:134
Foam::particle::TrackingData
Definition: particle.H:94
autoPtr.H