Figure: An early test of our new 3-D agent-based cell model, growing from 10 to 80,000 agents in about 25 days (24-threaded simulation required about 5 hours). Rendered in 3D using POVRAY (with a cutaway view). [Read more ...]
Showing posts with label diffusion. Show all posts
Showing posts with label diffusion. Show all posts

Friday, January 22, 2016

BioFVM warmup: 2D continuum simulation of tumor growth

Note: This is the fourth in a series of "how-to" blog posts to help new users and developers of BioFVM. See below for guides to setting up a C++ compiler in Windows or OSX. 

What you'll need

  1. A working C++ development environment with support for OpenMP. See these prior tutorials if you need help. 
  2. A download of BioFVM, available at http://BioFVM.MathCancer.org and http://BioFVM.sf.net. Use Version 1.0.3 or later. 
  3. Matlab or Octave for visualization. Matlab might be available for free at your university. Octave is open source and available from a variety of sources

Our modeling task

We will implement a basic 2-D model of tumor growth in a heterogeneous microenvironment, with inspiration by glioblastoma models by Kristin Swanson, Russell Rockne and others (e.g., this work), and continuum tumor growth models by Hermann Frieboes, John Lowengrub, and our own lab (e.g., this paper and this paper).

We will model tumor growth driven by a growth substrate, where cells die when the growth substrate is insufficient. The tumor cells will have motility. A continuum blood vasculature will supply the growth substrate, but tumor cells can degrade this existing vasculature. We will revisit and extend this model from time to time in future tutorials.

Mathematical model

Taking inspiration from the groups mentioned above, we'll model a live cell density ρ of a relatively low-adhesion tumor cell species (e.g., glioblastoma multiforme). We'll assume that tumor cells move randomly towards regions of low cell density (modeled as diffusion with motility μ). We'll assume that that the net birth rate rBis proportional to the concentration of growth substrate σ, which is released by blood vasculature with density b. Tumor cells can degrade the tissue and hence this existing vasculature. Tumor cells die at rate rDwhen the growth substrate level is too low. We assume that the tumor cell density cannot exceed a max level ρmax. A model that includes these effects is:
where for the birth and death rates, we'll use the constitutive relations:

Mapping the model onto BioFVM

BioFVM solves on a vector u of substrates. We'll set u = [ρ , b, σ ]. The code expects PDES of the general form:
So, we determine the decay rate (λ), source function (S), and uptake function (U) for the cell density ρ and the growth substrate σ.

Cell density

We first slightly rewrite the PDE:
and then try to match to the general form term-by-term. While BioFVM wasn't intended for solving nonlinear PDEs of this form, we can make it work by quasi-linearizing, with the following functions:
.
When implementing this, we'll evaluate σ and ρ at the previous time step. The diffusion coefficient is mu, and the decay rate is zero. The target or saturation density is ρmax.

Growth substrate

Similarly, by matching the PDE for σ term-by-term with the general form, we use:
The diffusion coefficient is D, the decay rate is λ1, and the saturation density is σmax.

Blood vessels

Lastly, a term-by-term matching of the blood vessel equation gives the following functions:
The diffusion coefficient, decay rate, and saturation density are all zero.

Implementation in BioFVM

