Reference

This page contains the full source code documentation and the literature referenced througout this website.

Literature

[1]
B. A. Madruga, C. C. Dorian, M. Sehgal, A. J. Silva, M. Shtrahman, D. Aharoni and P. Golshani. Open-source, high performance miniature multiphoton microscopy systems for freely behaving animals, bioRxiv (2024).
[2]
A. Manny. Simulation of the contrast mechanisms of a heterodyne two-color interferometer. Master's thesis, Faculty of Aerospace Engineering and Geodesy, University of Stuttgart (Dec 2024).
[3]
[4]
[5]
G. Fowles. Introduction to Modern Optics. Dover Books on Physics Series (Dover Publications, 1989).
[6]
M. Ware and J. Peatross. Physics of Light and Optics (Black & White) (Brigham Young University, Department of Physics, 2015).
[7]
B. Saleh and M. Teich. Fundamentals of Photonics. Wiley Series in Pure and Applied Optics (Wiley, 2019).
[8]
[9]
J. Ashcraft, poke v0.1.0, https://zenodo.org/10.5281/zenodo.7117214 (Sep 2022).
[10]
G. Wanner, E. Kochkina, C. Mahrdt, V. Müller, S. Schuster, G. Heinzel and K. Danzmann. Simulating laser interferometers for missions such as (E)Lisa, Lisa pathfinder and Grace follow-on. In: International Conference on Space Optics — ICSO 2014, Vol. 10563, edited by Z. Sodnik, B. Cugny and N. Karafolas (International Society for Optics and Photonics, SPIE, 2017); p. 105632G.
[11]
[12]

Index

BeamletOptics.NullableType
Nullable{T}

An alias which results in Union{T, Nothing} to provide a shorter notation for struct fields which can containing nothing.

source
BeamletOptics.NullableVectorType
NullableVector{T}

An alias which results in Union{Vector{T}, Nothing} to provide a shorter notation for struct fields which can containing nothing.

source
BeamletOptics.RefractiveIndexType
RefractiveIndex

Union type that represents valid means to pass a refractive index n to e.g. AbstractObjects. The core assumption is that:

  1. the refractive index is callable with a single Number argument λ to represent the wavelength in [m]
  2. the return value is a single Number value for the refractive index

Refer to e.g. DiscreteRefractiveIndex.

source
BeamletOptics.AbstractBeamType
AbstractBeam{T <: Real, R <: AbstractRay{T}}

A generic type for a container type which holds rays, beams etc.

Parametrization:

A subtype of AbstractBeam is parameterized by its main data type T <: Real, as well as the underlying ray representation R <: AbstractRay{T}. If a beam is to be compatible with different AbstractRay implementations, it must be parameterized by T and R. However, it can also be set to a fixed type for T and R, i.e. MyBeam <: AbstractBeam{Float32, MyRay}.

Implementation reqs.

Subtypes of AbstractBeam must implement the following:

Fields:

  • parent: a Nullable field that holds the same type as the subtype, used for tree navigation
  • children: a vector that holds the same type as the subtype, used for sub-beam tracking, i.e. beamsplitting

Functions:

  • _modify_beam_head!: modifies the beam path for retracing purposes
  • _last_beam_intersection: returns the last Beam intersection
source
BeamletOptics.AbstractBeamGroupType
AbstractBeamGroup

Provides a generic container type interface for bundles of Beams. This interface assumes that there exists a central beam around which the bundle propagates, e.g. akin to an optical axis.

Implementation reqs.

Subtypes of AbstractBeamGroup must implement the following:

Fields:

  • beams: a vector or tuple of Beams

Functions:

If the beams field does not exist, the following getters must be dispatched:

  • beams: getter for the beams field or equivalent return type
  • position: getter for the starting position of the central Beam
  • direction: getter for the starting direction of the central Beam
  • wavelength: getter for the common wavelength of the beam bundle
source
BeamletOptics.AbstractBeamsplitterType
AbstractBeamsplitter <: AbstractObject

A generic type to represent an AbstractObject that splits incoming beams by reflection and transmission.

Implementation reqs.

Subtypes of AbstractBeamsplitter should implement all supertype requirements.

Interaction logic

After intersection with the AbstractShape at which the beam splitting occurs, the interact3d function should appended the transmitted and reflected sub-beams to the parent beam via the children! function. The interact3d function should then return nothing in order to stop the tracing of the parent beam.

Appending convention

As a convention, when splitting an incoming beam, the order of children appended to the parent beam should be

  1. transmitted beam
  2. reflected beam

Functions

  • interact3d: see above
  • _beamsplitter_transmitted_beam: optional helper function
  • _beamsplitter_reflected_beam: optical helper function
source
BeamletOptics.AbstractCylindricalSurfaceSDFType
AbstractCylindricalSurfaceSDF <: AbstractLensSDF{T}

A class of AbstractSDFs that can be used to represent cylindric non-rotationally symmetric lens surfaces. It is implicity assumed that all surfaces are represented by closed volumes for ray-tracing correctness.

Implementation reqs.

Subtypes of AbstractCylindricalSurfaceSDF must implement the following additional methods.

Functions:

  • height: this function returns the height of the cylinder part of the surface
  • radius: this function returns the radius of the cylinder surface curvature
source
BeamletOptics.AbstractDetectorType
AbstractDetector <: AbstractObject

A generic representation of a detector that evaluates AbstractBeam data during interaction. Refer to Photodetector for more information.

Implementation reqs.

Subtypes of AbstractDetector should implement all supertype requirements as well as:

Functions

  • interact3d: see e.g. Photodetector for reference
  • empty!: resets data stored in the detector, see below

Additional information

The information provided below applies to the standard functional implementation of this type and may be overwritten by specialized subtypes.

Data mutability

In order to model field superposition effects, the concrete implementation of an AbstractDetector should be a mutable struct with a constant shape field. This is necessary, since (sub-)beams will interact sequentially with the detector during solve_system!. Only if the data can be accumulated sequentially, multiple beam interactions can be captured for a complex system, e.g. an interferometer.

Data reset

Since e.g. E-field data is supposed to be accumulated by mutability of the detector data, the burden of resetting the data for a new solver call is placed on the user. This function should be called empty!.

source
BeamletOptics.AbstractInteractionType
AbstractInteraction

Describes how an AbstractBeam and an AbstractObject interact with each other. This data type stores information from the interact3d function and provides it to the solver. The solver can use this data to extend the AbstractBeam.

Implementation reqs.

Subtypes of AbstractInteraction must implement the following:

Fields

  • hint: a nullable Hint for the solver (optional but recommended)

Beam data

It is required that concrete implementations of this type provide some form of data on how to extend the beam. For instance, refer to BeamInteraction and GaussianBeamletInteraction.

source
BeamletOptics.AbstractLensSDFType
AbstractLensSDF

A class of AbstractSDFs that can be used to represent rotationally symmetric lens surfaces. It is implicity assumed that all surfaces are represented by closed volumes for ray-tracing correctness.

Implementation reqs.

Subtypes of AbstractLensSDF must implement the following:

Functions:

  • thickness: this function returns the material thickness of the element along its symmetry axis
  • diameter: this function returns the outer diameter of the element
Shape orientation

For easy compatibility between subtypes, the follwing requirements should be fulfilled:

  1. Symmetry axis aligned onto the y-axis
  2. Surface contour aligned towards negative y-values
  3. Surface point with min(y) should satisfy min(y) = 0 on the symmetry axis
source
BeamletOptics.AbstractMeshType
AbstractMesh <: AbstractShape

A generic type for an shape whose volume can be described by a mesh. Must have a field mesh of type Mesh. See also Mesh{T}.

source
BeamletOptics.AbstractObjectType
AbstractObject

A generic type for 2D/3D objects that can be used to model optical elements. The geometry of the object is represented via a AbstractShape. The optical effect that occurs between the object and an incoming ray/beam of light is modeled via its interact3d method.

Implementation reqs.

Subtypes of AbstractObject must implement the following:

Shape trait

An AbstractObject can consist of a single AbstractShape, e.g. a lens element, or a collection of functionally dependant shapes, e.g. a cube beamsplitter. In order to model this, the API implementation of an AbstractObject requires the definition of the shape trait. This trait allows the dispatch onto specialized methods to handle the kinematic interface and tracing methods for objects consisting of one or more shapes.

  • shape_trait_of: defines the shape type of the AbstractObject, refer to AbstractShapeTrait for more information
Default shape trait

Unless specified otherwise, the shape_trait_of an AbstractObject is defined as SingleShape. This requires object.shape as a dedicated field. For MultiShapes the getter function shape(object) must return a tuple of all shapes that make up the object.

Getters/setters

All kinematic functions defined for the AbstractShape can also be called for a AbstractObject. In this case, the shape trait will define how the specific movement function is dispatched.

Functions:

source
BeamletOptics.AbstractObjectGroupType
AbstractObjectGroup <: AbstractObject

Container type for groups of optical elements, based on a tree-like data structure. Intended for easier kinematic handling of connected elements. See also ObjectGroup for a concrete implementation.

Implementation reqs.

Subtypes of AbstractObjectGroup must implement the following:

Fields:

  • objects: stores objects or additional subgroups of objects, allows for hierarchical structures

Functions:

  • for the kinematic API, all corresponding functions of AbstractObject must be implemented
source
BeamletOptics.AbstractPlateBeamsplitterType
AbstractPlateBeamsplitter <: AbstractBeamsplitter

A generic type to represent an AbstractBeamsplitter that consists of a substrate with a single coated face at which a beam splitting interaction occurs.

Implementation reqs.

Subtypes of AbstractPlateBeamsplitter should implement all supertype reqs. as well as:

Fields

  • coating: a ThinBeamsplitter that represents the splitter coating
  • substrate: a Prism that represents the substrate

Getters/setters

If the concrete implementation does not define the above fields, the following getters must be defined:

Additional information

Object orientation

This interact3d method of this type strongly assumes that the coating is positioned directly upon a single face of the substrate with a 100% fill factor.

Interaction logic

This type uses the Hint-API in order to ensure that the splitting interaction is correctly triggered at the coating.

source
BeamletOptics.AbstractRayType
AbstractRay{T<:Real}

An implementation for a geometrical optics ray in R³. In general, a AbstractRay is described by $\vec{p} + t\cdot\vec{d}$ with $t\in(0,\infty)$. AbstractRays are intended to model the propagation of light between optical interactions according to the laws of geometrical optics. To store the result of a ray tracing solution, refer to AbstractBeam.

Intersections:

Since the length of a ray can not be known before solving an optical system, the Intersection-type is used. This Nullable type can represent the intersection with an optical element, or lack thereof.

Implementation reqs.

