ComradeBase

Documentation for ComradeBase.

ComradeBase.AbstractModelType
AbstractModel

The Comrade abstract model type. To instantiate your own model type you should subtybe from this model. Additionally you need to implement the following methods to satify the interface:

Mandatory Methods

  • visanalytic: defines whether the model visibilities can be computed analytically. If yes then this should return IsAnalytic() and the user must to define visibility_point. If not analytic then visanalytic should return NotAnalytic().
  • imanalytic: defines whether the model intensities can be computed pointwise. If yes then this should return IsAnalytic() and the user must to define intensity_point. If not analytic then imanalytic should return NotAnalytic().
  • radialextent: Provides a estimate of the radial extent of the model in the image domain. This is used for estimating the size of the image, and for plotting.
  • flux: Returns the total flux of the model.
  • intensity_point: Defines how to compute model intensities pointwise. Note this is must be defined if imanalytic(::Type{YourModel})==IsAnalytic().
  • visibility_point: Defines how to compute model visibilties pointwise. Note this is must be defined if visanalytic(::Type{YourModel})==IsAnalytic().

Optional Methods:

  • ispolarized: Specified whether a model is intrinsically polarized (returns IsPolarized()) or is not (returns NotPolarized()), by default a model is NotPolarized()
  • visibilitymap_analytic: Vectorized version of visibility_point for models where visanalytic returns IsAnalytic()
  • visibilitymap_numeric: Vectorized version of visibility_point for models where visanalytic returns NotAnalytic() typically these are numerical FT's
  • intensitymap_analytic: Computes the entire image for models where imanalytic returns IsAnalytic()
  • intensitymap_numeric: Computes the entire image for models where imanalytic returns NotAnalytic()
  • intensitymap_analytic!: Inplace version of intensitymap
  • intensitymap_numeric!: Inplace version of intensitymap
ComradeBase.AbstractPolarizedModelType
abstract type AbstractPolarizedModel <: ComradeBase.AbstractModel

A generic polarized model. To implement the use needs to follow the AbstractModel implementation instructions. In addtion there is an optional method stokes(model, p::Symbol) which extracts the specific stokes parameter of the model. The default that the different stokes parameters are stored as fields of the model. To overwrite this behavior overload the function.

ComradeBase.DensityAnalyticType
DensityAnalytic

Internal type for specifying the nature of the model functions. Whether they can be easily evaluated pointwise analytic. This is an internal type that may change.

ComradeBase.IntensityMapType
struct IntensityMap{T, N, D, G<:(ComradeBase.AbstractRectiGrid{D}), A<:AbstractArray{T, N}, R<:Tuple, Na} <: DimensionalData.AbstractDimArray{T, N, D, A<:AbstractArray{T, N}}

This type is the basic array type for all images and models that obey the ComradeBase interface. The type is a subtype of DimensionalData.AbstractDimArray however, we make a few changes to support the Comrade API.

  1. The dimensions should be specified by an AbstractRectiGrid interface. Usually users just need the RectiGrid grid, for rectilinear grids.
  2. There are two ways to access the dimensions of the array. dims(img) will return the usual DimArray dimensions, i.e. a Tuple{DimensionalData.Dim, ...}. The other way to access the array dimensions is using the getproperty, e.g., img.X will return the RA/X grid locations but stripped of the usual DimensionalData.Dimension material. This getproperty behavior is *NOT CONSIDERED** part of the stable API and may be changed in the future.
  3. Metadata is stored in the AbstractRectiGrid type through the header property and can be accessed through metadata or header

The most common way to create a IntensityMap is to use the function definitions

julia> g = imagepixels(10.0, 10.0, 128, 128; header=NoHeader())
julia> X = g.X; Y = g.Y
julia> data = rand(128, 128)
julia> img1 = IntensityMap(data, g)
julia> img2 = IntensityMap(data, (;X, Y); header=header(g))
julia> img1 == img2
true
julia> img3 = IntensityMap(data, 10.0, 10.0; header=NoHeader())

Broadcasting, map, and reductions should all just obey the DimensionalData interface.

ComradeBase.IntensityMapMethod
IntensityMap(data::AbstractArray, g::AbstractRectiGrid; refdims=(), name=Symbol(""))

