def test_coalesce_decoalesce_with_iter(self):
     for rows in vis_timeslice_iter(self.blockvis):
         visslice = create_visibility_from_rows(self.blockvis, rows)
         cvisslice = convert_blockvisibility_to_visibility(visslice)
         assert numpy.min(cvisslice.frequency) == numpy.min(self.frequency)
         assert numpy.min(cvisslice.frequency) > 0.0
         dvisslice = decoalesce_visibility(cvisslice)
         assert dvisslice.nvis == visslice.nvis
Beispiel #2
0
 def test_vis_timeslice_iterator(self):
     self.actualSetUp()
     nchunks = vis_timeslices(self.vis, timeslice='auto')
     log.debug('Found %d chunks' % (nchunks))
     assert nchunks > 1
     total_rows = 0
     for chunk, rows in enumerate(vis_timeslice_iter(self.vis, nchunks)):
         visslice = create_visibility_from_rows(self.vis, rows)
         total_rows += visslice.nvis
         assert visslice.vis[0].real == visslice.time[0]
         assert len(rows)
         assert numpy.sum(rows) < self.vis.nvis
     assert total_rows == self.vis.nvis, "Total rows iterated %d, Original rows %d" % (total_rows, self.vis.nvis)
Beispiel #3
0
def create_gaintable_from_screen(vis,
                                 sc,
                                 screen,
                                 height=3e5,
                                 vis_slices=None,
                                 scale=1.0,
                                 **kwargs):
    """ Create gaintables from a screen calculated using ARatmospy

    Screen axes are ['XX', 'YY', 'TIME', 'FREQ']

    :param vis:
    :param sc: Sky components for which pierce points are needed
    :param screen:
    :param height: Height (in m) of screen above telescope e.g. 3e5
    :param scale: Multiply the screen by this factor
    :return:
    """
    assert isinstance(vis, BlockVisibility)

    station_locations = vis.configuration.xyz

    # Convert to TEC
    # phase = image[pixel] * -8.44797245e9 / frequency
    nant = station_locations.shape[0]
    t2r = numpy.pi / 43200.0
    gaintables = [
        create_gaintable_from_blockvisibility(vis, **kwargs) for i in sc
    ]

    # The time in the Visibility is hour angle in seconds!
    for iha, rows in enumerate(vis_timeslice_iter(vis, vis_slices=vis_slices)):
        v = create_visibility_from_rows(vis, rows)
        ha = numpy.average(v.time)
        tec2phase = -8.44797245e9 / numpy.array(vis.frequency)

        number_bad = 0
        number_good = 0

        for icomp, comp in enumerate(sc):
            pp = find_pierce_points(station_locations,
                                    (comp.direction.ra.rad + t2r * ha) * u.rad,
                                    comp.direction.dec,
                                    height=height,
                                    phasecentre=vis.phasecentre)
            scr = numpy.zeros([nant])
            for ant in range(nant):
                pp0 = pp[ant][0:2]
                worldloc = [pp0[0], pp0[1], ha, 1e8]
                try:
                    pixloc = screen.wcs.wcs_world2pix([worldloc],
                                                      0)[0].astype('int')
                    scr[ant] = scale * screen.data[pixloc[3], pixloc[2],
                                                   pixloc[1], pixloc[0]]
                    number_good += 1
                except:
                    number_bad += 1
                    scr[ant] = 0.0

            scr = scr[:, numpy.newaxis] * tec2phase[numpy.newaxis, :]
            # axes of gaintable.gain are time, ant, nchan, nrec
            gaintables[icomp].gain[iha, :, :, :] = numpy.exp(
                1j * scr[..., numpy.newaxis, numpy.newaxis])
            gaintables[icomp].phasecentre = comp.direction

        if number_bad > 0:
            log.warning(
                "create_gaintable_from_screen: %d pierce points are inside the screen image"
                % (number_good))
            log.warning(
                "create_gaintable_from_screen: %d pierce points are outside the screen image"
                % (number_bad))

    return gaintables