Subtypes of AbstractBeam must implement the following:

Fields:

  • pos: a R³-vector that stores the current position $\vec{p}$
  • dir: a R³-vector that stores the current direction $\vec{d}$
  • intersection: a Nullable field that stores the current [Intersection] or nothing
  • λ: wavelength in [m]
  • n: refractive index along the ray path

Additional information

Ray length

Base.length: this function is used to return the length of the AbstractRay (if no intersection exists, the ray length is Inf). The opl keyword can be used to obtain the optical path length instead.

Ray direction

Many functions assume that the direction vector has unit length (i.e. $|\vec{p}| = 1$). Violating this assumption might lead to spurious results.

source
BeamletOptics.AbstractReflectiveOpticType
AbstractReflectiveOptic <: AbstractObject

A generic type to represent an [AbstractObject] which reflects incoming rays.

Implementation reqs.

Subtypes of AbstractReflectiveOptic should implement all supertype reqs. as well as:

Fields

  • no specific fields required

Getters/setters

  • none required

Functions

  • interact3d: the interaction logic should be akin to reflection3d for each surface crossing

Additional information

The information provided below applies to the standard functional implementation of this type and may be overwritten by specialized subtypes.

Polarization ray tracing

Fresnel coefficients during reflection are set such that no reflection losses occur (i.e. |rₚ| = |rₛ| = 1).

source
BeamletOptics.AbstractRefractiveOpticType
AbstractRefractiveOptic <: AbstractObject

A generic type to represent an AbstractObject that refracts incoming rays.

Implementation reqs.

Subtypes of AbstractRefractiveOptic should implement all supertype reqs. as well as:

Fields

Getters/setters

  • refractive_index: gets the ref. index data of the optic

Functions

  • interact3d: the interaction logic should be akin to refraction3d for each surface crossing

Additional information

The information provided below applies to the standard functional implementation of this type and may be overwritten by specialized subtypes.

Uniform ref. index

It is assumed that the optic consists of a single transparent material with a homogeneous refractive index n. It does not consider coated surfaces.

Polarization ray tracing

Fresnel coefficients at the point of refraction are calculated via the fresnel_coefficients function with the refractive index data of the substrate and the previous medium.

source
BeamletOptics.AbstractRotationallySymmetricSurfaceType
AbstractRotationallySymmetricSurface{T} <: AbstractSurface{T}

A surface type which is rotationally symmetric around one axis.

Implementation reqs.

Subtypes of AbstractShape should implement the following:

Getters/setters

  • radius : Returns the radius of curvature of the AbstractRotationallySymmetricSurface
  • diameter : Returns the clear optical diameter of the AbstractRotationallySymmetricSurface
  • mechanical_diameter : Returns the mechanical diameter of the AbstractRotationallySymmetricSurface
  • edge_sag : Returns the edge sagitta of the AbstractRotationallySymmetricSurface

Functions:

source
BeamletOptics.AbstractSDFType
AbstractSDF <: AbstractShape

Provides a shape function based on signed distance functions. See https://iquilezles.org/articles/distfunctions/ for more information.

Implementation reqs.

Subtypes of AbstractSDF should implement all reqs. of AbstractShape as well as the following:

Functions

  • sdf(::AbstractSDF, point): a function that returns the signed distance for a point in 3D space
source
BeamletOptics.AbstractShapeType
AbstractShape{T<:Real}

A generic type for a shape that exists in 3D-space. Must have a position and orientation.
Types used to describe the geometry of a shape should be subtypes of Real.

Implementation reqs.

Subtypes of AbstractShape should implement the following:

Fields:

  • pos: a 3D-vector that stores the current position of the object-specific coordinate system
  • dir: a 3x3-matrix that represents the orthonormal basis of the object and therefore, the orientation

Getters/setters

  • position / position!: gets or sets the position vector of the AbstractShape
  • orientation / orientation!: gets or sets the orientation matrix of the AbstractShape

Kinematic:

Ray Tracing:

  • intersect3d: returns the intersection between an AbstractShape and AbstractRay, or lack thereof. See also Intersection

Rendering (with Makie):

Refer to the render! documentation.

source
BeamletOptics.AbstractSphericalSurfaceSDFType
AbstractSphericalSurfaceSDF{T} <: AbstractSDF{T}

An abstract type for SDF-based volumes which represent spherical lens surfaces, i.e. ConvexSphericalSurfaceSDF or ConcaveSphericalSurfaceSDF.

Implementation reqs.

Subtypes of AbstractSphericalSurfaceSDF should implement all supertype reqs. as well as the following:

Fields:

  • radius: the radius of curvature
  • diameter: the lens outer diameter
  • sag: the lens sagitta

Lens construction

It is intended that practical lens shapes are constructed from AbstractSphericalSurfaceSDFs using the UnionSDF type.

source
BeamletOptics.AbstractSurfaceType
AbstractSurface{T}

A generic type for a surface which is basically an information storage type in order to build shapes (volumes) from a combination of surfaces.

source
BeamletOptics.AbstractSystemType
AbstractSystem

A generic representation of a system of optical elements.

Implementation reqs.

Subtypes of AbstractSystem must implement the following:

Fields:

  • objects: a vector or tuple of AbstractObjects that make up the system
  • n: (optional) RefractiveIndex of the surrounding medium, default value is 1.0

Functions:

  • refractive_index: returns the RefractiveIndex n of the system medium, see above
source
BeamletOptics.AconcaveCylinderSDFMethod
AconcaveCylinderSDF(radius, diameter, height)

Constructs an aconcave cylinder cutout with radius r, diameter d and height h in [m]. The acylindric shape is defined by its conic_constant and the coefficients for the even aspheric polynomoial.

source
BeamletOptics.AconvexCylinderSDFMethod
AconvexCylinderSDF(radius, diameter, height)

Constructs an aconvex cylinder with radius r, diameter d and height h in [m]. The acylindric shape is defined by its conic_constant and the coefficients for the even aspheric polynomoial.

source
BeamletOptics.AcylindricalSurfaceType
AcylindricalSurface{T} <: AbstractAcylindricalSurface{T}

A type representing an acylindric optical surface defined by its radius of curvature, diameter, height, mechanical diameter, conic constant and even aspheric coefficients. It is therefore a cylindric surface with a deviation from the perfect cylindric shape.

Fields

  • radius::T: The radius of curvature of the curved surface.
  • diameter::T: The clear (optical) aperture of the surface.
  • height::T : The height/length of the uncurved surface direction
  • conic_constant::T : The conic_constant of the curved surface
  • coefficients::Vector{T} : The coefficients of the even aspherical equation for the curved surface.
  • mechanical_diameter::T: The overall mechanical diameter of the surface. In many cases, this is equal to the optical diameter, but it can be set independently if the mechanical mount requires a larger dimension.
source
BeamletOptics.AcylindricalSurfaceMethod
AcylindricalSurface(radius, diameter, height, conic_constant, coefficients)

Construct a CylindricalSurface given the radius of curvature, optical diameter and height. This constructor automatically sets the mechanical diameter equal to the optical diameter.

Arguments

  • radius: The radius of curvature of the curved surface.
  • diameter: The clear (optical) diameter of the surface.
  • height: The height/length of the uncurved surface direction.
  • conic_constant::T : The conic_constant of the curved surface
  • coefficients::Vector{T} : The coefficients of the even aspherical equation for the curved surface.
source
BeamletOptics.BeamType
Beam{T, R <: AbstractRay{T}} <: AbstractBeam{T, R}

Stores the rays that are calculated from geometric optics when propagating through an optical system. The Beam type is parametrically defined by the AbstractRay subtype that it stores.

Fields

  • rays: vector of AbstractRay objects, representing the rays that make up the beam
  • parent: reference to the parent beam, if any (Nullable to account for the root beam which has no parent)
  • children: vector of child beams, each child beam represents a branching or bifurcation of the original beam, i.e. beam-splitting
source
BeamletOptics.BoxSDFType
BoxSDF <: AbstractSDF

Implements the box SDF with edge lengths x, y, and z. Note that these values are stored in the dimensions field as:

  • dimensions::Point3 = ( leninx/2, leniny/2, leninz/2,

)

source
BeamletOptics.CircularFlatSurfaceType
CircularFlatSurface{T} <: AbstractRotationallySurface{T}

A type representing a planar circular surface, which is only parametrized by its diameter.

Fields

  • diameter::T: The diameter of the planar surface
source
BeamletOptics.CollimatedSourceType
CollimatedSource <: AbstractBeamGroup

Represents a parallel bundle of Beams being emitted from a disk in space.

Fields

  • beams: a vector of all Beams originating from the source
  • diameter: the diameter of the outermost beam ring

Functions

  • diameter: returns the diameter of the source
source
BeamletOptics.CollimatedSourceMethod
CollimatedSource(pos, dir, diameter, λ; num_rings, num_rays)

Spawns a bundle of collimated Beams at the specified position and direction. The source is modelled as a ring of concentric beam rings around the center beam. The amount of beam rings between the center ray and outer diameter can be specified via num_rings.

Info

Note that for correct sampling, the number of rays should be atleast 20x the number of rings.

Arguments

The following inputs and arguments can be used to configure the CollimatedSource:

Inputs

  • pos: center beam starting position
  • dir: center beam starting direction
  • diameter: outer beam bundle diameter in [m]
  • λ = 1e-6: wavelength in [m], default val. is 1000 nm

Keyword Arguments

  • num_rings: number of concentric beam rings, default is 10
  • num_rays: total number of rays in the source, default is 100x num_rings
Warning

The orthogonal basis vectors for the beam generation are generated randomly.

source
BeamletOptics.ConcaveAsphericalSurfaceSDFType
ConcaveAsphericalSurfaceSDF

Constructs an aspheric lens with a concave-like surface according to ISO10110.

Note

Currently, it is assumed that the aspheric surface is concave if the radius is negative. There might be unexpected effects for complex shapes which do not show a generally concave behavior.

  • coefficients : (even) coefficients of the asphere.
  • radius : radius of the lens (negative!)
  • conic_constant : conic constant of the lens surface
  • diameter: lens diameter
  • mechanical_diameter: mechanical lens diameter, defaults to be identical to the lens diameter, Otherwise an outer ring section will be added to the lens, if mechanical_diameter > diameter.
source
BeamletOptics.ConcaveSphericalMirrorMethod
ConcaveSphericalMirror(radius, thickness, diameter)

Constructor for a spherical mirror with a concave reflecting surface. The component is aligned with the positive y-axis. See also ConcaveSphericalMirror.

Inputs

  • radius: the spherical surface radius of curvature in [m]
  • thickness: substrate thickness in [m]
  • diameter: mirror outer diameter in [m]
