Feature Extraction Methods
Choose your preferred method for extracting radiomic features from Python:
Installation
Juliacall acts as a bridge between Python and Julia, enabling you to call Julia functions directly from Python code.
Installing Juliacall
Install juliacall via pip:
pip install juliacall
Note: If Julia is not already installed on your system, juliacall will download and use its own embedded version automatically.
Installing Radiomics.jl
To install the Radiomics.jl package, run the following in Python:
from juliacall import Main as jl
jl.seval('import Pkg; Pkg.add("Radiomics")')
Also, to update Radiomics.jl package, run:
from juliacall import Main as jl
jl.seval('import Pkg; Pkg.update("Radiomics")')
Multi Threading
If you want activate threads in Radiomics.jl from Python by using Juliacall, add also this line in .bashrc file:
export PYTHON_JULIACALL_HANDLE_SIGNALS=yes
Feature Extraction with Juliacall
Once the environment is set up, you can extract radiomic features as shown below:
import SimpleITK as sitk
import numpy as np
from juliacall import Main as jl
jl.seval("using Radiomics")
# Load medical images
ct_sitk = sitk.ReadImage('DATA_PATH/ct.nii.gz')
mask_sitk = sitk.ReadImage('DATA_PATH/mask.nii.gz')
# Convert to numpy arrays
ct = sitk.GetArrayFromImage(ct_sitk)
mask = sitk.GetArrayFromImage(mask_sitk)
# Extract spacing information
spacing = ct_sitk.GetSpacing()
# Extract radiomic features
radiomic_features = jl.Radiomics.extract_radiomic_features(ct, mask, spacing)
Specific cases
If you want to extract only one subset of features, you can use:
radiomic_features = jl.Radiomics.extract_radiomic_features(ct, mask, spacing, features="first_order")
To extract two or more subsets (e.g. first order and glcm), you can run:
radiomic_features = jl.Radiomics.extract_radiomic_features(ct, mask, spacing, features=["first_order","glcm"])
To extract features from a specific label (default is 1), execute:
radiomic_features = jl.Radiomics.extract_radiomic_features(ct, mask, spacing, labels=4)
To extract features just from a list of labels, execute:
radiomic_features = jl.Radiomics.extract_radiomic_features(ct, mask, spacing, labels=[3, 5, 11])
Generate C Shared Library
It is also possible to generate a C shared library (.dll, .so, or .dylib) and call it from C/C++ code or any language that provides a C shared library interface.
Building the Library
To generate the library, navigate to the Radiomics.jl source folder and run the following in Julia (this will take a few minutes):
using PackageCompiler
create_library(".", "radiomicsjl_build";
lib_name="libradiomicsjl",
force=true,
incremental=true,
filter_stdlibs=true)
Note: The folder radiomicsjl_build will contain the shared
libraries after compilation.
Using the Shared Library from Python
To extract features in Python using the shared library (in a Linux environment), run:
import ctypes
import os
import json
import numpy as np
import SimpleITK as sitk
# Load the shared library
lib_path = os.path.abspath("SHARED_LIB_PATH/radiomicsjl_build/lib/libradiomicsjl.so")
lib = ctypes.CDLL(lib_path)
lib.init_julia(0, None)
# Define function signatures
lib.c_extract_radiomic_features.argtypes = [
ctypes.POINTER(ctypes.c_float), # img_ptr (Float32)
ctypes.c_int64, # size_x (Int64)
ctypes.c_int64, # size_y (Int64)
ctypes.c_int64, # size_z (Int64)
ctypes.POINTER(ctypes.c_float), # mask_ptr (Float32)
ctypes.c_double, # spacing_x (Float64)
ctypes.c_double, # spacing_y (Float64)
ctypes.c_double, # spacing_z (Float64)
ctypes.c_double # bin_width (Float64)
]
lib.c_extract_radiomic_features.restype = ctypes.c_char_p
# Load medical images
ct_sitk = sitk.ReadImage('DATA_PATH/CT.nrrd')
mask_sitk = sitk.ReadImage('DATA_PATH/left_parotid.nrrd')
# Convert to numpy arrays (Fortran-contiguous)
ct = np.asfortranarray(sitk.GetArrayFromImage(ct_sitk), dtype=np.float32)
mask = np.asfortranarray(sitk.GetArrayFromImage(mask_sitk), dtype=np.float32)
# Get data pointers
ptr_ct = ct.ctypes.data_as(ctypes.POINTER(ctypes.c_float))
ptr_mask = mask.ctypes.data_as(ctypes.POINTER(ctypes.c_float))
# Extract parameters
sx, sy, sz = ct.shape
spacing = ct_sitk.GetSpacing()
bin_width = 25.0
# Extract features
raw_json = lib.c_extract_radiomic_features(
ptr_ct,
sx, sy, sz,
ptr_mask,
float(spacing[0]), float(spacing[1]), float(spacing[2]),
bin_width)
# Parse results
radiomic_features = json.loads(raw_json.decode('utf-8'))
Generate C Shared Library
It is also possible to generate a C shared library (.dll, .so, or .dylib) and call it from C/C++ code or any language that provides a C shared library interface.
Building the Library
To generate the library, navigate to the Radiomics.jl source folder and run the following in Julia (this will take a few minutes):
using PackageCompiler
create_library(".", "radiomicsjl_build";
lib_name="libradiomicsjl",
force=true,
incremental=true,
filter_stdlibs=true)
Note: The folder radiomicsjl_build will contain the shared
libraries after compilation.
Call the C Shared Library from C++
Here is a C++ code snippet to extract radiomic features using the shared library on Linux:
CMakeLists.txt
cmake_minimum_required(VERSION 3.10)
project(Radiomicsjl)
# Find ITK
find_package(ITK REQUIRED)
include(${ITK_USE_FILE})
add_executable(Radiomicsjl main.cpp)
# Link ITK and system library 'dl' to load .so file
target_link_libraries(Radiomicsjl
${ITK_LIBRARIES}
dl
)
main.cpp
#include "itkImage.h"
#include "itkImageFileReader.h"
#include <iostream>
#include <vector>
#include <dlfcn.h> // to load .so file in Linux/macOS
// Define function signature of the function in the .so file
typedef const char* (*ExtractFeaturesFunc)(
float*, // img_ptr (c_float)
int64_t, int64_t, int64_t, // size_x, size_y, size_z (c_int64)
float*, // mask_ptr (c_float)
double, double, double, // spacing_x, spacing_y, spacing_z (c_double)
double // bin_width (c_double)
);
int main(int argc, char* argv[]) {
if (argc < 4) {
std::cerr << "Usage: " << argv[0] << " anatomical_image binary_mask bin_width" << std::endl;
return 1;
}
// Get parameters from command line
std::string imagePath = argv[1];
std::string maskPath = argv[2];
double binWidth = std::stod(argv[3]);
// Read images
using PixelType = float;
using ImageType = itk::Image<PixelType, 3>;
using ReaderType = itk::ImageFileReader<ImageType>;
auto readerImg = ReaderType::New();
readerImg->SetFileName(imagePath);
auto readerMask = ReaderType::New();
readerMask->SetFileName(maskPath);
try {
readerImg->Update();
readerMask->Update();
} catch (itk::ExceptionObject &ex) {
std::cerr << "ITK error: " << ex.GetDescription() << std::endl;
return 1;
}
ImageType::Pointer img = readerImg->GetOutput();
ImageType::Pointer mask = readerMask->GetOutput();
// Prepare arguments for the shared library call
auto size = img->GetLargestPossibleRegion().GetSize();
int64_t nx = size[0];
int64_t ny = size[1];
int64_t nz = size[2];
auto spacingITK = img->GetSpacing();
double sx = spacingITK[0];
double sy = spacingITK[1];
double sz = spacingITK[2];
float* imgPtr = img->GetBufferPointer();
float* maskPtr = mask->GetBufferPointer();
// Load and initialize shared library
void* handle = dlopen("SHARED_LIB_PATH/radiomicsjl_build/lib/libradiomicsjl.so", RTLD_LAZY);
if (!handle) {
std::cerr << "Error during loading the library: " << dlerror() << std::endl;
return 1;
}
typedef void (*InitJuliaFunc)(int, char**);
auto init_julia = (InitJuliaFunc) dlsym(handle, "init_julia");
if (init_julia) {
init_julia(0, NULL);
}
auto extract_features = (ExtractFeaturesFunc) dlsym(handle, "c_extract_radiomic_features");
const char* dlsym_error = dlerror();
if (dlsym_error) {
std::cerr << "Failed to find function symbol: " << dlsym_error << std::endl;
dlclose(handle);
return 1;
}
// Run feature extraction
std::cout << "Running feature extraction with bin_width = " << binWidth << " ..." << std::endl;
const char* json_result = extract_features(
imgPtr,
nx, ny, nz,
maskPtr,
sx, sy, sz,
binWidth);
if (json_result) {
std::cout << "Features (JSON):\n" << json_result << std::endl;
} else {
std::cerr << "The function returned a null pointer!" << std::endl;
}
dlclose(handle);
return 0;
}