Beispiel #4
0
def apply_gaintable(vis: BlockVisibility,
                    gt: GainTable,
                    inverse=False,
                    vis_slices=None,
                    **kwargs) -> BlockVisibility:
    """Apply a gain table to a block visibility
    
    The corrected visibility is::
    
        V_corrected = {g_i * g_j^*}^-1 V_obs
        
    If the visibility data are polarised e.g. polarisation_frame("linear") then the inverse operator
    represents an actual inverse of the gains.
    
    :param vis: Visibility to have gains applied
    :param gt: Gaintable to be applied
    :param inverse: Apply the inverse (default=False)
    :return: input vis with gains applied
    
    """
    assert isinstance(
        vis, BlockVisibility), "vis is not a BlockVisibility: %r" % vis
    assert isinstance(gt, GainTable), "gt is not a GainTable: %r" % gt

    assert_vis_gt_compatible(vis, gt)

    if inverse:
        log.debug('apply_gaintable: Apply inverse gaintable')
    else:
        log.debug('apply_gaintable: Apply gaintable')

    is_scalar = gt.gain.shape[-2:] == (1, 1)
    if is_scalar:
        log.debug('apply_gaintable: scalar gains')

    for chunk, rows in enumerate(vis_timeslice_iter(vis,
                                                    vis_slices=vis_slices)):
        if numpy.sum(rows) > 0:
            vistime = numpy.average(vis.time[rows])
            gaintable_rows = abs(gt.time - vistime) < gt.interval / 2.0

            # Lookup the gain for this set of visibilities
            gain = gt.data['gain'][gaintable_rows]
            gainwt = gt.data['weight'][gaintable_rows]

            # The shape of the mueller matrix is
            ntimes, nant, nchan, nrec, _ = gain.shape

            original = vis.vis[rows]
            originalwt = vis.weight[rows]
            applied = copy.deepcopy(original)
            appliedwt = copy.deepcopy(originalwt)
            for time in range(ntimes):
                antantwt = numpy.outer(gainwt, gainwt)
                if is_scalar:
                    if inverse:
                        lgain = numpy.ones_like(gain)
                        lgain[numpy.abs(gain) >
                              0.0] = 1.0 / gain[numpy.abs(gain) > 0.0]
                    else:
                        lgain = gain
                    clgain = numpy.conjugate(lgain)
                    for chan in range(nchan):
                        smueller = numpy.ma.outer(
                            lgain[time, :, chan, 0],
                            clgain[time, :, chan, 0]).reshape([nant, nant])
                        applied[time, :, :, chan,
                                0] = original[time, :, :, chan, 0] * smueller
                        antantwt = numpy.outer(gainwt[time, :, chan, 0, 0],
                                               gainwt[time, :, chan, 0, 0])
                        applied[time, :, :, chan, 0][antantwt == 0.0] = 0.0
                        appliedwt[time, :, :, chan, 0][antantwt == 0.0] = 0.0
                else:
                    for a1 in range(vis.nants - 1):
                        for a2 in range(a1 + 1, vis.nants):
                            for chan in range(nchan):
                                mueller = numpy.kron(
                                    gain[time, a1, chan, :, :],
                                    numpy.conjugate(gain[time, a2,
                                                         chan, :, :]))
                                if inverse:
                                    # If the Mueller is singular, ignore it
                                    try:
                                        mueller = numpy.linalg.inv(mueller)
                                        applied[time, a2, a1, chan, :] = \
                                            numpy.matmul(mueller, original[time, a2, a1, chan, :])
                                    except numpy.linalg.linalg.LinAlgError:
                                        applied[time, a2, a1, chan, :] = \
                                            original[time, a2, a1, chan, :]
                                else:
                                    applied[time, a2, a1, chan, :] = \
                                        numpy.matmul(mueller, original[time, a2, a1, chan, :])
                                if (gainwt[time, a1, chan, 0, 0] <= 0.0) or (
                                        gainwt[time, a1, chan, 0, 0] <= 0.0):
                                    applied[time, a2, a1, chan, 0] = 0.0
                                    appliedwt[time, a2, a1, chan, 0] = 0.0

            vis.data['vis'][rows] = applied
    return vis
