distribute_vol_util

Short description

@author: Moritz F P Becker

pyrid.system.distribute_vol_util.add_molecule_single(Compartment_id, System, RBs, Particles, pos, quaternion, mol_type)[source]

Places a single new molecule into a compartment, given its position, orientation and type.

Parameters
Compartment_idint64

id of the compartment

Systemobject

Instance of System class

RBsobject

Instance of RBs class

Particlesobject

Instance of Particles class

posarray_like

Position of the molecul center

quaternionarray_like

Rotation quaternion of the molecule

mol_typestring

Molecule name

pyrid.system.distribute_vol_util.add_molecules(Compartment_id, System, RBs, Particles, pos, quaternion, mol_type_ids, mol_types)[source]

Places new molecules into a compartment, given their positions, orientations and types.

Parameters
Compartment_idint64

id of the compartment

Systemobject

Instance of System class

RBsobject

Instance of RBs class

Particlesobject

Instance of Particles class

posarray_like

Positions of the molecule centers

quaternionarray_like

Rotation quaternions of the molecules

mol_type_idsarray_like

Type ids of the molecules

mol_typesarray_like

Array of strings of the molecule names

pyrid.system.distribute_vol_util.mc(System, Compartment, mol_types, N)[source]

Does Monte Carlo sampling for N molecules inside a compartment returning vectors for the molecule positions, types and orientations.

Parameters
Compartmentobject

Instance of Compartment class

mol_typeslist of strings

List of molecule types to distribute.

Nint64[:]

Total number of molecules to distribute per molecule type.

Returns
tuple(float64[:,3], int64[:], float64[:,4])

molecule positions, molecule types, molecule orientations in quaternion representation

pyrid.system.distribute_vol_util.monte_carlo_distribution_3D(Compartment, System, N_total)[source]

Does Monte Carlo sampling for N points inside a compartment.

Parameters
Compartmentobject

Instance of Compartment class

Systemobject

Instance of System class

N_totalint64

Total number of points to distribute

Returns
float64[:,3]

Array of positions

pyrid.system.distribute_vol_util.normal(origin, jitter, number, mol_types, System)[source]

Distributes molecules spherically around a point by a gaussian distribution.

Parameters
originfloat64[3]

Origin of the distribution

jitterfloat64

Standard deviation of the normal distribution

numberint64[:]

Number of molecules to distribute per type

mol_typesarray_like

Array of strings of molecule names

Systemobject

Instance of System class

Returns
tuple(float64[N,3], float64[N,4], int64[N])

Positions, Quaternions, and Molecule ids

pyrid.system.distribute_vol_util.pds(Compartment, System, mol_types, number, clustering_factor=1, max_trials=100)[source]

Distributes a number of molecules of different types inside the volume of a compartment using a poisson disc sampling algorithm.

Parameters
Compartmentobject

Instance of Compartment class

Systemobject

Instance of System class

mol_typesarray_like

Array of strings of molecule names

numberfloat64[:]

Number of molecules to distribute per given type

clustering factorint64

Determines how much the molecules are spread inside the compartment. Higher values are necessary for uniform distributions at low density. Default = 1

max_trialsint64

Maximum number of trials to find a new position around a given active point in the poisson disc sampling algorithm.

Returns
tuple(float64[:,3], int64[:], float64[:,4])

Returns arrays of the positions, molecule types, and rotation quaternions of the distributed moelcules.

pyrid.system.distribute_vol_util.pds_uniform(Compartment, System, mol_types, N, multiplier=100)[source]

Uniformly distributes molecules of different types inside a compartment using a blue noise/poisson disc sampling algorithm originaly introduced by Corsini et al. [25] to distribute points on triangulated mesh surfaces but that has been adapted here for volume distributions. PyRID also uses a variant of the original algorithm to distribute molecules on mesh surfaces.

Parameters
Compartmentobject

Instance of Compartment class

Systemobject

Instance of System class

mol_typesarray_like

Array of strings of molecule names

Nfloat64[:]

Number of molecules to distribute per given type

multiplierint64

Specifies the multiplication factor by which the algorithm upscales the target number for the Monte Carlo sampling. Default = 100

Returns
tuple(float64[:,3], int64[:], float64[:,4])

Returns arrays of the positions, molecule types, and rotation quaternions of the distributed moelcules.

pyrid.system.distribute_vol_util.point_in_sphere_simple(radius)[source]

Returns a random point uniformly distributed in the volume of a sphere with radius r.

Parameters
radiusfloat
Returns
tuple(float, float, float)

Uniformly distributed random point in sphere volume

pyrid.system.distribute_vol_util.point_on_sphere(radius)[source]

Returns a random point uniformly distributed on the surface of a sphere with radius r.

Parameters
radiusfloat
Returns
tuple(float, float, float)

Uniformly distributed random point on sphere surface

pyrid.system.distribute_vol_util.poisson_disc_sampling(Compartment, System, radii, mol_type_ids, N, weights, clustering_factor, max_trials)[source]

