Anime API

Contents

Index

Public API

I/O

Anime.readobsconfigFunction
readobsconfig(configfile::String)

Read YAML configuration file with observation settings in a Dictionary.

source
Anime.readstationinfoFunction
readstationinfo(stationfile::String; delim::String=",", ignorerepeated::Bool=false)

Read additional information about stations from CSV file into a DataFrame.

source
Anime.readbandpassinfoFunction
readbandpassinfo(bandpassfile::String; delim::String=",", ignorerepeated::Bool=false)

Read bandpass information from CSV file into a DataFrame.

source
Anime.createmsfromconfigFunction
createmsfromconfig(msname::String, mscreationmode::String, casaanttemplate::String, stationinfo::DataFrame, spw_centrefreq::Array{Float64, 1}, 
spw_bw::Array{Float64, 1}, spw_channels::Array{Int64, 1}, sourcedict::Dict{String, Any}, starttime::String, exposure::Float64, scans::Int64,
scanlengths::Array{Float64, 1}, scanlag::Float64; autocorr::Bool=false, telescopename::String="VLBA", feed::String="perfect R L",
shadowlimit::Float64=1e-6, elevationlimit::String="10deg", stokes::String="RR RL LR LL", delim::String=",", ignorerepeated::Bool=false)

Create Measurement Set from observation parameters and DataFrame containing station information. Requires casatools.

source
createmsfromconfig(obsconfig, stationinfo::DataFrame)

Alias with fewer arguments for use in pipelines.

source
Anime.createmsfromuvfitsFunction
createmsfromuvfits(uvfits::String, msname::String, mscreationmode::String; overwrite::Bool=true)

Convert uvfits to MS. Requires casatools and casatasks.

source
Anime.createuvfitsfrommsFunction
createuvfitsfromms(msname::String, uvfits::String, datacolumn::String; field::String="", spw::String="", antenna::String="",
timerange::String="", overwrite::Bool=true)

Create UVFITS from existing MS. Requires casatasks.

source
Anime.writemsFunction
writems(ms::MeasurementSet; h5file::String="")

Write data to MS. Optionally, add weight columns from h5file.

source
Anime.makecasaantennatableFunction
makecasaantennatable(df::DataFrame, casaanttemplate::String; delim::String=",", ignorerepeated::Bool=false)

Create CASA antenna table from the DataFrame df containing station information using casaanttemplate as template. Requires casatools.

source
Anime.addweightcolumnsFunction
addweightcolumns(msname::String, mode::String, sigmaspec::Bool, weightspec::Bool)

Add WEIGHT_SPECTRUM and SIGMA_SPECTRUM columns to existing MS. Requires casatools to be installed.

source

Coherency

Anime.run_wscleanFunction
run_wsclean(msname::String, fitsdir::String, polarized::Bool, channelgroups::Int64, osfactor::Int64)

Compute source coherency matrix using WSClean and populate MS. fitsdir is the directory containing FITS source models, polarized indicates if source model is polarized, and channelgroups is the number of images in frequency. osfactor is the oversampling factor which is a WSClean argument that controls prediction/imaging accuracy.

source

Instrument Models

Anime.troposphere!Function
troposphere!(obs::CjlObservation, h5file::String; absorptionfile::String="", dispersivefile::String="", elevfile::String="")

Umbrella function for tropospheric model. Model is written to HDF5 file.

source
Anime.instrumentalpolarization!Function
instrumentalpolarization!(data::Array{Complex{Float32},4}, scanno::Vector{Int32}, times::Vector{Float64},
stationinfo::DataFrame, phasedir::Array{Float64,2}, pos::Array{Float64, 2}, chanfreqvec::Vector{Float64}, polframe::String,
polmode::String, antenna1::Vector{Int32}, antenna2::Vector{Int32}, exposure::Float64, corruptseed::Int; h5file::String="",
elevfile::String="", parangfile::String="")

Compute frequency-varying instrumental polarization (leakage, or "D-Jones" terms) and apply to data; write model to HDF5 file.

source
instrumentalpolarization!(ms::MeasurementSet, stationinfo::DataFrame, obsconfig::Dict; h5file::String="",
elevfile::String="", parangfile::String="")

Alias for use in pipelines.

source
Anime.pointing!Function
pointing!(data::Array{Complex{Float32},4}, stationinfo::DataFrame, scanno::Vector{Int32}, chanfreqvec::Vector{Float64}, 
ptgint::Float64, ptgscale::Float64, ptgmode::String, exposure::Float64, times::Vector{Float64}, corruptseed::Int, antenna1::Vector{Int32}, 
antenna2::Vector{Int32}, numchan::Int64; α=1.0, h5file::String="")

Compute pointing model and apply to data; write model to HDF5 file.

source
pointing!(ms::MeasurementSet, stationinfo::DataFrame, obsconfig::Dict; h5file::String="")

Alias with fewer arguments for use in pipelines.

source
Anime.stationgains!Function
stationgains!(data::Array{Complex{Float32},4}, scanno::Vector{Int32}, times::Vector{Float64}, exposure::Float64, 
stationinfo::DataFrame, mode::String, corruptseed::Int, antenna1::Vector{Int32}, antenna2::Vector{Int32}, numchan::Int64; h5file::String="")

