VTK  9.3.20240329
vtkFastSplatter.h
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2 // SPDX-FileCopyrightText: Copyright (c) Sandia Corporation
3 // SPDX-License-Identifier: BSD-3-Clause
52 #ifndef vtkFastSplatter_h
53 #define vtkFastSplatter_h
54 
55 #include "vtkImageAlgorithm.h"
56 #include "vtkImagingHybridModule.h" // For export macro
57 
58 VTK_ABI_NAMESPACE_BEGIN
59 class VTKIMAGINGHYBRID_EXPORT vtkFastSplatter : public vtkImageAlgorithm
60 {
61 public:
63  static vtkFastSplatter* New();
64  void PrintSelf(ostream& os, vtkIndent indent) override;
65 
67 
73  vtkSetVector6Macro(ModelBounds, double);
74  vtkGetVectorMacro(ModelBounds, double, 6);
76 
78 
81  vtkSetVector3Macro(OutputDimensions, int);
82  vtkGetVector3Macro(OutputDimensions, int);
84 
85  enum
86  {
90  FreezeScaleLimit
91  };
92 
94 
100  vtkSetMacro(LimitMode, int);
101  vtkGetMacro(LimitMode, int);
102  void SetLimitModeToNone() { this->SetLimitMode(NoneLimit); }
103  void SetLimitModeToClamp() { this->SetLimitMode(ClampLimit); }
104  void SetLimitModeToScale() { this->SetLimitMode(ScaleLimit); }
105  void SetLimitModeToFreezeScale() { this->SetLimitMode(FreezeScaleLimit); }
107 
109 
112  vtkSetMacro(MinValue, double);
113  vtkGetMacro(MinValue, double);
114  vtkSetMacro(MaxValue, double);
115  vtkGetMacro(MaxValue, double);
117 
119 
123  vtkGetMacro(NumberOfPointsSplatted, int);
125 
132 
133 protected:
135  ~vtkFastSplatter() override;
136 
137  double ModelBounds[6];
138  int OutputDimensions[3];
139 
141  double MinValue;
142  double MaxValue;
143  double FrozenScale;
144 
146 
151 
152  // Used internally for converting points in world space to indices in
153  // the output image.
154  double Origin[3];
155  double Spacing[3];
156 
157  // This is updated every time the filter executes
159 
160  // Used internally to track the data range. When the limit mode is
161  // set to FreezeScale, the data will be scaled as if this were the
162  // range regardless of what it actually is.
165 
166 private:
167  vtkFastSplatter(const vtkFastSplatter&) = delete;
168  void operator=(const vtkFastSplatter&) = delete;
169 };
170 
171 //-----------------------------------------------------------------------------
172 
173 template <class T>
174 void vtkFastSplatterClamp(T* array, vtkIdType arraySize, T minValue, T maxValue)
175 {
176  for (vtkIdType i = 0; i < arraySize; i++)
177  {
178  if (array[i] < minValue)
179  array[i] = minValue;
180  if (array[i] > maxValue)
181  array[i] = maxValue;
182  }
183 }
184 
185 //-----------------------------------------------------------------------------
186 
187 template <class T>
188 void vtkFastSplatterScale(T* array, int numComponents, vtkIdType numTuples, T minValue, T maxValue,
189  double* dataMinValue, double* dataMaxValue)
190 {
191  T* a;
192  T min, max;
193  *dataMinValue = 0;
194  *dataMaxValue = 0;
195  vtkIdType t;
196  for (int c = 0; c < numComponents; c++)
197  {
198  // Find the min and max values in the array.
199  a = array + c;
200  min = max = *a;
201  a += numComponents;
202  for (t = 1; t < numTuples; t++, a += numComponents)
203  {
204  if (min > *a)
205  min = *a;
206  if (max < *a)
207  max = *a;
208  }
209 
210  // Bias everything so that 0 is really the minimum.
211  if (min != 0)
212  {
213  for (t = 0, a = array + c; t < numTuples; t++, a += numComponents)
214  {
215  *a -= min;
216  }
217  }
218 
219  // Scale the values.
220  if (max != min)
221  {
222  for (t = 0, a = array + c; t < numTuples; t++, a += numComponents)
223  {
224  *a = ((maxValue - minValue) * (*a)) / (max - min);
225  }
226  }
227 
228  // Bias everything again so that it lies in the correct range.
229  if (minValue != 0)
230  {
231  for (t = 0, a = array + c; t < numTuples; t++, a += numComponents)
232  {
233  *a += minValue;
234  }
235  }
236  if (c == 0)
237  {
238  *dataMinValue = min;
239  *dataMaxValue = max;
240  }
241  }
242 }
243 
244 //-----------------------------------------------------------------------------
245 
246 template <class T>
248  T* array, int numComponents, vtkIdType numTuples, T minValue, T maxValue, double min, double max)
249 {
250  T* a;
251 
252  vtkIdType t;
253  for (int c = 0; c < numComponents; c++)
254  {
255  // Bias everything so that 0 is really the minimum.
256  if (min != 0)
257  {
258  for (t = 0, a = array + c; t < numTuples; t++, a += numComponents)
259  {
260  *a -= static_cast<T>(min);
261  }
262  }
263 
264  // Scale the values.
265  if (max != min)
266  {
267  for (t = 0, a = array + c; t < numTuples; t++, a += numComponents)
268  {
269  *a = static_cast<T>(((maxValue - minValue) * (*a)) / (max - min));
270  }
271  }
272 
273  // Bias everything again so that it lies in the correct range.
274  if (minValue != 0)
275  {
276  for (t = 0, a = array + c; t < numTuples; t++, a += numComponents)
277  {
278  *a += minValue;
279  }
280  }
281  }
282 }
283 
284 VTK_ABI_NAMESPACE_END
285 #endif // vtkFastSplatter_h
Proxy object to connect input/output ports.
A splatter optimized for splatting single kernels.
int FillInputPortInformation(int port, vtkInformation *info) override
These method should be reimplemented by subclasses that have more than a single input or single outpu...
int RequestUpdateExtent(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
Subclasses can reimplement this method to translate the update extent requests from each output port ...
vtkImageData * Buckets
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void SetLimitModeToFreezeScale()
Set/get the way voxel values will be limited.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called in response to a REQUEST_DATA request from the executive.
void SetLimitModeToNone()
Set/get the way voxel values will be limited.
void SetLimitModeToClamp()
Set/get the way voxel values will be limited.
static vtkFastSplatter * New()
void SetSplatConnection(vtkAlgorithmOutput *)
Convenience function for connecting the splat algorithm source.
void SetLimitModeToScale()
Set/get the way voxel values will be limited.
~vtkFastSplatter() override
int RequestInformation(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
Subclasses can reimplement this method to collect information from their inputs and set information f...
Generic algorithm superclass for image algs.
topologically and geometrically regular array of data
Definition: vtkImageData.h:156
a simple class to control print indentation
Definition: vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
@ info
Definition: vtkX3D.h:376
@ port
Definition: vtkX3D.h:447
void vtkFastSplatterClamp(T *array, vtkIdType arraySize, T minValue, T maxValue)
void vtkFastSplatterFrozenScale(T *array, int numComponents, vtkIdType numTuples, T minValue, T maxValue, double min, double max)
void vtkFastSplatterScale(T *array, int numComponents, vtkIdType numTuples, T minValue, T maxValue, double *dataMinValue, double *dataMaxValue)
int vtkIdType
Definition: vtkType.h:315
#define max(a, b)