pyLM  v1.1
A Python Problem Solving Environment for the simulation of stochastic biological systems
 All Classes Namespaces Files Functions Variables Pages
Public Member Functions | Public Attributes | List of all members
pyLM.RDME.RDMESimulation Class Reference

A class that contains all regions, reactions, diffusions and rules for a RDME simulation. More...

Public Member Functions

def __init__
 Specify a cuboid region that represents the extents to the reaction region as well as the lattice spacing. More...
 
def addRegion
 Add a region to the simulation. More...
 
def addCuboidRegion
 Add a cuboid to the builder. More...
 
def addShape
 
def modifyRegion
 Return a pointer to a region so that it may be modified. More...
 
def defineSpecies
 Define a specie/s of particles that exist in the simulation. More...
 
def siteVolume
 Get the actual volume of a specific site in L. More...
 
def buildCapsidCell
 Build a capsule based shell in this RDMESimulation centered within the simulation domain that includes a membrane and cytoplasm. More...
 
def buildSphericalCell
 Build a spherical based shell in this RDMESimulation centered within the simulation domain that includes a membrane and cytoplasm. More...
 
def addParticles
 Add a specified number of particles of the specified type to the specified region. More...
 
def packRegion
 Add nonmoving obstacles to a particular region. More...
 
def setTransitionRate
 Specify the diffusion rate between species; this is a one directional rate e.g. More...
 
def setTwoWayTransitionRate
 Specify the diffusion rate between species; this is a two directional rate. More...
 
def buildDiffusionModel
 Return the Lattice Microbes DiffusionModel object for fine-tuning. More...
 
def buildReactionModel
 Return the Lattice Microbes ReactionModel object for fine-tuning. More...
 
def buildSpatialModel
 Return the Lattice Microbes SpatialModel object for fine-tuning. More...
 
def setWriteInterval
 
def setLatticeWriteInterval
 
def setSimulationTime
 
def setTimestep
 
def setRandomSeed
 
def getLattice
 Get a discretized version of the simulation domain. More...
 
def setLatticeSite
 Set a particular lattice site type. More...
 
def getLatticeSite
 Get a particular lattice site. More...
 
def addParticleAt
 Add a particle/ to a particular site. More...
 
def addParticleAtIdx
 Create an HDF5 version of the simulation amenable for later running or stand-alone running. More...
 
def save
 Create an HDF5 version of the simulation amenable for later running or stand-alone running. More...
 
def run
 
def runMPI
 Run the simulation using a call to mpirun with the given options. More...
 
def runSolver
 
def buildVector
 
def buildPoint
 
def buildSphere
 
def buildEllipse
 
def buildDifference
 

Public Attributes

 siteTypes
 
 particleMap
 
 customAddedParticleList
 
 species_id
 
 regions
 
 transitionRates
 
 initial_counts
 
 parameters
 
 continousDimensions
 
 latticeSpacing
 
 lm_builder
 
 lattice
 
 hasBeenDiscretized
 
 name
 
 replicates
 
 filename
 

Detailed Description

A class that contains all regions, reactions, diffusions and rules for a RDME simulation.

Constructor & Destructor Documentation

def pyLM.RDME.RDMESimulation.__init__ (   self,
  dimensions,
  spacing,
  name = "unnamed",
  defaultRegion = "default" 
)

Specify a cuboid region that represents the extents to the reaction region as well as the lattice spacing.

Parameters
self
dimensionsA list of [x,y,z]
spacingLattice spacing
nameThe name of the RDME simulation; default: "unnamed"
defaultRegionThe name of the region that is associated with the lattice sites before any other regions are added; default:"default"
Returns
self

Member Function Documentation

def pyLM.RDME.RDMESimulation.addCuboidRegion (   self,
  name,
  a,
  b 
)

Add a cuboid to the builder.

Parameters
self
nameName of the site type for this region
atuple for the first corner in continous space
btuple for the second corner in continous space
def pyLM.RDME.RDMESimulation.addParticleAt (   self,
  index,
  particleType 
)

Add a particle/ to a particular site.

