Skip to content

API Reference

Crystal Modelling

citrine.Magnitude

Bases: Enum

Enum to define the magnitude prefixes for units such as pico, nano, micro, etc.

citrine.AngularFrequency dataclass

AngularFrequency(value: float | ndarray)

Class to represent angular frequency.

Attributes:

  • value (Union[float, ndarray]) –

    The angular frequency.

citrine.Wavelength dataclass

Wavelength(value: float | ndarray, unit: Magnitude = Magnitude.nano)

Class to represent wavelength with unit conversion.

Attributes:

  • value (Union[float, ndarray]) –

    The wavelength value.

  • unit (Magnitude) –

    The unit of the wavelength (default: nano).

Methods:

Functions

as_angular_frequency

as_angular_frequency() -> AngularFrequency

Convert the wavelength to angular frequency.

Returns:

  • AngularFrequency ( AngularFrequency ) –

    Angular frequency corresponding to the wavelength.

as_wavevector

as_wavevector(refractive_index: float = 1.0) -> float | NDArray[np.floating]

Convert the wavelength to a wavevector.

Returns:

  • float ( float | NDArray[floating] ) –

    The wavevector.

to_absolute

to_absolute() -> Wavelength

Convert the wavelength to base units (meters).

Returns:

  • Wavelength ( Wavelength ) –

    Wavelength in meters.

to_unit

to_unit(new_unit: Magnitude) -> Wavelength

Convert the wavelength to a new unit.

Parameters:

  • new_unit (Magnitude) –

    The desired unit for the wavelength.

Returns:

  • Wavelength ( Wavelength ) –

    New Wavelength object with converted units.

citrine.spectral_window

spectral_window(central_wavelength: Wavelength, spectral_width: Wavelength, steps: int, reverse: bool = False) -> Wavelength

Generate an array of wavelengths within a specified spectral window.

Parameters:

  • central_wavelength (Wavelength) –

    The central wavelength.

  • spectral_width (Wavelength) –

    The total spectral width.

  • steps (int) –

    The number of steps for the wavelength range.

Returns:

  • Wavelength ( Wavelength ) –

    An array of wavelengths within the specified window.

citrine.SellmeierCoefficientsSimple dataclass

SellmeierCoefficientsSimple(coefficients: list[float] | NDArray[floating], temperature: float, dn_dt: float | None = None)

citrine.SellmeierCoefficientsTemperatureDependent dataclass

SellmeierCoefficientsTemperatureDependent(first_order: list[float] | NDArray | None, second_order: list[float] | NDArray | None, temperature: float, zeroth_order: list[float] | NDArray | None = None)

citrine.SellmeierCoefficientsJundt dataclass

SellmeierCoefficientsJundt(a_terms: list[float] | NDArray[floating], b_terms: list[float] | NDArray[floating], temperature: float)

citrine._permittivity

_permittivity(sellmeier: list[float] | NDArray, wavelength_um: float | NDArray[floating]) -> float | NDArray[np.floating]

Compute the permittivity using the Sellmeier equation.

This modification adds comments to explain each step of the Sellmeier equation calculation. The code is implementing the Sellmeier equation, which is used to determine the refractive index of a material as a function of wavelength.

The equation typically takes the form:

\[ n^2 = A + \frac{B * λ^2}{λ^2 - C} + \frac{D * λ^2}{λ^2 - E} + \ldots \]

Where n is the refractive index, λ is the wavelength, and A, B, C, D, E, etc. are the Sellmeier coefficients specific to the material.

Parameters:

  • sellmeier (Union[List[float], NDArray]) –

    Sellmeier coefficients.

  • wavelength_um (float) –

    Wavelength in micrometers.

Returns:

  • float ( float | NDArray[floating] ) –

    The permittivity value.

Source code in src/citrine/citrine.py
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
def _permittivity(
    sellmeier: Union[List[float], NDArray],
    wavelength_um: Union[float, NDArray[np.floating]],
) -> Union[float, NDArray[np.floating]]:
    """
    Compute the permittivity using the Sellmeier equation.

    This modification adds comments to explain each step of the Sellmeier
    equation calculation. The code is implementing the Sellmeier equation,
    which is used to determine the refractive index of a material as a
    function of wavelength.

    The equation typically takes the form:

    $$
    n^2 = A + \\frac{B * λ^2}{λ^2 - C} + \\frac{D * λ^2}{λ^2 - E} + \\ldots
    $$

    Where n is the refractive index, λ is the wavelength, and A, B, C, D, E,
    etc. are the Sellmeier coefficients specific to the material.

    Args:
        sellmeier (Union[List[float], NDArray]): Sellmeier coefficients.
        wavelength_um (float): Wavelength in micrometers.

    Returns:
        float: The permittivity value.
    """

    # Take array of sellemeier coefficients of the form [A, B, C, D, E, ...]
    # and split into A, [B, C, ...], E
    first, *coeffs, last = sellmeier
    assert len(coeffs) % 2 == 0

    # Square the wavelength (in micrometers)
    wl = wavelength_um**2

    # Initialize the permittivity with the first coefficient
    p = first

    # Iterate through the Sellmeier coefficients in pairs
    # coeffs[0::2] correspond to B, D, ... -> numerator
    # coeffs[1::2] correspond to C, E, ... -> denominator
    for number, denom in zip(coeffs[0::2], coeffs[1::2]):
        # Apply the Sellmeier equation term:
        # Add (numerator) / (1 - denominator / wavelength^2)
        p += number / (1 - denom / wl)

    # Apply the final term of the Sellmeier equation
    # Subtract the last coefficient multiplied by the squared wavelength
    p -= last * wl
    return p

citrine._n_0

_n_0(sellmeier: list[float] | NDArray, wavelength_um: float | NDArray[floating]) -> float | NDArray[np.floating]

Calculate the refractive index (n_0) using the zeroth-order Sellmeier coefficients.

Parameters:

  • sellmeier (Union[List[float], NDArray]) –

    Zeroth-order Sellmeier coefficients.

  • wavelength_um (float) –

    Wavelength in micrometers.

