Complete Reference
This page contains what should be a complete list of all docstrings in the OpticSim module, and its submodule.
Index
OpticSim.Chebyshev
OpticSim.Examples
OpticSim.QType
OpticSim.Zernike
OpticSim.AbstractRayDirectionDistribution
OpticSim.AbstractSpectrum
OpticSim.BSplineCurve
OpticSim.BeamState
OpticSim.BezierCurve
OpticSim.CircularStopShape
OpticSim.CurveType
OpticSim.GridSagInterpolation
OpticSim.KnotVector
OpticSim.MeasuredSpectrum
OpticSim.Primitive
OpticSim.RayState
OpticSim.RectangularStopShape
OpticSim.Source
OpticSim.Spline
OpticSim.SplineSurface
OpticSim.StopSurface
OpticSim.UniformSpectrum
OpticSim.ZernikeIndexType
OpticSim.Annulus
OpticSim.ArizonaEye
OpticSim.AsphericLens
OpticSim.BoundedCylinder
OpticSim.Chebyshev.T
OpticSim.Chebyshev.U
OpticSim.Chebyshev.dTdq
OpticSim.Circle
OpticSim.CircularAperture
OpticSim.CircularAperture
OpticSim.ConicLens
OpticSim.Cuboid
OpticSim.FresnelLens
OpticSim.GridField
OpticSim.GridField
OpticSim.HexagonalPrism
OpticSim.HexapolarField
OpticSim.HexapolarField
OpticSim.ModelEye
OpticSim.QType.S
OpticSim.QType.S0
OpticSim.QType.dS0dx
OpticSim.QType.dSdx
OpticSim.RectangularAperture
OpticSim.RectangularAperture
OpticSim.RectangularPrism
OpticSim.SphericalLens
OpticSim.Spider
OpticSim.TriangularPrism
OpticSim.Zernike.NolltoNM
OpticSim.Zernike.OSAtoNM
OpticSim.Zernike.R
OpticSim.Zernike.normalisation
OpticSim.Zernike.δζ
OpticSim.Zernike.ζ
OpticSim.assembly
OpticSim.closestintersection
OpticSim.closestpointonray
OpticSim.csgdifference
OpticSim.csgintersection
OpticSim.csgunion
OpticSim.curvedimension
OpticSim.curveorder
OpticSim.detectorimage
OpticSim.direction
OpticSim.distance
OpticSim.doesintersect
OpticSim.evaluatecurve
OpticSim.extractmovinglines
OpticSim.fresnel
OpticSim.gendirection
OpticSim.generateray
OpticSim.generateray
OpticSim.generateray
OpticSim.generateray
OpticSim.generateray
OpticSim.generateray
OpticSim.genorigin
OpticSim.halfspaceintersection
OpticSim.identitytransform
OpticSim.interface
OpticSim.intersections
OpticSim.isemptyinterval
OpticSim.isinfinity
OpticSim.ispositivehalfspace
OpticSim.israyorigin
OpticSim.israyorigininterval
OpticSim.jacobian
OpticSim.leaf
OpticSim.leaf
OpticSim.linedimension
OpticSim.makemesh
OpticSim.matricesforeigen
OpticSim.newton
OpticSim.normal
OpticSim.numberoflines
OpticSim.orthogonalitymatrix
OpticSim.point
OpticSim.pressure
OpticSim.processintersection
OpticSim.randinsolidangle
OpticSim.raysample
OpticSim.reset!
OpticSim.resetdetector!
OpticSim.reversenormal
OpticSim.rotation
OpticSim.rotationd
OpticSim.rotmat
OpticSim.rotmatbetween
OpticSim.rotmatd
OpticSim.samplesurface
OpticSim.semidiameter
OpticSim.snell
OpticSim.spectrumpower
OpticSim.spectrumsample
OpticSim.sum!
OpticSim.surfaceintersection
OpticSim.temperature
OpticSim.trace
OpticSim.trace
OpticSim.trace
OpticSim.traceMT
OpticSim.tracehits
OpticSim.tracehitsMT
OpticSim.translation
OpticSim.triangulate
OpticSim.triangulatedintersection
OpticSim.uv
OpticSim.uvrange
OpticSim.uvtopix
OpticSim.α
OpticSim
OpticSim.AbstractRayDirectionDistribution
— TypeEvery ray direction distribution must implement the functions direction(a::Source)
OpticSim.AbstractSpectrum
— TypeEach AbstractSpectrum type defines a spectrumSample function which returns a uniformly sampled point from the spectrum of the light source and the power at that wavelength.
OpticSim.BSplineCurve
— TypeBSplineCurve{P,S,N,M} <: Spline{P,S,N,M}
N
is the spatial dimension of the curve. M
is the curve order, i.e., the highest power of the parameterizing variable, u
. All curve segments are assumed to be of the same order.
BSplineCurve{P,S,N,M}(knots::KnotVector{S}, controlpoints::AbstractArray{MVector{N,S},1})
OpticSim.BeamState
— TypeConvergingBeam
, DivergingBeam
or CollimatedBeam
, defines the behavior of a beam in a HologramInterface
.
OpticSim.BezierCurve
— TypeBezierCurve{P,S,N,M} <: Spline{P,S,N,M}
N
is the dimension of the curve, M
is the curve order
BezierCurve{P,S,N,M}(controlpoints::AbstractArray{<:AbstractArray{S,1},1})
OpticSim.CircularStopShape
— TypeCircularStopShape <: StopShape
OpticSim.CurveType
— TypeEither Rational
or Euclidean
, used for Spline
s and SplineSurface
s.
OpticSim.GridSagInterpolation
— TypeEither GridSagLinear
or GridSagBicubic
- determines the interpolation between sample points in the grid for a GridSagSurface
.
OpticSim.KnotVector
— TypeKnotVector{T<:Number}
Vector to define knots used for BSplineCurve
and BSplineSurface
.
OpticSim.MeasuredSpectrum
— TypeUse measured spectrum to compute emitter power. Create spectrum by reading CSV files.
Evaluate spectrum at arbitrary wavelength with spectrumpower
OpticSim.Primitive
— TypePrimitive{T<:Real}
T
is the number type used to represent the primitive, e.g., Float64
. Primitives are the basic elements which can be stored in bounding volume hierarchies and include surfaces and CSG objects
Must implement the following:
boundingbox(a::Primitive{T})::BoundingBox{T}
centroid(a::Primitive{T})::SVector{3,T}
OpticSim.RayState
— Typestores the state of ray generation. Keeps track of ray number so you can have a PointOrigin which emits multiple rays. This can be used in the future to do phase summations on the image plane, since only point sources are perfectly spatially coherent.
OpticSim.RectangularStopShape
— TypeRectangularStopShape <: StopShape
OpticSim.Source
— TypeCreate an emitter by choosing a type for the spectral distribution, ray origin distribution, ray direction distribution, and angular light power distribution. For example, to create an emitter with uniform spectrum, rays originating on a regular grid, ray directions distributed around a vector v up to a maximum angle θmax, with Lambertian power distribution call the Source construction with these types and arguments:
a = Source{UniformSpectrum,GridOrigin, ConeDistribution, Lambertian}(UniformSpectrum(),GridOrigin(),ConeDistribution(θmax),Lambertian())
OpticSim.Spline
— TypeSpline{P<:CurveType,S<:Number,N,M}
M
is the curve order, i.e., the highest power of the parameterizing variable, u. P
determines the CurveType
.
All Spline types must implement:
point(curve,u)
and have field controlpolygon
OpticSim.SplineSurface
— TypeSplineSurface{P,S,N,M} <: ParametricSurface{S,N}
Curve order, M
, is the same in the u and v direction and fixed over all spans. P
determines the CurveType
.
OpticSim.StopSurface
— TypeStopSurface{T} <: Surface{T}
Abstract type to encapsulate any surfaces acting as a stop.
OpticSim.UniformSpectrum
— Typeflat spectrum from 450nm to 680nm
OpticSim.ZernikeIndexType
— TypeEither ZernikeIndexingOSA
or ZernikeIndexingNoll
, see Zernike polynomials wikipedia entry for details.
OpticSim.Annulus
— MethodAnnulus(innerradius::T, outerradius::T, surfacenormal::SVector{3,T}, centrepoint::SVector{3,T})
Creates a circular aperture in a circle i.e. FiniteStop{T,CircularStopShape,CircularStopShape}
.
OpticSim.ArizonaEye
— MethodArizonaEye(::Type{T} = Float64; accommodation::T = 0.0)
The popular Arizona eye model taken from this definition. The accommodation
of the eye can be varied in this model. Returns a DataFrame
specifying the prescription of the eye model.
OpticSim.AsphericLens
— MethodAsphericLens(insidematerial, frontvertex, frontradius, frontconic, frontaspherics, backradius, backconic, backaspherics, thickness, semidiameter; lastmaterial = OpticSim.GlassCat.Air, nextmaterial = OpticSim.GlassCat.Air, frontsurfacereflectance = 0.0, backsurfacereflectance = 0.0, frontdecenter = (0, 0), backdecenter = (0, 0))
Cosntructs a simple cylindrical lens with front and back surfaces with a radius, conic and apsheric terms. The side walls of the lens are absorbing.
OpticSim.BoundedCylinder
— MethodBoundedCylinder(radius::T, height::T; interface::NullOrFresnel{T} = nullinterface(T)) -> CSGGenerator{T}
Create a cylinder with planar caps on both ends centred at (0, 0, 0)
with axis (0, 0, 1)
.
OpticSim.Circle
— MethodCircle(radius, [surfacenormal, centrepoint]; interface = nullinterface(T))
Shortcut method to create a circle. The minimal case returns a circle centred at the origin with normal = [0, 0, 1]
.
OpticSim.CircularAperture
— MethodCircularAperture(radius::T, surfacenormal::SVector{3,T}, centrepoint::SVector{3,T})
Creates a circular aperture in a plane i.e. InfiniteStop{T,CircularStopShape}
.
OpticSim.CircularAperture
— MethodCircularAperture(radius::T, outerhalfsizeu::T, outerhalfsizev::T, surfacenormal::SVector{3,T}, centrepoint::SVector{3,T}; rotationvec::SVector{3,T} = [0.0, 1.0, 0.0])
Creates a circular aperture in a rectangle i.e. FiniteStop{T,CircularStopShape,RectangularStopShape}
. The rotation of the rectangle around its normal is defined by rotationvec
. rotationvec×surfacenormal
is taken as the vector along the u axis.
OpticSim.ConicLens
— MethodConicLens(insidematerial, frontvertex, frontradius, frontconic, backradius, backconic, thickness, semidiameter; lastmaterial = OpticSim.GlassCat.Air, nextmaterial = OpticSim.GlassCat.Air, frontsurfacereflectance = 0.0, backsurfacereflectance = 0.0, frontdecenter = (0, 0), backdecenter = (0, 0))
Constructs a simple cylindrical lens with front and back surfaces with a radius and conic term. The side walls of the lens are absorbing.
OpticSim.Cuboid
— MethodCuboid(halfsizex::T, halfsizey::T, halfsizez::T; interface::NullOrFresnel{T} = nullinterface(T)) -> CSGGenerator{T}
Create a cuboid centred at (0, 0, 0)
.
OpticSim.FresnelLens
— MethodFresnelLens(insidematerial, frontvertex, radius, thickness, semidiameter, groovedepth; conic = 0.0, aspherics = nothing, outsidematerial = OpticSim.GlassCat.Air)
Create a Fresnel lens as a CSG object, can be concave or convex. Groove positions are found iteratively based on groovedepth
. For negative radii the vertex on the central surface is at frontvertex
, so the total thickness of the lens is thickness
+ groovedepth
. Aspherics currently not supported.
OpticSim.GridField
— MethodGridField(sys::AxisymmetricOpticalSystem; collimated = true, samples = 20, wavelength = 0.55, sourcepos = (0.0, 0.0, 3.0), sourceangle = 0.0, sourcenum = 0)
Distributes rays over the entrance pupil of the system in a rectangular grid pattern.
OpticSim.GridField
— MethodGridField(semidiameter, pupilpos; collimated = true, samples = 20, wavelength = 0.55, sourcepos = (0.0, 0.0, 3.0), sourceangle = 0.0, sourcenum = 0)
Distributes rays over a circular pupil with half-diameter defined by semidiameter
, centred at pupilpos
in a rectangular grid pattern. samples
is the number of rays on each side of the grid, so there are samples×samples
rays in total.
OpticSim.HexagonalPrism
— MethodHexagonalPrism(side_length::T, visheight::T = 2.0; interface::NullOrFresnel{T} = nullinterface(T)) -> CSGGenerator{T}
Create an infinitely tall hexagonal prism with axis (0, 0, 1)
, the longer hexagon diameter is along the x axis. For visualization visheight
is used, note that this does not fully represent the surface.
OpticSim.HexapolarField
— MethodHexapolarField(sys::AxisymmetricOpticalSystem; collimated = true, samples = 8, wavelength = 0.55, sourcepos = (0.0, 0.0, 3.0), sourceangle = 0.0, sourcenum = 0)
Distributes rays over the entrance pupil of the system in a hexapolar pattern.
OpticSim.HexapolarField
— MethodHexapolarField(semidiameter, pupilpos; collimated = true, samples = 8, wavelength = 0.55, sourcepos = (0.0, 0.0, 3.0), sourceangle = 0.0, sourcenum = 0)
Distributes rays over a circular pupil with half-diameter defined by semidiameter
, centred at pupilpos
in a hexapolar pattern. samples
is the number of rings in the hexapolar pattern, so the number of rays in total is samples * (samples + 1) / 2) * 6 + 1
.
OpticSim.ModelEye
— MethodModelEye(assembly::LensAssembly{T}, nsamples::Int = 17; pupil_radius::T = 3.0, detpixels::Int = 1000, transform::RigidBodyTransform{T} = identitytransform(T))
Geometrically accurate model of the human eye focussed at infinity with variable pupil_radius
. The eye is added to the provided assembly
to create a CSGOpticalSystem
with the retina of the eye as the detector.
The eye can be positioned in the scene using the transform
argument and the resolution of the detector specified with detpixels
. By default the eye is directed along the positive z-axis with the vertex of the cornea at the origin.
nsamples
determines the resolution at which accelerated surfaces within the eye are triangulated.
OpticSim.RectangularAperture
— MethodRectangularAperture(aphalfsizeu::T, aphalfsizev::T, surfacenormal::SVector{3,T}, centrepoint::SVector{3,T}; rotationvec::SVector{3,T} = [0.0, 1.0, 0.0])
Creates a rectangular aperture in a plane i.e. InfiniteStop{T,RectangularStopShape}
. The rotation of the rectangle around its normal is defined by rotationvec
. rotationvec×surfacenormal
is taken as the vector along the u axis.
OpticSim.RectangularAperture
— MethodRectangularAperture(innerhalfsizeu::T, innerhalfsizev::T, outerhalfsizeu::T, outerhalfsizev::T, surfacenormal::SVector{3,T}, centrepoint::SVector{3,T}; rotationvec::SVector{3,T} = [0.0, 1.0, 0.0])
Creates a rectangular aperture in a rectangle i.e. FiniteStop{T,RectangularStopShape,RectangularStopShape}
. The rotation of the rectangle around its normal is defined by rotationvec
. rotationvec×surfacenormal
is taken as the vector along the u axis.
OpticSim.RectangularPrism
— MethodRectangularPrism(halfsizex::T, halfsizey::T, visheight::T=2.0; interface::NullOrFresnel{T} = nullinterface(T)) -> CSGGenerator{T}
Create an infinitely tall rectangular prism with axis (0, 0, 1)
. For visualization visheight
is used, note that this does not fully represent the surface.
OpticSim.SphericalLens
— MethodSphericalLens(insidematerial, frontvertex, frontradius, backradius, thickness, semidiameter; lastmaterial = OpticSim.GlassCat.Air, nextmaterial = OpticSim.GlassCat.Air, frontsurfacereflectance = 0.0, backsurfacereflectance = 0.0, frontdecenter = (0, 0), backdecenter = (0, 0))
Constructs a simple cylindrical lens with spherical front and back surfaces. The side walls of the lens are absorbing.
OpticSim.Spider
— MethodSpider(narms::Int, armwidth::T, radius::T, origin::SVector{3,T} = SVector{3,T}(0.0, 0.0, 0.0), normal::SVector{3,T} = SVector{3,T}(0.0, 0.0, 1.0)) -> Vector{Rectangle{T}}
Creates a 'spider' obscuration with narms
rectangular arms evenly spaced around a circle defined by origin
and normal
. Each arm is a rectangle armwidth
×radius
.
e.g. for 3 and 4 arms we get:
| _|_
/ \ |
OpticSim.TriangularPrism
— MethodTriangularPrism(side_length::T, visheight::T = 2.0; interface::NullOrFresnel{T} = nullinterface(T)) -> CSGGenerator{T}
Create an infinitely tall triangular prism with axis (0, 0, 1)
. For visualization visheight
is used, note that this does not fully represent the surface.
OpticSim.assembly
— Methodassembly(system::OpticalSystem{T}) -> LensAssembly{T}
Get the LensAssembly
of system
.
OpticSim.closestintersection
— Functionclosestintersection(a::Union{EmptyInterval{T},Interval{T},DisjointUnion{T}}, ignorenull::Bool = true) -> Union{Nothing,Intersection{T,3}}
Returns the closest Intersection
from an Interval
or DisjointUnion
. Ignores intersection with null interfaces if ignorenull
is true. Will return nothing
if there is no valid intersection.
OpticSim.closestpointonray
— Methodclosestpointonray(r::Ray{T,N}, point::SVector{N,T}) -> SVector{T,N
Returns the point on the ray closest to point.
OpticSim.csgdifference
— Methodcsgdifference(a::CSGGenerator{T}, b::CSGGenerator{T}, transform::RigidBodyTransform{T} = identitytransform(T)) -> CSGGenerator{T}
Create a binary node in the CSG tree representing the difference of a and b, essentially a - b. A shortcut method for a
and b
as ParametricSurface
s is also available.
OpticSim.csgintersection
— Methodcsgintersection(a::CSGGenerator{T} b::CSGGenerator{T}, transform::RigidBodyTransform{T} = identitytransform(T)) -> CSGGenerator{T}
Create a binary node in the CSG tree representing an intersection between a and b. A shortcut method for a
and b
as ParametricSurface
s is also available.
OpticSim.csgunion
— Methodcsgunion(a::CSGGenerator{T}, b::CSGGenerator{T}, transform::RigidBodyTransform{T} = identitytransform(T)) -> CSGGenerator{T}
Create a binary node in the CSG tree representing a union between a and b. A shortcut method for a
and b
as ParametricSurface
s is also available.
OpticSim.curvedimension
— Methodspatial dimension of curve represented as an array of coefficients x[i] = ∑Bj(θ)*x[i,j]
where Bj(θ)
is the curve basis
OpticSim.curveorder
— Methodhighest polynomial power of the curve represented as an array of coefficients x[i] = ∑Bj(θ)*x[i,j]
where Bj(θ)
is the curve basis
OpticSim.detectorimage
— Methoddetectorimage(system::OpticalSystem{T}) -> HierarchicalImage{D}
Get the detector image of system
. D
is the datatype of the detector image and is not necessarily the same as the datatype of the system T
.
OpticSim.direction
— MethodGenerates a unit vector pointing somewhere within the cone with half angle θmax
around direction
which is the normal to the emitter surface.
OpticSim.distance
— Methoddistance(r::Ray{T,N}, point::SVector{N,T}) -> Union{T,Nothing}
Returns distance to the position on the ray closest to point. If t < 0 returns nothing.
OpticSim.doesintersect
— Methoddoesintersect(bbox::BoundingBox{T}, r::AbstractRay{T,3}) -> Bool
Tests whether r
intersects an axis-aligned BoundingBox
, bbox
.
OpticSim.evaluatecurve
— MethodEvaluates a curve defined in the power basis. Curves and moving lines accessed like this: [xi,ci]
where xi
is the dimension index, and ci
is the coefficient index.
OpticSim.extractmovinglines
— Methodreturns 3D array indexed like this: x[line curve order,spatial dimension, line number]
`
OpticSim.fresnel
— Methodfresnel(nᵢ::T, nₜ::T, sinθᵢ::T, sinθₜ::T) -> Tuple{T,T}
Returns reflectance and tranmission power coefficients according to the Fresnel equations. For geometric ray tracing this coefficient can be used directly to compute intensity on the detector plane. For Huygens phase optics need to take the square root to compute the amplitude. The power of the transmitted and refracted rays may not sum to one because of the area correction applied to the transmitted component. The intensity per area can increase or decrease depending on the indices of refraction.
nᵢ
is the RI of the material which the incident ray travels in, nₜ
is the RI of the material the transmitted ray travels in. sinθᵢ
and sinθₜ
are the sin of the angles of incidence and transmission respectively.
OpticSim.gendirection
— Methodgendirection(o::GeometricRayGenerator{T}, n::Int) -> SVector{3,T}
Generate directions for rays based on the type of the generator, e.g., randomly within a cone or collimated. n
is the index of the point being generated, starting from 0. This has little meaning for random generators, but is important for GridSource
.
OpticSim.generateray
— Methodgenerateray(a::PixelSource{T}, n::Int) -> OpticalRay{T,3}
Generates optical rays from all subpixels in the pixel. One ray is generated from each subpixel sequentially before looping back to the start.
OpticSim.generateray
— Methodgenerateray(o::GeometricRayGenerator{T}, n::Int) -> Ray{T,3}
Generate geometric rays distributed according to the type of the generator. n
is the index of the point being generated, starting from 0. This has little meaning for random generators, but is important for GridSource
, for example.
OpticSim.generateray
— Methodgenerateray(a::BasicDisplayPanel{T}, n::Int) -> OpticalRay{T,3}
Generates optical rays from all pixels in the display. One ray is generated from each pixel sequentially before looping back to the start of the display.
OpticSim.generateray
— Methodgenerateray(o::OpticalRayGenerator{T}, n::Int) -> OpticalRay{T,3}
Generate optical rays distributed according to the type of the generator. n
is the index of the point being generated, starting from 0. This has little meaning for random generators, but is important for generators using GridSource
or GridRectOriginPoints
, for example.
OpticSim.generateray
— Methodgenerateray(a::OpticalSourceArray{T}, n::Int) -> OpticalRay{T,3}
Generates optical rays from all generators in the array. One ray is generated from each element sequentially before looping back to the start of the array.
OpticSim.generateray
— Methodgenerateray(a::OpticalSourceGroup{T}, n::Int) -> OpticalRay{T,3}
Generate optical rays for each source in the group. All rays are generated for the first source, then all for the second source and so on as n
increases.
OpticSim.genorigin
— Methodgenorigin(o::RayOriginGenerator{T}, n::Int) -> SVector{3,T}
Generate origin positions for rays based on the type of the generator, e.g., randomly within a rectangle or ellipse. n
is the index of the point being generated, starting from 0. This has little meaning for random generators, but is important for HexapolarOriginPoints
and GridRectOriginPoints
.
OpticSim.halfspaceintersection
— Methodhalfspaceintersection(a::Interval{T}) -> Intersection{T,3}
Returns the Intersection
from a half space Interval
, throws an error if not a half space.
OpticSim.identitytransform
— Methodidentitytransform([S::Type]) -> RigidBodyTransform{S}
Returns the RigidBodyTransform
of type S
(default Float64
) representing the identity transform.
OpticSim.interface
— Methodinterface(surf::Surface{T}) -> OpticalInterface{T}
Return the OpticalInterface
associated with surf
.
OpticSim.intersections
— Methodreturns an array of intersection points. Each element in the array is ([x,y,...],alpha,theta)
where [x,y,...]
is the n-dimensional intersection point, alpha is the line parameter value at the intersection point, and theta is the curve parameter value at the intersection point
OpticSim.isemptyinterval
— Methodisemptyinterval(a) -> Bool
Returns true if a
is an EmptyInterval
. In performance critical contexts use a isa EmptyInterval{T}
.
OpticSim.isinfinity
— Methodisinfinity(a) -> Bool
Returns true if a
is Infinity
. In performance critical contexts use a isa Infinity{T}
.
OpticSim.ispositivehalfspace
— Methodispositivehalfspace(a) -> Bool
Returns true if upper(a)
is Infinity
. In performance critical contexts check directly i.e. upper(a) isa Infinity{T}
.
OpticSim.israyorigin
— Methodisrayorigin(a) -> Bool
Returns true if a
is RayOrigin
. In performance critical contexts use a isa RayOrigin{T}
.
OpticSim.israyorigininterval
— Methodisrayorigininterval(a) -> Bool
Returns true if lower(a)
is RayOrigin
. In performance critical contexts check directly i.e. lower(a) isa RayOrigin{T}
.
OpticSim.jacobian
— Methodjacobian(surf::ParametricSurface{T,N}, u::T, v::T, P1::SVector{M,T}, P2::SVector{M,T})
Computes Jacobian of f(t,u,v) = ( dot(P1,[surf(u,v),1],P2,[surf(u,v),1]) )
. P1
, P2
are orthogonal planes that pass through the ray. J = [ ∂f1/∂u ∂f1/∂v ; ∂f2/∂u ∂f2/∂v]
OpticSim.leaf
— Methodleaf(surf::CSGGenerator{T}, transform::RigidBodyTransform{T} = identitytransform(T)) -> CSGGenerator{T}
Create a (pseudo) leaf node from another CSGGenerator, this is useful if you want multiple copies of a premade CSG structure with different transforms, for example in an MLA.
OpticSim.leaf
— Methodleaf(surf::ParametricSurface{T}, transform::RigidBodyTransform{T} = identitytransform(T)) -> CSGGenerator{T}
Create a leaf node from a parametric surface with a given transform.
OpticSim.linedimension
— Methodspatial dimension of the moving line represented as an array of coefficients g[i] = ∑Bl(θ)*gl[i,j]
where Bl(θ)
is the polynomial basis
OpticSim.makemesh
— Methodmakemesh(object, subdivisions::Int = 30) -> TriangleMesh
Creates a TriangleMesh
from an object, either a ParametricSurface
, CSGTree
or certain surfaces (e.g. Circle
, Rectangle
). This is used for visualization purposes only.
OpticSim.matricesforeigen
— Methodmovinglines[:,i]
is the ith moving line. For li = movinglines[:,i] (dimension+1,lineorder) = size(li)
. rline[:,1] = pt1
and rline[:,2] = pt2
. The line equation is pt1 + alpha*pt2
.
OpticSim.newton
— Methodnewton(surf::ParametricSurface{T,N}, r::AbstractRay{T,N}, startingpoint::SVector{2,T})
Newton iteration to find the precise intersection of a parametric surface with a ray given a starting point (in uv space) on the surface.
OpticSim.normal
— Methodnormal(surf::ParametricSurface{T}, u::T, v::T) -> SVector{3,T}
normal(surf::ParametricSurface{T}, uv::SVector{2,T}) -> SVector{3,T}
Returns the normal to surf
at the given uv coordinate.
OpticSim.numberoflines
— Methodnumber of lines in moving line array
OpticSim.orthogonalitymatrix
— Methodreturns a matrix expressing the relationship [x(θ) 1]⋅g(θ) = 0
. The vectors in the right nullspace of this matrix contain the coefficients of the moving lines gᵢ(θ)
.
OpticSim.point
— Methodpoint(ray::AbstractRay{T,N}, alpha::T) -> SVector{T, N}
Returns a point on the ray at origin + alpha * direction. Alpha must be >= 0.
OpticSim.pressure
— Methodpressure(system::OpticalSystem{T}) -> T
Get the pressure of system
in Atm.
OpticSim.processintersection
— Methodprocessintersection(opticalinterface::OpticalInterface{T}, point::SVector{N,T}, normal::SVector{N,T}, incidentray::OpticalRay{T,N}, temperature::T, pressure::T, ::Bool, firstray::Bool = false) -> Tuple{SVector{N,T}, T, T}
Processes an intersection of an OpticalRay
with an OpticalInterface
, distinct behaviors must be implemented for each subclass of OpticalInterface
.
point
is the 3D intersection point in global space, normal
is the surface normal at the intersection point.
If test
is true then the behavior of the ray should be deterministic. firstray
indicates that this ray is the first segment of the trace and therefore the origin is not offset.
The values returned are the normalized direction of the ray after the intersection, the instantaneous power of the ray after the intersection and the optical path length of the ray up to the intersection.
nothing
is returned if the ray should stop here, in order to obtain the correct intensity on the detector through monte carlo integration nothing
should be returned proportionally to create the correct power distribution. i.e. If the interface should modulate power to 76% then 24% of calls to this function should return nothing
.
OpticSim.randinsolidangle
— Methodrandinsolidangle(direction::SVector{3,T}, uvec::SVector{3,T}, vvec::SVector{3,T}, θmax::T)
Generates a unit vector pointing somewhere within the cone with half angle θmax
around direction
. uvec
and vvec
should be orthogonal to each other and direction
.
OpticSim.raysample
— Methodraysample is the function that can couple origin and direction generation if necessary. The default function couples them in a simple way but more complex coupling should be possible. For each origin numsamples(a.raydirection) direction samples are taken with identical origin. Then the origin number is incremented. This repeats till all rays have been generated. The origin and direction functions receive an integer indicating the origin or direction number so regular patterns such as rectangular and hexapolar grids can be generated properly.
OpticSim.reset!
— Methodreset!(a::HierarchicalImage{T})
Resets the pixels in the image to zero(T)
. Do this rather than image .= zero(T)
because that will cause every pixel to be accessed, and therefore allocated. For large images this can cause huge memory traffic.
OpticSim.resetdetector!
— Methodresetdetector!(system::OpticalSystem{T})
Reset the deterctor image of system
to zero.
OpticSim.reversenormal
— Methodreversenormal(a::Intersection{T,N})
Used by the CSG complement operator (i.e. csgdifference
) to reverse the inside outside sense of the object.
OpticSim.rotation
— Methodrotation([S::Type], θ::T, ϕ::T, ψ::T) -> RigidBodyTransform{S}
Returns the RigidBodyTransform
of type S
(default Float64
) representing the rotation by θ
, ϕ
and ψ
around the x, y and z axes respectively in radians.
OpticSim.rotationd
— Methodrotationd([S::Type], θ::T, ϕ::T, ψ::T) -> RigidBodyTransform{S}
Returns the RigidBodyTransform
of type S
(default Float64
) representing the rotation by θ
, ϕ
and ψ
around the x, y and z axes respectively in degrees.
OpticSim.rotmat
— Methodrotmat([S::Type], θ::T, ϕ::T, ψ::T) -> SMatrix{3,3,S}
Returns the rotation matrix of type S
(default Float64
) representing the rotation by θ
, ϕ
and ψ
around the x, y and z axes respectively in radians.
OpticSim.rotmatbetween
— Methodrotmatbetween([S::Type], a::SVector{3,T}, b::SVector{3,T}) -> SMatrix{3,3,S}
Returns the rotation matrix of type S
(default Float64
) representing the rotation between vetors a
and b
.
OpticSim.rotmatd
— Methodrotmatd([S::Type], θ::T, ϕ::T, ψ::T) -> SMatrix{3,3,S}
Returns the rotation matrix of type S
(default Float64
) representing the rotation by θ
, ϕ
and ψ
around the x, y and z axes respectively in degrees.
OpticSim.samplesurface
— Methodsamplesurface(surf::ParametricSurface{T,N}, samplefunction::Function, numsamples::Int = 30)
Sample a parametric surface on an even numsamples
×numsamples
grid in UV space with provided function
OpticSim.semidiameter
— Methodsemidiameter(system::AxisymmetricOpticalSystem{T}) -> T
Get the semidiameter of system
, that is the semidiameter of the entrance pupil (i.e. first surface) of the system.
OpticSim.snell
— Methodsnell(surfacenormal::AbstractVector{T}, raydirection::AbstractVector{T}, nᵢ::T, nₜ::T) -> Tuple{T,T}
nᵢ
is the index of refraction on the incidence side of the interface. nₜ
is the index of refraction on the transmission side.
Returns sinθᵢ
and sinθₜ
according to Snell's law.
OpticSim.spectrumpower
— Methodexpects wavelength in nm not um
OpticSim.spectrumsample
— Methodreturns a tuple (power,wavelength)
OpticSim.sum!
— Methodsum!(a::HierarchicalImage{T}, b::HierarchicalImage{T})
Add the contents of b
to a
in an efficient way.
OpticSim.surfaceintersection
— Methodsurfaceintersection(surf::Surface{T}, r::AbstractRay{T}) where {T}
Calculates the intersection of r
with a surface of any type, surf
. Note that some surfaces cannot be intersected analytically so must be wrapped in an AcceleratedParametricSurface
in order to be intersected.
Returns an EmptyInterval
if there is no Intersection
, an Interval
if there is one or two intersections and a DisjointUnion
if there are more than two intersections.
OpticSim.temperature
— Methodtemperature(system::OpticalSystem{T}) -> T
Get the temperature of system
in °C.
OpticSim.trace
— Methodtrace(system::OpticalSystem{T}, ray::OpticalRay{T}; trackrays = nothing, test = false)
Traces system
with ray
, if test
is enabled then fresnel reflections are disabled and the power distribution will not be correct. Returns either a LensTrace
if the ray hits the detector or nothing
otherwise.
trackrays
can be passed an empty vector to accumulate the LensTrace
objects at each intersection of ray
with a surface in the system.
OpticSim.trace
— Methodtrace(assembly::LensAssembly{T}, r::OpticalRay{T}, temperature::T = 20.0, pressure::T = 1.0; trackrays = nothing, test = false)
Returns the ray as it exits the assembly in the form of a LensTrace
object if it hits any element in the assembly, otherwise nothing
. Recursive rays are offset by a small amount (RAY_OFFSET
) to prevent it from immediately reintersecting the same lens element.
trackrays
can be passed an empty vector to accumulate the LensTrace
objects at each intersection of ray
with a surface in the assembly.
OpticSim.trace
— Methodtrace(system::OpticalSystem{T}, raygenerator::OpticalRayGenerator{T}; printprog = true, test = false)
Traces system
with rays generated by raygenerator
on a single thread. Optionally the progress can be printed to the REPL. If test
is enabled then fresnel reflections are disabled and the power distribution will not be correct. If outpath
is specified then the result will be saved to this path.
Returns the detector image of the system.
OpticSim.traceMT
— MethodtraceMT(system::OpticalSystem{T}, raygenerator::OpticalRayGenerator{T}; printprog = true, test = false)
Traces system
with rays generated by raygenerator
using as many threads as possible. Optionally the progress can be printed to the REPL. If test
is enabled then fresnel reflections are disabled and the power distribution will not be correct. If outpath
is specified then the result will be saved to this path.
Returns the accumulated detector image from all threads.
OpticSim.tracehits
— Methodtracehits(system::OpticalSystem{T}, raygenerator::OpticalRayGenerator{T}; printprog = true, test = false)
Traces system
with rays generated by raygenerator
on a single thread. Optionally the progress can be printed to the REPL. If test
is enabled then fresnel reflections are disabled and the power distribution will not be correct.
Returns a list of LensTrace
s which hit the detector.
OpticSim.tracehitsMT
— MethodtracehitsMT(system::OpticalSystem{T}, raygenerator::OpticalRayGenerator{T}; printprog = true, test = false)
Traces system
with rays generated by raygenerator
using as many threads as possible. Optionally the progress can be printed to the REPL. If test
is enabled then fresnel reflections are disabled and the power distribution will not be correct.
Returns a list of LensTrace
s which hit the detector, accumulated from all threads.
OpticSim.translation
— Methodtranslation([S::Type], x::T, y::T, z::T) -> RigidBodyTransform{S}
Returns the RigidBodyTransform
of type S
(default Float64
) representing the translation by x
, y
and z
.
OpticSim.triangulate
— Methodtriangulate(surf::ParametricSurface{S,N}, quads_per_row::Int, extensionu::Bool = false, extensionv::Bool = false, radialu::Bool = false, radialv::Bool = false)
Create an array of triangles representing the parametric surface where vertices are sampled on an even grid in UV space. The surface can be extended by 1% in u and v separately, and specifying either u or v as being radial - i.e. detemining the radius on the surface e.g. rho for zernike - will result in that dimension being sampled using sqwrt so that area of triangles is uniform. The extension will also only apply to the maximum in this case.
OpticSim.triangulatedintersection
— Methodtriangulatedintersection(surf::AcceleratedParametricSurface{T,N,S}, r::AbstractRay{T,N})
Intersection of a ray, r
, with a triangulated surface, surf
, no concept of inside so never returns a RayOrigin
Interval
.
OpticSim.uv
— Methoduv(surf::ParametricSurface{T}, p::SVector{3,T}) -> SVector{2,T}
uv(surf::ParametricSurface{T}, x::T, y::T, z::T) -> SVector{2,T}
Returns the uv coordinate on surf
of a point, p
, in 3D space. If onsurface(surf, p)
is false then the behavior is undefined, it may return an inorrect uv, an invalid uv, NaN or crash.
OpticSim.uvrange
— Methoduvrange(s::ParametricSurface)
uvrange(::Type{S}) where {S<:ParametricSurface}
Returns a tuple of the form: ((umin, umax), (vmin, vmax))
specifying the limits of the parameterisation for this surface type. Also implemented for some Surface
s which are not ParametricSurface
s (e.g. Rectangle
).
OpticSim.uvtopix
— Methoduvtopix(surf::Surface{T}, uv::SVector{2,T}, imsize::Tuple{Int,Int}) -> Tuple{Int,Int}
Converts a uvcoordinate on surf
to an integer index to a pixel in an image of size imsize
. Not implemented on all Surface
objects. Used to determine where in the detector image a ray has hit when in intersects the detector surface of an OpticalSystem
.
OpticSim.α
— Methodα(ray::AbstractRay{T,N}, point::SVector{N,T}) -> T
Computes the alpha corresponding to the closest position on the ray to point
Zernike
OpticSim.Zernike
— ModuleModule to enclose Zernike polynomial specific functionality.
OpticSim.Zernike.NolltoNM
— MethodNolltoNM(j::Int) -> Tuple{Int, Int}
Convert Noll zernike index j
to (N,M)
form.
OpticSim.Zernike.OSAtoNM
— MethodOSAtoNM(j::Int) -> Tuple{Int, Int}
Convert OSA zernike index j
to (N,M)
form according to formula J = N * (N + 2) + M
.
OpticSim.Zernike.R
— MethodR(N::Int, M::Int, ρ::T) -> T
Evaluate radial polynomial $R_{n}^{m}(\rho)$.
OpticSim.Zernike.normalisation
— Methodnormalisation(::Type{T}, N::Int, M::Int) -> T
Normalisation coefficient for Zernike polynomial term $Z_{n}^{m}$.
OpticSim.Zernike.δζ
— Methodδζ(N::Int, M::Int, ρ::T, ϕ::T) -> Tuple{T,T}
Evaluate partial derivatives of Zernike polynomial term $Z_{n}^{m}(\rho, \phi)$.
OpticSim.Zernike.ζ
— Methodζ(N::Int, M::Int, ρ::T, ϕ::T) -> Tuple{T,T}
Evaluate Zernike polynomial term $Z_{n}^{m}(\rho, \phi)$.
QType
OpticSim.QType
— ModuleModule to enclose QType polynomial specific functionality. For reference see:
OpticSim.QType.S
— MethodS(coeffs::SVector{NP1,T}, m::Int x::T) -> T
Evaluates $\sum_{n=0}^{N}c_n^mQ_n^m(x)$ where $c_n^m$ is either an $\alpha$ or $\beta$ QType coefficient and $m \gt 0$.
OpticSim.QType.S0
— MethodS0(coeffs::SVector{NP1,T}, x::T) -> T
Evaluates $\sum_{n=0}^{N}\alpha_n^0Q_n^0(x)$.
OpticSim.QType.dS0dx
— MethoddS0dx(coeffs::SVector{NP1,T}, x::T) -> T
Evaluates $\frac{\partial}{\partial x}\sum_{n=0}^{N}\alpha_n^0Q_n^0(x)$.
OpticSim.QType.dSdx
— MethoddSdx(coeffs::SVector{NP1,T}, x::T) -> T
Evaluates $\frac{\partial}{\partial x}\sum_{n=0}^{N}c_n^mQ_n^m(x)$ where $c_n^m$ is either an $\alpha$ or $\beta$ QType coefficient and $m \gt 0$.
Chebyshev
OpticSim.Chebyshev
— ModuleModule to enclose Chebyshev polynomial specific functionality.
OpticSim.Chebyshev.T
— MethodT(n::Int, q::R, fast::Bool = true) -> R
Evaluate Chebyshev polynomial of the first kind $T_n(q)$.
fast
will use trigonometric definition, rather than the recursive definition which is much faster but slightly less precise.
OpticSim.Chebyshev.U
— MethodU(n::Int, q::R, fast::Bool = true) -> R
Evaluate Chebyshev polynomial of the second kind $U_n(q)$.
fast
will use trigonometric definition, rather than the recursive definition which is much faster but slightly less precise.
OpticSim.Chebyshev.dTdq
— MethoddTdq(n::Int, q::R, fast::Bool = true) -> R
Evaluate derivative of Chebyshev polynomial of the first kind $\frac{dT_n}{dq}(q)$.
fast
will use trigonometric definition, rather than the recursive definition which is much faster but slightly less precise.
Examples
OpticSim.Examples
— ModuleContains example usage of the features in the OpticSim.jl package.