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