Returns:

  • float ( float | NDArray[floating] ) –

    The refractive index (n_0).

Source code in src/citrine/citrine.py
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
def _n_0(
    sellmeier: Union[List[float], NDArray],
    wavelength_um: Union[float, NDArray[np.floating]],
) -> Union[float, NDArray[np.floating]]:
    """
    Calculate the refractive index (n_0) using the zeroth-order Sellmeier
        coefficients.

    Args:
        sellmeier (Union[List[float], NDArray]): Zeroth-order Sellmeier
            coefficients.
        wavelength_um (float): Wavelength in micrometers.

    Returns:
        float: The refractive index (n_0).
    """
    return np.sqrt(_permittivity(sellmeier, wavelength_um))

citrine._n_i

_n_i(sellmeier: list[float] | NDArray[floating], wavelength_um: float | NDArray[floating]) -> float | NDArray[np.floating]

Calculate the refractive index (n_i) using the higher-order Sellmeier coefficients.

Parameters:

  • sellmeier (Union[List[float], NDArray]) –

    Higher-order Sellmeier coefficients.

  • wavelength_um (float) –

    Wavelength in micrometers.

Returns:

  • float ( float | NDArray[floating] ) –

    The refractive index (n_i).

Source code in src/citrine/citrine.py
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
def _n_i(
    sellmeier: Union[List[float], NDArray[np.floating]],
    wavelength_um: Union[float, NDArray[np.floating]],
) -> Union[float, NDArray[np.floating]]:
    """
    Calculate the refractive index (n_i) using the higher-order Sellmeier
        coefficients.

    Args:
        sellmeier (Union[List[float], NDArray]): Higher-order Sellmeier
            coefficients.
        wavelength_um (float): Wavelength in micrometers.

    Returns:
        float: The refractive index (n_i).
    """
    n = 0
    for i, s in enumerate(sellmeier):
        n += s / (wavelength_um**i)
    return n

citrine.Orientation

Bases: Enum

Enum to represent the orientation: ordinary or extraordinary.

citrine.PhaseMatchingCondition

Bases: Enum

Phase-matching configurations for nonlinear optical processes.

Each configuration defines the polarization orientations for (pump, signal, idler) where:

- o/O: ordinary polarization (horizontal, H)

- e/E: extraordinary polarization (vertical, V)
The configurations are
  • Type 0 (o): All waves ordinary polarized (H,H,H), type0_o
  • Type 0 (e): All waves extraordinary polarized (V,V,V), type0_e
  • Type 1: Pump extraordinary, signal/idler ordinary (V,H,H), type1
  • Type 2 (o): Mixed polarizations (V,H,V), type2_o
  • Type 2 (e): Mixed polarizations (V,V,H), type2_e
Note
  • Ordinary (o) corresponds to horizontal (H) polarization
  • Extraordinary (e) corresponds to vertical (V) polarization
  • The tuple order is always (pump, signal, idler)

citrine.Crystal

Crystal(name: str, sellmeier_o: SellmeierCoefficientsSimple | SellmeierCoefficientsTemperatureDependent | SellmeierCoefficientsJundt | SellmeierCoefficientsBornAndWolf, sellmeier_e: SellmeierCoefficientsSimple | SellmeierCoefficientsTemperatureDependent | SellmeierCoefficientsJundt | SellmeierCoefficientsBornAndWolf, phase_matching: PhaseMatchingCondition, doi: str = None)

Class to represent a nonlinear crystal for refractive index and phase matching calculations.

Attributes:

  • name (str) –

    Name of the crystal.

  • sellmeier_o (SellmeierCoefficients) –

    Ordinary Sellmeier coefficients.

  • sellmeier_e (SellmeierCoefficients) –

    Extraordinary Sellmeier coefficients.

  • phase_matching (PhaseMatchingCondition) –

    Phase matching conditions.

Methods:

  • refractive_index

    Calculate the refractive index for a given wavelength and polarization

  • refractive_indices

    Calculate the refractive indices for the pump, signal, and idler

Attributes:

Attributes

phase_matching property writable

phase_matching: PhaseMatchingCondition

Get the current phase matching condition.

Returns:

Functions

refractive_index

refractive_index(wavelength: Wavelength, photon: Photon, temperature: float | None = None) -> float | NDArray[np.floating]

Calculate the refractive index for a given wavelength and polarization using the Sellmeier equation.

Parameters:

  • wavelength (Wavelength) –

    Wavelength.

  • photon (Photon) –

    Photon type.

  • temperature (Optional[float], default: None ) –

    Temperature in Celsius (defaults to reference temperature).

Returns:

  • float ( float | NDArray[floating] ) –

    Refractive index.

refractive_indices

refractive_indices(pump_wavelength: Wavelength, signal_wavelength: Wavelength, idler_wavelength: Wavelength, temperature: float | None = None) -> tuple[float | NDArray[np.floating], float | NDArray[np.floating], float | NDArray[np.floating]]

Calculate the refractive indices for the pump, signal, and idler photons.

Parameters:

  • pump_wavelength (Wavelength) –

    Pump photon wavelength.

  • signal_wavelength (Wavelength) –

    Signal photon wavelength.

  • idler_wavelength (Wavelength) –

    Idler photon wavelength.

  • temperature (Optional[float], default: None ) –

    Temperature in Celsius (defaults to reference temperature).

Returns:

  • tuple[float | NDArray[floating], float | NDArray[floating], float | NDArray[floating]]

    Tuple[float, float, float]: Refractive indices for the pump, signal, and idler photons.

citrine.calculate_grating_period

calculate_grating_period(lambda_p_central: Wavelength, lambda_s_central: Wavelength, lambda_i_central: Wavelength, crystal: Crystal, temperature: float = None) -> float | NDArray[np.floating]

Calculate the grating period (Λ) for the phase matching condition.

