VTK  9.3.0
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
33 #ifndef vtkFastSplatter_h
34 #define vtkFastSplatter_h
35 
36 #include "vtkImageAlgorithm.h"
37 #include "vtkImagingHybridModule.h" // For export macro
38 
39 VTK_ABI_NAMESPACE_BEGIN
40 class VTKIMAGINGHYBRID_EXPORT vtkFastSplatter : public vtkImageAlgorithm
41 {
42 public:
44  static vtkFastSplatter* New();
45  void PrintSelf(ostream& os, vtkIndent indent) override;
46 
48 
54  vtkSetVector6Macro(ModelBounds, double);
55  vtkGetVectorMacro(ModelBounds, double, 6);
57 
59 
62  vtkSetVector3Macro(OutputDimensions, int);
63  vtkGetVector3Macro(OutputDimensions, int);
65 
66  enum
67  {
71  FreezeScaleLimit
72  };
73 
75 
81  vtkSetMacro(LimitMode, int);
82  vtkGetMacro(LimitMode, int);
83  void SetLimitModeToNone() { this->SetLimitMode(NoneLimit); }
84  void SetLimitModeToClamp() { this->SetLimitMode(ClampLimit); }
85  void SetLimitModeToScale() { this->SetLimitMode(ScaleLimit); }
86  void SetLimitModeToFreezeScale() { this->SetLimitMode(FreezeScaleLimit); }
88 
90 
93  vtkSetMacro(MinValue, double);
94  vtkGetMacro(MinValue, double);
95  vtkSetMacro(MaxValue, double);
96  vtkGetMacro(MaxValue, double);
98 
100 
104  vtkGetMacro(NumberOfPointsSplatted, int);
106 
113 
114 protected:
116  ~vtkFastSplatter() override;
117 
118  double ModelBounds[6];
119  int OutputDimensions[3];
120 
122  double MinValue;
123  double MaxValue;
124  double FrozenScale;
125 
127 
132 
133  // Used internally for converting points in world space to indices in
134  // the output image.
135  double Origin[3];
136  double Spacing[3];
137 
138  // This is updated every time the filter executes
140 
141  // Used internally to track the data range. When the limit mode is
142  // set to FreezeScale, the data will be scaled as if this were the
143  // range regardless of what it actually is.
146 
147 private:
148  vtkFastSplatter(const vtkFastSplatter&) = delete;
149  void operator=(const vtkFastSplatter&) = delete;
150 };
151 
152 //-----------------------------------------------------------------------------
153 
154 template <class T>
155 void vtkFastSplatterClamp(T* array, vtkIdType arraySize, T minValue, T maxValue)
156 {
157  for (vtkIdType i = 0; i < arraySize; i++)
158  {
159  if (array[i] < minValue)
160  array[i] = minValue;
161  if (array[i] > maxValue)
162  array[i] = maxValue;
163  }
164 }
165 
166 //-----------------------------------------------------------------------------
167 
168 template <class T>
169 void vtkFastSplatterScale(T* array, int numComponents, vtkIdType numTuples, T minValue, T maxValue,
170  double* dataMinValue, double* dataMaxValue)
171 {
172  T* a;
173  T min, max;
174  *dataMinValue = 0;
175  *dataMaxValue = 0;
176  vtkIdType t;
177  for (int c = 0; c < numComponents; c++)
178  {
179  // Find the min and max values in the array.
180  a = array + c;
181  min = max = *a;
182  a += numComponents;
183  for (t = 1; t < numTuples; t++, a += numComponents)
184  {
185  if (min > *a)
186  min = *a;
187  if (max < *a)
188  max = *a;
189  }
190 
191  // Bias everything so that 0 is really the minimum.
192  if (min != 0)
193  {
194  for (t = 0, a = array + c; t < numTuples; t++, a += numComponents)
195  {
196  *a -= min;
197  }
198  }
199 
200  // Scale the values.
201  if (max != min)
202  {
203  for (t = 0, a = array + c; t < numTuples; t++, a += numComponents)
204  {
205  *a = ((maxValue - minValue) * (*a)) / (max - min);
206  }
207  }
208 
209  // Bias everything again so that it lies in the correct range.
210  if (minValue != 0)
211  {
212  for (t = 0, a = array + c; t < numTuples; t++, a += numComponents)
213  {
214  *a += minValue;
215  }
216  }
217  if (c == 0)
218  {
219  *dataMinValue = min;
220  *dataMaxValue = max;
221  }
222  }
223 }
224 
225 //-----------------------------------------------------------------------------
226 
227 template <class T>
229  T* array, int numComponents, vtkIdType numTuples, T minValue, T maxValue, double min, double max)
230 {
231  T* a;
232 
233  vtkIdType t;
234  for (int c = 0; c < numComponents; c++)
235  {
236  // Bias everything so that 0 is really the minimum.
237  if (min != 0)
238  {
239  for (t = 0, a = array + c; t < numTuples; t++, a += numComponents)
240  {
241  *a -= static_cast<T>(min);
242  }
243  }
244 
245  // Scale the values.
246  if (max != min)
247  {
248  for (t = 0, a = array + c; t < numTuples; t++, a += numComponents)
249  {
250  *a = static_cast<T>(((maxValue - minValue) * (*a)) / (max - min));
251  }
252  }
253 
254  // Bias everything again so that it lies in the correct range.
255  if (minValue != 0)
256  {
257  for (t = 0, a = array + c; t < numTuples; t++, a += numComponents)
258  {
259  *a += minValue;
260  }
261  }
262  }
263 }
264 
265 VTK_ABI_NAMESPACE_END
266 #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:43
a simple class to control print indentation
Definition: vtkIndent.h:29
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)