Creates a IntensityMap with the pixel fluxes data on the grid g. Optionally, you can specify a set of reference dimensions refdims as a tuple and a name for array name.

ComradeBase.IntensityMapMethod
IntensityMap(data::AbstractArray, fovx::Real, fovy::Real, x0::Real=0, y0::Real=0; header=NoHeader())

Creates a IntensityMap with the pixel fluxes data and a spatial grid with field of view (fovx, fovy) and center pixel offset (x0, y0) and header header.

ComradeBase.IsAnalyticType
struct IsAnalytic <: ComradeBase.DensityAnalytic

Defines a trait that a states that a model is analytic. This is usually used with an abstract model where we use it to specify whether a model has a analytic fourier transform and/or image.

ComradeBase.MinimalHeaderType
MinimalHeader{T}

A minimal header type for ancillary image information.

Fields

  • source: Common source name
  • ra: Right ascension of the image in degrees (J2000)
  • dec: Declination of the image in degrees (J2000)
  • mjd: Modified Julian Date in days
  • frequency: Frequency of the image in Hz
ComradeBase.NotAnalyticType
struct NotAnalytic <: ComradeBase.DensityAnalytic

Defines a trait that a states that a model is analytic. This is usually used with an abstract model where we use it to specify whether a model has does not have a easy analytic fourier transform and/or intensity function.

ComradeBase.RectiGridMethod
RectiGrid(dims::NamedTuple{Na}; executor=Serial(), header=ComradeBase.NoHeader())
RectiGrid(dims::NTuple{N, <:DimensionalData.Dimension}; executor=Serial(), header=ComradeBase.NoHeader())

Creates a rectilinear grid of pixels with the dimensions dims. The dims can either be a named tuple of dimensions or a tuple of dimensions. The dimensions can be in any order however the standard orders are:

  • (:X, :Y, :Ti, :Fr)
  • (:X, :Y, :Fr, :Ti)
  • (:X, :Y) # spatial only

where X,Y are the RA and DEC spatial dimensions respectively, Ti is the time dimension and Fr is the frequency dimension.

Note that the majority of the time users should just call imagepixels to create a spatial grid.

Optional Arguments

  • executor: specifies how different models are executed. The default is Serial which mean serial CPU computations. For threaded computations use ThreadsEx() or load OhMyThreads.jl to uses their schedulers.
  • header: specified underlying header information for the grid. This is used to store information about the image such as the source, RA and DEC, MJD.

Examples

dims = RectiGrid((X(-5.0:0.1:5.0), Y(-4.0:0.1:4.0), Ti([1.0, 1.5, 1.75]), Fr([230, 345])); executor=ThreadsEx())
dims = RectiGrid((X = -5.0:0.1:5.0, Y = -4.0:0.1:4.0, Ti = [1.0, 1.5, 1.75], Fr = [230, 345]); executor=ThreadsEx()))
ComradeBase.SerialType
Serial()

Uses serial execution when computing the intensitymap or visibilitymap

ComradeBase.ThreadsExType
ThreadsEx(;scheduler::Symbol = :dynamic)

Uses Julia's Threads @threads macro when computing the intensitymap or visibilitymap. You can choose from Julia's various schedulers by passing the scheduler as a parameter. The default is :dynamic, but it isn't considered part of the stable API and may change at any moment.

ComradeBase.UnstructuredDomainMethod
UnstructuredDomain(dims::AbstractArray; executor=Serial(), header=ComradeBase.NoHeader)

Builds an unstructured grid (really a vector of points) from the dimensions dims. The executor is used controls how the grid is computed when calling visibilitymap or intensitymap. The default is Serial which mean regular CPU computations. For threaded execution use ThreadsEx() or load OhMyThreads.jl to uses their schedulers.

Note that unlike RectiGrid which assigns dimensions to the grid points, UnstructuredDomain does not. This is becuase the grid is unstructured the points are a cloud in a space

ComradeBase.UnstructuredMapType
UnstructuredMap(data::AbstractVector, dims::UnstructuredDomain)

A map that is defined on an unstructured domain. This is typically just a vector of values. The vector of locations of the visibilities are stored in dims. Otherwise this behaves very similarly to IntensityMap, except that is isn't a grid.

