| Lattice Microbes 2.5
    This is for whole cell modeling | 
| Public Member Functions | |
| simulationTime (self) | |
| simulationTime (self, val) | |
| speciesWriteInterval (self) | |
| speciesWriteInterval (self, val) | |
| latticeWriteInterval (self) | |
| latticeWriteInterval (self, val) | |
| perfPrintInterval (self) | |
| perfPrintInterval (self, val) | |
| hookInterval (self) | |
| hookInterval (self, val) | |
| __init__ (self, name, filename, dimensions, latticeSpacing, latticeType=None) | |
| resizeLattice (self, dimensions, latticeSpacing, latticeType=None) | |
| placeNumber (self, sp, x, y, z, n) | |
| distributeNumber (self, sp, reg, count) | |
| distributeConcentration (self, sp, reg, conc) | |
| transitionRate (self, sp, rFrom, rTo, rate, value=None) | |
| assignReaction (self, reaction, region) | |
| species (self, name, **kwargs) | |
| region (self, name, **kwargs) | |
| reaction (self, reactants, products, rate, value=None, regions=None, **kwargs) | |
| rateConst (self, rate, value, order, **kwargs) | |
| diffusionConst (self, rate, value, **kwargs) | |
| diffusionZero (self) | |
| maxDiffusionRate (self, latticeSpacing=None, dt=None) | |
| diffusionFast (self) | |
| setMaximumTimestep (self) | |
| run (self, solver=None, replicate=1, seed=None, cudaDevices=None, checkpointInterval=0, sample_frame=False, max_frames=100) | |
| Public Attributes | |
| name = name | |
| filename = filename | |
| latticeType = latticeType | |
| float | NA = 6.02214085774e23 | 
| int | siteV = 1000 * latticeSpacing**3 | 
| int | siteNAV = self.siteV * self.NA | 
| shape = nz,ny,nx | |
| latticeSpacing = latticeSpacing | |
| pps = lm.getCompiledLatticeMaxOccupancy() | |
| lattice = lm.IntLattice(nz,ny,nx, latticeSpacing, self.pps) | |
| siteLattice = self.lattice.getSiteLatticeView() | |
| particleLattice = self.lattice.getParticleLatticeView() | |
| int | maxConcentration = self.pps/self.siteNAV | 
| diffRateList | |
| Protected Attributes | |
| int | _perfPrintInterval = 60 | 
| list | _particlePlacement = [] | 
| list | _particleDistCount = [] | 
| list | _particleDistConc = [] | 
| list | _transitionRates = [] | 
| list | _reactionLocations = [] | 
Base RDME model
| jLM.RDME.SpatialModel.__init__ | ( | self, | |
| name, | |||
| filename, | |||
| dimensions, | |||
| latticeSpacing, | |||
| latticeType = None ) | 
| jLM.RDME.SpatialModel.assignReaction | ( | self, | |
| reaction, | |||
| region ) | 
Assign a reaction to a region
Args:
    reaction (:py:class:`~jLM.Types.Reaction`):
        Reaction
    region (:py:class:`~jLM.Types.Region`):
        Target region
 | jLM.RDME.SpatialModel.diffusionConst | ( | self, | |
| rate, | |||
| value, | |||
| ** | kwargs ) | 
Lookup/create diffusion constant instance
Args:
    rate (str):
        Name of diffusion constant
    value (float):
        Value of diffusion constant
    texRepr (str):
        TeX math mode representation of the rate constant
    annotation (str):
        Description text
Returns:
    :py:class:`~jLM.Types.DiffusionConst`:
        DiffusionConst instance
 | jLM.RDME.SpatialModel.diffusionFast | ( | self | ) | 
Lookup/create a maximum-valued diffusion constant instance
Returns:
    :py:class:`~jLM.Types.DiffusionConst`: 
        DiffusionConst instance
 | jLM.RDME.SpatialModel.diffusionZero | ( | self | ) | 
Lookup/create a zero-valued diffusion constant instance
Returns:
    :py:class:`~jLM.Types.DiffusionConst`:
        DiffusionConst instance
 
| jLM.RDME.SpatialModel.distributeConcentration | ( | self, | |
| sp, | |||
| reg, | |||
| conc ) | 
Distribute a concentration of particles uniformly through a region
Args:
    sp (:py:class:`~jLM.Types.Species`]):
        Species type
    reg (:py:class:`~jLM.Types.Region`]):
        Region type
    conc (float):
        Concentration of particles
 
