Index of all fuctions in AGFFileReader

AGFFileReader.Air โ€” Constant

Special glass to represent air. Refractive index is defined to always be 1.0 for any temperature and pressure (other indices are relative to this).

source
AGFFileReader.Glass โ€” Type

Stores all attributes relating to a glass type specified in an .AGF glass catalog.

Never used directly, instead created using catalog glasses, e.g. AGFFileReader.SCHOTT.N_BK7.

In order to prevent type ambiguities in OpticSim.jl we can't have this type paramaterized.

source
AGFFileReader.GlassID โ€” Type

Object identifying a glass, containing a type (e.g. MODEL, MIL, OTHER or AGF) depending on how the glass is defined, and an integer ID. Air is AIR:0, others are on the form AGF:N, for example.

source
AGFFileReader.absairindex โ€” Method
absairindex(wavelength; temperature=20ยฐC, pressure=1Atm)

Compute the absolute refractive index of air at wavelength, optionally at specified temperature and pressure. If unitless, arguments are interpretted as ฮผm, ยฐC and Atm respectively.

Examples

julia> absairindex(700u"nm")
1.000271074905147

julia> absairindex(0.7, temperature=27.0)
1.000264738846504

julia> absairindex(532u"nm", temperature = 25u"ยฐC", pressure = 1.3)
1.0003494991178161
source
AGFFileReader.absorption โ€” Method
absorption(glass::AbstractGlass, wavelength; temperature=20ยฐC, pressure=1Atm)

Compute the intensity absorption per mm of glass at wavelength, optionally at specified temperature and pressure. Transmission values are linearly interpolated from the adjacent values in the data table of glass, if wavelength is below the minimum or above the maximum in the table then the nearest value is taken.

Absorption is defined as $\frac{-\log(t)}{\tau}$ where $t$ is the transmission value and $\tau$ is the thickness, both of which are provided in the data table.

If unitless, arguments are interpretted as ฮผm, ยฐC and Atm respectively.

Examples

julia> absorption(AGFFileReader.SUMITA.LAK7, 700u"nm")
0.0006018072325563021

julia> absorption(AGFFileReader.SCHOTT.N_BK7, 0.55, temperature = 22.0)
0.00016504471175660636

julia> absorption(AGFFileReader.SCHOTT.PSK3, 532u"nm", temperature = 25u"ยฐC", pressure = 1.3)
0.00020855284788532435
source
AGFFileReader.add_agf โ€” Method
add_agf(agffile; agfdir = AGF_DIR, sourcefile = SOURCES_PATH, name = nothing, rebuild = true)

Copies a file at agffile (this can be either a download link or local path) to agfdir and appends a corresponding entry to the source list at sourcefile.

If a name is not provided for the catalog, an implicit name is derived from agffile.

If rebuild is true, Pkg.build is called at the end to install the new catalog.

source
AGFFileReader.agffile_to_catalog โ€” Method

Parse a agffile (.agf) into a native Dict, catalogdict, where each kvp = (glassname, glassinfo) is a glass.

1234567891011
NMraw namedispform?NdVd[exclude_substatusmeltfreq]
EDTCE?pฮ”PgF[ignorethermalexp]
CDC1C2C3C4C5C6C7C8C9C10
TDDโ‚€Dโ‚Dโ‚‚Eโ‚€Eโ‚ฮปโ‚œโ‚–temp
ODrelcostCRFRSRARPR
LDฮปminฮปmax
ITT1T2T3
source
AGFFileReader.download_source โ€” Function
download_source(dest::AbstractString, url::AbstractString, POST_data::Union{Nothing,AbstractString} = nothing)

Download and unzip an AGF glass catalog from a publicly available source. Supports POST requests.

source
AGFFileReader.drawglassmap โ€” Method
drawglassmap(glasscatalog::Module; ฮป::Length = 550nm, glassfontsize::Integer = 3, showprefixglasses::Bool = false)

Draw a scatter plot of index vs dispersion (the derivative of index with respect to wavelength). Both index and dispersion are computed at wavelength ฮป.

Choose glasses to graph using the glassfilterprediate argument. This is a function that receives a Glass object and returns true if the glass should be graphed.

If showprefixglasses is true then glasses with names like F_BAK7 will be displayed. Otherwise glasses that have a leading letter prefix followed by an underscore, such as F_, will not be displayed.

The index formulas for some glasses may give incorrect results if ฮป is outside the valid range for that glass. This can give anomalous results, such as indices less than zero or greater than 6. To filter out these glasses set maximumindex to a reasonable value such as 3.0.

example: plot only glasses that do not contain the strings "E" and "J"

drawglassmap(NIKON,showprefixglasses = true,glassfilterpredicate = (x) -> !occursin("J",string(x)) && !occursin("E",string(x)))