For instance the locations of the visibilities can be accessed with axisdims, as well as the usual getproperty and propertynames functions. Like with IntensityMap during execution the executor is used to determine the execution context.

ComradeBase._visibilitymapFunction
_visibilitymap(model::AbstractModel, p)

Internal method used for trait dispatch and unpacking of args arguments in visibilities

Warn

Not part of the public API so it may change at any moment.

ComradeBase._visibilitymap!Function
_visibilitymap!(model::AbstractModel, p)

Internal method used for trait dispatch and unpacking of args arguments in visibilities!

Warn

Not part of the public API so it may change at any moment.

ComradeBase.allocate_imgmapFunction
allocate_imgmap(m::AbstractModel, g::AbstractSingleDomain)

Allocate the default map specialized by the grid g

ComradeBase.allocate_vismapFunction
allocate_vismap(m::AbstractModel, g::AbstractSingleDomain)

Allocate the default map specialized by the grid g

ComradeBase.amplitudeMethod
amplitude(model, p)

Computes the visibility amplitude of model m at the coordinate p. The coordinate p is expected to have the properties U, V, and sometimes Ti and Fr.

If you want to compute the amplitudemap at a large number of positions consider using the amplitudemap function.

ComradeBase.amplitudemapMethod
amplitudemap(m::AbstractModel, p)

Computes the visibility amplitudemap of the model m at the coordinates p. The coordinates p are expected to have the properties U, V, and sometimes Ti and Fr.

ComradeBase.axisdimsMethod
axisdims(img::IntensityMap)
axisdims(img::IntensityMap, p::Symbol)

Returns the keys of the IntensityMap as the actual internal AbstractRectiGrid object. Optionall the user can ask for a specific dimension with p

ComradeBase.bispectrumMethod
bispectrum(model, p1, p2, p3)

Computes the complex bispectrum of model m at the uv-triangle p1 -> p2 -> p3

If you want to compute the bispectrum over a number of triangles consider using the bispectrummap function.

ComradeBase.bispectrummapMethod
bispectrummap(m, p1, p2, p3)

Computes the closure phases of the model m at the triangles p1, p2, p3, where pi are coordinates.

ComradeBase.centroidMethod
centroid(im::AbstractIntensityMap)

Computes the image centroid aka the center of light of the image.

For polarized maps we return the centroid for Stokes I only.

ComradeBase.closure_phaseMethod
closure_phase(model, p1, p2, p3, p4)

Computes the closure phase of model m at the uv-triangle u1,v1 -> u2,v2 -> u3,v3

If you want to compute closure phases over a number of triangles consider using the closure_phasemap function.

ComradeBase.closure_phasemapMethod
closure_phasemap(m,
               p1::AbstractArray
               p2::AbstractArray
               p3::AbstractArray
               )

Computes the closure phases of the model m at the triangles p1, p2, p3, where pi are coordinates.

ComradeBase.create_imgmapMethod
create_imgmap(array, g::AbstractSingleDomain)

Create a map of values specialized by the grid g in the image domain. The default is to call create_map with the same arguments.

ComradeBase.create_mapFunction
create_map(array, g::AbstractSingleDomain)

Create a map of values specialized by the grid g.

ComradeBase.create_vismapMethod
create_vismap(array, g::AbstractSingleDomain)

Create a map of values specialized by the grid g in the visibility domain. The default is to call create_map with the same arguments.

ComradeBase.domainpointsFunction
domainpoints(g::AbstractSingleDomain)

Create a grid iterator that can be used to iterate through different points. All grid methods must implement this method.

ComradeBase.domainpointsMethod
domainpoints(k::IntensityMap)

Returns the grid the IntensityMap is defined as. Note that this is unallocating since it lazily computes the grid. The grid is an example of a DimArray and works similarly. This is useful for broadcasting a model across an abritrary grid.

ComradeBase.executorMethod
executor(g::AbstractSingleDomain)

Returns the executor used to compute the intensitymap or visibilitymap

ComradeBase.fieldofviewMethod
fieldofview(img::IntensityMap)
fieldofview(img::IntensityMap)

Returns a named tuple with the field of view of the image.

ComradeBase.fluxMethod
flux(im::IntensityMap)

Computes the flux of a intensity map

