1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130
|
// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
// SPDX-License-Identifier: BSD-3-Clause
/**
* @class vtkPLagrangianParticleTracker
* @brief parallel Lagrangian particle tracker
*
* This class implements parallel Lagrangian particle tracker.
* The implementation is as follows:
* First seeds input is parsed to create particle in each rank
* Particles which are not contained by the flow in a rank are sent to other ranks
* which can potentially contain it and will grab only if they actually contain it
* Then each rank begin integrating.
* When a particle goes out of domain, the particle will be sent to other ranks
* the same way.
* When a rank runs out of particle, it waits for other potential particles
* from other ranks.
* When all ranks run out of particles, integration is over.
* The master rank takes care of communications between rank regarding integration termination
* particles are directly streamed rank to rank, without going through the master
*
* @sa
* vtkStreamTracer
*/
#ifndef vtkPLagrangianParticleTracker_h
#define vtkPLagrangianParticleTracker_h
#include "vtkFiltersParallelFlowPathsModule.h" // For export macro
#include "vtkLagrangianParticleTracker.h"
#include "vtkNew.h" // for ivars
#include <map> // for std::map
VTK_ABI_NAMESPACE_BEGIN
class ParticleFeedManager;
class ParticleIdManager;
class ParticleStreamManager;
class vtkMPIController;
class vtkMultiBlockDataSet;
class vtkUnstructuredGrid;
class VTKFILTERSPARALLELFLOWPATHS_EXPORT vtkPLagrangianParticleTracker
: public vtkLagrangianParticleTracker
{
public:
vtkTypeMacro(vtkPLagrangianParticleTracker, vtkLagrangianParticleTracker);
void PrintSelf(ostream& os, vtkIndent indent) override;
static vtkPLagrangianParticleTracker* New();
protected:
vtkPLagrangianParticleTracker();
~vtkPLagrangianParticleTracker() override;
int RequestUpdateExtent(vtkInformation*, vtkInformationVector**, vtkInformationVector*) override;
void GenerateParticles(const vtkBoundingBox* bounds, vtkDataSet* seeds,
vtkDataArray* initialVelocities, vtkDataArray* initialIntegrationTimes, vtkPointData* seedData,
int nVar, std::queue<vtkLagrangianParticle*>& particles) override;
/**
* Flags description :
* Worker flag working : the worker has at least one particle in it's queue and
is currently integrating it.
* Worker flag empty : the worker has no more particle in it's queue and is
actively waiting for more particle to integrate from other ranks.
* Worker flag finished : the worker has received a master empty flag and after
checking one last time, still doesn't have any particle to integrate. It is
now just waiting for master to send the master finished flag.
* Master flag working : there is at least one worker or the master that have one
or more particle to integrate.
* Master flag empty : all ranks, including master, have no more particles to integrate
* Master flag finished : all workers ranks have sent the worker flag finished
*/
void GetParticleFeed(std::queue<vtkLagrangianParticle*>& particleQueue) override;
int Integrate(vtkInitialValueProblemSolver* integrator, vtkLagrangianParticle*,
std::queue<vtkLagrangianParticle*>& particleQueue, vtkPolyData* particlePathsOutput,
vtkPolyLine* particlePath, vtkDataObject* interactionOutput) override;
/**
* Non threadsafe methods to receive particles
*/
void ReceiveParticles(std::queue<vtkLagrangianParticle*>& particleQueue);
/**
* Non threadsafe methods to receive transferred particle ids
*/
void ReceiveTransferredParticleIds();
bool FinalizeOutputs(vtkPolyData* particlePathsOutput, vtkDataObject* interactionOutput) override;
bool UpdateSurfaceCacheIfNeeded(vtkDataObject*& surfaces) override;
/**
* Get an unique id for a particle
* This method is thread safe
*/
vtkIdType GetNewParticleId() override;
/**
* Delete a particle if not out of domain
* If out of domain, it will be stored and deleted later
* in case it needs to be registered as a transferred particle
*/
void DeleteParticle(vtkLagrangianParticle* particle) override;
/**
* Get the complete number of created particles
*/
vtkGetMacro(ParticleCounter, vtkIdType);
void SetController(vtkMPIController*);
vtkNew<vtkUnstructuredGrid> TmpSurfaceInput;
vtkNew<vtkMultiBlockDataSet> TmpSurfaceInputMB;
vtkMPIController* Controller;
ParticleStreamManager* StreamManager;
ParticleIdManager* TransferredParticleIdManager;
ParticleFeedManager* FeedManager;
std::mutex StreamManagerMutex;
std::mutex OutOfDomainParticleMapMutex;
std::map<vtkIdType, vtkLagrangianParticle*> OutOfDomainParticleMap;
private:
vtkPLagrangianParticleTracker(const vtkPLagrangianParticleTracker&) = delete;
void operator=(const vtkPLagrangianParticleTracker&) = delete;
};
VTK_ABI_NAMESPACE_END
#endif
|