SimpleRTK

From Openrtk
Jump to navigation Jump to search

Introduction

SimpleRTK provides a simple interface to RTK in a variety of languages. For now only Python and C# binding have been tested and are supported but other languages should also work, e.g., R, Ruby and Java. SimpleRTK is provided as part of RTK making it easier to configure, build and use. Note that SimpleRTK is derived from the SimpleITK project for the Insight Toolkit, however it's under its own project.

Building RTK with SimpleRTK

SimpleRTK uses Swig to generate the wrapping for different languages. As the version of Swig should be controlled, SimpleRTK will automatically download the right version of Swig and its dependencies for the different platforms. However you can always point to the installed version of Swig if you'd like to, but this is not recommended.

Here are the simple steps to build and configure wrapping for RTK:

  1. Run CMake on RTK as you normally do
  2. Check the option BUILD_SIMPLERTK (default is OFF)
  3. Ename WRAP_PYTHON and/or WRAP_CSHARP (On Windows)
  4. Configure
  5. Build as usual. Note that building will download the require libraries so you should have an internet connection enabled.

Note that the configuration and generation step can take some time as CMake as to generate the necessary classes.

Python installation

To install the built python package into the system Python, as root run:

cd SimpleRTK-build/Wrapping
python PythonPackage/setup.py install

Alternatively, a Python virtual environment can be created and the distribution installed there. If you build the "dist" target a Python egg will be created in the "Wrapping/dist" directory. Building Python wheels can be enabled with a CMake flag.

Testing SimpleRTK

A simple example is located in Utilities/SimpleRTK/Examples/RTKFirstReconstruction.py and shows how to use SimpleRTK

#!/usr/bin/env python
from __future__ import print_function
import SimpleRTK as srtk
import sys
import os
import matplotlib.pyplot as plt
import matplotlib.cm as cm

if len ( sys.argv ) < 2:
    print( "Usage: RTKFirstReconstruction <output>" )
    sys.exit ( 1 )

# Defines the RTK geometry object
geometry = srtk.ThreeDimCircularProjectionGeometry()
numberOfProjections = 360
firstAngle = 0
angularArc = 360
sid = 600 # source to isocenter distance in mm
sdd = 1200 # source to detector distance in mm
isox = 0 # X coordinate on the projection image of isocenter
isoy = 0 # Y coordinate on the projection image of isocenter
for x in range(0,numberOfProjections):
  angle = firstAngle + x * angularArc / numberOfProjections
  geometry.AddProjection(sid,sdd,angle,isox,isoy)

constantImageSource = srtk.ConstantImageSource()
origin = [ -127.5, -127.5, 0. ]
sizeOutput = [ 256, 256,  numberOfProjections ]
spacing = [ 1.0, 1.0, 1.0 ]
constantImageSource.SetOrigin( origin )
constantImageSource.SetSpacing( spacing )
constantImageSource.SetSize( sizeOutput )
constantImageSource.SetConstant(0.0)
source = constantImageSource.Execute()

rei = srtk.RayEllipsoidIntersectionImageFilter()
semiprincipalaxis = [ 50, 50, 50]
center = [ 0, 0, 0]
# Set GrayScale value, axes, center...
rei.SetDensity(20)
rei.SetAngle(0)
rei.SetCenter(center)
rei.SetAxis(semiprincipalaxis)
rei.SetGeometry( geometry )
reiImage = rei.Execute(source)
  
  
# Create reconstructed image
constantImageSource2 = srtk.ConstantImageSource()
origin = [ -63.5, -63.5, -63.5 ]
sizeOutput = [ 128, 128, 128 ]
constantImageSource2.SetOrigin( origin )
constantImageSource2.SetSpacing( spacing )
constantImageSource2.SetSize( sizeOutput )
constantImageSource2.SetConstant(0.0)
source2 = constantImageSource2.Execute()

print("Performing reconstruction")
feldkamp = srtk.FDKConeBeamReconstructionFilter()
feldkamp.SetGeometry( geometry );
feldkamp.SetTruncationCorrection(0.0);
feldkamp.SetHannCutFrequency(0.0);
image = feldkamp.Execute(source2,reiImage) 
   