ComradeBase.headerMethod
header(g::AbstractSingleDomain)

Returns the headerinformation of the dimensions g

ComradeBase.headerMethod
header(img::IntensityMap)

Retrieves the header of an IntensityMap

ComradeBase.imagepixelsFunction
imagepixels(fovx, fovy, nx, ny; x0=0, y0=0, executor=Serial(), header=NoHeader())

Construct a grid of pixels with a field of view fovx and fovy and nx and ny pixels. This points are the pixel centers and the field of view goes from the edge of the first pixel to the edge of the last pixel. The x0, y0 offsets shift the image origin over by (x0, y0) in the image plane.

Arguments:

  • fovx::Real: The field of view in the x-direction
  • fovy::Real: The field of view in the y-direction
  • nx::Integer: The number of pixels in the x-direction
  • ny::Integer: The number of pixels in the y-direction

Keyword Arguments:

  • x0::Real=0: The x-offset of the image
  • y0::Real=0: The y-offset of the image
  • executor=Serial(): The executor to use for the grid, default is serial execution
  • header=NoHeader(): The header to use for the grid
ComradeBase.imanalyticMethod
imanalytic(::Type{<:AbstractModel})

Determines whether the model is pointwise analytic in the image domain, i.e. we can evaluate its intensity at an arbritrary point.

If IsAnalytic() then it will try to call intensity_point to calculate the intensity.

ComradeBase.intensity_pointFunction
intensity_point(model::AbstractModel, p)

Function that computes the pointwise intensity if the model has the trait in the image domain IsAnalytic(). Otherwise it will use construct the image in visibility space and invert it.

ComradeBase.intensitymap!Function
intensitymap!(buffer::AbstractDimArray, model::AbstractModel)

Computes the intensity map of model by modifying the buffer

ComradeBase.intensitymap!Method
intensitymap!(img::AbstractIntensityMap, mode;, executor = SequentialEx())

Computes the intensity map or image of the model. This updates the IntensityMap object img.

Optionally the user can specify the executor that uses FLoops.jl to specify how the loop is done. By default we use the SequentialEx which uses a single-core to construct the image.

ComradeBase.intensitymapMethod
intensitymap(model::AbstractModel, dims::AbstractSingleDomain)

Computes the intensity map or image of the model. This returns an IntensityMap which is a IntensityMap with dims an AbstractSingleDomain as dimensions.

ComradeBase.intensitymap_numericFunction
intensitymap_numeric(m::AbstractModel, p::AbstractSingleDomain)

Computes the IntensityMap of a model m at the image positions p using a numerical method. This has to be specified uniquely for every model m if imanalytic(typeof(m)) === NotAnalytic(). See Comrade.jl for example implementations.

ComradeBase.intensitymap_numeric!Function
intensitymap_numeric!(img::IntensityMap, m::AbstractModel)

Updates the img using the model m using a numerical method. This has to be specified uniquely for every model m if imanalytic(typeof(m)) === NotAnalytic(). See Comrade.jl for example implementations.

ComradeBase.ispolarizedMethod
ispolarized(::Type)

Trait function that defines whether a model is polarized or not.

ComradeBase.loadMethod
ComradeBase.load(fitsfile::String, IntensityMap)

This loads in a fits file that is more robust to the various imaging algorithms in the EHT, i.e. is works with clean, smili, eht-imaging. The function returns an tuple with an intensitymap and a second named tuple with ancillary information about the image, like the source name, location, mjd, and radio frequency.

ComradeBase.logclosure_amplitudeMethod
logclosure_amplitude(model, p1, p2, p3, p4)

Computes the log-closure amplitude of model m at the uv-quadrangle u1,v1 -> u2,v2 -> u3,v3 -> u4,v4 using the formula

\[C = \log\left|\frac{V(u1,v1)V(u2,v2)}{V(u3,v3)V(u4,v4)}\right|\]

If you want to compute log closure amplitudemap over a number of triangles consider using the logclosure_amplitudemap function.

ComradeBase.logclosure_amplitudemapMethod
logclosure_amplitudemap(m::AbstractModel,
                      p1,
                      p2,
                      p3,
                      p4
                     )

Computes the log closure amplitudemap of the model m at the quadrangles p1, p2, p3, p4.

