VTK  9.4.20241218
vtkTemporalInterpolatedVelocityField.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
35#ifndef vtkTemporalInterpolatedVelocityField_h
36#define vtkTemporalInterpolatedVelocityField_h
37
38#include "vtkFiltersFlowPathsModule.h" // For export macro
39#include "vtkFunctionSet.h"
40#include "vtkSmartPointer.h" // For vtkSmartPointer
41
42#include <vector> // For internal structures
43
44VTK_ABI_NAMESPACE_BEGIN
48class vtkDataArray;
49class vtkDataSet;
50class vtkDoubleArray;
52class vtkGenericCell;
53class vtkLocator;
54class vtkPointData;
55
56class VTKFILTERSFLOWPATHS_EXPORT vtkTemporalInterpolatedVelocityField : public vtkFunctionSet
57{
58public:
60 void PrintSelf(ostream& os, vtkIndent indent) override;
61
67
72 {
73 INSIDE_ALL = 0,
74 OUTSIDE_ALL = 1,
75 OUTSIDE_T0 = 2,
76 OUTSIDE_T1 = 3
77 };
78
83 {
84 DIFFERENT = 0,
85 STATIC = 1,
86 LINEAR_TRANSFORMATION = 2,
87 SAME_TOPOLOGY = 3
88 };
89
91 /*
92 * Set/Get the type of variance of the mesh over time.
93 *
94 * DIFFERENT = 0,
95 * STATIC = 1,
96 * LINEAR_TRANSFORMATION = 2
97 * SAME_TOPOLOGY = 3
98 */
99 vtkSetClampMacro(MeshOverTime, int, DIFFERENT, SAME_TOPOLOGY);
100 void SetMeshOverTimeToDifferent() { this->SetMeshOverTime(DIFFERENT); }
101 void SetMeshOverTimeToStatic() { this->SetMeshOverTime(STATIC); }
102 void SetMeshOverTimeToLinearTransformation() { this->SetMeshOverTime(LINEAR_TRANSFORMATION); }
103 void SetMeshOverTimeToSameTopology() { this->SetMeshOverTime(SAME_TOPOLOGY); }
104 vtkGetMacro(MeshOverTime, int);
106
115
123
124 using Superclass::FunctionValues;
126
130 int FunctionValues(double* x, double* u) override;
131 int FunctionValuesAtT(int T, double* x, double* u);
133
139 void SelectVectors(const char* fieldName) { this->SetVectorsSelection(fieldName); }
140
142
147 void AddDataSetAtTime(int N, double T, vtkDataSet* dataset);
149
151
156 bool GetCachedCellIds(vtkIdType id[2], int ds[2]);
157 void SetCachedCellIds(vtkIdType id[2], int ds[2]);
159
165
167
171 int TestPoint(double* x);
172 int QuickTestPoint(double* x);
174
176
180 vtkGetVector3Macro(LastGoodVelocity, double);
182
184
187 vtkGetMacro(CurrentWeight, double);
189
190 bool InterpolatePoint(vtkPointData* outPD1, vtkPointData* outPD2, vtkIdType outIndex);
191
192 bool InterpolatePoint(int T, vtkPointData* outPD1, vtkIdType outIndex);
193
195 int T, double pcoords[3], double* weights, vtkGenericCell*& cell, vtkDoubleArray* cellVectors);
196
198
200
202
209 vtkGetObjectMacro(FindCellStrategy, vtkFindCellStrategy);
211
212protected:
215
216 virtual void SetVectorsSelection(const char* v);
217
218 int MeshOverTime = MeshOverTimeTypes::DIFFERENT;
219
221 const std::vector<vtkDataSet*>& datasets, vtkFindCellStrategy* strategy,
222 const std::vector<vtkSmartPointer<vtkLocator>>& locators,
223 const std::vector<vtkSmartPointer<vtkAbstractCellLinks>>& links);
224
225 void CreateLocators(const std::vector<vtkDataSet*>& datasets, vtkFindCellStrategy* strategy,
226 std::vector<vtkSmartPointer<vtkLocator>>& locators);
227 void CreateLinks(const std::vector<vtkDataSet*>& datasets,
228 std::vector<vtkSmartPointer<vtkAbstractCellLinks>>& links);
230 std::vector<vtkSmartPointer<vtkLocator>>& linearCellLocators);
231
232 double Vals1[3];
233 double Vals2[3];
234 double Times[2];
235 double LastGoodVelocity[3];
236
237 static const double WEIGHT_TO_TOLERANCE;
238
239 // The weight (0.0->1.0) of the value of T between the two available
240 // time values for the current computation
242 // One minus the CurrentWeight
244 // A scaling factor used when calculating the CurrentWeight { 1.0/(T2-T1) }
246
248 std::vector<vtkSmartPointer<vtkLocator>> Locators[2];
249 std::vector<vtkSmartPointer<vtkLocator>> InitialCellLocators;
250 std::vector<vtkSmartPointer<vtkAbstractCellLinks>> Links[2];
251 std::vector<size_t> MaxCellSizes[2];
252
254
255private:
256 // Hide this since we need multiple time steps and are using a different
257 // function prototype
258 virtual void AddDataSet(vtkDataSet*) {}
259
261 void operator=(const vtkTemporalInterpolatedVelocityField&) = delete;
262};
263
264VTK_ABI_NAMESPACE_END
265#endif
abstract superclass for composite (multi-block or AMR) datasets
An abstract class for obtaining the interpolated velocity values at a point.
abstract superclass for arrays of numeric data
abstract class to specify dataset behavior
Definition vtkDataSet.h:165
dynamic, self-adjusting array of double
helper class to manage the vtkPointSet::FindCell() method
Abstract interface for sets of functions.
provides thread-safe access to cells
a simple class to control print indentation
Definition vtkIndent.h:108
abstract base class for objects that accelerate spatial searches
Definition vtkLocator.h:78
represent and manipulate point attribute data
Hold a reference to a vtkObjectBase instance.
A helper class for interpolating between times during particle tracing.
int QuickTestPoint(double *x)
A utility function which evaluates the point at T1, T2 to see if it is inside the data at both times ...
int TestPoint(double *x)
A utility function which evaluates the point at T1, T2 to see if it is inside the data at both times ...
IDStates
States that define where the cell id are.
void AddDataSetAtTime(int N, double T, vtkDataSet *dataset)
In order to use this class, two sets of data must be supplied, corresponding to times T1 and T2.
MeshOverTimeTypes
Types of Variance of Mesh over time.
virtual void CopyParameters(vtkTemporalInterpolatedVelocityField *from)
Copy essential parameters between instances of this class.
bool InterpolatePoint(vtkPointData *outPD1, vtkPointData *outPD2, vtkIdType outIndex)
void SelectVectors(const char *fieldName)
If you want to work with an arbitrary vector array, then set its name here.
virtual void SetVectorsSelection(const char *v)
void CreateLinearTransformCellLocators(const std::vector< vtkSmartPointer< vtkLocator > > &locators, std::vector< vtkSmartPointer< vtkLocator > > &linearCellLocators)
bool GetVorticityData(int T, double pcoords[3], double *weights, vtkGenericCell *&cell, vtkDoubleArray *cellVectors)
std::vector< vtkSmartPointer< vtkLocator > > InitialCellLocators
bool InterpolatePoint(int T, vtkPointData *outPD1, vtkIdType outIndex)
virtual void SetFindCellStrategy(vtkFindCellStrategy *)
Set / get the strategy used to perform the FindCell() operation.
void SetCachedCellIds(vtkIdType id[2], int ds[2])
Between iterations of the Particle Tracer, Id's of the Cell are stored and then at the start of the n...
bool GetCachedCellIds(vtkIdType id[2], int ds[2])
Between iterations of the Particle Tracer, Id's of the Cell are stored and then at the start of the n...
void CreateLinks(const std::vector< vtkDataSet * > &datasets, std::vector< vtkSmartPointer< vtkAbstractCellLinks > > &links)
void ClearCache()
Set the last cell id to -1 so that the next search does not start from the previous cell.
int FunctionValues(double *x, double *u) override
Evaluate the velocity field, f, at (x, y, z, t).
void CreateLocators(const std::vector< vtkDataSet * > &datasets, vtkFindCellStrategy *strategy, std::vector< vtkSmartPointer< vtkLocator > > &locators)
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void InitializeWithLocators(vtkCompositeInterpolatedVelocityField *ivf, const std::vector< vtkDataSet * > &datasets, vtkFindCellStrategy *strategy, const std::vector< vtkSmartPointer< vtkLocator > > &locators, const std::vector< vtkSmartPointer< vtkAbstractCellLinks > > &links)
int FunctionValuesAtT(int T, double *x, double *u)
Evaluate the velocity field, f, at (x, y, z, t).
void Initialize(vtkCompositeDataSet *t0, vtkCompositeDataSet *t1)
The Initialize() method is used to build and cache supporting structures (such as locators) which are...
static vtkTemporalInterpolatedVelocityField * New()
Construct a vtkTemporalInterpolatedVelocityField with no initial data set.
int vtkIdType
Definition vtkType.h:315