Compute time-variable complex station gains and apply to data; write model to HDF5 file.

source
stationgains!(ms::MeasurementSet, stationinfo::DataFrame, obsconfig::Dict; h5file::String="")

Alias for use in pipelines.

source
Anime.bandpass!Function
bandpass!(data::Array{Complex{Float32},4}, stationinfo::DataFrame, bandpassinfo::DataFrame, corruptseed::Int,
antenna1::Vector{Int32}, antenna2::Vector{Int32}, numchan::Int64, chanfreqvec::Vector{Float64}; h5file::String="")

Compute the complex bandpass model and apply to data; write model to HDF5 file.

source
bandpass!(ms::MeasurementSet, stationinfo::DataFrame, bandpassinfo::DataFrame, obsconfig::Dict; h5file::String="")

Alias for use in pipelines.

source
Anime.thermalnoise!Function
thermalnoise!(data::Array{Complex{Float32},4}, times::Vector{Float64}, antenna1::Vector{Int32}, antenna2::Vector{Int32},
correff::Float64, exposure::Float64, chanwidth::Float64, corruptseed::Int, sefd::Vector{Float64}; h5file::String="",
noisefile::String="")

Compute per-baseline thermal noise using radiometer equation and apply to data. The actual numerical values are serialized in HDF5 format.

source
thermalnoise!(ms::MeasurementSet, stationinfo::DataFrame, obsconfig::Dict; h5file::String="", noisefile::String="")

Alias for use in pipelines.

source

Statistics

Anime.genseries1d!Function
genseries1d!(series::Vector{ComplexF32}, mode::String, location::ComplexF32, scale::Float32, driftrate::Float32, nsamples::Int64, rng::AbstractRNG)

Generate a complex-valued Gaussian process 1-D series of length nsamples with the given location, scale, and driftrate parameters. mode determines if "normal" distribution or "gaussian processes" is used to generate the samples.

source
genseries1d!(series::Vector{Float32}, mode::String, location::Float32, scale::Float32, driftrate::Float32, nsamples::Int64, rng::AbstractRNG)

Generate a Float32-valued Gaussian process 1-D series of length nsamples with the given location, scale, and driftrate parameters. mode determines if "normal" distribution or "gaussian processes" is used to generate the samples.

source
genseries1d!(series::Vector{Float64}, mode::String, location::Float64, scale::Float64, driftrate::Float64, nsamples::Int64, rng::AbstractRNG)

Generate a Float64-valued Gaussian process 1-D series of length nsamples with the given location, scale, and driftrate parameters. mode determines if "normal" distribution or "gaussian processes" is used to generate the samples.

source
genseries1d!(series, times, rng::AbstractRNG; μ=0.0, σ=1.0, ρ=1.0)

Generate 1-D series using SE kernel.

source
genseries1d!(series, times, rng::AbstractRNG, α::Float64; μ=0.0, σ=1.0, ρ=2.0)

Generate 1-D series using RQ kernel.

source

Utils

Anime.elevationangleFunction
elevationangle(times::Vector{Float64}, phasedir::Array{Float64,2}, stationinfo::DataFrame, pos::Array{Float64, 2})

Compute elevation angle for all stations in stationinfo for all times using phasedir and pos information.

source
Anime.parallacticangleFunction
parallacticangle(times::Vector{Float64}, phasedir::Array{Float64,2}, stationinfo::DataFrame, pos::Array{Float64,2})

Compute parallactic angle for all stations in stationinfo for all times using phasedir and pos information.

source

Plotting

Anime.plotvisFunction
plotvis(uvw::Matrix{Float64}, chanfreqvec::Array{Float64,1}, flag::Array{Bool,4}, data::Array{Complex{Float32},4},
numchan::Int64, times::Vector{Float64}; saveprefix="data_")

Plot complex visibilities against time and projected baseline length.

source
Anime.plotuvcovFunction
plotuvcov(uvw::Matrix{Float64}, flagrow::Vector{Bool}, chanfreqvec::Vector{Float64}; saveprefix="test_")

Plot uv-coverage of observation.

source
Anime.plotstationgainsFunction
plotstationgains(h5file::String, scanno::Vector{Int32}, times::Vector{Float64}, exposure::Float64, stationnames::Vector{String3})

Plot complex station gains against time.

source
Anime.plotbandpassFunction
plotbandpass(h5file::String, stationnames::Vector{String3}, chanfreqvec::Vector{Float64})

Plot bandpass gains against time.

source
Anime.plotpointingerrorsFunction
plotpointingerrors(h5file::String, scanno::Vector{Int32}, stationnames::Vector{String3})

Plot pointing errors.

source
Anime.plotelevationangleFunction
plotelevationangle(h5file::String, scanno::Vector{Int32}, times::Vector{Float64}, stationnames::Vector{String3})

Plot evolution of station elevation angles during the course of the observation.

source
Anime.plotparallacticangleFunction
plotparallacticangle(h5file::String, scanno::Vector{Int32}, times::Vector{Float64}, stationnames::Vector{String3})

