streamLineBase.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2015 OpenFOAM Foundation
9  Copyright (C) 2016-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::functionObjects::streamLineBase
29 
30 SeeAlso
31  Foam::functionObjects::streamLine
32  Foam::functionObjects::wallBoundedStreamLine
33 
34 SourceFiles
35  streamLineBase.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef Foam_functionObjects_streamLineBase_H
40 #define Foam_functionObjects_streamLineBase_H
41 
42 #include "fvMeshFunctionObject.H"
43 #include "DynamicList.H"
44 #include "scalarList.H"
45 #include "vectorList.H"
46 #include "coordSetWriter.H"
47 #include "indirectPrimitivePatch.H"
48 #include "interpolation.H"
49 #include "Enum.H"
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 namespace Foam
54 {
55 
56 // Forward Declarations
57 class meshSearch;
58 class sampledSet;
59 
60 namespace functionObjects
61 {
62 
63 /*---------------------------------------------------------------------------*\
64  Class streamLineBase Declaration
65 \*---------------------------------------------------------------------------*/
66 
67 class streamLineBase
68 :
70 {
71 public:
72 
73  // Public data types
74 
75  //- Enumeration defining the track direction
76  enum trackDirType : char
77  {
78  FORWARD,
79  BACKWARD,
81  };
82 
83  //- Names for the trackDir
85 
86 
87 protected:
88 
89  //- Seed set engine
91 
92  //- Axis of the sampled points to output
93  mutable word sampledSetAxis_;
94 
95  //- Input dictionary
97 
98  //- List of fields to sample
101  //- Field to transport particle with
102  word UName_;
103 
104  //- Interpolation scheme to use
106 
107  //- Whether to use +U, -U or both
109 
110  //- Maximum lifetime (= number of cells) of particle
111  label lifeTime_;
112 
113  //- Track length
114  scalar trackLength_;
116  //- Optional trimming of tracks
118 
119  //- Optional specified name of particles
121 
122  //- Type of seed
123  word seedSet_;
124 
125  //- Names of scalar fields
127 
128  //- Names of vector fields
131 
132  // Demand Driven
133 
134  //- File writer for tracks data
136 
137 
138  // Generated Data
139 
140  //- All tracks. Per track the points it passed through
142 
143  //- Per scalarField, per track, the sampled values
146  //- Per vectorField, per track, the sampled values
148 
149 
150  // Protected Member Functions
151 
152  //- The axis of the sampledSet. Creates sampledSet if required.
153  const word& sampledSetAxis() const;
154 
155  //- Demand driven construction of the sampledSet.
156  // Also updates sampledSetAxis_
157  const sampledSet& sampledSetPoints() const;
158 
159  //- Construct patch out of all wall patch faces
161 
162  //- Initialise interpolators and track storage
163  // Return velocity interpolator: standalone or part of vector
164  // interpolators
166  (
167  const label nSeeds,
168  PtrList<interpolation<scalar>>& vsInterp,
169  PtrList<interpolation<vector>>& vvInterp
170  );
172  //- Generate point and values by interpolating from existing values
173  void storePoint
174  (
175  const label tracki,
177  const scalar w,
178  const label lefti,
179  const label righti,
180 
182  DynamicList<List<scalar>>& newScalars,
183  DynamicList<List<vector>>& newVectors
184  ) const;
185 
186  //- Trim and possibly split a track
187  void trimToBox
188  (
189  const treeBoundBox& bb,
190  const label tracki,
191  PtrList<DynamicList<point>>& newTracks,
192  PtrList<DynamicList<scalarList>>& newScalars,
193  PtrList<DynamicList<vectorList>>& newVectors
194  ) const;
195 
196  //- Trim tracks to bounding box
197  void trimToBox(const treeBoundBox& bb);
198 
199  //- Do the actual tracking to fill the track data
200  virtual void track() = 0;
201 
202  //- Write tracks to file
203  virtual bool writeToFile();
204 
205  //- Reset the field names
206  virtual void resetFieldNames
207  (
208  const word& newUName,
209  const wordList& newFieldNames
210  );
211 
212 
213 public:
214 
215  //- Runtime type information
216  TypeName("streamLineBase");
217 
218 
219  // Constructors
220 
221  //- Construct for given objectRegistry and dictionary.
222  // Allow the possibility to load fields from files
224  (
225  const word& name,
226  const Time& runTime,
227  const dictionary& dict
228  );
229 
230  //- Construct from Time and dictionary and list of fields to sample
232  (
233  const word& name,
234  const Time& runTime,
235  const dictionary& dict,
236  const wordList& fieldNames
237  );
238 
239 
240  //- Destructor
241  virtual ~streamLineBase();
242 
243 
244  // Member Functions
245 
246  //- Read the field average data
247  virtual bool read(const dictionary&);
248 
249  //- Execute the averaging
250  virtual bool execute();
251 
252  //- Track and write
253  virtual bool write();
254 
255  //- Update for changes of mesh
256  virtual void updateMesh(const mapPolyMesh&);
257 
258  //- Update for mesh point-motion
259  virtual void movePoints(const polyMesh&);
260 };
261 
262 
263 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
264 
265 } // End namespace functionObjects
266 } // End namespace Foam
267 
268 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
269 
270 #endif
271 
272 // ************************************************************************* //
word UName_
Field to transport particle with.
streamLineBase(const word &name, const Time &runTime, const dictionary &dict)
Construct for given objectRegistry and dictionary.
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
dictionary dict
const sampledSet & sampledSetPoints() const
Demand driven construction of the sampledSet.
word sampledSetAxis_
Axis of the sampled points to output.
virtual bool writeToFile()
Write tracks to file.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
DynamicList< List< point > > allTracks_
All tracks. Per track the points it passed through.
trackDirType
Enumeration defining the track direction.
engineTime & runTime
virtual void track()=0
Do the actual tracking to fill the track data.
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
wordList scalarNames_
Names of scalar fields.
wordList fields_
List of fields to sample.
refPtr< interpolation< vector > > initInterpolations(const label nSeeds, PtrList< interpolation< scalar >> &vsInterp, PtrList< interpolation< vector >> &vvInterp)
Initialise interpolators and track storage.
List< DynamicList< scalarList > > allScalars_
Per scalarField, per track, the sampled values.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
A class for managing references or pointers (no reference counting)
Definition: HashPtrTable.H:49
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:157
const word & name() const noexcept
Return the name of this functionObject.
wordList vectorNames_
Names of vector fields.
List< DynamicList< vectorList > > allVectors_
Per vectorField, per track, the sampled values.
boundBox bounds_
Optional trimming of tracks.
trackDirType trackDir_
Whether to use +U, -U or both.
virtual bool write()
Track and write.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:51
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:79
virtual bool read(const dictionary &)
Read the field average data.
A class for handling words, derived from Foam::string.
Definition: word.H:63
autoPtr< coordSetWriter > trackWriterPtr_
File writer for tracks data.
void storePoint(const label tracki, const scalar w, const label lefti, const label righti, DynamicList< point > &newTrack, DynamicList< List< scalar >> &newScalars, DynamicList< List< vector >> &newVectors) const
Generate point and values by interpolating from existing values.
word cloudName_
Optional specified name of particles.
static const Enum< trackDirType > trackDirTypeNames
Names for the trackDir.
dictionary dict_
Input dictionary.
void trimToBox(const treeBoundBox &bb, const label tracki, PtrList< DynamicList< point >> &newTracks, PtrList< DynamicList< scalarList >> &newScalars, PtrList< DynamicList< vectorList >> &newVectors) const
Trim and possibly split a track.
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
virtual bool execute()
Execute the averaging.
word interpolationScheme_
Interpolation scheme to use.
label lifeTime_
Maximum lifetime (= number of cells) of particle.
TypeName("streamLineBase")
Runtime type information.
autoPtr< sampledSet > sampledSetPtr_
Seed set engine.
virtual void resetFieldNames(const word &newUName, const wordList &newFieldNames)
Reset the field names.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:90
const word & sampledSetAxis() const
The axis of the sampledSet. Creates sampledSet if required.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
autoPtr< indirectPrimitivePatch > wallPatch() const
Construct patch out of all wall patch faces.
Namespace for OpenFOAM.