1: Start a project: Create a new directory for your project (I'd recommend "BioFVM_2D_tumor"), and enter the directory. Place a copy of BioFVM (the zip file) into your directory. Unzip BioFVM, and copy BioFVM*.h, BioFVM*.cpp, and pugixml* files into that directory.

2: Copy the matlab visualization files: To help read and plot BioFVM data, we have provided matlab files. Copy all the *.m files from the matlab subdirectory to your project.

3: Copy the empty project: BioFVM Version 1.0.3 or later includes a template project and Makefile to make it easier to get started. Copy the Makefile and template_project.cpp file to your project. Rename template_project.cpp to something useful, like 2D_tumor_example.cpp.

4: Edit the makefile: Open a terminal window and browse to your project. Tailor the makefile to your new project:

> notepad++ Makefile

Change the PROGRAM_NAME to 2Dtumor.

Also, rename main to 2D_tumor_example throughout the Makefile.

Lastly, note that if you are using OSX, you'll probably need to change from "g++" to your installed compiler. See these tutorials.

5: Start adapting 2D_tumor_example.cpp: First, open 2D_tumor_example.cpp:

> notepad++ 2D_tumor_example.cpp

Just after the "using namespace BioFVM" section of the code, define useful globals. Here and throughout, new and/or modified code is in blue:

using namespace BioFVM:

// helpful -- have indices for each "species" 
int live_cells    = 0; 
int blood_vessels = 1; 
int oxygen        = 2; 

// some globals 
double prolif_rate = 1.0 /24.0; 
double death_rate =  1.0 / 6; // 
double cell_motility = 50.0 / 365.25 / 24.0 ;
// 50 mm^2 / year  --> mm^2 / hour 
double o2_uptake_rate = 3.673 * 60.0; // 165 micron length scale 
double vessel_degradation_rate = 1.0 / 2.0 / 24.0 ; 
// 2 days to disrupt tissue
double max_cell_density = 1.0; 

double o2_supply_rate = 10.0; 
double o2_normoxic    = 1.0; 
double o2_hypoxic     = 0.2; 

6: Set up the microenvironment: Within main(), make sure we have the right number of substrates, and set them up:

// create a microenvironment, and set units 

Microenvironment M; 
M.name = "Tumor microenvironment"; 
M.time_units = "hr"; 
M.spatial_units = "mm"; 
M.mesh.units = M.spatial_units;

// set up and add all the densities you plan

M.set_density( 0 , "live cells" , "cells" );
M.add_density( "blood vessels" , "vessels/mm^2" );
M.add_density( "oxygen" , "cells" );
// set the properties of the diffusing substrates

M.diffusion_coefficients[live_cells] = cell_motility
M.diffusion_coefficients[blood_vessels] = 0; 
M.diffusion_coefficients[oxygen] = 6.0; 
// 1e5 microns^2/min in units mm^2 / hr 
M.decay_rates[live_cells] = 0;
M.decay_rates[blood_vessels] = 0;
M.decay_rates[oxygen] = 0.01 * o2_uptake_rate; 
// 1650 micron length scale 

Notice how our earlier global definitions of "live_cells", "blood_vessels", and "oxygen" makes it easier to make sure we're referencing the correct substrates in lines like these.

7: Resize the domain and test: For this example (and so the code runs very quickly), we'll work in 2D in a 2 cm × 2 cm domain:

// set the mesh size 

double dx = 0.05; // 50 microns

M.resize_space( 0.0 , 20.0 , 0, 20.0 , -dx/2.0, dx/2.0 , dx, dx, dx );  

Notice that we use a tissue thickness of dx/2 to use the 3D code for a 2D simulation. Now, let's test: 

> make 
> 2Dtumor

Go ahead and cancel the simulation [Control]+C after a few seconds. You should see something like this:

Starting program ... 

Microenvironment summary: Tumor microenvironment: 

Mesh information: 
type: uniform Cartesian
Domain: [0,20] mm x [0,20] mm x [-0.025,0.025] mm
   resolution: dx = 0.05 mm
   voxels: 160000
   voxel faces: 0
   volume: 20 cubic mm
Densities: (3 total)
   live cells:
     units: cells
     diffusion coefficient: 0.00570386 mm^2 / hr
     decay rate: 0 hr^-1
     diffusion length scale: 75523.9 mm

   blood vessels:
     units: vessels/mm^2
     diffusion coefficient: 0 mm^2 / hr
     decay rate: 0 hr^-1
     diffusion length scale: 0 mm

   oxygen:
     units: cells
     diffusion coefficient: 6 mm^2 / hr
     decay rate: 2.2038 hr^-1
     diffusion length scale: 1.65002 mm


simulation time: 0 hr (100 hr max)

Using method diffusion_decay_solver__constant_coefficients_LOD_3D (implicit 3-D LOD with Thomas Algorithm) ... 

simulation time: 10 hr (100 hr max)
simulation time: 20 hr (100 hr max)

8: Set up initial conditions: We're going to make a small central focus of tumor cells, and a "bumpy"field of blood vessels. 

// set initial conditions 
// use this syntax to create a zero vector of length 3
// std::vector<double> zero(3,0.0); 

std::vector<double> center(3);
center[0] = M.mesh.x_coordinates[M.mesh.x_coordinates.size()-1] /2.0; 
center[1] = M.mesh.y_coordinates[M.mesh.y_coordinates.size()-1] /2.0; 
center[2] = 0;
double radius = 1.0; 
std::vector<double> one( M.density_vector(0).size() , 1.0 );
double pi = 2.0 * asin( 1.0 ); 
// use this syntax for a parallelized loop over all the 
// voxels in your mesh:
#pragma omp parallel for
for( int i=0; i < M.number_of_voxels() ; i++ )
{
std::vector<double> displacement = M.voxels(i).center - center;
double distance = norm( displacement );

if( distance < radius )
{
M.density_vector(i)[live_cells] = 0.1; 
}
M.density_vector(i)[blood_vessels]= 0.5 + 0.5*cos(0.4* pi * M.voxels(i).center[0])*cos(0.3*pi *M.voxels(i).center[1]);  
M.density_vector(i)[oxygen] = o2_normoxic; 
}

9: Change to a 2D diffusion solver: 

// set up the diffusion solver, sources and sinks 
M.diffusion_decay_solver = diffusion_decay_solver__constant_coefficients_LOD_2D;

10: Set the simulation times: We'll simulate 10 days, with output every 12 hours. 

double t     = 0.0; 
double t_max = 10.0 * 24.0; // 10 days
double dt    = 0.1; 
double output_interval  = 12.0;  // how often you save data 
double next_output_time = t;     // next time you save data 

11: Set up the source function: 

void supply_function( Microenvironment* microenvironment, int voxel_index, std::vector<double>* write_here )
{
// use this syntax to access the jth substrate write_here
// (*write_here)[j]
// use this syntax to access the jth substrate in voxel voxel_index of microenvironment: 
// microenvironment->density_vector(voxel_index)[j]
static double temp1 = prolif_rate / ( o2_normoxic - o2_hypoxic ); 

(*write_here)[live_cells] = 
microenvironment->density_vector(voxel_index)[oxygen];
(*write_here)[live_cells] -= o2_hypoxic; 
if( (*write_here)[live_cells] < 0.0 )
{
(*write_here)[live_cells] = 0.0; 
else
{
(*write_here)[live_cells] = temp1; 
(*write_here)[live_cells] *= 
microenvironment->density_vector(voxel_index)[live_cells]; 
}
(*write_here)[blood_vessels] = 0.0; 
(*write_here)[oxygen]  = o2_supply_rate; 
(*write_here)[oxygen]  *=  
microenvironment->density_vector(voxel_index)[blood_vessels]; 

return; 
}

Notice the use of the static internal variable temp1: the first time this function is called, it declares this helper variable (to save some multiplication operations down the road). The static variable is available to all subsequent calls of this function. 

12: Set up the target function (substrate saturation densities): 

void supply_target_function( Microenvironment* microenvironment, int voxel_index, std::vector<double>* write_here )
{
// use this syntax to access the jth substrate write_here
// (*write_here)[j]
// use this syntax to access the jth substrate in voxel voxel_index of microenvironment: 
// microenvironment->density_vector(voxel_index)[j]
(*write_here)[live_cells] = max_cell_density;
(*write_here)[blood_vessels] =  1.0; 
(*write_here)[oxygen] = o2_normoxic;

return; 
}

13: Set up the uptake function: 

void uptake_function( Microenvironment* microenvironment, int voxel_index, std::vector<double>* write_here )
{
// use this syntax to access the jth substrate write_here
// (*write_here)[j]
// use this syntax to access the jth substrate in voxel voxel_index of microenvironment: 
// microenvironment->density_vector(voxel_index)[j]

(*write_here)[live_cells] = o2_hypoxic; 
(*write_here)[live_cells] -= 
microenvironment->density_vector(voxel_index)[oxygen]; 
if( (*write_here)[live_cells] < 0.0 ) 
{
(*write_here)[live_cells] = 0.0;
}
else
{
(*write_here)[live_cells] *= death_rate; 
}

(*write_here)[oxygen] = o2_uptake_rate ; 
(*write_here)[oxygen] *= 
microenvironment->density_vector(voxel_index)[live_cells]; 
(*write_here)[blood_vessels] = vessel_degradation_rate ; 
(*write_here)[blood_vessels] *= 
microenvironment->density_vector(voxel_index)[live_cells];
return; 
}

And that's it. The source should be ready to go!

Source files

You can download completed source for this example here:

  1. 2D_tumor_example.cpp
  2. Makefile

Using the code

Running the code

First, compile and run the code:

> make
> 2Dtumor

The output should look like this.

Starting program ... 

Microenvironment summary: Tumor microenvironment: 

Mesh information: 
type: uniform Cartesian
Domain: [0,20] mm x [0,20] mm x [-0.025,0.025] mm
   resolution: dx = 0.05 mm
   voxels: 160000
   voxel faces: 0
   volume: 20 cubic mm
Densities: (3 total)
   live cells:
     units: cells
     diffusion coefficient: 0.00570386 mm^2 / hr
     decay rate: 0 hr^-1
     diffusion length scale: 75523.9 mm

   blood vessels:
     units: vessels/mm^2
     diffusion coefficient: 0 mm^2 / hr
     decay rate: 0 hr^-1
     diffusion length scale: 0 mm

   oxygen:
     units: cells
     diffusion coefficient: 6 mm^2 / hr
     decay rate: 2.2038 hr^-1
     diffusion length scale: 1.65002 mm


simulation time: 0 hr (240 hr max)

Using method diffusion_decay_solver__constant_coefficients_LOD_2D (2D LOD with Thomas Algorithm) ... 

simulation time: 12 hr (240 hr max)
simulation time: 24 hr (240 hr max)
simulation time: 36 hr (240 hr max)
simulation time: 48 hr (240 hr max)
simulation time: 60 hr (240 hr max)
simulation time: 72 hr (240 hr max)
simulation time: 84 hr (240 hr max)
simulation time: 96 hr (240 hr max)
simulation time: 108 hr (240 hr max)
simulation time: 120 hr (240 hr max)
simulation time: 132 hr (240 hr max)
simulation time: 144 hr (240 hr max)
simulation time: 156 hr (240 hr max)
simulation time: 168 hr (240 hr max)
simulation time: 180 hr (240 hr max)
simulation time: 192 hr (240 hr max)
simulation time: 204 hr (240 hr max)
simulation time: 216 hr (240 hr max)
simulation time: 228 hr (240 hr max)
simulation time: 240 hr (240 hr max)
Done!

Looking at the data

Now, let's pop it open in matlab (or octave): 

> matlab

To load and plot a single time (e.g., the last tim)

!ls *.mat 
M = read_microenvironment( 'output_240.000000.mat' ); 
plot_microenvironment( M ); 

To add some labels: 

labels{1} = 'tumor cells'; 
labels{2} = 'blood vessel density'; 
labels{3} = 'growth substrate'; 
plot_microenvironment( M ,labels ); 

Your output should look a bit like this: 
Lastly, you might want to script the code to create and save plots of all the times. 

labels{1} = 'tumor cells'; 
labels{2} = 'blood vessel density'; 
labels{3} = 'growth substrate'; 
for i=0:20
t = i*12;
input_file = sprintf( 'output_%3.6f.mat', t ); 
output_file = sprintf( 'output_%3.6f.png', t ); 
M = read_microenvironment( input_file ); 
plot_microenvironment( M , labels ); 
print( gcf , '-dpng' , output_file );
end

Tutorial series for BioFVM

Setting up your development environment

Using BioFVM

  1. BioFVM Warmup: a 2D continuum simulation of tumor growth

Return to News   Return to MathCancer  

Monday, December 14, 2015

BioFVM: an efficient, parallelized diffusive transport solver for 3-D biological simulations

I'm very excited to announce that our 3-D diffusion solver has been accepted for publication and is now online at Bioinformatics. Click here to check out the open access preprint!
A. Ghaffarizadeh, S.H. Friedman, and P. Macklin. BioFVM: an efficient, parallelized diffusive transport solver for 3-D biological simulations. Bioinformatics, 2015.
DOI: 10.1093/bioinformatics/btv730 (free; open access)
BioFVM (stands for "Finite Volume Method for biological problems) is an open source package to solve for 3-D diffusion of several substrates with desktop workstations, single supercomputer nodes, or even laptops (for smaller problems). We built it from the ground up for biological problems, with optimizations in C++ and OpenMP to take advantage of all those cores on your CPU. The code is available at SourceForge and BioFVM.MathCancer.org

The main idea here is to make it easier to simulate big, cool problems in 3-D multicellular biology. We'll take care of secretion, diffusion, and uptake of things like oxygen, glucose, metabolic waste products, signaling factors, and drugs, so you can focus on the rest of your model. 

Design philosophy and main capabilities

Solving diffusion equations efficiently and accurately is hard, especially in 3D. Almost all biological simulations deal with this, many by using explicit finite differences (easy to code and accurate, but very slow!) or implicit methods like ADI (accurate and relatively fast, but difficult to code with complex linking to libraries). While real biological systems often depend upon many diffusing things (lots of signaling factors for cell-cell communication, growth substrates, drugs, etc.), most solvers only scale well to simulating two or three. We solve a system of PDEs of the following form: 


Above, all vector-vector products are term-by-term.

Solving for many diffusing substrates

We set out to write a package that could simulate many diffusing substrates using algorithms that were fast but simple enough to optimize. To do this, we wrote the entire solver to work on vectors of substrates, rather than on individual PDEs. In performance testing, we found that simulating 10 diffusing things only takes about 2.6 times longer than simulating one. (In traditional codes, simulating ten things takes ten times as long as simulating one.) We tried our hardest to break the code in our testing, but we failed. We simulated all the way from 1 diffusing substrate up to 128 without any problems. Adding new substrates increases the computational cost linearly.

Combining simple but tailored solvers

We used an approach called operator splitting: breaking a complicated PDE into a series of simpler PDEs and ODEs, which can be solved one at a time with implicit methods.  This allowed us to write a very fast diffusion/decay solver, a bulk supply/uptake solver, and a cell-based secretion/uptake solver. Each of these individual solvers was individually optimized. Theory tells us that if each individual solver is first-order accurate in time and stable, then the overall approach is first-order accurate in time and stable.  

The beauty of the approach is that each solver can individually be improved over time. For example, in BioFVM 1.0.2, we doubled the performance of the cell-based secretion/uptake solver. The operator splitting approach also lets us add new terms to the "main" PDE by writing new solvers, rather than rewriting a large, monolithic solver. We will take advantage of this to add advective terms (critical for interstitial flow) in future releases. 

Optimizing the diffusion solver for large 3-D domains

For the first main release of BioFVM, we restricted ourselves to Cartesian meshes, which allowed us to write very tailored mesh data structures and diffusion solvers. (Note: the finite volume method reduces to finite differences on Cartesian meshes with trivial Neumann boundary conditions.) We intend to work on more general Voronoi meshes in a future release. (This will be particularly helpful for sources/sinks along blood vessels.)

By using constant diffusion and decay coefficients, we were able to write very fast solvers for Cartesian meshes. We use the locally one-dimensional (LOD) method--a specialized form of operator splitting--to break the 3-D diffusion problem into a series of 1-D diffusion problems. For each (y,z) in our mesh, we have a 1-D diffusion problem along x. This yields a tridiagonal linear system which we can solve efficiently with the Thomas algorithm. Moreover, because the forward-sweep steps only depend upon the coefficient matrix (which is unchanging over time), we can pre-compute and store the results in memory for all the x-diffusion problems. In fact, the structure of the matrix allows us to pre-compute part of the back-substitution steps as well. Same for y- and z-diffusion. This gives a big speedup.

Next, we can use all those CPU cores to speed up our work. While the back-substitution steps of the Thomas algorithm can't be easily parallelized (it's a serial operation), we can solve many x-diffusion problems at the same time, using independent copies (instances) of the Thomas solver. So, we break up all the x-diffusion problems up across a big OpenMP loop, and repeat for y- and z-diffusion.

Lastly, we used overloaded +=, axpy and similar operations on the vector of substrates, to avoid unnecessary (and very expensive) memory allocation and copy operations wherever we could. This was a really fun code to write!

The work seems to have payed off: we have found that solving on 1 million voxel meshes (about 8 mm3 at 20 μm resolution) is easy even for laptops.

Simulating many cells

We tailored the solver to allow both lattice- and off-lattice cell sources and sinks. Desktop workstations should have no trouble with 1,000,000 cells secreting and uptaking a few substrates. 

Simplifying the non-science

We worked to minimize external dependencies, because few things are more frustrating than tracking down a bunch of libraries that may not work together on your platform. The first release BioFVM only has one external dependency: pugixml (an XML parser). We didn't link an entire linear algebra library just to get axpy and a Thomas solver--it wouldn't have been optimized for our system anyway. We implemented what we needed of the freely available .mat file specification, rather than requiring a separate library for that. (We have used these matlab read/write routines in house for several years.)

Similarly, we stuck to a very simple mesh data structure so we wouldn't have to maintain compatibility with general mesh libraries (which can tend to favor feature sets and generality over performance and simplicity).  Rather than use general-purpose ODE solvers (with yet more library dependencies, and more work for maintaining compatibility), we wrote simple solvers tailored specifically to our equations.

The upshot of this is that you don't have to do anything fancy to replicate results with BioFVM. Just grab a copy of the source, drop it into your project directory, include it in your project (e.g., your makefile), and you're good to go. 

All the juicy details

The Bioinformatics paper is just 2 pages long, using the standard "Applications Note" format. It's a fantastic format for announcing and disseminating a piece of code, and we're grateful to be published there. But you should pop open the supplementary materials, because all the fun mathematics are there: 
  • The full details of the numerical algorithm, including information on our optimizations. 
  • Convergence tests: For several examples, we showed:  
    • First-order convergence in time (with respect to Δt), and stability
    • Second-order convergence in space (with respect to Δx)
  • Accuracy tests: For each convergence test, we looked at how small Δt has to be to ensure 5% relative accuracy at Δx = 20 μm resolution. For oxygen-like problems with cell-based sources and sinks, Δt = 0.01 min will do the trick. This is about 15 times larger than the stability-restricted time step for explicit methods.   
  • Performance tests:
    • Computational cost (wall time to simulate a fixed problem on a fixed domain size with fixed time/spatial resolution) increases linearly with the number of substrates. 5-10 substrates are very feasible on desktop workstations.
    • Computational cost increases linearly with the number of voxels
    • Computational cost increases linearly in the number of cell-based source/sinks 
And of course because this code is open sourced, you can dig through the implementation details all you like! (And improvements are welcome!)

What's next? 

  • As MultiCellDS (multicellular data standard) matures, we will implement read/write support for  <microenvironment> data in digital snapshots. 
  • We have a few ideas to improve the speed of the cell-based sources and sinks. In particular, switching to a higher-order accurate solver may allow larger time step sizes, so long as the method is still stable. For the specific form of the sources/sinks, the trapezoid rule could work well here. 
  • I'd like to allow a spatially-varying diffusion coefficient. We could probably do this (at very great memory cost) by writing separate Thomas solvers for each strip in x, y, and z, or by giving up the pre-computation part of the optimization. I'm still mulling this one over. 
  • I'd also like to implement non-Cartesian meshes. The data structure isn't a big deal, but we lose the LOD optimization and Thomas solvers. In this case, we'd either use explicit methods (very slow!), use an iterative matrix solver (trickier to parallelize nicely, except in matrix-vector multiplication operations), or start with quasi-steady problems that let us use Gauss-Seidel iterative type methods, like this old paper
  • Since advective flow (particularly interstitial flow) is so important for many problems, I'd like to add an advective solver. This will require some sort of upwinding to maintain stability. 
  • At some point, we'd like to port this to GPUs. However, I don't currently have time / resources to maintain a separate CUDA or OpenCL branch. (Perhaps this will be an excuse to learn Julia on GPUs.)
Well, we hope you find BioFVM useful. If you give it a shot, I'd love to hear back from you!

Very best -- Paul