VTK  9.1.0
vtkPLagrangianParticleTracker.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkPLagrangianParticleTracker.h
5 
6  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7  All rights reserved.
8  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9 
10  This software is distributed WITHOUT ANY WARRANTY; without even
11  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12  PURPOSE. See the above copyright notice for more information.
13 
14 =========================================================================*/
37 #ifndef vtkPLagrangianParticleTracker_h
38 #define vtkPLagrangianParticleTracker_h
39 
40 #include "vtkFiltersParallelFlowPathsModule.h" // For export macro
42 #include "vtkNew.h" // for ivars
43 
44 #include <map> // for std::map
45 
46 class ParticleFeedManager;
47 class ParticleIdManager;
48 class ParticleStreamManager;
49 class vtkMPIController;
52 
53 class VTKFILTERSPARALLELFLOWPATHS_EXPORT vtkPLagrangianParticleTracker
55 {
56 public:
58  void PrintSelf(ostream& os, vtkIndent indent) override;
60 
61 protected:
64 
66 
67  void GenerateParticles(const vtkBoundingBox* bounds, vtkDataSet* seeds,
68  vtkDataArray* initialVelocities, vtkDataArray* initialIntegrationTimes, vtkPointData* seedData,
69  int nVar, std::queue<vtkLagrangianParticle*>& particles) override;
70 
85  void GetParticleFeed(std::queue<vtkLagrangianParticle*>& particleQueue) override;
87  std::queue<vtkLagrangianParticle*>& particleQueue, vtkPolyData* particlePathsOutput,
88  vtkPolyLine* particlePath, vtkDataObject* interactionOutput) override;
89 
93  void ReceiveParticles(std::queue<vtkLagrangianParticle*>& particleQueue);
94 
99 
100  bool FinalizeOutputs(vtkPolyData* particlePathsOutput, vtkDataObject* interactionOutput) override;
101 
102  bool UpdateSurfaceCacheIfNeeded(vtkDataObject*& surfaces) override;
103 
109 
115  void DeleteParticle(vtkLagrangianParticle* particle) override;
116 
120  vtkGetMacro(ParticleCounter, vtkIdType);
121 
125  ParticleStreamManager* StreamManager;
126  ParticleIdManager* TransferredParticleIdManager;
127  ParticleFeedManager* FeedManager;
128 
129  std::mutex StreamManagerMutex;
131 
132  std::map<vtkIdType, vtkLagrangianParticle*> OutOfDomainParticleMap;
133 
134 private:
136  void operator=(const vtkPLagrangianParticleTracker&) = delete;
137 };
138 #endif
vtkPLagrangianParticleTracker::ReceiveTransferredParticleIds
void ReceiveTransferredParticleIds()
Non threadsafe methods to receive transferred particle ids.
vtkPLagrangianParticleTracker::GetNewParticleId
vtkIdType GetNewParticleId() override
Get an unique id for a particle This method is thread safe.
vtkPLagrangianParticleTracker::Controller
vtkMPIController * Controller
Definition: vtkPLagrangianParticleTracker.h:124
vtkLagrangianParticleTracker
Filter to inject and track particles in a flow.
Definition: vtkLagrangianParticleTracker.h:118
vtkPLagrangianParticleTracker::StreamManager
ParticleStreamManager * StreamManager
Definition: vtkPLagrangianParticleTracker.h:125
vtkPLagrangianParticleTracker::TransferredParticleIdManager
ParticleIdManager * TransferredParticleIdManager
Definition: vtkPLagrangianParticleTracker.h:126
vtkPointData
represent and manipulate point attribute data
Definition: vtkPointData.h:142
vtkIdType
int vtkIdType
Definition: vtkType.h:332
vtkInformationVector
Store zero or more vtkInformation instances.
Definition: vtkInformationVector.h:145
vtkPLagrangianParticleTracker::TmpSurfaceInputMB
vtkNew< vtkMultiBlockDataSet > TmpSurfaceInputMB
Definition: vtkPLagrangianParticleTracker.h:123
vtkPLagrangianParticleTracker::RequestUpdateExtent
int RequestUpdateExtent(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
vtkPLagrangianParticleTracker::New
static vtkPLagrangianParticleTracker * New()
vtkPLagrangianParticleTracker::vtkPLagrangianParticleTracker
vtkPLagrangianParticleTracker()
vtkPLagrangianParticleTracker::StreamManagerMutex
std::mutex StreamManagerMutex
Definition: vtkPLagrangianParticleTracker.h:129
vtkPLagrangianParticleTracker::UpdateSurfaceCacheIfNeeded
bool UpdateSurfaceCacheIfNeeded(vtkDataObject *&surfaces) override
vtkDataArray
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:159
vtkPLagrangianParticleTracker::Integrate
int Integrate(vtkInitialValueProblemSolver *integrator, vtkLagrangianParticle *, std::queue< vtkLagrangianParticle * > &particleQueue, vtkPolyData *particlePathsOutput, vtkPolyLine *particlePath, vtkDataObject *interactionOutput) override
This method is thread safe.
vtkMultiBlockDataSet
Composite dataset that organizes datasets into blocks.
Definition: vtkMultiBlockDataSet.h:155
vtkPLagrangianParticleTracker::GetParticleFeed
void GetParticleFeed(std::queue< vtkLagrangianParticle * > &particleQueue) override
Flags description : Worker flag working : the worker has at least one particle in it's queue and is c...
vtkLagrangianParticleTracker.h
vtkPLagrangianParticleTracker::TmpSurfaceInput
vtkNew< vtkUnstructuredGrid > TmpSurfaceInput
Definition: vtkPLagrangianParticleTracker.h:120
vtkBoundingBox
Fast, simple class for representing and operating on 3D bounds.
Definition: vtkBoundingBox.h:66
vtkPolyLine
cell represents a set of 1D lines
Definition: vtkPolyLine.h:146
vtkPLagrangianParticleTracker::~vtkPLagrangianParticleTracker
~vtkPLagrangianParticleTracker() override
vtkIndent
a simple class to control print indentation
Definition: vtkIndent.h:113
vtkNew< vtkUnstructuredGrid >
vtkPLagrangianParticleTracker
parallel Lagrangian particle tracker
Definition: vtkPLagrangianParticleTracker.h:55
vtkLagrangianParticle
Basis class for Lagrangian particles.
Definition: vtkLagrangianParticle.h:47
vtkPLagrangianParticleTracker::FinalizeOutputs
bool FinalizeOutputs(vtkPolyData *particlePathsOutput, vtkDataObject *interactionOutput) override
vtkDataSet
abstract class to specify dataset behavior
Definition: vtkDataSet.h:166
vtkPLagrangianParticleTracker::PrintSelf
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkInformation
Store vtkAlgorithm input/output information.
Definition: vtkInformation.h:183
vtkPLagrangianParticleTracker::DeleteParticle
void DeleteParticle(vtkLagrangianParticle *particle) override
Delete a particle if not out of domain If out of domain, it will be stored and deleted later in case ...
vtkPLagrangianParticleTracker::OutOfDomainParticleMap
std::map< vtkIdType, vtkLagrangianParticle * > OutOfDomainParticleMap
Definition: vtkPLagrangianParticleTracker.h:132
vtkNew.h
vtkPLagrangianParticleTracker::GenerateParticles
void GenerateParticles(const vtkBoundingBox *bounds, vtkDataSet *seeds, vtkDataArray *initialVelocities, vtkDataArray *initialIntegrationTimes, vtkPointData *seedData, int nVar, std::queue< vtkLagrangianParticle * > &particles) override
vtkPLagrangianParticleTracker::ReceiveParticles
void ReceiveParticles(std::queue< vtkLagrangianParticle * > &particleQueue)
Non threadsafe methods to receive particles.
vtkPolyData
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:195
vtkUnstructuredGrid
dataset represents arbitrary combinations of all possible cell types
Definition: vtkUnstructuredGrid.h:204
vtkInitialValueProblemSolver
Integrate a set of ordinary differential equations (initial value problem) in time.
Definition: vtkInitialValueProblemSolver.h:41
vtkMPIController
Process communication using MPI.
Definition: vtkMPIController.h:56
vtkDataObject
general representation of visualization data
Definition: vtkDataObject.h:169
vtkPLagrangianParticleTracker::FeedManager
ParticleFeedManager * FeedManager
Definition: vtkPLagrangianParticleTracker.h:127
vtkPLagrangianParticleTracker::OutOfDomainParticleMapMutex
std::mutex OutOfDomainParticleMapMutex
Definition: vtkPLagrangianParticleTracker.h:130