VTK  9.3.20240424
vtkLagrangianBasicIntegrationModel.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
42#ifndef vtkLagrangianBasicIntegrationModel_h
43#define vtkLagrangianBasicIntegrationModel_h
44
45#include "vtkFiltersFlowPathsModule.h" // For export macro
46#include "vtkFunctionSet.h"
47#include "vtkNew.h" // For arrays
48#include "vtkWeakPointer.h" // For weak pointer
49
50#include <map> // for array indexes
51#include <mutex> // for mutexes
52#include <queue> // for new particles
53
54VTK_ABI_NAMESPACE_BEGIN
57class vtkCell;
58class vtkCellData;
59class vtkDataArray;
60class vtkDataObject;
61class vtkDataSet;
62class vtkDataSetsType;
63class vtkDoubleArray;
64class vtkFieldData;
65class vtkGenericCell;
67class vtkIntArray;
70class vtkLocatorsType;
73class vtkPointData;
74class vtkPolyData;
75class vtkStringArray;
76class vtkSurfaceType;
78
79class VTKFILTERSFLOWPATHS_EXPORT vtkLagrangianBasicIntegrationModel : public vtkFunctionSet
80{
81public:
83 void PrintSelf(ostream& os, vtkIndent indent) override;
84
85 typedef enum SurfaceType
86 {
87 SURFACE_TYPE_MODEL = 0,
88 SURFACE_TYPE_TERM = 1,
89 SURFACE_TYPE_BOUNCE = 2,
90 SURFACE_TYPE_BREAK = 3,
91 SURFACE_TYPE_PASS = 4
92 } SurfaceType;
93
94 typedef enum VariableStep
95 {
96 VARIABLE_STEP_PREV = -1,
97 VARIABLE_STEP_CURRENT = 0,
98 VARIABLE_STEP_NEXT = 1,
99 } VariableStep;
100
101 typedef std::pair<unsigned int, vtkLagrangianParticle*> PassThroughParticlesItem;
102 typedef std::queue<PassThroughParticlesItem> PassThroughParticlesType;
103
104 using Superclass::FunctionValues;
112 int FunctionValues(double* x, double* f, void* userData) override;
113
115
121 virtual void SetLocator(vtkAbstractCellLocator* locator);
122 vtkGetObjectMacro(Locator, vtkAbstractCellLocator);
124
126
129 vtkGetMacro(LocatorsBuilt, bool);
130 vtkSetMacro(LocatorsBuilt, bool);
132
137
139
148 virtual void AddDataSet(
149 vtkDataSet* dataset, bool surface = false, unsigned int surfaceFlatIndex = 0);
150 virtual void ClearDataSets(bool surface = false);
152
154
157 vtkSetMacro(UseInitialIntegrationTime, bool);
158 vtkGetMacro(UseInitialIntegrationTime, bool);
159 vtkBooleanMacro(UseInitialIntegrationTime, bool);
161
163
166 vtkGetMacro(Tolerance, double);
168
170
173 vtkGetMacro(LocatorTolerance, double);
175
194 std::queue<vtkLagrangianParticle*>& particles, unsigned int& interactedSurfaceFlatIndex,
195 PassThroughParticlesType& passThroughParticles);
196
203 int idx, int port, int connection, int fieldAssociation, const char* name);
204
214 virtual bool FindInLocators(double* x, vtkLagrangianParticle* particle, vtkDataSet*& dataset,
215 vtkIdType& cellId, vtkAbstractCellLocator*& loc, double*& weights);
216
218
222 virtual bool FindInLocators(
223 double* x, vtkLagrangianParticle* particle, vtkDataSet*& dataset, vtkIdType& cellId);
224 virtual bool FindInLocators(double* x, vtkLagrangianParticle* particle);
226
231 virtual void InitializeParticle(vtkLagrangianParticle* vtkNotUsed(particle)) {}
232
241 virtual bool CheckAdaptiveStepReintegration(vtkLagrangianParticle* vtkNotUsed(particle))
242 {
243 return true;
244 }
245
254 virtual bool CheckFreeFlightTermination(vtkLagrangianParticle* vtkNotUsed(particle))
255 {
256 return false;
257 }
258
260
263 vtkSetMacro(NonPlanarQuadSupport, bool);
264 vtkGetMacro(NonPlanarQuadSupport, bool);
265 vtkBooleanMacro(NonPlanarQuadSupport, bool);
267
273
279
285
291
297
303
309
315
317
321 virtual int GetWeightsSize();
323
347 virtual bool ManualIntegration(vtkInitialValueProblemSolver* integrator, double* xcur,
348 double* xnext, double t, double& delT, double& delTActual, double minStep, double maxStep,
349 double maxError, double cellLength, double& error, int& integrationResult,
350 vtkLagrangianParticle* particle);
351
357 virtual void ParallelManualShift(vtkLagrangianParticle* vtkNotUsed(particle)) {}
358
364
370
376 virtual bool FinalizeOutputs(
377 vtkPolyData* vtkNotUsed(particlePathsOutput), vtkDataObject* vtkNotUsed(interractionOutput))
378 {
379 return true;
380 }
381
385 virtual void PreParticleInitalization() {}
386
390 virtual void PreIntegrate(std::queue<vtkLagrangianParticle*>& vtkNotUsed(particles)) {}
391
396 virtual vtkAbstractArray* GetSeedArray(int idx, vtkPointData* pointData);
397
405 vtkSetMacro(NumberOfTrackedUserData, int);
406 vtkGetMacro(NumberOfTrackedUserData, int);
407
413 virtual void InitializePathData(vtkFieldData* data);
414
421
427 virtual void InitializeParticleData(vtkFieldData* particleData, int maxTuples = 0);
428
434 virtual void InsertPathData(vtkLagrangianParticle* particle, vtkFieldData* data);
435
442
449 virtual void InsertParticleData(
450 vtkLagrangianParticle* particle, vtkFieldData* data, int stepEnum);
451
458
465 virtual void ParticleAboutToBeDeleted(vtkLagrangianParticle* vtkNotUsed(particle)) {}
466
472 vtkLagrangianParticle* vtkNotUsed(particle), vtkPointData* vtkNotUsed(particleData))
473 {
474 }
475
476protected:
479
485 virtual int FunctionValues(vtkLagrangianParticle* particle, vtkDataSet* dataSet, vtkIdType cellId,
486 double* weights, double* x, double* f) = 0;
487
494 virtual vtkIdType FindInLocator(vtkDataSet* dataSet, vtkAbstractCellLocator* locator, double* x,
495 vtkGenericCell* cell, double* weights);
496
503
509 virtual bool BounceParticle(
510 vtkLagrangianParticle* particle, vtkDataSet* surface, vtkIdType cellId);
511
519 virtual bool BreakParticle(vtkLagrangianParticle* particle, vtkDataSet* surface, vtkIdType cellId,
520 std::queue<vtkLagrangianParticle*>& particles);
521
531 virtual bool InteractWithSurface(int surfaceType, vtkLagrangianParticle* particle,
532 vtkDataSet* surface, vtkIdType cellId, std::queue<vtkLagrangianParticle*>& particles);
533
540 virtual bool IntersectWithLine(vtkLagrangianParticle* particle, vtkCell* cell, double p1[3],
541 double p2[3], double tol, double& t, double x[3]);
542
548 vtkLagrangianParticle* particle, double interpolationFactor, bool forceInside = false);
549
556 vtkLagrangianParticle* particle, vtkDataSet* surface, vtkIdType cellId);
557
565
573 virtual bool GetFlowOrSurfaceData(vtkLagrangianParticle* particle, int idx,
574 vtkDataSet* flowDataSet, vtkIdType tupleId, double* weights, double* data);
575
584
591
601 const char* arrayName, vtkDataSet* dataset, int nComponent, double* defaultValues);
602
605 vtkLocatorsType* Locators;
606 vtkDataSetsType* DataSets;
607 int WeightsSize = 0;
608
609 struct ArrayVal
610 {
611 int val[3];
612 };
613 typedef std::pair<ArrayVal, std::string> ArrayMapVal;
614 std::map<int, ArrayMapVal> InputArrays;
615
617 {
618 int nComp;
619 int type;
620 std::vector<std::pair<int, std::string>> enumValues;
621 } SurfaceArrayDescription;
622 std::map<std::string, SurfaceArrayDescription> SurfaceArrayDescriptions;
623
624 vtkSurfaceType* Surfaces;
625 vtkLocatorsType* SurfaceLocators;
626
627 double Tolerance;
628 double LocatorTolerance = 0.001;
631 int NumberOfTrackedUserData = 0;
632
641
644
645private:
647 void operator=(const vtkLagrangianBasicIntegrationModel&) = delete;
648};
649
650VTK_ABI_NAMESPACE_END
651#endif
Abstract superclass for all arrays.
an abstract base class for locators which find cells
represent and manipulate cell attribute data
abstract class to specify cell behavior
Definition vtkCell.h:130
abstract superclass for arrays of numeric data
general representation of visualization data
abstract class to specify dataset behavior
Definition vtkDataSet.h:166
dynamic, self-adjusting array of double
represent and manipulate fields of data
Abstract interface for sets of functions.
provides thread-safe access to cells
a simple class to control print indentation
Definition vtkIndent.h:108
Integrate a set of ordinary differential equations (initial value problem) in time.
dynamic, self-adjusting array of int
vtkFunctionSet abstract implementation to be used in the vtkLagrangianParticleTracker integrator.
virtual void InsertInteractionData(vtkLagrangianParticle *particle, vtkFieldData *data)
Method used by the LPT to insert data from the particle into the provided vtkFieldData.
virtual void InitializeParticle(vtkLagrangianParticle *vtkNotUsed(particle))
Initialize a particle by setting user variables and perform any user model specific operation.
virtual void SetInputArrayToProcess(int idx, int port, int connection, int fieldAssociation, const char *name)
Set a input array to process at a specific index, identified by a port, connection,...
virtual bool FinalizeOutputs(vtkPolyData *vtkNotUsed(particlePathsOutput), vtkDataObject *vtkNotUsed(interractionOutput))
Enable model post process on output Return true if successful, false otherwise Empty and Always retur...
virtual int GetFlowOrSurfaceDataFieldAssociation(int idx)
Recover a field association for a specified array index if it has been set using SetInputArrayToProce...
virtual vtkStringArray * GetSurfaceArrayEnumValues()
Get the surface arrays expected values and associated enums Used Only be the vtkLagrangianSurfaceHelp...
virtual vtkStringArray * GetSurfaceArrayNames()
Get the surface arrays expected name Used Only be the vtkLagrangianSurfaceHelper in ParaView plugins.
std::queue< PassThroughParticlesItem > PassThroughParticlesType
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
virtual void InitializeParticleData(vtkFieldData *particleData, int maxTuples=0)
Method used by the LPT to initialize data insertion in the provided vtkFieldData.
virtual bool TerminateParticle(vtkLagrangianParticle *particle)
Terminate a particle, by positioning flags.
virtual bool FindInLocators(double *x, vtkLagrangianParticle *particle, vtkDataSet *&dataset, vtkIdType &cellId)
Convenience methods to call FindInLocators with less arguments THESE METHODS ARE NOT THREAD-SAFE.
virtual void InterpolateNextParticleVariables(vtkLagrangianParticle *particle, double interpolationFactor, bool forceInside=false)
compute all particle variables using interpolation factor This method is thread-safe.
virtual void PreIntegrate(std::queue< vtkLagrangianParticle * > &vtkNotUsed(particles))
Enable model to modify particle before integration.
virtual void AddDataSet(vtkDataSet *dataset, bool surface=false, unsigned int surfaceFlatIndex=0)
Add a dataset to locate cells in This create a specific locator for the provided dataset using the Lo...
virtual void InitializePathData(vtkFieldData *data)
Method used by the LPT to initialize data insertion in the provided vtkFieldData.
virtual vtkLagrangianThreadedData * InitializeThreadedData()
Let the model allocate and initialize a threaded data.
virtual void SetTracker(vtkLagrangianParticleTracker *Tracker)
Set the parent tracker.
virtual bool CheckSurfacePerforation(vtkLagrangianParticle *particle, vtkDataSet *surface, vtkIdType cellId)
Given a particle, check if it perforate a surface cell ie : interact with next step after interacting...
virtual bool ManualIntegration(vtkInitialValueProblemSolver *integrator, double *xcur, double *xnext, double t, double &delT, double &delTActual, double minStep, double maxStep, double maxError, double cellLength, double &error, int &integrationResult, vtkLagrangianParticle *particle)
Let the model define it's own way to integrate Signature is very close to the integrator method signa...
virtual void SetLocator(vtkAbstractCellLocator *locator)
Set/Get the locator used to locate cells in the datasets.
virtual void InsertPathData(vtkLagrangianParticle *particle, vtkFieldData *data)
Method used by the LPT to insert data from the particle into the provided vtkFieldData.
virtual bool IntersectWithLine(vtkLagrangianParticle *particle, vtkCell *cell, double p1[3], double p2[3], double tol, double &t, double x[3])
Call vtkCell::IntersectWithLine This method is to be reimplemented in inherited classes willing to im...
std::pair< unsigned int, vtkLagrangianParticle * > PassThroughParticlesItem
virtual void InitializeInteractionData(vtkFieldData *data)
Method used by the LPT to initialize data insertion in the provided vtkFieldData.
virtual bool BounceParticle(vtkLagrangianParticle *particle, vtkDataSet *surface, vtkIdType cellId)
Bounce a particle, using the normal of the cell it bounces on.
virtual void ParticleAboutToBeDeleted(vtkLagrangianParticle *vtkNotUsed(particle))
Method to be reimplemented if needed in inherited classes.
virtual void FinalizeThreadedData(vtkLagrangianThreadedData *&data)
Let the model finalize and deallocate a user data at thread level This method is called serially for ...
int FunctionValues(double *x, double *f, void *userData) override
Evaluate integration model velocity f at position x.
virtual vtkIntArray * GetSurfaceArrayTypes()
Get the surface arrays expected type Used Only be the vtkLagrangianSurfaceHelper in ParaView plugins.
virtual bool CheckFreeFlightTermination(vtkLagrangianParticle *vtkNotUsed(particle))
Method to be reimplemented if needed in inherited classes.
virtual int GetFlowOrSurfaceDataNumberOfComponents(int idx, vtkDataSet *dataSet)
Recover the number of components for a specified array index if it has been set using SetInputArrayTo...
virtual bool CheckAdaptiveStepReintegration(vtkLagrangianParticle *vtkNotUsed(particle))
Method to be reimplemented if needed in inherited classes.
virtual void ClearDataSets(bool surface=false)
Add a dataset to locate cells in This create a specific locator for the provided dataset using the Lo...
virtual void InsertSurfaceInteractionData(vtkLagrangianParticle *vtkNotUsed(particle), vtkPointData *vtkNotUsed(particleData))
Method to be reimplemented if needed in inherited classes.
virtual vtkIntArray * GetSeedArrayComps()
Get the seed arrays expected number of components Used Only be the vtkLagrangianSeedHelper in ParaVie...
std::map< std::string, SurfaceArrayDescription > SurfaceArrayDescriptions
virtual vtkStringArray * GetSeedArrayNames()
Get the seed arrays expected name Used Only be the vtkLagrangianSeedHelper in ParaView plugins.
virtual void ComputeSurfaceDefaultValues(const char *arrayName, vtkDataSet *dataset, int nComponent, double *defaultValues)
Method used by ParaView surface helper to get default values for each leaf of each dataset of surface...
virtual bool FindInLocators(double *x, vtkLagrangianParticle *particle, vtkDataSet *&dataset, vtkIdType &cellId, vtkAbstractCellLocator *&loc, double *&weights)
Look for a dataset in this integration model containing the point x.
virtual vtkLagrangianParticle * ComputeSurfaceInteraction(vtkLagrangianParticle *particle, std::queue< vtkLagrangianParticle * > &particles, unsigned int &interactedSurfaceFlatIndex, PassThroughParticlesType &passThroughParticles)
Interact the current particle with a surfaces Return a particle to record as interaction point if not...
virtual bool FindInLocators(double *x, vtkLagrangianParticle *particle)
Convenience methods to call FindInLocators with less arguments THESE METHODS ARE NOT THREAD-SAFE.
virtual void PreParticleInitalization()
Allow for model setup prior to Particle Initialization.
virtual vtkAbstractArray * GetSeedArray(int idx, vtkPointData *pointData)
Get a seed array, as set in setInputArrayToProcess from the provided seed point data.
virtual vtkAbstractArray * GetSeedArray(int idx, vtkLagrangianParticle *particle)
Get a seed array, as set in setInputArrayToProcess from the provided particle seed data Access then t...
virtual void ParallelManualShift(vtkLagrangianParticle *vtkNotUsed(particle))
Method called by parallel algorithm after receiving a particle from stream if PManualShift flag has b...
virtual bool InteractWithSurface(int surfaceType, vtkLagrangianParticle *particle, vtkDataSet *surface, vtkIdType cellId, std::queue< vtkLagrangianParticle * > &particles)
Call vtkLagrangianBasicIntegrationModel::Terminate This method is to be reimplemented in inherited cl...
virtual void InsertParticleSeedData(vtkLagrangianParticle *particle, vtkFieldData *data)
Method used by the LPT to insert data from the particle into the provided vtkFieldData.
virtual vtkIntArray * GetSurfaceArrayComps()
Get the seed array expected number of components Used Only be the vtkLagrangianSurfaceHelper in ParaV...
virtual int FunctionValues(vtkLagrangianParticle *particle, vtkDataSet *dataSet, vtkIdType cellId, double *weights, double *x, double *f)=0
Actually compute the integration model velocity field pure abstract, to be implemented in inherited c...
virtual vtkIntArray * GetSeedArrayTypes()
Get the seed arrays expected type Used Only be the vtkLagrangianSeedHelper in ParaView plugins.
virtual vtkDoubleArray * GetSurfaceArrayDefaultValues()
Get the surface arrays default values for each leaf Used Only be the vtkLagrangianSurfaceHelper in Pa...
vtkWeakPointer< vtkLagrangianParticleTracker > Tracker
virtual vtkIdType FindInLocator(vtkDataSet *dataSet, vtkAbstractCellLocator *locator, double *x, vtkGenericCell *cell, double *weights)
Look in the given dataset and associated locator to see if it contains the point x,...
virtual int GetWeightsSize()
Get the maximum weights size necessary for calling FindInLocators with weights.
virtual void InsertParticleData(vtkLagrangianParticle *particle, vtkFieldData *data, int stepEnum)
Method used by the LPT to insert data from the particle into the provided vtkFieldData.
virtual bool BreakParticle(vtkLagrangianParticle *particle, vtkDataSet *surface, vtkIdType cellId, std::queue< vtkLagrangianParticle * > &particles)
Breakup a particle at intersection point, by terminating it and creating two new particle using the i...
virtual bool GetFlowOrSurfaceData(vtkLagrangianParticle *particle, int idx, vtkDataSet *flowDataSet, vtkIdType tupleId, double *weights, double *data)
Directly get a double value from flow or surface data as defined in SetInputArrayToProcess.
Filter to inject and track particles in a flow.
Basis class for Lagrangian particles.
Composite dataset that organizes datasets into blocks.
composite dataset to encapsulates pieces of dataset.
Allocate and hold a VTK object.
Definition vtkNew.h:160
represent and manipulate point attribute data
concrete dataset represents vertices, lines, polygons, and triangle strips
a vtkAbstractArray subclass for strings
a weak reference to a vtkObject.
struct to hold a user data
int vtkIdType
Definition vtkType.h:315