Full API Documentation¶
Functions of Trace Module
-
class
orsaytrace.trace.Simu(x, y, z, res)¶ Bases:
objectA simulation class that will allow one to create sereral objects and perform spacial transformations, such as rotations. This class also contains a few convenient functions in order to display data easily to the user.
- Parameters
x (float) – A float representing cell size in x
y (float) – A float representing cell size in y
z (float) – A float representing cell size in z
res (float) – A float for simulation resolution.
Notes
Other attributes of simulation available to user and used during simulation run:
- ss: int
A subsampling factor. Default is 2.01 and it says meshing will be done using sub pixels resolution.
- grid: array_like
A 3 dimensional array for the number of grid points.
- index: numpy array
A numpy array initialized as numpy.ones(grid). This means initial cell is under vacuum for all wavelengths.
- normal: numpy array
A numpy array initialized as numpy.asarray([numpy.zeros(grid), numpy.zeros(grid), numpy.zeros(grid)]) to all surface normals. Initial cell has no valid surface normals because there is no surfaces.
- photons: list
A list containing photon objects. Those photons will propagate and be saved along the simulation.
- photon_lists: numpy array
A numpy array cointaining photon_lists objects. photons will be saved here. Array length is given by the number of analyzes plans.
-
assign_block_n(min_index, max_index, value)¶ Assigns index of refraction for multiple points at once
- Parameters
min_index (array_like) – A three dimensional array with grid position (in pixels) for the botton left point.
max_index (array_like) – A three dimensional array with grid position (in pixels) for the top right point.
normal (array_like) – A three dimensiona array with the correspondent surface normal.
-
assign_block_normal(min_index, max_index, normal)¶ Assigns normal for multiple points at once.
- Parameters
min_index (array_like) – A three dimensional array with grid position (in pixels) for the botton left point.
max_index (array_like) – A three dimensional array with grid position (in pixels) for the top right point.
normal (array_like) – A three dimensiona array with the correspondent surface normal.
-
assign_n(index, value)¶ Assigns index of refraction for a given point
- Parameters
index (array_like) – A three dimensional array with grid position (in pixels).
value (float) – The desired index of refraction
-
assign_normal(index, normal)¶ Assigns a surface normal in a given index
- Parameters
index (array_like) – A three dimensional array with grid position (in pixels).
normal (array_like) – A three dimensional array with the correspondent surface normal.
-
check_analysis_plan(photon)¶ Check if photon is in a created analysis plan.
- Parameters
photon (class.photon)
rindex (int) – Run index. Used for multi process simulations.
-
create_analysis_plan(normal, value, **kargs)¶ Create an analysis plan.
- Parameters
normal (array_like) – A 3 dimensional array of plan normal vector.
value (float) – Indenpendent value of a plan equation.
**kargs (dict) – A dictionary that will be used as a condition dictionary. Photons will be added the plan only if satisfies this condition.
Notes
Condition must be an attribute of photon. It can take two forms: a single number or a tuple. Single number represents minimal values, while tuple represents a possible interval
Examples
>>> create_analysis_plan([0, 1, 0], a, reflection_count = 1)
Example above creates a analyses plan with equation
\[y=a\]and with a condition that reflection_count is higher than 1 (at least one reflection).
>>> create_analysis_plan([1, 2, 3], 3, reflection_count = (1, 1))
Example above creates a analyses plan with equation
\[(x + y + z) / \sqrt{14} = 3\]and with a condition that reflection_count is exactly 1. Factor square root of 14 comes from vector normalization
>>> create_analysis_plan([1, 0, 0], 0.4, refraction_count = (2, 30))
Example above creates a analyses plan with equation
\[x=0.4\]and with a condition that refraction_count is higher than 2 but smaller than 30.
-
create_chessboard(ts, center, mesh, normal=[0, 0, 1], na=0.0, angles=11)¶
-
create_cylinder_element(center, radius, length, n, normal)¶ Creates a cilindrical element.
- Parameters
center (array_like) – A 3 dimensional of center position.
radius (float) – Cilinder radius.
length (float) – Cilinder length.
n (float) – Index of refraction.
normal (array_like) – A 3 dimensional array of cilinder axis.
-
create_fiberbundle(tr, center, mesh, normal=[0, 0, 1], na=0.0, angles=11)¶
-
create_parabolic_surface_element(c, n, th, wid, pp)¶ Creates a parabolic element surface.
- Parameters
c (array_like) – A 3 dimensional position.
n (float) – The index of refraction. -1 creates a reflective material.
th (float) – Thickness of the element. Y direction total size is given by this value.
wid (float) – Width of the element. X direction total size is given by this value.
pp (float) – P parameter of the parabola.
Notes
Parabolic curve is defined as:
\[2pp * (Z- Z0) = (Y-Y0)^{2} - (X-X0)^{2}\]
-
create_rectangle_element(val, n, normal)¶ Creates a rectangular element.
- Parameters
val (array_like) – A 6 dimensional array containing [xmin, xmax, ymin, ymax, zmin, zmax]
n (float) – The index of refraction
normal (array_like) – A 3 dimensional array containing the surface normal.
Notes
You can use a rectangular element to destroy some of your active elements by setting index of refraction as 1.0 and normal vector as [0, 0, 0].
-
create_sphere_element(c, r, n)¶ Creates a spherical element.
- Parameters
c (array_like) – A 3 dimensional sphere center.
r (float) – Sphere radius.
n (float) – Element index of refraction.
-
create_thin_lens(cplane, focus, aperture, index_refr, lens_type='plane-convex')¶ Creates a thin lens based on other geometrical elements. Allow one to create either a plane-convex or a convex-plane lens.
- Parameters
cplane (array_like) – A 3 dimensional array representing the center of the lens. This point is always in the plane side of the lens.
focus (float) – Lens focus.
aperture (float) – Lens aperture.
index_refr (float) – Lens index of refraction.
lens_type (str) – Must be either plane-convex or convex-plane
- Raises
AssertionError – If aperture is bigger than focus or if type is not ‘plane-convex’ or ‘convex-plane’, raises an exception error.
-
d2_flex_source(r, c=[0, 0, 0], normal=[0, 0, 1], na=0.0, angles=11, plan='xy')¶ Creates a 2d source along axis z with a numerical aperture
- Parameters
r (float) – Source radius. If 0.0, creates a point source.
c (array_like:) – Source center.
normal (array_like) – Source normal. Photon propagation direction
na (float) – Source numerical aperture along axis [1, 1, 0].
angles (int) – Numerical aperture discretization
Examples
>>> d2_source(0.1, [0, 0, 0], [0, 0, 1], na=0.0, angles=1)
Creates a collimated source with radius 0.1 centered at [0, 0, 0].
>>> d2_source(0.0, [0, 0, 0], [0, 0, 1], na=0.22, angles=11)
Creates a point source centered at [0, 0, 0] with numerical aperture 0.22. Angles are discretized by 0.04 (-0.22 to 0.22) in 11 points.
>>> d2_source(0.3, [0, 0, 0], [0, 0, 1], na=0.22, angles=11)
Creates a diverging source centered at [0, 0, 0] with numerical aperture 0.22 and radius 0.3.
- Raises
AssertionError – If numerical aperture is higher than 0, propagation vector must be aligned with an axis. This function does not create an angle cone for an arbitrary propagation vector.
-
d2_source(r, c=[0, 0, 0], normal=[0, 0, 1], na=0.0, angles=11)¶ Creates a 2d source along axis z with a numerical aperture
- Parameters
r (float) – Source radius. If 0.0, creates a point source.
c (array_like:) – Source center.
normal (array_like) – Source normal. Photon propagation direction
na (float) – Source numerical aperture along axis [1, 1, 0].
angles (int) – Numerical aperture discretization
Examples
>>> d2_source(0.1, [0, 0, 0], [0, 0, 1], na=0.0, angles=1)
Creates a collimated source with radius 0.1 centered at [0, 0, 0].
>>> d2_source(0.0, [0, 0, 0], [0, 0, 1], na=0.22, angles=11)
Creates a point source centered at [0, 0, 0] with numerical aperture 0.22. Angles are discretized by 0.04 (-0.22 to 0.22) in 11 points.
>>> d2_source(0.3, [0, 0, 0], [0, 0, 1], na=0.22, angles=11)
Creates a diverging source centered at [0, 0, 0] with numerical aperture 0.22 and radius 0.3.
- Raises
AssertionError – If numerical aperture is higher than 0, propagation vector must be aligned with an axis. This function does not create an angle cone for an arbitrary propagation vector.
-
d2_source_rectangle(size, c=[0, 0, 0], normal=[0, 0, 1], na=0.0, angles=11)¶ Creates a 2d source along axis z with a numerical aperture
- Parameters
r (float) – Source radius. If 0.0, creates a point source.
c (array_like:) – Source center.
normal (array_like) – Source normal. Photon propagation direction
na (float) – Source numerical aperture along axis [1, 1, 0].
angles (int) – Numerical aperture discretization
- Raises
AssertionError – If numerical aperture is higher than 0, propagation vector must be aligned with an axis. This function does not create an angle cone for an arbitrary propagation vector.
-
distance(vec1, vec2)¶ Distance between two points.
- Parameters
vec1 (array_like) – A 3 dimensional array representing first position.
vec2 (array_like) – A 3 dimensional array representing second position.
-
grid_to_pos(grid)¶ Get position from grid index.
- Parameters
grid (array_like) – Grid Position.
- Returns
Position.
- Return type
array_like
- Raises
Exeception – If grid is outside cell boundaries.
-
index_from_pos(pos)¶ Get index of refraction from position.
- Parameters
pos (array_like) – Desired Position.
- Returns
The correspondent index of refraction.
- Return type
float
- Raises
Exception – If value is outside cell boundaries
-
is_photon_in_cell(photon)¶ Check if photon is in the cell.
- Parameters
photon (class.photon)
- Returns
- Return type
bool
-
merge_photon_lists(my_photon_lists)¶ Merge photon lists in a single array. This function is used for multi processing simulations.
- Parameters
my_photon_lists (array_like) – An array containing objects class.photon_list
- Returns
Return an array like containing photon_lists.
- Return type
array_like
-
normal_and_index_from_pos(pos)¶ Get both normal surface and index of refraction from pos
- Parameters
pos (array_like) – Desired position.
- Returns
A tuple in which first element is the normal and the second is the refractive index.
- Return type
tuple
-
normal_from_pos(pos)¶ Get normal direction from position.
- Parameters
pos (array_like) – Desired Position.
- Returns
The correspondent normal direction.
- Return type
array_like
-
pos_to_grid(pos)¶ Get grid index from position in the 3D space.
- Parameters
pos (float) – Desired position
- Returns
Grid position
- Return type
array_like
- Raises
Exception – If value is outside cell boundaries
-
pos_to_grid_1D(value, index)¶ Gets grid index from position for a single dimensional axis
- Parameters
value (float) – Desired position
index – Desired axis (0 for ‘x’, 1 for ‘y’ and 2 for ‘z’)
- Returns
Grid position.
- Return type
int
- Raises
Exception – If value is outside cell boundaries.
-
prepare_acquisition(split=1)¶
-
reset_photons()¶ Reset initial photon lists in order to run simulation a second time.
-
reset_structures()¶ Reset initial structures in order to run simulation a second time.
-
rotate(ang, axis, origin, ROI=None)¶ Rotate a selected ROI in a given direction and origin. This function only rotates grid points that have index of refraction different of 1.
- Parameters
ang (float) – Angle of rotation (in radians).
axis (array_like) – A 3 dimensional axis of rotation.
origin (array_like) – A 3 dimensional origin position.
ROI (array_like) – A 6 dimensional ROI in the form of [xmin, xmax, ymin, ymax, zmin, zmax].
- Raises
Exception – If there is no element to rotate.
-
run(run_index=0, multiprocessing=False, return_dict={}, xsym=False, ysym=False)¶ Simulation run main function.
- Parameters
run_index (int) – Run index. Used for multi process simulation
split (int) – Split initial photons in equal parts of ‘split’. Used for multi process simulation
xsym (bool) – If simulation is symmetric with respect to x axis, you can reduce number of running photons by two.
ysym (bool) – If simulation is symmetric with respect to y axis, you can reduce number of running photons by two.
- Returns
Returns an array of class.photon_list. Lenght is given by the name of planes created. Each element in the list has class.photon objects. Length for each element is the number of photons saved for the correspondent plane.
- Return type
array_like
- Raises
AssertionError – Raises an assertionError if run_index is smaller than split. You cannot run subset 1 if you have not splitted your initial photon array. Maximum run_index is always split-1. Raises an AssertionError also if there is no photons in the initial set.
-
run_photon(initial_photons, rindex)¶ Simulatiom run for the initial set (or subset) of photons. It terminates when all photons leave the cell.
- Parameters
initial_photon_list (array_like) – The complete set (or subset) of the intial photons.
rindex (int) – Run index. Used for multi process simulations.
See also
run()Starts the simulation for each set or subset of photon. Can be used with multi processing.
-
show_created_elements(mode)¶ Function shows in 3D elements created. Those are the initial conditions of the problem
- Parameters
mode (str) – Defines what is displayed. ‘all’ displays all, ‘all-noplan’ or ‘-noplan’ excludes analyzes plans and anything else display only photons.
-
show_elements(my_photon_lists, mode='all-noplan')¶ Shows in 3D elements at the end of simulation.
- Parameters
my_photon_lists (array_like) – Array_like containing class.photon_list objects.
mode (str) – Defines what is displayed. ‘all’ displays all, ‘all-noplan’ or ‘-noplan’ excludes analyzes plans and anything else display only photons.
-
show_photons2D(my_photon_lists, mode='all', binning=21)¶ See also
Deprecated()
-
show_single_photons2D(my_photon_lists, mode='all', binning=21)¶ See also
Deprecated()
-
class
orsaytrace.trace.photon(pos, normal, intensity=1.0)¶ Bases:
objectA photon is the basic structure of this simulation. It contains the attributes pos, normal, intensity, n, last_surface, refraction_count and reflection_count
- Parameters
pos (list) – A 3 dimensional list stating each photon position
normal (list) – A 3 dimensional list stating each photon normal (or wavevector)
intensity (float) – A float for the initial photon intensity
Notes
Other attributes of photon available to user and used during simulation run:
- n: float
Index of refraction
- last_surface: array_like
As photons crosses the grid, it overwrites the last valid surface normal
- refraction_count: int
Number os refractions photon has done
- reflection_count: int
Number os reflections photon has done
-
get_attr()¶ This gets all photon attributes. Can be used to clone photons.
- Returns
An array_like object as following: [pos, normal, intensity, n, last_surface, refraction_count, reflection_count]
- Return type
array_like
-
move(value: float)¶ Increment photon position.
- Parameters
value (float) – Increment photon position by value and given by normal (photon direction)
- Returns
True if is sucessfull.
- Return type
boolean
-
reflection()¶ Reflection photon object using last_surface attribute.
- Returns
True if reflection is sucessfull and false if you attempted a reflection in a grid which there was no surface normal vector.
- Return type
boolean
-
refraction(n2)¶ Refraction photon object using last_surface attribute and a given index of refraction.
- Parameters
n2 (float) – The index of refraction.
- Returns
True if refraction is sucessfull and false if you attempted a refraction in a grid which there was no surface normal vector
- Return type
boolean
-
set_attr(values)¶ This sets all photon attributes. Can be used to clone photons.
- Parameters
values (array_like) – An array_like object as following: [pos, normal, intensity, n, last_surface, refraction_count, reflection_count]
-
update(value)¶ Update is a core function that updates photon last_surface and photon normal by means of refractions of reflections.
- Parameters
value (array_like) – A 2 dimensional array where first index is the surface normal and the second is the index of refraction.
- Returns
True if photon suffered a refraction or reflection. False if nothing happens with photon normal.
- Return type
boolean
Notes
If return is True, photon moves a single step without refraction/reflection using move. This avoids non intentional double reflections/refractions
-
class
orsaytrace.trace.photon_list(normal: list, value: float)¶ Bases:
objectphoton_list is a class that bounds photons together by a plane and a condition. Plane is characterized by a normal vector and its value.
- Parameters
normal (list) – A 3 dimensional list stating the normal vector from plane
value (float) – The plane independent variable in the form of x + y + z = value
-
add_photon(sphoton)¶ Instanteates a new photon before calling append_photon. Append photon checks if photon object satifies conditions.
- Parameters
sphoton (class.photon) – A photon coming from photon class.
-
add_symmetric_xphoton(sphoton)¶ Instanteates a new X symmetric photon before callind append_photon. Append photon checks if photon object satifies conditions.
- Parameters
sphoton (class.photon) – A photon coming from photon class
-
add_symmetric_yphoton(sphoton)¶ Instanteates a new Y symmetric photon before callind append_photon. Append photon checks if photon object satifies conditions.
- Parameters
sphoton (class.photon) – A photon coming from photon class
-
append_photon(photon)¶ Append a photon tp the photon_list given that the condition is True.
- Parameters
photon (class.photon) – A photon coming from photon class.
-
avg_distance_axis_y(c=[0, 0])¶ Calculates the distance for all photons in photon_list relative to a (x, z) coordinate in a Y nornal plane.
- Parameters
c (The position in x, z plane)
- Returns
- Return type
float
-
avg_distance_axis_z(c=[0, 0])¶ Calculates the distance for all photons in photon_list relative to a (x, y) coordinate in a Z normal plane.
- Parameters
c (The position in x, y plane)
- Returns
- Return type
float
-
avg_divergence(vec_ref)¶ Calculates the average photon divergence from given a vector reference.
- Parameters
vec (array_like) – A 3 dimensional array of given direction.
- Returns
A float provenient from a numpy.average based on vector direction and propagation direction for each photon
- Return type
float
-
avg_position()¶ Calculates the average position of all photons in the photon_list.
- Returns
A 3 dimensional numpy array
- Return type
numpy.asarray
-
distance(pos)¶
-
distance_point_to_plane(pos: list)¶ Calculates the distance from a given point to plane.
- Parameters
pos (A 3 dimensional list of a given position)
- Returns
A float provenient from a numpy.abs based on point position and plane equation
- Return type
float
-
get_average_weighted_inverse()¶ Calculates the average of get_weighted_inverse.
- Returns
- Return type
float
-
get_average_weighted_inverse_axis_y(c=[0, 0])¶ Calculates the average of get_weighted_inverse relative to a (x, z) coordinate in Y normal plane.
- Parameters
c (The position in x, z plane)
- Returns
- Return type
float
-
get_intensities()¶ Gets the intensity of all photons in the photon_list
- Returns
A photons-dimensional numpy array
- Return type
numpy.asarray
-
get_positions()¶ Get positions of all photons in the photon_list.
- Returns
A numpy array
- Return type
numpy.asarray
-
get_relative_centroid_distances()¶ Calculates the distances from a given point to plane.
- Returns
- Return type
numpy.asarray
-
get_relative_centroid_positions()¶ Get the relative position of all photons in the photon_list relative from centroid (calculted using avg_position)
- Returns
- Return type
numpy.asarray
-
get_weighted_inverse()¶ Calculates the division of intensity and photon relative distance from centroid fro all photons in photon_list
- Returns
A n dimensional numpy array
- Return type
numpy array
-
max_position()¶ Calculates the max position in each axis for all photons in photon_list.
- Returns
A 3 dimensional numpy array
- Return type
numpy.asarray
-
min_position()¶ Calculates the min position in each axis for all photons in photon_list.
- Parameters
pos (A 3 dimensional list of a given position)
- Returns
A 3 dimensional numpy array
- Return type
numpy.asarray
-
std_position()¶ Calculates the standard deviation of the position for all photons in photon_list.
- Returns
A 3 dimensional numpy array
- Return type
numpy.asarray