VTK  9.2.6
vtkFastSplatter.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkFastSplatter.h
5 
6  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7  All rights reserved.
8  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9 
10  This software is distributed WITHOUT ANY WARRANTY; without even
11  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12  PURPOSE. See the above copyright notice for more information.
13 
14 =========================================================================*/
15 /*----------------------------------------------------------------------------
16  Copyright (c) Sandia Corporation
17  See Copyright.txt or http://www.paraview.org/HTML/Copyright.html for details.
18 ----------------------------------------------------------------------------*/
67 #ifndef vtkFastSplatter_h
68 #define vtkFastSplatter_h
69 
70 #include "vtkImageAlgorithm.h"
71 #include "vtkImagingHybridModule.h" // For export macro
72 
73 class VTKIMAGINGHYBRID_EXPORT vtkFastSplatter : public vtkImageAlgorithm
74 {
75 public:
77  static vtkFastSplatter* New();
78  void PrintSelf(ostream& os, vtkIndent indent) override;
79 
81 
87  vtkSetVector6Macro(ModelBounds, double);
88  vtkGetVectorMacro(ModelBounds, double, 6);
90 
92 
95  vtkSetVector3Macro(OutputDimensions, int);
96  vtkGetVector3Macro(OutputDimensions, int);
98 
99  enum
100  {
104  FreezeScaleLimit
105  };
106 
108 
114  vtkSetMacro(LimitMode, int);
115  vtkGetMacro(LimitMode, int);
116  void SetLimitModeToNone() { this->SetLimitMode(NoneLimit); }
117  void SetLimitModeToClamp() { this->SetLimitMode(ClampLimit); }
118  void SetLimitModeToScale() { this->SetLimitMode(ScaleLimit); }
119  void SetLimitModeToFreezeScale() { this->SetLimitMode(FreezeScaleLimit); }
121 
123 
126  vtkSetMacro(MinValue, double);
127  vtkGetMacro(MinValue, double);
128  vtkSetMacro(MaxValue, double);
129  vtkGetMacro(MaxValue, double);
131 
133 
137  vtkGetMacro(NumberOfPointsSplatted, int);
139 
146 
147 protected:
149  ~vtkFastSplatter() override;
150 
151  double ModelBounds[6];
152  int OutputDimensions[3];
153 
155  double MinValue;
156  double MaxValue;
157  double FrozenScale;
158 
160 
165 
166  // Used internally for converting points in world space to indices in
167  // the output image.
168  double Origin[3];
169  double Spacing[3];
170 
171  // This is updated every time the filter executes
173 
174  // Used internally to track the data range. When the limit mode is
175  // set to FreezeScale, the data will be scaled as if this were the
176  // range regardless of what it actually is.
179 
180 private:
181  vtkFastSplatter(const vtkFastSplatter&) = delete;
182  void operator=(const vtkFastSplatter&) = delete;
183 };
184 
185 //-----------------------------------------------------------------------------
186 
187 template <class T>
188 void vtkFastSplatterClamp(T* array, vtkIdType arraySize, T minValue, T maxValue)
189 {
190  for (vtkIdType i = 0; i < arraySize; i++)
191  {
192  if (array[i] < minValue)
193  array[i] = minValue;
194  if (array[i] > maxValue)
195  array[i] = maxValue;
196  }
197 }
198 
199 //-----------------------------------------------------------------------------
200 
201 template <class T>
202 void vtkFastSplatterScale(T* array, int numComponents, vtkIdType numTuples, T minValue, T maxValue,
203  double* dataMinValue, double* dataMaxValue)
204 {
205  T* a;
206  T min, max;
207  *dataMinValue = 0;
208  *dataMaxValue = 0;
209  vtkIdType t;
210  for (int c = 0; c < numComponents; c++)
211  {
212  // Find the min and max values in the array.
213  a = array + c;
214  min = max = *a;
215  a += numComponents;
216  for (t = 1; t < numTuples; t++, a += numComponents)
217  {
218  if (min > *a)
219  min = *a;
220  if (max < *a)
221  max = *a;
222  }
223 
224  // Bias everything so that 0 is really the minimum.
225  if (min != 0)
226  {
227  for (t = 0, a = array + c; t < numTuples; t++, a += numComponents)
228  {
229  *a -= min;
230  }
231  }
232 
233  // Scale the values.
234  if (max != min)
235  {
236  for (t = 0, a = array + c; t < numTuples; t++, a += numComponents)
237  {
238  *a = ((maxValue - minValue) * (*a)) / (max - min);
239  }
240  }
241 
242  // Bias everything again so that it lies in the correct range.
243  if (minValue != 0)
244  {
245  for (t = 0, a = array + c; t < numTuples; t++, a += numComponents)
246  {
247  *a += minValue;
248  }
249  }
250  if (c == 0)
251  {
252  *dataMinValue = min;
253  *dataMaxValue = max;
254  }
255  }
256 }
257 
258 //-----------------------------------------------------------------------------
259 
260 template <class T>
262  T* array, int numComponents, vtkIdType numTuples, T minValue, T maxValue, double min, double max)
263 {
264  T* a;
265 
266  vtkIdType t;
267  for (int c = 0; c < numComponents; c++)
268  {
269  // Bias everything so that 0 is really the minimum.
270  if (min != 0)
271  {
272  for (t = 0, a = array + c; t < numTuples; t++, a += numComponents)
273  {
274  *a -= static_cast<T>(min);
275  }
276  }
277 
278  // Scale the values.
279  if (max != min)
280  {
281  for (t = 0, a = array + c; t < numTuples; t++, a += numComponents)
282  {
283  *a = static_cast<T>(((maxValue - minValue) * (*a)) / (max - min));
284  }
285  }
286 
287  // Bias everything again so that it lies in the correct range.
288  if (minValue != 0)
289  {
290  for (t = 0, a = array + c; t < numTuples; t++, a += numComponents)
291  {
292  *a += minValue;
293  }
294  }
295  }
296 }
297 
298 #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:163
a simple class to control print indentation
Definition: vtkIndent.h:119
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
@ info
Definition: vtkX3D.h:382
@ port
Definition: vtkX3D.h:453
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:332
#define max(a, b)