RectangularRulings#
- class optika.rulings.RectangularRulings(spacing, depth, ratio_duty, diffraction_order)[source]#
Bases:
AbstractRulingsA ruling profile described by a rectangular wave.
Examples
Compute the 1st-order groove efficiency of rectangular rulings with a groove density of 2500 grooves/mm, a groove depth of 15 nm, and a duty cycle of 30 percent.
import numpy as np import matplotlib.pyplot as plt import astropy.units as u import named_arrays as na import optika # Define the groove density density = 2500 / u.mm # Define the groove depth depth = 15 * u.nm # Define ruling model rulings = optika.rulings.RectangularRulings( spacing=1 / density, depth=depth, ratio_duty=0.3, diffraction_order=1, ) # Define the wavelengths at which to sample the groove efficiency wavelength = na.geomspace(100, 1000, axis="wavelength", num=1001) * u.AA # Define the incidence angles at which to sample the groove efficiency angle = na.linspace(0, 30, num=3, axis="angle") * u.deg # Define the light rays incident on the grooves rays = optika.rays.RayVectorArray( wavelength=wavelength, direction=na.Cartesian3dVectorArray( x=np.sin(angle), y=0, z=np.cos(angle), ), ) # Compute the efficiency of the grooves for the given wavelength efficiency = rulings.efficiency( rays=rays, normal=na.Cartesian3dVectorArray(0, 0, -1), ) # Plot the groove efficiency as a function of wavelength fig, ax = plt.subplots() angle_str = angle.value.astype(str).astype(object) na.plt.plot( wavelength, efficiency, ax=ax, axis="wavelength", label=r"$\theta$ = " + angle_str + f"{angle.unit:latex_inline}", ); ax.set_xlabel(f"wavelength ({wavelength.unit:latex_inline})"); ax.set_ylabel(f"efficiency"); ax.legend();
Attributes
Depth of the ruling pattern.
The diffraction order to simulate.
The duty cycle of the ruling pattern.
The array shape of this object.
Spacing between adjacent rulings at the given position.
A normalized version of
spacingthat is guaranteed to be an instance ofoptika.rulings.AbstractRulingSpacing.Methods
__init__(spacing, depth, ratio_duty, ...)efficiency(rays, normal)The fraction of light diffracted into a given order.
incident_effective(rays, normal)Compute the effective propagation direction of the given rays using
incident_effective().to_string([prefix])Public-facing version of the
__repr__method that allows for defining a prefix string, which can be used to calculate how much whitespace to add to the beginning of each line of the result.Inheritance Diagram

- Parameters:
spacing (Quantity | AbstractScalar | AbstractRulingSpacing)
depth (Quantity | AbstractScalar)
ratio_duty (Quantity | AbstractScalar)
diffraction_order (int | AbstractScalar)
- efficiency(rays, normal)[source]#
The fraction of light diffracted into a given order.
Calculated using the expression given in Table 1 of Magnusson and Gaylord [1978].
- Parameters:
rays (RayVectorArray) – The light rays incident on the rulings
normal (AbstractCartesian3dVectorArray) – The vector normal to the surface on which the rulings are placed.
- Return type:
Notes
The theoretical efficiency of thin (wavelength much smaller than the groove spacing), rectangular rulings is given by Table 1 of Magnusson and Gaylord [1978],
\[\begin{split}\eta_i = \begin{cases} 1 - [(2 a / \pi) - (a / \pi)^2] \sin^2 \left\{ \pi \gamma / [2 (1 - \cos a)]^{1/2} \right\} & i = 0 \\ [2 / (i \pi)^2](1 - \cos i a) \sin^2 \left\{ \pi \gamma / [2 (1 - \cos a)]^{1/2} \right\} & i \ne 0, \\ \end{cases}\end{split}\]where \(\eta_i\) is the groove efficiency for diffraction order \(i\), \(a\) \((0 < a < 2 \pi)\) is the duty cycle of the rectangular wave, \(\gamma = \pi d n_1 / \lambda \cos \theta\) is the normalized amplitude of the fundamental grating, \(d\) is the thickness of the grating, \(n_1\) is the amplitude of the fundamental grating, \(\lambda\) is the free-space wavelength of the incident light, and \(\theta\) is the angle of incidence inside the medium.
- incident_effective(rays, normal)#
Compute the effective propagation direction of the given rays using
incident_effective().- Parameters:
rays (RayVectorArray) – The light rays incident on the rulings
normal (AbstractCartesian3dVectorArray) – The vector normal to the surface on which the rulings are placed.
- Return type:
- to_string(prefix=None)#
Public-facing version of the
__repr__method that allows for defining a prefix string, which can be used to calculate how much whitespace to add to the beginning of each line of the result.
- depth: Quantity | AbstractScalar = <dataclasses._MISSING_TYPE object>#
Depth of the ruling pattern.
- diffraction_order: int | AbstractScalar = <dataclasses._MISSING_TYPE object>#
The diffraction order to simulate.
- ratio_duty: Quantity | AbstractScalar = <dataclasses._MISSING_TYPE object>#
The duty cycle of the ruling pattern.
- spacing: Quantity | AbstractScalar | AbstractRulingSpacing = <dataclasses._MISSING_TYPE object>#
Spacing between adjacent rulings at the given position.
- property spacing_: AbstractRulingSpacing#
A normalized version of
spacingthat is guaranteed to be an instance ofoptika.rulings.AbstractRulingSpacing.