VTK  9.0.3
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 class MasterFlagManager;
45 class ParticleStreamManager;
46 class RankFlagManager;
47 class vtkMPIController;
50 
51 class VTKFILTERSPARALLELFLOWPATHS_EXPORT vtkPLagrangianParticleTracker
53 {
54 public:
56  virtual void PrintSelf(ostream& os, vtkIndent indent) override;
58 
59 protected:
62 
63  virtual int RequestUpdateExtent(
65 
66  void GenerateParticles(const vtkBoundingBox* bounds, vtkDataSet* seeds,
67  vtkDataArray* initialVelocities, vtkDataArray* initialIntegrationTimes, vtkPointData* seedData,
68  int nVar, std::queue<vtkLagrangianParticle*>& particles) override;
69 
84  virtual void GetParticleFeed(std::queue<vtkLagrangianParticle*>& particleQueue) override;
86  std::queue<vtkLagrangianParticle*>& particleQueue, vtkPolyData* particlePathsOutput,
87  vtkPolyLine* particlePath, vtkDataObject* interactionOutput) override;
88 
90 
94  void ReceiveParticles(std::queue<vtkLagrangianParticle*>& particleQueue);
96 
97  bool FinalizeOutputs(vtkPolyData* particlePathsOutput, vtkDataObject* interactionOutput) override;
98 
99  bool UpdateSurfaceCacheIfNeeded(vtkDataObject*& surfaces) override;
100 
105  virtual vtkIdType GetNewParticleId() override;
106 
110  vtkGetMacro(ParticleCounter, vtkIdType);
111 
115  ParticleStreamManager* StreamManager;
116  MasterFlagManager* MFlagManager;
117  RankFlagManager* RFlagManager;
118 
119  std::mutex StreamManagerMutex;
120 
121 private:
123  void operator=(const vtkPLagrangianParticleTracker&) = delete;
124 };
125 #endif
Fast, simple class for dealing with 3D bounds.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:50
general representation of visualization data
Definition: vtkDataObject.h:60
abstract class to specify dataset behavior
Definition: vtkDataSet.h:57
a simple class to control print indentation
Definition: vtkIndent.h:34
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Integrate a set of ordinary differential equations (initial value problem) in time.
Filter to inject and track particles in a flow.
Basis class for Lagrangian particles.
Process communication using MPI.
Composite dataset that organizes datasets into blocks.
parallel Lagrangian particle tracker
void SendParticle(vtkLagrangianParticle *particle)
Non threadsafe methods to send and receive particles.
virtual 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...
virtual int RequestUpdateExtent(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
virtual int Integrate(vtkInitialValueProblemSolver *integrator, vtkLagrangianParticle *, std::queue< vtkLagrangianParticle * > &particleQueue, vtkPolyData *particlePathsOutput, vtkPolyLine *particlePath, vtkDataObject *interactionOutput) override
This method is thread safe.
void GenerateParticles(const vtkBoundingBox *bounds, vtkDataSet *seeds, vtkDataArray *initialVelocities, vtkDataArray *initialIntegrationTimes, vtkPointData *seedData, int nVar, std::queue< vtkLagrangianParticle * > &particles) override
~vtkPLagrangianParticleTracker() override
static vtkPLagrangianParticleTracker * New()
vtkNew< vtkUnstructuredGrid > TmpSurfaceInput
virtual void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
bool UpdateSurfaceCacheIfNeeded(vtkDataObject *&surfaces) override
bool FinalizeOutputs(vtkPolyData *particlePathsOutput, vtkDataObject *interactionOutput) override
void ReceiveParticles(std::queue< vtkLagrangianParticle * > &particleQueue)
virtual vtkIdType GetNewParticleId() override
Get an unique id for a particle This method is thread safe.
vtkNew< vtkMultiBlockDataSet > TmpSurfaceInputMB
represent and manipulate point attribute data
Definition: vtkPointData.h:32
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:85
cell represents a set of 1D lines
Definition: vtkPolyLine.h:37
dataset represents arbitrary combinations of all possible cell types
int vtkIdType
Definition: vtkType.h:338