Complete Reference
This page contains what should be a complete list of all docstrings in the OpticSim module, and its submodule.
Index
OpticSim.ChebyshevOpticSim.QTypeOpticSim.ZernikeOpticSim.AbstractOpticalSystemOpticSim.AcceleratedParametricSurfaceOpticSim.AsphericSurfaceOpticSim.AxisymmetricOpticalSystemOpticSim.BSplineCurveOpticSim.BSplineSurfaceOpticSim.BeamStateOpticSim.BezierCurveOpticSim.BezierSurfaceOpticSim.BoundingBoxOpticSim.CSGGeneratorOpticSim.CSGOpticalSystemOpticSim.CSGTreeOpticSim.ChebyshevSurfaceOpticSim.ConvexPolygonOpticSim.CurveTypeOpticSim.CylinderOpticSim.DisjointUnionOpticSim.EllipseOpticSim.Emitters.AngularPower.CosineOpticSim.Emitters.AngularPower.GaussianOpticSim.Emitters.AngularPower.LambertianOpticSim.Emitters.Directions.ConstantOpticSim.Emitters.Directions.HexapolarConeOpticSim.Emitters.Directions.RectGridOpticSim.Emitters.Directions.UniformConeOpticSim.Emitters.Origins.HexapolarOpticSim.Emitters.Origins.PointOpticSim.Emitters.Origins.RectGridOpticSim.Emitters.Origins.RectJitterGridOpticSim.Emitters.Origins.RectUniformOpticSim.Emitters.Sources.CompositeSourceOpticSim.Emitters.Sources.SourceOpticSim.Emitters.Spectrum.DeltaFunctionOpticSim.Emitters.Spectrum.MeasuredOpticSim.Emitters.Spectrum.UniformOpticSim.EmptyIntervalOpticSim.FiniteStopOpticSim.FresnelInterfaceOpticSim.GeometricRayGeneratorOpticSim.Geometry.TransformOpticSim.Geometry.TransformOpticSim.Geometry.TransformOpticSim.Geometry.TransformOpticSim.Geometry.TransformOpticSim.Geometry.TransformOpticSim.Geometry.TransformOpticSim.Geometry.Vec3OpticSim.Geometry.Vec4OpticSim.Geometry.Vec4OpticSim.Geometry.Vec4OpticSim.GridSagInterpolationOpticSim.GridSagSurfaceOpticSim.HexagonOpticSim.HierarchicalImageOpticSim.HologramInterfaceOpticSim.HologramSurfaceOpticSim.InfiniteStopOpticSim.InfinityOpticSim.InterfaceModeOpticSim.IntersectionOpticSim.IntervalOpticSim.KnotVectorOpticSim.LensAssemblyOpticSim.LensTraceOpticSim.NullInterfaceOpticSim.OpticalInterfaceOpticSim.OpticalRayOpticSim.OpticalRayGeneratorOpticSim.ParametricSurfaceOpticSim.ParaxialInterfaceOpticSim.ParaxialLensOpticSim.PlaneOpticSim.QTypeSurfaceOpticSim.RayOpticSim.RayOriginOpticSim.RectangleOpticSim.SphereOpticSim.SphericalCapOpticSim.SphericalPolygonOpticSim.SplineOpticSim.SplineSurfaceOpticSim.SurfaceOpticSim.ThinGratingInterfaceOpticSim.ThinGratingSurfaceOpticSim.TriangleOpticSim.TriangleMeshOpticSim.WrapperSurfaceOpticSim.ZernikeIndexTypeOpticSim.ZernikeSurfaceOpticSim.AnnulusOpticSim.AsphericLensOpticSim.BoundedCylinderOpticSim.CircleOpticSim.CircularApertureOpticSim.CircularApertureOpticSim.ConicLensOpticSim.CuboidOpticSim.Data.ArizonaEyeOpticSim.Data.ModelEyeOpticSim.Data.comfortable_entrance_pupil_translationOpticSim.Data.cornea_to_eyecenterOpticSim.Data.entrancepupil_to_eyecenterOpticSim.Data.eyefocallengthOpticSim.Data.𝐃sdOpticSim.Emitters.applyOpticSim.Emitters.collimatedemitterOpticSim.Emitters.generateOpticSim.Emitters.pointemitterOpticSim.Emitters.visual_sizeOpticSim.EvenAsphericSurfaceOpticSim.FresnelLensOpticSim.Geometry.decomposeRTSOpticSim.Geometry.forwardOpticSim.Geometry.identitytransformOpticSim.Geometry.local2worldOpticSim.Geometry.rightOpticSim.Geometry.rotationOpticSim.Geometry.rotationOpticSim.Geometry.rotationXOpticSim.Geometry.rotationYOpticSim.Geometry.rotationZOpticSim.Geometry.rotationdOpticSim.Geometry.rotmatOpticSim.Geometry.rotmatdOpticSim.Geometry.scaleOpticSim.Geometry.scaleOpticSim.Geometry.scaleOpticSim.Geometry.translationOpticSim.Geometry.translationOpticSim.Geometry.unitW4OpticSim.Geometry.unitX3OpticSim.Geometry.unitX4OpticSim.Geometry.unitY3OpticSim.Geometry.unitY4OpticSim.Geometry.unitZ3OpticSim.Geometry.unitZ4OpticSim.Geometry.upOpticSim.Geometry.world2localOpticSim.HexagonalPrismOpticSim.OddAsphericSurfaceOpticSim.OddEvenAsphericSurfaceOpticSim.RectangularApertureOpticSim.RectangularApertureOpticSim.RectangularPrismOpticSim.SphericalLensOpticSim.SpiderOpticSim.TriangularPrismOpticSim.asphericTypeOpticSim.assemblyOpticSim.closestintersectionOpticSim.closestpointonrayOpticSim.detectorimageOpticSim.distancefromplaneOpticSim.doesintersectOpticSim.insideOpticSim.insideOpticSim.interfaceOpticSim.intersectionsOpticSim.leafOpticSim.makemeshOpticSim.makemeshOpticSim.normalOpticSim.onsurfaceOpticSim.onsurfaceOpticSim.originOpticSim.partialsOpticSim.pointOpticSim.pointOpticSim.pointOpticSim.pointOpticSim.pointOpticSim.pressureOpticSim.reset!OpticSim.resetdetector!OpticSim.samplesurfaceOpticSim.semidiameterOpticSim.sphericalcircleOpticSim.sum!OpticSim.surfaceintersectionOpticSim.surfaceintersectionOpticSim.surfaceintersectionOpticSim.temperatureOpticSim.traceOpticSim.traceOpticSim.traceOpticSim.traceMTOpticSim.tracehitsOpticSim.tracehitsMTOpticSim.triangulateOpticSim.uvOpticSim.uvrangeOpticSim.verticesOpticSim.verticesOpticSim.verticesOpticSim.vertices3dOpticSim.vertices3dOpticSim.virtualpoint
OpticSim
OpticSim.AbstractOpticalSystem — Type
AbstractOpticalSystem{T<:Real}Abstract type for any optical system, must parameterized by the datatype of entities within the system T.
OpticSim.AcceleratedParametricSurface — Type
AcceleratedParametricSurface{T,N,S} <: ParametricSurface{T,N}Wrapper class for ParametricSurfaces where analytical intersection isn't feasible (e.g. ZernikeSurface, ChebyshevSurface). The surface is instead triangulated and an iterative (newton raphson) process carried out to determine precise ray intersection points. S is the type of the ParametricSurface being wrapped.
AcceleratedParametricSurface(surf::ParametricSurface{T,N}, numsamples::Int = 17; interface::NullOrFresnel{T} = nullinterface(T))OpticSim.AsphericSurface — Type
AsphericSurface{T,N,Q,M} <: ParametricSurface{T,N}Surface incorporating an aspheric polynomial - radius, conic and aspherics are defined relative to absolute semi-diameter,. T is the datatype, N is the dimensionality, Q is the number of aspheric terms, and M is the type of aspheric polynomial.
AsphericSurface(semidiameter; radius, conic, aspherics=nothing, normradius = semidiameter)The surface is centered at the origin and treated as being the cap of an infinite cylinder, thus creating a true half-space. Outside of 0 <= ρ <= 1 the height of the surface is not necessarily well defined, so NaN may be returned.
aspherics aspherics should be vectors containing tuples of the form (i, v) where i is the polynomial power of the aspheric term. An empty vector is not permitted. Use nothing instead.
The sag is defined by the equation
\[z(r,\phi) = \frac{cr^2}{1 + \sqrt{1 - (1+k)c^2r^2}} + \sum_{i}^{Q}\alpha_ir^{i}\]
where $\rho = \frac{r}{\texttt{normradius}}$, $c = \frac{1}{\texttt{radius}}$, and $k = \texttt{conic}$ .
The function checks if the aspheric terms are missing, even, odd or both and uses an efficient polynomial evaluation strategy.
OpticSim.AxisymmetricOpticalSystem — Type
AxisymmetricOpticalSystem{T,C<:CSGOpticalSystem{T}} <: AbstractOpticalSystem{T}Optical system which has lens elements and an image detector, created from a DataFrame containing prescription data.
These tags are supported for columns: :Radius, :SemiDiameter, :SurfaceType, :Thickness, :Conic, :Parameters, :Reflectance, :Material.
These tags are supported for entries in a SurfaceType column: Object, Image, Stop. Assumes the Image row will be the last row in the DataFrame.
In practice a CSGOpticalSystem is generated automatically and stored within this system.
AxisymmetricOpticalSystem{T}(
prescription::DataFrame,
detectorpixelsx = 1000,
detectorpixelsy:: = 1000,
::Type{D} = Float32;
temperature = AGFFileReader.TEMP_REF,
pressure = AGFFileReader.PRESSURE_REF
)OpticSim.BSplineCurve — Type
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})OpticSim.BSplineSurface — Type
BSplineSurface{P,S,N,M} <: SplineSurface{P,S,N,M}Curve order is the same in the u and v direction and fixed over all spans. u and v knot vectors are allowed to be different - may change this to make them both the same.
Control points in the u direction correspond to columns, with the lowest value of u corresponding to row 1. Control points in the v direction correspond to rows, with the lowest value of v corresponding to col 1.
BSplineSurface{P,S,N,M}(knots::KnotVector{S}, controlpoints::AbstractArray{<:AbstractArray{S,1},2})
BSplineSurface{P,S,N,M}(uknots::KnotVector{S}, vknots::KnotVector{S}, controlpoints::AbstractArray{<:AbstractArray{S,1},2})OpticSim.BeamState — Type
ConvergingBeam, DivergingBeam or CollimatedBeam, defines the behavior of a beam in a HologramInterface.
OpticSim.BezierCurve — Type
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})OpticSim.BezierSurface — Type
BezierSurface{P,S,N,M} <: SplineSurface{P,S,N,M}Bezier surface defined by grid of control points.
BezierSurface{P,S,N,M}(controlpoints::AbstractArray{<:AbstractArray{S,1},2})OpticSim.BoundingBox — Type
BoundingBox{T<:Real}Axis-aligned three-dimensional bounding box.
BoundingBox(xmin::T, xmax::T, ymin::T, ymax::T, zmin::T, zmax::T)
BoundingBox(s::Surface{T})
BoundingBox(s::ParametricSurface{T,3}, transform::Transform{T} = identitytransform(T))
BoundingBox(c::CSGTree{T})
BoundingBox(tri::Triangle{T})
BoundingBox(triangles::AbstractVector{Triangle{T}})
BoundingBox(points::AbstractArray{SVector{3,T}})
BoundingBox(la::LensAssembly{T})OpticSim.CSGGenerator — Type
CSGGenerator{T<:Real}This is the type you should use when making CSG objects. This type allows for the construction of CSGTree objects with different transforms. When the generator is evaluated, all transforms are propagated down to the LeafNodes and stored there.
Example
a = Cylinder(1.0,1.0)
b = Plane([0.0,0.0,1.0], [0.0,0.0,0.0])
generator = a ∩ b
# now make a csg object that can be ray traced
csgobj = generator(Transform(1.0,1.0,2.0))OpticSim.CSGOpticalSystem — Type
CSGOpticalSystem{T,D<:Real,S<:Surface{T},L<:LensAssembly{T}} <: AbstractOpticalSystem{T}An optical system containing a lens assembly with all optical elements and a detector surface with associated image. The system can be at a specified temperature and pressure.
There are two number types in the type signature. The T type parameter is the numeric type for geometry in the optical system, the D type parameter is the numeric type of the pixels in the detector image. This way you can have Float64 geometry, where high precision is essential, but the pixels in the detector can be Float32 since precision is much less critical for image data, or Complex if doing wave optic simulations.
The detector can be any Surface which implements uv, uvtopix and onsurface, typically this is one of Rectangle, Ellipse or SphericalCap.
CSGOpticalSystem(
assembly::LensAssembly,
detector::Surface,
detectorpixelsx = 1000,
detectorpixelsy = 1000, ::Type{D} = Float32;
temperature = AGFFileReader.TEMP_REF,
pressure = AGFFileReader.PRESSURE_REF
)OpticSim.CSGTree — Type
Abstract type representing any evaluated CSG structure.
OpticSim.ChebyshevSurface — Type
ChebyshevSurface{T,N,P,Q} <: ParametricSurface{T,N}Rectangular surface incorporating Chebyshev polynomials as well as radius and conic terms. T is the datatype, N is the dimensionality, P is the number of Chebyshev terms in u and Q is the number of Chebyshev terms in v.
The surface is centered at the origin and treated as being the cap of an infinite rectangular prism, thus creating a true half-space. Note that the surface is vertically offset so that the center (i.e., (u,v) == (0,0)) lies at 0 on the z-axis.
ChebyshevSurface(halfsizeu, halfsizev, chebycoeff; radius = Inf, conic = 0)chebycoeff is a vector containing tuples of the form (i, j, v) where v is the value of the coefficient $c_{ij}$.
The sag is defined by the equation
\[z(u,v) = \frac{c(u^2 + v^2)^2}{1 + \sqrt{1 - (1+k)c^2(u^2 + v^2)}} + \sum_{i}^{P}\sum_{j}^{Q}c_{ij}T_i(u)T_j(v)\]
where $c = \frac{1}{\texttt{radius}}$, $k = \texttt{conic}$ and $T_n$ is the nᵗʰ Chebyshev polynomial of the first kind.
OpticSim.ConvexPolygon — Type
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.
OpticSim.CurveType — Type
Either Rational or Euclidean, used for Splines and SplineSurfaces.
OpticSim.Cylinder — Type
Cylinder{T,N} <: ParametricSurface{T,N}Cylinder of infinite height centered at the origin, oriented along the z-axis. visheight is used for visualization purposes only, note that this does not fully represent the surface.
Cylinder(radius::T, visheight::T = 2.0; interface::NullOrFresnel{T} = nullinterface(T))OpticSim.DisjointUnion — Type
Datatype representing an ordered series of disjoint intervals on a ray. An arbitrary array of Intervals can be input to the constructor and they will automatically be processed into a valid DisjointUnion (or a single Interval if appropriate).
DisjointUnion(intervals::AbstractVector{Interval{R}})OpticSim.Ellipse — Type
Ellipse{T} <: Surface{T}Elliptical surface, not a valid CSG object. The rotation of the rectangle around its normal is defined by rotationvec. rotationvec×surfacenormal is taken as the vector along the u axis.
Can be used as a detector in AbstractOpticalSystems.
Ellipse(halfsizeu::T, halfsizev::T, [surfacenormal::SVector{3,T}, centrepoint::SVector{3,T}]; interface::NullOrFresnel{T} = nullinterface(T))The minimal case returns an ellipse centered at the origin with surfacenormal = [0, 0, 1].
OpticSim.EmptyInterval — Type
EmptyInterval{T} <: AbstractRayInterval{T}An interval with no Intersections which is also not infinite.
EmptyInterval(T = Float64)
EmptyInterval{T}()OpticSim.FiniteStop — Type
FiniteStop{T,P<:StopShape,Q<:StopShape} <: Surface{T}Stop surface with finite extent. P refers to the shape of the aperture and Q represents the shape of the bounds of the stop surface.
OpticSim.FresnelInterface — Type
FresnelInterface{T} <: OpticalInterface{T}Interface between two materials with behavior defined according to the Fresnel equations, with a specified reflectance and transmission. Assumes unpolarized light.
FresnelInterface{T}(insidematerial, outsidematerial; reflectance = 0, transmission = 1, interfacemode = ReflectOrTransmit)The interfacemode can be used to trace rays deterministically. Valid values are defined in the InterfaceMode enum. Reflect means that all values are reflected, Transmit means that all values are transmitted. ReflectOrTransmit will randomly reflect and transmit rays with the distribution given by the reflection and transmission arguments. This is also the default. In all cases the power recorded with the ray is correctly updated. This can be used to fake sequential raytracing. For example a beamsplitter surface may be set to either Reflect or Transmit to switch between the two outgoing ray paths.
OpticSim.GeometricRayGenerator — Type
GeometricRayGenerator{T,O<:RayOriginGenerator{T}} <: AbstractRayGenerator{T}Generates geometric Rays according to the specific implementation of the subclass.
OpticSim.GridSagInterpolation — Type
Either GridSagLinear or GridSagBicubic - determines the interpolation between sample points in the grid for a GridSagSurface.
OpticSim.GridSagSurface — Type
GridSagSurface{T,N,S<:Union{ZernikeSurface{T,N},ChebyshevSurface{T,N}},Nu,Nv} <: ParametricSurface{T,N}Either a Zernike (circular) or Chebyshev (rectangular) surface with grid sag height added to the base sag. The surface shape is determined by either a linear or a bicubic spline interpolation of the Nu×Nv grid of sag values, set by the interpolation argument taking either GridSagLinear or GridSagBicubic.
Each entry in the grid is a vector of the form $[z, \frac{\partial z}{\partial x}, \frac{\partial z}{\partial y}, \frac{\partial^2 z}{\partial x \partial y}]$. The first data item corresponds to the lower left corner of the surface, that is, the corner defined by the -u and -v limit. Each point that follows is read across the face of the surface from left to right moving upwards. If zero is given for the partials (and using bicubic interpolation) then the partials will be approximated using finite differences.
The sag grid can be decentered from the surface in uv space, if so the surface may become wild outside of the area over which the grid is defined. It is advised to clip the surface to the valid area using CSG operations in this case.
A surface can also be generated from a .GRD file by passing in the filename as the first and only positional argument. In this case the surface will be rectangular with optional radius and conic.
See docs for ZernikeSurface and ChebyshevSurface for details of the base surface.
GridSagSurface(basesurface::Union{ZernikeSurface{T,N},ChebyshevSurface{T,N}}, sag_grid::AbstractArray{T,3}; interpolation = GridSagBicubic, decenteruv = (0, 0))
GridSagSurface{T}(filename::String; radius = Inf, conic = 0, interpolation = GridSagBicubic)OpticSim.Hexagon — Type
Hexagon{T} <: Surface{T}Hexagonal surface, not a valid CSG object. The rotation of the hexagon around its normal is defined by rotationvec. rotationvec×surfacenormal is taken as the vector along the u axis.
Hexagon(side_length::T, [surfacenormal::SVector{3,T}, centrepoint::SVector{3,T}]; rotationvec::SVector{3,T} = [0.0, 1.0, 0.0], interface::NullOrFresnel{T} = nullinterface(T))The minimal case returns a rectangle centered at the origin with surfacenormal = [0, 0, 1].
OpticSim.HierarchicalImage — Type
HierarchicalImage{T<:Number} <: AbstractArray{T,2}Image type which dynamically allocated memory for pixels when their value is set, the value of unset pixels is assumed to be zero.
This is used for the detector image of AbstractOpticalSystems which can typically be very high resolution, but often have a large proportion of the image blank.
OpticSim.HologramInterface — Type
HologramInterface{T} <: OpticalInterface{T}Interface representing a thick hologram (though geometrically thin). The efficiency, η, is calculated using Kogelnik's coupled wave theory so is only valid for the first order. If the zero order is included then it has efficiency 1 - η. Also assumes that the HOE was recorded under similar conditions to the playback conditions, thickness is in microns.
BeatState arguments can be one of ConvergingBeam, DivergingBeam and CollimatedBeam. In the first two cases signalpointordir and referencepointordir are 3D point in global coordinate space. For CollimatedBeam they are normalized direction vectors.
For reference, see:
- Coupled Wave Theory for Thick Hologram Gratings - H Kogelnik, 1995
- Sequential and non-sequential simulation of volume holographic gratings - M Kick et al, 2018
HologramInterface(signalpointordir::SVector{3,T}, signalbeamstate::BeamState, referencepointordir::SVector{3,T}, referencebeamstate::BeamState, recordingλ::T, thickness::T, beforematerial, substratematerial, aftermaterial, signalrecordingmaterial, referencerecordingmaterial, RImodulation::T, include0order = false)OpticSim.HologramSurface — Type
HologramSurface{T,S} <: WrapperSurface{T,S}Surface type for use with HologramInterface.
HologramSurface(surface::Surface{T}, interface::HologramInterface{T})OpticSim.InfiniteStop — Type
InfiniteStop{T,P<:StopShape} <: Surface{T}Stop surface with infinite extent (outside of the aperture). P refers to the shape of the aperture.
OpticSim.Infinity — Type
Infinity{T} <: IntervalPoint{T}Point representing ∞ within an Interval.
Infinity(T = Float64)
Infinity{T}()OpticSim.InterfaceMode — Type
Valid modes for deterministic raytracing
OpticSim.Intersection — Type
Intersection{T,N} <: IntervalPoint{T}Represents the point at which an Ray hits a Surface. This consists of the distance along the ray, the intersection point in world space, the normal in world space, the UV on the surface and the OpticalInterface hit.
Has the following accessor methods:
point(a::Intersection{T,N}) -> SVector{N,T}
normal(a::Intersection{T,N}) -> SVector{N,T}
uv(a::Intersection{T,N}) -> SVector{2,T}
u(a::Intersection{T,N}) -> T
v(a::Intersection{T,N}) -> T
α(a::Intersection{T,N}) -> T
interface(a::Intersection{T,N}) -> OpticalInterface{T}
flippednormal(a::Intersection{T,N}) -> BoolOpticSim.Interval — Type
Interval{T} <: AbstractRayInterval{T}Datatype representing an interval between two IntervalPoints on a ray.
The lower element can either be RayOrigin or an Intersection. The upper element can either be an Intersection or Infinity.
positivehalfspace(int::Intersection) -> Interval with lower = int, upper = Infinity
rayorigininterval(int::Intersection) -> Interval with lower = RayOrigin, upper = int
Interval(low, high)Has the following accessor methods:
lower(a::Interval{T}) -> Union{RayOrigin{T},Intersection{T,3}}
upper(a::Interval{T}) -> Union{Intersection{T,3},Infinity{T}}OpticSim.KnotVector — Type
KnotVector{T<:Number}Vector to define knots used for BSplineCurve and BSplineSurface.
OpticSim.LensAssembly — Type
LensAssembly{T<:Real}Structure which contains the elements of the optical system, these can be CSGTree or Surface objects.
In order to prevent type ambiguities bespoke structs are created for each possible number of elements e.g. LensAssembly3. These are parameterized by the types of the elements to prevent ambiguities. Basic surface types such as Rectangle (which can occur in large numbers) are stored independently in Vectors, so type parameters are only needed for CSG objects.
Each struct looks like this:
struct LensAssemblyN{T,T1,T2,...,TN} <: LensAssembly{T}
axis::SVector{3,T}
rectangles::Vector{Rectangle{T}}
ellipses::Vector{Ellipse{T}}
hexagons::Vector{Hexagon{T}}
paraxials::Vector{ParaxialLens{T}}
E1::T1
E2::T2
...
EN::TN
endWhere Ti <: Union{Surface{T},CSGTree{T}}.
To create a LensAssembly object the following functions can be used:
LensAssembly(elements::Vararg{Union{Surface{T},CSGTree{T},LensAssembly{T}}}; axis = SVector(0.0, 0.0, 1.0)) where {T<:Real}OpticSim.LensTrace — Type
LensTrace{T<:Real,N}Contains an intersection point and the ray segment leading to it from within an optical trace. The ray carries the path length, power, wavelength, number of intersections and source number, all of which are accessible directly on this class too.
Has the following accessor methods:
ray(a::LensTrace{T,N}) -> OpticalRay{T,N}
intersection(a::LensTrace{T,N}) -> Intersection{T,N}
power(a::LensTrace{T,N}) -> T
wavelength(a::LensTrace{T,N}) -> T
pathlength(a::LensTrace{T,N}) -> T
point(a::LensTrace{T,N}) -> SVector{N,T}
uv(a::LensTrace{T,N}) -> SVector{2,T}
sourcenum(a::LensTrace{T,N}) -> Int
nhits(a::LensTrace{T,N}) -> IntOpticSim.NullInterface — Type
NullInterface{T} <: OpticalInterface{T}Interface which will be ignored totally by any rays, used only in construction of CSG objects.
NullInterface(T = Float64)
NullInterface{T}()OpticSim.OpticalInterface — Type
OpticalInterface{T<:Real}Any subclass of OpticalInterface must implement the following:
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}
insidematerial(i::OpticalInterface{T}) -> AGFFileReader.AbstractGlass
outsidematerial(i::OpticalInterface{T}) -> AGFFileReader.AbstractGlass
reflectance(i::OpticalInterface{T}) -> T
transmission(i::OpticalInterface{T}) -> TSee documentation for processintersection for details.
OpticSim.OpticalRay — Type
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}) -> IntOpticSim.OpticalRayGenerator — Type
OpticalRayGenerator{T} <: AbstractRayGenerator{T}Generates OpticalRays according to the specific implementation of the subclass.
OpticSim.ParametricSurface — Type
ParametricSurface{T,N} <: Surface{T}T is the number type used to represent the surface, e.g., Float64. N is the dimension of the space the surface is embedded in. ParametricSurfaces are valid CSG objects, in some cases (where analytic intersection isn't possible) they must be wrapped in an AcceleratedParametricSurface for use.
Must implement the following:
uv(surface::ParametricSurface{T,N}, p::SVector{N,T}) -> SVector{2,T}
uvrange(surface::ParametricSurface{T,N}) -> Tuple{Tuple{T,T},Tuple{T,T}}
point(surface::ParametricSurface{T,N}, u::T, v::T) -> SVector{N,T}
partials(surface::ParametricSurface{T,N}, u::T, v::T) -> Tuple{SVector{N,T}, SVector{N,T}}
normal(surface::ParametricSurface{T,N}, u::T, v::T) -> SVector{N,T}
inside(surface::ParametricSurface{T,N}, p: :SVector{N,T}) -> Bool
onsurface(surface::ParametricSurface{T,N}, p::SVector{N,T}) -> Bool
surfaceintersection(surface::ParametricSurface{T,N}, AbstractRay::Ray{T,N}) -> Union{EmptyInterval{T},Interval{T},DisjointUnion{T}}
interface(surface::ParametricSurface{T,N}) -> OpticalInterface{T}OpticSim.ParaxialInterface — Type
ParaxialInterface{T} <: OpticalInterface{T}Interface describing an idealized planar lens, i.e. one that is thin and with no aberrations.
In general this interface should not be constructed directly, the ParaxialLensEllipse and ParaxialLensRect functions should be used to create a ParaxialLens object directly.
ParaxialInterface(focallength::T, centroid::SVector{3,T}, outsidematerial::Y)OpticSim.ParaxialLens — Type
ParaxialLens{T} <: Surface{T}surfacenormal is the output direction of the lens. Paraxial lens cannot act as the interface between two materials, hence only a single outside material is specified, by default Air.
Create with the following functions
ParaxialLensEllipse(focaldistance, halfsizeu, halfsizev, surfacenormal, centrepoint; rotationvec = [0.0, 1.0, 0.0], outsidematerial = AGFFileReader.Air, decenteruv = (0.0, 0.0))
ParaxialLensRect(focaldistance, halfsizeu, halfsizev, surfacenormal, centrepoint; rotationvec = [0.0, 1.0, 0.0], outsidematerial = AGFFileReader.Air, decenteruv = (0.0, 0.0))
ParaxialLensHex(focaldistance, side_length, surfacenormal, centrepoint; rotationvec = [0.0, 1.0, 0.0], outsidematerial = AGFFileReader.Air, decenteruv = (0.0, 0.0))
ParaxialLensConvexPoly(focaldistance, local_frame, local_polygon_points, local_center_point; outsidematerial = AGFFileReader.Air)OpticSim.Plane — Type
Plane{T,N} <: ParametricSurface{T,N}Infinite planar surface where the positive normal side is outside the surface.
By default this will not create any geometry for visualization, the optional vishalfsizeu and vishalfsizev arguments can be used to draw the plane as a rectangle for visualization note that this does not fully represent the surface. In this case, the rotation of the rectangle around the normal to the plane is defined by visvec - surfacenormal×visvec is taken as the vector along the u axis.
Plane(surfacenormal::SVector{N,T}, pointonplane::SVector{N,T}; interface::NullOrFresnel{T} = nullinterface(T), vishalfsizeu::T = 0.0, vishalfsizev::T = 0.0, visvec::SVector{N,T} = [0.0, 1.0, 0.0])
Plane(nx::T, ny::T, nz::T, x::T, y::T, z::T; interface::NullOrFresnel{T} = nullinterface(T), vishalfsizeu::T = 0.0, vishalfsizev::T = 0.0, visvec::SVector{N,T} = [0.0, 1.0, 0.0])OpticSim.QTypeSurface — Type
QTypeSurface{T,D,M,N} <: ParametricSurface{T,D}Surface incorporating the QType polynomials - radius and conic are defined relative to absolute semi-diameter, QType terms are normalized according to the normradius parameter. T is the datatype, D is the dimensionality, M and N are the maximum QType terms used.
The surface is centered at the origin and treated as being the cap of an infinite cylinder, thus creating a true half-space. Outside of 0 <= ρ <= 1 the height of the surface is not necessarily well defined, so NaN may be returned.
QTypeSurface(semidiameter; radius = Inf, conic = 0.0, αcoeffs = nothing, βcoeffs = nothing, normradius = semidiameter)αcoeffs and βcoeffs should be a vector of tuples of the form (m, n, v) where v is the value of the coefficient $α_n^m$ or $β_n^m$ respectively.
The sag is defined by the equation
\[\begin{aligned} z(r,\phi) = & \frac{cr^2}{1 + \sqrt{1 - (1+k)c^2r^2}} + \frac{\sqrt{1 + kc^2r^2}}{\sqrt{1-(1+k)c^2r^2}} \cdot \\ & \left\{ \rho^2(1-\rho^2)\sum_{n=0}^{N}\alpha_n^0 Q_n^0 (\rho^2) + \sum_{m=1}^{M}\rho^m\sum_{n=0}^N \left[ \alpha_n^m\cos{m\phi} +\beta_n^m\sin{m\phi}\right]Q_n^m(\rho^2) \right\} \end{aligned}\]
where $\rho = \frac{r}{\texttt{normradius}}$, $c = \frac{1}{\texttt{radius}}$, $k = \texttt{conic}$ and $Q_n^m$ is the QType polynomial index $m$, $n$.
OpticSim.Ray — Type
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}OpticSim.RayOrigin — Type
RayOrigin{T} <: IntervalPoint{T}Point representing 0 within an Interval, i.e. the start of the ray.
RayOrigin(T = Float64)
RayOrigin{T}()OpticSim.Rectangle — Type
Rectangle{T} <: Surface{T}Rectangular surface, not a valid CSG object. The rotation of the rectangle around its normal is defined by rotationvec. rotationvec×surfacenormal is taken as the vector along the u axis.
Can be used as a detector in AbstractOpticalSystems.
Rectangle(halfsizeu::T, halfsizev::T, [surfacenormal::SVector{3,T}, centrepoint::SVector{3,T}]; rotationvec::SVector{3,T} = [0.0, 1.0, 0.0], interface::NullOrFresnel{T} = nullinterface(T))The minimal case returns a rectangle centered at the origin with surfacenormal = [0, 0, 1].
OpticSim.Sphere — Type
Sphere{T,N} <: ParametricSurface{T,N}Spherical surface centered at the origin.
Sphere(radius::T = 1.0; interface::NullOrFresnel{T} = nullinterface(T))OpticSim.SphericalCap — Type
SphericalCap{T} <: ParametricSurface{T}Spherical cap surface, creates a half-space which is essentially the subtraction of a sphere from an infinite plane. Only the spherical cap itself is visualized, not the plane. The positive normal side is outside the surface.
Can be used as a detector in AbstractOpticalSystems.
SphericalCap(radius::T, ϕmax::T, [surfacenormal::SVector{3,T}, centrepoint::SVector{3,T}]; interface::NullOrFresnel{T} = nullinterface(T))The minimal case returns a spherical cap centered at the origin with surfacenormal = [0, 0, 1].
OpticSim.SphericalPolygon — Type
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.
OpticSim.Spline — Type
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
OpticSim.SplineSurface — Type
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.
OpticSim.Surface — Type
Surface{T<:Real}T is the number type used to represent the surface, e.g., Float64. Basic Surfaces are not valid CSG objects, they function only in a stand-alone capacity.
Must implement the following:
surfaceintersection(surface::Surface{T}, ray::AbstractRay{T,3}) -> Union{EmptyInterval{T},Interval{T}}
normal(surface::Surface{T}) -> SVector{3,T}
interface(surface::Surface{T}) -> OpticalInterface{T}
makemesh(surface::Surface{T}) -> TriangleMesh{T}In a conventional ray tracer the surface intersection function would only return the first surface the ray intersects. Because our ray tracer does CSG operations the surface intersection function intersects the ray with all leaf surfaces which are part of the CSG tree.
Each leaf surface returns one or more 1D intervals along the ray. These intervals contain the part of the ray which is inside the surface. The intervals computed at the leaves are propagated upward through the CSG tree and the CSG operations of union, intersection, and difference are applied to generate new intervals which are themselves propagated upward.
The result is a union of 1D intervals, which may be disjoint, a single interval, or empty. The union of intervals represents the parts of the ray which are inside the CSG object.
Inside is well defined for halfspaces such as cylinders and spheres which divide space into two parts, but not for Bezier or NURBS patches which generally do not enclose a volume. For surfaces which are not halfspaces the notion of inside is defined locally by computing the angle between the incoming ray and the normal of the surface at the point of intersection. All surfaces must be defined so that the normal points to the outside of the surface.
A negative dot product between the incoming ray and the normal indicates the ray is coming from the outside of the surface and heading toward the inside. A positive dot product indicates the ray is coming from the inside of the surface and heading toward the outside.
Intervals are defined along the ray which is being intersected with the surface, so they are one dimensional. For example, assume we have a ray with origin o on the outside of a plane and an intersection with the plane at point int = o + td where t is a scalar and d is the unit direction of the ray. The inside interval will be (Intersection(t),Infinity). This interval begins at the intersection point on the plane and continues to positive infinity. The Intersection struct stores both the parametric value t and the 3D point of intersection to make various operations more efficient. But the interval operations only depend on the parametric value t.
If the origin o is on the inside of the plane then the inside interval will be (RayOrigin,Intersection(t)). Only the part of the ray from the ray origin to the intersection point is inside the plane.
It is the programmer's responsibility to return Interval results from surfaceintersection that maintain these properties.
The following must be implemented only if the surface is being used as a detector
uv(surface::Surface{T}, p::SVector{3,T}) -> SVector{2,T}
uvtopix(surface::Surface{T}, uv::SVector{2,T}, imsize::Tuple{Int,Int}) -> Tuple{Int,Int}
onsurface(surface::Surface{T}, p::SVector{3,T}) -> BoolOpticSim.ThinGratingInterface — Type
ThinGratingInterface{T} <: OpticalInterface{T}Interface representing an idealized thin grating. period is in microns, vector should lie in the plane of the surface. Transmission and reflectance can be specified for an arbitrary number of orders up to 10, selected using the maxorder and minorder parameters. If nothing then reflectance is assumed to be 0 and transmission is assumed to be 1.
ThinGratingInterface(vector, period, insidematerial, outsidematerial; maxorder = 1, minorder = -1, reflectance = nothing, transmission = nothing)OpticSim.ThinGratingSurface — Type
ThinGratingSurface{T,S} <: WrapperSurface{T,S}Surface type for use with ThinGratingInterface.
ThinGratingSurface(surface::Surface{T}, interface::ThinGratingInterface{T})OpticSim.Triangle — Type
Triangle{T} <: Surface{T}Triangular surface, not a valid CSG object. Primarily used as a component part of TriangleMesh or to enable intersection of AcceleratedParametricSurfaces. Can never be used directly as an optical surface as it doesn't have an OpticalInterface.
Triangle(v1::SVector{3,T}, v2::SVector{3,T}, v3::SVector{3,T}, [uv1::SVector{2,T}, uv2::SVector{2,T}, uv3::SVector{2,T}])OpticSim.TriangleMesh — Type
TriangleMesh{T} <: Surface{T}An array of Triangles forming a mesh. Used for visualization purposes only.
TriangleMesh(tris::Vector{Triangle{T}})OpticSim.WrapperSurface — Type
WrapperSurface{T,S<:Surface{T}} <: Surface{T}A generic surface type which serves as a basis for extension of Surfaces for custom OpticalInterface subclasses. Essentially just forwards all Surface and ParametricSurface methods to a field of the WrapperSurface named surface. Also provides a generic implementation of surfaceintersection which tests for an intersection with the underlying surface and returns either an EmptyInterval or a half space (never a closed interval).
OpticSim.ZernikeIndexType — Type
Either ZernikeIndexingOSA or ZernikeIndexingNoll, see Zernike polynomials wikipedia entry for details.
OpticSim.ZernikeSurface — Type
ZernikeSurface{T,N,P,Q,M} <: ParametricSurface{T,N}Surface incorporating the Zernike polynomials - radius, conic and aspherics are defined relative to absolute semi-diameter, Zernike terms are normalized according to the normradius parameter. T is the datatype, N is the dimensionality, P is the number of Zernike terms, Q is the number of aspheric terms and M is the Aspheric Type.
The surface is centered at the origin and treated as being the cap of an infinite cylinder, thus creating a true half-space. Outside of 0 <= ρ <= 1 the height of the surface is not necessarily well defined, so NaN may be returned.
For convenience the input zcoeff can be indexed using either OSA or Noll convention, indicated using the indexing argument as either ZernikeIndexingOSA or ZernikeIndexingNoll.
ZernikeSurface(semidiameter, radius = Inf, conic = 0, zcoeff = nothing, aspherics = nothing, normradius = semidiameter, indexing = ZernikeIndexingOSA)zcoeff and aspherics should be vectors containing tuples of the form (i, v) where i is either the index of the Zernike term for the corresponding indexing, or the polynomial power of the aspheric term (may be even or odd) and v is the corresponding coefficient $A_i$ or $\alpha_i$ respectively.. M will be determined from the terms entered to optimize the evaluation of the aspheric polynomial.
The sag is defined by the equation
\[z(r,\phi) = \frac{cr^2}{1 + \sqrt{1 - (1+k)c^2r^2}} + \sum_{i}^{Q}\alpha_ir^{2i} + \sum_{i}^PA_iZ_i(\rho, \phi)\]
where $\rho = \frac{r}{\texttt{normradius}}$, $c = \frac{1}{\texttt{radius}}$, $k = \texttt{conic}$ and $Z_n$ is the nᵗʰ Zernike polynomial.
OpticSim.Annulus — Method
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}.
OpticSim.AsphericLens — Method
AsphericLens(insidematerial, frontvertex, frontradius, frontconic, frontaspherics, backradius, backconic, backaspherics, thickness, semidiameter; lastmaterial = AGFFileReader.Air, nextmaterial = AGFFileReader.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.
OpticSim.BoundedCylinder — Method
BoundedCylinder(radius::T, height::T; interface::NullOrFresnel{T} = nullinterface(T)) -> CSGGenerator{T}Create a cylinder with planar caps on both ends centered at (0, 0, 0) with axis (0, 0, 1).
OpticSim.Circle — Method
Circle(radius, [surfacenormal, centrepoint]; interface = nullinterface(T))Shortcut method to create a circle. The minimal case returns a circle centered at the origin with normal = [0, 0, 1].
OpticSim.CircularAperture — Method
CircularAperture(radius::T, surfacenormal::SVector{3,T}, centrepoint::SVector{3,T})Creates a circular aperture in a plane i.e. InfiniteStop{T,CircularStopShape}.
OpticSim.CircularAperture — Method
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.
OpticSim.ConicLens — Method
ConicLens(insidematerial, frontvertex, frontradius, frontconic, backradius, backconic, thickness, semidiameter; lastmaterial = AGFFileReader.Air, nextmaterial = AGFFileReader.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.
OpticSim.Cuboid — Method
Cuboid(halfsizex::T, halfsizey::T, halfsizez::T; interface::NullOrFresnel{T} = nullinterface(T)) -> CSGGenerator{T}Create a cuboid centered at (0, 0, 0).
OpticSim.EvenAsphericSurface — Method
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
OpticSim.FresnelLens — Method
Code does not work.
FresnelLens(insidematerial, frontvertex, radius, thickness, semidiameter, groovedepth; conic = 0.0, aspherics = nothing, outsidematerial = AGFFileReader.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.HexagonalPrism — Method
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.
OpticSim.OddAsphericSurface — Method
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
OpticSim.OddEvenAsphericSurface — Method
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
OpticSim.RectangularAperture — Method
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.
OpticSim.RectangularAperture — Method
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.
OpticSim.RectangularPrism — Method
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.
OpticSim.SphericalLens — Method
SphericalLens(insidematerial, frontvertex, frontradius, backradius, thickness, semidiameter; lastmaterial = AGFFileReader.Air, nextmaterial = AGFFileReader.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.
OpticSim.Spider — Method
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:
| _|_
/ \ |OpticSim.TriangularPrism — Method
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.
OpticSim.asphericType — Method
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.
OpticSim.assembly — Method
assembly(system::AbstractOpticalSystem{T}) -> LensAssembly{T}Get the LensAssembly of system.
OpticSim.closestintersection — Function
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.
OpticSim.closestpointonray — Method
closestpointonray(r::Ray{T,N}, point::SVector{N,T}) -> SVector{T,NReturns the point on the ray closest to point.
OpticSim.detectorimage — Method
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.
OpticSim.distancefromplane — Method
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.
OpticSim.doesintersect — Method
doesintersect(bbox::BoundingBox{T}, r::AbstractRay{T,3}) -> BoolTests whether r intersects an axis-aligned BoundingBox, bbox.
OpticSim.inside — Method
inside(obj::CSGTree{T}, point::SVector{3,T}) -> Bool
inside(obj::CSGTree{T}, x::T, y::T, z::T) -> BoolTests whether a 3D point in world space is inside obj.
OpticSim.inside — Method
inside(surf::ParametricSurface{T}, p::SVector{3,T}) -> Bool
inside(surf::ParametricSurface{T}, x::T, y::T, z::T) -> BoolTests whether a 3D point in world space is inside surf.
OpticSim.interface — Method
interface(surf::Surface{T}) -> OpticalInterface{T}Return the OpticalInterface associated with surf.
OpticSim.intersections — Method
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
OpticSim.leaf — Method
leaf(surf::ParametricSurface{T}, transform::Transform{T} = identitytransform(T)) -> CSGGenerator{T}Create a leaf node from a parametric surface with a given transform.
OpticSim.makemesh — Method
makemesh(poly::ConvexPolygon{N, T}, ::Int = 0) where {N, T<:Real} -> TriangleMeshCreate 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.
OpticSim.makemesh — Method
makemesh(object, subdivisions::Int = 30) -> TriangleMeshCreates a TriangleMesh from an object, either a ParametricSurface, CSGTree or certain surfaces (e.g. Circle, Rectangle). This is used for visualization purposes only.
OpticSim.normal — Method
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.
OpticSim.onsurface — Method
onsurface(obj::CSGTree{T}, point::SVector{3,T}) -> Bool
onsurface(obj::CSGTree{T}, x::T, y::T, z::T) -> BoolTests whether a 3D point in world space is on the surface (i.e. shell) of obj.
OpticSim.onsurface — Method
onsurface(surf::ParametricSurface{T}, p::SVector{3,T}) -> Bool
onsurface(surf::ParametricSurface{T}, x::T, y::T, z::T) -> BoolTests whether a 3D point in world space is on surf.
OpticSim.partials — Method
partials(surf::ParametricSurface{T}, u::T, v::T) -> (SVector{3,T}, SVector{3,T})
partials(surf::ParametricSurface{T}, uv::SVector{2,T}) -> (SVector{3,T}, SVector{3,T})Returns a tuple of the 3D partial derivatives of surf with respect to u and v at the given uv coordinate.
OpticSim.point — Method
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
OpticSim.point — Method
point(ray::AbstractRay{T,N}, alpha::T) -> SVector{T, N}Returns a point on the ray at origin + alpha * direction. Alpha must be >= 0.
OpticSim.point — Method
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.
OpticSim.point — Method
point(surf::ParametricSurface{T}, u::T, v::T) -> SVector{3,T}
point(surf::ParametricSurface{T}, uv::SVector{2,T}) -> SVector{3,T}Returns the 3D point on surf at the given uv coordinate.
OpticSim.point — Method
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.
OpticSim.pressure — Method
pressure(system::AbstractOpticalSystem{T}) -> TGet the pressure of system in Atm.
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.
OpticSim.resetdetector! — Method
resetdetector!(system::AbstractOpticalSystem{T})Reset the deterctor image of system to zero.
OpticSim.samplesurface — Method
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
OpticSim.semidiameter — Method
semidiameter(system::AxisymmetricOpticalSystem{T}) -> TGet the semidiameter of system, that is the semidiameter of the entrance pupil (i.e. first surface) of the system.
OpticSim.sphericalcircle — Function
creates a circular polygon that subtends a half angle of θ
OpticSim.sum! — Method
sum!(a::HierarchicalImage{T}, b::HierarchicalImage{T})Add the contents of b to a in an efficient way.
OpticSim.surfaceintersection — Method
surfaceintersection(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.surfaceintersection — Method
surfaceintersection(obj::CSGTree{T}, r::AbstractRay{T,N})Calculates the intersection of r with CSG object, obj.
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.
The ray is intersected with the LeafNodes that make up the CSG object and the resulting Intervals and DisjointUnions are composed with the same boolean operations to give a final result. The ray is transformed by the inverse of the transform associated with the leaf node to put it in object space for that node before the intersection is carried out, typically this object space is centered at the origin, but may differ for each primitive.
Some intersections are culled without actually evaluating them by first checking if the ray intersects the BoundingBox of each node in the CSGTree, this can substantially improve performance in some cases.
OpticSim.surfaceintersection — Method
surfaceintersection(bbox::BoundingBox{T}, r::AbstractRay{T,3}) -> Union{EmptyInterval{T},Interval{T}}Calculates the intersection of r with an axis-aligned BoundingBox, bbox.
Returns an EmptyInterval if there is no intersection or an Interval if there is one or two intersections. Note that the uv of the returned intersection is always 0.
OpticSim.temperature — Method
temperature(system::AbstractOpticalSystem{T}) -> TGet the temperature of system in °C.
OpticSim.trace — Method
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.
OpticSim.trace — Method
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.
OpticSim.trace — Method
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.
OpticSim.traceMT — Method
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.
OpticSim.tracehits — Method
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.
OpticSim.tracehitsMT — Method
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.
OpticSim.triangulate — Method
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. determining 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.uv — Method
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.
OpticSim.uvrange — Method
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).
OpticSim.vertices — Method
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
OpticSim.vertices — Method
returns the 2 dimensional vertex points of the shape defining the lens aperture. These points lie in the plane of the shape
OpticSim.vertices — Method
returns the 2D vertices in the plane of the rectangle
OpticSim.vertices3d — Method
Returns the vertices of the Hexagon represented in the local coordinate frame. The vertices lie in the z = 0 plane and are 2D
OpticSim.vertices3d — Method
returns the vertices of the rectangle in 3D
OpticSim.virtualpoint — Method
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
Optical Emitters
OpticSim.Emitters.apply — Function
apply(???)[TODO] Returns ray power
OpticSim.Emitters.collimatedemitter — Method
collimatedemitter(origin::AbstractVector{T}, halfsquaresize; λ::Length = 500nm, numrays = 100) where {T<:Real}Creates a square collimated emitter, emitting rays in the -z direction. Rays are emitted on a square grid with sqrt(numrays) on a side. λ can be a unitful quantity, e.g., 550nm, or a number. In the latter case the units are implicitly microns.
OpticSim.Emitters.generate — Function
generate(???)[TODO]
OpticSim.Emitters.pointemitter — Method
pointemitter(origin::AbstractVector{T}, coneangle; λ::Length = 500nm, numrays = 100) where {T<:Real}Creates a point source with Lambertian emission power and cone distribution of rays, emitting in the -z direction. λ is a unitful Length quantity, e.g., 550nm.
OpticSim.Emitters.visual_size — Function
visual_size(???)[TODO]
OpticSim.Emitters.AngularPower.Cosine — Type
Cosine{T} <: AbstractAngularPowerDistribution{T}Cosine power distribution. Ray power is calculated by:
power = power * (cosine_angle ^ cosine_exp)
OpticSim.Emitters.AngularPower.Gaussian — Type
Gaussian{T} <: AbstractAngularPowerDistribution{T}GGaussian power distribution. Ray power is calculated by:
power = power * exp(-(gaussianu * l^2 + gaussianv * m^2)) where l and m are the cos_angles between the two axes respectively.
OpticSim.Emitters.AngularPower.Lambertian — Type
Lambertian{T} <: AbstractAngularPowerDistribution{T}Ray power is unaffected by angle.
OpticSim.Emitters.Directions.Constant — Type
Constant{T} <: AbstractDirectionDistribution{T}Encapsulates a single ray direction, where the default direction is unitZ3 [0, 0, 1].
Constant(direction::Vec3{T}) where {T<:Real}
Constant(::Type{T} = Float64) where {T<:Real}OpticSim.Emitters.Directions.HexapolarCone — Type
HexapolarCone{T} <: AbstractDirectionDistribution{T}Rays are generated by sampling a cone with θmax angle in an hexapolar fashion. The number of rays depends on the requested rings and is computed using the following formula: 1 + round(Int64, (nrings * (nrings + 1) / 2) * 6)
HexapolarCone(direction::Vec3{T}, θmax::T, nrings::Int64) where {T<:Real}
HexapolarCone(θmax::T, nrings::Int64 = 3) where {T<:Real}OpticSim.Emitters.Directions.RectGrid — Type
RectGrid{T} <: AbstractDirectionDistribution{T}Encapsulates a single ray direction, where the default direction is unitZ3 [0, 0, 1].
Constant(direction::Vec3{T}) where {T<:Real}
Constant(::Type{T} = Float64) where {T<:Real}OpticSim.Emitters.Directions.UniformCone — Type
UniformCone{T} <: AbstractDirectionDistribution{T}Encapsulates numsamples rays sampled uniformly from a cone with max angle θmax.
UniformCone(direction::Vec3{T}, θmax::T, numsamples::Int64) where {T<:Real}
UniformCone(θmax::T, numsamples::Int64) where {T<:Real}OpticSim.Emitters.Origins.Hexapolar — Type
Hexapolar{T} <: AbstractOriginDistribution{T}Encapsulates an ellipse (or a circle where halfsizeu=halfsizev) sampled in an hexapolar fashion (rings).
Hexapolar(nrings::Int64, halfsizeu::T, halfsizev::T) where {T<:Real} OpticSim.Emitters.Origins.Point — Type
Point{T} <: AbstractOriginDistribution{T}Encapsulates a single point origin.
Point(position::Vec3{T}) where {T<:Real}
Point(x::T, y::T, z::T) where {T<:Real}
Point(::Type{T} = Float64) where {T<:Real}OpticSim.Emitters.Origins.RectGrid — Type
RectGrid{T} <: AbstractOriginDistribution{T}Encapsulates a rectangle sampled in a grid fashion.
RectGrid(width::T, height::T, usamples::Int64, vsamples::Int64) where {T<:Real} OpticSim.Emitters.Origins.RectJitterGrid — Type
RectJitterGrid{T} <: AbstractOriginDistribution{T}Encapsulates a rectangle sampled in a grid fashion with jitter.
RectGrid(width::T, height::T, ures::Int64, vres::Int64, samplesPerRegion::Int64) where {T<:Real} OpticSim.Emitters.Origins.RectUniform — Type
RectUniform{T} <: AbstractOriginDistribution{T}Encapsulates a uniformly sampled rectangle with user defined number of samples.
RectUniform(width::T, height::T, samples_count::Int64) where {T<:Real}OpticSim.Emitters.Sources.CompositeSource — Type
CompositeSource{T} <: AbstractSource{T}This data-type represents the composite emitter (Source) which is constructed with a list of basic or composite emitters and a 3D Transform.
CompositeSource(transform::Transform{T}, sources::Vector{<:AbstractSource}) where {T<:Real} OpticSim.Emitters.Sources.Source — Type
Source{T<:Real, Tr<:Transform{T}, S<:Spectrum.AbstractSpectrum{T}, O<:Origins.AbstractOriginDistribution{T}, D<:Directions.AbstractDirectionDistribution{T}, P<:AngularPower.AbstractAngularPowerDistribution{T}} <: AbstractSource{T}This data-type represents the basic emitter (Source), which is a combination of a Spectrum, Angular Power Distribution, Origins and Directions distribution and a 3D Transform.
Source(::Type{T} = Float64;
transform::Tr = Transform(),
spectrum::S = Spectrum.Uniform(),
origins::O = Origins.Point(),
directions::D = Directions.Constant(),
power::P = AngularPower.Lambertian(),
sourcenum::Int64 = 0) where {
Tr<:Transform,
S<:Spectrum.AbstractSpectrum,
O<:Origins.AbstractOriginDistribution,
D<:Directions.AbstractDirectionDistribution,
P<:AngularPower.AbstractAngularPowerDistribution,
T<:Real}
Source(transform::Tr, spectrum::S, origins::O, directions::D, power::P, ::Type{T} = Float64; sourcenum::Int64 = 0) where {
Tr<:Transform,
S<:Spectrum.AbstractSpectrum,
O<:Origins.AbstractOriginDistribution,
D<:Directions.AbstractDirectionDistribution,
P<:AngularPower.AbstractAngularPowerDistribution,
T<:Real}OpticSim.Emitters.Spectrum.DeltaFunction — Type
DeltaFunction{T} <: AbstractSpectrum{T}Encapsulates a constant spectrum.
DeltaFunction{T<:Real}OpticSim.Emitters.Spectrum.Measured — Type
Measured{T} <: AbstractSpectrum{T}Encapsulates a measured spectrum to compute emitter power. Create spectrum by reading CSV files. Assumes spectrum samples are evenly spaced - exception otherwise. Evaluate spectrum at arbitrary wavelength with spectrumpower (more technical details coming soon)
Measured(samples::DataFrame)OpticSim.Emitters.Spectrum.Uniform — Type
Uniform{T} <: AbstractSpectrum{T}Encapsulates a flat spectrum range which is sampled uniformly. Unless stated diferrently, the range used will be 450nm to 680nm.
Uniform(low_end::T, high_end::T) where {T<:Real}
Uniform(::Type{T} = Float64) where {T<:Real}Geometry
OpticSim.Geometry.Transform — Type
Transform{S<:Real}Transform encapsulating rotation, translation and scale in 3D space. Translation happens after rotation.
Transform{S}(θ::T, ϕ::T, ψ::T, x::T, y::T, z::T)
Transform(rotation::SMatrix{3,3,S}, translation::SVector{3,S})
Transform(rotation::AbstractArray{S,2}, translation::AbstractArray{S,1})θ, ϕ and ψ in first constructor are in radians.
OpticSim.Geometry.Transform — Method
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 automatically.
OpticSim.Geometry.Transform — Method
Transform(colx::Vec3{T}, coly::Vec3{T},colz::Vec3{T}, colw::Vec3{T}, ::Type{T} = Float64) where {T<:Real}Construct a transform from the input columns.
OpticSim.Geometry.Transform — Method
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.
OpticSim.Geometry.Transform — Method
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.
OpticSim.Geometry.Transform — Method
Transform(colx::Vec3{T}, coly::Vec3{T},colz::Vec3{T}, colw::Vec3{T}, ::Type{T} = Float64) where {T<:Real}Construct a transform from the input columns.
OpticSim.Geometry.Transform — Method
Transform([S::Type]) -> Transform{S}Returns the Transform of type S (default Float64) representing the identity transform.
OpticSim.Geometry.Vec3 — Type
Vec3{T} provides an immutable vector of fixed length 3 and type T.
Vec3 defines a series of convenience constructors, so you can just type e.g. Vec3(1, 2, 3) or Vec3([1.0, 2.0, 3.0]). It also supports comprehensions, and the zeros(), ones(), fill(), rand() and randn() functions, such as Vec3(rand(3)).
OpticSim.Geometry.Vec4 — Type
Vec4{T} provides an immutable vector of fixed length 4 and type T.
Vec4 defines a series of convenience constructors, so you can just type e.g. Vec3(1, 2, 3, 4) or Vec3([1.0, 2.0, 3.0, 4.0]). It also supports comprehensions, and the zeros(), ones(), fill(), rand() and randn() functions, such as Vec4(rand(4)).
OpticSim.Geometry.Vec4 — Method
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.
OpticSim.Geometry.Vec4 — Method
Vec4(v::SVector{3, T}) where {T<:Real} -> Vec4{T}Accept SVector and create a Vec4 type [v[1], v[2], v[3], 1]
OpticSim.Geometry.decomposeRTS — Method
decomposeRTS(tr::Transform{T}) where {T<:Real}return a tuple containing the rotation matrix, the translation vector and the scale vecto represnting the transform.
OpticSim.Geometry.forward — Method
forward(t::Transform{<:Real}) -> Vec3Assuming t is a 3D rigid transform representing a local left-handed coordinate system, this function will return the third column, representing the "Z" axis.
OpticSim.Geometry.identitytransform — Method
identitytransform([S::Type]) -> Transform{S}
Returns the Transform of type S (default Float64) representing the identity transform.
OpticSim.Geometry.local2world — Method
local2world(t::Transform{T}) where {T<:Real}return the transform matrix that takes a point in the local coordinate system to the global one
OpticSim.Geometry.right — Method
right(t::Transform{<:Real}) -> Vec3Assuming t is a 3D rigid transform representing a local left-handed coordinate system, this function will return the first column, representing the "X" axis.
OpticSim.Geometry.rotation — Method
rotation(t::Transform{T}) where {T<:Real} -> SMatrix{3,3,T}returns the rotation part of the transform t - a 3x3 matrix.
OpticSim.Geometry.rotation — Method
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.
OpticSim.Geometry.rotationX — Method
rotationX(angle::T) where {T<:Real} -> TransformBuilds a rotation matrix for a rotation around the x-axis. Parameters: The counter-clockwise angle in radians.
OpticSim.Geometry.rotationY — Method
rotationY(angle::T) where {T<:Real} -> TransformBuilds a rotation matrix for a rotation around the y-axis. Parameters: The counter-clockwise angle in radians.
OpticSim.Geometry.rotationZ — Method
rotationZ(angle::T) where {T<:Real} -> TransformBuilds a rotation matrix for a rotation around the z-axis. Parameters: The counter-clockwise angle in radians.
OpticSim.Geometry.rotationd — Method
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.
OpticSim.Geometry.rotmat — Method
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.
OpticSim.Geometry.rotmatd — Method
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.
OpticSim.Geometry.scale — Method
scale(s::T) where {T<:Real}Creates a uniform scaling transform
OpticSim.Geometry.scale — Method
scale(t::Vec3{T}) where {T<:Real}Creates a scaling transform
OpticSim.Geometry.scale — Method
scale(x::T, y::T, z::T) where {T<:Real}Creates a scaling transform
OpticSim.Geometry.translation — Method
translation(x::T, y::T, z::T) where {T<:Real}Creates a translation transform
OpticSim.Geometry.translation — Method
translation(x::T, y::T, z::T) where {T<:Real}Creates a translation transform
OpticSim.Geometry.unitW4 — Method
returns the unit vector [0, 0, 0, 1]
OpticSim.Geometry.unitX3 — Method
returns the unit vector [1, 0, 0]
OpticSim.Geometry.unitX4 — Method
returns the unit vector [1, 0, 0, 0]
OpticSim.Geometry.unitY3 — Method
returns the unit vector [0, 1, 0]
OpticSim.Geometry.unitY4 — Method
returns the unit vector [0, 1, 0, 0]
OpticSim.Geometry.unitZ3 — Method
returns the unit vector [0, 0, 1]
OpticSim.Geometry.unitZ4 — Method
returns the unit vector [0, 0, 1, 0]
OpticSim.Geometry.up — Method
up(t::Transform{<:Real}) -> Vec3Assuming t is a 3D rigid transform representing a local left-handed coordinate system, this function will return the second column, representing the "Y" axis.
OpticSim.Geometry.world2local — Method
world2local(t::Transform{T}) where {T<:Real}return the transform matrix that takes a point in the global coordinate system to the local one
OpticSim.origin — Method
origin(t::Transform{<:Real}) -> Vec3Assuming 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.
Zernike
OpticSim.Zernike — Module
Module to enclose Zernike polynomial specific functionality.
QType
OpticSim.QType — Module
Module to enclose QType polynomial specific functionality. For reference see:
Chebyshev
OpticSim.Chebyshev — Module
Module to enclose Chebyshev polynomial specific functionality.
Data
OpticSim.Data.ArizonaEye — Method
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.
OpticSim.Data.ModelEye — Method
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.
OpticSim.Data.comfortable_entrance_pupil_translation — Method
Average one sided translation of the entrance pupil associated with comfortable eye rotation. If you are using this to define an eyebox multiply this value by 2
OpticSim.Data.cornea_to_eyecenter — Method
distance from vertex of cornea to center of rotation
OpticSim.Data.entrancepupil_to_eyecenter — Method
distance from entrance pupil to center of rotation.
OpticSim.Data.eyefocallength — Method
Posterior focal length, i.e., optical distance from entrance pupil to the retina. Focal length will change depending on accommodation. This value is for focus at ∞. When the eye is focused at 25cm focal length will be ≈ 22mm. Because the index of refraction of the vitreous humor is approximately 1.33 the physical distance from the entrance pupil to the retina will be 24mm/1.33 = 18mm.
OpticSim.Data.𝐃sd — Method
computes pupil diameter as a function of scene luminance L, in cd/m², and the angular area, a, over which this luminance is presented to the eye.