HolographicRulingSpacing#

class optika.rulings.HolographicRulingSpacing(x1, x2, wavelength, is_diverging_1=True, is_diverging_2=False, transformation=None)[source]#

Bases: AbstractRulingSpacing

Rulings created by interfering two beams.

Examples

Create some holographic rulings from two source points, launch rays from the first source point and confirm they are refocused onto the second source point.

import matplotlib.pyplot as plt
import astropy.units as u
import astropy.visualization
import named_arrays as na
import optika

# Define the origins of the two recording beams
x1 = na.Cartesian3dVectorArray(5, 0, -10) * u.mm
x2 = na.Cartesian3dVectorArray(10, 0, -10) * u.mm

# Define the wavelength of the recording beams
wavelength = 500 * u.nm

# Define the surface normal
normal = na.Cartesian3dVectorArray(0, 0, -1)

# Define input rays emanating from the origin of the first recording beam.
position = na.Cartesian3dVectorArray(
    x=na.linspace(-5, +5, axis="x", num=5) * u.mm,
    y=0 * u.mm,
    z=0 * u.mm,
)
direction_input = position - x1

# Initialize the holographic ruling spacing representation.
rulings = optika.rulings.HolographicRulingSpacing(
    x1=x1,
    x2=x2,
    wavelength=wavelength,
)

# Evaluate the ruling spacing where the rays strike the surface.
d = rulings(position, normal)

# Compute the effective propagation direction of the input rays
direction_input = optika.rulings.incident_effective(
    wavelength=wavelength,
    direction=direction_input.normalized,
    index_refraction=1,
    normal=normal,
    diffraction_order=1,
    spacing_rulings=d.length,
    normal_rulings=d.normalized,
)

# Compute the output direction of the diffracted rays.
direction_output = optika.materials.snells_law(
    direction=direction_input,
    index_refraction=1,
    index_refraction_new=1,
    normal=normal,
    is_mirror=True,
)
direction_output = direction_output * 20 * u.mm

# Plot the results
with astropy.visualization.quantity_support():
    fig, ax = plt.subplots()
    na.plt.plot(
        na.stack([x1, position], axis="t"),
        components=("z", "x"),
        axis="t",
        color="tab:blue",
    )
    na.plt.plot(
        na.stack([position, position + direction_output], axis="t"),
        components=("z", "x"),
        axis="t",
        color="tab:orange",
    )
    ax.axvline(0)
    ax.scatter(x1.z, x1.x)
    ax.scatter(x2.z, x2.x)
../_images/optika.rulings.HolographicRulingSpacing_0_0.png

Notes

From Welford [1975], the ruling spacing is given by

\[\mathbf{d} = \frac{\lambda}{a} \mathbf{q} \times \mathbf{n}\]

where \(\lambda\) is the wavelength of the recording beam, \(a\) is an arbitrary scalar to be solved for, \(\mathbf{n}\) is a unit vector perpendicular to the surface, and \(\mathbf{q}\) is a unit vector parallel to the rulings. The vector \(a \mathbf{q}\) is given by

\[a \mathbf{q} = \mathbf{n} \times (\pm \mathbf{r}_1 \mp \mathbf{r}_2),\]

where \(\mathbf{r}_1\) is a unit vector in the direction of the first recording beam, and \(\mathbf{r}_2\) is a unit vector in the direction of the second recording beam. If rays are diverging from the origin of the recording beams, the top branch is used, otherwise the bottom branch is used. Since \(\mathbf{q}\) is a unit vector,

\[a = | a \mathbf{q} |.\]

Attributes

is_diverging_1

A boolean flag indicating if rays are diverging from the origin of the first beam.

is_diverging_2

A boolean flag indicating if rays are diverging from the origin of the second beam.

shape

The array shape of this object.

transformation

A transformation from surface coordinates to ruling coordinates.

x1

The origin of the first recording beam in local coordinates.

x2

The origin of the second recording beam in local coordinates.

wavelength

The wavelength of the recording beams.

Methods

__init__(x1, x2, wavelength[, ...])

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

Inheritance diagram of optika.rulings.HolographicRulingSpacing
Parameters:
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.

Parameters:

prefix (None | str) – an optional string, the length of which is used to calculate how much whitespace to add to the result.

Return type:

str

is_diverging_1: bool | AbstractScalar = True#

A boolean flag indicating if rays are diverging from the origin of the first beam.

is_diverging_2: bool | AbstractScalar = False#

A boolean flag indicating if rays are diverging from the origin of the second beam.

property shape: dict[str, int]#

The array shape of this object.

transformation: None | AbstractTransformation = None#

A transformation from surface coordinates to ruling coordinates.

wavelength: Quantity | AbstractScalar#

The wavelength of the recording beams.

x1: AbstractCartesian3dVectorArray#

The origin of the first recording beam in local coordinates.

x2: AbstractCartesian3dVectorArray#

The origin of the second recording beam in local coordinates.