Parameters
index(x,y,z) a list of spatial location
specieThe specie type to add
def pyLM.RDME.RDMESimulation.addParticleAtIdx (   self,
  index,
  particleType 
)

Create an HDF5 version of the simulation amenable for later running or stand-alone running.

Parameters
self
filenameA file to write the simulation to Add a particle/ to a particular site
index(x,y,z) a list of lattice site indices
specieThe specie type to add
def pyLM.RDME.RDMESimulation.addParticles (   self,
  species = 'unknown',
  region = 'default',
  count = 1 
)

Add a specified number of particles of the specified type to the specified region.

Parameters
self
speciesThe species to add to the region
regionThe region to add particles to
countNumber of particles to add (default: 1)
Returns
The simulation object
def pyLM.RDME.RDMESimulation.addRegion (   self,
  region 
)

Add a region to the simulation.

Parameters
self
regionThe region to add to the simulation
Returns
The region just added
def pyLM.RDME.RDMESimulation.addShape (   self,
  shape 
)
def pyLM.RDME.RDMESimulation.buildCapsidCell (   self,
  length,
  diameter,
  membraneThickness,
  points = False 
)

Build a capsule based shell in this RDMESimulation centered within the simulation domain that includes a membrane and cytoplasm.

Parameters
self
lengthThe length of the capsule from one sphere origin to the other
diameterThe diameter of the cell
membraneThicknessThe thickness of the membrane
pointsOPTIONAL: List of lists containing the coordinates of the centers of the sphereoids that cap the capsid cell, e.g. [[x1, y1, z1], [x2, y2, z2]]. If not given, cell is centered in the volume and aligned in the z-direction
Returns
The simulation object
def pyLM.RDME.RDMESimulation.buildDifference (   self,
  arg1,
  arg2,
  arg3 
)
def pyLM.RDME.RDMESimulation.buildDiffusionModel (   self)

Return the Lattice Microbes DiffusionModel object for fine-tuning.

Parameters
selfThe diffusion modle for the simulation
def pyLM.RDME.RDMESimulation.buildEllipse (   self,
  arg1,
  arg2,
  arg3,
  arg4,
  arg5,
  arg6,
  arg7 
)
def pyLM.RDME.RDMESimulation.buildPoint (   self,
  arg1,
  arg2,
  arg3 
)
def pyLM.RDME.RDMESimulation.buildReactionModel (   self)

Return the Lattice Microbes ReactionModel object for fine-tuning.

Parameters
self
Returns
The reaction model for this simulation
def pyLM.RDME.RDMESimulation.buildSpatialModel (   self)

Return the Lattice Microbes SpatialModel object for fine-tuning.

Parameters
self
Returns
The spatial model (i.e. obstacles, sites, etc.) for this simulation
def pyLM.RDME.RDMESimulation.buildSphere (   self,
  arg1,
  arg2,
  arg3 
)
def pyLM.RDME.RDMESimulation.buildSphericalCell (   self,
  diameter,
  membraneThickness,
  point = False 
)

Build a spherical based shell in this RDMESimulation centered within the simulation domain that includes a membrane and cytoplasm.

Parameters
self
diameterThe diameter of the cell
membraneThicknessThe thickness of the membrane
pointThe center of the spherical cell
Returns
The simulation object
def pyLM.RDME.RDMESimulation.buildVector (   self,
  arg1,
  arg2,
  arg3 
)
def pyLM.RDME.RDMESimulation.defineSpecies (   self,
  species 
)

Define a specie/s of particles that exist in the simulation.

Parameters
self
speciesA list of species to add to the simulation
Returns
The simulation object
def pyLM.RDME.RDMESimulation.getLattice (   self)

Get a discretized version of the simulation domain.

Call this after building all spherical and capsule cells

Parameters
self
Returns
A lattice object. This function should only be called once.
def pyLM.RDME.RDMESimulation.getLatticeSite (   self,
  index 
)

Get a particular lattice site.

Parameters
index(x,y,z) a list of coordinates
Returns
The type of the lattice site (string)
def pyLM.RDME.RDMESimulation.modifyRegion (   self,
  region 
)

Return a pointer to a region so that it may be modified.

