Full API Documentation

Functions of Trace Module

class orsaytrace.trace.Simu(x, y, z, res)

Bases: object

A simulation class that will allow one to create a space grid. This class also contains a few convenient functions in order to easily display data using matplotlib.

Parameters
  • x (float) – Cell size in x.

  • y (float) – Cell size in y.

  • z (float) – Cell size in z.

  • res (float) – Simulation resolution.

Other Parameters
  • ss (int) – A subsampling factor. Default is 2.01 and it says meshing will be done using sub pixels resolution for object creation.

  • 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.

  • photons (list) – A list containing photon objects. Those photons will propagate and be saved in class.photon_lists along the simulation.

  • photon_lists (numpy array) – A numpy array containing photon_lists objects. Array length is given by the number of created plans.

assign_block_n(min_index, max_index, value)

Assigns index of refraction for multiple points at once.

Parameters
  • min_index (array_like) – A positional 3 dimensional array (bottom left point).

  • max_index (array_like) – A positional 3 dimensional array (top right point).

  • normal (array_like) – A vectorial 3 dimensional array.

assign_block_normal(min_index, max_index, normal)

Assigns normal for multiple points at once.

Parameters
  • min_index (array_like) – A positional 3 dimensional array (bottom left point).

  • max_index (array_like) – A tpositional 3 dimensional array (top right point).

  • normal (array_like) – A vectorial three dimensional array.

assign_n(index, value)

Assigns index of refraction for a given point.

Parameters
  • index (array_like) – A positional 3 dimensional array.

  • value (float) – The index of refraction.

assign_normal(index, normal)

Assigns a surface normal in a given index.

Parameters
  • index (array_like) – A positional 3 dimensional array.

  • normal (array_like) – A vectorial 3 dimensional array.

check_analysis_plan(photon)

Check if photon is in any created analysis plan.

Parameters

photon (class.photon)

create_analysis_plan(normal, value, **kargs)

Create an analysis plan.

Parameters
  • normal (array_like) – A 3 dimensional 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. For int type attribute, it can take two forms: a int or a tuple. int represents minimal values, while tuple represents a possible interval. It is also possible to create conditional plans using photon ‘pos’ attribute. In this case, see example above for how to create it.

Examples

>>> create_analysis_plan([0, 1, 0], a, reflection_count = 1)

Example above creates an 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 an 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 an analyses plan with equation

\[x=0.4\]

and with a condition that refraction_count is higher than 2 but smaller than 30.

>>> a.create_analysis_plan([0, 0, 1], 1.0, reflection_count = (1, 1), pos=([1.2, 0.0, 0.0], (0.1, 0.4)))

Example above creates an analyses plan with equation

\[z = 1.0\]

and with a condition that reflection_count is exactly 1 and photon position has a distance from the point [1.2, 0.0, 0.0] between 0.1 and 0.4.

create_chessboard(total_size, center, mesh, normal=[0, 0, 1], na=0.0, angles=11)

Create a chessboard like source.

Parameters
  • total_size (float) – Total size of the chessboard.

  • center (array_like) – Center of the chessboard.

  • mesh (int) – Chessboard discretization. Unit square has total_size / mesh size.

  • normal (array_like) – Propagation vector.

  • na (float) – Numerical aperture.

  • angles (int) – Numerical aperture angle discretization

Raises

AssertionError – If mesh is not an integer.

create_cylinder_element(center, radius, length, n_refr, normal)

Creates a cylindrical element.

Parameters
  • center (array_like) – A 3 dimensional of center position.

  • radius (float) – Cylinder radius.

  • length (float) – Cylinder length.

  • n_refr (float) – Index of refraction.

  • normal (array_like) – A 3 dimensional vector of cylinder axis.

create_fiberbundle(total_radius, center, mesh, normal=[0, 0, 1], na=0.0, angles=11)

Create a optical fiber bundle like source.

Parameters
  • total_radius (float) – Total radius in which the bundle in necessarily inside.

  • center (array_like) – The center of the bundle.

  • mesh (int) – Bundle discretization factor. Each fiber has total_radius / mesh radius.

  • normal (array_like) – Propagation direction.

  • na (float) – Numerical aperture.

  • angles (int) – Numerical aperture angle discretization.