Reimplemented in jLM.RDMEExtended.SimExtended.
| jLM.RDME.SpatialModel.distributeNumber | ( | self, | |
| sp, | |||
| reg, | |||
| count ) | 
Distribute a number of particles uniformly through a region
Args:
    sp (:py:class:`~jLM.Types.Species`]):
        Species type
    reg (:py:class:`~jLM.Types.Region`]):
        Region type
    count (int):
        Number of particles
 
Reimplemented in jLM.RDMEExtended.SimExtended.
| jLM.RDME.SpatialModel.hookInterval | ( | self | ) | 
Simulation time between calls to hookSimulation
| jLM.RDME.SpatialModel.hookInterval | ( | self, | |
| val ) | 
| jLM.RDME.SpatialModel.latticeWriteInterval | ( | self | ) | 
Simulation time between writing the lattice to the HDF5 file
| jLM.RDME.SpatialModel.latticeWriteInterval | ( | self, | |
| val ) | 
| jLM.RDME.SpatialModel.maxDiffusionRate | ( | self, | |
| latticeSpacing = None, | |||
| dt = None ) | 
Compute max allowed diffusion constant for the simulation 
Args:
    latticeSpacing (float):
        Lattice spacing in meters
    dt (float):
        Timestep in seconds
Returns:
    float:
        Maximum diffusion constant in m^/s
 | jLM.RDME.SpatialModel.perfPrintInterval | ( | self | ) | 
Real elapsed time between writing simulation progress to stdout
| jLM.RDME.SpatialModel.perfPrintInterval | ( | self, | |
| val ) | 
| jLM.RDME.SpatialModel.placeNumber | ( | self, | |
| sp, | |||
| x, | |||
| y, | |||
| z, | |||
| n ) | 
Place a number of particles at a specific subvolume
Args:
    sp (:py:class:`~jLM.Types.Species`):
        Species type
    x (int):
        x-coordinate
    y (int):
        y-coordinate
    z (int):
        z-coordinate
    n (int):
        Number of particles
 
| jLM.RDME.SpatialModel.rateConst | ( | self, | |
| rate, | |||
| value, | |||
| order, | |||
| ** | kwargs ) | 
Lookup/create reaction rate constant instance
Args:
    rate (str):  
        Name of rate constant
    value (float): 
        Value of rate constant
    order (int): 
        Order of reaction
    texRepr (str): 
        TeX math mode representation of the rate constant
    annotation (str): 
        Description text
Returns:
    :py:class:`~jLM.Types.RateConst`:
        RateConst instance
 | jLM.RDME.SpatialModel.reaction | ( | self, | |
| reactants, | |||
| products, | |||
| rate, | |||
| value = None, | |||
| regions = None, | |||
| ** | kwargs ) | 
Lookup/create reaction instance
The product/reactants can be specified as a string, Species 
instance, list of strings, list of Species instances, or the 
empty list, denoting no species. The reaction rate can be 
created in place if a string is given for `rate` and the value
of the rate constant is specified by `value`.  Rates must be 
provided in standard chemical kinetic units, e.g. for a reaction
.. math:: \mathrm{A} + \mathrm{B} \overset{k}{\to} \mathrm{C},
the rate coefficent, :math:`k` has units of 
:math:`\mathrm{M}^{-1}\mathrm{s}^{-1}`
Args:
    reactants ([:py:class:`~jLM.Types.Species`]):
        Species, or list of species acting as reactants
    products ([:py:class:`~jLM.Types.Species`]):
        Species, or list of species acting as products
    rate ([:py:class:`~jLM.Types.RateConst`]):
        Rate constant
    value (float):
        Optional value of rate constant
    regions (py:class:`~jLM.Types.Region`]):
        List of applicable regions
Returns:
    :py:class:`~jLM.Types.Reaction`:
       Reaction instance
 | jLM.RDME.SpatialModel.region | ( | self, | |
| name, | |||
| ** | kwargs ) | 
Lookup/create region instance
Args:
    name (str):
        Name of the region
    texRepr (str):
        TeX math mode representation of the region
    annotation (str):
        Description text