Beispiel #5
0
def create_gaintable_from_screen(vis,
                                 sc,
                                 screen,
                                 height=3e5,
                                 vis_slices=None,
                                 scale=1.0,
                                 r0=5e3,
                                 type_atmosphere='ionosphere',
                                 **kwargs):
    """ Create gaintables from a screen calculated using ARatmospy

    Screen axes are ['XX', 'YY', 'TIME', 'FREQ']

    :param vis:
    :param sc: Sky components for which pierce points are needed
    :param screen:
    :param height: Height (in m) of screen above telescope e.g. 3e5
    :param r0: r0 in meters
    :param type_atmosphere: 'ionosphere' or 'troposphere'
    :param scale: Multiply the screen by this factor
    :return:
    """
    assert isinstance(vis, BlockVisibility)

    station_locations = vis.configuration.xyz

    scale = numpy.power(r0 / 5000.0, -5.0 / 3.0)
    if type_atmosphere == 'troposphere':
        # In troposphere files, the units are phase in radians.
        screen_to_phase = scale
    else:
        # In the ionosphere file, the units are dTEC.
        screen_to_phase = -scale * 8.44797245e9 / numpy.array(vis.frequency)

    nant = station_locations.shape[0]
    t2r = numpy.pi / 43200.0
    gaintables = [
        create_gaintable_from_blockvisibility(vis, **kwargs) for i in sc
    ]

    number_bad = 0
    number_good = 0

    ha_zero = numpy.average(calculate_blockvisibility_hourangles(vis))
    for iha, rows in enumerate(vis_timeslice_iter(vis, vis_slices=vis_slices)):
        v = create_visibility_from_rows(vis, rows)
        ha = numpy.average(calculate_blockvisibility_hourangles(v) -
                           ha_zero).to('rad').value
        for icomp, comp in enumerate(sc):
            pp = find_pierce_points(station_locations,
                                    (comp.direction.ra.rad + t2r * ha) *
                                    units.rad,
                                    comp.direction.dec,
                                    height=height,
                                    phasecentre=vis.phasecentre)
            scr = numpy.zeros([nant])
            for ant in range(nant):
                pp0 = pp[ant][0:2]
                # Using narrow band approach - we should loop over frequency
                try:
                    worldloc = [
                        pp0[0], pp0[1], ha,
                        numpy.average(vis.frequency)
                    ]
                    pixloc = screen.wcs.wcs_world2pix([worldloc],
                                                      0)[0].astype('int')
                    if type_atmosphere == 'troposphere':
                        pixloc[3] = 0
                    scr[ant] = screen_to_phase * screen.data[
                        pixloc[3], pixloc[2], pixloc[1], pixloc[0]]
                    number_good += 1
                except (ValueError, IndexError):
                    number_bad += 1
                    scr[ant] = 0.0

            # axes of gaintable.gain are time, ant, nchan, nrec
            gaintables[icomp].gain[iha, :, :, :] = numpy.exp(
                1j * scr)[..., numpy.newaxis, numpy.newaxis, numpy.newaxis]
            gaintables[icomp].phasecentre = comp.direction

    assert number_good > 0, "create_gaintable_from_screen: There are no pierce points inside the atmospheric screen image"
    if number_bad > 0:
        log.warning(
            "create_gaintable_from_screen: %d pierce points are inside the atmospheric screen image"
            % (number_good))
        log.warning(
            "create_gaintable_from_screen: %d pierce points are outside the atmospheric screen image"
            % (number_bad))

    return gaintables
