streamlinesParticle.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-2025 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::streamlinesParticle
26 
27 Description
28  Particle class that samples fields as it passes through. Used in streamlines
29  calculation.
30 
31 SourceFiles
32  streamlinesParticle.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef streamlinesParticle_H
37 #define streamlinesParticle_H
38 
39 #include "particle.H"
40 #include "Cloud.H"
41 #include "autoPtr.H"
42 #include "interpolation.H"
43 #include "vectorList.H"
44 #include "DynamicField.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 // Forward declaration of friend functions and operators
52 
53 class streamlinesParticle;
54 class streamlinesCloud;
55 
56 Ostream& operator<<(Ostream&, const streamlinesParticle&);
57 
58 /*---------------------------------------------------------------------------*\
59  Class streamlinesParticle Declaration
60 \*---------------------------------------------------------------------------*/
61 
63 :
64  public particle
65 {
66 public:
67 
68  class trackingData
69  :
71  {
72  public:
73 
74  // Public data
75 
76  #define DeclareTypeInterpolator(Type, nullArg) \
77  const PtrList<interpolation<Type>>& Type##Interp_;
79  #undef DeclareTypeInterpolator
80 
82 
83  bool trackForward_;
84 
85  bool trackOutside_;
86 
87  const label nSubCycle_;
88 
89  const scalar trackLength_;
90 
91  label tracki;
92 
94 
96 
98 
100 
101  #define DeclareAllTypes(Type, nullArg) \
102  List<DynamicField<Type>>& all##Type##s_;
104  #undef DeclareAllTypes
105 
106 
107  // Constructors
108 
109  //- Construct from components
111  (
113  #define TypeInterpolatorArg(Type, nullArg) \
114  const PtrList<interpolation<Type>>& Type##Interp,
116  #undef TypeInterpolatorArg
117  const interpolation<vector>& UInterp,
118  const bool trackForward,
119  const bool trackOutside,
120  const label nSubCycle,
121  const scalar trackLength,
122  DynamicField<point>& allPositions,
123  DynamicField<label>& allTracks,
124  DynamicField<label>& allTrackParts,
125  DynamicField<scalar>& allAges
126  #define AllTypesArg(Type, nullArg) \
127  , List<DynamicField<Type>>& all##Type##s
129  #undef AllTypesArg
130  )
131  :
133  #define TypeInterpolatorInit(Type, nullArg) \
134  Type##Interp_(Type##Interp),
136  #undef TypeInterpolatorInit
137  UInterp_(UInterp),
138  trackForward_(trackForward),
139  trackOutside_(trackOutside),
140  nSubCycle_(nSubCycle),
141  trackLength_(trackLength),
142  allPositions_(allPositions),
143  allTracks_(allTracks),
144  allTrackParts_(allTrackParts),
145  allAges_(allAges)
146  #define AllTypesInit(Type, nullArg) \
147  , all##Type##s_(all##Type##s)
149  #undef AllTypesInit
150  {}
151  };
152 
153 
154 private:
155 
156  // Private Data
157 
158  //- Lifetime of particle. Particle dies when reaches 0.
159  label lifeTime_;
160 
161  //- Index of the track
162  label trackIndex_;
163 
164  //- Index of the part of the track
165  label trackPartIndex_;
166 
167  //- Age of the particle
168  scalar age_;
169 
170  //- Current compound transform
171  transformer transform_;
172 
173  //- Sampled positions
174  DynamicField<point> sampledPositions_;
175 
176  //- Sampled ages
177  DynamicField<scalar> sampledAges_;
178 
179  //- Sampled types
180  #define DeclareSampledTypes(Type, nullArg) \
181  List<DynamicField<Type>> sampled##Type##s_;
182  FOR_ALL_FIELD_TYPES(DeclareSampledTypes);
183  #undef DeclareSampledTypes
184 
185 
186  // Private Member Functions
187 
188  //- Interpolate all quantities; return interpolated velocity.
189  vector interpolateFields
190  (
191  const trackingData&,
192  const point&,
193  const label celli,
194  const label facei
195  );
196 
197  //- End the current track
198  void endTrack(trackingData&);
199 
200 
201 public:
202 
203  // Static Data Members
204 
205  //- Streamlines are computed at an instant
206  static const bool instantaneous = true;
207 
208 
209  // Constructors
210 
211  //- Construct from components
213  (
214  const polyMesh& mesh,
215  const vector& position,
216  const label celli,
217  label& nLocateBoundaryHits,
218  const label lifeTime,
219  const label trackIndex
220  );
221 
222  //- Construct from Istream
223  streamlinesParticle(Istream& is, bool readFields = true);
224 
225  //- Construct copy
227 
228  //- Construct and return a clone
229  autoPtr<particle> clone() const
230  {
231  return autoPtr<particle>(new streamlinesParticle(*this));
232  }
233 
234  //- Construct from Istream and return
236  {
238  }
239 
240 
241  // Member Functions
242 
243  // Tracking
244 
245  //- Track all particles to their end point
246  bool move(streamlinesCloud&, trackingData&);
247 
248  //- Overridable function to handle the particle hitting a wedge
249  void hitWedgePatch(streamlinesCloud&, trackingData&);
250 
251  //- Overridable function to handle the particle hitting a
252  // symmetry plane
253  void hitSymmetryPlanePatch(streamlinesCloud&, trackingData&);
254 
255  //- Overridable function to handle the particle hitting a
256  // symmetry patch
257  void hitSymmetryPatch(streamlinesCloud&, trackingData&);
258 
259  //- Overridable function to handle the particle hitting a cyclic
260  void hitCyclicPatch(streamlinesCloud&, trackingData&);
261 
262  //- Overridable function to handle the particle hitting an
263  // nonConformalCyclicPolyPatch
265  (
266  const vector& displacement,
267  const scalar fraction,
268  const label patchi,
270  trackingData& td
271  );
272 
273  //- Overridable function to handle the particle hitting a
274  // processorPatch
275  void hitProcessorPatch(streamlinesCloud&, trackingData&);
276 
277  //- Overridable function to handle the particle hitting a wallPatch
278  void hitWallPatch(streamlinesCloud&, trackingData&);
279 
280 
281  // Transformations
282 
283  //- Transform the physical properties of the particle
284  // according to the given transformation tensor
285  virtual void transformProperties(const transformer&);
286 
287 
288  // I-O
289 
290  //- Read
291  static void readFields
292  (
294  );
295 
296  //- Write
297  static void writeFields
298  (
300  );
301 
302 
303  // Ostream Operator
304 
305  friend Ostream& operator<<(Ostream&, const streamlinesParticle&);
306 };
307 
308 
309 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
310 
311 } // End namespace Foam
312 
313 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
314 
315 #endif
316 
317 // ************************************************************************* //
Dynamically sized Field.
Definition: DynamicField.H:72
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Base class for clouds. Provides a basic evolution algorithm, models, and a database for caching deriv...
Definition: cloud.H:63
Abstract base class for interpolation.
Definition: interpolation.H:55
Base particle class.
Definition: particle.H:83
vector position(const polyMesh &mesh) const
Return current particle position.
Definition: particleI.H:159
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
A Cloud of streamlines particles.
const interpolation< vector > & UInterp_
FOR_ALL_FIELD_TYPES(DeclareTypeInterpolator)
trackingData(streamlinesCloud &cloud, #define TypeInterpolatorArg(Type, nullArg) const interpolation< vector > &UInterp, const bool trackForward, const bool trackOutside, const label nSubCycle, const scalar trackLength, DynamicField< point > &allPositions, DynamicField< label > &allTracks, DynamicField< label > &allTrackParts, DynamicField< scalar > &allAges #define AllTypesArg(Type, nullArg))
Construct from components.
Particle class that samples fields as it passes through. Used in streamlines calculation.
bool move(streamlinesCloud &, trackingData &)
Track all particles to their end point.
static const bool instantaneous
Streamlines are computed at an instant.
void hitProcessorPatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a.
bool hitNonConformalCyclicPatch(const vector &displacement, const scalar fraction, const label patchi, streamlinesCloud &cloud, trackingData &td)
Overridable function to handle the particle hitting an.
static void writeFields(const lagrangian::Cloud< streamlinesParticle > &)
Write.
void hitSymmetryPatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a.
static autoPtr< streamlinesParticle > New(Istream &is)
Construct from Istream and return.
void hitCyclicPatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a cyclic.
void hitWallPatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a wallPatch.
void hitSymmetryPlanePatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a.
streamlinesParticle(const polyMesh &mesh, const vector &position, const label celli, label &nLocateBoundaryHits, const label lifeTime, const label trackIndex)
Construct from components.
virtual void transformProperties(const transformer &)
Transform the physical properties of the particle.
static void readFields(lagrangian::Cloud< streamlinesParticle > &)
Read.
friend Ostream & operator<<(Ostream &, const streamlinesParticle &)
autoPtr< particle > clone() const
Construct and return a clone.
void hitWedgePatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a wedge.
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space.
Definition: transformer.H:84
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
label patchi
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Namespace for OpenFOAM.
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
vector point
Point is a vector.
Definition: point.H:41
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
Ostream & operator<<(Ostream &os, const fvConstraints &constraints)
volScalarField & p
#define AllTypesInit(Type, nullArg)
#define DeclareAllTypes(Type, nullArg)
#define DeclareTypeInterpolator(Type, nullArg)
#define TypeInterpolatorInit(Type, nullArg)
#define DeclareSampledTypes(Type, nullArg)
Sampled types.
#define AllTypesArg(Type, nullArg)
#define TypeInterpolatorArg(Type, nullArg)