source
BeamletOptics.ConcaveSphericalSurfaceSDFType
ConcaveSphericalSurfaceSDF

AbstractSDF-based representation of a concave spherical lens surface. When constructed, it is assumed that the plano-surface lies at the origin and the optical axis is aligned with the y-axis. The concave surface is orientated towards negative y-values for R > 0 and vice versa.

Fields:

  • radius: the radius of curvature of the convex spherical surface.
  • diameter: the outer diameter of the lens surface
  • sag: the sagitta of the opposing convex shape
source
BeamletOptics.ConvexAsphericalSurfaceSDFType
ConvexAsphericalSurfaceSDF

Constructs an aspheric lens with a convex-like surface according to ISO10110.

Note

Currently, it is assumed that the aspheric surface is convex if the radius is positive. There might be unexpected effects for complex shapes which do not show a generally convex behavior.

  • coefficients : (even) coefficients of the asphere.
  • radius : radius of the lens
  • conic_constant : conic constant of the lens surface
  • diameter: lens diameter
source
BeamletOptics.ConvexSphericalSurfaceSDFType
ConvexSphericalSurfaceSDF

AbstractSDF-based representation of a convex spherical lens surface. When constructed, it is assumed that the plano-surface lies at the origin and the optical axis is aligned with the y-axis. The convex surface is orientated towards negative y-values for R > 0 and vice versa.

Fields:

  • radius: the radius of curvature of the concave spherical surface.
  • diameter: the outer diameter of the lens surface
  • sag: the sagitta of the convex shape
  • height: the sphere cutoff height, see also CutSphereSDF
source
BeamletOptics.CubeBeamsplitterType
CubeBeamsplitter <: AbstractBeamsplitter

A cuboid beamsplitter where the splitting interaction occurs between two RightAnglePrisms. For more information refer to the AbstractPlateBeamsplitter docs.

Fields

Additional information

Hints and interaction logic

In order to model gap-free beam propagation, the interact3d model relies heavily on the Hint-API. If the front or back substrate is hit, the Hint will ensure that the beam intersects the coating.

source
BeamletOptics.CubeBeamsplitterMethod
CubeBeamsplitter(leg_length, n; reflectance=0.5)

Creates a CubeBeamsplitter. The cuboid is centered at the origin. The splitter coating is orientated at a 45° angle with respect to the y-axis.

Inputs

  • leg_length: the x-, y- and z-edge length in [m]
  • n: the RefractiveIndex of the front and back prism

Keywords

  • reflectance: defines the splitting ratio in [-], i.e. R = 0 ... 1.0
source
BeamletOptics.CutSphereSDFType
CutSphereSDF <: AbstractSDF

Implements SDF of a sphere which is cut off in the x-z-plane at some point along the y-axis.

source
BeamletOptics.CylinderSDFType
CylinderSDF <: AbstractSDF

Implements cylinder SDF. Cylinder is initially orientated along the y-axis and symmetrical in x-z.

source
BeamletOptics.CylindricalSurfaceType
CylindricalSurface{T} <: AbstractCylindricalSurface{T}

A type representing a cylindric optical surface defined by its radius of curvature, diameter, height and mechanical diameter

Fields

  • radius::T: The radius of curvature of the curved surface.
  • diameter::T: The clear (optical) aperture of the surface.
  • height::T : The height/length of the uncurved surface direction
  • mechanical_diameter::T: The overall mechanical diameter of the surface. In many cases, this is equal to the optical diameter, but it can be set independently if the mechanical mount requires a larger dimension.
source
BeamletOptics.CylindricalSurfaceMethod
CylindricalSurface(radius::T, diameter::T, height::T) where T

Construct a CylindricalSurface given the radius of curvature, optical diameter and height. This constructor automatically sets the mechanical diameter equal to the optical diameter.

Arguments

  • radius::T: The radius of curvature of the curved surface.
  • diameter::T: The clear (optical) diameter of the surface.
  • height::T: The height/length of the uncurved surface direction.
source
BeamletOptics.DiscreteRefractiveIndexType
DiscreteRefractiveIndex{T}

Represents a incomplete set of dispersion data where for each exact wavelength one refractive index value is stored in the data field. Can be called like a function n = n(λ). Does not interpolate between data points. Refer to RefractiveIndex for more information.

source
BeamletOptics.DoubletLensType
DoubletLens

Represents a two-component cemented doublet lens with two respective refractive indices n = n(λ). See also SphericalDoubletLens.

Fields

  • front: front Lens component
  • back: back Lens component

Additional information

Air gap

This component type strongly assumes that both lenses are mounted fully flush with respect to each other. Gaps between the components might lead to incorrect results.

source
BeamletOptics.EvenAsphericalSurfaceType
EvenAsphericalSurface{T} <: AbstractRotationallySymmetricSurface{T}

A type representing an aspherical optical surface defined by its radius of curvature, clear (optical) diameter, conic constant, aspheric coefficients and mechanical diameter. This surface is rotationally symmetric about its optical axis.

Fields

  • spherical::SphericalSurface{T}: The base spherical surface portion of the asphere.
  • conic_constant::T: The conic constant defining the deviation from a spherical shape.
  • coefficients::AbstractVector{T}: A vector of even aspheric coefficients for higher-order corrections.
source
BeamletOptics.EvenAsphericalSurfaceMethod
EvenAsphericalSurface(radius, diameter, conic_constant, coefficients::AbstractVector, mechanical_diameter = diameter)

Construct a EvenAsphericalSurface given the radius of curvature, the optical diameter, conic constant, the aspheric coefficients and optionally the mechanical diameter. This constructor automatically sets the mechanical diameter equal to the optical diameter.

Arguments

  • radius: The radius of curvature of the base spherical surface.
  • diameter: The clear (optical) diameter of the surface.
  • conic_constant: The conic constant defining the deviation from a spherical shape.
  • coefficients::AbstractVector: A vector of even aspheric coefficients for higher-order corrections.
  • mechanical_diameter: The mechanical diameter of the surface; if not provided, it defaults to diameter.
source
BeamletOptics.GaussianBeamletType
GaussianBeamlet{T} <: AbstractBeam{T, Ray{T}}

Ray representation of the stigmatic Gaussian beam as per J. Arnaud (1985). The beam quality M2 is fully considered via the divergence angle. The formalism for the beam parameter calculation is based on the following publications:

Jacques Arnaud, "Representation of Gaussian beams by complex rays," Appl. Opt. 24, 538-543 (1985)

and

Donald DeJager and Mark Noethen, "Gaussian beam parameters that use Coddington-based Y-NU paraprincipal ray tracing," Appl. Opt. 31, 2199-2205 (1992)

Fields

  • chief: a Beam of Rays to store the chief ray
  • waist: a Beam of Rays to store the waist ray
  • divergence: a Beam of Rays to store the divergence ray
  • λ: beam wavelength in [m]
  • w0: local beam waist radius in [m]
  • E0: complex field value in [V/m]
  • parent: reference to the parent beam, if any (Nullable to account for the root beam which has no parent)
  • children: vector of child beams, each child beam represents a branching or bifurcation of the original beam, i.e. beam-splitting

Additional information

Beam parameters

Parameters of the beam, e.g. $w(z)$ or $R(z)$, can be obtained through the gauss_parameters function.

Astigmatism and abberations

It is assumed, but not forbidden, that the optical system contains non-flat or non-parabolic beam-surface-interactions that cause the beam to obtain astigmatism or higher-order abberations. These can not be represented by the GaussianBeamlet.

source
BeamletOptics.GaussianBeamletMethod
GaussianBeamlet(position, direction, λ, w0; kwargs...)

Constructs a Gaussian beamlet at its waist with the specified beam parameters.

Arguments

The following inputs and arguments can be used to configure the beamlet:

Inputs

  • position: origin of the beamlet
  • direction: direction of the beamlet
  • λ: wavelength of the beamlet in [m]. Default value is 1000 nm.
  • w0: beam waist (radius) in [m]. Default value is 1 mm.

Keyword Arguments

  • M2: beam quality factor. Default is 1
  • P0: beam total power in [W]. Default is 1 mW
  • z0: beam waist offset in [m]. Default is 0 m
  • support: Nullable support vector for the construction of the waist and div rays

Additional information

Waist offset

The z0 keyword arg. can be used in order to spawn a beam where the waist is not located at the specified position, but rather at an offset z0 in [m] along the chief ray axis.

Support vector

In order to calculate the basis vectors required for the beamlet construction, a random orthogonal vector is chosen. If results fluctuate due to the randomness of this vector, make sure to specify a fixed orthogonal support vector.

source
BeamletOptics.HintType
Hint

A Hint can be passed as part of an AbstractInteraction and will inform the tracing algorithm about which AbstractObject in the AbstractSystem will be hit next.

Info

The hint does not need to result in a guaranteed Intersection.

Fields

  • object: the object that might or will be intersected as next
  • shape: the underlying shape that will be intersected next, i.e. shape(object), relevant for multi-shape objects
source
BeamletOptics.IntersectionType
Intersection{T}

Stores data calculated by the intersect3d method. This information can be reused, i.e. for retracing.

Fields:

  • object: a Nullable reference to the AbstractObject that has been hit (optional but recommended)
  • shape: a Nullable reference to the AbstractShape of the object that has been hit (optional but recommended)
  • t: length of the ray parametrization in [m]
  • n: normal vector at the point of intersection
source
BeamletOptics.LensType
Lens{T, S <: AbstractShape{T}, N <: RefractiveIndex} <: AbstractRefractiveOptic{T, S, N}

Represents an uncoated Lens with a homogeneous RefractiveIndex n = n(λ). Refer to the Lens and SphericalLens constructors for more information on how to generate lenses.

Fields

Additional information

Refractive index

The chromatic dispersion of the lens is represented by a λ-dependent function for n and must be provided by the user. For testing purposes, an anonymous function, e.g. λ -> 1.5 can be passed such that the lens has the same refractive index for all wavelengths.

source
BeamletOptics.LensMethod
 Lens(front_surface::AbstractCylindricalSurface, back_surface::AbstractCylindricalSurface, center_thickness::Real, n::RefractiveIndex)

Constructs a new Lens object using the cylindric surface specifications front_surface and back_surface and the center_thickness. These inputs are used to construct a UnionSDF that consists of the appropriate sub-SDFs to represent the shape of the lens.

This method of Lens is specific for cylindric lenses and has some limitations: - The cylinder height of both surfaces has to be identical - No mixture with non-cylindric surfaces is supported at the moment

The material properties are supplied via the n parameter.

Additional information

Radius of curvature (ROC) sign definition

