Lattice Microbes 2.5
This is for whole cell modeling
Loading...
Searching...
No Matches
pyLM.RDME.RDMESimulation Class Reference
Inheritance diagram for pyLM.RDME.RDMESimulation:
[legend]

Public Member Functions

 __init__ (self, dimensions, spacing, name="unnamed", defaultRegion="default")
 addRegion (self, region)
 addCuboidRegion (self, name, a, b)
 addShape (self, shape)
 modifyRegion (self, region)
 defineSpecies (self, species)
 siteVolume (self)
 buildCapsidCell (self, length, diameter, membraneThickness, points=False)
 buildSphericalCell (self, diameter, membraneThickness, point=False)
 addParticles (self, species='unknown', region='default', count=1)
 packRegion (self, region, radius, percentage, obstacleID)
 setTransitionRate (self, species, via, to, rate)
 setTwoWayTransitionRate (self, species, one, two, rate)
 buildDiffusionModel (self)
 buildReactionModel (self)
 buildSpatialModel (self)
 setWriteInterval (self, time)
 setLatticeWriteInterval (self, time)
 setSimulationTime (self, time)
 setTimestep (self, time)
 setHookInterval (self, time)
 setRandomSeed (self, seed)
 setOverflowHandler (self, mode)
 getLattice (self, force=False)
 allocateLattice (self)
 setLatticeSite (self, index, siteType)
 getLatticeSite (self, index)
 addParticleAt (self, index, particleType)
 addParticleAtIdx (self, index, particleType)
 Create an HDF5 version of the simulation amenable for later running or stand-alone running.
 save (self, filename)
 run (self, filename, method, replicates=1, seed=None, cudaDevices=None, checkpointInterval=0)
 runMPI (self, filename, method, replicates=1, driver="mpirun", ppe=1, seed=None)
 runSolver (self, filename, solver, replicates=1, cudaDevices=None, checkpointInterval=0)
 buildVector (self, arg1, arg2, arg3)
 buildPoint (self, arg1, arg2, arg3)
 buildSphere (self, arg1, arg2, arg3)
 buildEllipse (self, arg1, arg2, arg3, arg4, arg5, arg6, arg7)
 buildDifference (self, arg1, arg2, arg3)

Public Attributes

dict siteTypes = {}
dict particleMap = {}
list customAddedParticleList = []
list species_id = []
dict regions = {}
list transitionRates = []
dict initial_counts = {}
dict parameters = {}
 continousDimensions = dimensions
 latticeSpacing = spacing
int bytesPerParticle = 1
 lm_builder = lm.LatticeBuilder(dimensions[0], dimensions[1], dimensions[2], spacing, 1, 0)
bool hasBeenDiscretized = False
 name = name
list replicates = []
str filename = ""
 lattice = self.allocateLattice()

Protected Member Functions

 _repr_html_ (self)

Detailed Description

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

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

Args:
    dimensions:
        A list of [x,y,z]
    spacing:
        Lattice spacing
    name:
        The name of the RDME simulation; default: "unnamed"
    defaultRegion:
        The name of the region that is associated with the lattice sites before any other regions are added; default:"default"

Constructor & Destructor Documentation

◆ __init__()

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

Member Function Documentation

◆ _repr_html_()

pyLM.RDME.RDMESimulation._repr_html_ ( self)
protected
Here is the call graph for this function:

◆ addCuboidRegion()

pyLM.RDME.RDMESimulation.addCuboidRegion ( self,
name,
a,
b )
Add a cuboid to the builder

Args:
    name:
        Name of the site type for this region
    a:
        tuple for the first corner in continous space
    b:
        tuple for the second corner in continous space
Here is the call graph for this function:

◆ addParticleAt()

pyLM.RDME.RDMESimulation.addParticleAt ( self,
index,
particleType )
Add a particle/ to a particular site 

Args:
    index (x,y,z):
        a list of spatial location
    specie:
        The specie type to add

◆ addParticleAtIdx()

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 

Args:
    index (x,y,z):
        a list of lattice site indices
    specie:
        The specie type to add

◆ addParticles()

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

Args:
    species:
        The species to add to the region
    region:
        The region to add particles to
    count:
        Number of particles to add (default: 1)

Returns:
    The simulation object
Here is the call graph for this function:
Here is the caller graph for this function:

◆ addRegion()

pyLM.RDME.RDMESimulation.addRegion ( self,
region )
Add a region to the simulation

Args:
    region:
        The region to add to the simulation

Returns:
    region just added
Here is the caller graph for this function:

◆ addShape()

pyLM.RDME.RDMESimulation.addShape ( self,
shape )
Add a Shape to the builder

Args:
    shape:
        The Shape to add to the builder
