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