Anime API
Contents
Index
Anime.MeasurementSet
Anime.addweightcolumns
Anime.attenuate!
Anime.bandpass!
Anime.compute_meandelays!
Anime.compute_skynoise!
Anime.compute_transmission!
Anime.compute_turbulence!
Anime.computeweights!
Anime.copymodeltodata
Anime.createmsfromconfig
Anime.createmsfromuvfits
Anime.createuvfitsfromms
Anime.elevationangle
Anime.genseries1d!
Anime.instrumentalpolarization!
Anime.makecasaantennatable
Anime.parallacticangle
Anime.plotbandpass
Anime.plotdterms
Anime.plotelevationangle
Anime.plotmeandelays
Anime.plotparallacticangle
Anime.plotpointingerrors
Anime.plotstationgains
Anime.plottransmission
Anime.plotuvcov
Anime.plotvis
Anime.pointing!
Anime.rationalquadratickernel
Anime.readalistv5
Anime.readalistv6
Anime.readbandpassinfo
Anime.readms
Anime.readobsconfig
Anime.readstationinfo
Anime.run_aatm
Anime.run_wsclean
Anime.squaredexponentialkernel
Anime.stationgains!
Anime.thermalnoise!
Anime.troposphere!
Anime.writems
Public API
I/O
Anime.readms
— Functionreadms(msname::String)
Read data from Measurement Set.
Anime.readobsconfig
— Functionreadobsconfig(configfile::String)
Read YAML configuration file with observation settings in a Dictionary.
Anime.readstationinfo
— Functionreadstationinfo(stationfile::String; delim::String=",", ignorerepeated::Bool=false)
Read additional information about stations from CSV file into a DataFrame.
Anime.readbandpassinfo
— Functionreadbandpassinfo(bandpassfile::String; delim::String=",", ignorerepeated::Bool=false)
Read bandpass information from CSV file into a DataFrame.
Anime.readalistv5
— Functionreadalistv5(alistfile::String)
Read in an alist v5 file generated by HOPS.
Anime.readalistv6
— Functionreadalistv6(alistfile::String)
Read in an alist v6 file generated by HOPS.
Anime.createmsfromconfig
— Functioncreatemsfromconfig(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
.
createmsfromconfig(obsconfig, stationinfo::DataFrame)
Alias with fewer arguments for use in pipelines.
Anime.createmsfromuvfits
— Functioncreatemsfromuvfits(uvfits::String, msname::String, mscreationmode::String; overwrite::Bool=true)
Convert uvfits
to MS. Requires casatools
and casatasks
.
Anime.createuvfitsfromms
— Functioncreateuvfitsfromms(msname::String, uvfits::String, datacolumn::String; field::String="", spw::String="", antenna::String="",
timerange::String="", overwrite::Bool=true)
Create UVFITS from existing MS. Requires casatasks
.
Anime.writems
— Functionwritems(ms::MeasurementSet; h5file::String="")
Write data to MS. Optionally, add weight columns from h5file
.
Anime.makecasaantennatable
— Functionmakecasaantennatable(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
.
Anime.addweightcolumns
— Functionaddweightcolumns(msname::String, mode::String, sigmaspec::Bool, weightspec::Bool)
Add WEIGHT_SPECTRUM
and SIGMA_SPECTRUM
columns to existing MS. Requires casatools
to be installed.
Coherency
Anime.run_wsclean
— Functionrun_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.
Instrument Models
Anime.troposphere!
— Functiontroposphere!(obs::CjlObservation, h5file::String; absorptionfile::String="", dispersivefile::String="", elevfile::String="")
Umbrella function for tropospheric model. Model is written to HDF5 file.
Anime.instrumentalpolarization!
— Functioninstrumentalpolarization!(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.
instrumentalpolarization!(ms::MeasurementSet, stationinfo::DataFrame, obsconfig::Dict; h5file::String="",
elevfile::String="", parangfile::String="")
Alias for use in pipelines.
Anime.pointing!
— Functionpointing!(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.
pointing!(ms::MeasurementSet, stationinfo::DataFrame, obsconfig::Dict; h5file::String="")
Alias with fewer arguments for use in pipelines.
Anime.stationgains!
— Functionstationgains!(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.
stationgains!(ms::MeasurementSet, stationinfo::DataFrame, obsconfig::Dict; h5file::String="")
Alias for use in pipelines.
Anime.bandpass!
— Functionbandpass!(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.
bandpass!(ms::MeasurementSet, stationinfo::DataFrame, bandpassinfo::DataFrame, obsconfig::Dict; h5file::String="")
Alias for use in pipelines.
Anime.thermalnoise!
— Functionthermalnoise!(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.
thermalnoise!(ms::MeasurementSet, stationinfo::DataFrame, obsconfig::Dict; h5file::String="", noisefile::String="")
Alias for use in pipelines.
Statistics
Anime.genseries1d!
— Functiongenseries1d!(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.
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.
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.
genseries1d!(series, times, rng::AbstractRNG; μ=0.0, σ=1.0, ρ=1.0)
Generate 1-D series using SE kernel.
genseries1d!(series, times, rng::AbstractRNG, α::Float64; μ=0.0, σ=1.0, ρ=2.0)
Generate 1-D series using RQ kernel.
Utils
Anime.elevationangle
— Functionelevationangle(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.
Anime.parallacticangle
— Functionparallacticangle(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.
Plotting
Anime.plotvis
— Functionplotvis(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.
Anime.plotuvcov
— Functionplotuvcov(uvw::Matrix{Float64}, flagrow::Vector{Bool}, chanfreqvec::Vector{Float64}; saveprefix="test_")
Plot uv-coverage of observation.
Anime.plotstationgains
— Functionplotstationgains(h5file::String, scanno::Vector{Int32}, times::Vector{Float64}, exposure::Float64, stationnames::Vector{String3})
Plot complex station gains against time.
Anime.plotbandpass
— Functionplotbandpass(h5file::String, stationnames::Vector{String3}, chanfreqvec::Vector{Float64})
Plot bandpass gains against time.
Anime.plotpointingerrors
— Functionplotpointingerrors(h5file::String, scanno::Vector{Int32}, stationnames::Vector{String3})
Plot pointing errors.
Anime.plotelevationangle
— Functionplotelevationangle(h5file::String, scanno::Vector{Int32}, times::Vector{Float64}, stationnames::Vector{String3})
Plot evolution of station elevation angles during the course of the observation.
Anime.plotparallacticangle
— Functionplotparallacticangle(h5file::String, scanno::Vector{Int32}, times::Vector{Float64}, stationnames::Vector{String3})
Plot evolution of station parallactic angles during the course of the observation.
Anime.plotdterms
— Functionplotdterms(h5file::String, stationnames::Vector{String3}, chanfreqvec::Vector{Float64})
Plot frequency-dependent complex instrumental polarization.
Anime.plottransmission
— Functionplottransmission(h5file::String, stationnames::Vector{String3}, times::Vector{Float64}, chanfreqvec::Vector{Float64})
Plot tropospheric transmission variation with frequency.
Anime.plotmeandelays
— Functionplotmeandelays(h5file::String, stationnames::Vector{String3}, times::Vector{Float64}, chanfreqvec::Vector{Float64})
Plot mean delays against time.
Data Types
Anime.MeasurementSet
— Typestruct MeasurementSet{T}
Composite type for storing data from a Measurement Set.
Fields
name
: Measurement Set name
data
: Complex visibility data (MSDATA
column)
flag
: Flag array of the same dimensions as data (MSFLAG
column)
flagrow
: Flag row Boolean vector (MSFLAG_ROW
column)
antenna1
: Antenna 1 in a baseline pair (MSANTENNA1
column)
antenna2
: Antenna 2 in a baseline pair (MSANTENNA2
column)
uvw
: uvw coordinates (MSUVW
column)
times
: Timestamps (MSTIME
column)
exposure
: Integration time
scanno
: Scan numbers (SCAN_NUMBER
column in MS)
numchan
: Number of frequency channels
chanfreqvec
: Channel frequencies (Hz) (MSCHAN_FREQ
column)
chanwidth
: Width of frequency channel
phasedir
: Direction of phase centre
pos
: Antenna positions in x,y,z coordinates
Internal
Anime.copymodeltodata
— Functioncopymodeltodata(msname::String)
Copy MODEL_DATA
to DATA
in MS msname
.
Anime.computeweights!
— Functioncomputeweights!(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
.
Anime.run_aatm
— Functionrun_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
.
Anime.compute_transmission!
— Functioncompute_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}}}\]
Anime.attenuate!
— Functionattenuate!(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}}}\]
Anime.compute_skynoise!
— Functioncompute_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.
Anime.compute_meandelays!
— Functioncompute_meandelays!(ms::MeasurementSet, stationinfo::DataFrame, obsconfig::Dict, atmdf::DataFrame, elevationmatrix::Array{Float64, 2}, g::HDF5.Group)
Compute delays due to (mean) troposphere.
Anime.compute_turbulence!
— Functioncompute_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).
Anime.squaredexponentialkernel
— Functionsquaredexponentialkernel(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.
Anime.rationalquadratickernel
— Functionrationalquadratickernel(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}\]