Beispiel #6
0
def grid_gaintable_to_screen(vis,
                             gaintables,
                             screen,
                             height=3e5,
                             gaintable_slices=None,
                             scale=1.0,
                             r0=5e3,
                             type_atmosphere='ionosphere',
                             vis_slices=None,
                             **kwargs):
    """ Grid a gaintable to a screen image

    Screen axes are ['XX', 'YY', 'TIME', 'FREQ']

    The phases are just averaged per grid cell, no phase unwrapping is performed.

    :param vis:
    :param gaintables: input gaintables
    :param screen:
    :param height: Height (in m) of screen above telescope e.g. 3e5
    :param r0: r0 in meters
    :param type_atmosphere: 'ionosphere' or 'troposphere'
    :param scale: Multiply the screen by this factor
    :return: gridded screen image, weights image
    """
    assert isinstance(vis, BlockVisibility)

    station_locations = vis.configuration.xyz

    nant = station_locations.shape[0]
    t2r = numpy.pi / 43200.0

    newscreen = create_empty_image_like(screen)
    weights = create_empty_image_like(screen)
    nchan, ntimes, ny, nx = screen.shape

    number_no_weight = 0

    for gaintable in gaintables:
        ha_zero = numpy.average(calculate_blockvisibility_hourangles(vis))
        for iha, rows in enumerate(
                vis_timeslice_iter(vis, vis_slices=vis_slices)):
            v = create_visibility_from_rows(vis, rows)
            ha = numpy.average(
                calculate_blockvisibility_hourangles(v) -
                ha_zero).to('rad').value
            pp = find_pierce_points(station_locations,
                                    (gaintable.phasecentre.ra.rad + t2r * ha) *
                                    units.rad,
                                    gaintable.phasecentre.dec,
                                    height=height,
                                    phasecentre=vis.phasecentre)

            scr = numpy.angle(gaintable.gain[0, :, 0, 0, 0])
            wt = gaintable.weight[0, :, 0, 0, 0]
            for ant in range(nant):
                pp0 = pp[ant][0:2]
                for freq in vis.frequency:
                    scale = numpy.power(r0 / 5000.0, -5.0 / 3.0)
                    if type_atmosphere == 'troposphere':
                        # In troposphere files, the units are phase in radians.
                        screen_to_phase = scale
                    else:
                        # In the ionosphere file, the units are dTEC.
                        screen_to_phase = -scale * 8.44797245e9 / freq
                    worldloc = [pp0[0], pp0[1], ha, freq]
                    pixloc = newscreen.wcs.wcs_world2pix([worldloc],
                                                         0)[0].astype('int')
                    assert pixloc[0] >= 0
                    assert pixloc[0] < nx
                    assert pixloc[1] >= 0
                    assert pixloc[1] < ny
                    pixloc[3] = 0
                    newscreen.data[
                        pixloc[3], pixloc[2], pixloc[1],
                        pixloc[0]] += wt[ant] * scr[ant] / screen_to_phase
                    weights.data[pixloc[3], pixloc[2], pixloc[1],
                                 pixloc[0]] += wt[ant]
                    if wt[ant] == 0.0:
                        number_no_weight += 1
    if number_no_weight > 0:
        log.warning(
            "grid_gaintable_to_screen: %d pierce points are have no weight" %
            (number_no_weight))

    assert numpy.max(weights.data) > 0.0, "No points were gridded"

    newscreen.data[weights.data > 0.0] = newscreen.data[
        weights.data > 0.0] / weights.data[weights.data > 0.0]

    return newscreen, weights