Parameters:

  • lambda_p_central (Wavelength) –

    Central wavelength of the pump.

  • lambda_s_central (Wavelength) –

    Central wavelength of the signal.

  • lambda_i_central (Wavelength) –

    Central wavelength of the idler.

Returns:

  • float ( float | NDArray[floating] ) –

    Grating period in microns.

citrine.delta_k_matrix

delta_k_matrix(lambda_p: Wavelength, lambda_s: Wavelength, lambda_i: Wavelength, crystal: Crystal, temperature: float = None) -> float | NDArray[np.floating]

Calculate the Δk matrix for the phase matching function using the wavevector

Parameters:

  • lambda_p (Wavelength) –

    Central wavelength of the pump.

  • lambda_s (Wavelength) –

    Signal wavelengths (Wavelength type).

  • lambda_i (Wavelength) –

    Idler wavelengths (Wavelength type).

  • grating_period (float) –

    Grating period in microns.

Returns:

  • float | NDArray[floating]

    np.ndarray: A 2D matrix of Δk values.

citrine.phase_mismatch

phase_mismatch(lambda_p: Wavelength, lambda_s: Wavelength, crystal: Crystal, temperature: float, grating_period: float)

citrine.phase_matching_function

phase_matching_function(delta_k: ndarray, grating_period, crystal_length) -> np.ndarray

Compute the phase matching function as a sinc function of Δk.

Parameters:

  • delta_k (ndarray) –

    The Δk matrix from the phase matching condition.

Returns:

  • ndarray

    np.ndarray: The phase matching function.

citrine.PhotonType

Bases: Enum

citrine.calculate_marginal_spectrum

calculate_marginal_spectrum(jsa_matrix: ndarray, lambda_s: Wavelength, lambda_i: Wavelength, photon_type: PhotonType) -> tuple[np.ndarray, np.ndarray]

Calculate the marginal spectrum for either signal or idler photons from a JSA matrix.

Parameters:

  • jsa_matrix (ndarray) –

    A 2D complex matrix representing the joint spectral amplitude. Shape (len(lambda_s), len(lambda_i)) with indexing='ij'.

  • lambda_s (Wavelength) –

    Signal wavelengths (Wavelength type).

  • lambda_i (Wavelength) –

    Idler wavelengths (Wavelength type).

  • photon_type (PhotonType) –

    Enum specifying whether to calculate for SIGNAL or IDLER photon.

Returns:

  • tuple[ndarray, ndarray]

    Tuple containing: - Array of wavelengths for the selected photon - Array of corresponding intensity values (complex)

citrine.joint_spectral_amplitude

joint_spectral_amplitude(phase_mismatch_matrix: ndarray, pump_envelope_matrix: ndarray, normalisation: bool = True) -> np.ndarray

citrine.bandwidth_conversion

bandwidth_conversion(delta_lambda_FWHM: Wavelength, pump_wl: Wavelength) -> float

Convert pump bandwidth from FWHM in wavelength to frequency.

Parameters:

  • delta_lambda_FWHM (Wavelength) –

    Bandwidth in wavelength units.

  • pump_wl (Wavelength) –

    Central pump wavelength.

Returns:

  • float ( float ) –

    Converted bandwidth.

citrine.Time dataclass

Time(value: float | ndarray, unit: Magnitude = Magnitude.base)

Class to represent time in chosen units

Attributes:

  • value (Union[float, ndarray]) –

    The time in chosen units

  • unit (Magnitude) –

    The unit of time

Methods:

  • to_absolute

    Convert the time to base units (seconds).

  • to_unit

    Convert the time to a new unit.

Functions

to_absolute

to_absolute() -> Time

Convert the time to base units (seconds).

Returns:

  • Time ( Time ) –

    Time in seconds.

to_unit

to_unit(new_unit: Magnitude) -> Time

Convert the time to a new unit.

Parameters:

  • new_unit (Magnitude) –

    The desired unit for the Time.

Returns:

  • time ( Time ) –

    New Time object with converted units.

citrine.Bunching

Bases: Enum

citrine.hom_interference_from_jsa

hom_interference_from_jsa(joint_spectral_amplitude: ndarray, wavelengths_signal: Wavelength, wavelengths_idler: Wavelength, time_delay: Time = Time(0.0, Magnitude.base), bunching: Bunching = Bunching.Bunching) -> tuple[float, np.ndarray]

Calculates the Hong-Ou-Mandel (HOM) interference pattern and total coincidence rate

This function computes the quantum interference that occurs when two photons enter a 50:50 beam splitter from different inputs. The interference pattern depends on the joint spectral amplitude (JSA) of the photon pair and the time delay between their arrivals at the beam splitter.

The calculation follows these steps: 1. Create the second JSA with swapped signal and idler (anti-transpose) 2. Apply a phase shift based on the time delay and wavelength difference 3. Calculate the interference between the two quantum pathways 4. Compute the joint spectral intensity (JSI) by taking the squared magnitude

Parameters:

  • jsa (ndarray) –

    Complex 2D numpy array of shape (M, N) representing the Joint Spectral Amplitude. The first dimension corresponds to signal wavelengths and the second dimension to idler wavelengths.

  • wavelengths_signal (Wavelength) –

    1D array of signal wavelengths in meters, must match the first dimension of jsa.

  • wavelengths_idler (Wavelength) –

    1D array of idler wavelengths in meters, must match the second dimension of jsa. Note: should be in descending order to match standard convention.

  • time_delay (float = 0.0, default: Time(0.0, base) ) –

    Time delay between signal and idler photons in seconds. At zero delay, maximum quantum interference occurs for indistinguishable photons.

  • bunching (Bunching = Bunching.Bunching, default: Bunching ) –

    Type of interference to be simulated

Returns:

  • tuple[float, ndarray]

    Tuple containing: - float: Total coincidence probability (sum of all JSI values). - np.ndarray: 2D array representing the joint spectral intensity after interference.

Notes

