TGSTK
0.0.1
The Tumour Growth Simulation ToolKit
|
Finite difference solver for reaction-diffusion tumour growth simulation over regular grids. More...
#include <tgstkFiniteDifferenceReactionDiffusionTumourGrowthFilter.h>
Public Member Functions | |
tgstkFiniteDifferenceReactionDiffusionTumourGrowthFilter () | |
~tgstkFiniteDifferenceReactionDiffusionTumourGrowthFilter () | |
bool | check () |
void | execute () |
vtkSmartPointer< vtkImageData > | getFinalCellDensityImage () |
vtkSmartPointer< vtkImageData > | getFinalCellDensityGradientImage () |
void | setBrainMapImage (vtkImageData *image) |
void | setDiffusionTensorImage (vtkSmartPointer< vtkImageData > image) |
void | setInitialCellDensityImage (vtkSmartPointer< vtkImageData > image) |
void | setProliferationRateImage (vtkSmartPointer< vtkImageData > image) |
void | setSimulatedTime (double time) |
void | setTimeStep (double step) |
Public Member Functions inherited from tgstkImageProcessorBase | |
virtual | ~tgstkImageProcessorBase () |
Public Member Functions inherited from tgstkAlgorithmBase | |
virtual | ~tgstkAlgorithmBase () |
bool | update () |
Public Member Functions inherited from tgstkObjectBase | |
virtual | ~tgstkObjectBase () |
std::string | getObjectName () |
Private Attributes | |
double | simulatedTime |
double | timeStep |
vtkSmartPointer< vtkImageData > | brainMapImage |
vtkSmartPointer< vtkImageData > | finalCellDensityImage |
vtkSmartPointer< vtkImageData > | finalCellDensityGradientImage |
vtkSmartPointer< vtkImageData > | initialCellDensityImage |
vtkSmartPointer< vtkImageData > | proliferationRateImage |
vtkSmartPointer< vtkImageData > | diffusionTensorImage |
Additional Inherited Members | |
Protected Member Functions inherited from tgstkImageProcessorBase | |
tgstkImageProcessorBase () | |
bool | _assertImageNumberOfScalarComponents (vtkSmartPointer< vtkImageData > image, std::vector< int > numberOfScalarComponents, std::string name) |
bool | _assertImageScalarType (vtkSmartPointer< vtkImageData > image, std::vector< int > scalarTypes, std::string name) |
bool | assertEqualImageDimensions (std::vector< vtkSmartPointer< vtkImageData >> images) |
bool | assertEqualImageSpacings (std::vector< vtkSmartPointer< vtkImageData >> images) |
Protected Member Functions inherited from tgstkAlgorithmBase | |
tgstkAlgorithmBase () | |
template<typename Type > | |
bool | _assertNotNullPtr (Type var, std::string name) |
template<typename Type > | |
bool | _assertValueInRange (Type var, Type min, Type max, std::string name) |
template<typename Type1 , typename Type2 > | |
bool | _assertValueIsEqual (Type1 var1, Type2 var2, std::string name1, std::string name2) |
Protected Member Functions inherited from tgstkObjectBase | |
tgstkObjectBase () | |
Static Protected Member Functions inherited from tgstkImageProcessorBase | |
template<typename Type > | |
static void | fillImage (vtkSmartPointer< vtkImageData > image, Type value) |
static vtkSmartPointer< vtkImageData > | getNewImageFromReferenceImage (vtkSmartPointer< vtkImageData > reference, int type, int numberOfComponents=1) |
Protected Attributes inherited from tgstkObjectBase | |
std::string | objectName |
Finite difference solver for reaction-diffusion tumour growth simulation over regular grids.
tgstkFiniteDifferenceReactionDiffusionTumourGrowthFilter solves the reaction-diffusion tumour growth problem introduced in [4] [2] [7] [6] [8] :
\begin{cases} &\frac{\partial c(\bar{r}, t)}{\partial t} = \bar{\nabla} \cdot \left( \bar{\bar{D}}(\bar{r}) \, \bar{\nabla} c(\bar{r}, t) \right) + \rho(\bar{r}) \, c(\bar{r}, t) \left( 1-c(\bar{r}, t) \right) & \forall \bar{r} \in \Omega, \; \forall t > 0 \\ &c(\bar{r}, 0) = c_0(\bar{r}) & \forall \bar{r} \in \Omega \\ &\left( \bar{\bar{D}}(\bar{r}) \, \bar{\nabla} c(\bar{r}, t) \right) \cdot \bar{n}_{\partial \Omega}(\bar{r}) = 0 & \forall \bar{r} \in \partial \Omega \end{cases}
where \(c(\bar{r}, t)\) is the nomarlised tumour cell density at location \(\bar{r}\) and time \(t\) with \(c(\bar{r}, t) \in [0, 1], \; \forall \bar{r}, t\); \(\bar{\bar{D}}(\bar{r})\) is the tumour diffusion tensor field; \(\rho(\bar{r})\) is the tumour proliferation rate field; \(\Omega\) is the solving domain; \(c_0(\bar{r})\) is the initial normalised cell density field at time \(t=0\); and \(\bar{n}_{\partial \Omega}(\bar{r})\) is the unit normal vector pointing outwards the domain boundary \(\partial \Omega\) at location \(\bar{r} \in \partial \Omega\).
The problem is solved using a forward Euler finite difference scheme implemented in CUDA:
\[ c_{i,j,k}^{n+1} = c_{i,j,k}^n + \Delta t \left(\text{div}_{i,j,k}^n + \rho_{i,j,k} \, c_{i,j,k}^n \, (1-c_{i,j,k}^n) \right) \]
where subscript \(i, j, k\) refers to voxel \((i,j,k)\); superscript \(n\) refers to simulation step \(n\); \(\Delta t\) is the time step; and \(\text{div}_{i,j,k}^n\) is the divergence of the tumour cell flux given by \(\bar{\nabla} \cdot \left(\bar{\bar{D}}(\bar{r}) \, \bar{\nabla} c(\bar{r}, t) \right)\). The divergence term \(\text{div}_{i,j,k}^n\) is computed using a 3D extension of the standard stencil presented in [5].
The initial tumour cell density \(c_0(\bar{r})\), diffusion tensor \(\bar{\bar{D}}(\bar{r})\), and proliferation rate \(\rho(\bar{r})\) fields as well as the brain map defining the domain \(\Omega\) can be specified by the user as vtkImageData objects. The simulated time \(T\) and time step \(\Delta t\) can also be specified.
tgstkFiniteDifferenceReactionDiffusionTumourGrowthFilter::tgstkFiniteDifferenceReactionDiffusionTumourGrowthFilter | ( | ) |
tgstkFiniteDifferenceReactionDiffusionTumourGrowthFilter::~tgstkFiniteDifferenceReactionDiffusionTumourGrowthFilter | ( | ) |
|
virtual |
Reimplemented from tgstkAlgorithmBase.
|
virtual |
Implements tgstkAlgorithmBase.
vtkSmartPointer< vtkImageData > tgstkFiniteDifferenceReactionDiffusionTumourGrowthFilter::getFinalCellDensityGradientImage | ( | ) |
Gets the final normalised tumour cell density gradient image \(\bar{\nabla}c(\bar{r}, T)\) in \(\text{mm}^{-1}\).
The image scalar type is VTK_DOUBLE. The image has 3 scalar components ( \(\!\partial_x c, \partial_y c, \partial_z c\)).
vtkSmartPointer< vtkImageData > tgstkFiniteDifferenceReactionDiffusionTumourGrowthFilter::getFinalCellDensityImage | ( | ) |
Gets the final normalised tumour cell density image \(c(\bar{r}, T)\).
The image scalar type is VTK_DOUBLE.
void tgstkFiniteDifferenceReactionDiffusionTumourGrowthFilter::setBrainMapImage | ( | vtkImageData * | image | ) |
Sets the bain map image defining the solving domain \(\Omega\).
The image scalar type must be VTK_UNSIGNED_SHORT. The different brain tissues must be referenced as specified by the tgstkBrainTissueType enum type.
void tgstkFiniteDifferenceReactionDiffusionTumourGrowthFilter::setDiffusionTensorImage | ( | vtkSmartPointer< vtkImageData > | image | ) |
Sets the tumour diffusion tensor image \(\bar{\bar{D}} (\bar{r})\) in \(\text{mm}^2\text{ d}^{-1}\).
The image scalar type must be VTK_DOUBLE. The image must have 6 scalar components ( \(\!D_{xx}, D_{xy}, D_{xz}, D_{yy}, D_{yz}, D_{zz}\) – implicitly symmetric).
void tgstkFiniteDifferenceReactionDiffusionTumourGrowthFilter::setInitialCellDensityImage | ( | vtkSmartPointer< vtkImageData > | image | ) |
Sets the initial normalised tumour cell density image \(c_0(\bar{r})\).
The image scalar type must be VTK_DOUBLE. The image values must be in range \([0; 1]\).
void tgstkFiniteDifferenceReactionDiffusionTumourGrowthFilter::setProliferationRateImage | ( | vtkSmartPointer< vtkImageData > | image | ) |
Sets the tumour proliferation rate image \(\rho(\bar{r})\) in \(\text{d}^{-1}\).
The image scalar type must be VTK_DOUBLE. The image values must be \(\geq 0\).
void tgstkFiniteDifferenceReactionDiffusionTumourGrowthFilter::setSimulatedTime | ( | double | time | ) |
Sets the time duration to simulate \(T\) in \(\text{d}\).
The value must be \(\geq 0\).
void tgstkFiniteDifferenceReactionDiffusionTumourGrowthFilter::setTimeStep | ( | double | step | ) |
Sets the smilation time step \(\Delta t\) in \(\text{d}\).
The value must be \(> 0\).
\[ \Delta t \leq \min_{\bar{r}} \frac{1}{2} \left(\frac{D_{xx} (\bar{r})}{\Delta x^2} + \frac{D_{yy}(\bar{r})}{\Delta y^2} + \frac{D_{zz}(\bar{r})}{\Delta z^2}\right)^{-1} \]
where \(\Delta x\), \(\Delta y\), and \(\Delta z\) are respectively the voxel size in \(x\), \(y\), and \(z\) in \(\text{mm}\).
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |