Spatiocyte A lattice-based particle simulator

Spatiocyte sources

  1. Serial version of Spatiocyte: Github

  2. Parallel version of Spatiocyte (pSpatiocyte): Github

Spatiocyte Quick Start

  1. Get and run the install script

    Spatiocyte runs stably on Ubuntu Linux, whereas on Mac OS X Yosemite and Windows systems, it is still experimental. On a fresh Ubuntu or a Mac OS X system, Spatiocyte requires several additional packages to run. To install these packages and Spatiocyte, open a terminal and execute the following instructions:

    Ubuntu Linux:
    $ wget https://raw.githubusercontent.com/ecell/spatiocyte/master/install-spatiocyte-ubuntu.sh
    $ sh -x install-spatiocyte-ubuntu.sh 
    
    Mac OS X:
    $ curl -O https://raw.githubusercontent.com/ecell/spatiocyte/master/install-spatiocyte-mac.sh
    $ sh -x install-spatiocyte-mac.sh 
    

    Enter your password when requested since some packages require the administrator privilege to install. If you have any issues during install, post the error messages to the Spatiocyte Users forum. On Mac, the installation script will take a longer time to finish executing because more packages need to be downloaded. We can use Blender to render simulation snapshots and VLC to view movies of the snapshots. Since the installation script does not install them on Mac, you can download and install them yourself from http://www.blender.org/ and http://www.videolan.org/ .

    Windows (with Docker container):
    1. Install the latest release of Docker toolbox for Windows and UltraVNC.
      If you are using Windows 10 Pro or Enterprise 64-bit, then you can install Docker for Windows and Kitematic instead of Docker toolbox for Windows. Kitematic can be downloaded from Dashboard menu after installing Docker for Windows.
    2. Run Kitematic from your desktop or start menu.
    3. Search Spatiocyte and create container for ecell/spatiocyte.
      ecell/spatiocyte
    4. Clic EXEC button and go into the Spatiocyte container. go into container
    5. Issue the following command in the Spatiocyte container.
      USER=root vncserver :1 -geometry 1360x768 -depth 24
      
    6. Connect Docker container ACCESS URL written in Kitematic IP & PORTS panel with UltraVNC. The password for this VNC login is password. ACCESS URL for VNC UltraVNC
  2. Test the installation

    Close and reopen the terminal for the installation to take effect. To test if the installation is successful, run the following command in the terminal:

    $ ecell3-session-monitor
    

    The above command will open the E-Cell Session Monitor, as shown below. You have now sucessfully installed Spatiocyte! If for some reason, the Session Monitor does not come up, please post the error message to the Spatiocyte Users forum.

    E-Cell Session Monitor
  3. Run a model

    Now try running a simple 1D diffusion Python model available in the examples/1D directory:

    $ cd $HOME/wrk/spatiocyte/examples/1D
    $ ecell3-session 1D.py
    

    The simulator will terminate itself once the simulation has ended. To view the simulated molecules, run the Spatiocyte Visualizer by issuing:

    $ spatiocyte
    

    More models can be found under the examples directory . To run them, follow the instructions in the README file in the directory.

  4. Any issues? Contact us!

    If you have any issues or questions about Spatiocyte, message us. Spatiocyte documentation can be found here .