The ROC is defined to be positive if the center is to the right of the surface. Otherwise it is negative.

source
BeamletOptics.LensMethod
 Lens(front_surface::AbstractRotationallySymmetricSurface, back_surface::AbstractRotationallySymmetricSurface, center_thickness::Real, n::RefractiveIndex)

Constructs a new Lens object using the surface specifications front_surface and back_surface and the center_thickness. These inputs are used to construct a UnionSDF that consists of the appropriate sub-SDFs to represent the shape of the lens.

The material properties are supplied via the n parameter.

Additional information

Radius of curvature (ROC) sign definition

The ROC is defined to be positive if the center is to the right of the surface. Otherwise it is negative.

Meniscus

If your specification results in a meniscus lens, only spherical meniscus lenses are supported at the moment.

source
BeamletOptics.MeniscusLensSDFType
MeniscusLensSDF

AbstractSDF-based representation of a positive or negative meniscus lens. When constructed, it is assumed that the lens lies at the origin and the optical axis is aligned with the y-axis. Parameters that lead to a sharp lens edge will cause an error.

Notes

Radius of curvature sign convention

The ROC is defined to be positive if the center is to the right of the surface. Otherwise it is negative.

Fields:

  • convex: the convex part of the lens composite SDF
  • cylinder: the cylindrical part of the lens composite SDF
  • concave: the concave part of the lens composite SDF
  • thickness: lens thickness on the optical axis
source
BeamletOptics.MeniscusLensSDFMethod
MeniscusLensSDF(r1::R1, r2::R2, l::L, d::D, md::MD)

Constructs a positive or negative MeniscusLensSDF with:

  • r1: front surface radis or curvature
  • r2: back surface radis or curvature
  • l: lens thickness
  • d: lens diameter, default value is one inch
  • md: mechanical lens diameter, must be > d
source
BeamletOptics.MeshType
Mesh <: AbstractMesh

Contains the STL mesh information for an arbitrary shape, that is the vertices that make up the mesh and a matrix of faces, i.e. the connectivity matrix of the mesh. The data is read in using the FileIO.jl and MeshIO.jl packages. Translations and rotations of the mesh are directly saved in absolute coordinates in the vertex matrix. For orientation and translation tracking, a positional and directional matrix are stored.

Fields

  • vertices: (m x 3)-matrix that stores the edge points of all triangles
  • faces: (n x 3)-matrix that stores the connectivity data for all faces
  • dir: (3 x 3)-matrix that represents the current orientation of the mesh
  • pos: 3-element vector that is used as the mesh location reference
  • scale: scalar value that represents the current scale of the original mesh
source
BeamletOptics.MeshMethod
Mesh(mesh)

Parametric type constructor for struct Mesh. Takes data of type GeometryBasics.Mesh and extracts the vertices and faces. The mesh is initialized at the global origin. Data type of Mesh is variably selected based on type of vertex data (i.e Float32). Mesh data is scaled by factor 1e-3, assuming m scale.

source
BeamletOptics.MirrorType
Mirror{S <: AbstractShape} <: AbstractReflectiveOptic

Concrete implementation of a perfect mirror (R = 1) with arbitrary shape.

Reflecting surfaces

It is important to consider that all surfaces of this mirror type are reflecting!

source
BeamletOptics.MultiShapeType
MultiShape <: AbstractShapeTrait

Represents that the AbstractObject consists of a two or more AbstractShapes.

AbstractObject implementation reqs.

If shape_trait_of(::Foo) = MultiShape() is defined, Foo must implement the following:

Functions

  • shape(::Foo): a getter function that returns a Tuple of all relevant shapes, e.g. (foo.front, foo.middle, foo.back)

Additional information

Kinematic center

Unless specified otherwise by dispatching position / position! and orientation / orientation! onto custom pos and dir data fields, the position and orientation of the first element returned by shape(object) will be used as the kinematic center for e.g. translate3d!.

source
BeamletOptics.NonInteractableObjectType
NonInteractableObject

A passive AbstractObject which does not interact with the ray tracing simulation but can be moved via the kinematic API.

Fields

Usage

This type is intended mainly for visualization purposes, e.g. kinematic mount Meshs, or similar applications. In essence, this object behaves fully transparent. The intersect3d and interact3d methods default to nothing.

source
BeamletOptics.ObjectGroupType
ObjectGroup <: AbstractObjectGroup

A tree-like storage container for groups of objects. Can store individual objects and subgroups. Main purpose is handling of, i.e., groups of lenses.

Fields

  • center: a point in 3D space which is regarded as the reference origin of the group
  • dir: a 3x3 matrix that describes the common orientation of the group
  • objects: stores AbstractObject, can also store subgroups of type AbstractObjectGroup

Kinematic

A ObjectGroup implements the kinematic functions of AbstractObject. The following logic is applied to

  • translate3d!: all objects in the group are translated by the offset vector
  • translate_to3d!: all objects are moved in parallel such that the group center is equal to the target position
  • rotate3d!: all objects are rotated around the center point with respect to their relative position
source
BeamletOptics.PSFDetectorType
PSFDetector{T} <: AbstractDetector{T, Mesh{T}}

Represents a flat quadratic surface in R³ for capturing the point spread function of a <: AbstractSystem The active surface is discretized in the local R² x-y-coordinate system. Ray hits are recorded by the corresponding interact3d method.

Fields

  • shape: geometry of the active surface, must represent 2D-field in x any y dimensions
  • data : A vector of PSFData capturing the information about each ray hit.

Additional information

Reset behavior

The PSFDetector must be reset between each call of solve_system! in order to overwrite previous results using the empty! function. Otherwise, the current result will be added onto the previous result.

Supported beams

Currently, only the Beam is supported.

source
BeamletOptics.PSFDetectorMethod
Photodetector(width)

Spawns a quadratic rectangular 2D PSFDetector that is aligned with the positive y-axis. Refer to the type docs for more information.

Inputs:

  • width: edge length in [m]
source
BeamletOptics.PhotodetectorType
Photodetector{T, S <: AbstractShape{T}} <: AbstractDetector{T, S}

Represents a flat rectangular or quadratic surface in R³ that is the active surface of a photodetector. The active surface is discretized in the local R² x-y-coordinate system. Field contributions Eᵢ are added by the corresponding interact3d method.

Fields

  • shape: geometry of the active surface, must represent 2D-field in x any y dimensions
  • x: linear range of local x-coordinates
  • y: linear range of local y-coordinates
  • field: size(x) by size(y) matrix of complex values to store superposition E₀

Additional information

Reset behavior

The Photodetector must be reset between each call of solve_system! in order to overwrite previous results using the empty! function. Otherwise, the current result will be added onto the previous result.

Supported beams

Currently, only the GaussianBeamlet is supported.

source
BeamletOptics.PhotodetectorMethod
Photodetector(width, n)

Spawns a quadratic rectangular 2D Photodetector that is aligned with the positive y-axis. Refer to the type docs for more information.

Inputs:

  • width: edge length in [m]
  • n: field discretization factor, higher results in more computational cost
source
BeamletOptics.PlanoSurfaceSDFType
PlanoSurfaceSDF

AbstractSDF-based representation of two flat optical surfaces, i.e. equivalent to the CylinderSDF. When constructed, it is assumed that the first flat surface lies at the origin and the optical axis is aligned with the positive y-axis.

Fields:

  • diameter: the outer diameter of the circular flat lens surface
  • thickness: the distance between the flat surfaces
source
BeamletOptics.PointSourceType
PointSource <: AbstractBeamGroup

Represents a cone of Beams being emitted from a single point in space.

Fields

  • beams: a vector of all Beams originating from the source
  • NA: the numerical_aperture of the point source spread angle

Functions

  • numerical_aperture: returns the NA of the source
source
BeamletOptics.PointSourceMethod
PointSource(pos, dir, θ, λ; num_rings, num_rays)

Spawns a point source of Beams at the specified position and direction. The point source is modelled as a collection of concentric beam fans centered around the center beam. The amount of beam rings between the center ray and half-spread-angle θ can be specified via num_rings.

Info

Note that for correct sampling, the number of rays should be atleast 20x the number of rings.

Arguments

The following inputs and arguments can be used to configure the PointSource:

Inputs

  • pos: center beam starting position
  • dir: center beam starting direction
  • θ: half spread angle in rad
  • λ = 1e-6: wavelength in [m], default val. is 1000 nm

Keyword Arguments

  • num_rings: number of concentric beam rings, default is 10
  • num_rays: total number of rays in the source, default is 100x num_rings
Warning

The orthogonal basis vectors for the beam generation are generated randomly.

source
BeamletOptics.PolarizedRayType
PolarizedRay{T} <: AbstractRay{T}

A ray type to model the propagation of an electric field vector based on the publication:

Yun, Garam, Karlton Crabtree, and Russell A. Chipman. "Three-dimensional polarization ray-tracing calculus I: definition and diattenuation." Applied optics 50.18 (2011): 2855-2865.

The geometrical ray description is identical to the standard Ray. The polarization interaction can be described in local s-p-coordinates but must be transformed into global coordinates using the method described in the publication above, see also _calculate_global_E0.

Fields

  • pos: a point in R³ that describes the Ray origin
  • dir: a normalized vector in R³ that describes the Ray direction
  • intersection: refer to Intersection
  • λ: wavelength in [m]
  • n: refractive index along the beam path
  • E0: complex-valued 3-tuple to represent the electric field in global coordinates

Jones matrices

In local coordinates the Jones matrices in the case of reflection/refraction are defined as

  • reflection: [-rₛ 0; 0 rₚ]
  • transmission: [tₛ 0; 0 tₚ]

where r and t are the complex-valued Fresnel coefficients (see also fresnel_coefficients).

Additional information

Field vector

It is assumed that the electric field vector $E_0$ stays orthogonal to the direction of propagation throughout the optical system.

Intensity

E0 can not be converted into an intensity value, since a single PolarizedRay can not directly model the change in intensity during imaging by an optical system.

source
BeamletOptics.PrismType
Prism{T, S <: AbstractShape{T}, N <: RefractiveIndex} <: AbstractRefractiveOptic{T, S, N}

Essentially represents the same functionality as Lens. Refer to its documentation.

source
BeamletOptics.RayType
Ray{T} <: AbstractRay{T}

Mutable struct to store ray information.

Fields

  • pos: a point in R³ that describes the Ray origin
  • dir: a normalized vector in R³ that describes the Ray direction
  • intersection: refer to Intersection
  • λ: wavelength in [m]
  • n: refractive index along the beam path
source
BeamletOptics.RayMethod
Ray(pos, dir, λ=1000e-9)

Constructs a Ray where:

  • pos: is the Ray origin
  • dir: is the Ray direction of propagation, normalized to unit length