Hong-Ou-Mandel interference results from the destructive interference between two indistinguishable two-photon quantum states. For perfectly indistinguishable photons at zero time delay, they will always exit the beam splitter together from the same port, resulting in zero coincidence count rate (the HOM dip).

The time delay introduces a wavelength-dependent phase shift to one pathway, reducing the interference and increasing the coincidence probability.

citrine.hong_ou_mandel_interference

hong_ou_mandel_interference(jsa: ndarray, signal_wavelengths: ndarray, idler_wavelengths: ndarray, time_delays: Time = Time(0.0, Magnitude.base), bunching: Bunching = Bunching.Bunching) -> np.ndarray

Calculates the complete Hong-Ou-Mandel dip by scanning the time delays.

This function generates the characteristic HOM dip by calculating the coincidence rate at each time delay value in the provided array.

Parameters:

  • jsa (ndarray) –

    Complex 2D numpy array representing the Joint Spectral Amplitude.

  • signal_wavelengths (ndarray) –

    1D array of signal wavelengths in meters.

  • idler_wavelengths (ndarray) –

    1D array of idler wavelengths in meters.

  • time_delays (Time, default: Time(0.0, base) ) –

    1D array of time delay values in seconds to scan across.

Returns:

  • ndarray

    np.ndarray: 1D array of coincidence rates corresponding to each time delay, showing the characteristic HOM dip.

Physics Notes

The width of the HOM dip is inversely proportional to the spectral bandwidth of the photons. Spectrally narrower photons produce wider dips, while broader bandwidth photons produce narrower dips, demonstrating the time-frequency uncertainty principle.

The shape and visibility of the dip reveals information about the spectral entanglement and distinguishability of the photon pairs.

citrine.spectral_purity

spectral_purity(jsa: ndarray) -> tuple[float, float, float]

Calculate the spectral purity of a biphoton state

This function calculate the spectral purity, schmidt number and entropy of a biphoton state. This is achieved by performing a singular value decomposition (SVD) a discrete analogue to the Schmidt decomposition (continuous space).

Parameters:

  • joint_spectral_amplitude

    2D array representing the joint spectral amplitude.

Returns:

  • tuple ( tuple[float, float, float] ) –

    A tuple containing: - probabilities (NDArray): Normalized HOM interference pattern - purity (float): Spectral purity of the state (0 to 1) - schmidt_number (float): Schmidt number indicating mode entanglement

Notes

The function performs these steps: 1. Singular value decomposition of the JSA 2. Computes quantum metrics (purity, Schmidt number, entropy)

Phase Matching Functions

citrine.pmf_gaussian

pmf_gaussian(delta_k: ndarray, poling_period, crystal_length) -> np.ndarray

Calculate the Gaussian phase-matching function.

Parameters:

  • delta_k (ndarray) –

    The delta k values.

  • poling_period (float) –

    The poling period of the crystal.

  • crystal_length (float) –

    The length of the crystal.

Returns:

  • ndarray

    np.ndarray: The Gaussian phase-matching function values.

Source code in src/citrine/phase_matching.py
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
def pmf_gaussian(
    delta_k: np.ndarray, poling_period, crystal_length
) -> np.ndarray:
    """
    Calculate the Gaussian phase-matching function.

    Args:
        delta_k (np.ndarray): The delta k values.
        poling_period (float): The poling period of the crystal.
        crystal_length (float): The length of the crystal.

    Returns:
        np.ndarray: The Gaussian phase-matching function values.
    """
    _delta_k = _adjusted_delta_k(delta_k, poling_period)
    return np.exp((1j * _delta_k * (crystal_length / 2)) ** 2)

citrine.pmf_antisymmetric

pmf_antisymmetric(delta_k: ndarray, poling_period, crystal_length) -> np.ndarray

Calculate the antisymmetric phase-matching function.

Parameters:

  • delta_k (ndarray) –

    The delta k values.

  • poling_period (float) –

    The poling period of the crystal.

  • crystal_length (float) –

    The length of the crystal.

Returns:

  • ndarray

    np.ndarray: The antisymmetric phase-matching function values.

Source code in src/citrine/phase_matching.py
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
def pmf_antisymmetric(
    delta_k: np.ndarray, poling_period, crystal_length
) -> np.ndarray:
    """
    Calculate the antisymmetric phase-matching function.

    Args:
        delta_k (np.ndarray): The delta k values.
        poling_period (float): The poling period of the crystal.
        crystal_length (float): The length of the crystal.

    Returns:
        np.ndarray: The antisymmetric phase-matching function values.
    """
    _delta_k = _adjusted_delta_k(delta_k, poling_period)
    return np.exp((1j * _delta_k * (crystal_length / 2)) ** 2) * _delta_k

Pump Envelope Functions

citrine.gaussian

gaussian(lambda_p: Wavelength, sigma_p: float, lambda_s: Wavelength, lambda_i: Wavelength, k: float = None) -> np.ndarray

Generate a pump envelope matrix with a Gaussian profile.

Implements a Gaussian pump envelope in frequency space according to: $$ \alpha(\omega_s + \omega_i) = \exp\left(-\frac{(\omega_s + \omega_i - \omega_p)^2}{\sigma_p^2}\right) $$

Can be chirped according to: $$ \chi(\omega_s + \omega_i) = \exp(-ik(\omega_s + \omega_i - \omega_p)^2), $$

returning \(\alpha\omega_s + \omega_i \cdot \chi\omega_s + \omega_i\)

Parameters:

  • lambda_p (Wavelength) –

    Central wavelength of the pump.

  • sigma_p (float) –

    Pump bandwidth.

  • lambda_s (Wavelength) –

    Signal wavelengths array.

  • lambda_i (Wavelength) –

    Idler wavelengths array.

  • k (float, default: None ) –

    Chirp parameter, controls strength and direction of chirp. Positive values create up-chirp, negative values create down-chirp.

Returns:

  • ndarray

    A 2D pump envelope matrix.

Notes

The equation uses the following variables: \(\omega_s\): signal frequency \(\omega_i\): idler frequency \(\omega_p\): pump central frequency \(\sigma_p\): pump bandwidth

