Refractive Index
NonlinearCrystals.SellmeierFunction — TypeSellmeierFunction(n_fun, lambda_range=nothing; temp_ref=293.15u"K", temp_range=nothing)Wraps a wavelength- and temperature-dependent refractive index function in a RefractiveIndex model. This is typically used to represent Sellmeier equations or any function that evaluates n(λ, T), possibly dependent on a reference temperature temp_ref. If lambda_range and temp_range are provided, they define the validity range and are checked optionally during evaluation.
NonlinearCrystals.default_lambda — Methoddefault_lambda(ri::RefractiveIndex)Returns a default wavelength for use with the refractive index model ri. If ri.lambda_range is defined, it returns the midpoint of the range; otherwise, it defaults to 633 nm.
NonlinearCrystals.default_temp — Methoddefault_temp(ri::RefractiveIndex)Returns the reference temperature temp_ref defined in the refractive index model ri. This is typically used as the default when no explicit temperature is provided.
NonlinearCrystals.dn_dtemp — Functiondn_dtemp(ri::RefractiveIndex, lambda::Length=default_lambda(ri), temp::Temperature=default_temp(ri))Computes the temperature derivative of the RefractiveIndex at the specified wavelength lambda and temperature temp. This is used to analyze thermal drift, tuning, or thermal lensing effects in nonlinear and resonant systems.
NonlinearCrystals.group_index — Functiongroup_index(ri, lambda=default_lambda(ri), temp=default_temp(ri))Computes the group index n_g, which determines the group velocity in a dispersive medium. This is calculated as β₁ · c, where β₁ is the first derivative of the propagation constant with respect to angular frequency. The result is unitless.
NonlinearCrystals.group_velocity — Functiongroup_velocity(ri::RefractiveIndex, lambda::Length=default_lambda(ri), temp::Temperature=default_temp(ri))Computes the group velocity of a wavepacket centered at lambda, using the first derivative of the wavevector with respect to frequency. This quantity reflects how the envelope of a pulse propagates through the medium.
NonlinearCrystals.is_lambda_valid — Methodis_lambda_valid(lambda::Length, sri::SellmeierFunction; warn_tol::Length=1u"nm")Checks whether the given wavelength is within the valid lambda_range of a SellmeierFunction. A tolerance warn_tol can be specified to allow for small numerical deviations near the boundaries. Returns true if the wavelength lies within the specified valid range, false otherwise.
NonlinearCrystals.phase_velocity — Functionphase_velocity(ri::RefractiveIndex, lambda::Length=default_lambda(ri), temp::Temperature=default_temp(ri))Returns the phase velocity v_p = c / n(λ) of a monochromatic wave in the material described by ri at the given conditions. The result has physical units of velocity.
NonlinearCrystals.plot_refractiveindex! — Methodplot_refractiveindex!(ri::RefractiveIndex; n_sample_pts=500, temp=nothing, lambda_min=nothing,
lambda_max=nothing, label="", colormap=:vik, digits::Integer=3)Adds a refractive index vs. wavelength curve to the current GLMakie axis, optionally for multiple temperatures. This is useful for overlaying plots from different models or conditions. The interactive inspector shows wavelength, index, and temperature.
NonlinearCrystals.plot_refractiveindex — Methodplot_refractiveindex(ri::RefractiveIndex; n_sample_pts=500, temp=nothing, lambda_min=nothing, lambda_max=nothing, label="", colormap=:vik)Creates and returns a new GLMakie plot of the refractive index as a function of wavelength, optionally at multiple temperatures if a vector of temperatures is provided. If temp is omitted, the reference temperature of ri is used.
NonlinearCrystals.refractive_index — Functionrefractive_index(sri::SellmeierFunction, lambda::Length, temp::Temperature=sri.temp_ref;
check_lambda_range::Symbol=:warn, check_temp_range::Symbol=:warn,
warn_tol::Length=1u"nm")Evaluates the refractive index for a given wavelength and temperature using the stored function inside a SellmeierFunction. If bounds are set and check_lambda_range or check_temp_range is not :none, warnings or errors are issued for out-of-range inputs. Returns a unitless number corresponding to the refractive index.
NonlinearCrystals.BidirectionalCrystal — TypeBidirectionalCrystal(metadata::Dict,
n_X_principal::RefractiveIndex,
n_Y_principal::RefractiveIndex,
n_Z_principal::RefractiveIndex,
d_XYZ_ref_full::AbstractArray{<:Number,3};
miller_delta=nothing,
warn_n_Z_smaller_n_X=true)Constructs a biaxial crystal using separate refractive index models for the principal X, Y, and Z axes, along with a full nonlinear tensor d_XYZ_ref_full. The metadata must include a valid :point_group. If warn_n_Z_smaller_n_X is true, a warning is issued if the crystal axes are not sorted as n_Z ≥ n_Y ≥ n_X, which is the expected convention in this package. If miller_delta is provided, it is used for Miller scaling of the effective nonlinear coefficient d_eff.
NonlinearCrystals.RefractionData — TypeRefractionData(hi_or_lo::Symbol, theta::Angle, phi::Angle, cr::NonlinearCrystal, lambda::Length=default_lambda(cr);
temp::Temperature=default_temp(cr))Convenience constructor that combines computation and extraction: builds a RefractionDataHiLo and returns the single-polarization RefractionData corresponding to :hi or :lo.
NonlinearCrystals.RefractionData — TypeRefractionDataA single-branch slice of RefractionDataHiLo corresponding to either the :hi or :lo polarization branch. Stores all relevant scalar and vector optical quantities for one mode: refractive index, group index, Poynting vector, polarization directions, walk-off, and β dispersion parameters.
NonlinearCrystals.RefractionDataHiLo — TypeRefractionDataHiLo(theta::Angle, phi::Angle, cr::NonlinearCrystal, lambda::Length=default_lambda(cr);
temp::Temperature=default_temp(cr), angle_tol::Angle=0.1u"°")Constructs a high/low index birefringent refraction object RefractionDataHiLo for a given propagation direction defined by spherical angles theta, phi. Internally calculates refractive indices, energy and phase directions, walk-off angles, group indices, and frequency derivatives for both polarization branches. Polarization classification is automatic based on principal planes or optical axis orientation.
NonlinearCrystals.RefractionDataHiLo — TypeRefractionDataHiLoContains all relevant optical data for the two birefringent branches (:hi and :lo) in a nonlinear crystal. Includes refractive indices, group indices, energy flow vectors, walk-off angles, and dispersion derivatives (β₀ through β₃). Computed for a specific propagation direction, wavelength, and temperature. RefractionDataHiLo can be split into RefractionData instances to represent only one of both polarization branches.
NonlinearCrystals.RefractionType — TypeRefractionTypeDescribes a single polarization state (:o or :e) in a specified principal plane.
NonlinearCrystals.RefractionType — MethodRefractionType(hi_or_lo::Symbol, rt::RefractionTypeHiLo)Extracts a single-polarization RefractionType from a two-branch RefractionTypeHiLo, corresponding to either the high or low index solution.
NonlinearCrystals.RefractionTypeHiLo — TypeRefractionTypeHiLoDescribes the polarization type of the two birefringent solutions (:hi and :lo) in a given principal plane. Each polarization is labeled as ordinary (:o) or extraordinary (:e) based on field orientation.
NonlinearCrystals.RefractionTypeHiLo — MethodRefractionTypeHiLo(principal_plane::Symbol, E_dir_hi_lo::NTuple{2,<:AbstractVector{<:Number}})Classifies a pair of electric field directions as ordinary or extraordinary relative to the specified principal plane. Used to track the polarization type of high and low index branches in birefringent crystals.
NonlinearCrystals.UnidirectionalCrystal — TypeUnidirectionalCrystal(metadata::Dict, n_o_principal::RefractiveIndex,
n_e_principal::RefractiveIndex,
d_XYZ_ref_full::AbstractArray{<:Number,3};
miller_delta=nothing)Constructs a uniaxial crystal with ordinary (n_o_principal) and extraordinary (n_e_principal) refractive index models, and a full 3×3×3 nonlinear tensor d_XYZ_ref_full. The metadata must include a recognized :point_group keyword. If miller_delta is provided, it is used for Miller scaling of the effective nonlinear coefficient d_eff. The Z axis is treated as the optical axis.
NonlinearCrystals.assign_o_or_e — Methodassign_o_or_e(principal_plane::Symbol, E_dir::AbstractVector{<:Number}; angle_tol_ud::Angle=0.2u"°")Classifies a polarization vector E_dir as ordinary (:o) or extraordinary (:e) with respect to the given principal_plane. In uniaxial crystals (:UD), the classification depends on the angle between E_dir and the optical axis. In biaxial crystals, the function assumes E_dir lies in or perpendicular to the specified principal plane (:XY, :XZ, or :YZ), and emits an assertion otherwise. Returns the symbol :o or :e depending on the geometric configuration.
NonlinearCrystals.calc_E_dir_S_dir — Methodcalc_E_dir_S_dir(k_dir::AbstractVector{<:Real},
ε_tensor::AbstractMatrix{<:Real},
D_dir_hi_lo::Tuple)From the D vectors and dielectric tensor, computes the electric field directions E = ε⁻¹·D and corresponding Poynting vectors S = E × H, where H = k × E. The resulting vectors are normalized and aligned with k_dir.
NonlinearCrystals.calc_k_dir_ε_tensor_n_XYZ — Functioncalc_k_dir_ε_tensor_n_XYZ(θ::Angle, ϕ::Angle, cr::NonlinearCrystal,
lambda::Length=default_lambda(cr),
temp::Temperature=default_temp(cr))Computes the wavevector direction from spherical angles and constructs the dielectric tensor ε for the given crystal cr and conditions. Also returns the principal refractive indices along the X, Y, and Z axes for the given wavelength lambda and temperature temp.
NonlinearCrystals.calc_n_hi_lo — Methodcalc_n_hi_lo(θ::Angle, ϕ::Angle, cr::NonlinearCrystal, lambda::Length; temp::Temperature)Computes the refractive indices of the two birefringent eigenmodes for a given direction in spherical angles. This is a minimal calculation used for derivative tracing (e.g., with ForwardDiff) and does not return field vectors.
NonlinearCrystals.calc_n_hi_lo_D_dir_hi_lo — Methodcalc_n_hi_lo_D_dir_hi_lo(k_dir::AbstractVector{<:Real}, ε_tensor::AbstractMatrix{<:Real})Given a propagation direction k_dir and dielectric tensor ε_tensor, computes the two effective refractive indices and corresponding displacement vectors (D) for the birefringent modes. Works by diagonalizing the projected dielectric tensor in the plane orthogonal to k_dir.
NonlinearCrystals.default_lambda — Methoddefault_lambda(cr::NonlinearCrystal)Returns a default wavelength for use with a crystal cr. If a wavelength range is defined for the X-axis principal index, the midpoint of that range is returned. Otherwise, the default value is 633 nm.
NonlinearCrystals.default_temp — Methoddefault_temp(cr::NonlinearCrystal)Returns the reference temperature from the X-axis principal refractive index of the crystal cr. This is typically used when no explicit temperature is provided.
NonlinearCrystals.hi_lo_to_o_e — Functionhi_lo_to_o_e(hi_or_lo::Symbol, cr::UnidirectionalCrystal, lambda::Length=default_lambda(cr); temp::Temperature=default_temp(cr))Maps a :hi or :lo refractive index label to the corresponding polarization :o or :e, based on the actual indices at the given wavelength and temperature. This is the inverse of o_e_to_hi_lo and assumes uniaxial behavior.
NonlinearCrystals.hi_lo_to_o_e — Methodhi_lo_to_o_e(hi_or_lo_rrb::NTuple{3,Symbol}, cr::UnidirectionalCrystal, lambda_rrb::NTuple{3,Length}; temp::Temperature=default_temp(cr))Applies hi_lo_to_o_e elementwise to a 3-tuple of :hi/:lo labels and corresponding wavelengths, for use in triple-wave processes. Returns a tuple of :o/:e polarization labels.
NonlinearCrystals.is_lambda_valid — Methodis_lambda_valid(lambda::Length, cr::NonlinearCrystal; warn_tol::Length=1u"nm")Checks whether the given wavelength lies within the valid wavelength range of all three principal refractive index models of the crystal cr. A tolerance warn_tol allows small deviations near the boundaries. Returns true if the wavelength is valid for all three axes, otherwise false.
NonlinearCrystals.o_e_to_hi_lo — Functiono_e_to_hi_lo(o_or_e::Symbol, cr::UnidirectionalCrystal, lambda::Length=default_lambda(cr); temp::Temperature=default_temp(cr))Maps a polarization label :o or :e to a high (:hi) or low (:lo) refractive index, based on the comparison between ordinary and extraordinary indices at the given wavelength and temperature. Returns :hi if the refractive index for the given polarization is larger than the other, and :lo otherwise.
NonlinearCrystals.o_e_to_hi_lo — Methodo_e_to_hi_lo(o_or_e_rrb::NTuple{3,Symbol}, cr::UnidirectionalCrystal, lambda_rrb::NTuple{3,Length}; temp::Temperature=default_temp(cr))Applies o_e_to_hi_lo elementwise to a 3-tuple of polarizations o_e_to_hi_lo and a 3-tuple of wavelengths lambda_rrb, corresponding to the first (typically the longest wavelength) and second red wavelength and the blue wavelength. Returns a 3-tuple of :hi/:lo symbols.
NonlinearCrystals.optical_axis_angle — Functionoptical_axis_angle(cr::UnidirectionalCrystal, lambda::Length=default_lambda(cr), temp::Temperature=default_temp(cr))Returns the optical axis angle for a uniaxial crystal, which is always zero by definition since the optical axis aligns with the principal Z-axis.
NonlinearCrystals.optical_axis_angle — Functionoptical_axis_angle(cr::BidirectionalCrystal, lambda::Length=default_lambda(cr), temp::Temperature=default_temp(cr))Computes the angle between the optical axis and the principal Z-axis for a biaxial crystal cr, based on the principal refractive indices at the specified wavelength and temperature. Returns an angle between 0 and 90°.
NonlinearCrystals.plot_birefringent_refraction — Functionplot_birefringent_refraction(theta, phi, cr::NonlinearCrystal,
lambda::Length=default_lambda(cr),
temp::Temperature=default_temp(cr))Creates a 3D GLMakie visualization of birefringent refraction for the given direction and crystal. The plot shows the index ellipsoid, wavevector direction k, and the Poynting vector (S), electric fields (E), and displacement (D) vectors for both high and low index solutions. Useful for verifying polarization, walk-off behavior, and eigenvector orientation.
NonlinearCrystals.plot_refractiveindex — Methodplot_refractiveindex(cr::BidirectionalCrystal; n_sample_pts=500,
temp=[default_temp(cr)], lambda_min=nothing, lambda_max=nothing)Plots the refractive indices n_X, n_Y, and n_Z as functions of wavelength for a biaxial crystal using GLMakie. Each principal axis is drawn with a different color and can optionally show temperature dependence. Returns a figure with legend and interactive inspection enabled.
NonlinearCrystals.plot_refractiveindex — Methodplot_refractiveindex(cr::UnidirectionalCrystal; n_sample_pts=500,
temp=[default_temp(cr)], lambda_min=nothing, lambda_max=nothing)Plots the ordinary (n_X = n_Y) and extraordinary (n_Z) refractive indices of a uniaxial crystal over wavelength using GLMakie. The ordinary index is shown once for both X and Y, and different colormaps distinguish the two polarizations. Returns a Makie figure with labeled axes, legend, and interactive inspection.
NonlinearCrystals.polarization_rrb_to_hi_lo — Methodpolarization_rrb_to_hi_lo(pol_rrb::NTuple{3,Symbol}, cr::NonlinearCrystal, lambda_rrb::NTuple{3,Length}; temp::Temperature=default_temp(cr))Converts a 3-tuple of polarization symbols (:o/:e or :hi/:lo) to :hi/:lo, depending on the type of crystal and the refractive indices at the given wavelengths. Only uniaxial crystals support conversion from :o/:e; other crystal types require the input to already use :hi/:lo. Returns a 3-tuple of :hi/:lo symbols or emits an error if the input is invalid.
NonlinearCrystals.valid_lambda_range — Methodvalid_lambda_range(cr::NonlinearCrystal)Returns the range of wavelengths that is valid across all three principal axes of the crystal cr. The result is the intersection of the individual wavelength ranges for n_X, n_Y, and n_Z.
NonlinearCrystals.crystal_system — Methodcrystal_system(point_group::AbstractString)Returns the name of the crystal system corresponding to the given crystallographic point group (e.g., "mm2" → "Orthorhombic"). Throws an error if the point group is unknown. Based on conventional assignments from nonlinear optics literature.
NonlinearCrystals.plot_symmetry — Methodplot_symmetry(point_group::AbstractString)Visualizes the symmetry relations between nonlinear tensor components d_ij for a given point group.
Black lines indicate symmetry relations that follow from the basic point group; dashed lines show those that appear only under Kleinman symmetry assumptions. Black circles represent same-sign components; white circles represent opposite signs. Square markers indicate coefficients that become zero under Kleinman conditions.
NonlinearCrystals.calc_d_eff — Methodcalc_d_eff(cr::NonlinearCrystal,
E_dir_r1, E_dir_r2, E_dir_b;
lambda_rrb=nothing, temp=default_temp(cr), use_miller_scaling=true)Computes the effective nonlinear coefficient d_eff for a given nonlinear crystal and polarization directions E_dir_r1, E_dir_r2, E_dir_b.
If lambda_rrb (a 3-tuple of wavelengths) is provided and Miller data is available, the tensor is scaled using Miller's rule. Otherwise, the reference tensor d_XYZ_ref_full is used directly without Miller scaling.
NonlinearCrystals.calc_miller_delta — Methodcalc_miller_delta(d_ref_XYZ_full, n_X_principal, n_Y_principal, n_Z_principal, temp_ref;
lambda_r1=nothing, lambda_r2=nothing, lambda_b=nothing)Computes the Miller scaling tensor Δᵢⱼₖ from a known nonlinear tensor d_ref_XYZ_full and refractive index models.
Each χ⁽²⁾ component is scaled by the product of χ⁽¹⁾ values (defined as n² − 1) for the corresponding directions and wavelengths. You may also call this function with n_o_principal and n_e_principal for uniaxial crystals.
NonlinearCrystals.contract_d_tensor — Methodcontract_d_tensor(d_full::AbstractArray{<:Number,3})Contracts a 3×3×3 nonlinear tensor into a 3×6 Voigt representation, assuming symmetry in the last two indices. This is the inverse of expand_d_contract.
NonlinearCrystals.contract_voigt_index — Methodcontract_voigt_index(i::Integer, j::Integer, k::Integer) -> (i::Int, l::Int)Given full tensor indices (i, j, k), return the contracted Voigt index (i, l). This assumes symmetry in (j, k): for example, (1,3) and (3,1) both map to l = 5.
NonlinearCrystals.expand_d_contract — Methodexpand_d_contract(d_contract::AbstractMatrix{<:Number})Expands a 3×6 Voigt-contracted nonlinear tensor into a full 3×3×3 tensor with explicit symmetry in the second and third indices.
NonlinearCrystals.expand_voigt_index — Methodexpand_voigt_index(i::Integer, l::Integer)Expands a contracted Voigt index pair (i, l) into the corresponding set of full tensor indices (i, j, k). If l refers to a symmetric index (like l = 5 for j,k = 1,3 and 3,1), both permutations are returned.
NonlinearCrystals.miller_rescale — Methodmiller_rescale(cr::NonlinearCrystal, lambda_rrb; temp=default_temp(cr))Applies Miller scaling to the reference nonlinear tensor stored in cr, using its miller_delta data (if available).
The result is a scaled tensor d_XYZ_full for the given wavelength triplet lambda_rrb = (λ₁, λ₂, λ₃). This models the variation of nonlinear response with wavelength using Miller's empirical rule.
NonlinearCrystals.plot_axes_assignment_crys_to_diel — Methodplot_axes_assignment_crys_to_diel(axes_assignment_crys_to_diel::NTuple{3,Symbol};
rotate_about::Union{Symbol,Nothing}=nothing,
phi::Angle=0.0u"°")Visualizes the transformation from crystallophysical axes (x, y, z) to dielectric axes (X, Y, Z) in 3D using GLMakie. The original axes are shown as dashed blue vectors, and the rotated dielectric frame is shown as solid red vectors. Used to verify the effect of rot_mat_crys_to_diel.
NonlinearCrystals.plot_miller_scaling_coeffs_shg — Methodplot_miller_scaling_coeffs_shg(cr::NonlinearCrystal;
temp=default_temp(cr),
size=(800, 600))Visualizes the Miller-scaled nonlinear tensor coefficients d_ij(λ) for second harmonic generation (SHG) as a function of blue wavelength λ_b. Assumes type-0 SHG, where the red wavelengths are 2λ_b. All nonzero contracted components are plotted with interactive tooltips showing wavelength and coefficient values.
NonlinearCrystals.rot_mat_crys_to_diel — Methodrot_mat_crys_to_diel(axes_assignment_crys_to_diel::NTuple{3,Symbol};
rotate_about::Union{Symbol,Nothing}=nothing,
phi::Angle=0.0u"°")Returns a 3×3 rotation matrix that maps crystallophysical coordinates (x, y, z) to dielectric coordinates (X, Y, Z), based on an axis assignment tuple like (:Z, :X, :Y).
An optional additional rotation can be applied around the assigned dielectric axis (:X, :Y, or :Z) by an angle phi. This is used to align the nonlinear tensor with the dielectric frame expected by the rest of the simulation.