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

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,
) -> 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: float,
) -> float

Compute the permittivity using the Sellmeier equation.

Parameters:

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

    Sellmeier coefficients.

  • wavelength_um (float) –

    Wavelength in micrometers.

Returns:

  • float ( float ) –

    The permittivity value.

Source code in citrine/citrine.py
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
def _permittivity(
    sellmeier: Union[List[float], NDArray], wavelength_um: float
) -> float:
    """
    Compute the permittivity using the Sellmeier equation.

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

    Returns:
        float: The permittivity value.
    """
    wl = wavelength_um**2
    # a = 1
    # b = a + 1
    # p = sellmeier[0]
    # lim: int = int(np.ceil((len(sellmeier) - 2) / 2))
    # for i in range(0, lim):
    #     p += sellmeier[a] / (1 - sellmeier[b] / wl)
    #     a += 2
    #     b += 2

    # p -= sellmeier[-1] * wl
    # return p

    first, *coeffs, last = sellmeier
    p = first
    assert len(coeffs) % 2 == 0
    for numer, denom in zip(coeffs[0::2], coeffs[1::2]):
        p += numer / (1 - denom / wl)
    p -= last * wl
    return p

citrine._n_0

_n_0(
    sellmeier: Union[List[float], NDArray],
    wavelength_um: float,
) -> float

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

    The refractive index (n_0).

Source code in citrine/citrine.py
188
189
190
191
192
193
194
195
196
197
198
199
def _n_0(sellmeier: Union[List[float], NDArray], wavelength_um: float) -> float:
    """
    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],
    wavelength_um: float,
) -> float

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

    The refractive index (n_i).

Source code in citrine/citrine.py
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
def _n_i(sellmeier: Union[List[float], NDArray], wavelength_um: float) -> float:
    """
    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 citrine/citrine.py
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
def refractive_index(
    sellmeier: SellmeierCoefficients,
    wavelength: Wavelength,
    temperature: Optional[float] = None,
) -> float:
    """
    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 = 0
    if sellmeier.zeroth_order is not None:
        n0 = _n_0(sellmeier.zeroth_order, wavelength_um)
    n1 = _n_i(sellmeier.first_order, wavelength_um)
    n2 = _n_i(sellmeier.second_order, wavelength_um)

    t_offset = temperature - sellmeier.temperature

    n = n0 + (n1 * t_offset) + (n2 * t_offset**2)
    return n

citrine.Orientation

Bases: Enum

Enum to represent the orientation: ordinary or extraordinary.

citrine.Crystal dataclass

Crystal(
    name,
    sellmeier_o,
    sellmeier_e,
    pump_orientation,
    signal_orientation,
    idler_orientation,
)

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

Attributes:

Functions

refractive_index

refractive_index(
    wavelength: Wavelength,
    polarization: Literal["ordinary", "extraordinary"],
    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.

  • polarization (Literal['ordinary', 'extraordinary']) –

    Polarization ('ordinary' or 'extraordinary').

  • photon (Literal['pump', 'signal', 'idler']) –

    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,
) -> 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,
) -> 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_matching_function

phase_matching_function(
    delta_k: np.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.pump_envelope_gaussian

pump_envelope_gaussian(
    lambda_p: Wavelength,
    sigma_p: float,
    lambda_s: Wavelength,
    lambda_i: Wavelength,
) -> np.ndarray

Generate a pump envelope matrix with a Gaussian profile.

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.

Returns:

  • ndarray

    np.ndarray: A 2D pump envelope matrix.

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.

Crystal Library

citrine.crystals.sellmeier_o_ppKTP module-attribute

sellmeier_o_ppKTP: Incomplete = 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,
)

citrine.crystals.sellmeier_e_ppKTP module-attribute

sellmeier_e_ppKTP: Incomplete = 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,
)

citrine.crystals.ppKTP module-attribute

ppKTP: Incomplete = Crystal(
    name="ppKTP",
    sellmeier_o=sellmeier_o_ppKTP,
    sellmeier_e=sellmeier_e_ppKTP,
    pump_orientation=extraordinary,
    signal_orientation=extraordinary,
    idler_orientation=ordinary,
)