3.1. Main Optical Classes¶
3.1.1. Optical Rays¶
- class raysect.optical.ray.Ray¶
Bases:
raysect.core.ray.Ray
Optical Ray class for optical applications, inherits from core Ray class.
Provides the trace(world) method.
- Parameters
origin (Point3D) – Point defining ray’s origin (default=Point3D(0, 0, 0))
direction (Vector3D) – Vector defining ray’s direction (default=Vector3D(0, 0, 1))
min_wavelength (float) – Lower wavelength bound for observed spectrum
max_wavelength (float) – Upper wavelength bound for observed spectrum
bins (int) – Number of samples to use over the spectral range
max_distance (float) – The terminating distance of the ray
extinction_prob (float) – Probability of path extinction at every material surface interaction (default=0.1)
extinction_min_depth (int) – Minimum number of paths before triggering extinction probability (default=3)
max_depth (int) – Maximum number of material interactions before terminating ray trajectory.
importance_sampling (bool) – Toggles use of importance sampling for important primitives. See help documentation on importance sampling, (default=True).
important_path_weight (float) – Weight to use for important paths when using importance sampling.
>>> from raysect.core import Point3D, Vector3D >>> from raysect.optical import World, Ray >>> >>> world = World() >>> >>> ray = Ray(origin=Point3D(0, 0, -5), >>> direction=Vector3D(0, 0, 1), >>> min_wavelength=375, >>> max_wavelength=785, >>> bins=100) >>> >>> spectrum = ray.trace(world) >>> spectrum <raysect.optical.spectrum.Spectrum at 0x7f5b08b6e048>
- bins¶
Number of spectral bins across wavelength range.
- Return type
int
- copy()¶
Obtain a new Ray object with the same configuration settings.
- Parameters
- Return type
>>> from raysect.core import Point3D, Vector3D >>> from raysect.optical import Ray >>> >>> ray = Ray(origin=Point3D(0, 0, -5), >>> direction=Vector3D(0, 0, 1), >>> min_wavelength=375, >>> max_wavelength=785, >>> bins=100) >>> >>> ray.copy() Ray(Point3D(0.0, 0.0, -5.0), Vector3D(0.0, 0.0, 1.0), inf)
- extinction_min_depth¶
Minimum number of paths before triggering extinction probability.
- Return type
int
- extinction_prob¶
Probability of path extinction at every material surface interaction.
- Return type
float
- important_path_weight¶
Weight to use for important paths when using importance sampling.
- Return type
float
- max_depth¶
Maximum number of material interactions before terminating ray trajectory.
- Return type
int
- max_wavelength¶
Upper bound on wavelength range.
- Return type
float
- min_wavelength¶
Lower bound on wavelength range.
- Return type
float
- new_spectrum()¶
Returns a new Spectrum compatible with the ray spectral settings.
- Return type
>>> from raysect.core import Point3D, Vector3D >>> from raysect.optical import Ray >>> >>> ray = Ray(origin=Point3D(0, 0, -5), >>> direction=Vector3D(0, 0, 1), >>> min_wavelength=375, >>> max_wavelength=785, >>> bins=100) >>> >>> ray.new_spectrum() <raysect.optical.spectrum.Spectrum at 0x7f5b08b6e1b0>
- sample()¶
Samples the radiance directed along the ray direction.
This methods calls trace repeatedly to obtain a statistical sample of the radiance directed along the ray direction from the world. The count parameter specifies the number of samples to obtain. The mean spectrum accumulated from these samples is returned.
- Parameters
world (World) – World object defining the scene.
count (int) – Number of samples to take.
- Returns
The accumulated spectrum collected by the ray.
- Return type
>>> from raysect.core import Point3D, Vector3D >>> from raysect.optical import World, Ray >>> >>> world = World() >>> >>> ray = Ray(origin=Point3D(0, 0, -5), >>> direction=Vector3D(0, 0, 1), >>> min_wavelength=375, >>> max_wavelength=785, >>> bins=100) >>> >>> ray.sample(world, 10) <raysect.optical.spectrum.Spectrum at 0x7f5b08b6e318>
- spawn_daughter()¶
Spawns a new daughter of the ray.
A daughter ray has the same spectral configuration as the source ray, however the ray depth is increased by 1.
- trace()¶
Traces a single ray path through the world.
- Parameters
world (World) – World object defining the scene.
keep_alive (bool) – If true, disables Russian roulette termination of the ray.
- Returns
The resulting Spectrum object collected by the ray.
- Return type
>>> from raysect.core import Point3D, Vector3D >>> from raysect.optical import World, Ray >>> >>> world = World() >>> >>> ray = Ray(origin=Point3D(0, 0, -5), >>> direction=Vector3D(0, 0, 1), >>> min_wavelength=375, >>> max_wavelength=785, >>> bins=100) >>> >>> spectrum = ray.trace(world) >>> spectrum <raysect.optical.spectrum.Spectrum at 0x7f5b08b6e048>
- wavelength_range¶
Upper and lower wavelength range.
- Return type
tuple
3.1.2. Spectral Functions¶
- class raysect.optical.spectralfunction.SpectralFunction¶
SpectralFunction abstract base class.
A common interface for representing optical properties that are a function of wavelength. It provides methods for sampling, integrating and averaging a spectral function over specified wavelength ranges. The optical package uses SpectralFunctions to represent a number of different wavelength dependent optical properties, for example emission spectra, refractive indices and attenuation curves.
Deriving classes must implement the integrate method.
It is also recommended that subclasses implement __call__(). This should accept a single argument - wavelength - and return a single sample of the function at that wavelength. The units of wavelength are nanometers.
A number of utility sub-classes exist to simplify SpectralFunction development.
see also: NumericallyIntegratedSF, InterpolatedSF, ConstantSF, Spectrum
- average()¶
Average radiance over the requested spectral range (W/m^2/sr/nm).
Virtual method, to be implemented in child classes.
- Parameters
min_wavelength (float) – lower wavelength for calculation
max_wavelength (float) – upper wavelength for calculation
- Return type
float
>>> spectrum = ray.trace(world) >>> spectrum.average(400, 700) 1.095030870970234
- evaluate()¶
Evaluate the spectral function f(wavelength)
- Parameters
wavelength (float) – Wavelength in nanometers.
- Return type
float
- integrate()¶
Calculates the integrated radiance over the specified spectral range.
Virtual method, to be implemented in child classes.
- Parameters
min_wavelength (float) – The minimum wavelength in nanometers
max_wavelength (float) – The maximum wavelength in nanometers
- Returns
Integrated radiance in W/m^2/str
- Return type
float
>>> spectrum = ray.trace(world) >>> spectrum.integrate(400, 700) 328.50926129107023
- sample()¶
Re-sample the spectral function with a new wavelength range and resolution.
- Parameters
min_wavelength (float) – lower wavelength for calculation
max_wavelength (float) – upper wavelength for calculation
bins (int) – The number of spectral bins
- Return type
ndarray
>>> spectrum <raysect.optical.spectrum.Spectrum at 0x7f56c22bd8b8> >>> spectrum.min_wavelength, spectrum.max_wavelength (375.0, 785.0) >>> sub_spectrum = spectrum.sample(450, 550, 100)
- class raysect.optical.spectrum.Spectrum¶
Bases:
raysect.optical.spectralfunction.SpectralFunction
A class for working with spectra.
Describes the distribution of light at each wavelength in units of radiance (W/m^2/str/nm). Spectral samples are regularly spaced over the wavelength range and lie in the centre of the wavelength bins.
- Parameters
min_wavelength (float) – Lower wavelength bound for this spectrum
max_wavelength (float) – Upper wavelength bound for this spectrum
bins (int) – Number of samples to use over the spectral range
>>> from raysect.optical import Spectrum >>> >>> spectrum = Spectrum(400, 720, 250)
- clear()¶
Resets the sample values in the spectrum to zero.
- is_compatible()¶
Returns True if the stored samples are consistent with the specified wavelength range and sample size.
- Parameters
min_wavelength (float) – The minimum wavelength in nanometers
max_wavelength (float) – The maximum wavelength in nanometers
bins (int) – The number of bins.
- Returns
True if the samples are compatible with the range/samples, False otherwise.
- Return type
boolean
- is_zero()¶
Can be used to determine if all the samples are zero.
True if the spectrum is zero, False otherwise.
- Return type
bool
>>> spectrum = ray.trace(world) >>> spectrum.is_zero() False
- new_spectrum()¶
Returns a new Spectrum compatible with the same spectral settings.
- Return type
- to_photons()¶
Converts the spectrum sample array from radiance W/m^2/str/nm to Photons/s/m^2/str/nm and returns the data in a numpy array.
- Return type
ndarray
>>> spectrum = ray.trace(world) >>> spectrum.to_photons() array([2.30744985e+17, 3.12842916e+17, ...])
- total()¶
Calculates the total radiance integrated over the whole spectral range.
Returns radiance in W/m^2/str
- Return type
float
>>> spectrum = ray.trace(world) >>> spectrum.total() 416.6978223103715
- wavelengths¶
Wavelength array in nm
- Return type
ndarray
- class raysect.optical.spectralfunction.NumericallyIntegratedSF¶
Bases:
raysect.optical.spectralfunction.SpectralFunction
Numerically integrates a supplied function.
This abstract class provides an implementation of the integrate method that numerically integrates a supplied function (typically a non-integrable analytical function). The function to numerically integrate is supplied by sub-classing this class and implementing the function() method.
The function is numerically sampled at regular intervals. A sampling resolution may be specified in the class constructor (default: 1 sample/nm).
- Parameters
sample_resolution (double) – The numerical sampling resolution in nanometers.
- class raysect.optical.spectralfunction.InterpolatedSF¶
Bases:
raysect.optical.spectralfunction.SpectralFunction
Linearly interpolated spectral function.
Spectral function defined by samples of regular or irregular spacing, ends are extrapolated. You must set the ends to zero if you want the function to go to zero at the edges!
wavelengths and samples will be sorted during initialisation.
If normalise is set to True the data is rescaled so the integrated area of the spectral function over the full range of the input data is normalised to 1.0.
- Parameters
wavelengths (object) – 1D array of wavelengths in nanometers.
samples (object) – 1D array of spectral samples.
normalise (bool) – True/false toggle for whether to normalise the spectral function so its integral equals 1.
>>> from raysect.optical import InterpolatedSF >>> >>> # defining a set of spectral filters >>> filter_red = InterpolatedSF([100, 650, 660, 670, 680, 800], [0, 0, 1, 1, 0, 0]) >>> filter_green = InterpolatedSF([100, 530, 540, 550, 560, 800], [0, 0, 1, 1, 0, 0]) >>> filter_blue = InterpolatedSF([100, 480, 490, 500, 510, 800], [0, 0, 1, 1, 0, 0])
- class raysect.optical.spectralfunction.ConstantSF¶
Bases:
raysect.optical.spectralfunction.SpectralFunction
Constant value spectral function
- Parameters
value (float) – Constant radiance value
>>> from raysect.optical import ConstantSF >>> >>> unity_spectral_function = ConstantSF(1.0)
- raysect.optical.spectrum.photon_energy()¶
Returns the energy of a photon with the specified wavelength.
- Parameters
wavelength (float) – Photon wavelength in nanometers.
- Returns
Photon energy in Joules.
- Return type
float
3.1.3. Colours¶
The CIE 1931 colour spaces define quantitatively the link between pure physical colours (i.e. wavelengths in the visible spectrum) and human perceivable colours. The mathematical relationships between the spectrum and perceivable colours are based on the sensitivity curves for the three different cone cells in the human eye. Raysect implements three X, Y, Z normalised spectral functions from the CIE 1931 Standard Colorimetric Observer. For more information see Wikipedia.
- raysect.optical.colour.ciexyz_x¶
X spectral function from CIE 1931 Standard Colorimetric Observer (normalised)
- raysect.optical.colour.ciexyz_y¶
Y spectral function from CIE 1931 Standard Colorimetric Observer (normalised)
- raysect.optical.colour.ciexyz_z¶
Z spectral function from CIE 1931 Standard Colorimetric Observer (normalised)
- raysect.optical.colour.d65_white¶
CIE D65 standard illuminant, normalised to 1.0 over visual range 375-785 nm
- raysect.optical.colour.resample_ciexyz()¶
Pre-calculates samples of XYZ sensitivity curves over desired spectral range.
Returns ndarray of shape [N, 3] where the last dimension (0, 1, 2) corresponds to (X, Y, Z).
- Parameters
min_wavelength (float) – Lower wavelength bound on spectrum
max_wavelength (float) – Upper wavelength bound on spectrum
bins (int) – Number of spectral bins in spectrum
- Return type
memoryview
- raysect.optical.colour.spectrum_to_ciexyz()¶
Calculates a tuple of CIE X, Y, Z values from an input spectrum
- Parameters
spectrum (Spectrum) – Spectrum to process
resampled_xyz (memoryview) – Pre-calculated XYZ sensitivity curves optimised for this spectral range (default=None).
- Return type
tuple
- raysect.optical.colour.ciexyy_to_ciexyz()¶
Performs conversion from CIE xyY to CIE XYZ colour space
Returns a tuple of (X, Y, Z)
- Parameters
cx (float) – chromaticity x
cy (float) – chromaticity y
y (float) – tristimulus Y
- Return type
tuple
- raysect.optical.colour.ciexyz_to_ciexyy()¶
Performs conversion from CIE XYZ to CIE xyY colour space
Returns a tuple of (cx, cy, Y)
- Parameters
x (float) – tristimulus X
y (float) – tristimulus y
z (float) – tristimulus Z
- Return type
tuple
Raysect also supports conversion of CIE colour space values to standard sRGB colour space as defined by HP and Microsoft in 1996 as per IEC 61966-2-1:1999. For more information see Wikipedia.
- raysect.optical.colour.ciexyz_to_srgb()¶
Convert CIE XYZ values to sRGB colour space.
x, y, z in range [0, 1] r, g, b in range [0, 1]
- Parameters
x (float) – tristimulus X
y (float) – tristimulus y
z (float) – tristimulus Z
- Return type
tuple
- raysect.optical.colour.srgb_to_ciexyz()¶
Convert sRGB values to CIE XYZ colour space.
r, g, b in range [0, 1] x, y, z in range [0, 1]
- Parameters
r (float) – Red value
g (float) – Green value
b (float) – Blue value
- Return type
tuple
3.1.4. Optical Scenegraph¶
- class raysect.optical.scenegraph.world.World¶
Bases:
raysect.core.scenegraph.world.World
The root node of the optical scene-graph.
Inherits a lot of functionality and attributes from the core World object.
The world node tracks all primitives and observers in the world. It maintains acceleration structures to speed up the ray-tracing calculations. The particular acceleration algorithm used is selectable. The default acceleration structure is a kd-tree.
- Parameters
name – A string defining the node name.
- build_importance()¶
This method manually triggers a rebuild of the importance manager object.
If the importance manager object is already in a consistent state this method will do nothing unless the force keyword option is set to True.
- Parameters
force (bint) – If set to True, forces rebuilding of acceleration structure.
- has_important_primitives()¶
Returns true if any primitives in this scene-graph have an importance weighting.
- Return type
bool
- important_direction_pdf()¶
Calculates the value of the PDF for the specified sample point and direction.