Plot evolution of station parallactic angles during the course of the observation.

source
Anime.plotdtermsFunction
plotdterms(h5file::String, stationnames::Vector{String3}, chanfreqvec::Vector{Float64})

Plot frequency-dependent complex instrumental polarization.

source
Anime.plottransmissionFunction
plottransmission(h5file::String, stationnames::Vector{String3}, times::Vector{Float64}, chanfreqvec::Vector{Float64})

Plot tropospheric transmission variation with frequency.

source
Anime.plotmeandelaysFunction
plotmeandelays(h5file::String, stationnames::Vector{String3}, times::Vector{Float64}, chanfreqvec::Vector{Float64})

Plot mean delays against time.

source

Data Types

Anime.MeasurementSetType
struct MeasurementSet{T}

Composite type for storing data from a Measurement Set.

Fields

  • name: Measurement Set name
  • data: Complex visibility data (MS DATA column)
  • flag: Flag array of the same dimensions as data (MS FLAG column)
  • flagrow: Flag row Boolean vector (MS FLAG_ROW column)
  • antenna1: Antenna 1 in a baseline pair (MS ANTENNA1 column)
  • antenna2: Antenna 2 in a baseline pair (MS ANTENNA2 column)
  • uvw: uvw coordinates (MS UVW column)
  • times: Timestamps (MS TIME column)
  • exposure: Integration time
  • scanno: Scan numbers (SCAN_NUMBER column in MS)
  • numchan: Number of frequency channels
  • chanfreqvec: Channel frequencies (Hz) (MS CHAN_FREQ column)
  • chanwidth: Width of frequency channel
  • phasedir: Direction of phase centre
  • pos: Antenna positions in x,y,z coordinates
source

Internal

Anime.computeweights!Function
computeweights!(totalrmsspec::Array{Float32, 4}, totalwtspec::Array{Float32, 4}; h5file::String="")

Compute total rms (sigma) values and inverse-squared visibility weights from thermal+sky noise terms stored in h5file.

source
Anime.run_aatmFunction
run_aatm(ms::MeasurementSet, stationinfo::DataFrame; absorptionfile::String="", dispersivefile::String="")::DataFrame

Run AATM (Bjona Nikolic; Pardo et al. 2001) to compute absorption by and dispersive delay in the troposphere. If AATM is not installed, this function can still accept input absorption and dispersion values in a specific CSV format and populate atm.csv.

source
Anime.compute_transmission!Function
compute_transmission!(transmission::Array{Float64, 3}, ms::MeasurementSet, stationinfo::DataFrame, obsconfig::Dict, atmdf::DataFrame,
elevationmatrix::Array{Float64, 2}, g::HDF5.Group)::Array{Float64, 3}

Compute elevation-dependent (mean) tropospheric transmission given opacity τ and elevation angle θ for each station.

\[e^{-τ/\sin{\theta_{\rm el}}}\]

source
Anime.attenuate!Function
attenuate!(ms::MeasurementSet, transmission::Array{Float64, 3})

Attenuate the signal as it passes through (mean) troposphere using precomputed transmission values.

\[I = I_0 e^{-τ/\sin{\theta_{\rm el}}}\]

source
Anime.compute_skynoise!Function
compute_skynoise!(ms::MeasurementSet, stationinfo::DataFrame, obsconfig::Dict, atmdf::DataFrame, transmission::Array{Float64, 3}, g::HDF5.Group)

Compute sky contribution to visibility noise using the radiometer equation.

source
Anime.compute_meandelays!Function
compute_meandelays!(ms::MeasurementSet, stationinfo::DataFrame, obsconfig::Dict, atmdf::DataFrame, elevationmatrix::Array{Float64, 2}, g::HDF5.Group)

Compute delays due to (mean) troposphere.

source
Anime.compute_turbulence!Function
compute_turbulence!(ms::MeasurementSet, stationinfo::DataFrame, obsconfig::Dict, elevationmatrix::Array{Float64, 2}, g::HDF5.Group)

Compute phase delays due to tropospheric turbulence. The time series of phase errors for station p is given by

\[\{{\delta \phi_p(t, \nu)}\} = \frac{1}{\sqrt{\sin({\theta_{\mathrm{el}}(t)})}} \{\delta\phi^{\prime}_p(t)\} \big(\frac{\nu}{\nu_0}\big)\]

where ν is the list of channel frequencies and ν_0 is the reference frequency (lowest in the band).

source
Anime.squaredexponentialkernelFunction
squaredexponentialkernel(x1, x2; σ=1.0, ρ=1.0)

Generate squared exponential kernel function of the form

\[k_{SE}(x-x') = \sigma^2 e^{-\frac{(x-x')^2}{2\rho^2}}\]

where σ^2 is the variance and ρ is the characteristic length.

source
Anime.rationalquadratickernelFunction
rationalquadratickernel(x1, x2; σ=1.0, α=1.0, ρ=2.0)

Generate rational quadratic kernel of the form

\[k_{RQ}(x-x') = \sigma^2 \big(1+\frac{(x-x')^2}{2\alpha\rho^2}\big)^{-\alpha}\]

source