Optionally, a wavelength λ can be specified. The start refractive index is assumed to be in vacuum (n = 1).

source
BeamletOptics.RectangularFlatSurfaceType
RectangularFlatSurface{T} <: AbstracCylindricalSurface{T}

A type representing a planar rectangular surface, which is only parametrized by its size.

Fields

  • size::T: The size of the planar surface
source
BeamletOptics.RectangularPlateBeamsplitterType
RectangularPlateBeamsplitter <: AbstractPlateBeamsplitter

A plate beamsplitter with rectangular substrate and a single coated face. For more information refer to the AbstractPlateBeamsplitter docs.

Fields

  • substrate: a rectangular Prism that acts as the substrate
  • coating: a ThinBeamsplitter that acts as the coating

Additional information

Kinematic center

The center of kinematics of this splitter lies at the center of the coating.

source
BeamletOptics.RectangularPlateBeamsplitterMethod
RectangularPlateBeamsplitter(width, height, thickness, n; reflectance=0.5)

Creates a RectangularPlateBeamsplitter. The splitter is aligned with the negative y-axis. The splitter coating is centered at the origin. See also RoundPlateBeamsplitter.

Inputs

  • width: substrate width along the x-axis in [m]
  • height: substrate height along the z-axis in [m]
  • thickness: substrate thickness along the y-axis in [m]
  • n: the RefractiveIndex of the substrate

Keywords

  • reflectance: defines the splitting ratio in [-], i.e. R = 0 ... 1.0
source
BeamletOptics.RetroreflectorType
Retroreflector

A Retroreflector reflects incoming rays back toward their source, independent of the incident angle. The shape is represented by a tetrahedral Mesh.

Fields

  • mesh: shape of the Retroreflector
source
BeamletOptics.RightAnglePrismMirrorMethod
RightAnglePrismMirror(leg_length, height)

Constructs a right angle prism mirror. The primary surface is aligned with the pos. y-axis.

Inputs

  • leg_length: edge length in x and y in [m]
  • height: in z-axis in [m]
source
BeamletOptics.RightAnglePrismSDFType
RightAnglePrismSDF <: AbstractSDF

Implements the SDF of a right angle prism with symmetric leg length l and height h. Note that these values are stored in the dimensions field as:

dimensions::Point3 = ( leglength, # dim in x leglength, # dim in y height, # dim in z )

Alignment

Note that the prism is not aligned with the positive y-axis!

source
BeamletOptics.RingSDFType
RingSDF <: AbstractSDF

Implements the SDF of a ring in the x-z-plane for some distance in the y axis. This allows to add planar outer sections to any SDF which fits inside of the ring.

source
BeamletOptics.RingSDFMethod
RingSDF(inner_radius, width, thickness)

Constructs a ring with inner_radius with a width and some thickness.

source
BeamletOptics.RoundPlateBeamsplitterType
RoundPlateBeamsplitter <: AbstractPlateBeamsplitter

A plate beamsplitter with cylindrical substrate and a single coated face. For more information refer to the AbstractPlateBeamsplitter docs.

Fields

Additional information

Kinematic center

The center of kinematics of this splitter lies at the center of the coating.

source
BeamletOptics.SphericalSurfaceType
SphericalSurface{T} <: AbstractRotationallySymmetricSurface{T}

A type representing a spherical optical surface defined by its radius of curvature, clear (optical) diameter, and mechanical diameter. This surface is rotationally symmetric about its optical axis.

Fields

  • radius::T: The radius of curvature of the spherical surface. A positive value indicates that the center of curvature lies to the right of the vertex (following ISO 10110).
  • diameter::T: The clear (optical) aperture of the surface.
  • mechanical_diameter::T: The overall mechanical diameter of the surface. In many cases, this is equal to the optical diameter, but it can be set independently if the mechanical mount requires a larger dimension.
source
BeamletOptics.SphericalSurfaceMethod
SphericalSurface(radius, diameter)

Construct a SphericalSurface given the radius of curvature and the optical diameter. This constructor automatically sets the mechanical diameter equal to the optical diameter.

Arguments

  • radius: The radius of curvature of the surface.
  • diameter: The clear (optical) diameter of the surface.
source
BeamletOptics.SpotdetectorType
Spotdetector <: AbstractDetector

Simple 2D screen that stores the intersection point of incoming Beams. The intersection points are stored in local coordinates of the detector with respect to the screen origin.

Fields

  • shape: a 2D QuadraticFlatMesh that is aligned with the negative y-axis
  • data: stores the intersection points as Point2
  • hw: half-width of the detector plane

Additional information

Normal vector

Check the normal vector orientation of the detector plane if the spot diagram looks mirrored.

Reset behavior

Spot diagram data must be manually reset between traces via empty!

source
BeamletOptics.SpotdetectorMethod
Spotdetector(width)

Generates a quadratic rectangular 2D Spotdetector that is aligned with the negative y-axis. Refer to the type docs for more information.

Inputs:

  • width: edge length in [m]
source
BeamletOptics.StaticSystemType
StaticSystem <: AbstractSystem

A static container storing the optical elements of, i.e. a camera lens or lab setup. Compared to System this way defining the system is less flexible, i.e. no elements can be added or removed after construction but it allows for more performant ray-tracing.

Warning

This type uses long tuples for storing the elements. This container should not be used for very large optical systems as it puts a lot of stress onto the compiler.

Fields

  • objects: vector containing the different objects that are part of the system (subtypes of AbstractObject)
source
BeamletOptics.SystemType
System <: AbstractSystem

A container storing the optical elements of, i.e. a camera lens or lab setup.

Fields

  • objects: vector containing the different objects that are part of the system (subtypes of AbstractObject)
source
BeamletOptics.ThinBeamsplitterType
ThinBeamsplitter <: AbstractBeamsplitter

Represents a 2D beam-splitting device.

Fields

  • shape: 2D AbstractShape at which the splitting process occurs (e.g. a 2D-Mesh)
  • reflectance: scalar reflection factor
  • transmittance: scalar transmission factor
Warning

Note that the transmittance should be calculated from an input reflectance in order to ensure that R² + T² = 1.

source
BeamletOptics.ThinBeamsplitterMethod
ThinBeamsplitter(width, height; reflectance=0.5)

Creates a zero-thickness, lossless, non-polarizing 2D rectangular ThinBeamsplitter where

  • width: is the x-dir. edge length in [m]
  • height: is the z-dir. edge length in [m]
  • reflectance: kw-arg that determines how much light is reflected, i.e. 0.7 for a 70:30 splitter

Additional information

Reflectance

The input value for the reflectance R is normed such that R² + T² = 1, where T is the transmittance. The transmittance is calculated via T = √(1 - R²).

Reflection phase jump

Note that the reflection phase jump θᵣ is implemented by the individual interact3d-methods. Refer to them for more information.

source
BeamletOptics.UnionSDFType
UnionSDF{T, TT <: Tuple} <: AbstractSDF{T}

This SDF represents the merging of two or more SDFs. If the constituent SDFs do not overlap (they can and should touch) the resulting SDF should be still exact if the constituent SDFs are exact.

The intended way to construct these is not explicitely but by just adding two AbstractSDFs using the regular + operator.

s1 = SphereSDF(1.0)
translate3d!(s1, Point3(0, 1.0, 0.0))

s2 = SphereSDF(1.0)

# will result in a SDF with two spheres touching each other.
s_merged = s1 + s2
source
Base.lengthMethod
Base.length(ray::AbstractRay)

Returns the geometric length of a ray between its start and intersection point. If no intersection exists, Inf is returned.

Tip

Use optical_path_length to get the optical path length instead.

source
Base.lengthMethod
Base.length(beam::Beam)

Calculate the length of a beam up to the point of the last intersection.

Tip

Use optical_path_length to get the optical path length instead.

source
Base.positionMethod
position(object) -> Point3

Returns the current position of the object in R³ as a Point3 where (x, y, z) in a right-hand coordinate system.

In general, position(object) returns position(shape(object)) unless specified otherwise.

source
Base.positionMethod

Enforces that shape has to have the field pos or implement position().

source
BeamletOptics.BiConcaveLensSDFMethod
BiConcaveLensSDF(r1, r2, l, d=1inch)

Constructs a bi-concave lens SDF with:

  • r1 > 0: radius of concave front
  • r2 > 0: radius of convex back
  • l: lens thickness
  • d: lens diameter, default value is one inch
  • md: mechanical lens diameter, adds an outer ring section to the lens, if md > d.

The spherical surfaces are constructed flush with the cylinder surface.

source
BeamletOptics.BiConvexLensSDFMethod
BiConvexLensSDF(r1, r2, l, d=1inch)

Constructs a cylindrical bi-convex lens SDF with:

  • r1 > 0: radius of convex front
  • r2 > 0: radius of convex back
  • l: lens thickness
  • d: lens diameter, default value is one inch

The spherical surfaces are constructed flush with the cylinder surface.

source
BeamletOptics.CircularFlatMeshMethod
CircularFlatMesh(radius, n)

Creates a 2D rectangular Mesh that is centered around the origin and aligned with respect to the negative y-axis.

Inputs

  • radius: of the mesh in [m]
  • n: slice discretization factor (higher equals better resolution)
source
BeamletOptics.CuboidMeshMethod
CuboidMesh(x, y, z, θ=π/2)

Constructs the Mesh of a rectangular cuboid as per the dimensions specified by x, y and z. In addition, one side of the mesh can be tilted by an angle θ in order to generate the mesh of a rhomb. The mesh is initialized such that one corner of the cuboid lies at the origin.

Arguments

  • x, y, z: dimensions for the cube in [m]
  • θ: parallel tilt angle
source
BeamletOptics.MoellerTrumboreAlgorithmMethod
MoellerTrumboreAlgorithm(face::Matrix, ray::Ray)

A culling implementation of the Möller-Trumbore algorithm for ray-triangle-intersection. This algorithm evaluates the possible intersection between a ray and a face that is defined by three vertices. If no intersection occurs, Inf is returned. is the abort threshold for backfacing and non-intersecting triangles. is the threshold for negative values of t. This algorithm is fast due to multiple breakout conditions.

source
BeamletOptics.PlanoConcaveAsphericalLensSDFMethod
PlanoConcaveAsphericalLensSDF(r, l, d=1inch)

Constructs a plano-concave aspheric lens SDF with:

  • r > 0: front radius
  • l: lens thickness
  • d: lens diameter
  • cz: aspheric surface chip zone
  • k : The conic constant of the surface
  • α_coeffs : The (even) aspheric coefficients, starting with A4.
  • md: lens mechanical diameter (default: md = d)

The spherical surface is constructed flush with the cylinder surface.

source
BeamletOptics.PlanoConcaveLensSDFMethod
PlanoConcaveLensSDF(r, l, d=1inch)

Constructs a plano-concave lens SDF with:

  • r > 0: front radius
  • l: lens thickness
  • d: lens diameter, default value is one inch
  • md: mechanical lens diameter, must be > d

The spherical surface is constructed flush with the cylinder surface.

source
BeamletOptics.PlanoConvexAsphericalLensSDFMethod
PlanoConvexAsphericalLensSDF(r, l, d=1inch)

Constructs a plano-convex aspheric lens SDF with:

  • r > 0: front radius
  • l: lens thickness
  • d: lens diameter
  • k : The conic constant of the surface
  • α_coeffs : The (even) aspheric coefficients, starting with A4.

The spherical surface is constructed flush with the cylinder surface.

source
BeamletOptics.PlanoConvexLensSDFMethod
PlanoConvexLensSDF(r, l, d=1inch)

Constructs a plano-convex lens SDF with:

  • r > 0: front radius
  • l: lens thickness
  • d: lens diameter, default value is one inch

The spherical surface is constructed flush with the cylinder surface.

source
BeamletOptics.RectangularCompensatorPlateMethod
RectangularCompensatorPlate(width, height, thickness, n)

Creates a compensator plate (modeled as a Prism) that can be used to remove parallel beam offsets created by e.g. the RectangularPlateBeamsplitter. The compensator is aligned with the positive y-axis. The first surface lies at the origin.

Inputs

  • width: compensator width along the x-axis in [m]
  • height: compensator height along the z-axis in [m]
  • thickness: compensator thickness along the y-axis in [m]
  • n: the RefractiveIndex of the substrate
source
BeamletOptics.RectangularFlatMeshMethod
RectangularFlatMesh(width, height)

Creates a 2D rectangular Mesh that is centered around the origin and aligned with respect to the y-axis. Vertex normals are parallel to the positive y-axis.

Inputs

  • width: width along the x-axis in [m]
  • height: height along the z-axis in [m]
source
BeamletOptics.RectangularPlanoMirrorMethod
RectangularPlanoMirror(width, height, thickness)

Constructs a rectangular plano Mirror based on the input dimensions. The front reflecting surface is normal to the y-axis and lies at the origin.

Inputs

  • width: of the mirror in x-direction [m]
  • height: of the mirror in z-direction [m]
  • thickness: of the mirror in y-direction [m]
source
BeamletOptics.RetroMeshMethod
RetroMesh(scale::Real; T = Float64)

Creates an open tetrahedral Mesh with edges derived from the vertices of a unit cube. Can be scaled with a scale factor. The data type for the vertices and internal computations can be adjusted using T (default: Float64).

source
BeamletOptics.SphericalDoubletLensMethod
SphericalDoubletLens(r1, r2, r3, l1, l2, d, n1, n2)

Generates a two-component "cemented" doublet lens consisting of two spherical lenses. For radii sign definition, refer to the SphericalLens constructor.

Arguments

  • r1: radius of curvature for first surface
  • r2: radius of curvature for second (cemented) surface
  • r3: radius of curvature for third surface
  • l1: first lens thickness
  • l2: second lens thickness
  • d: lens diameter
  • n1: first lens RefractiveIndex
  • n1: second lens RefractiveIndex
source
BeamletOptics.SphericalLensFunction
SphericalLens(r1, r2, l, d=1inch, n=λ->1.5)

Creates a spherical Lens based on:

  • r1: front radius
  • r2: back radius
  • l: lens thickness
  • d: lens diameter, default is one inch
  • n: RefractiveIndex as a function of λ, i.e. n = n(λ)

Notes

Radius of curvature (ROC) sign

The ROC is defined to be positive if the center is to the right of the surface. Otherwise it is negative.

Thin lenses

If l is set to zero, a ThinLens will be created. However, note that the actual lens thickness will be different from zero.

source
BeamletOptics.SquarePlanoMirrorMethod
SquarePlanoMirror(width, thickness)

Constructs a square plano Mirror with equal width and height. The front reflecting surface is normal to the y-axis and lies at the origin. See also RectangularPlanoMirror.

Inputs

  • width: the side length of the square mirror in x- and y-direction [m]
  • thickness: of the mirror in [m]
source
BeamletOptics.SquarePlanoMirror2DMethod
SquarePlanoMirror2D(edge_length)

Constructs a 2D square plano Mirror with a given edge_length. The reflecting surface is normal to the y-axis.

Inputs

  • edge_length: the edge length of the square mirror in [m]
source
BeamletOptics.ThinLensSDFMethod
ThinLensSDF(r1, r2, d=1inch)

Constructs a bi-convex thin lens SDF-based shape with:

  • r1 > 0: radius of convex front
  • r2 > 0: radius of convex back
  • d: lens diameter, default value is one inch

The spherical surfaces are constructed flush.

source
BeamletOptics.UniformDiscSourceMethod
UniformDiscSource(pos, dir, diameter, λ; num_rays=1_000)

Generates a ray fan with equal area per ray across a circular pupil using the deterministic sunflower (Fibonacci) pattern.

Note

This is merely a CollimatedSource constructor which uses Fibonacci sampling instead of a linear grid.

Arguments

The following inputs and arguments can be used to configure the underlying CollimatedSource:

Inputs

  • pos: center beam starting position
  • dir: center beam starting direction
  • diameter: outer beam bundle diameter in [m]
  • λ = 1e-6: wavelength in [m]

Keyword Arguments

  • num_rays=1000: total number of rays in the source
source
BeamletOptics._calculate_global_E0Method
_calculate_global_E0(in_dir, out_dir, J, E0)

Calculates the resulting polarization vector as per the publication by Yun et al. for each surface interaction. If the in- and out-directions of propagation are parallel, an arbitrary basis is chosen for the s- and p-components.

Arguments

  • in_dir: propagation direction before surface interaction
  • out_dir: propagation direction after surface interaction
  • J: Jones matrix extended to 3x3, e.g. [-rₛ 0 0; 0 rₚ 0; 0 0 1] for reflection
  • E0: Polarization vector before surface interaction
source
BeamletOptics._raymarch_insideMethod
_raymarch_inside(object::AbstractSDF, pos, dir; num_iter=1000, dl=0.1)

Perform the ray marching algorithm if the starting pos is inside of object.

source
BeamletOptics._raymarch_outsideMethod
_raymarch_outside(shape::AbstractSDF, pos, dir; num_iter=1000, eps=1e-10)

Perform the ray marching algorithm if the starting pos is outside of shape.

source
BeamletOptics._world_to_sdfMethod
_world_to_sdf(sdf, point)

Transforms the coordinates of point into a reference frame where the sdf lies at the origin. Useful to represent translation and rotation. If rotations are applied, the rotation is applied around the local sdf coordinate system.

source
BeamletOptics.align3d!Method
align3d!(shape, target_axis)

Rotates the shape such that its local y-axis aligns with the target_axis.

source
BeamletOptics.align3dMethod
align3d(start::AbstractVector, target::AbstractVector)

Returns the rotation matrix R that will align the start vector to be parallel to the target vector. Based on 'Avoiding Trigonometry' by Íñigo Quílez. The resulting matrix was transposed due to column/row major issues. Vector length is maintained. This function is very fast.

source
BeamletOptics.angle3dFunction
angle3d(ray::AbstractRay, intersect::Intersection=intersection(ray))

Calculates the angle between a ray and its or some other intersection.

source
BeamletOptics.angle3dMethod
angle3d(target::AbstractVector, reference::AbstractVector)

Returns the angle between the target and reference vector in rad.

source
BeamletOptics.aspheric_equationMethod

asphericequation(r, c, k, αcoeffs)

The aspheric surface equation. The asphere is defined by:

  • c : The curvature (1/radius) of the surface
  • k : The conic constant of the surface
  • α_coeffs : The (even) aspheric coefficients, starting with A4.

This function returns NaN if the square root argument becomes negative.

Note

Only even aspheres are implemented at the moment. This will change soon.

source
BeamletOptics.base_transformFunction
base_transform(base, base2=I(3))

Return the base transformation matrix for transforming from vectors given relative to base2 into base.

source
BeamletOptics.bounding_sphereMethod
bounding_sphere(sdf)

Returns nothing or a center point and the radius of a sphere which encloses the shape of the SDF. This function is currently only used for rendering SDFs but might be used in the future to optimize the raymarching algorithm, by tracing against the bounding sphere of and SDF first, instead of calling the more costly complex SDF.

source
BeamletOptics.calc_local_limsMethod
calc_local_lims(psf; crop_factor=1, center=:centroid)

Compute a symmetric [xmin,xmax]×[zmin,zmax] box around the PSF’s weighted centroid.

• If center==:centroid (the default), uses x0 = ∑ wᵢ·xᵢ / ∑ wᵢ, z0 = ∑ wᵢ·yᵢ / ∑ wᵢ with wᵢ = projection_factor. • If center==:bbox, falls back to the midpoint of [min,max].

Returns (x_min, x_max, z_min, z_max).

source
BeamletOptics.children!Method
children!(beam::B, child::B) where {B<:AbstractBeam}

Handles the inclusion of adding a single child to an existing beam. The function behaves as follows:

  1. If no previous children exist, add child
  2. If beam already has a single child, modify child beam starting ray (retracing)
  3. Else throw error
source
BeamletOptics.convex_aspheric_surface_distanceMethod
convex_aspheric_surface_distance(r, z, c, k, d, α_coeffs)

Calculates the 2D distance field for an aspheric surface at radius r away from the optical axis position z. The asphere is defined by:

  • c : The curvature (1/radius) of the surface
  • k : The conic constant of the surface
  • d : The diameter of the asphere
  • α_coeffs : The (even) aspheric coefficients, starting with A2.

Note that this is not just an infinite aspheric surface and also not a surface segment but a closed 2D perimeter.

It is intended to pair the SDF derived from this distance field with a cylinder SDF to build a real lens.

source
BeamletOptics.electric_fieldFunction
electric_field(I::Real, Z=Z_vacuum, ϕ=0)

Calculates the E-field phasor in [V/m] for a given intensity I and phase ϕ. Vacuum wave impedance is assumed.

source
BeamletOptics.electric_fieldMethod
electric_field(gauss::GaussianBeamlet, r, z)

Calculates the electric field phasor [V/m] of gauss at the radial and longitudinal positions r and z.

source
BeamletOptics.electric_fieldMethod
electric_field(r, z, E0, w0, w, k, ψ, R) -> ComplexF64

Computes the analytical complex electric field distribution of a stigmatic TEM₀₀ Gaussian beam which is described by:

\[E(r,z) = {E_0}\frac{{{w_0}}}{{w(z)}}\exp\left( { - \frac{{{r^2}}}{{w{{(z)}^2}}}} \right)\exp\left(i\left[ {kz + \psi + \frac{{k{r^2}}}{2 R(z)}} \right] \right)\]

Arguments

  • r: radial distance from beam origin
  • z: axial distance from beam origin
  • E0: peak electric field amplitude
  • w0: waist radius
  • w: local beam radius
  • k: wave number, equal to 2π/λ
  • ψ: Gouy phase shift (defined as $-\text{atan}\left(\frac{z}{z_r}\right)$ !)
  • R: wavefront curvature, i.e. 1/r (radius of curvature)
source
BeamletOptics.fresnel_coefficientsMethod

fresnel_coefficients(θ, n)

Calculates the complex Fresnel coefficients for reflection and transmission based on the incident angle θ in [rad] and the refractive index ratio n = n₂ / n₁. Returns rₛ, rₚ, tₛ and tₚ.

Signs

Info

The signs of rₛ, rₚ are based on the definition by Fowles (1975, 2nd Ed. p. 44) and Peatross (2015, 2023 Ed. p. 78)

source
BeamletOptics.gauss_parametersMethod
gauss_parameters(gauss::GaussianBeamlet, z; hint::Union{Nothing, Tuple{Int, Vector{<:Real}}}=nothing)

Calculate the local waist radius and Gouy phase of an unastigmatic Gaussian beamlet at a specific distance z based on the method of J. Arnaud (1985) and D. DeJager (1992).

Arguments

  • gauss: the GaussianBeamlet object for which parameters are to be calculated.
  • z: the position along the beam at which to calculate the parameters.
  • hint: an optional hint parameter for the relevant point/index of the appropriate beam segment. If not provided, the function will automatically select the ray.

Returns

  • w: local radius
  • R: curvature, i.e. 1/r where r is the radius of curvature
  • ψ: Gouy phase (note that -atan definition is used)
  • w0: local beam waist radius
source
BeamletOptics.gauss_parametersMethod
gauss_parameters(gauss::GaussianBeamlet, zs::AbstractArray)

Return the parameters of the GaussianBeamlet along the specified positions in zs.

source
BeamletOptics.intensityFunction

Calculates the intensity in [W/m²] for a given complex electric field phasor E. Vacuum wave impedance is assumed.

source
BeamletOptics.intensityMethod
intensity(psf::PSFDetector{T};
          n::Int=100,
          crop_factor::Real=1,
          center::Symbol=:centroid,
          x_min = Inf,
          x_max = Inf,
          z_min = Inf,
          z_max = Inf,
          x0_shift::Real=0,
          z0_shift::Real=0) where T

Compute the two‐dimensional point‐spread function (PSF) of an optical system as captured by a PSFDetector. The returned intensity map is sampled on a regular n×n grid in the detector’s local (x,z)-plane.

Keyword Arguments

  • n::Int=100 Number of sample points per axis.
  • crop_factor::Real=1 Scales the half‐width of the sampling window returned by calc_local_lims; values >1 expand, <1 shrink.
  • center::Symbol=:centroid How the sampling window is centred. :centroid uses the projection‑weighted centroid, :bbox uses the geometric mid‑point of the bounding box.
  • x_min, x_max, z_min, z_max Manually override the sampling bounds in the local x or z directions. If left as Inf, the bounds from calc_local_lims are used.
  • x0_shift::Real=0, z0_shift::Real=0 Apply a constant offset to the entire x or z coordinate arrays, useful for recentring or testing alignment.

Returns

A tuple (xs, zs, I) where

  • xs::LinRange{T} and zs::LinRange{T} are the sampled coordinates in the detector’s local x and z axes,
  • I::Matrix{T} is the corresponding raw/unscaled intensity map
Resetting detectors

Be sure to call empty!(psf) before each new measurement if reusing the same detector.

Scaling

The returned values are raw/unscaled and not a Strehl ratio. This feature is not yet added. In future versions a pupil finder along with a Strehl estimator will be added.

source
BeamletOptics.interact3dMethod
interact3d(::AbstractSystem, pd::Photodetector, gauss::GaussianBeamlet, ray_id::Int)

Implements the Photodetector interaction with a GaussianBeamlet. On hit, the scalar E-field of the gauss is added to the current PD field matrix. Tilt and tip between beam and PD surface are considered via projection factors.

source
BeamletOptics.interact3dMethod
interact3d(::AbstractSystem, bs::ThinBeamsplitter, gauss::GaussianBeamlet, ray_id::Int)

Models the interaction between a ThinBeamsplitter and a GaussianBeamlet.

Reflection phase jump

The reflection phase jump is modeled here as θᵣ = π for simplicity. This is since in practice it will have only a relative effect on the signal at the detector for interferometric setups. The phase jump is applied to the reflected portion of any incoming beam that faces the ThinBeamsplitter normal vector, which assumes that the splitter has an unambigous normal, i.e. a 2D mesh. This is intended to model the effect of the Fresnel equations without full polarization calculus.

source
BeamletOptics.interact3dMethod
interact3d(::AbstractSystem, object::AbstractObject, ::AbstractBeam)

Defines the optical interaction between an incoming/outgoing beam/ray of light and an optical element, must return an AbstractInteraction or nothing. The default behavior is that no interaction occurs, i.e. return of nothing, which should stop the system tracing procedure. Refer to the AbstractInteraction typedocs for more information on the return type value.

source
BeamletOptics.interact3dMethod
interact3d(system::AbstractSystem, object::AbstractObject, gauss::GaussianBeamlet{R}, ray_id::Int)

Generic dispatch for the interact3d method of a GaussianBeamlet with an AbstractObject. Unless a more concrete implementation exists, the interaction of the Gaussian is assumed to be the interaction of the chief, waist and divergence rays with an object.

Returns

The interact3d method for the GaussianBeamlet must return a GaussianBeamletInteraction.

source
BeamletOptics.interact3dMethod
interact3d(AbstractReflectiveOptic, PolarizedRay)

Implements the ideal reflection of a PolarizedRay via the normal at the intersection point on an optical surface. A Jones matrix of [-1 0 0; 0 1 0] is assumed as per Peatross (2015, 2023 Ed. p. 154) and Yun et al. (see PolarizedRay for more information).

source
BeamletOptics.interact3dMethod
interact3d(AbstractReflectiveOptic, Ray)

Implements the reflection of a Ray via the normal at the intersection point on an optical surface.

source
BeamletOptics.interact3dMethod
interact3d(AbstractSystem, AbstractRefractiveOptic, Beam, PolarizedRay)

Implements the refraction of a PolarizedRay at an uncoated optical surface. The "outside" ref. index is obtained from the system unless specified otherwise. Reflection and transmission values are calculated via the fresnel_coefficients. Stray light is not tracked. In the case of total internal reflection, only the reflected light is traced.

source
BeamletOptics.interact3dMethod
interact3d(AbstractSystem, AbstractRefractiveOptic, Beam, Ray)

Implements the refraction of a Ray at an optical surface. The "outside" ref. index is obtained from the system unless specified otherwise. At the critical angle, total internal reflection occurs (see refraction3d).

source
BeamletOptics.interact3dMethod
interact3d(::AbstractSystem, psf::PSFDetector, beam::Beam{T, Ray{T}}, ray::Ray{T}) where T

Implements the PSFDetector interaction with a Beam. On hit, the hit position, direction, optical path length, wavelength and projection factor are captured and stored with the data field of the detector.

source
BeamletOptics.intersect3dMethod
intersect3d(shape::AbstractShape, ::AbstractRay)

Defines the intersection between an AbstractShape and an AbstractRay, must return an Intersection or nothing. The default behavior for concrete shapes and rays is to indicate no intersection, that is nothing, which will inform the tracing algorithm to stop. Refer to the Intersection documentation for more information on the return type value.

source
BeamletOptics.intersect3dMethod
intersect3d(mesh::Mesh, ray::Ray)

This function is a generic implementation to check if a ray intersects the shape mesh.

source
BeamletOptics.intersect3dMethod
intersect3d(plane_position, plane_normal, ray)

Returns the intersection between a ray and an infinitely large plane which is characterized by its position and normal.

source
BeamletOptics.isenteringMethod
isentering(ray)

Tests whether the ray is entering a shape based on the orientation of the ray direction and surface normal. If no intersection is present, default behavior is to return false.

source
BeamletOptics.isinfrontofMethod
isinfrontof(point::AbstractVector, pos::AbstractVector, dir::AbstractVector)

Tests if a point is in front of the plane defined by the position and direction vectors.

source
BeamletOptics.isinfrontofMethod
isinfrontof(shape::AbstractShape, ray::AbstractRay)

A simple test to check if a shape lies "in front of" a ray. The forward direction is here defined as the ray orientation. Only works well if ray is outside of the volume of shape. Can be dispatched to return more accurate results for subtypes of AbstractShape.

source
BeamletOptics.isparaxialFunction
isparaxial(system, gb::GaussianBeamlet, threshold=π/4)

Tests the angle between the waist and divergence beams and refractive surfaces. A target threshold of π/4 or 45° is assumed before abberations become dominant.

source
BeamletOptics.istiltedMethod
istilted(system::System, gb::GaussianBeamlet)

Tests if refractive elements are tilted with respect to the beamlet optical axis, i.e. introduce simple astigmatism.

source
BeamletOptics.lensmakers_eqMethod
lensmakers_eq(R1, R2, n)

Calculates the thin lens focal length based on the radius of curvature R1/R2 and the lens refractive index n. If center of sphere is on left then R < 0. If center of sphere is on right then R > 0.

source
BeamletOptics.line_plane_distance3dMethod
line_plane_distance3d(plane_position, plane_normal, line_position, line_direction)

Returns the distance between a line and an infinitely large plane which are characterized by their position and normal/direction.

source
BeamletOptics.line_point_distance3dMethod
line_point_distance3d(pos, dir, point)

Computes the shortes distance between a line described by pos+t*dir and a point in 3D. This function is slow and should be used only for debugging purposes.

source
BeamletOptics.list_subtypesFunction
list_subtypes(T::Type; max_depth::Int=5)

Prints a tree of all subtypes, e.g. list_subtypes(AbstractObject). Maximum exploration depth can be limited by passing max_depth. Returns the total number of types encountered.

source
BeamletOptics.mechanical_diameterMethod
mechanical_diameter(s::AbstractRotationallySymmetricSurface)

Returns the mechanical diameter of the surface.

Note

It is assumed that mechanical_diameter(s) >= diameter(s) always holds.

source
BeamletOptics.normal3dMethod
normal3d(target, reference)

Returns a vector with unit length that is perpendicular to the target and an additional reference vector. Vector orientation is determined according to right-hand rule.

source
BeamletOptics.normal3dMethod
normal3d(mesh::AbstractMesh, fID::Int)

Returns a vector with unit length that is perpendicular to the target face` according to the right-hand rule. The vertices must be listed row-wise within the face matrix.

source
BeamletOptics.objectsMethod
objects(system::System)

Exposes all objects stored within the system. By exposing the Leaves of the tree only, it is ensured that AbstractObjectGroups are flattened into a regular vector.

source
BeamletOptics.op_extrude_xMethod
op_extrude_x(p, sdf2d::Function, height)

Calculates the SDF at point p for the given 2D-SDF function and extrudes the shape to height along the x-axis.

source
BeamletOptics.op_extrude_zMethod
op_extrude_z(p, sdf2d::Function, height)

Calculates the SDF at point p for the given 2D-SDF function and extrudes the shape to height along the z-axis.

source
BeamletOptics.op_revolve_yMethod
op_revolve_y(p, sdf2d::Function, offset)

Calculates the SDF at point p for the given 2D-SDF function with offset by revolving the 2D shape around the y-axis.

source
BeamletOptics.op_revolve_zMethod
op_revolve_z(p, sdf2d::Function, offset)

Calculates the SDF at point p for the given 2D-SDF function with offset by revolving the 2D shape around the z-axis.

source
BeamletOptics.orientationMethod
orientation(object) -> Matrix

Returns the current orientation of the object in R³ as a matrix. The matrix represents the local fixed-body coordinate system.

In general, orientation(object) returns orientation(shape(object)) unless specified otherwise.

source
BeamletOptics.parent!Method
parent!(beam::GaussianBeamlet, parent::GaussianBeamlet)

Ensures that the GaussianBeamlet knows about its parent beam. In addition, links the chief beams of child and parent. Important for correct functioning of point_on_beam and length.

source
BeamletOptics.point_on_beamMethod
point_on_beam(beam::Beam, t::Real)

Function to find a point given a specific distance t along the beam. Return the ray index aswell. For negative distances, assume first ray backwards.

source
BeamletOptics.radiusMethod
radius(s::AbstractRotationallySymmetricSurface)

Returns the radius of curvature of the surface. This might return Inf for planar surfaces or surfaces which cannot be described by just one curvature radius.

source
BeamletOptics.rayleigh_rangeMethod
rayleigh_range(g::GaussianBeamlet; M2=1)

Returns the Rayleigh range for the first beam section of the GaussianBeamlet g. Note: M2 is not stored in g during construction and must be specified by the user.

source
BeamletOptics.reflection3dMethod
reflection3d(dir, normal)

Calculates the reflection between an input vector dir and surface normal vector in R³. Vectors dir and normal must have unit length!

source
BeamletOptics.refraction3dMethod
refraction3d(dir, normal, n1, n2)

Calculates the refraction between an input vector dir and surface normal vector in R³. n1 is the "outside" refractive index and n2 is the "inside" refractive index. The function returns the new direction of propagation and a boolean flag to indicate if internal refraction has occured.

Vectors dir and normal must have unit length!

Total internal reflection

If the critical angle for n1, n2 and the incident angle is reached, the ray is reflected internally instead!

Arguments

  • dir: direction vector of incoming ray
  • normal: surface normal at point of intersection
  • n1: index of ref. before refraction
  • n2: index of ref. after refraction
source
BeamletOptics.render!Method
render!(axis, thing; kwargs...)

The render! function allows for the visualization of optical system and beams under the condition that a suitable backend is loaded. This means that either one of the following packages must be loaded in combination with BeamletOptics via using:

  1. GLMakie
    • preferred for 3D viewing
    • use LScene or Axis3 environments
  2. CairoMakie
    • preferred for the generation of high-quality .pngs
    • only Axis3 is supported

If no suitable backend is loaded, a MissingBackendError will be thrown.

Implementations reqs.

All concrete implementations of render! must adhere to the following minimal interface:

render!(axis, thing; kwargs...)

  • axis: an axis type of the union of LScene or Axis3
  • thing: an abstract or concrete object or beam type
  • kwargs: custom or Makie keyword arguments that are passed to the underlying backend

Refer to the BeamletOptics extension docs for Makie for more information.

source
BeamletOptics.reset_rotation3d!Method
reset_rotation3d!(::MultiShape, object)

Reset all applied rotations of the object, i.e. resets the local coordinate system to the standard base.

Parts within parts

Sub-part relative rotations are not reset!

source
BeamletOptics.reset_translation3d!Method
reset_translation3d!(::MultiShape, object)

Resets all applied translations of the object, i.e. moves the center back to the origin.

Parts within parts

Sub-part relative translations are not reset!

source
BeamletOptics.retrace_system!Method
retrace_system!(system, beam)

This function tries to reuse data from a previous solution in order to solve the system using a sequential approach.

Retracing

The retracing logic for an already solved beam loops over the rays and children and is as follows:

Begin

  1. Test if current ray has a valid intersection
    • If not, mark beam tail for cleanup and go to End
  2. Recalculate the intersection
    • If a hint was provided by a previous interaction, use hinted object
    • Else, test against previous intersection
  3. Test if the ray still has a valid intersection after recalculation
    • If no object is hit, mark beam tail for cleanup and go to End

Interact

  1. Recalculate the optical interaction
    • Catch hints provided for next ray
    • If no interaction occurs, mark beam tail for conditional cleanup and go to End
  2. Add the interaction to the current beam
    • If another ray follows, modify the next starting position - Go to Begin
    • Else mark children for cleanup, push new ray to beam tail - Go to End

End

  1. If cleanup is required, do conditionally
    • remove all beam tail rays after current ray
    • remove all beam children
    • reset beam tail ray intersection to nothing
source
BeamletOptics.retrace_system!Method
retrace_system!(system::System, gauss::GaussianBeamlet{T}) where {T <: Real}

Retrace the beam stored in GaussianBeamlet through the optical system. Chief, waist and divergence ray intersections and interactions are recalculated. All rays must hit the same object, or the retracing step is aborted. If retracing is stopped before the end of the beam is reached, further rays are dropped.

source
BeamletOptics.rotate3d!Method
rotate3d!(mesh, axis, θ)

Mutating function that rotates the mesh around the specified rotation axis by the angle θ.

source
BeamletOptics.rotate3d!Method
rotate3d!(shape::AbstractShape, axis, θ)

Rotates the dir-matrix of shape around the reference-axis by an angle of θ.

source
BeamletOptics.rotate3dMethod
rotate3d(reference::Vector, θ)

Returns the rotation matrix that will rotate a vector around the reference axis at an angle θ in radians. Vector length is maintained. Rotation in clockwise direction?

source
BeamletOptics.sagMethod
sag(r::Real, l::Real)

Calculates the sag of a cut circle with radius r and chord length l

source
BeamletOptics.sdfMethod
sdf(::AbstractRotationallySymmetricSurface, ::Union{Nothing, AbstractOrientationType})

Takes the surface specification and an optional AbstractOrientationType as trait parameter and returns a corresponding AbstractSDF type.

Surface vs. volume based tracing

This function is a mere convenience provider for users coming from other optic simulations frameworks which are surface oriented. The goal of this function is to return the best matching closed volume SDF which posesses a surface with the given specs on one side and most often a boundary and planar surface on the other side.

Info

Always keep in mind that this package performs closed-volume baced ray tracing using either SDFs or meshes.

source
BeamletOptics.set_new_origin3d!Method
set_new_origin3d!(mesh::AbstractMesh)

Resets the mesh directional matrix and position vector to their initial values.
Warning: this operation is non-reversible!

source
BeamletOptics.solve_system!Method
solve_system!(system::System, beam::AbstractBeam; r_max=100, retrace=true)

Manage the tracing of an AbstractBeam through an optical system. The function retraces the beam if possible and then proceeds to trace each leaf of the beam tree through the system. The condition to stop ray tracing is that the last beam intersection is nothing or the beam interaction is nothing. Then, the system is considered to be solved. A maximum number of rays per beam (r_max) can be specified in order to avoid infinite calculations under resonant conditions, i.e. two facing mirrors.

Arguments

  • system::System: The optical system in which the beam will be traced.
  • beam::AbstractBeam: The beam object to be traced through the system.
  • r_max::Int=100 (optional): Maximum number of tracing iterations for each leaf. Default is 100.
  • retrace::Bool=true (optional): Flag to indicate if the system should be retraced. Default is true.
source
BeamletOptics.trace_system!Method
trace_system!(system::AbstractSystem, beam::Beam{T}; r_max::Int = 20) where {T <: Real}

Trace a Beam through an optical system. Maximum number of tracing steps can be capped by r_max.

Tracing logic

The intersection of the last ray of the beam with any objects in the system is tested. If an object is hit, the optical interaction is analyzed and tracing continues. Else the tracing procedure is stopped.

Arguments

  • system:: The optical system through which the Beam is traced.
  • beam: The Beam object to be traced.
  • r_max: Maximum number of tracing iterations. Default is 20.
source
BeamletOptics.trace_system!Method
trace_system!(system::System, gauss::GaussianBeamlet{T}; r_max::Int = 20) where {T <: Real}

Trace a GaussianBeamlet through an optical system. Maximum number of tracing steps can be capped by r_max.

Tracing logic

The chief, waist and divergence beams are traced step-by-step through the system. For each intersection after a tracing_step!, the intersections are compared. If all rays hit the same target, the optical interaction is analyzed, else the tracing stops.

Arguments

  • system: The optical system through which the GaussianBeamlet is traced.
  • gauss: The GaussianBeamlet object to be traced.
  • r_max: Maximum number of tracing iterations. Default is 20.
source
BeamletOptics.tracing_step!Method
tracing_step!(system::AbstractSystem, ray::AbstractRay{R}, hint::Hint)

Tests if the ray intersects an object in the optical system. Returns the closest intersection.

Hint

An optional Hint can be provided to test against a specific object (and shape) in the system first.

Warning

If a hint is provided and the object intersection is valid, the intersection will be returned immediately. However, it is not guaranteed that this is the true closest intersection.

source
BeamletOptics.translate3d!Method
translate3d!(mesh::AbstractMesh, offset)

Mutating function that translates the vertices of an mesh in relation to the offset vector. In addition, the mesh position vector is overwritten to reflect the new "center of gravity".

source
BeamletOptics.translate_to3d!Method
translate_to3d!(::MultiShape, object, target)

Translates all parts of the MultiShape object in parallel to the specified target position. The object center point will be equal to the target.

source