VTK  9.3.20240418
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
73 #ifndef vtkFitToHeightMapFilter_h
74 #define vtkFitToHeightMapFilter_h
75 
76 #include "vtkFiltersModelingModule.h" // For export macro
77 #include "vtkPolyDataAlgorithm.h"
78 
79 VTK_ABI_NAMESPACE_BEGIN
80 class vtkImageData;
81 
82 class VTKFILTERSMODELING_EXPORT vtkFitToHeightMapFilter : public vtkPolyDataAlgorithm
83 {
84 public:
86 
91  void PrintSelf(ostream& os, vtkIndent indent) override;
93 
101 
103 
108 
110 
116 
117  // Strategies to fit the polydata.
119  {
120  POINT_PROJECTION = 0,
121  POINT_MINIMUM_HEIGHT = 1,
122  POINT_MAXIMUM_HEIGHT = 2,
123  POINT_AVERAGE_HEIGHT = 3,
124  CELL_MINIMUM_HEIGHT = 4,
125  CELL_MAXIMUM_HEIGHT = 5,
126  CELL_AVERAGE_HEIGHT = 6,
127  };
128 
130 
142  vtkSetMacro(FittingStrategy, int);
143  vtkGetMacro(FittingStrategy, int);
144  void SetFittingStrategyToPointProjection() { this->SetFittingStrategy(POINT_PROJECTION); }
145  void SetFittingStrategyToPointMinimumHeight() { this->SetFittingStrategy(POINT_MINIMUM_HEIGHT); }
146  void SetFittingStrategyToPointMaximumHeight() { this->SetFittingStrategy(POINT_MAXIMUM_HEIGHT); }
147  void SetFittingStrategyToAverageHeight() { this->SetFittingStrategy(POINT_AVERAGE_HEIGHT); }
148  void SetFittingStrategyToCellMinimumHeight() { this->SetFittingStrategy(CELL_MINIMUM_HEIGHT); }
149  void SetFittingStrategyToCellMaximumHeight() { this->SetFittingStrategy(CELL_MAXIMUM_HEIGHT); }
150  void SetFittingStrategyToCellAverageHeight() { this->SetFittingStrategy(CELL_AVERAGE_HEIGHT); }
152 
154 
160  vtkSetMacro(UseHeightMapOffset, vtkTypeBool);
161  vtkGetMacro(UseHeightMapOffset, vtkTypeBool);
162  vtkBooleanMacro(UseHeightMapOffset, vtkTypeBool);
164 
165 protected:
168 
171 
174  double Offset;
175 
176  void AdjustPoints(vtkPolyData* output, vtkIdType numCells, vtkPoints* newPts);
178  vtkPolyData* output, vtkIdType numCells, double* cellHts, vtkPoints* inPts, vtkPoints* newPts);
179 
180 private:
182  void operator=(const vtkFitToHeightMapFilter&) = delete;
183 };
184 
185 VTK_ABI_NAMESPACE_END
186 #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:156
a simple class to control print indentation
Definition: vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
represent and manipulate 3D points
Definition: vtkPoints.h:139
Superclass for algorithms that produce only polydata as output.
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:181
int vtkTypeBool
Definition: vtkABI.h:64
int vtkIdType
Definition: vtkType.h:315