Build a simple 3D reaction-diffusion model

  1. Specify the stepper and voxel radius

    To build a simple 3D model, open a text editor and add the Python scripts as provided in the code boxes below. Each Spatiocyte model requires the SpatiocyteStepper. It advances the simulation in an event-driven manner. The radius of the hexagonal close-packed lattice voxels can be set in the SpatiocyteStepper using the VoxelRadius property. In the code below, we have set the radius to 4.4 nm.

    sim = theSimulator.createStepper('SpatiocyteStepper', 'SS')
    sim.VoxelRadius = 4.4e-9
    theSimulator.rootSystem.StepperID = 'SS'
    

    The memory usage increases linearly with the number of voxels. In a 64-bit system, each voxel typically takes up 108 bytes of memory. The number of voxels with radius, r in a volume, V is given by V/(4r 3 2 0.5 ). In each compartment, the StepperID must be set to the SpatiocyteStepper ID. Here we only have the root compartment, rootSystem and we set its StepperID to 'SS' which is the SpatiocyteStepper ID given in the first line.

  2. Define compartment geometry

    We will use a cuboid compartment geometry for this model. The GEOMETRY variable of a volume compartment specifies one of the six supported geometric primitives: cuboid (‘0’), ellipsoid (‘1’), cylinder (‘2’), rod (‘3’), pyramid (‘4’) and erythrocyte (‘5’). More complex forms can be constructed using a combination of these primitives. Compartments without the GEOMETRY definition is set to the cuboid form since the default value is ‘0’.

    theSimulator.createEntity('Variable', 'Variable:/:GEOMETRY').Value = 0
    theSimulator.createEntity('Variable', 'Variable:/:LENGTHX').Value = 2.5e-7
    theSimulator.createEntity('Variable', 'Variable:/:LENGTHY').Value = 2.5e-7
    theSimulator.createEntity('Variable', 'Variable:/:LENGTHZ').Value = 2.5e-7
    

    The three variables LENGTH[X, Y, Z] can specify the compartment lengths in the directions of [x, y, z]-axes, respectively. In the example above, the lengths are set to 250 nm. More information about the LENGTH variable definition can be found in the documentation .

  3. Assign compartment boundary types

    When a volume compartment has the cuboid geometry, as in our example model, the boundary type can be specified using [XY, XZ, YZ]PLANE variables. The boundary type can be reflective (‘0’), periodic (‘1’) or semi-periodic (‘2’). A semi-periodic boundary allows nonHD molecules to move unidirectionally from one boundary to the other. We set all the boundaries to 1 to make them periodic:

    theSimulator.createEntity('Variable', 'Variable:/:XYPLANE').Value = 1
    theSimulator.createEntity('Variable', 'Variable:/:XZPLANE').Value = 1
    theSimulator.createEntity('Variable', 'Variable:/:YZPLANE').Value = 1
    
  4. Set species and their initial numbers

    Every compartment must have a VACANT variable that represents the 'species' of empty voxels within the compartment. We would like to add the reaction A + B -> C + C in the model, with A and B species each made up of 500 molecules initially. These three variables can be declared as below:

    theSimulator.createEntity('Variable', 'Variable:/:VACANT')
    theSimulator.createEntity('Variable', 'Variable:/:A').Value = 500
    theSimulator.createEntity('Variable', 'Variable:/:B').Value = 500
    theSimulator.createEntity('Variable', 'Variable:/:C').Value = 0
    
  5. Populate non-zero species

    We must populate the molecules of A and B species in the compartment before the simulation is started because they initially have non-zero numbers. We populate them using the MoleculePopulateProcess.

    populator = theSimulator.createEntity('MoleculePopulateProcess', 'Process:/:pop')
    populator.VariableReferenceList = [['_', 'Variable:/:A']]
    populator.VariableReferenceList = [['_', 'Variable:/:B']]
    
  6. Fix the diffusion property of each species

    We are going to diffuse all three species, A, B and C in the compartment with a diffusion coefficient of 0.1 um 2 /s. We specify the diffusion property of each species with the DiffusionProcess:

    diffuser = theSimulator.createEntity('DiffusionProcess', 'Process:/:diffuseA')
    diffuser.VariableReferenceList = [['_', 'Variable:/:A']]
    diffuser.D = 1e-13
    
    diffuser = theSimulator.createEntity('DiffusionProcess', 'Process:/:diffuseB')
    diffuser.VariableReferenceList = [['_', 'Variable:/:B']]
    diffuser.D = 1e-13
    
    diffuser = theSimulator.createEntity('DiffusionProcess', 'Process:/:diffuseC')
    diffuser.VariableReferenceList = [['_', 'Variable:/:C']]
    diffuser.D = 1e-13
    
  7. Specify the reaction

    Next, we specify the reaction A + B -> C + C using the DiffusionInfluencedReaction process. We set the reaction probability, p to unity, which means a C molecule will be created whenever an A and a B molecule collide.

    binder = theSimulator.createEntity('DiffusionInfluencedReactionProcess', 'Process:/:reaction1')
    binder.VariableReferenceList = [['_', 'Variable:/:A','-1']]
    binder.VariableReferenceList = [['_', 'Variable:/:B','-1']]
    binder.VariableReferenceList = [['_', 'Variable:/:C','1']]
    binder.VariableReferenceList = [['_', 'Variable:/:C','1']]
    binder.p = 1
    
  8. Set the logger for Spatiocyte visualizer

    We can use the VisualizationLogProcess to log the coordinates of the diffusing species at specified intervals or at the shortest stepper intervals (default). The coordinates will be saved in the file VisualLog.dat in binary format. The Spatiocyte Visualizer can load the log file to display the molecules in 3D while the simulation is running or after it has ended.

    logger = theSimulator.createEntity('VisualizationLogProcess', 'Process:/:logger')
    logger.VariableReferenceList = [['_', 'Variable:/:A']]
    logger.VariableReferenceList = [['_', 'Variable:/:B']]
    logger.VariableReferenceList = [['_', 'Variable:/:C']]
    
  9. Set molecule coordinates logger

    The coordinates of the molecules can also be saved in csv format with CoordinateLogProcess:

    coord = theSimulator.createEntity('CoordinateLogProcess', 'Process:/:coord')
    coord.VariableReferenceList = [['_', 'Variable:/:A']]
    coord.VariableReferenceList = [['_', 'Variable:/:B']]
    coord.VariableReferenceList = [['_', 'Variable:/:C']]
    
  10. Specify the model run time

    Finally, we need to tell the simulator how long the model should be run. Here, we will run it for 0.08 s.

    run(0.08)
    
  11. Run the model

    Save the file as 3D.py and run the model in a terminal by issuing,

    $ ecell3-session 3D.py
    

    We can visualize the simulation with

    $ spatiocyte
    

    The above instruction will start the Spatiocyte visualizer, which will load the VisualLog.dat file and display the molecule positions in 3D. The shortcut keys to manipulate the visualizer are provided here. The complete model file can be downloaded here.