citrine.sech2

sech2(lambda_p: Wavelength, sigma_p: float, lambda_s: Wavelength, lambda_i: Wavelength, k: float = None) -> np.ndarray

Generate a pump envelope matrix with a hyperbolic secant squared profile.

Implements a sech² pump envelope in frequency space according to: $$ \alpha(\omega_s + \omega_i) = \text{sech}^2(\pi\sigma_p(\omega_s + \omega_i - \omega_p)) $$

Can be chirped according to: $$ \chi(\omega_s + \omega_i) = \exp(-ik(\omega_s + \omega_i - \omega_p)^2), $$

returning \(\alpha\omega_s + \omega_i \cdot \chi\omega_s + \omega_i\)

Parameters:

  • lambda_p (Wavelength) –

    Central wavelength of the pump.

  • sigma_p (float) –

    Pump bandwidth parameter.

  • lambda_s (Wavelength) –

    Signal wavelengths array.

  • lambda_i (Wavelength) –

    Idler wavelengths array.

  • k (float, default: None ) –

    Chirp parameter, controls strength and direction of chirp. Positive values create up-chirp, negative values create down-chirp.

Returns:

  • ndarray

    A 2D pump envelope matrix.

Notes

The equation uses the following variables: \(\omega_s\): signal frequency \(\omega_i\): idler frequency \(\omega_p\): pump central frequency \(\sigma_p\): pump bandwidth parameter

citrine.skewed_gaussian

skewed_gaussian(lambda_p: Wavelength, sigma_p: float, lambda_s: Wavelength, lambda_i: Wavelength, alpha: float = 3.0, k: float = None) -> np.ndarray

Generate a pump envelope matrix with a skewed Gaussian profile.

Implements a skewed Gaussian pump envelope in frequency space according to: $$ \alpha(\omega_s + \omega_i) = \exp\left(-\frac{(\omega_s + \omega_i - \omega_p)^2}{2\sigma_p^2}\right) \cdot \left(1 + \text{erf}\left(\alpha \cdot \frac{\omega_s + \omega_i - \omega_p}{\sigma_p \sqrt{2}}\right)\right) $$

Can be chirped according to: $$ \chi(\omega_s + \omega_i) = \exp(-ik(\omega_s + \omega_i - \omega_p)^2), $$

returning \(\alpha(\omega_s + \omega_i) \cdot \chi(\omega_s + \omega_i)\)

Parameters:

  • lambda_p (Wavelength) –

    Central wavelength of the pump.

  • sigma_p (float) –

    Pump bandwidth.

  • lambda_s (Wavelength) –

    Signal wavelengths array.

  • lambda_i (Wavelength) –

    Idler wavelengths array.

  • alpha (float, default: 3.0 ) –

    Skewness parameter (default 3.0). Positive values skew toward lower frequencies (longer wavelengths), negative values skew toward higher frequencies (shorter wavelengths).

  • k (float, default: None ) –

    Chirp parameter, controls strength and direction of chirp. Positive values create up-chirp, negative values create down-chirp.

Returns:

  • ndarray

    A 2D pump envelope matrix.

Notes

The skewed Gaussian combines a standard Gaussian with the error function (erf) to create asymmetry while maintaining smoothness.

Crystal Library

citrine.KTiOPO4_Fradkin module-attribute

KTiOPO4_Fradkin: Incomplete = Crystal(name='Potassium Titanyl Phosphate', doi='10.1063/1.123408', sellmeier_o=SellmeierCoefficientsTemperatureDependent(zeroth_order=array([2.12725, 1.18431, 0.0514852, 0.6603, 100.00507, 0.00968956]), first_order=array([9.9587, 9.9228, -8.9603, 4.101]) * 1e-06, second_order=array([-1.1882, 10.459, -9.8136, 3.1481]) * 1e-08, temperature=25), sellmeier_e=SellmeierCoefficientsTemperatureDependent(zeroth_order=array([2.0993, 0.922683, 0.0467695, 0.0138408]), first_order=array([6.2897, 6.3061, -6.0629, 2.6486]) * 1e-06, second_order=array([-0.14445, 2.2244, -3.577, 1.347]) * 1e-08, temperature=25), phase_matching=type2_o)

citrine.KTiOPO4_Emanueli module-attribute

KTiOPO4_Emanueli: Incomplete = Crystal(name='Potassium Titanyl Phosphate', doi='10.1364/AO.42.006661', sellmeier_o=SellmeierCoefficientsTemperatureDependent(zeroth_order=None, first_order=array([6.2897, 6.3061, -6.0629, 2.6486]) * 1e-06, second_order=array([-0.14445, 2.2244, -3.577, 1.347]) * 1e-08, temperature=25), sellmeier_e=SellmeierCoefficientsTemperatureDependent(zeroth_order=None, first_order=array([9.9587, 9.9228, -8.9603, 4.101]) * 1e-06, second_order=array([-1.1882, 10.459, -9.8136, 3.1481]) * 1e-08, temperature=25), phase_matching=type2_o)

citrine.KTiOAsO4_Emanueli module-attribute

KTiOAsO4_Emanueli: Incomplete = Crystal(name='Potassium Titanyl Aresnate', doi='10.1364/AO.42.006661', sellmeier_o=SellmeierCoefficientsTemperatureDependent(zeroth_order=None, first_order=array([-4.1053, 44.261, -38.012, 11.302]) * 1e-06, second_order=array([0.5857, 3.9386, -4.0081, 1.4316]) * 1e-08, temperature=25), sellmeier_e=SellmeierCoefficientsTemperatureDependent(zeroth_order=None, first_order=array([-6.1537, 64.505, -56.447, 17.169]) * 1e-06, second_order=array([-0.96751, 13.192, -11.78, 3.6292]) * 1e-08, temperature=25), phase_matching=type2_o)

citrine.LiNbO3_Zelmon module-attribute

