InterpolateTerrain
VTKEx/Cxx/PolyData/InterpolateTerrain
Description¶
This example samples a "terrain map" using two approaches.
- Creates an image from the x,y points of the terrain and then uses vtkProbeFilter to interpolate the heights.
- Uses vtkCellLocator to directly interpolate the terrain polydata.
Note that the results differ when the point is not one of the original terrain points. This is because the image has quadrilateral elements and the polydata has triangles. As the resolution of the grid increases, the two results converge.
Question
If you have a simple question about this example contact us at VTKExProject If your question is more complex and may require extended discussion, please use the VTK Discourse Forum
Code¶
InterpolateTerrain.cxx
#include <vtkSmartPointer.h>
#include <vtkCellArray.h>
#include <vtkPoints.h>
#include <vtkTriangle.h>
#include <vtkPolyData.h>
#include <vtkPointData.h>
#include <vtkLine.h>
#include <vtkImageData.h>
#include <vtkProbeFilter.h>
#include <vtkDelaunay2D.h>
#include <vtkXMLPolyDataWriter.h>
#include <vtkDoubleArray.h>
#include <vtkMath.h>
#include <vtkCellLocator.h>
int main(int, char *[])
{
vtkSmartPointer<vtkImageData> image =
vtkSmartPointer<vtkImageData>::New();
image->SetExtent(0, 9, 0, 9, 0, 0);
image->AllocateScalars(VTK_DOUBLE,1);
//Create a random set of heights on a grid. This is often called a
//"terrain map"
vtkSmartPointer<vtkPoints> points =
vtkSmartPointer<vtkPoints>::New();
unsigned int GridSize = 10;
for ( unsigned int x = 0; x < GridSize; x++ )
{
for ( unsigned int y = 0; y < GridSize; y++ )
{
double val = vtkMath::Random(-1.0, 1.0);
points->InsertNextPoint ( x, y, val);
image->SetScalarComponentFromDouble(x,y,0,0,val);
}
}
//add the grid points to a polydata object
vtkSmartPointer<vtkPolyData> polydata =
vtkSmartPointer<vtkPolyData>::New();
polydata->SetPoints ( points );
//triangulate the grid points
vtkSmartPointer<vtkDelaunay2D> delaunay =
vtkSmartPointer<vtkDelaunay2D>::New();
delaunay->SetInputData ( polydata );
delaunay->Update();
vtkSmartPointer<vtkXMLPolyDataWriter> writer =
vtkSmartPointer<vtkXMLPolyDataWriter>::New();
writer->SetFileName ( "surface.vtp" );
writer->SetInputConnection ( delaunay->GetOutputPort() );
writer->Write();
// Add some points to interpolate
vtkSmartPointer<vtkPoints> probePoints =
vtkSmartPointer<vtkPoints>::New();
probePoints->InsertNextPoint(5.2, 3.2, 0);
probePoints->InsertNextPoint(5.0, 3.0, 0);
probePoints->InsertNextPoint(0.0, 0.0, 0);
vtkSmartPointer<vtkPolyData> probePolyData =
vtkSmartPointer<vtkPolyData>::New();
probePolyData->SetPoints(probePoints);
vtkSmartPointer<vtkProbeFilter> probe =
vtkSmartPointer<vtkProbeFilter>::New();
probe->SetSourceData(image);
probe->SetInputData(probePolyData);
probe->Update();
vtkDataArray* data = probe->GetOutput()->GetPointData()->GetScalars();
vtkDoubleArray* doubleData = dynamic_cast<vtkDoubleArray*> (data);
for(int i = 0; i < doubleData->GetNumberOfTuples(); i++)
{
double val = doubleData->GetValue(i);
cout << "Interpolation using ProbeFilter ";
cout << "doubleData->GetValue(" << i << "): " << val << endl;
}
// Now find the elevation with a CellLocator
vtkSmartPointer<vtkCellLocator> cellLocator =
vtkSmartPointer<vtkCellLocator>::New();
cellLocator->SetDataSet(delaunay->GetOutput());
cellLocator->BuildLocator();
for(int i = 0; i < doubleData->GetNumberOfTuples(); i++)
{
int subId;
double t, xyz[3], pcoords[3];
double rayStart[3], rayEnd[3];
probePoints->GetPoint(i, rayStart);
rayStart[2] += 1000.0;
probePoints->GetPoint(i, rayEnd);
rayEnd[2] -= 1000.0;
if (cellLocator->IntersectWithLine(
rayStart,
rayEnd,
0.0001,
t,
xyz,
pcoords,
subId))
{
cout << "Interpolation using CellLocator ";
cout << "Elevation at " << rayStart[0] << ", " << rayStart[1] << " is " << xyz[2] << endl;
}
}
return EXIT_SUCCESS;
}
CMakeLists.txt¶
cmake_minimum_required(VERSION 3.3 FATAL_ERROR)
project(InterpolateTerrain)
find_package(VTK COMPONENTS
vtkvtkCommonCore
vtkvtkCommonDataModel
vtkvtkFiltersCore
vtkvtkIOXML QUIET)
if (NOT VTK_FOUND)
message("Skipping InterpolateTerrain: ${VTK_NOT_FOUND_MESSAGE}")
return ()
endif()
message (STATUS "VTK_VERSION: ${VTK_VERSION}")
if (VTK_VERSION VERSION_LESS "8.90.0")
# old system
include(${VTK_USE_FILE})
add_executable(InterpolateTerrain MACOSX_BUNDLE InterpolateTerrain.cxx )
target_link_libraries(InterpolateTerrain PRIVATE ${VTK_LIBRARIES})
else ()
# include all components
add_executable(InterpolateTerrain MACOSX_BUNDLE InterpolateTerrain.cxx )
target_link_libraries(InterpolateTerrain PRIVATE ${VTK_LIBRARIES})
# vtk_module_autoinit is needed
vtk_module_autoinit(
TARGETS InterpolateTerrain
MODULES ${VTK_LIBRARIES}
)
endif ()
Download and Build InterpolateTerrain¶
Click here to download InterpolateTerrain and its CMakeLists.txt file. Once the tarball InterpolateTerrain.tar has been downloaded and extracted,
cd InterpolateTerrain/build
If VTK is installed:
cmake ..
If VTK is not installed but compiled on your system, you will need to specify the path to your VTK build:
cmake -DVTK_DIR:PATH=/home/me/vtk_build ..
Build the project:
make
and run it:
./InterpolateTerrain
WINDOWS USERS
Be sure to add the VTK bin directory to your path. This will resolve the VTK dll's at run time.