ShepardMethod
Repository source: ShepardMethod
Description¶
This example samples unstructured points onto structured points using the Shepard method. The example starts with two points which have associated scalars (0 (black) and 1(white)). The results are displayed by coloring planes between the two points with the corresponding interpolated values. The values are reflected by black (0) to white (1).
Question
If you have a question about this example, please use the VTK Discourse Forum
Code¶
ShepardMethod.cxx
#include <vtkActor.h>
#include <vtkCamera.h>
#include <vtkCellArray.h>
#include <vtkColorTransferFunction.h>
#include <vtkContourFilter.h>
#include <vtkFloatArray.h>
#include <vtkNamedColors.h>
#include <vtkNew.h>
#include <vtkPointData.h>
#include <vtkPolyDataMapper.h>
#include <vtkProperty.h>
#include <vtkRenderWindow.h>
#include <vtkRenderWindowInteractor.h>
#include <vtkRenderer.h>
#include <vtkShepardMethod.h>
#include <vtkUnsignedCharArray.h>
#include <vtkVertexGlyphFilter.h>
// For compatibility with new VTK generic data arrays
#ifdef vtkGenericDataArray_h
#define InsertNextTupleValue InsertNextTypedTuple
#endif
#include <iostream>
#include <string>
int main(int, char*[])
{
vtkNew<vtkNamedColors> colors;
// Create a set of vertices (polydata)
vtkNew<vtkPoints> points;
points->InsertNextPoint(100.0, 0.0, 0.0);
points->InsertNextPoint(300.0, 0.0, 0.0);
// Setup colors
unsigned char white[3] = {255, 255, 255};
unsigned char black[3] = {0, 0, 0};
vtkNew<vtkUnsignedCharArray> vertexColors;
vertexColors->SetNumberOfComponents(3);
vertexColors->SetName("Colors");
vertexColors->InsertNextTupleValue(black);
vertexColors->InsertNextTupleValue(white);
// Create a scalar array for the pointdata, each value represents the distance
// of the vertices from the first vertex.
vtkNew<vtkFloatArray> values;
values->SetNumberOfComponents(1);
values->SetName("Values");
values->InsertNextValue(0.0);
values->InsertNextValue(1.0);
// We must make two objects, because the ShepardMethod uses the ActiveScalars,
// as does the renderer!
vtkNew<vtkPolyData> polydataToProcess;
polydataToProcess->SetPoints(points);
polydataToProcess->GetPointData()->SetScalars(values);
vtkNew<vtkPolyData> polydataToVisualize;
polydataToVisualize->SetPoints(points);
polydataToVisualize->GetPointData()->SetScalars(vertexColors);
vtkNew<vtkVertexGlyphFilter> vertexGlyphFilter;
vertexGlyphFilter->AddInputData(polydataToVisualize);
vertexGlyphFilter->Update();
// Create a mapper and actor.
vtkNew<vtkPolyDataMapper> vertsMapper;
// vertsMapper->ScalarVisibilityOff();
vertsMapper->SetInputConnection(vertexGlyphFilter->GetOutputPort());
vtkNew<vtkActor> vertsActor;
vertsActor->SetMapper(vertsMapper);
vertsActor->GetProperty()->SetColor(1, 0, 0);
vertsActor->GetProperty()->SetPointSize(3);
// Create a shepard filter to interpolate the vertices over a regularized
// image grid.
vtkNew<vtkShepardMethod> shepard;
shepard->SetInputData(polydataToProcess);
shepard->SetSampleDimensions(2, 2, 2);
shepard->SetModelBounds(100, 300, -10, 10, -10, 10);
shepard->SetMaximumDistance(1);
// Contour the shepard generated image at 3 isovalues.
// The accuracy of the results are highly dependent on how the shepard filter
// is set up.
vtkNew<vtkContourFilter> contourFilter;
contourFilter->SetNumberOfContours(3);
contourFilter->SetValue(0, 0.25);
contourFilter->SetValue(1, 0.50);
contourFilter->SetValue(2, 0.75);
contourFilter->SetInputConnection(shepard->GetOutputPort());
contourFilter->Update();
// Create a mapper and actor for the resulting isosurfaces.
vtkNew<vtkPolyDataMapper> contourMapper;
contourMapper->SetInputConnection(contourFilter->GetOutputPort());
contourMapper->ScalarVisibilityOn();
contourMapper->SetColorModeToMapScalars();
vtkNew<vtkActor> contourActor;
contourActor->SetMapper(contourMapper);
contourActor->GetProperty()->SetAmbient(1);
contourActor->GetProperty()->SetSpecular(0);
contourActor->GetProperty()->SetDiffuse(0);
// Report the results of the interpolation.
double* range = contourFilter->GetOutput()->GetScalarRange();
std::cout << "Shepard interpolation:" << std::endl;
std::cout << "Contour output scalar range: " << range[0] << ", " << range[1]
<< std::endl;
vtkIdType nCells = contourFilter->GetOutput()->GetNumberOfCells();
double bounds[6];
for (vtkIdType i = 0; i < nCells; ++i)
{
if (i %
2) // each isosurface value only has 2 cells so report on the odd ones.
{
contourFilter->GetOutput()->GetCellBounds(i, bounds);
std::cout << "Cell " << i << ", x position: " << bounds[0] << std::endl;
}
}
// Create a transfer function to color the isosurfaces.
vtkNew<vtkColorTransferFunction> lut;
lut->SetColorSpaceToRGB();
lut->AddRGBPoint(range[0], 0, 0, 0); // black
lut->AddRGBPoint(range[1], 1, 1, 1); // white
lut->SetScaleToLinear();
contourMapper->SetLookupTable(lut);
// Create a renderer, render window and interactor.
vtkNew<vtkRenderer> renderer;
renderer->GradientBackgroundOn();
renderer->AddActor(contourActor);
renderer->AddActor(vertsActor);
renderer->SetBackground(colors->GetColor3d("Blue").GetData());
renderer->SetBackground2(colors->GetColor3d("Magenta").GetData());
vtkNew<vtkRenderWindow> renderWindow;
renderWindow->AddRenderer(renderer);
renderWindow->SetWindowName("ShepardMethod");
vtkNew<vtkRenderWindowInteractor> renderWindowInteractor;
renderWindowInteractor->SetRenderWindow(renderWindow);
// Position the camera so that the image produced is viewable
vtkCamera* camera = renderer->GetActiveCamera();
camera->SetPosition(450, 100, 100);
camera->SetFocalPoint(200, 0, 0);
camera->SetViewUp(0, 0, 1);
renderWindow->Render();
renderWindowInteractor->Start();
return EXIT_SUCCESS;
}
CMakeLists.txt¶
cmake_minimum_required(VERSION 3.12 FATAL_ERROR)
project(ShepardMethod)
find_package(VTK COMPONENTS
CommonColor
CommonCore
CommonDataModel
FiltersCore
FiltersGeneral
ImagingHybrid
InteractionStyle
RenderingContextOpenGL2
RenderingCore
RenderingFreeType
RenderingGL2PSOpenGL2
RenderingOpenGL2
)
if (NOT VTK_FOUND)
message(FATAL_ERROR "ShepardMethod: Unable to find the VTK build folder.")
endif()
# Prevent a "command line is too long" failure in Windows.
set(CMAKE_NINJA_FORCE_RESPONSE_FILE "ON" CACHE BOOL "Force Ninja to use response files.")
add_executable(ShepardMethod MACOSX_BUNDLE ShepardMethod.cxx )
target_link_libraries(ShepardMethod PRIVATE ${VTK_LIBRARIES}
)
# vtk_module_autoinit is needed
vtk_module_autoinit(
TARGETS ShepardMethod
MODULES ${VTK_LIBRARIES}
)
Download and Build ShepardMethod¶
Click here to download ShepardMethod and its CMakeLists.txt file. Once the tarball ShepardMethod.tar has been downloaded and extracted,
cd ShepardMethod/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:
./ShepardMethod
WINDOWS USERS
Be sure to add the VTK bin directory to your path. This will resolve the VTK dll's at run time.