Sky Model¶
-
class
oskar.
Sky
(precision=None, settings=None)¶ This class provides a Python interface to an OSKAR sky model.
The
oskar.Sky
class can be used to create an OSKAR sky model from Python. The sky model can be thought of as simply a table of source data, where each row of the table contains parameters for a single source. The sky model table format is described in the sky model documentation and a summary of the columns is also included below for convenience. Only the first three are manadatory (others can be omitted, and will default to zero).Column Parameter Unit 1 Right Ascension deg 2 Declination deg 3 Stokes I flux Jy 4 Stokes Q flux Jy 5 Stokes U flux Jy 6 Stokes V flux Jy 7 Reference frequency Hz 8 Spectral index N/A 9 Rotation measure rad / m^2 10 Major axis FWHM arcsec 11 Minor axis FWHM arcsec 12 Position angle deg A sky model can be loaded from a text file using the
load()
method, or created directly from a numpy array using thefrom_array()
method. Sky models can be concatenated usingappend()
, which is useful if building up a sky model in multiple stages.Sky models can be filtered to exclude sources outside specified radii from a given point using the method
filter_by_radius()
, and filtered by Stokes-I flux using the methodfilter_by_flux()
.The number of sources in a sky model is available from the
num_sources
property. A copy of the sky model data can be returned as a numpy array using theto_array()
method, and written as a text file using thesave()
method.Note that the sky model used by OSKAR exists independently from the simulated observation parameters: make sure the phase centre is set appropriately when running simulations, otherwise you may not see what you were expecting!
Some examples of setting up a sky model from Python are shown below.
Examples
To load a sky model text file
my_sky.txt
:>>> sky = oskar.Sky.load('my_sky.txt')
To create a sky model from a numpy array:
>>> # Specifying only RA, Dec and Stokes I (other columns default to 0). >>> data = numpy.array([[20.0, -30.0, 1], [20.0, -30.5, 3], [20.5, -30.5, 3]]) >>> sky = oskar.Sky.from_array(data) >>> print(sky.num_sources) 3
To create a sky model from columns in a FITS binary table
GLEAM_EGC.fits
:>>> from astropy.io import fits >>> hdulist = fits.open('GLEAM_EGC.fits') >>> cols = hdulist[1].data[0].array >>> data = numpy.column_stack( (cols['RAJ2000'], cols['DEJ2000'], cols['peak_flux_wide'])) >>> sky = oskar.Sky.from_array(data) >>> print(sky.num_sources) 307455
To create an all-sky model containing 100000 sources with fluxes between 1 mJy and 100 mJy, and a power law index for the luminosity function of -2:
>>> sky = oskar.Sky.generate_random_power_law(100000, 1e-3, 100e-3, -2) >>> print(sky.num_sources) 100000
To filter the sky model to contain sources only between 5 degrees and 15 degrees from the point (RA, Dec) = (0, 80) degrees:
>>> # (continued from previous section) >>> ra0 = 0 >>> dec0 = 80 >>> sky.filter_by_radius(5, 15, ra0, dec0) >>> print(sky.num_sources) 1521
To plot the sky model using matplotlib:
>>> # (continued from previous section) >>> import matplotlib.pyplot as plt >>> data = sky.to_array() # Get sky model data as numpy array. >>> ra = numpy.radians(data[:, 0] - ra0) >>> dec = numpy.radians(data[:, 1]) >>> log_flux = numpy.log10(data[:, 2]) >>> x = numpy.cos(dec) * numpy.sin(ra) >>> y = numpy.cos(numpy.radians(dec0)) * numpy.sin(dec) - \ numpy.sin(numpy.radians(dec0)) * numpy.cos(dec) * numpy.cos(ra) >>> sc = plt.scatter(x, y, s=5, c=log_flux, cmap='plasma', vmin=numpy.min(log_flux), vmax=numpy.max(log_flux)) >>> plt.axis('equal') >>> plt.xlabel('x direction cosine') >>> plt.ylabel('y direction cosine') >>> plt.colorbar(sc, label='Log10(Stokes I flux [Jy])') >>> plt.show()
-
__init__
(precision=None, settings=None)¶ Creates an OSKAR sky model.
Parameters: - precision (Optional[str]) – Either ‘double’ or ‘single’ to specify the numerical precision of the sky model. Default ‘double’.
- settings (Optional[oskar.SettingsTree]) – Optional settings to use to set up the sky model.
-
append
(other)¶ Appends data from another sky model.
Parameters: other (oskar.Sky) – Another sky model.
-
append_sources
(ra_deg, dec_deg, I, Q=None, U=None, V=None, ref_freq_hz=None, spectral_index=None, rotation_measure=None, major_axis_arcsec=None, minor_axis_arcsec=None, position_angle_deg=None)¶ Appends source data to an OSKAR sky model from arrays in memory.
Parameters: - ra_deg (float, array-like) – Source Right Ascension values, in degrees.
- dec_deg (float, array-like) – Source Declination values, in degrees.
- I (float, array-like) – Source Stokes I fluxes, in Jy.
- Q (Optional[float, array-like]) – Source Stokes Q fluxes, in Jy.
- U (Optional[float, array-like]) – Source Stokes U fluxes, in Jy.
- V (Optional[float, array-like]) – Source Stokes V fluxes, in Jy.
- ref_freq_hz (Optional[float, array-like]) – Source reference frequency values, in Hz.
- spectral_index (Optional[float, array-like]) – Source spectral index values.
- rotation_measure (Optional[float, array-like]) – Source rotation measure values, in rad/m^2.
- major_axis_arcsec (Optional[float, array-like]) – Source Gaussian major axis values, in arcsec.
- minor_axis_arcsec (Optional[float, array-like]) – Source Gaussian minor axis values, in arcsec.
- position_angle_deg (Optional[float, array-like]) – Source Gaussian position angle values, in degrees.
-
append_file
(filename)¶ Appends data to the sky model from a text file.
Parameters: filename (str) – Name of file to load.
-
create_copy
()¶ Creates a copy of the sky model.
Returns: A copy of the sky model. Return type: oskar.Sky
-
filter_by_flux
(min_flux_jy, max_flux_jy)¶ Filters the sky model according to Stokes-I flux.
Sources with flux values outside the range are deleted.
Parameters: - min_flux_jy (float) – Minimum allowed flux, in Jy.
- max_flux_jy (float) – Maximum allowed flux, in Jy.
-
filter_by_radius
(inner_radius_deg, outer_radius_deg, ra0_deg, dec0_deg)¶ Filters the sky model according to source radius from phase centre.
Sources outside the range are deleted.
Parameters: - inner_radius_deg (float) – Minimum allowed radius, in degrees.
- outer_radius_deg (float) – Maximum allowed radius, in degrees.
- ra0_deg (float) – Right Ascension of phase centre, in degrees.
- dec0_deg (float) – Declination of phase centre, in degrees.
-
classmethod
from_array
(array, precision='double')¶ Creates a new sky model from a 2D numpy array.
The format of the array is the same as that in sky model text files. Each column specifies a different source parameter, and each row specifies data for a different source.
The array could be created, for example, like: array = numpy.zeros((num_sources, num_parameters))
There must be at least 3 columns present (RA, Dec, Stokes I). Parameters for missing columns will take default values.
If the array is 1-dimensional, it will be treated as specifying parameters only for a single source.
Parameters: - array (float, array-like) – Input array.
- precision (Optional[str]) – Either ‘double’ or ‘single’ to specify the numerical precision of the sky model.
Returns: The new sky model.
Return type:
-
classmethod
from_fits_file
(filename, min_peak_fraction=0.0, min_abs_val=0.0, default_map_units='K', override_units=False, frequency_hz=0.0, spectral_index=-0.7, precision='double')¶ Loads data from a FITS file and returns it as a new sky model.
The file can be either a regular FITS image or a HEALPix FITS file in RING format.
Parameters: - filename (str) – Name of FITS file to load.
- min_peak_fraction (Optional[float]) – Minimum pixel value loaded, as a fraction of the image peak.
- min_abs_val (Optional[float]) – Minimum pixel value loaded.
- default_map_units (Optional[str]) – Default map units, if not found in the file. Can be ‘Jy/beam’, ‘Jy/pixel’, ‘K’ or ‘mK’.
- override_units (Optional[boolean]) – If true, override image units with the default.
- frequency_hz (Optional[float]) – Frequency of the image data in Hz, if not found in the file.
- spectral_index (Optional[float]) – Spectral index value to give to each pixel.
- precision (Optional[str]) – Either ‘double’ or ‘single’ to specify the numerical precision of the sky model.
Returns: The new sky model.
Return type:
-
classmethod
generate_grid
(ra0_deg, dec0_deg, side_length, fov_deg, mean_flux_jy=1.0, std_flux_jy=0.0, seed=1, precision='double')¶ Generates a grid of sources and returns it as a new sky model.
Parameters: - ra0_deg (float) – Right Ascension of grid centre, in degrees.
- dec0_deg (float) – Declination of grid centre, in degrees.
- side_length (int) – Side length of generated grid.
- fov_deg (float) – Grid field-of-view, in degrees.
- mean_flux_jy (float) – Mean Stokes-I source flux, in Jy.
- std_flux_jy (float) – Standard deviation Stokes-I flux, in Jy.
- seed (int) – Random generator seed.
- precision (Optional[str]) – Either ‘double’ or ‘single’ to specify the numerical precision of the sky model.
Returns: The new sky model.
Return type:
-
classmethod
generate_random_power_law
(num_sources, min_flux_jy, max_flux_jy, power_law_index, seed=1, precision='double')¶ Generates sources scattered randomly over the celestial sphere.
Parameters: - num_sources (int) – The number of sources to generate.
- min_flux_jy (float) – Minimum Stokes-I source flux, in Jy.
- max_flux_jy (float) – Maximum Stokes-I source flux, in Jy.
- power_law_index (float) – Power law index/exponent.
- seed (int) – Random generator seed.
- precision (Optional[str]) – Either ‘double’ or ‘single’ to specify the numerical precision of the sky model.
Returns: The new sky model.
Return type:
-
classmethod
load
(filename, precision='double')¶ Loads data from a text file and returns it as a new sky model.
Parameters: - filename (str) – Name of file to load.
- precision (Optional[str]) – Either ‘double’ or ‘single’ to specify the numerical precision of the sky model.
Returns: The new sky model.
Return type:
-
save
(filename)¶ Saves data to a sky model text file.
Parameters: filename (str) – Name of file to write.
-
to_array
()¶ Returns a copy of the sky model as a 2D numpy array.
Returns: A copy of the sky model. Return type: numpy.ndarray
-
to_ds9_regions
(filename, colour='green', width=1)¶ Writes a ds9 region file from coordinates in the sky model.
Parameters: - filename (str) – Name of file to write.
- colour (Optional[str]) – Name of region colour. Default ‘green’.
- width (Optional[int]) – Line width. Default 1.
-
num_sources
¶ Returns the number of sources in the sky model.
- Type
- int
-