Complete Reference

This page contains what should be a complete list of all docstrings in the OpticSim module, and its submodule.

Index

OpticSim

OpticSim.BSplineCurveType
BSplineCurve{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})
source
OpticSim.BezierCurveType
BezierCurve{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})
source
OpticSim.ConvexPolygonType
ConvexPolygon{N, T<:Real} <: PlanarShape{T}

General Convex Polygon surface, not a valid CSG object. The rotation of the polygon around its normal is defined by rotationvec. rotationvec×surfacenormal is taken as the vector along the u axis.

ConvexPolygon(local_frame::Transform{T}, local_polygon_points::Vector{SVector{2, T}}, interface::NullOrFresnel{T} = nullinterface(T))

The local frame defines the plane (spans by the right and up vectors) with the plane normal given by the forward vector. the localpolygonpoints are given with respect to the local frame and are 2D points. NOTE: This class uses static vectors to hold the points which will lead to more efficient performance, but should not be used with polygons with more than 20-30 points.

source
OpticSim.GeometricRayGeneratorType
GeometricRayGenerator{T,O<:RayOriginGenerator{T}} <: AbstractRayGenerator{T}

Generates geometric Rays according to the specific implementation of the subclass.

source
OpticSim.OpticalRayType
OpticalRay{T,N} <: AbstractRay{T,N}

Ray with power, wavelength and optical path length.

NOTE: we use monte carlo integration to get accurate results on the detector, this means that all rays essentially hit the detector with power = 1 and some rays are thrown away at any interface to correctly match the reflection/transmission at that interface. For inspection purposes we also track the 'instantaneous' power of the ray in the power field of the OpticalRay.

OpticalRay(ray::Ray{T,N}, power::T, wavelength::T, opl=zero(T))
OpticalRay(origin::SVector{N,T}, direction::SVector{N,T}, power::T, wavelength::T, opl=zero(T))

Has the following accessor methods:

ray(r::OpticalRay{T,N}) -> Ray{T,N}
direction(r::OpticalRay{T,N}) -> SVector{N,T}
origin(r::OpticalRay{T,N}) -> SVector{N,T}
power(r::OpticalRay{T,N}) -> T
wavelength(r::OpticalRay{T,N}) -> T
pathlength(r::OpticalRay{T,N}) -> T
sourcepower(r::OpticalRay{T,N}) -> T
nhits(r::OpticalRay{T,N}) -> Int
sourcenum(r::OpticalRay{T,N}) -> Int
source
OpticSim.PlanarShapeType

The PlanarShape interface:

distancefromplane(p::PlanarShape,point) returns distance of the point from the plane the planar shape lies within normal(p::PlanarShape) returns normal of plane interface(p::PlanarShape) returns optical interface of plane vertices(p::PlanarShape) returns vertices of shape. For Ellipse this is an approximation.

There are default functions for plane,normal,interface,vertices which assume each PlanarShape type has a field of the same name plane(a::PlanarShape) = a.plane normal(a::PlanaShape) = a.plane.normal etc.

If your type doesn't have these fields then you should define a more specialized method to handle this.

source
OpticSim.PrimitiveType
Primitive{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}
source
OpticSim.RayType
Ray{T,N} <: AbstractRay{T,N}

Purely geometric ray, defined as origin + alpha * direction.

Ray(origin::SVector{N,T}, direction::SVector{N,T})

Has the following accessor methods:

direction(ray::Ray{T,N}) -> SVector{N,T}
origin(ray::Ray{T,N}) -> SVector{N,T}
source
OpticSim.SphericalPolygonType

SphericalPolygon uses StaticArrays to represent vertices. Expect performance degradation for polygons with large numbers of vertices. Performance appears to be good up to perhaps 100 vertices, perhaps as much as 1000 vertices. By 10,000 vertices performance is terrible.