LiNbO3_Zelmon: Incomplete = Crystal(name='Lithium Niobate', doi='10.1364/JOSAB.14.003319', sellmeier_o=SellmeierCoefficientsBornAndWolf(coefficients=array([2.9804, 0.02047, 0.5981, 0.0666, 8.9543, 416.08]), temperature=21), sellmeier_e=SellmeierCoefficientsBornAndWolf(coefficients=array([2.6734, 0.01764, 1.229, 0.05914, 12.614, 474.6]), temperature=21), phase_matching=type0_e)

citrine.LiNbO3_5molMgO_Zelmon module-attribute

LiNbO3_5molMgO_Zelmon: Incomplete = Crystal(name='Lithium Niobate 5Mol Magnesium Doped', doi='10.1364/JOSAB.14.003319', sellmeier_o=SellmeierCoefficientsBornAndWolf(coefficients=array([2.2454, 0.01242, 1.3005, 0.0513, 6.8972, 331.33]) * 1e-06, temperature=21), sellmeier_e=SellmeierCoefficientsBornAndWolf(coefficients=array([2.4272, 0.01478, 1.4617, 0.005612, 9.6536, 371.216]) * 1e-06, temperature=21), phase_matching=type0_e)

citrine.LiNbO3_Newlight module-attribute

LiNbO3_Newlight: Incomplete = Crystal(name='Lithium Niobate', doi='', sellmeier_o=SellmeierCoefficientsSimple(coefficients=array([4.9048, 0.11768, -0.0475, -0.027169]), temperature=20, dn_dt=-8.74e-07), sellmeier_e=SellmeierCoefficientsSimple(coefficients=array([4.582, 0.099169, -0.04443, -0.02195]), temperature=20, dn_dt=3.9073e-05), phase_matching=type0_e)

citrine.LiNbO3_5molMgO_Gayer module-attribute

LiNbO3_5molMgO_Gayer: Incomplete = Crystal(name='Lithium Niobate 5Mol Magnesium Doped', doi=' 10.1007/s00340-008-2998-2', sellmeier_e=SellmeierCoefficientsJundt(a_terms=array([5.756, 0.0983, 0.202, 189.32, 12.52, 0.0132]), b_terms=array([2.86e-06, 4.7e-08, 6.113e-08, 0.0001516]), temperature=24.5), sellmeier_o=SellmeierCoefficientsJundt(a_terms=array([5.653, 0.1185, 0.2091, 89.61, 10.85, 0.0197]), b_terms=array([7.941e-07, 3.134e-08, -4.641e-09, -2.188e-06]), temperature=24.5), phase_matching=type0_e)

Brightness and Heralding

citrine.photon_pair_coupling_filter

photon_pair_coupling_filter(xi: float, alpha_sq: float, phi_0: float, np_over_ns: float, np_over_ni: float, rho_bounds: tuple[float, float] = (0, 2), theta_bounds: tuple[float, float] = (0, 2 * np.pi)) -> float

Calculate the photon pair coupling efficiency K2 for SPDC into single spatial modes.

This function implements Equation (26) from Smirr et al., which describes the spatial dependence of the two-photon coupling probability P2. The function K2 accounts for spatial interferences caused by both phase matching and coupling to the target mode.

The integrand Q2 contains three main physical terms: 1. Pump-target mode overlap (spatial filtering) 2. Signal-idler correlation from SPDC process 3. Phase matching with longitudinal and transverse mismatch

Parameters:

  • xi (float) –

    Focusing parameter ξ = L/(2zR), where L is crystal length and zR is Rayleigh range.

  • alpha_sq (float) –

    Squared normalized target mode waist parameter α² = (a0/w0)², where a0 is target mode waist and w0 is pump beam waist.

  • phi_0 (float) –

    Longitudinal phase mismatch φ0 = Δk0·L (Eq. 21).

  • np_over_ns (float) –

    Ratio of refractive indices np/ns (pump to signal).

  • np_over_ni (float) –

    Ratio of refractive indices np/ni (pump to idler).

  • rho_bounds (Tuple[float, float], default: (0, 2) ) –

    Integration bounds for normalized radial coordinates ρs, ρi. Defaults to (0, 2).

  • theta_bounds (Tuple[float, float], default: (0, 2 * pi) ) –

    Integration bounds for angular coordinate θ = θs - θi. Defaults to (0, 2π).

Returns:

  • float ( float ) –

    Dimensionless photon pair coupling factor K2/(kp0·L).

Note

Reference: Smirr et al., Eq. (26), page 5: K2 = |∫∫ d²κs/(2π)² ∫∫ d²κi/(2π)² S̃p(κs + κi, 0) × Õ0(κs, z0)Õ0(κi, z0) × e^(-iz0(np/n's |κs|²/kp0 + np/n'i |κi|²/kp0)) × L sinc(ΔK(κs, κi)L/2)|²

The normalized coordinates are defined as: - ρs = w0·κs/2, ρi = w0·κi/2 (Eq. 50) - The phase mismatch ΔK includes both longitudinal (φ0) and transverse terms (Eq. 21)

citrine.single_photon_coupling_filter

single_photon_coupling_filter(xi: float, alpha_sq: float, phi_0: float, np_over_ns: float, np_over_ni: float, rho_bounds: tuple[float, float] = (0, 2), theta_bounds: tuple[float, float] = (-np.pi, np.pi)) -> float

Calculate the single photon coupling efficiency K1 for heralding ratio calculation.

This function implements Equation (39) from Smirr et al., which describes the single-photon coupling probability P1. This is used to calculate the heralding ratio Γ2|1 = K2/K1 (Eq. 41).

The calculation involves: 1. Inner integral over signal photon modes (with spatial filtering) 2. Outer integral over idler photon modes (free space) 3. Phase matching for the combined signal-idler system