ComradeBase.named_dimsMethod
named_dims(g::AbstractSingleDomain)

Returns a named tuple containing the dimensions of g. For a unnamed version see dims

ComradeBase.phasecenterMethod
phasecenter(img::IntensityMap)

Computes the phase center of an intensity map. Note this is the pixels that is in the middle of the image.

ComradeBase.pixelsizesMethod
pixelsizes(img::IntensityMap)
pixelsizes(img::AbstractRectiGrid)

Returns a named tuple with the spatial pixel sizes of the image.

ComradeBase.radialextentFunction
radialextent(model::AbstractModel)

Provides an estimate of the radial size/extent of the model. This is used internally to estimate image size when plotting and using modelimage

ComradeBase.saveMethod
ComradeBase.save(file::String, img::IntensityMap, obs)

Saves an image to a fits file. You can optionally pass an EHTObservation so that ancillary information will be added.

ComradeBase.second_momentMethod
second_moment(im::AbstractIntensityMap; center=true)

Computes the image second moment tensor of the image. By default we really return the second cumulant or centered second moment, which is specified by the center argument.

ComradeBase.second_momentMethod
second_moment(im::AbstractIntensityMap; center=true)

Computes the image second moment tensor of the image. By default we really return the second cumulant or centered second moment, which is specified by the center argument.

For polarized maps we return the second moment for Stokes I only.

ComradeBase.stokesMethod
stokes(m::AbstractArray{<:StokesParams}, p::Symbol)

Extract the specific stokes component p from the polarized image m.

ComradeBase.stokesMethod
stokes(m::AbstractPolarizedModel, p::Symbol)

Extract the specific stokes component p from the polarized model m

ComradeBase.visanalyticMethod
visanalytic(::Type{<:AbstractModel})

Determines whether the model is pointwise analytic in Fourier domain, i.e. we can evaluate its fourier transform at an arbritrary point.

If IsAnalytic() then it will try to call visibility_point to calculate the complex visibilities. Otherwise it fallback to using the FFT that works for all models that can compute an image.

ComradeBase.visibilityMethod
visibility(mimg, p)

Computes the complex visibility of model m at coordinates p. p corresponds to the coordinates of the model. These need to have the properties U, V and sometimes Ti for time and Fr for frequency.

Notes

If you want to compute the visibilities at a large number of positions consider using the visibilitymap.

Warn

This is only defined for analytic models. If you want to compute the visibility for a single point for a non-analytic model, please use the visibilitymap function and create an UnstructuredDomain with a single point.

ComradeBase.visibility_pointFunction
visibility_point(model::AbstractModel, p)

Function that computes the pointwise visibility. This must be implemented in the model interface if visanalytic(::Type{MyModel}) == IsAnalytic()

ComradeBase.visibilitymapFunction
visibilitymap(model::AbstractModel, p)

Computes the complex visibilities at the locations p.

ComradeBase.visibilitymap!Function
visibilitymap!(vis::AbstractArray, model::AbstractModel, p)

Computes the complex visibilities vis in place at the locations p

ComradeBase.visibilitymap!Method
visibilitymap!(vis, m, p)

Computes the visibilities vis in place of the model m using the coordinates p. The coordinates p are expected to have the properties U, V, and sometimes T and F.

ComradeBase.visibilitymapMethod
visibilitymap(m, p)

Computes the visibilities of the model m using the coordinates p. The coordinates p are expected to have the properties U, V, and sometimes T and F.

ComradeBase.visibilitymap_analyticFunction
visibilties_analytic(model, p)

Computes the visibilties of a model using using the analytic visibility expression given by visibility_point.

ComradeBase.visibilitymap_analytic!Function
visibilties_analytic!(vis, model)

Computes the visibilties of a model in-place, using using the analytic visibility expression given by visibility_point.

ComradeBase.visibilitymap_numericFunction
visibilties_numeric(model, p)

Computes the visibilties of a model using a numerical fourier transform. Note that none of these are implemented in ComradeBase. For implementations please see Comrade.

ComradeBase.visibilitymap_numeric!Function
visibilties_numeric!(vis, model)

Computes the visibilties of a model in-place using a numerical fourier transform. Note that none of these are implemented in ComradeBase. For implementations please see Comrade.