VTK  9.3.0
vtkFitToHeightMapFilter.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
63 #ifndef vtkFitToHeightMapFilter_h
64 #define vtkFitToHeightMapFilter_h
65 
66 #include "vtkFiltersModelingModule.h" // For export macro
67 #include "vtkPolyDataAlgorithm.h"
68 
69 VTK_ABI_NAMESPACE_BEGIN
70 class vtkImageData;
71 
72 class VTKFILTERSMODELING_EXPORT vtkFitToHeightMapFilter : public vtkPolyDataAlgorithm
73 {
74 public:
76 
81  void PrintSelf(ostream& os, vtkIndent indent) override;
83 
91 
93 
98 
100 
106 
107  // Strategies to fit the polydata.
109  {
110  POINT_PROJECTION = 0,
111  POINT_MINIMUM_HEIGHT = 1,
112  POINT_MAXIMUM_HEIGHT = 2,
113  POINT_AVERAGE_HEIGHT = 3,
114  CELL_MINIMUM_HEIGHT = 4,
115  CELL_MAXIMUM_HEIGHT = 5,
116  CELL_AVERAGE_HEIGHT = 6,
117  };
118 
120 
132  vtkSetMacro(FittingStrategy, int);
133  vtkGetMacro(FittingStrategy, int);
134  void SetFittingStrategyToPointProjection() { this->SetFittingStrategy(POINT_PROJECTION); }
135  void SetFittingStrategyToPointMinimumHeight() { this->SetFittingStrategy(POINT_MINIMUM_HEIGHT); }
136  void SetFittingStrategyToPointMaximumHeight() { this->SetFittingStrategy(POINT_MAXIMUM_HEIGHT); }
137  void SetFittingStrategyToAverageHeight() { this->SetFittingStrategy(POINT_AVERAGE_HEIGHT); }
138  void SetFittingStrategyToCellMinimumHeight() { this->SetFittingStrategy(CELL_MINIMUM_HEIGHT); }
139  void SetFittingStrategyToCellMaximumHeight() { this->SetFittingStrategy(CELL_MAXIMUM_HEIGHT); }
140  void SetFittingStrategyToCellAverageHeight() { this->SetFittingStrategy(CELL_AVERAGE_HEIGHT); }
142 
144 
150  vtkSetMacro(UseHeightMapOffset, vtkTypeBool);
151  vtkGetMacro(UseHeightMapOffset, vtkTypeBool);
152  vtkBooleanMacro(UseHeightMapOffset, vtkTypeBool);
154 
155 protected:
158 
161 
164  double Offset;
165 
166  void AdjustPoints(vtkPolyData* output, vtkIdType numCells, vtkPoints* newPts);
168  vtkPolyData* output, vtkIdType numCells, double* cellHts, vtkPoints* inPts, vtkPoints* newPts);
169 
170 private:
172  void operator=(const vtkFitToHeightMapFilter&) = delete;
173 };
174 
175 VTK_ABI_NAMESPACE_END
176 #endif
Proxy object to connect input/output ports.
adjust polydata to fit image height map
void SetFittingStrategyToPointProjection()
Specify a strategy for fitting, or projecting, the polydata to the height field.
int FillInputPortInformation(int, vtkInformation *) override
Fill the input port information objects for this algorithm.
vtkImageData * GetHeightMap()
Get a pointer to the height map.
void SetFittingStrategyToCellMinimumHeight()
Specify a strategy for fitting, or projecting, the polydata to the height field.
~vtkFitToHeightMapFilter() override
void SetFittingStrategyToCellAverageHeight()
Specify a strategy for fitting, or projecting, the polydata to the height field.
void AdjustPoints(vtkPolyData *output, vtkIdType numCells, vtkPoints *newPts)
void SetFittingStrategyToCellMaximumHeight()
Specify a strategy for fitting, or projecting, the polydata to the height field.
static vtkFitToHeightMapFilter * New()
Standard methods for construction, type and printing.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods for construction, type and printing.
void SetHeightMapConnection(vtkAlgorithmOutput *algOutput)
Specify the pipeline connection to the height map.
void SetFittingStrategyToPointMinimumHeight()
Specify a strategy for fitting, or projecting, the polydata to the height field.
void AdjustCells(vtkPolyData *output, vtkIdType numCells, double *cellHts, vtkPoints *inPts, vtkPoints *newPts)
vtkImageData * GetHeightMap(vtkInformationVector *sourceInfo)
Get a pointer to the height map.
void SetHeightMapData(vtkImageData *idata)
Set the height map for the filter.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
void SetFittingStrategyToPointMaximumHeight()
Specify a strategy for fitting, or projecting, the polydata to the height field.
void SetFittingStrategyToAverageHeight()
Specify a strategy for fitting, or projecting, the polydata to the height field.
topologically and geometrically regular array of data
Definition: vtkImageData.h:52
a simple class to control print indentation
Definition: vtkIndent.h:38
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
represent and manipulate 3D points
Definition: vtkPoints.h:38
Superclass for algorithms that produce only polydata as output.
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:89
int vtkTypeBool
Definition: vtkABI.h:64
int vtkIdType
Definition: vtkType.h:315