Parameters
self
regionGet a region that is attached to the simulation for modification
Returns
The region to modify
def pyLM.RDME.RDMESimulation.packRegion (   self,
  region,
  radius,
  percentage,
  obstacleID 
)

Add nonmoving obstacles to a particular region.

Parameters
self
regionThe name of the region in which to add particles to
radiusThe radius of the particles
percentageThe percentage of the total region volume that should be packed
obstacleIDan identifier for the obstacle
Returns
The simulation object
def pyLM.RDME.RDMESimulation.run (   self,
  filename,
  method,
  replicates = 1,
  seed = None 
)
def pyLM.RDME.RDMESimulation.runMPI (   self,
  filename,
  method,
  replicates = 1,
  driver = "mpirun",
  ppe = 1,
  seed = None 
)

Run the simulation using a call to mpirun with the given options.

Parameters
self
filenameThe HDF file to write to
methodThe class name for the solver to use (e.g., lm::cme::GillespieDSolver")
replicatesThe number of replicates to serially run
driverThe program to execute the parallel run, e.g. "mpirun", "aprun", "ibrun", etc.
ppeThe number of processing elements to use
seedA seed for the random number generator to use when running the simulation; None indicates default
def pyLM.RDME.RDMESimulation.runSolver (   self,
  filename,
  solver,
  replicates = 1 
)
def pyLM.RDME.RDMESimulation.save (   self,
  filename 
)

Create an HDF5 version of the simulation amenable for later running or stand-alone running.

Parameters
self
filenameA file to write the simulation to
def pyLM.RDME.RDMESimulation.setLatticeSite (   self,
  index,
  siteType 
)

Set a particular lattice site type.

Parameters
index(x,y,z) a list of coordinates
siteTypeThe type to set the lattice point to. This would be the name of a region that has previously been performed"
def pyLM.RDME.RDMESimulation.setLatticeWriteInterval (   self,
  time 
)
def pyLM.RDME.RDMESimulation.setRandomSeed (   self,
  seed 
)
def pyLM.RDME.RDMESimulation.setSimulationTime (   self,
  time 
)
def pyLM.RDME.RDMESimulation.setTimestep (   self,
  time 
)
def pyLM.RDME.RDMESimulation.setTransitionRate (   self,
  species,
  via,
  to,
  rate 
)

Specify the diffusion rate between species; this is a one directional rate e.g.

membrane->cytosol or extracellular->membrane

Parameters
self
speciesThe specie that can transition between regions
viaFrom this region
toTo this region
rateDiffusion rate between regions
Returns
The simulation object
def pyLM.RDME.RDMESimulation.setTwoWayTransitionRate (   self,
  species,
  one,
  two,
  rate 
)

Specify the diffusion rate between species; this is a two directional rate.

Parameters
self
speciesThe specie that can transition between regions
oneA region
twoThe other region
rateDiffusion rate between regions
Returns
The simulation object
def pyLM.RDME.RDMESimulation.setWriteInterval (   self,
  time 
)
def pyLM.RDME.RDMESimulation.siteVolume (   self)

Get the actual volume of a specific site in L.

Returns
the volume of the site in L

Member Data Documentation

pyLM.RDME.RDMESimulation.continousDimensions
pyLM.RDME.RDMESimulation.customAddedParticleList
pyLM.RDME.RDMESimulation.filename
pyLM.RDME.RDMESimulation.hasBeenDiscretized
pyLM.RDME.RDMESimulation.initial_counts
pyLM.RDME.RDMESimulation.lattice
pyLM.RDME.RDMESimulation.latticeSpacing
pyLM.RDME.RDMESimulation.lm_builder
pyLM.RDME.RDMESimulation.name
pyLM.RDME.RDMESimulation.parameters
pyLM.RDME.RDMESimulation.particleMap
pyLM.RDME.RDMESimulation.regions
pyLM.RDME.RDMESimulation.replicates
pyLM.RDME.RDMESimulation.siteTypes
pyLM.RDME.RDMESimulation.species_id
pyLM.RDME.RDMESimulation.transitionRates

The documentation for this class was generated from the following file: