streamLine.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-2017 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::functionObjects::streamLine
26 
27 Group
28  grpFieldFunctionObjects
29 
30 Description
31  Generates streamline data by sampling a set of user-specified fields along a
32  particle track, transported by a user-specified velocity field.
33 
34  Example of function object specification:
35  \verbatim
36  streamLine1
37  {
38  type streamLine;
39  libs ("libfieldFunctionObjects.so");
40 
41  writeControl writeTime;
42 
43  setFormat vtk;
44  U U;
45  trackForward yes;
46 
47  fields
48  (
49  U
50  p
51  );
52 
53  lifeTime 10000;
54  trackLength 1e-3;
55  nSubCycle 5;
56  cloudName particleTracks;
57 
58  seedSampleSet
59  {
60  type uniform;
61  axis x; //distance;
62  start (-0.0205 0.0001 0.00001);
63  end (-0.0205 0.0005 0.00001);
64  nPoints 100;
65  }
66  }
67  \endverbatim
68 
69 Usage
70  \table
71  Property | Description | Required | Default value
72  type | Type name: streamLine | yes |
73  setFormat | Output data type | yes |
74  U | Tracking velocity field name | no | U
75  fields | Fields to sample | yes |
76  lifetime | Maximum number of particle tracking steps | yes |
77  trackLength | Tracking segment length | no |
78  nSubCycle | Number of tracking steps per cell | no|
79  cloudName | Cloud name to use | yes |
80  seedSampleSet| Seeding method (see below)| yes |
81  \endtable
82 
83  Where \c seedSampleSet \c type is typically one of
84  \plaintable
85  uniform | uniform particle seeding
86  cloud | cloud of points
87  triSurfaceMeshPointSet | points according to a tri-surface mesh
88  \endplaintable
89 
90 Note
91  When specifying the track resolution, the \c trackLength OR \c nSubCycle
92  option should be used
93 
94 See also
95  Foam::functionObject
96  Foam::functionObjects::timeControl
97  Foam::sampledSet
98  Foam::wallBoundedStreamLine
99 
100 SourceFiles
101  streamLine.C
102 
103 \*---------------------------------------------------------------------------*/
104 
105 #ifndef functionObjects_streamLine_H
106 #define functionObjects_streamLine_H
107 
108 #include "fvMeshFunctionObject.H"
109 #include "volFieldsFwd.H"
110 #include "DynamicList.H"
111 #include "scalarList.H"
112 #include "vectorList.H"
113 #include "writer.H"
114 #include "indirectPrimitivePatch.H"
115 
116 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
117 
118 namespace Foam
119 {
120 
121 // Forward declaration of classes
122 class meshSearch;
123 class sampledSet;
124 
125 namespace functionObjects
126 {
127 
128 /*---------------------------------------------------------------------------*\
129  Class streamLine Declaration
130 \*---------------------------------------------------------------------------*/
131 
132 class streamLine
133 :
134  public fvMeshFunctionObject
135 {
136  // Private data
137 
138  //- Input dictionary
139  dictionary dict_;
140 
141  //- List of fields to sample
142  wordList fields_;
143 
144  //- Field to transport particle with
145  word UName_;
146 
147  //- Interpolation scheme to use
148  word interpolationScheme_;
149 
150  //- Whether to use +u or -u
151  bool trackForward_;
152 
153  //- Maximum lifetime (= number of cells) of particle
154  label lifeTime_;
155 
156  //- Number of subcycling steps
157  label nSubCycle_;
158 
159  //- Track length
160  scalar trackLength_;
161 
162  //- Optional specified name of particles
163  word cloudName_;
164 
165  //- Type of seed
166  word seedSet_;
167 
168  //- Names of scalar fields
169  wordList scalarNames_;
170 
171  //- Names of vector fields
172  wordList vectorNames_;
173 
174 
175  // Demand driven
176 
177  //- Mesh searching enigne
178  autoPtr<meshSearch> meshSearchPtr_;
179 
180  //- Seed set engine
181  autoPtr<sampledSet> sampledSetPtr_;
182 
183  //- Axis of the sampled points to output
184  word sampledSetAxis_;
185 
186  //- File writer for scalar data
187  autoPtr<writer<scalar>> scalarFormatterPtr_;
188 
189  //- File writer for vector data
190  autoPtr<writer<vector>> vectorFormatterPtr_;
191 
192 
193  // Generated data
194 
195  //- All tracks. Per particle the points it passed through
196  DynamicList<List<point>> allTracks_;
197 
198  //- Per scalarField, per particle, the sampled value.
199  List<DynamicList<scalarList>> allScalars_;
200 
201  //- Per scalarField, per particle, the sampled value.
202  List<DynamicList<vectorList>> allVectors_;
203 
204 
205  //- Construct patch out of all wall patch faces
206  autoPtr<indirectPrimitivePatch> wallPatch() const;
207 
208  //- Do all seeding and tracking
209  void track();
210 
211  //- Disallow default bitwise copy construct
212  streamLine(const streamLine&);
213 
214  //- Disallow default bitwise assignment
215  void operator=(const streamLine&);
216 
217 
218 public:
219 
220  //- Runtime type information
221  TypeName("streamLine");
222 
223 
224  // Constructors
225 
226  //- Construct from Time and dictionary
227  streamLine
228  (
229  const word& name,
230  const Time& runTime,
231  const dictionary& dict
232  );
233 
234 
235  //- Destructor
236  virtual ~streamLine();
237 
238 
239  // Member Functions
240 
241  //- Read the field average data
242  virtual bool read(const dictionary&);
243 
244  //- Do nothing
245  virtual bool execute();
246 
247  //- Calculate and write the steamlines
248  virtual bool write();
249 
250  //- Update for changes of mesh
251  virtual void updateMesh(const mapPolyMesh&);
252 
253  //- Update for mesh point-motion
254  virtual void movePoints(const polyMesh&);
255 };
256 
257 
258 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
259 
260 } // End namespace functionObjects
261 } // End namespace Foam
262 
263 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
264 
265 #endif
266 
267 // ************************************************************************* //
dictionary dict
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
const word & name() const
Return the name of this functionObject.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:60
virtual ~streamLine()
Destructor.
Definition: streamLine.C:298
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
virtual bool read(const dictionary &)
Read the field average data.
Definition: streamLine.C:304
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
A class for handling words, derived from string.
Definition: word.H:59
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
Definition: streamLine.C:641
List< word > wordList
A List of words.
Definition: fileName.H:54
TypeName("streamLine")
Runtime type information.
Generates streamline data by sampling a set of user-specified fields along a particle track...
Definition: streamLine.H:193
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
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
Definition: streamLine.C:650
virtual bool execute()
Do nothing.
Definition: streamLine.C:394
virtual bool write()
Calculate and write the steamlines.
Definition: streamLine.C:400
Namespace for OpenFOAM.