source
AGFFileReader.findglass โ€” Method
findglass(condition::Function) -> Vector{Glass}

Returns the list of glasses which satisfy condition where condition::(Glass -> Bool).

Example

julia> findglass(x -> (x.Nd > 2.3 && x.ฮปmin < 0.5 && x.ฮปmax > 0.9))
8-element Array{AGFFileReader.Glass,1}:
 BIREFRINGENT.TEO2_E
 BIREFRINGENT.PBMOO4
 BIREFRINGENT.LINBO3
 INFRARED.CLEARTRAN_OLD
 INFRARED.CLEARTRAN
 INFRARED.SRTIO3
 INFRARED.ZNS_BROAD
 INFRARED.ZNS_VIS
source
AGFFileReader.generate_jls โ€” Method
generate_jls(sourcenames::Vector{<:AbstractString}, mainfile::AbstractString, jldir::AbstractString, agfdir::AbstractString; test::Bool = false)

Generates .jl files in jldir: a mainfile and several catalog files.

Each catalog file is a module representing a distinct glass catalog (e.g. NIKON, SCHOTT), generated from corresponding AGF files in agfdir. These are then included and exported in mainfile.

In order to avoid re-definition of constants AGF_GLASS_NAMES and AGF_GLASSES during testing, we have an optional test argument. If true, we generate a .jl file that defines glasses with GlassType.TEST to avoid namespace/ID clashes.

source
AGFFileReader.glasscatalogs โ€” Method
glasscatalogs()

Returns the complete list of glass catalogs available from AGFFileReader.

Example

julia> glasscatalogs()
41-element Array{Any,1}:
 OpticSim.AGFFileReader.AMTIR
 OpticSim.AGFFileReader.ANGSTROMLINK
 OpticSim.AGFFileReader.APEL
 OpticSim.AGFFileReader.ARCHER
 OpticSim.AGFFileReader.ARTON
 OpticSim.AGFFileReader.AUER_LIGHTING
 OpticSim.AGFFileReader.BIREFRINGENT
 โ‹ฎ
source
AGFFileReader.glassfromMIL โ€” Method
glassfromMIL(glasscode::Union{Float64,Int}) -> Glass

Generates a glass object for the given glass code based on U.S. military standard MIL-G-174, see the MIL specification for further details.

The glass code is a six-digit number specifying the glass according to its refractive index Nd at d-light (587.5618nm), and its Abbe number Vd also taken at d-light. The resulting glass code is the value of Nd - 1 rounded to three digits, followed by Vd rounded to three digits, with all decimal points ignored. For example, N_BK7 has Nd = 1.5168 and Vd = 64.17, giving a six-digit glass code of 517642.

For Nd > 1.999 the format 1.123642 can be used representing Nd = 2.123 and Vd = 64.2.

Accuracy is poor given the low precision of the input parameters, the mean error to measured data may be significant. Behavior may differ from other optical simulation tools when using MIL glasses. The approximate dispersion calculation used these glasses is generally only valid for visible wavelengths, in this case a limit of 360nm to 750nm is imposed.

Examples

julia> index(glassfromMIL(517642), 0.5875618)
1.5170003960064509

julia> index(glassfromMIL(1.134642), 0.5875618)
2.1340008686098946
source
AGFFileReader.glassname โ€” Method
glassname(g::Union{AbstractGlass,GlassID})

Get the name (including catalog) of the glass, or glass with this ID.

source
AGFFileReader.glassnames โ€” Method
glassnames(catalog::Module)

Returns the glass names available from a given catalog.

Example

julia> glassnames(AGFFileReader.CARGILLE)
3-element Array{Any,1}:
 "OG0607"
 "OG0608"
 "OG081160"
source
AGFFileReader.glassnames โ€” Method
glassnames()

Returns the glass names available from all catalogs.

Example

julia> glassnames()
6-element Array{Pair{Module,Array{Any,1}},1}:
 OpticSim.AGFFileReader.CARGILLE => ["OG0607", "OG0608", "OG081160"]
     OpticSim.AGFFileReader.HOYA => ["BAC4", "BACD11"  โ€ฆ  "TAFD65"]
    OpticSim.AGFFileReader.NIKON => ["BAF10", "BAF11"  โ€ฆ  "_7054"]
    OpticSim.AGFFileReader.OHARA => ["L_BAL35", "L_BAL35P"  โ€ฆ  "S_TIM8"]
   OpticSim.AGFFileReader.SCHOTT => ["AF32ECO", "BAFN6"  โ€ฆ  "SFL6"]
   OpticSim.AGFFileReader.SUMITA => ["BAF1", "BAF10"  โ€ฆ  "ZNSF8"]
