Loading [MathJax]/extensions/tex2jax.js
VTK  9.4.20250323
vtkParticleTracerBase.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
16#ifndef vtkParticleTracerBase_h
17#define vtkParticleTracerBase_h
18
19#include "vtkDeprecation.h" // For VTK_DEPRECATED_IN_9_4_0
20#include "vtkFiltersFlowPathsModule.h" // For export macro
22#include "vtkSmartPointer.h" // For vtkSmartPointer
23#include "vtkTemporalAlgorithm.h" // For vtkTemporalAlgorithm
24#include "vtkWeakPointer.h" // For vtkWeakPointer
25
26#include <list> // STL Header
27#include <mutex> // STL Header
28#include <numeric> // STL Header
29#include <unordered_map> // STL Header
30#include <vector> // STL Header
31
32#ifndef __VTK_WRAP__
33#define vtkPolyDataAlgorithm vtkTemporalAlgorithm<vtkPolyDataAlgorithm>
34#endif
35
36VTK_ABI_NAMESPACE_BEGIN
39class vtkDataArray;
40class vtkDataSet;
41class vtkDoubleArray;
42class vtkFloatArray;
43class vtkGenericCell;
45class vtkIntArray;
49class vtkPointData;
50class vtkPoints;
51class vtkPolyData;
54VTK_ABI_NAMESPACE_END
55
57{
58VTK_ABI_NAMESPACE_BEGIN
60{
61 double x[4];
62};
63using Position = struct Position_t;
64
66{
67 // These are used during iteration
72 // These are computed scalars we might display
74 int TimeStepAge; // amount of time steps the particle has advanced
76 int InjectedStepId; // time step the particle was injected
78 // These are useful to track for debugging etc
80 float age;
81 // these are needed across time steps to compute vorticity
82 float rotation;
84 float time;
85 float speed;
86 // once the particle is added, PointId is valid and is the tuple location
87 // in ProtoPD.
89
90 double velocity[3];
91};
93
94typedef std::vector<ParticleInformation> ParticleVector;
95typedef ParticleVector::iterator ParticleIterator;
96typedef std::list<ParticleInformation> ParticleDataList;
97typedef ParticleDataList::iterator ParticleListIterator;
98struct ParticleTracerFunctor;
99VTK_ABI_NAMESPACE_END
100}
101
102VTK_ABI_NAMESPACE_BEGIN
103class VTKFILTERSFLOWPATHS_EXPORT vtkParticleTracerBase : public vtkPolyDataAlgorithm
104{
105public:
106 friend struct vtkParticleTracerBaseNamespace::ParticleTracerFunctor;
108 {
113 UNKNOWN
114 };
115
118
119#ifndef __VTK_WRAP__
120#undef vtkPolyDataAlgorithm
121#endif
123 void PrintSelf(ostream& os, vtkIndent indent) override;
125
126#if defined(__VTK_WRAP__) || defined(__WRAP_GCCXML)
128#endif
129
131
138
140
145 vtkGetMacro(ComputeVorticity, bool);
148
150
153 vtkGetMacro(TerminalSpeed, double);
154 void SetTerminalSpeed(double);
156
158
162 vtkGetMacro(RotationScale, double);
163 void SetRotationScale(double);
165
167
171 vtkSetMacro(IgnorePipelineTime, vtkTypeBool);
172 vtkGetMacro(IgnorePipelineTime, vtkTypeBool);
173 vtkBooleanMacro(IgnorePipelineTime, vtkTypeBool);
175
177
186 vtkGetMacro(ForceReinjectionEveryNSteps, int);
189
191
198 "Please edit the TIME_STEPS information key in vtkAlgorithm::RequestInformation() instead")
199 virtual void SetTerminationTime(double) {}
201 "Please edit the TIME_STEPS information key in vtkAlgorithm::RequestInformation() instead")
202 virtual double GetTerminationTime() { return std::numeric_limits<double>::quiet_NaN(); }
204 "Please edit the TIME_STEPS information key in vtkAlgorithm::RequestInformation() instead")
205 virtual void SetStartTime(double) {}
207 "Please edit the TIME_STEPS information key in vtkAlgorithm::RequestInformation() instead")
208 virtual double GetStartTime() { return std::numeric_limits<double>::quiet_NaN(); }
210
212
215 VTK_DEPRECATED_IN_9_4_0("Caching is now automated")
216 virtual void SetDisableResetCache(bool) {}
217 VTK_DEPRECATED_IN_9_4_0("Caching is now automated")
218 virtual bool GetDisableResetCache() { return false; }
219 VTK_DEPRECATED_IN_9_4_0("Caching is now automated")
220 virtual void DisableResetCacheOn() {}
221 VTK_DEPRECATED_IN_9_4_0("Caching is now automated")
222 virtual void DisableResetCacheOff() {}
224
226 vtkGetObjectMacro(Integrator, vtkInitialValueProblemSolver);
227
228 void SetIntegratorType(int type);
230
232
241 vtkSetMacro(StaticSeeds, vtkTypeBool);
242 vtkGetMacro(StaticSeeds, vtkTypeBool);
244
249 {
250 DIFFERENT = 0,
251 STATIC = 1,
252 LINEAR_TRANSFORMATION = 2,
253 SAME_TOPOLOGY = 3
254 };
255
257 /*
258 * Set/Get the type of variance of the mesh over time.
259 *
260 * DIFFERENT = 0,
261 * STATIC = 1,
262 * LINEAR_TRANSFORMATION = 2
263 * SAME_TOPOLOGY = 3
264 */
265 virtual void SetMeshOverTime(int meshOverTime);
266 virtual int GetMeshOverTimeMinValue() { return DIFFERENT; }
267 virtual int GetMeshOverTimeMaxValue() { return SAME_TOPOLOGY; }
268 void SetMeshOverTimeToDifferent() { this->SetMeshOverTime(DIFFERENT); }
269 void SetMeshOverTimeToStatic() { this->SetMeshOverTime(STATIC); }
270 void SetMeshOverTimeToLinearTransformation() { this->SetMeshOverTime(LINEAR_TRANSFORMATION); }
271 void SetMeshOverTimeToSameTopology() { this->SetMeshOverTime(SAME_TOPOLOGY); }
272 vtkGetMacro(MeshOverTime, int);
274
275 enum
276 {
278 INTERPOLATOR_WITH_CELL_LOCATOR
279 };
280
292 void SetInterpolatorType(int interpolatorType);
293
302
310
312
319 vtkGetObjectMacro(ParticleWriter, vtkAbstractParticleWriter);
321
323
327 vtkSetFilePathMacro(ParticleFileName);
328 vtkGetFilePathMacro(ParticleFileName);
330
332
336 vtkSetMacro(EnableParticleWriting, vtkTypeBool);
337 vtkGetMacro(EnableParticleWriting, vtkTypeBool);
338 vtkBooleanMacro(EnableParticleWriting, vtkTypeBool);
340
342
348
350
354 vtkGetMacro(ForceSerialExecution, bool);
355 vtkSetMacro(ForceSerialExecution, bool);
356 vtkBooleanMacro(ForceSerialExecution, bool);
358protected:
360
366 vtkTypeBool IgnorePipelineTime; // whether to use the pipeline time for termination
368
369 // Control execution as serial or threaded
371
373
376
377 //
378 // Make sure the pipeline knows what type we expect as input
379 //
380 int FillInputPortInformation(int port, vtkInformation* info) override;
381
385 int Initialize(vtkInformation* request, vtkInformationVector** inputVector,
386 vtkInformationVector* outputVector) override;
387
395 int Execute(vtkInformation* request, vtkInformationVector** inputVector,
396 vtkInformationVector* outputVector) override;
397
402
407 virtual std::vector<vtkDataSet*> GetSeedSources(vtkInformationVector* inputVector);
408
409 // Initialization of input (vector-field) geometry
411
417
424
426 vtkParticleTracerBaseNamespace::ParticleVector& candidates, std::vector<int>& passed);
427
434 virtual void AssignSeedsToProcessors(double time, vtkDataSet* source,
436
442 bool SendReceiveParticles(std::vector<vtkIdType>&);
443
449
456
461 double currentTime, double targetTime, vtkInitialValueProblemSolver* integrator,
462 vtkTemporalInterpolatedVelocityField* interpolator, vtkDoubleArray* cellVectors,
463 std::atomic<vtkIdType>& particleCount, std::mutex& eraseMutex, bool sequential);
464
472 double pos[4], double p2[4], double intersection[4], vtkGenericCell* cell);
473
474 //
475 // Scalar arrays that are generated as each particle is updated
476 //
478
488
489 // utility function we use to test if a point is inside any of our local datasets
490 bool InsideBounds(double point[]);
491
493 vtkGenericCell* cell, double pcoords[3], vtkDoubleArray* cellVectors, double vorticity[3]);
494
495 //------------------------------------------------------
496
498 vtkTemporalInterpolatedVelocityField* interpolator, vtkDoubleArray* cellVectors);
499
500 std::unordered_map<vtkIdType, int> InjectedPointIdToProcessId;
501
503
508 virtual bool IsPointDataValid(vtkDataObject* input);
509 bool IsPointDataValid(vtkCompositeDataSet* input, std::vector<std::string>& arrayNames);
510 void GetPointDataArrayNames(vtkDataSet* input, std::vector<std::string>& names);
512
513 vtkGetMacro(ReinjectionCounter, int);
514
515 void ResizeArrays(vtkIdType numTuples);
516
521 virtual void InitializeExtraPointDataArrays(vtkPointData* vtkNotUsed(outputPD)) {}
522
525 {
526 }
527
529
534 virtual void AddRestartSeeds(vtkInformationVector** /*inputVector*/) {}
535
546 double delT, int subSteps, vtkTemporalInterpolatedVelocityField* interpolator);
547
549
551
552 // Parameters of tracing
559
560 // A counter to keep track of how many times we reinjected
562
563 // Important for Caching of Cells/Ids/Weights etc
567
568 // Innjection parameters
572
573 // Particle writing to disk
577
579
580 // The velocity interpolator
582
583 // MPI controller needed when running in parallel
585
591
599 std::unordered_map<vtkIdType, vtkParticleTracerBaseNamespace::ParticleInformation> MPIRecvList;
600
601 // Cache bounds info for each dataset we will use repeatedly
602 struct bounds_t
603 {
604 double b[6];
605 };
606 using bounds = struct bounds_t;
607 std::vector<bounds> CachedBounds[2];
608
609 // variables used by Execute() to produce output
610
612
615
626
627 // temp array
629
631 void operator=(const vtkParticleTracerBase&) = delete;
633
634 unsigned int NumberOfParticles();
635
636 friend class ParticlePathFilterInternal;
637 friend class StreaklineFilterInternal;
638
639 static const double Epsilon;
640
641private:
642 // Internal method to initialize CachedData[1] using the provided input
643 void InitializeNextCachedData(vtkDataObject* input);
644
645 // Data for time step CurrentTimeStep-1 and CurrentTimeStep
647};
648
649VTK_ABI_NAMESPACE_END
650#endif
abstract class to write particle data to file
Proxy object to connect input/output ports.
abstract superclass for composite (multi-block or AMR) datasets
abstract superclass for arrays of numeric data
general representation of visualization data
abstract class to specify dataset behavior
Definition vtkDataSet.h:165
dynamic, self-adjusting array of double
dynamic, self-adjusting array of float
provides thread-safe access to cells
a simple class to control print indentation
Definition vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Integrate a set of ordinary differential equations (initial value problem) in time.
dynamic, self-adjusting array of int
Composite dataset that organizes datasets into blocks.
Multiprocessing communication superclass.
Allocate and hold a VTK object.
Definition vtkNew.h:167
A particle tracer for vector fields.
vtkFloatArray * GetParticleVorticity(vtkPointData *)
vtkInitialValueProblemSolver * Integrator
vtkSmartPointer< vtkPointData > ProtoPD
ProtoPD is used just to keep track of the input array names and number of components for copy allocat...
vtkNew< vtkPoints > OutputCoordinates
vtkNew< vtkPointData > OutputPointData
void SetComputeVorticity(bool)
Turn on/off vorticity computation at streamline points (necessary for generating proper stream-ribbon...
vtkIntArray * GetInjectedPointIds(vtkPointData *)
vtkFloatArray * GetParticleAge(vtkPointData *)
void TestParticles(vtkParticleTracerBaseNamespace::ParticleVector &candidates, std::vector< int > &passed)
virtual void AddRestartSeeds(vtkInformationVector **)
For restarts of particle paths, we add in the ability to add in particles from a previous computation...
bool IsPointDataValid(vtkCompositeDataSet *input, std::vector< std::string > &arrayNames)
Methods that check that the input arrays are ordered the same on all data sets.
void SetInterpolatorTypeToCellLocator()
Set the velocity field interpolator type to one that uses a cell locator to perform spatial searching...
vtkNew< vtkFloatArray > ParticleRotation
vtkNew< vtkDoubleArray > CellVectors
vtkNew< vtkIdTypeArray > ParticleCellsConnectivity
vtkIntArray * GetInjectedStepIds(vtkPointData *)
virtual bool UpdateParticleListFromOtherProcesses()
this is used during classification of seed points and also between iterations of the main loop as par...
vtkIntArray * GetErrorCodeArr(vtkPointData *)
bool RetryWithPush(vtkParticleTracerBaseNamespace::ParticleInformation &info, double *point1, double delT, int subSteps, vtkTemporalInterpolatedVelocityField *interpolator)
When particles leave the domain, they must be collected and sent to the other processes for possible ...
vtkSmartPointer< vtkTemporalInterpolatedVelocityField > Interpolator
int Initialize(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
Resets internal cache for a clean start.
vtkNew< vtkIdTypeArray > InjectedPointIds
void operator=(const vtkParticleTracerBase &)=delete
void UpdateParticleList(vtkParticleTracerBaseNamespace::ParticleVector &candidates)
and sending between processors, into a list, which is used as the master list on this processor
vtkGetFilePathMacro(ParticleFileName)
Set/Get the filename to be used with the particle writer when dumping particles to disk.
virtual void SetToExtraPointDataArrays(vtkIdType, vtkParticleTracerBaseNamespace::ParticleInformation &)
vtkFloatArray * GetParticleAngularVel(vtkPointData *)
bool SendReceiveParticles(std::vector< vtkIdType > &)
this is used during classification of seed points and also between iterations of the main loop as par...
std::unordered_map< vtkIdType, int > InjectedPointIdToProcessId
void AddSourceConnection(vtkAlgorithmOutput *input)
Provide support for multiple seed sources.
vtkNew< vtkIntArray > ErrorCodeArray
void RemoveAllSources()
Provide support for multiple seed sources.
void CalculateVorticity(vtkGenericCell *cell, double pcoords[3], vtkDoubleArray *cellVectors, double vorticity[3])
void CreateProtoPD(vtkDataObject *input)
std::unordered_map< vtkIdType, vtkParticleTracerBaseNamespace::ParticleInformation > MPIRecvList
Storage of the particles we received.
virtual void SetController(vtkMultiProcessController *)
Get/Set the controller to use.
vtkIntArray * GetParticleIds(vtkPointData *)
vtkNew< vtkIntArray > InjectedStepIds
MeshOverTimeTypes
Types of Variance of Mesh over time.
void SetForceReinjectionEveryNSteps(int)
When animating particles, it is nice to inject new ones every Nth step to produce a continuous flow.
void SetTerminalSpeed(double)
Specify the terminal speed value, below which integration is terminated.
virtual void SetParticleWriter(vtkAbstractParticleWriter *pw)
Set/Get the Writer associated with this Particle Tracer Ideally a parallel IO capable vtkH5PartWriter...
virtual std::vector< vtkDataSet * > GetSeedSources(vtkInformationVector *inputVector)
Method to get the data set seed sources.
void SetRotationScale(double)
This can be used to scale the rate with which the streamribbons twist.
vtkParticleTracerBaseNamespace::ParticleVector MPISendList
Storage of the particles we want to send to another rank.
vtkParticleTracerBaseNamespace::ParticleDataList ParticleHistories
ProtoPD is used just to keep track of the input array names and number of components for copy allocat...
vtkParticleTracerBase(const vtkParticleTracerBase &)=delete
vtkFloatArray * GetParticleRotation(vtkPointData *)
void IntegrateParticle(vtkParticleTracerBaseNamespace::ParticleListIterator &it, double currentTime, double targetTime, vtkInitialValueProblemSolver *integrator, vtkTemporalInterpolatedVelocityField *interpolator, vtkDoubleArray *cellVectors, std::atomic< vtkIdType > &particleCount, std::mutex &eraseMutex, bool sequential)
particle between the two times supplied.
void GetPointDataArrayNames(vtkDataSet *input, std::vector< std::string > &names)
Methods that check that the input arrays are ordered the same on all data sets.
vtkNew< vtkFloatArray > ParticleVorticity
void SetInterpolatorType(int interpolatorType)
Set the type of the velocity field interpolator to determine whether INTERPOLATOR_WITH_DATASET_POINT_...
virtual void SetMeshOverTime(int meshOverTime)
bool SetTerminationTimeNoModify(double t)
int Finalize(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
Generates an output using the data provided after Execute was ran.
vtkTypeBool IgnorePipelineTime
ProtoPD is used just to keep track of the input array names and number of components for copy allocat...
bool InsideBounds(double point[])
void TestParticles(vtkParticleTracerBaseNamespace::ParticleVector &candidates, vtkParticleTracerBaseNamespace::ParticleVector &passed)
inside our data.
virtual bool IsPointDataValid(vtkDataObject *input)
Methods that check that the input arrays are ordered the same on all data sets.
vtkTemporalInterpolatedVelocityField * GetInterpolator()
void SetIntegratorType(int type)
vtkNew< vtkIntArray > ParticleIds
vtkParticleTracerBaseNamespace::ParticleVector LocalSeeds
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
unsigned int NumberOfParticles()
vtkNew< vtkSignedCharArray > ParticleSourceIds
virtual void InitializeExtraPointDataArrays(vtkPointData *vtkNotUsed(outputPD))
Methods to append values to existing point data arrays that may only be desired on specific concrete ...
void SetParticle(vtkParticleTracerBaseNamespace::ParticleInformation &info, vtkTemporalInterpolatedVelocityField *interpolator, vtkDoubleArray *cellVectors)
void EnqueueParticleToAnotherProcess(vtkParticleTracerBaseNamespace::ParticleInformation &)
vtkNew< vtkFloatArray > ParticleAngularVel
~vtkParticleTracerBase() override
static const double Epsilon
virtual void AssignSeedsToProcessors(double time, vtkDataSet *source, vtkParticleTracerBaseNamespace::ParticleVector &localSeedPoints)
all the injection/seed points according to which processor they belong to.
void SetIntegrator(vtkInitialValueProblemSolver *)
bool ComputeDomainExitLocation(double pos[4], double p2[4], double intersection[4], vtkGenericCell *cell)
This is an old routine kept for possible future use.
void ResizeArrays(vtkIdType numTuples)
vtkNew< vtkFloatArray > ParticleAge
void SetInterpolatorTypeToDataSetPointLocator()
Set the velocity field interpolator type to one that uses a point locator to perform local spatial se...
vtkSmartPointer< vtkMultiProcessController > Controller
vtkSignedCharArray * GetParticleSourceIds(vtkPointData *)
vtkSetFilePathMacro(ParticleFileName)
Set/Get the filename to be used with the particle writer when dumping particles to disk.
virtual vtkMultiProcessController * GetController()
Get/Set the controller to use.
vtkAbstractParticleWriter * ParticleWriter
vtkSmartPointer< vtkDataSet > Seeds
All ranks have the same representation of the seeds.
int Execute(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
Moves the particles one time step further.
composite dataset to encapsulates a dataset consisting of partitions.
represent and manipulate point attribute data
represent and manipulate 3D points
Definition vtkPoints.h:139
Superclass for algorithms that produce only polydata as output.
concrete dataset represents vertices, lines, polygons, and triangle strips
dynamic, self-adjusting array of signed char
Hold a reference to a vtkObjectBase instance.
A helper class for interpolating between times during particle tracing.
record modification and/or execution time
ParticleVector::iterator ParticleIterator
std::list< ParticleInformation > ParticleDataList
ParticleDataList::iterator ParticleListIterator
std::vector< ParticleInformation > ParticleVector
int vtkTypeBool
Definition vtkABI.h:64
boost::graph_traits< vtkGraph * >::vertex_descriptor source(boost::graph_traits< vtkGraph * >::edge_descriptor e, vtkGraph *)
#define VTK_DEPRECATED_IN_9_4_0(reason)
#define vtkCreateWrappedTemporalAlgorithmInterface()
int vtkIdType
Definition vtkType.h:332