Here is the call graph for this function:

◆ allocateLattice()

pyLM.RDME.RDMESimulation.allocateLattice ( self)

Reimplemented in pyLM.RDME.IntRDMESimulation.

Here is the caller graph for this function:

◆ buildCapsidCell()

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

Args:
    length:
        The length of the capsule from one sphere origin to the other
        diameter The diameter of the cell
    membraneThickness:
        The thickness of the membrane
    points:
        OPTIONAL: 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
Here is the call graph for this function:

◆ buildDifference()

pyLM.RDME.RDMESimulation.buildDifference ( self,
arg1,
arg2,
arg3 )

◆ buildDiffusionModel()

pyLM.RDME.RDMESimulation.buildDiffusionModel ( self)
Return the Lattice Microbes DiffusionModel object for fine-tuning

Returns:
    Simulation diffusion model
Here is the caller graph for this function:

◆ buildEllipse()

pyLM.RDME.RDMESimulation.buildEllipse ( self,
arg1,
arg2,
arg3,
arg4,
arg5,
arg6,
arg7 )

◆ buildPoint()

pyLM.RDME.RDMESimulation.buildPoint ( self,
arg1,
arg2,
arg3 )

◆ buildReactionModel()

pyLM.RDME.RDMESimulation.buildReactionModel ( self)
Return the Lattice Microbes ReactionModel object for fine-tuning

Returns:
    The reaction model for this simulation
Here is the caller graph for this function:

◆ buildSpatialModel()

pyLM.RDME.RDMESimulation.buildSpatialModel ( self)
Return the Lattice Microbes SpatialModel object for fine-tuning

Returns:
    The spatial model (i.e. obstacles, sites, etc.) for this simulation
Here is the caller graph for this function:

◆ buildSphere()

pyLM.RDME.RDMESimulation.buildSphere ( self,
arg1,
arg2,
arg3 )

◆ buildSphericalCell()

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

Args:
    diameter:
        The diameter of the cell
    membraneThickness:
        The thickness of the membrane
    point:
        The center of the spherical cell

Returns:
    The simulation object
Here is the call graph for this function:

◆ buildVector()

pyLM.RDME.RDMESimulation.buildVector ( self,
arg1,
arg2,
arg3 )

◆ defineSpecies()

pyLM.RDME.RDMESimulation.defineSpecies ( self,
species )
Define a specie/s of particles that exist in the simulation

Args:
    species:
        A list of species to add to the simulation

Returns:
    The simulation object

◆ getLattice()

pyLM.RDME.RDMESimulation.getLattice ( self,
force = False )
Get a discretized version of the simulation domain.  Call this after building all spherical and capsule cells

Args:
    force:
        Force a rediscretization of the lattice. This will throw away any additional changes made after the first call.

Returns:
    A lattice object.  This function should usually only be called once.
Here is the caller graph for this function:

◆ getLatticeSite()

pyLM.RDME.RDMESimulation.getLatticeSite ( self,
index )
Get a particular lattice site

Args:
    index (x,y,z):
        a list of coordinates

Returns:
    The type of the lattice site (string)

◆ modifyRegion()

pyLM.RDME.RDMESimulation.modifyRegion ( self,
region )
Return a pointer to a region so that it may be modified

Args:
    region:
        Get a region that is attached to the simulation for modification

Returns:
    The region to modify

◆ packRegion()

pyLM.RDME.RDMESimulation.packRegion ( self,
region,
radius,
percentage,
obstacleID )
Add nonmoving obstacles to a particular region

Args:
    region:
        The name of the region in which to add particles to
    radius:
        The radius of the particles
    percentage:
        The percentage of the total region volume that should be packed
    obstacleID:
        an identifier for the obstacle

Returns:
    The simulation object

◆ run()

pyLM.RDME.RDMESimulation.run ( self,
filename,
method,
replicates = 1,
seed = None,
cudaDevices = None,
checkpointInterval = 0 )
Run the simulation using the specified solver (e.g.  NextSubvolume, MultiParticleDiffusion, etc.) for the specified amount of time