pixelID = image.GetPixelIDValue()

plt.imshow(srtk.GetArrayFromImage(image[:,64,:]), cmap = cm.Greys_r)

caster = srtk.CastImageFilter()
caster.SetOutputPixelType( pixelID )
image = caster.Execute( image )

writer = srtk.ImageFileWriter()
writer.SetFileName ( sys.argv[1] )
writer.Execute ( image );

The example displays the central slice of a reconstructed sphere:

Error creating thumbnail: File missing

SimpleRTK with SimpleITK

Even if SimpleRTK is based on SimpleITK, the basic object types (images, transforms) are different and there is no direct compatibility between the two toolkits. However it is fairly straightforward to pass the data from RTK to ITK using the following example:

 sitk.GetImageFromArray( srtk.GetArrayFromImage( image ) )

Extending SimpleRTK

Adding new classes to SimpleRTK should be fairly straightforward. In this section we cover how to add common types as well as filters. Note that common types should only be added when necessary.

The file locate in Wrapping/SimpleRTK.i is the main SWIG file and new types and manually wrapped files should be added there.

Common Type

One example is for the ThreeDCircularProjectionGeometry object. The code is located in utilities/SimpleRTK/Code/Common. Note that the class is named srtkThreeDimCircularProjectionGeometry to avoid some conflicting type resolution by Swig (this could be tricky). This class in the namespace rtk::simple provides a PIMPL class implementation of the object as well as a conversion from SimpleRTK to RTK itself to expose the necessary functions.

Filters

Filters are implemented using a json description for simplicity. The best way is to look into already wrapped filters in the: utilities/SimpleRTK/Code/BasicFilters/json directory.

Here we provide a description of the current JSON format.

{
  # Name of the class to be create. This would be prefixed by srtk. 
  "name" : "FDKConeBeamReconstructionFilter",

  # File that indicates the template to use. This should be RTKImageFilter for image filters
  "template_code_filename" : "RTKImageFilter",

  # File that indicates the template to use for the test.
  "template_test_filename" : "ImageFilter",

  # Number of inputs to the filter.
  "number_of_inputs" : 2,

  # Documentation for the class, if any.
  "doc" : "",

  # Type of output image. This creates a typedef OutputImageType. Here we set it to the input image.
  "output_image_type" : "TImageType",

  # List of supported Pixel types
  "pixel_types" : "RealPixelIDTypeList",

  # Definition of the RTK filter
  "filter_type" : "rtk::FDKConeBeamReconstructionFilter<InputImageType>",

  # Array of include files.
  "include_files" : [
    "srtkThreeDimCircularProjectionGeometry.h"
  ],

The next section describes the member functions. Each member function has a name and a type. The Set...() and Get...() functions are automatically generated.

  "members" : [
    {
      # Name of the variable to be Set/Get
      "name" : "Geometry",

      # Type of the variable
      "type" : "ThreeDimCircularProjectionGeometry*",

      # Default value of the variable
      "default" : "0",
      
      # Internal Simple RTK type of the function
      "itk_type" : "ThreeDimCircularProjectionGeometry",

      # Casting from the itk_type to the internal type
      "custom_itk_cast" : "filter->SetGeometry( this->m_Geometry->GetRTKBase() );\n",

      # Documentation for the member function
      "doc" : "",
      "briefdescriptionSet" : "",
      "detaileddescriptionSet" : "Set the object pointer to projection geometry.",
      "briefdescriptionGet" : "",
      "detaileddescriptionGet" : "Get the object pointer to projection geometry."
    }
  ],

At the end of the file we add the description of the class for documentation purposes.

  "briefdescription" : "Brief description of the filter",
  "detaileddescription" : "Detailed description of the filter"
  }

Note that you might need to rebuild the solution in order for CMake to take the changes into account. The code generated is located in the "SimpleRTK-build/Code/BasicFilters" of your build tree.