def create_test_image( canonical=True, cellsize=None, frequency=None, channel_bandwidth=None, phasecentre=None, polarisation_frame=PolarisationFrame("stokesI")) -> Image: """Create a useful test image This is the test image M31 widely used in ALMA and other simulations. It is actually part of an Halpha region in M31. :param canonical: Make the image into a 4 dimensional image :param cellsize: :param frequency: Frequency (array) in Hz :param channel_bandwidth: Channel bandwidth (array) in Hz :param phasecentre: Phase centre of image (SkyCoord) :param polarisation_frame: Polarisation frame :return: Image """ check_data_directory() if frequency is None: frequency = [1e8] im = import_image_from_fits(rascil_path("data/models/M31.MOD")) if canonical: if polarisation_frame is None: im.polarisation_frame = PolarisationFrame("stokesI") elif isinstance(polarisation_frame, PolarisationFrame): im.polarisation_frame = polarisation_frame else: raise ValueError("polarisation_frame is not valid") im = replicate_image(im, frequency=frequency, polarisation_frame=im.polarisation_frame) if cellsize is not None: im.wcs.wcs.cdelt[0] = -180.0 * cellsize / numpy.pi im.wcs.wcs.cdelt[1] = +180.0 * cellsize / numpy.pi if frequency is not None: im.wcs.wcs.crval[3] = frequency[0] if channel_bandwidth is not None: im.wcs.wcs.cdelt[3] = channel_bandwidth[0] else: if len(frequency) > 1: im.wcs.wcs.cdelt[3] = frequency[1] - frequency[0] else: im.wcs.wcs.cdelt[3] = 0.001 * frequency[0] im.wcs.wcs.radesys = 'ICRS' im.wcs.wcs.equinox = 2000.00 if phasecentre is not None: im.wcs.wcs.crval[0] = phasecentre.ra.deg im.wcs.wcs.crval[1] = phasecentre.dec.deg # WCS is 1 relative im.wcs.wcs.crpix[0] = im.data.shape[3] // 2 + 1 im.wcs.wcs.crpix[1] = im.data.shape[2] // 2 + 1 return im
def create_configuration_from_MIDfile(antfile: str, location=None, mount: str = 'azel', rmax=None, name='') -> Configuration: """ Define configuration from a SKA MID format file :param antfile: Antenna file name :param mount: mount type: 'azel', 'xy' :param rmax: Maximum distance from array centre (m) :param name: Name of array :return: Configuration """ check_data_directory() # X Y Z Diam Station # 9.36976 35.48262 1052.99987 13.50 M001 antxyz = numpy.genfromtxt(antfile, skip_header=5, usecols=[0, 1, 2], delimiter=" ") antxyz = xyz_at_latitude(antxyz, location.lat.rad) antxyz += [ location.geocentric[0].to(u.m).value, location.geocentric[1].to(u.m).value, location.geocentric[2].to(u.m).value ] nants = antxyz.shape[0] assert antxyz.shape[ 1] == 3, "Antenna array has wrong shape %s" % antxyz.shape anames = numpy.genfromtxt(antfile, dtype='str', skip_header=5, usecols=[4], delimiter=" ") mounts = numpy.repeat(mount, nants) diameters = numpy.genfromtxt(antfile, dtype='str', skip_header=5, usecols=[3], delimiter=" ") antxyz, diameters, anames, mounts = limit_rmax(antxyz, diameters, anames, mounts, rmax) fc = Configuration(location=location, names=anames, mount=mounts, xyz=antxyz, diameter=diameters, name=name) return fc
def create_configuration_from_SKAfile(antfile: str, mount: str = 'azel', names: str = "%d", rmax=None, name='', location=None) -> Configuration: """ Define configuration from a SKA format file :param antfile: Antenna file name :param location: Earthlocation of array :param mount: mount type: 'azel', 'xy', 'equatorial' :param names: Antenna names e.g. "VLA%d" :param rmax: Maximum distance from array centre (m) :param name: Name of array :return: Configuration """ check_data_directory() antdiamlonglat = numpy.genfromtxt(antfile, usecols=[0, 1, 2], delimiter="\t") assert antdiamlonglat.shape[1] == 3, ("Antenna array has wrong shape %s" % antdiamlonglat.shape) antxyz = numpy.zeros([antdiamlonglat.shape[0] - 1, 3]) diameters = numpy.zeros([antdiamlonglat.shape[0] - 1]) for ant in range(antdiamlonglat.shape[0] - 1): loc = EarthLocation(lon=antdiamlonglat[ant, 1], lat=antdiamlonglat[ant, 2], height=0.0).geocentric antxyz[ant] = [ loc[0].to(u.m).value, loc[1].to(u.m).value, loc[2].to(u.m).value ] diameters[ant] = antdiamlonglat[ant, 0] nants = antxyz.shape[0] anames = [names % ant for ant in range(nants)] mounts = numpy.repeat(mount, nants) antxyz, diameters, anames, mounts = limit_rmax(antxyz, diameters, anames, mounts, rmax) fc = Configuration(location=location, names=anames, mount=mounts, xyz=antxyz, diameter=diameters, name=name) return fc
def create_configuration_from_file(antfile: str, location: EarthLocation = None, mount: str = 'azel', names: str = "%d", diameter=35.0, rmax=None, name='') -> Configuration: """ Define configuration from a text file :param antfile: Antenna file name :param location: Earthlocation of array :param mount: mount type: 'azel', 'xy', 'equatorial' :param names: Antenna names e.g. "VLA%d" :param diameter: Effective diameter of station or antenna :param rmax: Maximum distance from array centre (m) :param name: Name of array :return: Configuration """ check_data_directory() antxyz = numpy.genfromtxt(antfile, delimiter=",") assert antxyz.shape[1] == 3, ("Antenna array has wrong shape %s" % antxyz.shape) latitude = location.geodetic[1].to(u.rad).value antxyz = xyz_at_latitude(antxyz, latitude) antxyz += [ location.geocentric[0].to(u.m).value, location.geocentric[1].to(u.m).value, location.geocentric[2].to(u.m).value ] nants = antxyz.shape[0] diameters = diameter * numpy.ones(nants) anames = [names % ant for ant in range(nants)] mounts = numpy.repeat(mount, nants) antxyz, diameters, anames, mounts = limit_rmax(antxyz, diameters, anames, mounts, rmax) fc = Configuration(location=location, names=anames, mount=mounts, xyz=antxyz, diameter=diameters, name=name) return fc
def create_LOFAR_configuration(antfile: str, location, rmax=1e6) -> Configuration: """ Define configuration from the LOFAR configuration file :param antfile: :param location: EarthLocation :param rmax: Maximum distance from array centre (m) :return: Configuration """ check_data_directory() antxyz = numpy.genfromtxt(antfile, skip_header=2, usecols=[1, 2, 3], delimiter=",") nants = antxyz.shape[0] assert antxyz.shape[ 1] == 3, "Antenna array has wrong shape %s" % antxyz.shape anames = numpy.genfromtxt(antfile, dtype='str', skip_header=2, usecols=[0], delimiter=",") mounts = numpy.repeat('XY', nants) diameters = numpy.repeat(35.0, nants) antxyz, diameters, mounts, anames = limit_rmax(antxyz, diameters, anames, mounts, rmax) fc = Configuration(location=location, names=anames, mount=mounts, xyz=antxyz, diameter=diameters, name='LOFAR') return fc
def create_low_test_skycomponents_from_gleam(flux_limit=0.1, polarisation_frame=PolarisationFrame("stokesI"), frequency=numpy.array([1e8]), kind='cubic', phasecentre=None, radius=1.0) \ -> List[Skycomponent]: """Create sky components from the GLEAM survey Stokes I is estimated from a cubic spline fit to the measured fluxes. The polarised flux is always zero. See http://www.mwatelescope.org/science/gleam-survey The catalog is available from Vizier. VIII/100 GaLactic and Extragalactic All-sky MWA survey (Hurley-Walker+, 2016) GaLactic and Extragalactic All-sky Murchison Wide Field Array (GLEAM) survey. I: A low-frequency extragalactic catalogue. Hurley-Walker N., et al., Mon. Not. R. Astron. Soc., 464, 1146-1167 (2017), 2017MNRAS.464.1146H :param flux_limit: Only write components brighter than this (Jy) :param polarisation_frame: Polarisation frame (default PolarisationFrame("stokesI")) :param frequency: Frequencies at which the flux will be estimated :param kind: Kind of interpolation (see scipy.interpolate.interp1d) Default: linear :param phasecentre: Desired phase centre (SkyCoord) default None implies all sources :param radius: Radius of sources selected around phasecentre (default 1.0 rad) :return: List of Skycomponents """ check_data_directory() fitsfile = rascil_path("data/models/GLEAM_EGC.fits") rad2deg = 180.0 / numpy.pi decmin = phasecentre.dec.to('deg').value - rad2deg * radius / 2.0 decmax = phasecentre.dec.to('deg').value + rad2deg * radius / 2.0 hdulist = fits.open(fitsfile, lazy_load_hdus=False) recs = hdulist[1].data[0].array fluxes = recs['peak_flux_wide'] mask = fluxes > flux_limit filtered_recs = recs[mask] decs = filtered_recs['DEJ2000'] mask = decs > decmin filtered_recs = filtered_recs[mask] decs = filtered_recs['DEJ2000'] mask = decs < decmax filtered_recs = filtered_recs[mask] ras = filtered_recs['RAJ2000'] decs = filtered_recs['DEJ2000'] names = filtered_recs['Name'] if polarisation_frame is None: polarisation_frame = PolarisationFrame("stokesI") npol = polarisation_frame.npol nchan = len(frequency) # For every source, we read all measured fluxes and interpolate to the # required frequencies gleam_freqs = numpy.array([ 76, 84, 92, 99, 107, 115, 122, 130, 143, 151, 158, 166, 174, 181, 189, 197, 204, 212, 220, 227 ]) gleam_flux_freq = numpy.zeros([len(names), len(gleam_freqs)]) for i, f in enumerate(gleam_freqs): gleam_flux_freq[:, i] = filtered_recs['int_flux_%03d' % (f)][:] skycomps = [] directions = SkyCoord(ra=ras * u.deg, dec=decs * u.deg) if phasecentre is not None: separations = directions.separation(phasecentre).to('rad').value else: separations = numpy.zeros(len(names)) for isource, name in enumerate(names): direction = directions[isource] if separations[isource] < radius: fint = interpolate.interp1d(gleam_freqs * 1.0e6, gleam_flux_freq[isource, :], kind=kind) flux = numpy.zeros([nchan, npol]) flux[:, 0] = fint(frequency) if not numpy.isnan(flux).any(): skycomps.append( Skycomponent(direction=direction, flux=flux, frequency=frequency, name=name, shape='Point', polarisation_frame=polarisation_frame)) log.info( 'create_low_test_skycomponents_from_gleam: %d sources above flux limit %.3f' % (len(skycomps), flux_limit)) hdulist.close() return skycomps
def create_low_test_skymodel_from_gleam(npixel=512, polarisation_frame=PolarisationFrame( "stokesI"), cellsize=0.000015, frequency=numpy.array([1e8]), channel_bandwidth=numpy.array([1e6]), phasecentre=None, kind='cubic', applybeam=True, flux_limit=0.1, flux_max=numpy.inf, flux_threshold=1.0, insert_method='Nearest', telescope='LOW') -> SkyModel: """Create LOW test skymodel from the GLEAM survey Stokes I is estimated from a cubic spline fit to the measured fluxes. The polarised flux is always zero. See http://www.mwatelescope.org/science/gleam-survey The catalog is available from Vizier. VIII/100 GaLactic and Extragalactic All-sky MWA survey (Hurley-Walker+, 2016) GaLactic and Extragalactic All-sky Murchison Wide Field Array (GLEAM) survey. I: A low-frequency extragalactic catalogue. Hurley-Walker N., et al., Mon. Not. R. Astron. Soc., 464, 1146-1167 (2017), 2017MNRAS.464.1146H :param telescope: :param npixel: Number of pixels :param polarisation_frame: Polarisation frame (default PolarisationFrame("stokesI")) :param cellsize: cellsize in radians :param frequency: :param channel_bandwidth: Channel width (Hz) :param phasecentre: phasecentre (SkyCoord) :param kind: Kind of interpolation (see scipy.interpolate.interp1d) Default: cubic :param applybeam: Apply the primary beam? :param flux_limit: Weakest component :param flux_max: Maximum strength component to be included in components :param flux_threshold: Split between components (brighter) and image (weaker) :param insert_method: Nearest | PSWF | Lanczos :return: :return: SkyModel """ check_data_directory() if phasecentre is None: phasecentre = SkyCoord(ra=+15.0 * u.deg, dec=-35.0 * u.deg, frame='icrs', equinox='J2000') radius = npixel * cellsize sc = create_low_test_skycomponents_from_gleam( flux_limit=flux_limit, polarisation_frame=polarisation_frame, frequency=frequency, phasecentre=phasecentre, kind=kind, radius=radius) sc = filter_skycomponents_by_flux(sc, flux_max=flux_max) if polarisation_frame is None: polarisation_frame = PolarisationFrame("stokesI") npol = polarisation_frame.npol nchan = len(frequency) shape = [nchan, npol, npixel, npixel] w = WCS(naxis=4) # The negation in the longitude is needed by definition of RA, DEC w.wcs.cdelt = [ -cellsize * 180.0 / numpy.pi, cellsize * 180.0 / numpy.pi, 1.0, channel_bandwidth[0] ] w.wcs.crpix = [npixel // 2 + 1, npixel // 2 + 1, 1.0, 1.0] w.wcs.ctype = ["RA---SIN", "DEC--SIN", 'STOKES', 'FREQ'] w.wcs.crval = [phasecentre.ra.deg, phasecentre.dec.deg, 1.0, frequency[0]] w.naxis = 4 w.wcs.radesys = 'ICRS' w.wcs.equinox = 2000.0 model = create_image_from_array(numpy.zeros(shape), w, polarisation_frame=polarisation_frame) if applybeam: beam = create_pb(model, telescope=telescope, use_local=False) sc = apply_beam_to_skycomponent(sc, beam) weaksc = filter_skycomponents_by_flux(sc, flux_max=flux_threshold) brightsc = filter_skycomponents_by_flux(sc, flux_min=flux_threshold) model = insert_skycomponent(model, weaksc, insert_method=insert_method) log.info( 'create_low_test_skymodel_from_gleam: %d bright sources above flux threshold %.3f, %d weak sources below ' % (len(brightsc), flux_threshold, len(weaksc))) return SkyModel(components=brightsc, image=model, mask=None, gaintable=None)
def create_low_test_image_from_gleam(npixel=512, polarisation_frame=PolarisationFrame( "stokesI"), cellsize=0.000015, frequency=numpy.array([1e8]), channel_bandwidth=numpy.array([1e6]), phasecentre=None, kind='cubic', applybeam=False, flux_limit=0.1, flux_max=numpy.inf, flux_min=-numpy.inf, radius=None, insert_method='Nearest') -> Image: """Create LOW test image from the GLEAM survey Stokes I is estimated from a cubic spline fit to the measured fluxes. The polarised flux is always zero. See http://www.mwatelescope.org/science/gleam-survey The catalog is available from Vizier. VIII/100 GaLactic and Extragalactic All-sky MWA survey (Hurley-Walker+, 2016) GaLactic and Extragalactic All-sky Murchison Wide Field Array (GLEAM) survey. I: A low-frequency extragalactic catalogue. Hurley-Walker N., et al., Mon. Not. R. Astron. Soc., 464, 1146-1167 (2017), 2017MNRAS.464.1146H :param npixel: Number of pixels :param polarisation_frame: Polarisation frame (default PolarisationFrame("stokesI")) :param cellsize: cellsize in radians :param frequency: :param channel_bandwidth: Channel width (Hz) :param phasecentre: phasecentre (SkyCoord) :param kind: Kind of interpolation (see scipy.interpolate.interp1d) Default: linear :return: Image """ check_data_directory() if phasecentre is None: phasecentre = SkyCoord(ra=+15.0 * u.deg, dec=-35.0 * u.deg, frame='icrs', equinox='J2000') if radius is None: radius = npixel * cellsize / numpy.sqrt(2.0) sc = create_low_test_skycomponents_from_gleam( flux_limit=flux_limit, polarisation_frame=polarisation_frame, frequency=frequency, phasecentre=phasecentre, kind=kind, radius=radius) sc = filter_skycomponents_by_flux(sc, flux_min=flux_min, flux_max=flux_max) if polarisation_frame is None: polarisation_frame = PolarisationFrame("stokesI") npol = polarisation_frame.npol nchan = len(frequency) shape = [nchan, npol, npixel, npixel] w = WCS(naxis=4) # The negation in the longitude is needed by definition of RA, DEC w.wcs.cdelt = [ -cellsize * 180.0 / numpy.pi, cellsize * 180.0 / numpy.pi, 1.0, channel_bandwidth[0] ] w.wcs.crpix = [npixel // 2 + 1, npixel // 2 + 1, 1.0, 1.0] w.wcs.ctype = ["RA---SIN", "DEC--SIN", 'STOKES', 'FREQ'] w.wcs.crval = [phasecentre.ra.deg, phasecentre.dec.deg, 1.0, frequency[0]] w.naxis = 4 w.wcs.radesys = 'ICRS' w.wcs.equinox = 2000.0 model = create_image_from_array(numpy.zeros(shape), w, polarisation_frame=polarisation_frame) model = insert_skycomponent(model, sc, insert_method=insert_method) if applybeam: beam = create_pb(model, telescope='LOW', use_local=False) model.data[...] *= beam.data[...] return model
def create_test_skycomponents_from_s3(polarisation_frame=PolarisationFrame( "stokesI"), frequency=numpy.array([1e8]), channel_bandwidth=numpy.array([1e6]), phasecentre=None, fov=20, flux_limit=1e-3, radius=None): """Create test image from S3 The input catalog was generated at http://s-cubed.physics.ox.ac.uk/s3_sex using the following query:: Database: s3_sex SQL: select * from Galaxies where (pow(10,itot_151)*1000 > 1.0) and (right_ascension between -5 and 5) and (declination between -5 and 5);; Number of rows returned: 29966 For frequencies < 610MHz, there are three tables to use:: data/models/S3_151MHz_10deg.csv, use fov=10 data/models/S3_151MHz_20deg.csv, use fov=20 data/models/S3_151MHz_40deg.csv, use fov=40 For frequencies > 610MHz, there are three tables: data/models/S3_1400MHz_1mJy_10deg.csv, use flux_limit>= 1e-3 data/models/S3_1400MHz_100uJy_10deg.csv, use flux_limit < 1e-3 data/models/S3_1400MHz_1mJy_18deg.csv, use flux_limit>= 1e-3 data/models/S3_1400MHz_100uJy_18deg.csv, use flux_limit < 1e-3 The component spectral index is calculated from the 610MHz and 151MHz or 1400MHz and 610MHz, and then calculated for the specified frequencies. If polarisation_frame is not stokesI then the image will a polarised axis but the values will be zero. :param polarisation_frame: Polarisation frame (default PolarisationFrame("stokesI")) :param frequency: :param channel_bandwidth: Channel width (Hz) :param phasecentre: phasecentre (SkyCoord) :param fov: fov 10 | 20 | 40 :param flux_limit: Minimum flux (Jy) :return: Image """ check_data_directory() ras = [] decs = [] fluxes = [] names = [] if phasecentre is None: phasecentre = SkyCoord(ra=+180.0 * u.deg, dec=-60.0 * u.deg, frame='icrs', equinox='J2000') if polarisation_frame is None: polarisation_frame = PolarisationFrame("stokesI") if numpy.max(frequency) > 6.1E8: if fov > 10: fovstr = '18' else: fovstr = '10' if flux_limit >= 1e-3: csvfilename = rascil_data_path('models/S3_1400MHz_1mJy_%sdeg.csv' % fovstr) else: csvfilename = rascil_data_path( 'models/S3_1400MHz_100uJy_%sdeg.csv' % fovstr) log.info( 'create_test_skycomponents_from_s3: Reading S3-SEX sources from %s ' % csvfilename) else: assert fov in [ 10, 20, 40 ], "Field of view invalid: use one of %s" % ([10, 20, 40]) csvfilename = rascil_data_path('models/S3_151MHz_%ddeg.csv' % (fov)) log.info( 'create_test_skycomponents_from_s3: Reading S3-SEX sources from %s ' % csvfilename) skycomps = list() with open(csvfilename) as csvfile: readCSV = csv.reader(csvfile, delimiter=',') r = 0 for row in readCSV: # Skip first row if r > 0: ra = float(row[4]) / numpy.cos( phasecentre.dec.rad) + phasecentre.ra.deg dec = float(row[5]) + phasecentre.dec.deg if numpy.max(frequency) > 6.1E8: alpha = (float(row[11]) - float(row[10])) / numpy.log10( 1400.0 / 610.0) flux = numpy.power(10, float(row[10])) * numpy.power( frequency / 1.4e9, alpha) else: alpha = (float(row[10]) - float(row[9])) / numpy.log10( 610.0 / 151.0) flux = numpy.power(10, float(row[9])) * numpy.power( frequency / 1.51e8, alpha) if numpy.max(flux) > flux_limit: ras.append(ra) decs.append(dec) if polarisation_frame == PolarisationFrame("stokesIQUV"): polscale = numpy.array([1.0, 0.0, 0.0, 0.0]) fluxes.append(numpy.outer(flux, polscale)) else: fluxes.append([[f] for f in flux]) names.append("S3_%s" % row[0]) r += 1 csvfile.close() assert len(fluxes) > 0, "No sources found above flux limit %s" % flux_limit directions = SkyCoord(ra=ras * u.deg, dec=decs * u.deg) if phasecentre is not None: separations = directions.separation(phasecentre).to('rad').value else: separations = numpy.zeros(len(names)) for isource, name in enumerate(names): direction = directions[isource] if separations[isource] < radius: if not numpy.isnan(flux).any(): skycomps.append( Skycomponent(direction=direction, flux=fluxes[isource], frequency=frequency, name=names[isource], shape='Point', polarisation_frame=polarisation_frame)) log.info( 'create_test_skycomponents_from_s3: %d sources found above fluxlimit inside search radius' % len(skycomps)) return skycomps
def create_test_image_from_s3(npixel=16384, polarisation_frame=PolarisationFrame("stokesI"), cellsize=0.000015, frequency=numpy.array([1e8]), channel_bandwidth=numpy.array([1e6]), phasecentre=None, fov=20, flux_limit=1e-3) -> Image: """Create MID test image from S3 The input catalog was generated at http://s-cubed.physics.ox.ac.uk/s3_sex using the following query:: Database: s3_sex SQL: select * from Galaxies where (pow(10,itot_151)*1000 > 1.0) and (right_ascension between -5 and 5) and (declination between -5 and 5);; Number of rows returned: 29966 For frequencies < 610MHz, there are three tables to use:: data/models/S3_151MHz_10deg.csv, use fov=10 data/models/S3_151MHz_20deg.csv, use fov=20 data/models/S3_151MHz_40deg.csv, use fov=40 For frequencies > 610MHz, there are three tables: data/models/S3_1400MHz_1mJy_10deg.csv, use flux_limit>= 1e-3 data/models/S3_1400MHz_100uJy_10deg.csv, use flux_limit < 1e-3 data/models/S3_1400MHz_1mJy_18deg.csv, use flux_limit>= 1e-3 data/models/S3_1400MHz_100uJy_18deg.csv, use flux_limit < 1e-3 The component spectral index is calculated from the 610MHz and 151MHz or 1400MHz and 610MHz, and then calculated for the specified frequencies. If polarisation_frame is not stokesI then the image will a polarised axis but the values will be zero. :param npixel: Number of pixels :param polarisation_frame: Polarisation frame (default PolarisationFrame("stokesI")) :param cellsize: cellsize in radians :param frequency: :param channel_bandwidth: Channel width (Hz) :param phasecentre: phasecentre (SkyCoord) :param fov: fov 10 | 20 | 40 :param flux_limit: Minimum flux (Jy) :return: Image """ check_data_directory() ras = [] decs = [] fluxes = [] if phasecentre is None: phasecentre = SkyCoord(ra=+180.0 * u.deg, dec=-60.0 * u.deg, frame='icrs', equinox='J2000') if polarisation_frame is None: polarisation_frame = PolarisationFrame("stokesI") npol = polarisation_frame.npol nchan = len(frequency) shape = [nchan, npol, npixel, npixel] w = WCS(naxis=4) # The negation in the longitude is needed by definition of RA, DEC w.wcs.cdelt = [ -cellsize * 180.0 / numpy.pi, cellsize * 180.0 / numpy.pi, 1.0, channel_bandwidth[0] ] w.wcs.crpix = [npixel // 2 + 1, npixel // 2 + 1, 1.0, 1.0] w.wcs.ctype = ["RA---SIN", "DEC--SIN", 'STOKES', 'FREQ'] w.wcs.crval = [phasecentre.ra.deg, phasecentre.dec.deg, 1.0, frequency[0]] w.naxis = 4 w.wcs.radesys = 'ICRS' w.wcs.equinox = 2000.0 model = create_image_from_array(numpy.zeros(shape), w, polarisation_frame=polarisation_frame) if numpy.max(frequency) > 6.1E8: if fov > 10: fovstr = '18' else: fovstr = '10' if flux_limit >= 1e-3: csvfilename = rascil_data_path('models/S3_1400MHz_1mJy_%sdeg.csv' % fovstr) else: csvfilename = rascil_data_path( 'models/S3_1400MHz_100uJy_%sdeg.csv' % fovstr) log.info('create_test_image_from_s3: Reading S3 sources from %s ' % csvfilename) else: assert fov in [ 10, 20, 40 ], "Field of view invalid: use one of %s" % ([10, 20, 40]) csvfilename = rascil_data_path('models/S3_151MHz_%ddeg.csv' % (fov)) log.info('create_test_image_from_s3: Reading S3 sources from %s ' % csvfilename) with open(csvfilename) as csvfile: readCSV = csv.reader(csvfile, delimiter=',') r = 0 for row in readCSV: # Skip first row if r > 0: ra = float(row[4]) + phasecentre.ra.deg dec = float(row[5]) + phasecentre.dec.deg if numpy.max(frequency) > 6.1E8: alpha = (float(row[11]) - float(row[10])) / numpy.log10( 1400.0 / 610.0) flux = numpy.power(10, float(row[10])) * numpy.power( frequency / 1.4e9, alpha) else: alpha = (float(row[10]) - float(row[9])) / numpy.log10( 610.0 / 151.0) flux = numpy.power(10, float(row[9])) * numpy.power( frequency / 1.51e8, alpha) if numpy.max(flux) > flux_limit: ras.append(ra) decs.append(dec) fluxes.append(flux) r += 1 csvfile.close() assert len(fluxes) > 0, "No sources found above flux limit %s" % flux_limit log.info('create_test_image_from_s3: %d sources read' % (len(fluxes))) p = w.sub(2).wcs_world2pix(numpy.array(ras), numpy.array(decs), 1) fluxes = numpy.array(fluxes) total_flux = numpy.sum(fluxes) ip = numpy.round(p).astype('int') ok = numpy.where((0 <= ip[0, :]) & (npixel > ip[0, :]) & (0 <= ip[1, :]) & (npixel > ip[1, :]))[0] ps = ip[:, ok] fluxes = fluxes[ok] actual_flux = numpy.sum(fluxes) log.info('create_test_image_from_s3: %d sources inside the image' % (ps.shape[1])) log.info( 'create_test_image_from_s3: average channel flux in S3 model = %.3f, actual average channel flux in ' 'image = %.3f' % (total_flux / float(nchan), actual_flux / float(nchan))) for chan in range(nchan): for iflux, flux in enumerate(fluxes): model.data[chan, 0, ps[1, iflux], ps[0, iflux]] = flux[chan] return model
from rascil.data_models.polarisation import PolarisationFrame from rascil.processing_components.calibration.chain_calibration import create_calibration_controls from rascil.processing_components.calibration.operations import create_gaintable_from_blockvisibility, apply_gaintable from rascil.processing_components.image.operations import create_image_from_array from rascil.processing_components.image.operations import import_image_from_fits from rascil.processing_components.imaging import predict_2d, dft_skycomponent_visibility, \ create_image_from_visibility, advise_wide_field from rascil.processing_components.imaging.primary_beams import create_pb from rascil.processing_components.skycomponent.operations import create_skycomponent, insert_skycomponent, \ apply_beam_to_skycomponent, filter_skycomponents_by_flux from rascil.processing_components.visibility.base import create_blockvisibility, create_visibility from rascil.processing_components.visibility.coalesce import convert_blockvisibility_to_visibility, \ convert_visibility_to_blockvisibility from rascil.processing_components.util.installation_checks import check_data_directory check_data_directory() log = logging.getLogger('logger') def create_test_image( canonical=True, cellsize=None, frequency=None, channel_bandwidth=None, phasecentre=None, polarisation_frame=PolarisationFrame("stokesI")) -> Image: """Create a useful test image This is the test image M31 widely used in ALMA and other simulations. It is actually part of an Halpha region in M31.
def create_named_configuration(name: str = 'LOWBD2', **kwargs) -> Configuration: """ Create standard configurations e.g. LOWBD2, MIDBD2 Possible configurations are:: LOWBD1 LOWBD2 LOWBD2-core LOW == LOWR3 MID == MIDR5 ASKAP LOFAR VLAA VLAA_north :param name: name of Configuration LOWBD2, LOWBD1, LOFAR, VLAA, ASKAP :param rmax: Maximum distance of station from the average (m) :return: For LOWBD2, setting rmax gives the following number of stations 100.0 13 300.0 94 1000.0 251 3000.0 314 10000.0 398 30000.0 476 100000.0 512 """ check_data_directory() if name == 'LOWBD2': location = EarthLocation(lon="116.76444824", lat="-26.824722084", height=300.0) log.info("create_named_configuration: %s\n\t%s\n\t%s" % (name, location.geocentric, location.geodetic)) fc = create_configuration_from_file( antfile=rascil_path("data/configurations/LOWBD2.csv"), location=location, mount='xy', names='LOWBD2_%d', diameter=35.0, name=name, **kwargs) elif name == 'LOWBD1': location = EarthLocation(lon="116.76444824", lat="-26.824722084", height=300.0) log.info("create_named_configuration: %s\n\t%s\n\t%s" % (name, location.geocentric, location.geodetic)) fc = create_configuration_from_file( antfile=rascil_path("data/configurations/LOWBD1.csv"), location=location, mount='xy', names='LOWBD1_%d', diameter=35.0, name=name, **kwargs) elif name == 'LOWBD2-CORE': location = EarthLocation(lon="116.76444824", lat="-26.824722084", height=300.0) log.info("create_named_configuration: %s\n\t%s\n\t%s" % (name, location.geocentric, location.geodetic)) fc = create_configuration_from_file( antfile=rascil_path("data/configurations/LOWBD2-CORE.csv"), location=location, mount='xy', names='LOWBD2_%d', diameter=35.0, name=name, **kwargs) elif (name == 'LOW') or (name == 'LOWR3'): location = EarthLocation(lon="116.76444824", lat="-26.824722084", height=300.0) log.info("create_named_configuration: %s\n\t%s\n\t%s" % (name, location.geocentric, location.geodetic)) fc = create_configuration_from_MIDfile( antfile=rascil_path("data/configurations/ska1low_local.cfg"), mount='xy', name=name, location=location, **kwargs) elif (name == 'MID') or (name == "MIDR5"): location = EarthLocation(lon="21.443803", lat="-30.712925", height=0.0) log.info("create_named_configuration: %s\n\t%s\n\t%s" % (name, location.geocentric, location.geodetic)) fc = create_configuration_from_MIDfile( antfile=rascil_path("data/configurations/ska1mid_local.cfg"), mount='azel', name=name, location=location, **kwargs) elif name == 'ASKAP': location = EarthLocation(lon="+116.6356824", lat="-26.7013006", height=377.0) log.info("create_named_configuration: %s\n\t%s\n\t%s" % (name, location.geocentric, location.geodetic)) fc = create_configuration_from_file( antfile=rascil_path("data/configurations/A27CR3P6B.in.csv"), mount='equatorial', names='ASKAP_%d', diameter=12.0, name=name, location=location, **kwargs) elif name == 'LOFAR': location = EarthLocation(x=3826923.9 * u.m, y=460915.1 * u.m, z=5064643.2 * u.m) log.info("create_named_configuration: %s\n\t%s\n\t%s" % (name, location.geocentric, location.geodetic)) assert get_parameter(kwargs, "meta", False) is False fc = create_LOFAR_configuration( antfile=rascil_path("data/configurations/LOFAR.csv"), location=location) elif name == 'VLAA': location = EarthLocation(lon="-107.6184", lat="34.0784", height=2124.0) log.info("create_named_configuration: %s\n\t%s\n\t%s" % (name, location.geocentric, location.geodetic)) fc = create_configuration_from_file( antfile=rascil_path("data/configurations/VLA_A_hor_xyz.csv"), location=location, mount='azel', names='VLA_%d', diameter=25.0, name=name, **kwargs) elif name == 'VLAA_north': location = EarthLocation(lon="-107.6184", lat="90.000", height=0.0) log.info("create_named_configuration: %s\n\t%s\n\t%s" % (name, location.geocentric, location.geodetic)) fc = create_configuration_from_file( antfile=rascil_path("data/configurations/VLA_A_hor_xyz.csv"), location=location, mount='azel', names='VLA_%d', diameter=25.0, name=name, **kwargs) else: raise ValueError("No such Configuration %s" % name) return fc