Beispiel #7
0
def simulate_gaintable_from_voltage_pattern(vis,
                                            sc,
                                            vp,
                                            vis_slices=None,
                                            scale=1.0,
                                            order=3,
                                            use_radec=False,
                                            elevation_limit=15.0 * numpy.pi /
                                            180.0,
                                            **kwargs):
    """ Create gaintables from a list of components and voltagr patterns

    :param elevation_limit:
    :param use_radec:
    :param vis_slices:
    :param vis:
    :param sc: Sky components for which pierce points are needed
    :param vp: Voltage pattern in AZELGEO frame
    :param scale: Multiply the screen by this factor
    :param order: order of spline (default is 3)
    :return:
    """

    nant = vis.vis.shape[1]
    gaintables = [
        create_gaintable_from_blockvisibility(vis, **kwargs) for i in sc
    ]

    nchan, npol, ny, nx = vp.data.shape

    real_spline = [
        RectBivariateSpline(range(ny),
                            range(nx),
                            vp.data[0, pol, ...].real,
                            kx=order,
                            ky=order) for pol in range(npol)
    ]
    imag_spline = [
        RectBivariateSpline(range(ny),
                            range(nx),
                            vp.data[0, pol, ...].imag,
                            kx=order,
                            ky=order) for pol in range(npol)
    ]

    if not use_radec:
        assert isinstance(vis, BlockVisibility)
        assert vp.wcs.wcs.ctype[0] == 'AZELGEO long', vp.wcs.wcs.ctype[0]
        assert vp.wcs.wcs.ctype[1] == 'AZELGEO lati', vp.wcs.wcs.ctype[1]

        assert vis.configuration.mount[
            0] == 'azel', "Mount %s not supported yet" % vis.configuration.mount[
                0]

        number_bad = 0
        number_good = 0

        # For each hourangle, we need to calculate the location of a component
        # in AZELGEO. With that we can then look up the relevant gain from the
        # voltage pattern
        for iha, rows in enumerate(
                vis_timeslice_iter(vis, vis_slices=vis_slices)):
            v = create_visibility_from_rows(vis, rows)
            if v is not None:
                utc_time = Time([numpy.average(v.time) / 86400.0],
                                format='mjd',
                                scale='utc')
                azimuth_centre, elevation_centre = calculate_azel(
                    v.configuration.location, utc_time, vis.phasecentre)
                azimuth_centre = azimuth_centre[0].to('deg').value
                elevation_centre = elevation_centre[0].to('deg').value

                # Calculate the az el for this time
                wcs_azel = vp.wcs.sub(2).deepcopy()

                for icomp, comp in enumerate(sc):

                    if elevation_centre >= elevation_limit:

                        antgain = numpy.zeros([nant, npol], dtype='complex')
                        antwt = numpy.zeros([nant, npol])

                        # Calculate the azel of this component
                        azimuth_comp, elevation_comp = calculate_azel(
                            v.configuration.location, utc_time, comp.direction)
                        cosel = numpy.cos(elevation_comp[0]).value
                        azimuth_comp = azimuth_comp[0].to('deg').value
                        elevation_comp = elevation_comp[0].to('deg').value
                        if azimuth_comp - azimuth_centre > 180.0:
                            azimuth_centre += 360.0
                        elif azimuth_comp - azimuth_centre < -180.0:
                            azimuth_centre -= 360.0

                        try:
                            worldloc = [[
                                (azimuth_comp - azimuth_centre) * cosel,
                                elevation_comp - elevation_centre
                            ]]
                            # radius = numpy.sqrt(((azimuth_comp-azimuth_centre)*cosel)**2 +
                            #                     (elevation_comp-elevation_centre)**2)
                            pixloc = wcs_azel.wcs_world2pix(worldloc, 1)[0]
                            assert pixloc[0] > 2
                            assert pixloc[0] < nx - 3
                            assert pixloc[1] > 2
                            assert pixloc[1] < ny - 3
                            for pol in range(npol):
                                gain = real_spline[pol].ev(pixloc[1], pixloc[0]) \
                                       + 1j * imag_spline[pol].ev(pixloc[1], pixloc[0])
                                antgain[:, pol] = gain
                            for ant in range(nant):
                                ag = antgain[ant, :].reshape([2, 2])
                                ag = numpy.linalg.inv(ag)
                                antgain[ant, :] = ag.reshape([4])
                                number_good += 1
                        except (ValueError, AssertionError):
                            number_bad += 1
                            antgain[...] = 0.0
                            antwt[...] = 0.0

                        gaintables[icomp].gain[
                            iha, :, :, :] = antgain[:,
                                                    numpy.newaxis, :].reshape(
                                                        [nant, nchan, 2, 2])
                        gaintables[icomp].weight[
                            iha, :, :, :] = antwt[:, numpy.newaxis, :].reshape(
                                [nant, nchan, 2, 2])
                        gaintables[icomp].phasecentre = comp.direction
                    else:
                        gaintables[icomp].gain[...] = 1.0 + 0.0j
                        gaintables[icomp].weight[iha, :, :, :] = 0.0
                        gaintables[icomp].phasecentre = comp.direction
                        number_bad += nant

    else:
        assert isinstance(vis, BlockVisibility)
        assert vp.wcs.wcs.ctype[0] == 'RA---SIN', vp.wcs.wcs.ctype[0]
        assert vp.wcs.wcs.ctype[1] == 'DEC--SIN', vp.wcs.wcs.ctype[1]

        # The time in the Visibility is hour angle in seconds!
        number_bad = 0
        number_good = 0

        d2r = numpy.pi / 180.0
        ra_centre = vp.wcs.wcs.crval[0] * d2r
        dec_centre = vp.wcs.wcs.crval[1] * d2r

        r2d = 180.0 / numpy.pi
        s2r = numpy.pi / 43200.0
        # For each hourangle, we need to calculate the location of a component
        # in AZELGEO. With that we can then look up the relevant gain from the
        # voltage pattern
        for iha, rows in enumerate(
                vis_timeslice_iter(vis, vis_slices=vis_slices)):
            v = create_visibility_from_rows(vis, rows)
            ha = numpy.average(v.time)

            for icomp, comp in enumerate(sc):
                antgain = numpy.zeros([nant, npol], dtype='complex')
                antwt = numpy.zeros([nant, npol])
                # Calculate the location of the component in AZELGEO, then add the pointing offset
                # for each antenna
                ra_comp = comp.direction.ra.rad
                dec_comp = comp.direction.dec.rad
                for ant in range(nant):
                    wcs_azel = vp.wcs.deepcopy()
                    ra_pointing = ra_centre * r2d
                    dec_pointing = dec_centre * r2d

                    # We use WCS sensible coordinate handling by labelling the axes misleadingly
                    wcs_azel.wcs.crval[0] = ra_pointing
                    wcs_azel.wcs.crval[1] = dec_pointing
                    wcs_azel.wcs.ctype[0] = 'RA---SIN'
                    wcs_azel.wcs.ctype[1] = 'DEC--SIN'

                    for pol in range(npol):
                        worldloc = [
                            ra_comp * r2d, dec_comp * r2d, vp.wcs.wcs.crval[2],
                            vp.wcs.wcs.crval[3]
                        ]
                        try:
                            pixloc = wcs_azel.sub(2).wcs_world2pix(
                                [worldloc[:2]], 1)[0]
                            assert pixloc[0] > 2
                            assert pixloc[0] < nx - 3
                            assert pixloc[1] > 2
                            assert pixloc[1] < ny - 3
                            gain = real_spline[pol].ev(
                                pixloc[1], pixloc[0]
                            ) + 1j * imag_spline[pol].ev(pixloc[1], pixloc[0])
                            if numpy.abs(gain) > 0.0:
                                antgain[ant, pol] = 1.0 / (scale * gain)
                                antwt[ant, pol] = 1.0
                            else:
                                antgain[ant, pol] = 0.0
                                antwt[ant, pol] = 0.0
                            antwt[ant, pol] = 1.0
                            number_good += 1
                        except (ValueError, AssertionError):
                            number_bad += 1
                            antgain[ant, pol] = 1e15
                            antwt[ant, pol] = 0.0

                gaintables[icomp].gain[
                    iha, :, :, :] = antgain[:, numpy.newaxis, :].reshape(
                        [nant, nchan, 2, 2])
                gaintables[icomp].weight[
                    iha, :, :, :] = antwt[:, numpy.newaxis, :].reshape(
                        [nant, nchan, 2, 2])
                gaintables[icomp].phasecentre = comp.direction

    assert number_good > 0, "simulate_gaintable_from_voltage_pattern: No points inside the voltage pattern image"
    if number_bad > 0:
        log.warning(
            "simulate_gaintable_from_voltage_pattern: %d points are inside the voltage pattern image"
            % (number_good))
        log.warning(
            "simulate_gaintable_from_voltage_pattern: %d points are outside the voltage pattern image"
            % (number_bad))

    return gaintables
