distribute_surface_util
Short description
@author: Moritz F P Becker
- pyrid.system.distribute_surface_util.add_molecule_single(Compartment_id, System, RBs, Particles, pos, quaternion, triangle_id, mol_type)[source]
Places a single new molecule into a compartment, given its position, orientation, triangle id 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
- triangle_idsarray_like
Triangle ids of the molecules
- mol_typestring
Molecule name
- pyrid.system.distribute_surface_util.add_molecules(comp_id, System, RBs, Particles, pos, quaternion, triangle_ids, mol_type_ids, Types)[source]
Places new molecules on teh surface of a compartment, given their positions, orientations, triangle ids 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
- triangle_idsarray_like
Triangle ids of the molecules
- mol_type_idsarray_like
Type ids of the molecules
- Typesarray_like
Array of strings of the molecule names
- pyrid.system.distribute_surface_util.evenly_on_sphere(n, r)[source]
Approximately distributes n points evenly on the sphere surface.
- Parameters
- nint64
Number of points.
- `r`float64
Sphere radius.
- Returns
- float64[n,3]
Array of n points that are approximatels evenly distributed on a sphere surface.
Notes
To distribute an arbitrary number of points evenly on a sphere surface, maximizing the minimum distance between all points, is not easily solved and only for a few cases exact theoretical solutions exist. For the cases n=1 to n=6 the known solutions to the Thompson problem are used, which is closesly related to the problem of distributing n points evenly on a sphere surface. N=1 and N=2 are trivial, for the other cases we find
N = 3: points reside at vertices of equilateral triangle
N = 4: points reside at the vertices of a regular tetrahedron.
N = 5: points reside at vertices of a triangular dipyramid.
N = 6: points reside at vertices of a regular octahedron.
For n>6, points an approximation using the Fibonacci lattice is used (see e.g. http://extremelearning.com.au/how-to-evenly-distribute-points-on-a-sphere-more-effectively-than-the-canonical-fibonacci-lattice/#more-3069, https://stackoverflow.com/questions/9600801/evenly-distributing-n-points-on-a-sphere). Thereby
\[\begin{split}\begin{align*} x_i = R (\cos(\theta_i) \sin(\phi_i)) \\ y_i = R (\sin(\theta_i) \sin(\phi_i)) \\ z_i = R \cos(\phi_i) \\ \end{align*}\end{split}\]where
\[\begin{split}\begin{align*} G = \frac{1 + \sqrt{5}}{2} \\ \phi = \arccos \Big(1 - \frac{2 i+1}{n} \Big) \\ \theta = \frac{2 \pi i}{G} \end{align*}\end{split}\]Here, \(G\) is the so called golden ratio, \(i = 1,2,...,n\), and \(R\) is the sphere radius.
- pyrid.system.distribute_surface_util.mc(System, Compartment, mol_types, N, face_group=None)[source]
Does Monte Carlo sampling for N molecules on the surface of a triangulated mesh compartment returning vectors for the molecule positions, types and orientations and triangle indices.
- 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.
- face_groupstring
Name of the face group on whose triangles to distribute the molecules.
- Returns
- tuple(float64[:,3], int64[:], float64[:,4], int64[:])
molecule positions, molecule types, molecule orientations in quaternion representation
- pyrid.system.distribute_surface_util.monte_carlo_distribution_2D(Compartment, System, N_total, facegroup=None)[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
- facegroupstring
Name of the facegroup to which to restrict the molecule distribution. Default = None
- Returns
- float64[:,3]
Array of positions
- pyrid.system.distribute_surface_util.normal(triangle_id, jitter, number, mol_types, System, center=None)[source]
Distributes molecules spherically around a point on a mesh surface by a gaussian distribution and using a ray marching algorithm.
- Parameters
- triangle_idint64
triangle id around whose centroid the molecules are distributed
- 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_surface_util.pds(Compartment, System, mol_types, number, facegroup=None, multiplier=50)[source]
Uniformly distributes molecules of different types on the mesh surface of a compartment using a blue noise/poisson disc sampling algorithm originaly introduced by Corsini et al. [25].
- 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
- facegroupstring
Name of the facegroup to which to restrict the molecule distribution. Default = None
- multiplierint64
Specifies the multiplication factor by which the algorithm upscales the target number for the Monte Carlo sampling. Default = 50
- Returns
- tuple(float64[:,3], int64[:], float64[:,4])
Returns arrays of the positions, molecule types, and rotation quaternions of the distributed moelcules.
- pyrid.system.distribute_surface_util.point_in_disc(triangle_id, radius, System)[source]
Returns a random point uniformly distributed inside a disc with radius r that lies within the plane of given mesh triangle.
- Parameters
- triangle_idint64
Index of the triangle from where the point is distributed
- radiusfloat
Radius of the disc
- Systemobject
Instance of System class
- Returns
- float64[3]
Random point uniformly distributed in a disc.
- pyrid.system.distribute_surface_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_surface_util.poisson_disc_sampling_2D(Compartment, radii, mol_type_ids, weights, N, System, facegroup=None, multiplier=50)[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
- facegroupstring
Name of the facegroup to which to restrict the molecule distribution. Default = None
- 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:
Generate a sample pool S using monte_carlo_distribution_3D().
Divide space into cells and count the number of samples in each cell.
Randomly select a cell weighted by the number of active samples in each cell (active sample: sample that is not yet occupied or deleted).
Randomly select a sample from the selected cell.
Randomly choose a particle type of radius Ri (weighted by the relative number of each type that we want to distribute).
Check whether the distance of the selected sample to the neighboring samples that are already occupied is larger or equal to Ri+Rj.
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.
Update the number count of samples for the current cell.
While the desired number of molecules is not reached, return to 3. However, set a maximum number of trials.
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_surface_util.random_point_in_triangle(p0, p1, p2)[source]
Returns a random point uniformly distributed in a triangle.
- Parameters
- p0float64[3]
Triangle vertex 1
- p1float64[3]
Triangle vertex 2
- p2float64[3]
Triangle vertex 3
- Returns
- float64[3]
Uniformly, random point in triangle
Notes
Random points, uniformly distributed inside a triangle can be sampled using an algorithm introduced by Osada et al. [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.
- pyrid.system.distribute_surface_util.random_point_on_edge(p0, p1)[source]
Returns a random point uniformly distributed on a line segment/edge between two vertices.
- Parameters
- p0float64[3]
Vertex 1
- p1float64[3]
Vertex 2
- Returns
- float64[3]
Random point on edge
- pyrid.system.distribute_surface_util.release_molecules_boundary_2d(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 edge length \(L\) from one side within a time step \(\Delta t\) can be calculated from the molecule surface concentration \(C\) and the average distance a diffusing molecule travels normal to an line \(l_{n}\) within \(\Delta t\) [10]:
\[N_{hits} = \frac{L l_{n}}{2 C},\]where
\[l_{n} = \sqrt{\frac{4D\Delta t}{\pi}}.\]Here \(D = Tr(\boldsymbol{D}_{xy}^{tt,b})/2\) 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 edges’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}}\).
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 edge:
\[P(\boldsymbol{r}) = p_0+\mu*(p_1-p_0) ,\]where \(\mu\) is a random number between 0 and 1. \(\boldsymbol{p}_0, \boldsymbol{p}_1\) are the two vertices of the edge.
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_surface_util.trace_direction_vector(dX, origin, triangle_id, 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) triangle id -1 is returned.
- Parameters
- dXfloat64[3]
Direction vector
- originfloat64[3]
Origin of the molecule
- triangle_idint64
Triangle id of the starting position/origin.
- Systemobject
Instance of System class
- Returns
- tuple(float64[3], float64[4], int64)
Returns position, quaternion, and triangle_id of the new position. If path hit an absorptive boundary triangle_id = -1 is returned.