source
AGFFileReader.index โ€” Method
index(glass::AbstractGlass, wavelength; temperature=20ยฐC, pressure=1Atm)

Compute the refractive index of glass at wavelength, optionally at specified temperature and pressure. Result is relative to the refractive index of air at given temperature and pressure.

If unitless, arguments are interpretted as ฮผm, ยฐC and Atm respectively.

This is defined to always equal 1.0 for Air at any temperature and pressure, use absairindex for the absolute refractive index of air at a given temperature and pressure.

Examples

julia> index(AGFFileReader.SUMITA.LAK7, 700u"nm")
1.646494204478318

julia> index(AGFFileReader.SCHOTT.N_BK7, 0.55, temperature = 22.0)
1.51852824383283

julia> index(AGFFileReader.HOYA.FF1, 532u"nm", temperature = 25u"ยฐC", pressure = 1.3)
1.5144848290944655
source
AGFFileReader.info โ€” Method
info([io::IO], glass::AbstractGlass)

Print out all data associated with glass in an easily readable format.

Examples

julia> info(AGFFileReader.RPO.IG4)
ID:                                                AGF:52
Dispersion formula:                                Schott (1)
Dispersion formula coefficients:
     aโ‚€:                                           6.91189161
     aโ‚:                                           -0.000787956404
     aโ‚‚:                                           -4.22296071
     aโ‚ƒ:                                           142.900646
     aโ‚„:                                           -1812.32748
     aโ‚…:                                           7766.33028
Valid wavelengths:                                 3.0ฮผm to 12.0ฮผm
Reference temperature:                              20.0ยฐC
Thermal ฮ”RI coefficients:
     Dโ‚€:                                           3.24e-5
     Dโ‚:                                           0.0
     Dโ‚‚:                                           0.0
     Eโ‚€:                                           0.0
     Eโ‚:                                           0.0
     ฮปโ‚œโ‚–:                                          0.0
TCE (รท1e-6):                                       20.4
Ignore thermal expansion:                          false
Density (p):                                       4.47g/mยณ
ฮ”PgF:                                              0.0
RI at sodium D-Line (587nm):                       1.0
Abbe Number:                                       0.0
Cost relative to N_BK7:                              ?
Status:                                            Standard (0)
Melt frequency:                                    0
Exclude substitution:                              false
source
AGFFileReader.modelglass โ€” Method
modelglass(Nd::Float64, Vd::Float64, ฮ”PgF::Float64) -> Glass

Generates a glass object for the given refractive index at d-light (587.5618nm), Nd, the Abbe number also at d-light, Vd, and partial dispersion, ฮ”PgF. The mean error to measured data for these models is typically small - usually < 0.0001. Behavior may differ from other optical simulation tools when using model glasses.

The approximate dispersion calculation used for these glasses is generally only valid for visible wavelengths, in this case a limit of 360nm to 750nm is imposed.

Examples

julia> index(modelglass(1.5168, 64.17, 0.0), 0.5875618)
1.5168003970108495

julia> index(modelglass(1.2344, 61.57, 0.003), 0.678)
1.2329425902693352
source
AGFFileReader.plot_indices โ€” Method
plot_indices(glass::AbstractGlass; polyfit=false, fiterror=false, degree=5, temperature=20ยฐC, pressure=1Atm, nsamples=300, sampling_domain="wavelength")

Plot the refractive index for glass for nsamples within its valid range of wavelengths, optionally at temperature and pressure. polyfit will show a polynomial of optionally specified degree fitted to the data, fiterror will also show the fitting error of the result. sampling_domain specifies whether the samples will be spaced uniformly in "wavelength" or "wavenumber".

source
AGFFileReader.polyfit_indices โ€” Method
polyfit_indices(wavelengths, n_rel; degree=5)

Fit a polynomial to indices at wavelengths, optionally specifying the degree of the polynomial. Returns tuple of array of fitted indices at wavelengths and the polynomial.

source
AGFFileReader.verify_source โ€” Method
verify_source(agffile::AbstractString, expected_sha256sum::AbstractString)

Verify a source file using SHA256, returning true if successful. Otherwise, remove the file and return false.

source
AGFFileReader.verify_sources! โ€” Method
verify_sources!(sources::AbstractVector{<:AbstractVector{<:AbstractString}}, agfdir::AbstractString)

Verify a list of sources located in agfdir. If AGF files are missing or invalid, try to download them using the information provided in sources.

Each source โˆˆ sources is a collection of strings in the format name, sha256sum, url [, POST_data], where the last optional string is used to specify data to be sent in a POST request. This allows us to download a greater range of sources (e.g. SUMITA).

Modifies sources in-place such that only verified sources remain.

source