Poisson disc sampling algorithm for different sized spheres, accounting for different boundary conditions and mesh compartment boundaries. The algorithm is based on Bridson [27].

Parameters
Compartmentobject

Instance of Compartment class

Systemobject

Instance of System class

radiifloat64[:]

Radii of the different sphere types

mol_type_idsarray_like

Type ids of the different spheres

Nfloat64[:]

Number of spheres to distribute per given type

clustering factorint64

Determines how much the spheres are spread inside the compartment. Higher values are necessary for uniform distributions at low density. Default = 1

max_trialsint64

Maximum number of trials to find a new position around a given active point in the poisson disc sampling algorithm.

Returns
tuple(float64[:,3], int64[:], float64[:,4], int64)

Returns arrays of the positions, types and rotations quaternions. Also return the total count of distributed spheres.

Raises
ValueError(‘There is no overlap of the Compartment with the Simulation box! Unable to distribute points in Compartment!’)

Raised if the selected compartment does not overlap with the simulation box.

pyrid.system.distribute_vol_util.poisson_disc_sampling_uniform(Compartment, radii, mol_type_ids, weights, N, System, multiplier=100)[source]

Uniformly distributes spheres of different radii inside a compartment using a blue noise/poisson disc sampling algorithm originaly introduced by Corsini et al. [25] to distribute points on triangulated mesh surfaces but that has been adapted here for volume distributions. PyRID also uses a variant of the original algorithm to distribute molecules on mesh surfaces.

Parameters
Compartmentobject

Instance of Compartment class

radiiint64[:]

radii of the different sphere types

mol_type_idsarray_like

Array of molecule type ids

weightsfloat64[:]

Weights specifying how often the different sphere types are selected for distribution. Weights should be higher for large spheres since the probability of a successfull trial is lower than for small spheres.

Nfloat64[:]

Number of molecules to distribute per given type

Systemobject

Instance of System class

multiplierint64

Specifies the multiplication factor by which the algorithm upscales the target number for the Monte Carlo sampling. Default = 100

Returns
tuple(float64[:,3], int64[:], float64[:,4], int64)

Returns arrays of the positions, types and rotations quaternions. Also return the total count of distributed spheres.

Raises
ValueError(‘There is no overlap of the Compartment with the Simulation box! Unable to distribute points in Compartment!’)

Raised if the selected compartment does not overlap with the simulation box.

Notes

The algorithm is based on Corsini et al. 2012, Efficient and Flexible Sampling with Blue Noise Properties of Triangular Meshes.

The algorithm consists of 10 steps:

  1. Generate a sample pool S using monte_carlo_distribution_3D().

  2. Divide space into cells and count the number of samples in each cell.

  3. Randomly select a cell weighted by the number of active samples in each cell (active sample: sample that is not yet occupied or deleted).

  4. Randomly select a sample from the selected cell.

  5. Randomly choose a particle type of radius Ri (weighted by the relative number of each type that we want to distribute).

  6. Check whether the distance of the selected sample to the neighboring samples that are already occupied is larger or equal to Ri+Rj.

  7. If True, accept the sample and add the molecule type and position to an occupied sample list. Next, delete all other samples within radius Ri, as these wont ever become occupied anyway.

  8. Update the number count of samples for the current cell.

  9. While the desired number of molecules is not reached, return to 3. However, set a maximum number of trials.

  10. If their are no active samples left before we reach the desired molecule number and the maximum number of trials, generate a new sample pool.

pyrid.system.distribute_vol_util.random_direction_Halfsphere(sphere_radius, normal)[source]

Returns a random direction vector uniformly distributed in the volume of a half-sphere with radius r, splitted in half in direction of a given plane normal vector.

Parameters
sphere_radiusfloat64

Radius of the sphere

normalfloat64[3]

Normal vector of the plane splitting the sphere in half

Returns
float64[3]

Direction vector

Notes

First a uniformly distributed random vector \(d\boldsymbol{X}\) that sites somewhere within the full sphere is drawn. Next, we test whether the vector points into the same direction as the plane normal vector \(\hat{\boldsymbol{n}}\) via their dot product. If the dot product is negative, the direction vector is reflected at the plane by \(d\boldsymbol{X}_{refl} = d\boldsymbol{X} -2 d\boldsymbol{X} \cdot \hat{\boldsymbol{n}}\).

pyrid.system.distribute_vol_util.random_direction_sphere(sphere_radius)[source]

Returns a random direction vector uniformly distributed in the volume of a sphere with radius r.

Parameters
sphere_radiusfloat64

Radius of the sphere

Returns
float64[3]

Direction vector

pyrid.system.distribute_vol_util.random_quaternion()[source]

Returns a unit quaternion representing a uniformly distributed random rotation.

Returns
float64[4]

Uniformly distributed, random rotation quaternion

Notes

The method can be found in K. Shoemake, Graphics Gems III by D. Kirk, pages 124-132 [28] and on http://planning.cs.uiuc.edu/node198.html:

\[\boldsymbol{q} = [\sqrt{1-u_1} \, \sin(2 \pi u_2), \sqrt{1-u_1} \, \cos(2 \pi u_2), \sqrt{u_1} \, \sin(2 \pi u_3), \sqrt{u_1} \, \cos(2 \pi u_3)] ,\]

where \(u_1,u_2,u_3 \in [0,1]\) are uniformly distributed random numbers. \(\boldsymbol{q}\) is a uniformly distributed, random rotation quaternion.

pyrid.system.distribute_vol_util.random_quaternion_tuple()[source]

Returns a unit quaternion representing a uniformly distributed random rotation.

Returns
tuple[4]

Uniformly distributed, random rotation quaternion

Notes

The method can be found in K. Shoemake, Graphics Gems III by D. Kirk, pages 124-132 [28] and on http://planning.cs.uiuc.edu/node198.html:

\[\boldsymbol{q} = [\sqrt{1-u_1} \, \sin(2 \pi u_2), \sqrt{1-u_1} \, \cos(2 \pi u_2), \sqrt{u_1} \, \sin(2 \pi u_3), \sqrt{u_1} \, \cos(2 \pi u_3)] ,\]

where \(u_1,u_2,u_3 \in [0,1]\) are uniformly distributed random numbers. \(\boldsymbol{q}\) is a uniformly distributed, random rotation quaternion.

pyrid.system.distribute_vol_util.release_molecules_boundary(System, RBs, Particles)[source]

Releases molecules into the volume at the simulation box boundary given an outside concentration of the different molecule types.

Parameters
Systemobject

Instance of System class

RBsobject

Instance of RBs class

Particlesobject

Notes

The average number of molecules that hit a boundary of area \(A\) from one side within a time step \(\Delta t\) can be calculated from the molecule concentration \(C\) and the average distance a diffusing molecule travels normal to a plane \(l_{n}\) within \(\Delta t\) [10]:

\[N_{hits} = \frac{A l_{n}}{2 C},\]

where

\[l_{n} = \sqrt{\frac{4D\Delta t}{\pi}}.\]

Here \(D = Tr(\boldsymbol{D}^{tt,b})/3\) is the scalar translational diffusion constant. The boundary crossing of molecules can be described as a poisson process. As such, the number of molecules that cross the boundary each time step is drawn from a poisson distribution with a rate \(N_{hits}\).

The normalized distance that a crossing molecule ends up away from the plane/boundary follows the following distribution [10]:

\[P(d\tilde{x}) = 1-e^{-d\tilde{x}^2}+\sqrt{\pi}*dx*\text{erfc}(d\tilde{x})\]

The distance vector normal to the plane after the crossing can then be calculated from the diffusion length constant \(\lambda\) and the plane’s normal vector \(\hat{\boldsymbol{n}}\) by \(d\boldsymbol{x} = \lambda \, d\tilde{x} \, \hat{\boldsymbol{n}} = \sqrt{4Dt} \, d\tilde{x} \, \hat{\boldsymbol{n}}\).

In the case that a molecule enters the simulation box near to another boundary, e.g. of a mesh compartment, we may also want to account for the distance traveled parallel to the plane in order to correctly resolve collision with the mesh. However, currently PyRID does not account for this. For small intregration time steps and meshes that are further than \(\sqrt{4Dt}\) away from the simulation box border, the error introduced should however be negligable.

Now that the number of molecules and their distance away from the plane are determined, the molecules are distributed in the simualtion box. Since the diffusion along each dimension is independent we can simply pick a random point uniformly distributed on the respective plane. For triangulated mesh surfaces, triangles are picked randomly, weighted by their area. Sampling a uniformly distributed random point in a triangle is done by [26]

\[P(\boldsymbol{r}) = (1-\sqrt{\mu_1})*\boldsymbol{p}_0+(\sqrt{\mu_1}*(1-\mu_2))*\boldsymbol{p}_1+(\mu_2*\sqrt{\mu_1})*\boldsymbol{p}_2 ,\]

where \(\mu_1, \mu_2\) are random numbers between 0 and 1. \(\boldsymbol{p}_0, \boldsymbol{p}_1, \boldsymbol{p}_2\) are the three vertices of the triangle.

Note that, in general, here we neglect any interactions between the virtual molecules. Therefore, fixed concentration boundary conditions only result in the same inside and outside concentrations if no molecular in teractions are simulated.

pyrid.system.distribute_vol_util.trace_direction_vector(dX, origin, pos, System)[source]

Traces the path along a direction vector, and, depending on the boundary conditions and the presence of mesh compartments updates the direction vector. If the direction vector hits an absorptive boundary (fixed concentration boundary conditions), returns False, indicating that the corresponding molecule needs to be deleted.

Parameters
dXfloat64[3]

Direction vector

originfloat64[3]

Origin of the molecule

posfloat64[3]

New position/origin of the moelcule is saved in this array.

Systemobject

Instance of System class

Returns
bool

Indicates whether the position tracing was successful or the molecule hit an absorptive boundary.