Raises

AssertionError – If mesh is not an integer or if mesh is an even number

Notes

Mesh must be an odd number. This comes from the selection rule of fibers inside bundle to do not overlap.

create_parabolic_surface_element(center, n_refr, thickness, width, pp)

Creates a parabolic element surface.

Parameters
  • center (array_like) – A 3 dimensional position.

  • n_refr (float) – The index of refraction. -1 creates a reflective material.

  • thickness (float) – Thickness of the element. Y direction total size is given by this value.

  • width (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(ROI, n_refr, normal)

Creates a rectangular element.

Parameters
  • ROI (array_like) – A 6 dimensional array in the form [xmin, xmax, ymin, ymax, zmin, zmax].

  • n_refr (float) – The index of refraction.

  • normal (array_like) – A 3 dimensional array containing the surface normal.

Raises

AssertionError – If ROI length in any dimension is negative.

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(center, radius, n_refr)

Creates a spherical element.

Parameters
  • center (array_like) – A 3 dimensional sphere center.

  • radius (float) – Sphere radius.

  • n_refr (float) – Index of refraction.

create_thin_lens(cplane, focus, aperture, index_refr, lens_type='plane-convex')

Creates a thin lens based on other basical geometrical elements. Allows 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 flat side of the lens.

  • focus (float) – Lens focus.

  • aperture (float) – Lens aperture.

  • index_refr (float) – Lens index of refraction.

  • lens_type (str) – The type of the lens.

Raises

AssertionError – If aperture is bigger than focus or if lens_type is not ‘plane-convex’ or ‘convex-plane’.

d2_flex_source(radius, center=[0, 0, 0], normal=[0, 0, 1], na=0.0, angles=11, plan='xy')

A circular flexible source method. It can create planar sources in X, Y or Z with a given normal vector and a given numerical aperture referenced by the normal. Angle discretization is controllable. A point source is also available.

Parameters
  • radius (float) – Source radius. If 0.0, creates a point source. Use ‘plan’ to create a point source.

  • center (array_like:) – Source center.

  • normal (array_like) – Source normal. Photon propagation direction.

  • na (float) – Source numerical aperture referenced by the normal vector axis.

  • angles (int) – Numerical aperture discretization.

Notes

When using the point source, radius value is not meaningful. You can create a point source by inputing 0 in radius value. This is however not 100% correct because you can have points inside your resolution grid volume cell.

Examples

>>> d2_flex_source(0.1, [0, 0, 0], [0, 0, 1], na=0.0, angles=1, plan='xy')

Creates a collimated source with radius 0.1 centered at [0, 0, 0] propagating in Z direction at X-Y plan.

>>> d2_flex_source(0.1, [0, 0, 0], [0, 1, 0], na=0.0, angles=1, plan='xy')

Creates a collimated source with radius 0.1 centered at [0, 0, 0] propagating in Y direction at X-Y plan.

>>> d2_flex_source(0.25, [0, 0, 2.0], [0, -1, 0], na=0.0, angles=1, plan='xz')

Creates a collimated source with radius 0.25 centered at [0, 0, 2.0] propagating in -Y direction at X-Z plan.

>>> d2_flex_source(10.0, [0, 0, 0], [0, 0, 1], na=0.22, angles=11, plan='point')

Creates a point source centered at [0, 0, 0] with numerical aperture 0.22 propagating in Z axis. Angles are discretized by 0.044 (-0.22 to 0.22) in 11 points. Note that 10.0 in radius changes nothing.

Raises

AssertionError – Negative numerical aperture and - if numerical aperture is higher than 0 - normal not in an axis (can be [0, 0, 1]; [0, 1, 0] or [1, 0, 0]).

d2_source(radius, center=[0, 0, 0], normal=[0, 0, 1], na=0.0, angles=11)

Creates a 2d source along a given normal and an numerical aperture. Source can only be created in X-Y plan.

See also

d2_flex_source()

A more complete option. d2_source is a subset of d2_flex_source in terms of functionality. d2_flex_source can create a numerical aperture cone for any given normal vector and can create sources in plans other than X-Y.

Parameters
  • radius (float) – Source radius. If 0.0, creates a point source.

  • center (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] in X-Y plane.

>>> 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.044 (-0.22 to 0.22) in 11 points. Source is in X-Y plane.

>>> 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. Source is in X-Y plane.

Raises

AssertionError – If numerical aperture is higher than 0, propagation vector must be aligned with an axis.

d2_source_rectangle(size, center=[0, 0, 0], normal=[0, 0, 1], na=0.0, angles=11)

Creates a rectangular 2D source in X-Y plane for any given normal vector and numerical aperture.

Parameters
  • size (array_like) – Source size in X-Y.

  • center (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 – len(size) is not 2 or any of the values are negative or equal to zero. Also if numerical aperture is a negative value.

distance(vec1, vec2)

Calculates the distance between two points.

Parameters
  • vec1 (array_like) – A 3 dimensional array for the first position.

  • vec2 (array_like) – A 3 dimensional array for the second position.

Returns

The distance.

Return type

float

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) – Position.

Returns

The correspondent index of refraction.

Return type

float

Raises

AssertionError – If value is outside cell boundaries.

is_photon_in_cell(photon)

Check if photon is in the simulation 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 at the end of multi processing simulations to bind together an array of class.photon_list.

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 position.

Parameters

pos (array_like) – Position.

Returns

First element is the normal and the second is the refractive index.

Return type

tuple

Raises

AssertionError – If value is outside cell boundaries.

normal_from_pos(pos)

Get normal direction from position.

Parameters

pos (array_like) – 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) – Position.