source
OpticSim.SplineType
Spline{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

source
OpticSim.SplineSurfaceType
SplineSurface{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.

source
Base.:*Method

Apply a Transform to an Intersection object

source
Base.:*Method

Apply a Transform to a TriangleMesh object

source
OpticSim.AnnulusMethod
Annulus(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}.

source
OpticSim.AsphericLensMethod
AsphericLens(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), interfacemode = ReflectOrTransmit)

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.

source
OpticSim.BoundedCylinderMethod
BoundedCylinder(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).

source
OpticSim.CircleMethod
Circle(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].

source
OpticSim.CircularApertureMethod
CircularAperture(radius::T, surfacenormal::SVector{3,T}, centrepoint::SVector{3,T})

Creates a circular aperture in a plane i.e. InfiniteStop{T,CircularStopShape}.

source
OpticSim.CircularApertureMethod
CircularAperture(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.

source
OpticSim.ConicLensMethod
ConicLens(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), interfacemode = ReflectOrTransmit)

Constructs a simple cylindrical lens with front and back surfaces with a radius and conic term. The side walls of the lens are absorbing.

source
OpticSim.CuboidMethod
Cuboid(halfsizex::T, halfsizey::T, halfsizez::T; interface::NullOrFresnel{T} = nullinterface(T)) -> CSGGenerator{T}

Create a cuboid centred at (0, 0, 0).

source
OpticSim.EvenAsphericSurfaceMethod
EvenAsphericSurface(semidiameter, curvature::T, conic::T, aspherics::Vector{T}; normradius::T=semidiameter)

Surface incorporating an aspheric polynomial - radius, conic and aspherics are defined relative to absolute semi-diameter.

aspherics should be an array of the even coefficients of the aspheric polynomial starting with A2

source
OpticSim.FresnelLensMethod
FresnelLens(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.

source
OpticSim.HexagonalPrismMethod
HexagonalPrism(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.

source
OpticSim.OddAsphericSurfaceMethod
OddAsphericSurface(semidiameter, curvature::T, conic::T, aspherics::Vector{T}; normradius::T=semidiameter)

Surface incorporating an aspheric polynomial - radius, conic and aspherics are defined relative to absolute semi-diameter.

aspherics should be an array of the odd coefficients of the aspheric polynomial starting with A1

source
OpticSim.OddEvenAsphericSurfaceMethod
OddEvenAsphericSurface(semidiameter, curvature::T, conic::T, aspherics::Vector{T}; normradius::T=semidiameter)

Surface incorporating an aspheric polynomial - radius, conic and aspherics are defined relative to absolute semi-diameter.

aspherics should be an array of the both odd and even coefficients of the aspheric polynomial starting with A1

source
OpticSim.RectangularApertureMethod
RectangularAperture(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.

source
OpticSim.RectangularApertureMethod
RectangularAperture(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.

source
OpticSim.RectangularPrismMethod
RectangularPrism(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.

source
OpticSim.SphericalLensMethod
SphericalLens(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), interfacemode = ReflectOrTransmit)

Constructs a simple cylindrical lens with spherical front and back surfaces. The side walls of the lens are absorbing.

source
OpticSim.SpiderMethod
Spider(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:

   |         _|_
  / \         |
source
OpticSim.TriangularPrismMethod
TriangularPrism(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.

source
OpticSim.areaMethod

Conceptually breaks the convex spherical polygon into spherical triangles and computes the sum of the angles of all the triangles. The sum of all the angles around the centroid is 2π. Have to subtract π for each of the N triangles. Rather than compute the angles of triangles formed by taking edges from the centroid to each vertex, can instead just compute the internal angle of neighboring edges. Total polygon area is 2π -Nπ + ∑(interior angles).

source
OpticSim.asphericTypeMethod
asphericType(surf::AsphericSurface)

Query the polynomial type of `asp. Returns CONIC, ODD, EVEN, or ODDEVEN. CONIC corresponds to no aspheric terms, ODD means it only has odd aspheric terms, EVEN means only even aspheric terms and ODDEVEN means both even and odd terms.

This function is to enable proper interpretation of surf.aspherics by any optimization routines that directly query the aspheric coefficients.

source
OpticSim.closestintersectionFunction
closestintersection(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.

source
OpticSim.curvedimensionMethod

spatial dimension of curve represented as an array of coefficients x[i] = ∑Bj(θ)*x[i,j] where Bj(θ) is the curve basis

source
OpticSim.curveorderMethod

highest polynomial power of the curve represented as an array of coefficients x[i] = ∑Bj(θ)*x[i,j] where Bj(θ) is the curve basis

source
OpticSim.detectorimageMethod
detectorimage(system::AbstractOpticalSystem{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.

source
OpticSim.distanceMethod
distance(r::Ray{T,N}, point::SVector{N,T}) -> Union{Nothing,T}

Returns distance to the position on the ray closest to point. If t < 0 returns nothing.

source
OpticSim.distancefromplaneMethod

All planar shapes lie on a plane. This function computes the distance from a point to that plane. This is a signed distance. If the point is on the positive side of the plane (the side the normal points toward) the distance will be positive, otherwise negative or 0 if the point lies in the plane.

source
OpticSim.evalcsgFunction
evalcsg(
    a::Union{UnionNode{T},IntersectionNode{T},ComplementNode{T},LeafNode{T}},
    ray::AbstractRay{T,N},
    normalreverse::Bool = false
)::Union{EmptyInterval{T},DisjointUnion{T},Interval{T}}

[TODO]

source
OpticSim.evaluatecurveMethod

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

source
OpticSim.fresnelMethod
fresnel(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.

source
OpticSim.intersectionsMethod

returns 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

source
OpticSim.jacobianMethod
jacobian(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]

source
OpticSim.leafMethod
leaf(surf::ParametricSurface{T}, transform::Transform{T} = identitytransform(T)) -> CSGGenerator{T}

Create a leaf node from a parametric surface with a given transform.

source
OpticSim.linedimensionMethod

spatial dimension of the moving line represented as an array of coefficients g[i] = ∑Bl(θ)*gl[i,j] where Bl(θ) is the polynomial basis

source
OpticSim.makemeshMethod
makemesh(poly::ConvexPolygon{N, T}, ::Int = 0) where {N, T<:Real} -> TriangleMesh

Create a triangle mesh that can be rendered by iterating on the polygon's edges and for each edge use the centroid as the third vertex of the triangle.

source
OpticSim.matricesforeigenMethod

movinglines[:,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.

source
OpticSim.newtonMethod
newton(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.

source
OpticSim.normalMethod
normal(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.

source
OpticSim.orthogonalitymatrixMethod

returns 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ᵢ(θ).

source
OpticSim.plane_from_pointsMethod
plane_from_points(points::SMatrix{D, N, P}}) ->  centroid, normal, local_to_world transform

Points to be fitted are assumed to be stored by column in the points matrix. Estimate the best fitting plane for a set of points in 3D. D is the dimension of the plane. N is the number of points to fit. P is the number type used to represent points.

source
OpticSim.pointMethod

This will return (Inf,Inf,Inf) if the point is at infinity. In this case you probably should be using the direction of the VirtualPoint rather than its position

source
OpticSim.pointMethod
point(ray::AbstractRay{T,N}, alpha::T) -> SVector{T, N}

Returns a point on the ray at origin + alpha * direction. Alpha must be >= 0.

source
OpticSim.pointMethod

returns a 3D point. This takes into account the offset of centerpoint and the rotation vector used to construct the Rectangle. u and v are scaled by the size of the rectangle so that u=0,v=0 is one corner and u=v=1 is the diagonal corner. This function should go away once we have a sensible object transform hierarchy system.

source
OpticSim.pointMethod

returns a 3D point in the plane of the rectangle. This takes into account the offset of centerpoint and the rotation vector used to construct the Rectangle. u and v are scaled by the size of the rectangle so that u=0,v=0 is one corner and u=v=1 is the diagonal corner. This function should go away once we have a sensible object transform hierarchy system.

source
OpticSim.pressureMethod
pressure(system::AbstractOpticalSystem{T}) -> T

Get the pressure of system in Atm.

source
OpticSim.processintersectionMethod
processintersection(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.

source
OpticSim.reset!Method
reset!(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.

source
OpticSim.reversenormalMethod
reversenormal(a::Intersection{T,N})

Used by the CSG complement operator (i.e. -) to reverse the inside outside sense of the object.

source
OpticSim.samplesurfaceMethod
samplesurface(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

source
OpticSim.semidiameterMethod
semidiameter(system::AxisymmetricOpticalSystem{T}) -> T

Get the semidiameter of system, that is the semidiameter of the entrance pupil (i.e. first surface) of the system.

source
OpticSim.snellMethod
snell(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.

source
OpticSim.sphericalangleMethod

returns the spherical angle formed by the cone with centervector at its center with neighbor1,neighbor2 the edges

source
OpticSim.sum!Method
sum!(a::HierarchicalImage{T}, b::HierarchicalImage{T})

Add the contents of b to a in an efficient way.

source
OpticSim.traceMethod
trace(system::AbstractOpticalSystem{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.

source
OpticSim.traceMethod
trace(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.

source
OpticSim.traceMethod
trace(system::AbstractOpticalSystem{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.

source
OpticSim.traceMTMethod
traceMT(system::AbstractOpticalSystem{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.

source
OpticSim.tracehitsMethod
tracehits(system::AbstractOpticalSystem{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 LensTraces which hit the detector.

source
OpticSim.tracehitsMTMethod
tracehitsMT(system::AbstractOpticalSystem{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 LensTraces which hit the detector, accumulated from all threads.

source
OpticSim.transformMethod
transform(surf::CSGGenerator{T}, transform::Transform{T} = identitytransform(T)) -> CSGGenerator{T}

Returns a new CSGGenerator with another transform applied. This is useful if you want multiple copies of a premade CSG structure with different transforms, for example in an MLA.

source
OpticSim.triangulateMethod
triangulate(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.

source
OpticSim.uvMethod
uv(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.

source
OpticSim.uvrangeMethod
uvrange(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 Surfaces which are not ParametricSurfaces (e.g. Rectangle).

source
OpticSim.uvtopixMethod
uvtopix(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 AbstractOpticalSystem.

source
OpticSim.verticesMethod

The vertices of planar shapes are defined in a plane so they are two dimensional. In the local coordinate frame this is the x,y plane, so the implied z coordinate is 0

source
OpticSim.verticesMethod

returns the 2 dimensional vertex points of the shape defining the lens aperture. These points lie in the plane of the shape

source
OpticSim.vertices3dMethod

Returns the vertices of the Hexagon represented in the local coordinate frame. The vertices lie in the z = 0 plane and are 2D

source
OpticSim.virtualdistanceMethod

returns the virtual distance of the point from the lens plane. When |distance| == focallength then virtualdistance = ∞

source
OpticSim.virtualpointMethod

computes the virtual point position corresponding to the input point, or returns nothing for points at infinity. point is specified in the lens coordinate frame

source
OpticSim.αMethod
α(ray::AbstractRay{T,N}, point::SVector{N,T}) -> T

Computes the alpha corresponding to the closest position on the ray to point

source

Geometry

OpticSim.Geometry.TransformMethod
Transform(origin, forward) -> Transform{S}

Returns the Transform of type S (default Float64) representing the local frame with origin and forward direction. the other 2 axes are computed automaticlly.

source
OpticSim.Geometry.TransformMethod
Transform(colx::Vec3{T}, coly::Vec3{T},colz::Vec3{T}, colw::Vec3{T}, ::Type{T} = Float64) where {T<:Real}

Costruct a transform from the input columns.

source
OpticSim.Geometry.TransformMethod
Transform(rotation::AbstractArray{T,2}, translation::AbstractArray{T,1}) where {T<:Real} -> Transform{S}

Returns the Transform of type S (default Float64) created by a rotation matrix (3x3) and translation vector of length 3.

source
OpticSim.Geometry.TransformMethod
Transform(rotation::SMatrix{3,3,T}, translation::SVector{3,T}) where {T<:Real} -> Transform{S}

Returns the Transform of type S (default Float64) created by a rotation matrix and translation vector.

source
OpticSim.Geometry.TransformMethod
Transform(colx::Vec3{T}, coly::Vec3{T},colz::Vec3{T}, colw::Vec3{T}, ::Type{T} = Float64) where {T<:Real}

Costruct a transform from the input columns.

source
OpticSim.Geometry.Vec4Method
Vec4(m::SMatrix{3,N,T} where{N,T<:Real} -> SMatrix{3,N,T})

Input is matrix of 3d points, each column is one point. Returns matrix of 3d points with 1 appended in the last row.

source
OpticSim.Geometry.Vec4Method
Vec4(v::SVector{3, T}) where {T<:Real} -> Vec4{T}

Accept SVector and create a Vec4 type [v[1], v[2], v[3], 1]

source
Base.:*Method

The t and m matrices are allowed to be of different element type. This allows transforming a Unitful matrix for example:

id = identitytransform()
m = fill(1mm,3,4)
id*m #returns a matrix filled with Unitful quantities. If both matrices had to be the same type this would not work
source
Base.:*Method

The t and m matrices are allowed to be of different element type. This allows transforming a Unitful matrix for example: WARNING: this doesn't work. The translation component of the transform matrix has to be in Unitful units but the rotation part has to be in unitless units for this to work. Only works if one assumes that the translation part of the transform implicitly has the same units as the Unitful vectors being transformed. Brittle and likely to cause obscure bugs.

id = identitytransform()
m = fill(1mm,3,4)
id*m #returns a matrix filled with Unitful quantities. If both matrices had to be the same type this would not work
source
OpticSim.Geometry.decomposeRTSMethod
decomposeRTS(tr::Transform{T}) where {T<:Real}

return a touple containing the rotation matrix, the translation vector and the scale vecto represnting the transform.

source
OpticSim.Geometry.forwardMethod
forward(t::Transform{<:Real}) -> Vec3

Assuming t is a 3D rigid transform representing a local left-handed coordinate system, this function will return the third column, representing the "Z" axis.

source
OpticSim.Geometry.local2worldMethod
local2world(t::Transform{T}) where {T<:Real}

return the transform matrix that takes a point in the local coordinate system to the global one

source
OpticSim.Geometry.rightMethod
right(t::Transform{<:Real}) -> Vec3

Assuming t is a 3D rigid transform representing a local left-handed coordinate system, this function will return the first column, representing the "X" axis.

source
OpticSim.Geometry.rotateMethod
rotate(a::Transform{T}, vector::Union{Vec3{T}, SVector{3,T}}) where {T<:Real} -> Vec3{T}

apply the rotation part of the transform a to the vector vector - this operation is usually used to rotate direction vectors.

source
OpticSim.Geometry.rotationMethod
rotation(t::Transform{T}) where {T<:Real} -> SMatrix{3,3,T}

returns the rotation part of the transform t - a 3x3 matrix.

source
OpticSim.Geometry.rotationMethod
rotation([S::Type], θ::T, ϕ::T, ψ::T) -> Transform{S}

Returns the Transform of type S (default Float64) representing the rotation by θ, ϕ and ψ around the x, y and z axes respectively in radians.

source
OpticSim.Geometry.rotationXMethod
rotationX(angle::T) where {T<:Real} -> Transform

Builds a rotation matrix for a rotation around the x-axis. Parameters: The counter-clockwise angle in radians.

source
OpticSim.Geometry.rotationYMethod
rotationY(angle::T) where {T<:Real} -> Transform

Builds a rotation matrix for a rotation around the y-axis. Parameters: The counter-clockwise angle in radians.

source
OpticSim.Geometry.rotationZMethod
rotationZ(angle::T) where {T<:Real} -> Transform

Builds a rotation matrix for a rotation around the z-axis. Parameters: The counter-clockwise angle in radians.

source
OpticSim.Geometry.rotationdMethod
rotationd([S::Type], θ::T, ϕ::T, ψ::T) -> Transform{S}

Returns the Transform of type S (default Float64) representing the rotation by θ, ϕ and ψ around the x, y and z axes respectively in degrees.

source
OpticSim.Geometry.rotmatMethod
rotmat([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.

source
OpticSim.Geometry.rotmatbetweenMethod
rotmatbetween([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, i.e. rotation(a,b) * a = b.

source
OpticSim.Geometry.rotmatdMethod
rotmatd([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.

source
OpticSim.Geometry.upMethod
up(t::Transform{<:Real}) -> Vec3

Assuming t is a 3D rigid transform representing a local left-handed coordinate system, this function will return the second column, representing the "Y" axis.

source
OpticSim.Geometry.world2localMethod
world2local(t::Transform{T}) where {T<:Real}

return the transform matrix that takes a point in the global coordinate system to the local one

source
OpticSim.originMethod
origin(t::Transform{<:Real}) -> Vec3

Assuming t is a 3D rigid transform representing a local left-handed coordinate system, this function will return the fourth column, containing the translation part of the transform in 3D.

source

Zernike

OpticSim.Zernike.OSAtoNMMethod
OSAtoNM(j::Int) -> Tuple{Int, Int}

Convert OSA zernike index j to (N,M) form according to formula J = N * (N + 2) + M.

source
OpticSim.Zernike.δζMethod
δζ(N::Int, M::Int, ρ::T, ϕ::T) -> Tuple{T,T}

Evaluate partial derivatives of Zernike polynomial term $Z_{n}^{m}(\rho, \phi)$.

source
OpticSim.Zernike.ζMethod
ζ(N::Int, M::Int, ρ::T, ϕ::T) -> Tuple{T,T}

Evaluate Zernike polynomial term $Z_{n}^{m}(\rho, \phi)$.

source

QType

OpticSim.QType.SMethod
S(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$.

source
OpticSim.QType.S0Method
S0(coeffs::SVector{NP1,T}, x::T) -> T

Evaluates $\sum_{n=0}^{N}\alpha_n^0Q_n^0(x)$.

source
OpticSim.QType.dS0dxMethod
dS0dx(coeffs::SVector{NP1,T}, x::T) -> T

Evaluates $\frac{\partial}{\partial x}\sum_{n=0}^{N}\alpha_n^0Q_n^0(x)$.

source
OpticSim.QType.dSdxMethod
dSdx(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$.

source

Chebyshev

OpticSim.Chebyshev.TMethod
T(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.

source
OpticSim.Chebyshev.UMethod
U(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.

source
OpticSim.Chebyshev.dTdqMethod
dTdq(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.

source

Examples

OpticSim.Examples.ArizonaEyeMethod
ArizonaEye(::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.

source
OpticSim.Examples.ModelEyeMethod
ModelEye(assembly::LensAssembly{T}, nsamples::Int = 17; pupil_radius::T = 3.0, detpixels::Int = 1000, transform::Transform{T} = identitytransform(T))

Geometrically accurate model of the human eye focused 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.

source
OpticSim.Examples.opticalhemisphereMethod
opticalhemisphere()

Create an optical hemisphere that has optical material properties so it will reflect and refract light. In the previous example the hemisphere object had optical properties of Air, which is the default optical interface, so it won't refract or reflect light.

source