CarotidFlow
vtk-examples/Cxx/VisualizationAlgorithms/CarotidFlow
Description¶
This example generates streamtubes of blood velocity. an isosurface of speed provides context. The starting positions for the streamtubes were determined by experimenting with the data. Because of the way the data was measured and the resolution of the velocity field, many streamers travel outside the artery. This is because the boundary layer of the blood flow is not captured due to limitations in data resolution. Consequently, as the blood flows around curves, there is a component of the velocity field that directs the streamtube outside the artery. As a result it is hard to find starting positions for the streamtubes that yield interesting results. The examples uses the source object vtkPointSource in combination with vtkThresholdPoints to work around this problem. vtkPointSource generates random points centered around a sphere of a specified radius. We need only find an approximate position for the starting points of the streamtubes and then generate a cloud of random seed points. vtkThresholdPoints is used to cull points that may be generated outside the regions of high flow velocity.
Cite
See 3D Phase Contrast MRI of Cerebral Blood Flow and Surface Anatomy for background.
Info
See Figure 6-44 in Chapter 6 the VTK Textbook.
Other languages
See (Python)
Question
If you have a question about this example, please use the VTK Discourse Forum
Code¶
CarotidFlow.cxx
#include <vtkActor.h>
#include <vtkCamera.h>
#include <vtkContourFilter.h>
#include <vtkLookupTable.h>
#include <vtkNamedColors.h>
#include <vtkNew.h>
#include <vtkOutlineFilter.h>
#include <vtkPointData.h>
#include <vtkPointSource.h>
#include <vtkPolyDataMapper.h>
#include <vtkProperty.h>
#include <vtkRenderWindow.h>
#include <vtkRenderWindowInteractor.h>
#include <vtkRenderer.h>
#include <vtkStreamTracer.h>
#include <vtkStructuredPointsReader.h>
#include <vtkThresholdPoints.h>
#include <vtkTubeFilter.h>
#include <iostream>
#include <string>
int main(int argc, char* argv[])
{
if (argc < 2)
{
std::cout << "Usage: " << argv[0] << " carotid.vtk" << std::endl;
return EXIT_FAILURE;
}
vtkNew<vtkNamedColors> colors;
vtkNew<vtkRenderer> ren1;
vtkNew<vtkRenderWindow> renWin;
renWin->AddRenderer(ren1);
vtkNew<vtkRenderWindowInteractor> iren;
iren->SetRenderWindow(renWin);
// create pipeline
//
vtkNew<vtkStructuredPointsReader> reader;
reader->SetFileName(argv[1]);
vtkNew<vtkPointSource> psource;
psource->SetNumberOfPoints(25);
psource->SetCenter(133.1, 116.3, 5.0);
psource->SetRadius(2.0);
vtkNew<vtkThresholdPoints> threshold;
threshold->SetInputConnection(reader->GetOutputPort());
threshold->ThresholdByUpper(275);
vtkNew<vtkStreamTracer> streamers;
streamers->SetInputConnection(reader->GetOutputPort());
streamers->SetSourceConnection(psource->GetOutputPort());
// streamers->SetMaximumPropagationUnitToTimeUnit();
streamers->SetMaximumPropagation(100.0);
// streamers->SetInitialIntegrationStepUnitToCellLengthUnit();
streamers->SetInitialIntegrationStep(0.2);
streamers->SetTerminalSpeed(.01);
streamers->Update();
double range[2];
range[0] =
streamers->GetOutput()->GetPointData()->GetScalars()->GetRange()[0];
range[1] =
streamers->GetOutput()->GetPointData()->GetScalars()->GetRange()[1];
vtkNew<vtkTubeFilter> tubes;
tubes->SetInputConnection(streamers->GetOutputPort());
tubes->SetRadius(0.3);
tubes->SetNumberOfSides(6);
tubes->SetVaryRadius(0);
vtkNew<vtkLookupTable> lut;
lut->SetHueRange(.667, 0.0);
lut->Build();
vtkNew<vtkPolyDataMapper> streamerMapper;
streamerMapper->SetInputConnection(tubes->GetOutputPort());
streamerMapper->SetScalarRange(range[0], range[1]);
streamerMapper->SetLookupTable(lut);
vtkNew<vtkActor> streamerActor;
streamerActor->SetMapper(streamerMapper);
// contours of speed
vtkNew<vtkContourFilter> iso;
iso->SetInputConnection(reader->GetOutputPort());
iso->SetValue(0, 175);
vtkNew<vtkPolyDataMapper> isoMapper;
isoMapper->SetInputConnection(iso->GetOutputPort());
isoMapper->ScalarVisibilityOff();
vtkNew<vtkActor> isoActor;
isoActor->SetMapper(isoMapper);
isoActor->GetProperty()->SetRepresentationToWireframe();
isoActor->GetProperty()->SetOpacity(0.25);
// outline
vtkNew<vtkOutlineFilter> outline;
outline->SetInputConnection(reader->GetOutputPort());
vtkNew<vtkPolyDataMapper> outlineMapper;
outlineMapper->SetInputConnection(outline->GetOutputPort());
vtkNew<vtkActor> outlineActor;
outlineActor->SetMapper(outlineMapper);
outlineActor->GetProperty()->SetColor(colors->GetColor3d("Black").GetData());
// Add the actors to the renderer, set the background and size.
//
ren1->AddActor(outlineActor);
ren1->AddActor(streamerActor);
ren1->AddActor(isoActor);
ren1->SetBackground(colors->GetColor3d("Wheat").GetData());
renWin->SetSize(640, 480);
renWin->SetWindowName("CarotidFlow");
vtkNew<vtkCamera> cam1;
cam1->SetClippingRange(17.4043, 870.216);
cam1->SetFocalPoint(136.71, 104.025, 23);
cam1->SetPosition(204.747, 258.939, 63.7925);
cam1->SetViewUp(-0.102647, -0.210897, 0.972104);
cam1->Zoom(1.2);
ren1->SetActiveCamera(cam1);
// Render the image.
//
renWin->Render();
iren->Start();
return EXIT_SUCCESS;
}
CMakeLists.txt¶
cmake_minimum_required(VERSION 3.12 FATAL_ERROR)
project(CarotidFlow)
find_package(VTK COMPONENTS
)
if (NOT VTK_FOUND)
message(FATAL_ERROR "CarotidFlow: 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(CarotidFlow MACOSX_BUNDLE CarotidFlow.cxx )
target_link_libraries(CarotidFlow PRIVATE ${VTK_LIBRARIES}
)
# vtk_module_autoinit is needed
vtk_module_autoinit(
TARGETS CarotidFlow
MODULES ${VTK_LIBRARIES}
)
Download and Build CarotidFlow¶
Click here to download CarotidFlow and its CMakeLists.txt file. Once the tarball CarotidFlow.tar has been downloaded and extracted,
cd CarotidFlow/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:
./CarotidFlow
WINDOWS USERS
Be sure to add the VTK bin directory to your path. This will resolve the VTK dll's at run time.