def phaserotate_visibility(vis: Visibility, newphasecentre: SkyCoord, tangent=True, inverse=False) -> Visibility: """ Phase rotate from the current phase centre to a new phase centre If tangent is False the uvw are recomputed and the visibility phasecentre is updated. Otherwise only the visibility phases are adjusted :param vis: Visibility to be rotated :param newphasecentre: :param tangent: Stay on the same tangent plane? (True) :param inverse: Actually do the opposite :return: Visibility """ assert isinstance(vis, Visibility), "vis is not a Visibility: %r" % vis l, m, n = skycoord_to_lmn(newphasecentre, vis.phasecentre) # No significant change? if numpy.abs(n) > 1e-15: # Make a new copy newvis = copy_visibility(vis) phasor = simulate_point(newvis.uvw, l, m) nvis, npol = vis.vis.shape # TODO: Speed up (broadcast rules not obvious to me) if inverse: for i in range(nvis): for pol in range(npol): newvis.data['vis'][i, pol] *= phasor[i] else: for i in range(nvis): for pol in range(npol): newvis.data['vis'][i, pol] *= numpy.conj(phasor[i]) # To rotate UVW, rotate into the global XYZ coordinate system and back. We have the option of # staying on the tangent plane or not. If we stay on the tangent then the raster will # join smoothly at the edges. If we change the tangent then we will have to reproject to get # the results on the same image, in which case overlaps or gaps are difficult to deal with. if not tangent: if inverse: xyz = uvw_to_xyz(vis.data['uvw'], ha=-newvis.phasecentre.ra.rad, dec=newvis.phasecentre.dec.rad) newvis.data['uvw'][...] = \ xyz_to_uvw(xyz, ha=-newphasecentre.ra.rad, dec=newphasecentre.dec.rad)[...] else: # This is the original (non-inverse) code xyz = uvw_to_xyz(newvis.data['uvw'], ha=-newvis.phasecentre.ra.rad, dec=newvis.phasecentre.dec.rad) newvis.data['uvw'][...] = xyz_to_uvw(xyz, ha=-newphasecentre.ra.rad, dec=newphasecentre.dec.rad)[ ...] newvis.phasecentre = newphasecentre return newvis else: return vis
def test_phase_rotate(self): uvw = numpy.array([(1, 0, 0), (0, 1, 0), (0, 0, 1)]) pos = [ SkyCoord(17, 35, unit=u.deg), SkyCoord(17, 30, unit=u.deg), SkyCoord(12, 30, unit=u.deg), SkyCoord(11, 35, unit=u.deg), SkyCoord(51, 35, unit=u.deg), SkyCoord(15, 70, unit=u.deg) ] # Sky coordinates to reproject to for phasecentre in pos: for newphasecentre in pos: # Rotate UVW xyz = uvw_to_xyz(uvw, -phasecentre.ra.rad, phasecentre.dec.rad) uvw_rotated = xyz_to_uvw(xyz, -newphasecentre.ra.rad, newphasecentre.dec.rad) # Determine phasor l_p, m_p, n_p = skycoord_to_lmn(phasecentre, newphasecentre) phasor = simulate_point(uvw_rotated, l_p, m_p) for sourcepos in pos: # Simulate visibility at old and new phase centre l, m, _ = skycoord_to_lmn(sourcepos, phasecentre) vis = simulate_point(uvw, l, m) l_r, m_r, _ = skycoord_to_lmn(sourcepos, newphasecentre) vis_rotated = simulate_point(uvw_rotated, l_r, m_r) # Difference should be given by phasor assert_allclose(vis * phasor, vis_rotated, atol=1e-10)
def transform(x, y, z, ha, dec): """ :param x: :param y: :param z: :param ha: :param dec: :return: """ res = xyz_to_uvw(numpy.array([x, y, z]), numpy.radians(ha), numpy.radians(dec)) assert_allclose(numpy.linalg.norm(res), numpy.linalg.norm([x, y, z])) assert_allclose(uvw_to_xyz(res, numpy.radians(ha), numpy.radians(dec)), [x, y, z]) return res
def create_visibility(config: Configuration, times: numpy.array, frequency: numpy.array, channel_bandwidth, phasecentre: SkyCoord, weight: float, polarisation_frame=PolarisationFrame('stokesI'), integration_time=1.0) -> Visibility: """ Create a Visibility from Configuration, hour angles, and direction of source Note that we keep track of the integration time for BDA purposes :param config: Configuration of antennas :param times: hour angles in radians :param frequency: frequencies (Hz] [nchan] :param weight: weight of a single sample :param phasecentre: phasecentre of observation :param channel_bandwidth: channel bandwidths: (Hz] [nchan] :param integration_time: Integration time ('auto' or value in s) :param polarisation_frame: PolarisationFrame('stokesI') :return: Visibility """ assert phasecentre is not None, "Must specify phase centre" if polarisation_frame is None: polarisation_frame = correlate_polarisation(config.receptor_frame) nch = len(frequency) ants_xyz = config.data['xyz'] nants = len(config.data['names']) nbaselines = int(nants * (nants - 1) / 2) ntimes = len(times) npol = polarisation_frame.npol nrows = nbaselines * ntimes * nch nrowsperintegration = nbaselines * nch row = 0 rvis = numpy.zeros([nrows, npol], dtype='complex') rweight = weight * numpy.ones([nrows, npol]) rtimes = numpy.zeros([nrows]) rfrequency = numpy.zeros([nrows]) rchannel_bandwidth = numpy.zeros([nrows]) rantenna1 = numpy.zeros([nrows], dtype='int') rantenna2 = numpy.zeros([nrows], dtype='int') ruvw = numpy.zeros([nrows, 3]) # Do each hour angle in turn for iha, ha in enumerate(times): # Calculate the positions of the antennas as seen for this hour angle # and declination ant_pos = xyz_to_uvw(ants_xyz, ha, phasecentre.dec.rad) rtimes[row:row + nrowsperintegration] = ha * 43200.0 / numpy.pi # Loop over all pairs of antennas. Note that a2>a1 for a1 in range(nants): for a2 in range(a1 + 1, nants): rantenna1[row:row + nch] = a1 rantenna2[row:row + nch] = a2 # Loop over all frequencies and polarisations for ch in range(nch): # noinspection PyUnresolvedReferences k = frequency[ch] / constants.c.value ruvw[row, :] = (ant_pos[a2, :] - ant_pos[a1, :]) * k rfrequency[row] = frequency[ch] rchannel_bandwidth[row] = channel_bandwidth[ch] row += 1 assert row == nrows rintegration_time = numpy.full_like(rtimes, integration_time) vis = Visibility(uvw=ruvw, time=rtimes, antenna1=rantenna1, antenna2=rantenna2, frequency=rfrequency, vis=rvis, weight=rweight, imaging_weight=rweight, integration_time=rintegration_time, channel_bandwidth=rchannel_bandwidth, polarisation_frame=polarisation_frame) vis.phasecentre = phasecentre vis.configuration = config log.info("create_visibility: %s" % (vis_summary(vis))) assert isinstance(vis, Visibility), "vis is not a Visibility: %r" % vis return vis
def create_blockvisibility(config: Configuration, times: numpy.array, frequency: numpy.array, phasecentre: SkyCoord, weight: float = 1.0, polarisation_frame: PolarisationFrame = None, integration_time=1.0, channel_bandwidth=1e6, zerow=False, **kwargs) -> BlockVisibility: """ Create a BlockVisibility from Configuration, hour angles, and direction of source Note that we keep track of the integration time for BDA purposes :param config: Configuration of antennas :param times: hour angles in radians :param frequency: frequencies (Hz] [nchan] :param weight: weight of a single sample :param phasecentre: phasecentre of observation :param channel_bandwidth: channel bandwidths: (Hz] [nchan] :param integration_time: Integration time ('auto' or value in s) :param polarisation_frame: :return: BlockVisibility """ assert phasecentre is not None, "Must specify phase centre" if polarisation_frame is None: polarisation_frame = correlate_polarisation(config.receptor_frame) nch = len(frequency) ants_xyz = config.data['xyz'] nants = len(config.data['names']) ntimes = len(times) npol = polarisation_frame.npol visshape = [ntimes, nants, nants, nch, npol] rvis = numpy.zeros(visshape, dtype='complex') rweight = weight * numpy.ones(visshape) rtimes = numpy.zeros([ntimes]) ruvw = numpy.zeros([ntimes, nants, nants, 3]) # Do each hour angle in turn for iha, ha in enumerate(times): # Calculate the positions of the antennas as seen for this hour angle # and declination ant_pos = xyz_to_uvw(ants_xyz, ha, phasecentre.dec.rad) rtimes[iha] = ha * 43200.0 / numpy.pi # Loop over all pairs of antennas. Note that a2>a1 for a1 in range(nants): for a2 in range(a1 + 1, nants): ruvw[iha, a2, a1, :] = (ant_pos[a2, :] - ant_pos[a1, :]) ruvw[iha, a1, a2, :] = (ant_pos[a1, :] - ant_pos[a2, :]) rintegration_time = numpy.full_like(rtimes, integration_time) rchannel_bandwidth = numpy.full_like(frequency, channel_bandwidth) if zerow: ruvw[..., 2] = 0.0 vis = BlockVisibility(uvw=ruvw, time=rtimes, frequency=frequency, vis=rvis, weight=rweight, integration_time=rintegration_time, channel_bandwidth=rchannel_bandwidth, polarisation_frame=polarisation_frame) vis.phasecentre = phasecentre vis.configuration = config log.info("create_blockvisibility: %s" % (vis_summary(vis))) assert isinstance( vis, BlockVisibility), "vis is not a BlockVisibility: %r" % vis return vis