Returns

Grid position.

Return type

array_like

Raises

AssertionError – 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) – The position.

  • index – Desired axis (0 for ‘x’, 1 for ‘y’ and 2 for ‘z’).

Returns

Grid position.

Return type

int

Raises

AssertionError – If value is outside cell boundaries.

prepare_acquisition(split=1)

Internal function to prepare for multiprocess simulation. If no multiprocess is used, this function is called nonetheless with split = 1.

Parameters

split (int) – The number of processes. Slits initial photon list in ‘split’ smaller lists.

Raises

AssertionError – If there is no initial photon list or if split is less or equal to zero.

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. This value is 0 for single process simulation.

  • multiprocessing (bool) – Explicit if simulation must run using multi processes.

  • return_dict (dict) – You must pass this argument from multiprocessing.manager in order to use the multiprocessing capability.

  • 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

Array containing objects from class.photon_list. len(array) is given by the name of planes created. Length of each element is the number of photons saved for the correspondent plane.

Return type

array_like

Raises

AssertionError – If run_index is smaller than split (maximum run_index is always split-1) or if there is no photons in the initial set or if it is a multiprocessing simulation but return dict is not an instance of multiprocessing.managers.DictProxy.

See also

run_photon()

Starts the simulation for each set or subset of photon. This is an internal function and in normal conditions must not be called by user.

prepare_acquisition()

Must be called by user for multi process simulation. If a single process simulation is used, this function is called internaly and user must not bother.

merge_photon_list()

Must be called by user if using multi process simulation. For single process simulation, this function must not be used.

run_photon(initial_photons, rindex)

Simulation run for the initial set (or subset) of photons. It terminates when all photons leave the cell.

Parameters
  • initial_photons (array_like) – The complete set (or subset) of the initial photons.

  • rindex (int) – Run index. Used for multi process simulations.

show_created_elements(mode)

Shows 3D projection of the initial condition of the problem.

Parameters

mode (str) –

‘all’ displays everything. If ‘-noplan’ is appended, it excludes analyzes plans. Appending ‘-verbose’

prints the number of initial photons.

show_elements(my_photon_lists, mode='all-noplan')

Shows 3D projection at the end of simulation.

Parameters
  • my_photon_lists (array_like) – Must contain class.photon_list objects.

  • mode (str) –

    ‘all’ displays everything. If ‘-noplan’ is appended, it excludes analyzes plans. Appending ‘-verbose’

    prints the number of initial photons.

Raises

AssertionError: – If my_photon_list is not a list of class.photon_list objects.

show_photons2D(my_photon_lists, mode='all', binning=21)

