VTK  9.1.0
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 ----------------------------------------------------------------------------*/
64 #ifndef vtkFastSplatter_h
65 #define vtkFastSplatter_h
66 
67 #include "vtkImageAlgorithm.h"
68 #include "vtkImagingHybridModule.h" // For export macro
69 
70 class VTKIMAGINGHYBRID_EXPORT vtkFastSplatter : public vtkImageAlgorithm
71 {
72 public:
74  static vtkFastSplatter* New();
75  void PrintSelf(ostream& os, vtkIndent indent) override;
76 
78 
84  vtkSetVector6Macro(ModelBounds, double);
85  vtkGetVectorMacro(ModelBounds, double, 6);
87 
89 
92  vtkSetVector3Macro(OutputDimensions, int);
93  vtkGetVector3Macro(OutputDimensions, int);
95 
96  enum
97  {
101  FreezeScaleLimit
102  };
103 
105 
111  vtkSetMacro(LimitMode, int);
112  vtkGetMacro(LimitMode, int);
113  void SetLimitModeToNone() { this->SetLimitMode(NoneLimit); }
114  void SetLimitModeToClamp() { this->SetLimitMode(ClampLimit); }
115  void SetLimitModeToScale() { this->SetLimitMode(ScaleLimit); }
116  void SetLimitModeToFreezeScale() { this->SetLimitMode(FreezeScaleLimit); }
118 
120 
123  vtkSetMacro(MinValue, double);
124  vtkGetMacro(MinValue, double);
125  vtkSetMacro(MaxValue, double);
126  vtkGetMacro(MaxValue, double);
128 
130 
134  vtkGetMacro(NumberOfPointsSplatted, int);
136 
143 
144 protected:
146  ~vtkFastSplatter() override;
147 
148  double ModelBounds[6];
149  int OutputDimensions[3];
150 
152  double MinValue;
153  double MaxValue;
154  double FrozenScale;
155 
157 
162 
163  // Used internally for converting points in world space to indices in
164  // the output image.
165  double Origin[3];
166  double Spacing[3];
167 
168  // This is updated every time the filter executes
170 
171  // Used internally to track the data range. When the limit mode is
172  // set to FreezeScale, the data will be scaled as if this were the
173  // range regardless of what it actually is.
176 
177 private:
178  vtkFastSplatter(const vtkFastSplatter&) = delete;
179  void operator=(const vtkFastSplatter&) = delete;
180 };
181 
182 //-----------------------------------------------------------------------------
183 
184 template <class T>
185 void vtkFastSplatterClamp(T* array, vtkIdType arraySize, T minValue, T maxValue)
186 {
187  for (vtkIdType i = 0; i < arraySize; i++)
188  {
189  if (array[i] < minValue)
190  array[i] = minValue;
191  if (array[i] > maxValue)
192  array[i] = maxValue;
193  }
194 }
195 
196 //-----------------------------------------------------------------------------
197 
198 template <class T>
199 void vtkFastSplatterScale(T* array, int numComponents, vtkIdType numTuples, T minValue, T maxValue,
200  double* dataMinValue, double* dataMaxValue)
201 {
202  T* a;
203  T min, max;
204  *dataMinValue = 0;
205  *dataMaxValue = 0;
206  vtkIdType t;
207  for (int c = 0; c < numComponents; c++)
208  {
209  // Find the min and max values in the array.
210  a = array + c;
211  min = max = *a;
212  a += numComponents;
213  for (t = 1; t < numTuples; t++, a += numComponents)
214  {
215  if (min > *a)
216  min = *a;
217  if (max < *a)
218  max = *a;
219  }
220 
221  // Bias everything so that 0 is really the minimum.
222  if (min != 0)
223  {
224  for (t = 0, a = array + c; t < numTuples; t++, a += numComponents)
225  {
226  *a -= min;
227  }
228  }
229 
230  // Scale the values.
231  if (max != min)
232  {
233  for (t = 0, a = array + c; t < numTuples; t++, a += numComponents)
234  {
235  *a = ((maxValue - minValue) * (*a)) / (max - min);
236  }
237  }
238 
239  // Bias everything again so that it lies in the correct range.
240  if (minValue != 0)
241  {
242  for (t = 0, a = array + c; t < numTuples; t++, a += numComponents)
243  {
244  *a += minValue;
245  }
246  }
247  if (c == 0)
248  {
249  *dataMinValue = min;
250  *dataMaxValue = max;
251  }
252  }
253 }
254 
255 //-----------------------------------------------------------------------------
256 
257 template <class T>
259  T* array, int numComponents, vtkIdType numTuples, T minValue, T maxValue, double min, double max)
260 {
261  T* a;
262 
263  vtkIdType t;
264  for (int c = 0; c < numComponents; c++)
265  {
266  // Bias everything so that 0 is really the minimum.
267  if (min != 0)
268  {
269  for (t = 0, a = array + c; t < numTuples; t++, a += numComponents)
270  {
271  *a -= static_cast<T>(min);
272  }
273  }
274 
275  // Scale the values.
276  if (max != min)
277  {
278  for (t = 0, a = array + c; t < numTuples; t++, a += numComponents)
279  {
280  *a = static_cast<T>(((maxValue - minValue) * (*a)) / (max - min));
281  }
282  }
283 
284  // Bias everything again so that it lies in the correct range.
285  if (minValue != 0)
286  {
287  for (t = 0, a = array + c; t < numTuples; t++, a += numComponents)
288  {
289  *a += minValue;
290  }
291  }
292  }
293 }
294 
295 #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:157
a simple class to control print indentation
Definition: vtkIndent.h:113
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)