Beispiel #8
0
def simulate_gaintable_from_zernikes(vis,
                                     sc,
                                     vp_list,
                                     vp_coeffs,
                                     vis_slices=None,
                                     order=3,
                                     use_radec=False,
                                     elevation_limit=15.0 * numpy.pi / 180.0,
                                     **kwargs):
    """ Create gaintables for a set of zernikes

    :param vis:
    :param sc: Sky components for which pierce points are needed
    :param vp: List of Voltage patterns in AZELGEO frame
    :param vp_coeffs: Fractional contribution [nants, nvp]
    :param order: order of spline (default is 3)
    :return:
    """

    ntimes, nant = vis.vis.shape[0:2]
    vp_coeffs = numpy.array(vp_coeffs)
    gaintables = [
        create_gaintable_from_blockvisibility(vis, **kwargs) for i in sc
    ]

    if not use_radec:
        assert isinstance(vis, BlockVisibility)
        assert vis.configuration.mount[
            0] == 'azel', "Mount %s not supported yet" % vis.configuration.mount[
                0]

        # The time in the Visibility is hour angle in seconds!
        number_bad = 0
        number_good = 0

        # Cache the splines, one per voltage pattern
        real_splines = list()
        imag_splines = list()
        for ivp, vp in enumerate(vp_list):
            assert vp.wcs.wcs.ctype[0] == 'AZELGEO long', vp.wcs.wcs.ctype[0]
            assert vp.wcs.wcs.ctype[1] == 'AZELGEO lati', vp.wcs.wcs.ctype[1]

            nchan, npol, ny, nx = vp.data.shape
            real_splines.append(
                RectBivariateSpline(range(ny),
                                    range(nx),
                                    vp.data[0, 0, ...].real,
                                    kx=order,
                                    ky=order))
            imag_splines.append(
                RectBivariateSpline(range(ny),
                                    range(nx),
                                    vp.data[0, 0, ...].imag,
                                    kx=order,
                                    ky=order))

        latitude = vis.configuration.location.lat.rad

        r2d = 180.0 / numpy.pi
        s2r = numpy.pi / 43200.0
        # For each hourangle, we need to calculate the location of a component
        # in AZELGEO. With that we can then look up the relevant gain from the
        # voltage pattern
        for iha, rows in enumerate(
                vis_timeslice_iter(vis, vis_slices=vis_slices)):
            v = create_visibility_from_rows(vis, rows)
            ha = numpy.average(
                calculate_blockvisibility_hourangles(v).to('rad').value)

            # Calculate the az el for this hourangle and the phasecentre declination
            utc_time = Time([numpy.average(v.time) / 86400.0],
                            format='mjd',
                            scale='utc')
            azimuth_centre, elevation_centre = calculate_azel(
                v.configuration.location, utc_time, vis.phasecentre)
            azimuth_centre = azimuth_centre[0].to('deg').value
            elevation_centre = elevation_centre[0].to('deg').value

            for icomp, comp in enumerate(sc):

                if elevation_centre >= elevation_limit:

                    antgain = numpy.zeros([nant], dtype='complex')
                    # Calculate the location of the component in AZELGEO, then add the pointing offset
                    # for each antenna
                    hacomp = comp.direction.ra.rad - vis.phasecentre.ra.rad + ha
                    deccomp = comp.direction.dec.rad
                    azimuth_comp, elevation_comp = hadec_to_azel(
                        hacomp, deccomp, latitude)

                    for ant in range(nant):
                        for ivp, vp in enumerate(vp_list):
                            nchan, npol, ny, nx = vp.data.shape
                            wcs_azel = vp.wcs.deepcopy()

                            # We use WCS sensible coordinate handling by labelling the axes misleadingly
                            wcs_azel.wcs.crval[0] = azimuth_centre
                            wcs_azel.wcs.crval[1] = elevation_centre
                            wcs_azel.wcs.ctype[0] = 'RA---SIN'
                            wcs_azel.wcs.ctype[1] = 'DEC--SIN'

                            worldloc = [
                                azimuth_comp * r2d, elevation_comp * r2d,
                                vp.wcs.wcs.crval[2], vp.wcs.wcs.crval[3]
                            ]
                            try:
                                pixloc = wcs_azel.sub(2).wcs_world2pix(
                                    [worldloc[:2]], 1)[0]
                                assert pixloc[0] > 2
                                assert pixloc[0] < nx - 3
                                assert pixloc[1] > 2
                                assert pixloc[1] < ny - 3
                                gain = real_splines[ivp].ev(pixloc[1], pixloc[0]) \
                                       + 1j * imag_splines[ivp](pixloc[1], pixloc[0])
                                antgain[ant] += vp_coeffs[ant, ivp] * gain
                                number_good += 1
                            except (ValueError, AssertionError):
                                number_bad += 1
                                antgain[ant] = 1.0

                        antgain[ant] = 1.0 / antgain[ant]

                    gaintables[icomp].gain[
                        iha, :, :, :] = antgain[:, numpy.newaxis,
                                                numpy.newaxis, numpy.newaxis]
                    gaintables[icomp].phasecentre = comp.direction
            else:
                gaintables[icomp].gain[...] = 1.0 + 0.0j
                gaintables[icomp].phasecentre = comp.direction
                number_bad += nant

    else:
        assert isinstance(vis, BlockVisibility)
        number_bad = 0
        number_good = 0

        # Cache the splines, one per voltage pattern
        real_splines = list()
        imag_splines = list()
        for ivp, vp in enumerate(vp_list):
            nchan, npol, ny, nx = vp.data.shape
            real_splines.append(
                RectBivariateSpline(range(ny),
                                    range(nx),
                                    vp.data[0, 0, ...].real,
                                    kx=order,
                                    ky=order))
            imag_splines.append(
                RectBivariateSpline(range(ny),
                                    range(nx),
                                    vp.data[0, 0, ...].imag,
                                    kx=order,
                                    ky=order))

        for iha, rows in enumerate(
                vis_timeslice_iter(vis, vis_slices=vis_slices)):

            # The time in the Visibility is hour angle in seconds!
            r2d = 180.0 / numpy.pi
            # For each hourangle, we need to calculate the location of a component
            # in AZELGEO. With that we can then look up the relevant gain from the
            # voltage pattern
            v = create_visibility_from_rows(vis, rows)
            ha = numpy.average(calculate_blockvisibility_hourangles(v))

            for icomp, comp in enumerate(sc):
                antgain = numpy.zeros([nant], dtype='complex')
                antwt = numpy.zeros([nant])
                ra_comp = comp.direction.ra.rad
                dec_comp = comp.direction.dec.rad
                for ant in range(nant):
                    for ivp, vp in enumerate(vp_list):

                        assert vp.wcs.wcs.ctype[
                            0] == 'RA---SIN', vp.wcs.wcs.ctype[0]
                        assert vp.wcs.wcs.ctype[
                            1] == 'DEC--SIN', vp.wcs.wcs.ctype[1]

                        worldloc = [
                            ra_comp * r2d, dec_comp * r2d, vp.wcs.wcs.crval[2],
                            vp.wcs.wcs.crval[3]
                        ]
                        nchan, npol, ny, nx = vp.data.shape

                        try:
                            pixloc = vp.wcs.sub(2).wcs_world2pix(
                                [worldloc[:2]], 1)[0]
                            assert pixloc[0] > 2
                            assert pixloc[0] < nx - 3
                            assert pixloc[1] > 2
                            assert pixloc[1] < ny - 3
                            gain = real_splines[ivp].ev(pixloc[1], pixloc[0]) \
                                   + 1j * imag_splines[ivp](pixloc[1], pixloc[0])
                            antgain[ant] += vp_coeffs[ant, ivp] * gain
                            antwt[ant] = 1.0
                            number_good += 1
                        except (ValueError, AssertionError):
                            number_bad += 1
                            antgain[ant] = 1e15
                            antwt[ant] = 0.0

                        antgain[ant] = 1.0 / antgain[ant]

                    gaintables[icomp].gain[
                        iha, :, :, :] = antgain[:, numpy.newaxis,
                                                numpy.newaxis, numpy.newaxis]
                    gaintables[icomp].weight[
                        iha, :, :, :] = antwt[:, numpy.newaxis, numpy.newaxis,
                                              numpy.newaxis]
                    gaintables[icomp].phasecentre = comp.direction

    if number_bad > 0:
        log.warning(
            "simulate_gaintable_from_zernikes: %d points are inside the voltage pattern image"
            % (number_good))
        log.warning(
            "simulate_gaintable_from_zernikes: %d points are outside the voltage pattern image"
            % (number_bad))

    return gaintables
Beispiel #9
0
 def test_vis_timeslice_iterator_single(self):
     self.actualSetUp(times=numpy.zeros([1]))
     nchunks = vis_timeslices(self.vis, timeslice='auto')
     log.debug('Found %d chunks' % (nchunks))
     for chunk, rows in enumerate(vis_timeslice_iter(self.vis)):
         assert len(rows)