References

  1. Watabe, M., Yoshimura, H., Arjunan, S. N. V., Kaizu, K. & Takahashi, K. (2020). Signaling activations through G-protein-coupled-receptor aggregations. Physical Review E 102: 032413.
  2. Arjunan, S. N. V., Miyauchi, A., Iwamoto, K. & Takahashi, K. (2020). pSpatiocyte: a high-performance simulator for intracellular reaction-diffusion systems. BMC Bioinformatics 21, 33.
  3. Watabe, M., Arjunan, S. N. V., Chew, W.-X., Kaizu, K. & Takahashi, K. (2019). Cooperativity transitions driven by higher-order oligomer formations in ligand-induced receptor dimerization. Physical Review E, 100: 062407.
  4. Watabe, M., Arjunan, S. N. V., Chew, W.-X., Kaizu, K. & Takahashi, K. (2019). Simulation of live-cell imaging system reveals hidden uncertainties in cooperative binding measurements. Physical Review E, 100: 010402(R).
  5. Chew, W.-X., Kaizu, K., Watabe, M., Muniandy, S. V., Takahashi, K. & Arjunan, S. N. V. (2019). Surface reaction-diffusion kinetics on lattice at the microscopic scale. Physical Review E, 99: 042411.
  6. Chew, W.-X., Kaizu, K., Watabe, M., Muniandy, S. V., Takahashi, K. & Arjunan, S. N. V. (2018). Reaction-diffusion kinetics on lattice at the microscopic scale. Physical Review E, 98: 032418.
  7. Arjunan, S. N. V. & Takahashi, K. (2017). Multi-algorithm particle simulations with Spatiocyte. Methods in Molecular Biology, 1611: 219-236.
  8. Watabe, M., Arjunan, S. N. V., Fukushima, S., Iwamoto, K., Kozuka, J., Matsuoka, S., Shindo, Y., Ueda, M. & Takahashi, K. (2015). A Computational Framework for Bioimaging Simulation. PLoS ONE 10(7): e0130089.
  9. Shimo, H., Arjunan, S.N.V., Machiyama, H., Nishino, T., Suematsu, M., Fujita, H., Tomita, M. & Takahashi, K. (2015). Particle Simulation of Oxidation Induced Band 3 Clustering in Human Erythrocytes. PLoS Computational Biology, 11(6): e1004210.
  10. Arjunan, S. N. V. (2013). A guide to modeling reaction–diffusion of molecules with the E-Cell system, in: E-Cell System: basic concepts and applications, Springer, pp. 43–62. The latest version of the user guide can be found at http://spatiocyte.readthedocs.org/en/latest/intro.html
  11. Arjunan, S. N. V., & Tomita, M. (2010). A new multicompartmental reaction-diffusion modeling method links transient membrane attachment of E. coli MinE to E-ring formation. Systems and Synthetic Biology, 4(1), 35–53.
  12. Arjunan, S. N. V., & Tomita, M. (2009). Modeling reaction-diffusion of molecules on surface and in volume spaces with the E-Cell System. International Journal of Computer Science and Information Security, 3(1), 211–216.