Show a 2 dimensional histogram for an array of photon lists. Can be used for plans perpendicular to x, y and z axis.

Parameters
  • my_photon_lists (array_like) – An array of length 1. Element is an class.photon_list

  • mode (str) – If ‘-verbose’, prints the number of photons plotted.

  • binning (int) – Number of bins in the histogram. X and Y bins are equal.

Raises
  • AssertionError: – If my_photon_lists must have len() less or equal to 1.

  • Exception: – If elements in my_photon_list are not a plan perpendicular to x, y or z.

show_single_photons2D(my_photon_lists, mode='all', binning=21)

Show a 2 dimensional histogram for a single photon list. Can be used for plans perpendicular to x, y and z axis.

Parameters
  • my_photon_lists (array_like) – An array of length 1. Element is an class.photon_list

  • mode (str) – If ‘-verbose’, prints the number of photons plotted.

  • binning (int) – Number of bins in the histogram. X and Y bins are equal.

Raises
  • AssertionError: – If my_photon_lists have len() different from 1.

  • Exception: – If my_photon_list[0] is not a plan perpendicular to x, y or z.

class orsaytrace.trace.photon(pos, normal, intensity=1.0)

Bases: object

A photon is the basic structure of this simulation. It contains the attributes pos, normal, intensity, n, last_surface, refraction_count, reflection_count and a init dictionary.

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.

Other Parameters
  • 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.

  • init (dict) – A dictionary that does not change during simulation. Can be used to track photons.

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 successful.

Return type

boolean

reflection()

Reflection photon object using last_surface.

Returns

True if reflection is successful and false if you attempted a reflection in a grid which there was no surface normal vector.

Return type

boolean

refraction(n_refr)

Refraction photon object using last_surface and a given index of refraction.

Parameters

n_refr (float) – The index of refraction.

Returns

True if refraction is successful 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 not.

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: object

photon_list is a class that bounds photons together by a plane and a condition. Plane is characterized

by a normal vector and its independent 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)
Instantiates a new photon before calling append_photon. Append photon checks if photon

object satisfies all conditions.

Parameters

sphoton (class.photon) – A photon coming from photon class.

add_symmetric_xphoton(sphoton)

Instantiates a new X symmetric photon before calling append_photon. Append photon checks if photon object satisfies conditions.

Parameters

sphoton (class.photon) – A photon coming from photon class

add_symmetric_yphoton(sphoton)

Instantiates a new Y symmetric photon before calling append_photon. Append photon checks if photon object satisfies conditions.

Parameters

sphoton (class.photon) – A photon coming from photon class

append_photon(photon)

Append a photon to the photon_list given that the condition is satisfied.

Parameters

photon (class.photon) – A photon coming from photon class.

avg_distance_axis_y(center=[0, 0])

Calculates the distance for all photons in photon_list relative to a (x, z) coordinate in a Y normal plane.

Parameters

center (The position in x, z plane)

Returns

Return type

float

avg_distance_axis_z(center=[0, 0])

Calculates the distance for all photons in photon_list relative to a (x, y) coordinate in a Z normal plane.

Parameters

center (array_like) – 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 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 positional 3 dimensional numpy array.

Return type

numpy.asarray

distance_point_to_plane(pos)

Calculates the distance from a given point to plane.

Parameters

pos (array_like) – 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(center=[0, 0])

Calculates the average of get_weighted_inverse relative to a (x, z) coordinate in Y normal plane.

Parameters

center (The position in x, z plane)

Returns

Return type

float

Raises

AssertionError – If normal is not in Y direction.

get_intensities()

Gets the intensity of all photons in the photon_list.

Returns

A photon-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

float

get_relative_centroid_positions()

Get the relative position of all photons in the photon_list relative from plane center of mass (calculted using avg_position).

Returns

Return type

numpy.asarray

get_weighted_inverse()

Calculates the division of intensity and photon relative distance from center of mass from all photons in photon_list.

Returns

A photon-dimensional numpy array

Return type

numpy array

max_position()

Calculates the max position in each axis for all photons in photon_list.

Returns

A positional 3 dimensional numpy array.

Return type

numpy.asarray

min_position()

Calculates the min position in each axis for all photons in photon_list.

Returns

A positional 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