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)

Class to represent angular frequency.

Attributes:

  • value (Union[float, ndarray]) –

    The angular frequency.

citrine.Wavelength dataclass

Wavelength(value, unit=...)

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

Convert the wavelength to a wavevector.

Returns:

  • float ( float ) –

    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.SellmeierCoefficients dataclass

SellmeierCoefficients(first_order, second_order, temperature, zeroth_order=...)

Sellmeier coefficients for calculating refractive indices.

Attributes:

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

    Zeroth-order Sellmeier coefficients.

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

    First-order Sellmeier coefficients.

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

    Second-order Sellmeier coefficients.

  • temperature (float) –

    The reference temperature for the coefficients (in Celsius).

citrine._permittivity

_permittivity(sellmeier: Union[List[float], NDArray], wavelength_um: Union[float, NDArray[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.

Parameters:

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

    Sellmeier coefficients.

  • wavelength_um (float) –

    Wavelength in micrometers.

Returns:

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

    The permittivity value.

Source code in src/citrine/citrine.py
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
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: Union[List[float], NDArray], wavelength_um: Union[float, NDArray[floating]]) -> Union[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 ( Union[float, NDArray[floating]] ) –

    The refractive index (n_0).

Source code in src/citrine/citrine.py
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
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: Union[List[float], NDArray[floating]], wavelength_um: Union[float, NDArray[floating]]) -> Union[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 ( Union[float, NDArray[floating]] ) –

    The refractive index (n_i).

Source code in src/citrine/citrine.py
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
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.refractive_index

refractive_index(sellmeier: SellmeierCoefficients, wavelength: Wavelength, temperature: float | None = None) -> float

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

Parameters:

  • sellmeier (SellmeierCoefficients) –

    The Sellmeier coefficients for the material.

  • wavelength (Wavelength) –

    The wavelength at which to calculate the refractive index.

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

    The temperature (in Celsius), defaults to the reference temperature.

Returns:

  • float ( float ) –

    The refractive index at the given wavelength and temperature.

Source code in src/citrine/citrine.py
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
def refractive_index(
    sellmeier: SellmeierCoefficients,
    wavelength: Wavelength,
    temperature: Optional[float] = None,
) -> Union[float, NDArray[np.floating]]:
    """
    Calculate the refractive index for a given wavelength and temperature using the Sellmeier equation.

    Args:
        sellmeier (SellmeierCoefficients): The Sellmeier coefficients for the material.
        wavelength (Wavelength): The wavelength at which to calculate the refractive index.
        temperature (Optional[float]): The temperature (in Celsius), defaults to the reference temperature.

    Returns:
        float: The refractive index at the given wavelength and temperature.
    """
    if temperature is None:
        temperature = sellmeier.temperature

    wavelength_um = wavelength.to_unit(Magnitude.micro).value
    n0: Union[float, NDArray[np.floating]] = 0.0
    if sellmeier.zeroth_order is not None:
        n0 = _n_0(sellmeier.zeroth_order, wavelength_um)

    n1 = None
    if sellmeier.first_order is not None:
        n1 = _n_i(sellmeier.first_order, wavelength_um)

    n2 = None
    if sellmeier.second_order is not None:
        n2 = _n_i(sellmeier.second_order, wavelength_um)

    t_offset = temperature - sellmeier.temperature

    if n1 is None:
        n1 = 0.0

    if n2 is None:
        n2 = 0.0

    n = n0 + (n1 * t_offset) + (n2 * t_offset**2)
    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, sellmeier_o, sellmeier_e, phase_matching: PhaseMatchingCondition = PhaseMatchingCondition.type0_e, doi: str = None)

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

Attributes:

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

Attributes:

Methods:

  • refractive_index

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

  • refractive_indices

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

Attributes

phase_matching property writable

phase_matching: PhaseMatchingCondition

Get the current phase matching condition.

Returns:

Functions

refractive_index

refractive_index(wavelength: Wavelength, photon: Literal['pump', 'signal', 'idler'], temperature: float | None = None) -> float

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

    Refractive index.

refractive_indices

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

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, float, float]

    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

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

    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) -> np.ndarray

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:

  • ndarray

    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: Union[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. SVD decomposition of the JSA 2. Computes quantum metrics (purity, Schmidt number, entropy)

Phase Matching Functions

citrine.phase_matching.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.phase_matching.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.pump_envelope.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.pump_envelope.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.pump_envelope.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.crystals.KTiOPO4_Fradkin module-attribute

KTiOPO4_Fradkin = Crystal(name='Potassium Titanyl Phosphate', doi='10.1063/1.123408', sellmeier_o=SellmeierCoefficients(zeroth_order=[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=SellmeierCoefficients(zeroth_order=[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.crystals.KTiOPO4_Emanueli module-attribute

KTiOPO4_Emanueli = Crystal(name='Potassium Titanyl Phosphate', doi='10.1364/AO.42.006661', sellmeier_o=SellmeierCoefficients(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=SellmeierCoefficients(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.crystals.KTiOAsO4_Emanueli module-attribute

KTiOAsO4_Emanueli = Crystal(name='Potassium Titanyl Aresnate', doi='10.1364/AO.42.006661', sellmeier_o=SellmeierCoefficients(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=SellmeierCoefficients(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.crystals.LiNbO3_Zelmon module-attribute

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

citrine.crystals.LiNbO3_5molMgO_Zelmon module-attribute

LiNbO3_5molMgO_Zelmon = Crystal(name='Lithium Niobate 5Mol Magnesium Doped', doi='10.1364/JOSAB.14.003319', sellmeier_o=SellmeierCoefficients(zeroth_order=None, first_order=array([2.2454, 0.01242, 1.3005, 0.0513, 6.8972, 331.33]), second_order=None, temperature=21), sellmeier_e=SellmeierCoefficients(zeroth_order=None, first_order=array([2.4272, 0.01478, 1.4617, 0.005612, 9.6536, 371.216]), second_order=None, temperature=21), phase_matching=type0_e)