VTK  9.3.20240423
vtkProjectedTerrainPath.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
59#ifndef vtkProjectedTerrainPath_h
60#define vtkProjectedTerrainPath_h
61
62#include "vtkFiltersHybridModule.h" // For export macro
64
65VTK_ABI_NAMESPACE_BEGIN
67class vtkImageData;
68class vtkEdgeList;
69class vtkPoints;
70
71class VTKFILTERSHYBRID_EXPORT vtkProjectedTerrainPath : public vtkPolyDataAlgorithm
72{
73public:
75
79 void PrintSelf(ostream& os, vtkIndent indent) override;
81
86
88
97
103
104 enum
105 {
106 SIMPLE_PROJECTION = 0,
108 HUG_PROJECTION
109 };
110
112
120 vtkSetClampMacro(ProjectionMode, int, SIMPLE_PROJECTION, HUG_PROJECTION);
121 vtkGetMacro(ProjectionMode, int);
122 void SetProjectionModeToSimple() { this->SetProjectionMode(SIMPLE_PROJECTION); }
123 void SetProjectionModeToNonOccluded() { this->SetProjectionMode(NONOCCLUDED_PROJECTION); }
124 void SetProjectionModeToHug() { this->SetProjectionMode(HUG_PROJECTION); }
126
128
133 vtkSetMacro(HeightOffset, double);
134 vtkGetMacro(HeightOffset, double);
136
138
143 vtkSetClampMacro(HeightTolerance, double, 0.0, VTK_FLOAT_MAX);
144 vtkGetMacro(HeightTolerance, double);
146
148
153 vtkSetClampMacro(MaximumNumberOfLines, vtkIdType, 1, VTK_ID_MAX);
154 vtkGetMacro(MaximumNumberOfLines, vtkIdType);
156
157protected:
160
162 int FillInputPortInformation(int port, vtkInformation* info) override;
163
164 // Supporting methods
165 void GetImageIndex(double x[3], double loc[2], int ij[2]);
166 double GetHeight(double loc[2], int ij[2]);
170 void SplitEdge(vtkIdType eId, double t);
171
172 // ivars that the API addresses
177
178 // Bookkeeping arrays
179 int Dimensions[3];
180 int Extent[6];
181 double Origin[3];
182 double Spacing[3];
186
187 // Errors above/below terrain. In both instances, negative values are
188 // inserted because the priority queue puts smallest values on top.
189 vtkPriorityQueue* PositiveLineError; // errors above terrain
190 vtkPriorityQueue* NegativeLineError; // errors below terrain
191
192 // This is a PIMPL'd vector representing edges
193 vtkEdgeList* EdgeList;
194
195private:
197 void operator=(const vtkProjectedTerrainPath&) = delete;
198};
199
200VTK_ABI_NAMESPACE_END
201#endif
Proxy object to connect input/output ports.
abstract superclass for arrays of numeric data
topologically and geometrically regular array of data
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.
a list of ids arranged in priority order
project a polyline onto a terrain
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
~vtkProjectedTerrainPath() override
void SetSourceData(vtkImageData *source)
Specify the second input (the terrain) onto which the polyline(s) should be projected.
void SetProjectionModeToNonOccluded()
Determine how to control the projection process.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
double GetHeight(double loc[2], int ij[2])
void ComputeError(vtkIdType edgeId)
void SetProjectionModeToHug()
Determine how to control the projection process.
void SplitEdge(vtkIdType eId, double t)
static vtkProjectedTerrainPath * New()
Instantiate the class.
void SetSourceConnection(vtkAlgorithmOutput *algOutput)
Specify the second input (the terrain) onto which the polyline(s) should be projected.
void GetImageIndex(double x[3], double loc[2], int ij[2])
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods for printing and determining type information.
void SetProjectionModeToSimple()
Determine how to control the projection process.
vtkImageData * GetSource()
Specify the second input (the terrain) onto which the polyline(s) should be projected.
boost::graph_traits< vtkGraph * >::vertex_descriptor source(boost::graph_traits< vtkGraph * >::edge_descriptor e, vtkGraph *)
int vtkIdType
Definition vtkType.h:315
#define VTK_ID_MAX
Definition vtkType.h:319
#define VTK_FLOAT_MAX
Definition vtkType.h:152