Args:
    filename:
        The HDF file to write to
    method:
        The class name for the solver to use (e.g., lm::rdme::MpdRdmeSolver")
    replicates:
        The number of replicates to serially run
    seed:
        A seed for the random number generator to use when running the simulation; None indicates default
Here is the call graph for this function:

◆ runMPI()

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

Args:
    filename:
        The HDF file to write to
    method:
        The class name for the solver to use (e.g., lm::cme::GillespieDSolver")
    replicates:
        The number of replicates to serially run
    driver:
        The program to execute the parallel run, e.g. "mpirun", "aprun", "ibrun", etc.
    ppe:
        The number of processing elements to use  (e.g. number of nodes)
    seed:
        A seed for the random number generator to use when running the simulation; None indicates default

◆ runSolver()

pyLM.RDME.RDMESimulation.runSolver ( self,
filename,
solver,
replicates = 1,
cudaDevices = None,
checkpointInterval = 0 )
Run a simulation with a given solver

Args:
    filename:
        The HDF file to write to
    solver:
        An MESolver object
    replicates:
        The number of replicates to serially run
Here is the call graph for this function:

◆ save()

pyLM.RDME.RDMESimulation.save ( self,
filename )
Create an HDF5 version of the simulation amenable for later running or stand-alone running

Args:
    filename:
        A file to write the simulation to
Here is the call graph for this function:

◆ setHookInterval()

pyLM.RDME.RDMESimulation.setHookInterval ( self,
time )
Set the hook interval

Args:
    time:
        The interval to call hookSimulation

◆ setLatticeSite()

pyLM.RDME.RDMESimulation.setLatticeSite ( self,
index,
siteType )
Set a particular lattice site type

Args:
    index (x,y,z):
        A list of coordinates
    siteType:
        The type to set the lattice point to. This would be the name of a region that
        has previously been performed

◆ setLatticeWriteInterval()

pyLM.RDME.RDMESimulation.setLatticeWriteInterval ( self,
time )
Set the time interval to write latticedata at during simulation

Args:
     time:
        The length of time between lattice data writes

◆ setOverflowHandler()

pyLM.RDME.RDMESimulation.setOverflowHandler ( self,
mode )

◆ setRandomSeed()

pyLM.RDME.RDMESimulation.setRandomSeed ( self,
seed )
Set a known random seed

Args:
    seed:
        The seed value

◆ setSimulationTime()

pyLM.RDME.RDMESimulation.setSimulationTime ( self,
time )
Set the simulation time length

Args:
    time:
        The length of simulation time 

◆ setTimestep()

pyLM.RDME.RDMESimulation.setTimestep ( self,
time )
Set the simulation time step

Args:
    time:
        The length of simulation timestep for RDME

◆ setTransitionRate()

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

Args:
    species:
        The specie that can transition between regions
    via:
        From this region
    to:
        To this region
    rate:
        Diffusion rate between regions

Returns:
    The simulation object
Here is the caller graph for this function:

◆ setTwoWayTransitionRate()

pyLM.RDME.RDMESimulation.setTwoWayTransitionRate ( self,
species,
one,
two,
rate )
Specify the diffusion rate between species; this is a two directional rate 

Args:
    species:
        The specie that can transition between regions
    one:
        A region
    two:
        The other region
    rate:
        Diffusion rate between regions

Returns:
    The simulation object
Here is the call graph for this function:

◆ setWriteInterval()

pyLM.RDME.RDMESimulation.setWriteInterval ( self,
time )
Set the time interval to write data at during simulation

Args:
    time
        The length of time between data writes

◆ siteVolume()

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

◆ bytesPerParticle

pyLM.RDME.RDMESimulation.bytesPerParticle = 1

◆ continousDimensions

pyLM.RDME.RDMESimulation.continousDimensions = dimensions

◆ customAddedParticleList

pyLM.RDME.RDMESimulation.customAddedParticleList = []

◆ filename

str pyLM.RDME.RDMESimulation.filename = ""

◆ hasBeenDiscretized

pyLM.RDME.RDMESimulation.hasBeenDiscretized = False

◆ initial_counts

dict pyLM.RDME.RDMESimulation.initial_counts = {}

◆ lattice

pyLM.RDME.RDMESimulation.lattice = self.allocateLattice()

◆ latticeSpacing

pyLM.RDME.RDMESimulation.latticeSpacing = spacing

◆ lm_builder

pyLM.RDME.RDMESimulation.lm_builder = lm.LatticeBuilder(dimensions[0], dimensions[1], dimensions[2], spacing, 1, 0)

◆ name

pyLM.RDME.RDMESimulation.name = name

◆ parameters

dict pyLM.RDME.RDMESimulation.parameters = {}

◆ particleMap

dict pyLM.RDME.RDMESimulation.particleMap = {}

◆ regions

pyLM.RDME.RDMESimulation.regions = {}

◆ replicates

pyLM.RDME.RDMESimulation.replicates = []

◆ siteTypes

dict pyLM.RDME.RDMESimulation.siteTypes = {}

◆ species_id

pyLM.RDME.RDMESimulation.species_id = []

◆ transitionRates

pyLM.RDME.RDMESimulation.transitionRates = []

The documentation for this class was generated from the following file:
  • /data2/LM_zls_github/Lattice-Microbes/src/pylm/pyLM/RDME.py