Returns:
    :py:class:`~jLM.Types.Region`:
       Region instance
 | jLM.RDME.SpatialModel.resizeLattice | ( | self, | |
| dimensions, | |||
| latticeSpacing, | |||
| latticeType = None ) | 
| jLM.RDME.SpatialModel.run | ( | self, | |
| solver = None, | |||
| replicate = 1, | |||
| seed = None, | |||
| cudaDevices = None, | |||
| checkpointInterval = 0, | |||
| sample_frame = False, | |||
| max_frames = 100 ) | 
Run the RDME simulation
Args:
    solver (:py:class:`lm.RDMESolver`):
        Rdme solver
    replicate (int):
        Index of replicate
    seed (int):
        RNG seed
    cudaDevices ([int]):
        List of CUDA device indexes
    checkpointInterval (int):
        Number of seconds between checkpoints
    sample_frame (bool):
        Sample frames or not
    max_frames (int):
        Number of frames after sampling
Returns:
    :py:class:`jLM.File`:
        Simulation result
 
Reimplemented in jLM.RDMEExtended.SimExtended.
| jLM.RDME.SpatialModel.setMaximumTimestep | ( | self | ) | 
Set the simulation timestep using the fastest diffusion rate
| jLM.RDME.SpatialModel.simulationTime | ( | self | ) | 
Total simulated time
| jLM.RDME.SpatialModel.simulationTime | ( | self, | |
| val ) | 
| jLM.RDME.SpatialModel.species | ( | self, | |
| name, | |||
| ** | kwargs ) | 
Lookup/create species instance
Args:
    name (str):
        Name of the species
    texRepr (str):
        TeX math mode representation of the species
    annotation (str):
        Description text
Returns:
    :py:class:`~jLM.Types.Species`: 
        Species instance
 
Reimplemented in jLM.RDMEExtended.SimExtended.
| jLM.RDME.SpatialModel.speciesWriteInterval | ( | self | ) | 
Simulation time between writing the species counts to the HDF5 file
| jLM.RDME.SpatialModel.speciesWriteInterval | ( | self, | |
| val ) | 
| jLM.RDME.SpatialModel.transitionRate | ( | self, | |
| sp, | |||
| rFrom, | |||
| rTo, | |||
| rate, | |||
| value = None ) | 
Define the diffusive transition rate between regions
This allows for fine grained detail over the diffusion matrix.
The transition is defined for the species `sp` from region
`rFrom` to `rTo`. If `None` is given for any of these
parameters, then the entire axis of the matrix is affected,
i.e. `None` is a wildcard. The rate is specified with the
`rate` parameter, either as a DiffusionConst instance, or as
a string to lookup the previously created DiffusionConst. If a
string is provided and a DiffusionConst instance has not been
created, providing the parameter `value` will create the new
constant.
Args:
    sp (:py:class:`~jLM.Types.Species`]): 
        Species type or None
    rFrom (:py:class:`~jLM.Types.Region`]):
        Region type or None
    rTo (:py:class:`~jLM.Types.Region`]):
        Region type or None
    rate (:py:class:`~jLM.Types.DiffusionConst`]):
        Diffusion rate
    value (float):
        Value of diffusion constant if new
 | 
 | protected | 
| 
 | protected | 
| 
 | protected | 
| 
 | protected | 
| 
 | protected | 
| 
 | protected | 
| jLM.RDME.SpatialModel.diffRateList | 
| jLM.RDME.SpatialModel.filename = filename | 
| jLM.RDME.SpatialModel.lattice = lm.IntLattice(nz,ny,nx, latticeSpacing, self.pps) | 
| jLM.RDME.SpatialModel.latticeSpacing = latticeSpacing | 
| jLM.RDME.SpatialModel.latticeType = latticeType | 
| int jLM.RDME.SpatialModel.maxConcentration = self.pps/self.siteNAV | 
| float jLM.RDME.SpatialModel.NA = 6.02214085774e23 | 
| jLM.RDME.SpatialModel.name = name | 
| jLM.RDME.SpatialModel.particleLattice = self.lattice.getParticleLatticeView() | 
| jLM.RDME.SpatialModel.pps = lm.getCompiledLatticeMaxOccupancy() | 
| jLM.RDME.SpatialModel.shape = nz,ny,nx | 
| jLM.RDME.SpatialModel.siteLattice = self.lattice.getSiteLatticeView() | 
| int jLM.RDME.SpatialModel.siteNAV = self.siteV * self.NA | 
| int jLM.RDME.SpatialModel.siteV = 1000 * latticeSpacing**3 |