Parameters:

  • xi (float) –

    Focusing parameter ξ = L/(2zR).

  • alpha_sq (float) –

    Squared normalized target mode waist parameter α².

  • phi_0 (float) –

    Longitudinal phase mismatch φ0 = Δk0·L.

  • np_over_ns (float) –

    Ratio of refractive indices np/ns.

  • np_over_ni (float) –

    Ratio of refractive indices np/ni.

  • rho_bounds (Tuple[float, float], default: (0, 2) ) –

    Integration bounds for normalized radial coordinates. Defaults to (0, 2).

  • theta_bounds (Tuple[float, float], default: (-pi, pi) ) –

    Integration bounds for angular coordinate. Defaults to (-π, π).

Returns:

  • float ( float ) –

    Dimensionless single photon coupling factor K1/(kp0·L).

Note

Reference: Smirr et al., Eq. (39), page 6: K1 = ∫∫ d²κi/(2π)² |∫∫ d²κs/(2π)² S̃p(κs + κi, 0)Õ0(κs, z0) × e^(-iz0 np/n's |κs|²/kp0) L sinc(ΔK L/2)|²

K1 represents the probability of detecting at least one photon in the target spatial mode, which is needed to calculate the conditional probability (heralding ratio) of detecting the second photon given that the first was detected.

citrine.brightness_and_heralding

brightness_and_heralding(xi: float, alpha: float, phi_0: float, n_p: float, n_i: float, n_s: float, rho_max: float = 2.0) -> tuple[float, float, float]

Calculate brightness factors K1, K2 and heralding ratio for SPDC source optimization.

This function combines the calculations from Equations (26), (39), and (41) in Smirr et al. to provide the key figures of merit for SPDC source performance.

Parameters:

  • xi (float) –

    Focusing parameter ξ = L/(2zR), where L is crystal length and zR is pump beam Rayleigh range. Typical range: 0.1 to 10.

  • alpha (float) –

    Normalized target mode waist α = a0/w0, where a0 is target mode waist and w0 is pump beam waist. Optimal value around 1-2.

  • phi_0 (float) –

    Longitudinal phase mismatch φ0 = Δk0·L. Optimization parameter for temperature tuning.

  • n_p (float) –

    Pump refractive index at pump wavelength.

  • n_i (float) –

    Idler refractive index at idler wavelength.

  • n_s (float) –

    Signal refractive index at signal wavelength.

  • rho_max (float, default: 2.0 ) –

    Maximum integration bound for normalized radial coordinates. Defaults to 2.0, which provides good convergence.

Returns:

  • tuple[float, float, float]

    Tuple[float, float, float]: A tuple containing: - K1 (float): Single photon coupling factor K1/(kp0·L). Related to heralded single photon brightness.

    • K2 (float): Photon pair coupling factor K2/(kp0·L). Related to coincidence brightness.

    • heralding (float): Heralding ratio Γ2|1 = K2/K1. Measures the conditional probability of detecting the idler photon given detection of the signal photon.

Note

References: - K2: Smirr et al., Eq. (26) - spatial filtering factor for pair collection - K1: Smirr et al., Eq. (39) - spatial filtering factor for single photon - Γ2|1: Smirr et al., Eq. (41) - heralding ratio for single photon sources

The heralding ratio Γ2|1 is a key figure of merit for single photon sources, as it determines the quality of heralded single photons. Values closer to 1 indicate better performance, with typical experimental values of 0.5-0.9.

For entangled pair sources, both K1 and K2 should be maximized to optimize the coincidence detection rate while maintaining good spatial mode purity.

citrine.calculate_optimisation_grid

calculate_optimisation_grid(xi: float, alpha_range: float | NDArray[floating], phi0_range: float | NDArray[floating], n_p: float, n_s: float, n_i: float, rho_max: float = 2.0, progress_bar: bool = True)

Calculate optimization grids for SPDC source design across parameter space.

This function systematically evaluates the brightness factors K1, K2 and heralding ratio across ranges of the normalized target mode waist (α) and longitudinal phase mismatch (φ0) parameters. This enables optimization of SPDC sources as described in Section III.C of Smirr et al.

The results can be used to generate contour plots like Figures 3, 4, and 6 in the paper, showing the optimal parameter combinations for different focusing conditions.

Parameters:

  • xi (float) –

    Focusing parameter ξ = L/(2zR). Different values correspond to different pump beam focusing conditions (see Fig. 4 in paper).

  • alpha_range (Union[float, NDArray[floating]]) –

    Range of normalized target mode waist values α = a0/w0 to evaluate. Typical range: 0.5 to 3.5 (see Fig. 3 in paper).

  • phi0_range (Union[float, NDArray[floating]]) –

    Range of longitudinal phase mismatch values φ0 = Δk0·L to evaluate. This parameter is typically optimized via crystal temperature control.

  • n_p (float) –

    Pump refractive index.

  • n_s (float) –

    Signal refractive index.

  • n_i (float) –

    Idler refractive index.

  • rho_max (float, default: 2.0 ) –

    Maximum integration bound for radial coordinates. Defaults to 2.0.

  • progress_bar (bool, default: True ) –

    Whether to display progress bar during calculation. Defaults to True.

Returns:

  • Tuple[NDArray, NDArray, NDArray]: A tuple containing: - K1_grid (NDArray): 2D array with shape (n_phi, n_alpha) of single photon coupling efficiency values for heralding calculations. - K2_grid (NDArray): 2D array with shape (n_phi, n_alpha) of photon pair coupling efficiency values for brightness optimization (Fig. 3). - ratio_grid (NDArray): 2D array with shape (n_phi, n_alpha) of heralding ratio Γ2|1 values for heralding optimization (Fig. 6).

Note

This function enables the systematic optimization described in Section III.C: 1. For brightness optimization: find (α, φ0) that maximizes K2_grid 2. For heralding optimization: find (α, φ0) that maximizes ratio_grid 3. For balanced optimization: find compromise between brightness and heralding

The optimal parameters typically follow the Boyd-Kleinman conditions for second harmonic generation, adapted for the down-conversion geometry.

citrine.make_lookup_table

make_lookup_table(xi_values: float | NDArray[floating], alpha: float | NDArray[floating], phi_0: float | NDArray[floating], crystal: Crystal, wavelength_pump: Wavelength, wavelength_signal: Wavelength, wavelength_idler: Wavelength, name: str | None = None, temp_dir: str | None = None, rho_max: float = 2.0, temperature: float = None, overwrite: bool = False) -> dict

Generate lookup tables for SPDC optimization across focusing parameter range.

This function creates comprehensive lookup tables that enable fast interpolation of brightness and heralding calculations across different focusing conditions. This is essential for the systematic optimization described in Section III of Smirr et al., where multiple ξ values must be evaluated to find global optima.

Parameters:

  • xi_values (Union[float, NDArray[floating]]) –

    Range of focusing parameter values ξ = L/(2zR) to calculate. Typical range: 0.1 to 10 (from collimated to tightly focused).

  • alpha (Union[float, NDArray[floating]]) –

    Range of normalized target mode waist values α = a0/w0.

  • phi_0 (Union[float, NDArray[floating]]) –

    Range of longitudinal phase mismatch values φ0 = Δk0·L.

  • crystal (Crystal) –

    Crystal object containing material properties and Sellmeier equations.

  • wavelength_pump (Wavelength) –

    Wavelength object for pump beam.

  • wavelength_signal (Wavelength) –

    Wavelength object for signal beam.

  • wavelength_idler (Wavelength) –

    Wavelength object for idler beam.

  • name (Optional[str], default: None ) –

    Output filename for the lookup table pickle file. Defaults to None (auto-generated).

  • temp_dir (Optional[str], default: None ) –

    Directory for temporary storage of intermediate results. Defaults to None.

  • rho_max (float, default: 2.0 ) –

    Maximum radial integration bound. Defaults to 2.0.

  • temperature (float, default: None ) –

    Crystal temperature for refractive index calculations. Defaults to None.

  • overwrite (bool, default: False ) –

    Whether to overwrite existing intermediate files. Defaults to False.

Returns:

  • dict

    Tuple[BrightnessHeraldingLUT, str]: A tuple containing: - lookup_table (BrightnessHeraldingLUT): Lookup table object containing all calculated grids and parameters. - filename (str): Path to the saved lookup table file.

Note

The lookup table enables: 1. Fast parameter optimization without recalculating integrals 2. Interpolation between calculated points for smooth optimization 3. Systematic comparison across different focusing regimes

The lookup table structure enables the optimization workflow from Fig. 4: 1. Calculate grids for multiple ξ values 2. Find optimal (α, φ0) for each ξ 3. Determine global optimum across all ξ values 4. Use interpolation for fine-tuning around optima

citrine.BrightnessHeraldingLUT dataclass

BrightnessHeraldingLUT(refractive_index: dict[str, float], wavelength: dict[str, Wavelength], temperature: float, xi: float | NDArray[floating], alpha: float | NDArray[floating], phi: float | NDArray[floating], K1: NDArray[floating], K2: NDArray[floating], Heralding: NDArray[floating])

Lookup table for SPDC brightness and heralding optimization.

This class encapsulates the complete parameter space calculations needed for SPDC source optimization as described in Smirr et al. It provides methods for interpolating between calculated points and extracting optimal parameters.

The lookup table enables rapid optimization across the multidimensional parameter space (ξ, α, φ0) without recalculating expensive integrals.

Attributes:

  • refractive_index (Dict[str, float]) –

    Refractive indices for pump, signal, and idler wavelengths.

  • wavelength (Dict[str, Wavelength]) –

    Wavelength specifications for the three interacting beams.

  • temperature (float) –

    Crystal temperature used for calculations.

  • xi (Union[float, NDArray[floating]]) –

    Focusing parameter values ξ = L/(2zR).

  • alpha (Union[float, NDArray[floating]]) –

    Normalized target mode waist values α = a0/w0.

  • phi (Union[float, NDArray[floating]]) –

    Longitudinal phase mismatch values φ0 = Δk0·L.

  • K1 (NDArray[floating]) –

    Single photon coupling efficiency grids with shape (n_xi, n_phi, n_alpha).

  • K2 (NDArray[floating]) –

    Photon pair coupling efficiency grids with shape (n_xi, n_phi, n_alpha).

  • Heralding (NDArray[floating]) –

    Heralding ratio grids Γ2|1 = K2/K1 with shape (n_xi, n_phi, n_alpha).

Methods:

  • from_file

    Load lookup table from pickled file.

  • interpolate

    Create interpolated grids for specified focusing parameter values.

Functions

from_file classmethod

from_file(path: str)

Load lookup table from pickled file.

Parameters:

  • path (str) –

    Path to the pickled lookup table file.

Returns:

  • BrightnessHeraldingLUT

    Loaded lookup table object.

interpolate

interpolate(target_xi_values: float | list[float] | NDArray[floating]) -> tuple[NDArray[np.floating], NDArray[np.floating], NDArray[np.floating]]

Create interpolated grids for specified focusing parameter values.

This method enables fine-grained optimization by interpolating between the calculated lookup table points. This is particularly useful for finding precise optimal parameters around the global maxima identified from the coarse grid calculations.

Parameters:

  • target_xi_values (Union[float, List[float], NDArray[floating]]) –

    Focusing parameter values for which to create interpolated grids. These can be between the original calculated ξ values.

Returns:

  • tuple[NDArray[floating], NDArray[floating], NDArray[floating]]

    Tuple[NDArray[np.floating], NDArray[np.floating], NDArray[np.floating]]: A tuple containing: - interpolated_grids (dict): Dictionary with ξ values as keys, each containing 'K2_grid' and 'heralding_grid' for interpolated brightness and heralding optimization grids. - alpha_range (NDArray): Alpha values used in the grids (from original lookup table). - phi0_range (NDArray): Phi0 values used in the grids (from original lookup table).

Note

The interpolation uses RegularGridInterpolator for smooth, continuous parameter optimization. This enables finding optimal parameters with higher precision than the original grid spacing.

Warning will be issued if target ξ values are outside the lookup table range.