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
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
def predict_list_serial_workflow(vis_list, model_imagelist, context, vis_slices=1, facets=1, gcfcf=None, **kwargs): """Predict, iterating over both the scattered vis_list and image The visibility and image are scattered, the visibility is predicted on each part, and then the parts are assembled. :param vis_list: list of vis :param model_imagelist: Model used to determine image parameters :param vis_slices: Number of vis slices (w stack or timeslice) :param facets: Number of facets (per axis) :param context: Type of processing e.g. 2d, wstack, timeslice or facets :param gcfcg: tuple containing grid correction and convolution function :param kwargs: Parameters for functions in components :return: List of vis_lists """ assert len(vis_list) == len( model_imagelist), "Model must be the same length as the vis_list" # Predict_2d does not clear the vis so we have to do it here. vis_list = zero_list_serial_workflow(vis_list) c = imaging_context(context) vis_iter = c['vis_iterator'] predict = c['predict'] if facets % 2 == 0 or facets == 1: actual_number_facets = facets else: actual_number_facets = facets - 1 def predict_ignore_none(vis, model, g): if vis is not None: assert isinstance(vis, Visibility) or isinstance( vis, BlockVisibility), vis assert isinstance(model, Image), model return predict(vis, model, context=context, gcfcf=g, **kwargs) else: return None if gcfcf is None: gcfcf = [create_pswf_convolutionfunction(m) for m in model_imagelist] # Loop over all frequency windows if facets == 1: image_results_list = list() for ivis, sub_vis_list in enumerate(vis_list): if len(gcfcf) > 1: g = gcfcf[ivis] else: g = gcfcf[0] # Loop over sub visibility vis_predicted = copy_visibility(sub_vis_list, zero=True) for rows in vis_iter(sub_vis_list, vis_slices): row_vis = create_visibility_from_rows(sub_vis_list, rows) row_vis_predicted = predict_ignore_none( row_vis, model_imagelist[ivis], g) if row_vis_predicted is not None: vis_predicted.data['vis'][ rows, ...] = row_vis_predicted.data['vis'] image_results_list.append(vis_predicted) return image_results_list else: image_results_list = list() for ivis, sub_vis_list in enumerate(vis_list): # Create the graph to divide an image into facets. This is by reference. facet_lists = image_scatter_facets(model_imagelist[ivis], facets=facets) facet_vis_lists = list() sub_vis_lists = visibility_scatter(sub_vis_list, vis_iter, vis_slices) # Loop over sub visibility for sub_sub_vis_list in sub_vis_lists: facet_vis_results = list() # Loop over facets for facet_list in facet_lists: # Predict visibility for this subvisibility from this facet facet_vis_list = predict_ignore_none( sub_sub_vis_list, facet_list, None) facet_vis_results.append(facet_vis_list) # Sum the current sub-visibility over all facets facet_vis_lists.append(sum_predict_results(facet_vis_results)) # Sum all sub-visibilties image_results_list.append( visibility_gather(facet_vis_lists, sub_vis_list, vis_iter)) return image_results_list
def invert_list_serial_workflow(vis_list, template_model_imagelist, dopsf=False, normalize=True, facets=1, vis_slices=1, context='2d', gcfcf=None, **kwargs): """ Sum results from invert, iterating over the scattered image and vis_list :param vis_list: list of vis :param template_model_imagelist: list of template models :param dopsf: Make the PSF instead of the dirty image :param facets: Number of facets :param normalize: Normalize by sumwt :param vis_slices: Number of slices :param context: Imaging context :param gcfcg: tuple containing grid correction and convolution function :param kwargs: Parameters for functions in components :return: List of (image, sumwt) tuples, one per vis in vis_list For example:: model_list = [create_image_from_visibility (v, npixel=npixel, cellsize=cellsize, polarisation_frame=pol_frame) for v in vis_list] dirty_list = invert_list_serial_workflow(vis_list, template_model_imagelist=model_list, context='wstack', vis_slices=51) dirty, sumwt = dirty_list[centre] """ if not isinstance(template_model_imagelist, collections.abc.Iterable): template_model_imagelist = [template_model_imagelist] c = imaging_context(context) vis_iter = c['vis_iterator'] invert = c['invert'] if facets % 2 == 0 or facets == 1: actual_number_facets = facets else: actual_number_facets = max(1, (facets - 1)) def gather_image_iteration_results(results, template_model): result = create_empty_image_like(template_model) i = 0 sumwt = numpy.zeros([template_model.nchan, template_model.npol]) for dpatch in image_scatter_facets(result, facets=facets): assert i < len( results), "Too few results in gather_image_iteration_results" if results[i] is not None: assert len(results[i]) == 2, results[i] dpatch.data[...] = results[i][0].data[...] sumwt += results[i][1] i += 1 return result, sumwt def invert_ignore_none(vis, model, gg): if vis is not None: return invert(vis, model, context=context, dopsf=dopsf, normalize=normalize, gcfcf=gg, **kwargs) else: return create_empty_image_like(model), numpy.zeros( [model.nchan, model.npol]) # If we are doing facets, we need to create the gcf for each image if gcfcf is None and facets == 1: gcfcf = [create_pswf_convolutionfunction(template_model_imagelist[0])] # Loop over all vis_lists independently results_vislist = list() if facets == 1: for ivis, sub_vis_list in enumerate(vis_list): if len(gcfcf) > 1: g = gcfcf[ivis] else: g = gcfcf[0] # Iterate within each vis_list result_image = create_empty_image_like( template_model_imagelist[ivis]) result_sumwt = numpy.zeros([ template_model_imagelist[ivis].nchan, template_model_imagelist[ivis].npol ]) for rows in vis_iter(sub_vis_list, vis_slices): row_vis = create_visibility_from_rows(sub_vis_list, rows) result = invert_ignore_none(row_vis, template_model_imagelist[ivis], g) if result is not None: result_image.data += result[1][:, :, numpy.newaxis, numpy. newaxis] * result[0].data result_sumwt += result[1] result_image = normalize_sumwt(result_image, result_sumwt) results_vislist.append((result_image, result_sumwt)) else: for ivis, sub_vis_list in enumerate(vis_list): # Create the graph to divide an image into facets. This is by reference. facet_lists = image_scatter_facets(template_model_imagelist[ivis], facets=facets) # Create the graph to divide the visibility into slices. This is by copy. sub_sub_vis_lists = visibility_scatter(sub_vis_list, vis_iter, vis_slices=vis_slices) # Iterate within each vis_list vis_results = list() for sub_sub_vis_list in sub_sub_vis_lists: facet_vis_results = list() for facet_list in facet_lists: facet_vis_results.append( invert_ignore_none(sub_sub_vis_list, facet_list, None)) vis_results.append( gather_image_iteration_results( facet_vis_results, template_model_imagelist[ivis])) results_vislist.append(sum_invert_results(vis_results)) return results_vislist
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
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
def test_apply_voltage_pattern_image_pointsource(self): self.createVis(rmax=1e3) telescope = 'MID_FEKO_B2' vpol = PolarisationFrame("linear") self.times = numpy.linspace(-4, +4, 8) * numpy.pi / 12.0 bvis = create_blockvisibility(self.config, self.times, self.frequency, channel_bandwidth=self.channel_bandwidth, phasecentre=self.phasecentre, weight=1.0, polarisation_frame=vpol, zerow=True) cellsize = advise_wide_field(bvis)['cellsize'] pbmodel = create_image_from_visibility( bvis, cellsize=self.cellsize, npixel=self.npixel, override_cellsize=False, polarisation_frame=PolarisationFrame("stokesIQUV")) vpbeam = create_vp(pbmodel, telescope=telescope, use_local=False) vpbeam.wcs.wcs.ctype[0] = 'RA---SIN' vpbeam.wcs.wcs.ctype[1] = 'DEC--SIN' vpbeam.wcs.wcs.crval[0] = pbmodel.wcs.wcs.crval[0] vpbeam.wcs.wcs.crval[1] = pbmodel.wcs.wcs.crval[1] s3_components = create_test_skycomponents_from_s3( flux_limit=0.1, phasecentre=self.phasecentre, frequency=self.frequency, polarisation_frame=PolarisationFrame('stokesI'), radius=1.5 * numpy.pi / 180.0) for comp in s3_components: comp.polarisation_frame = PolarisationFrame('stokesIQUV') comp.flux = numpy.array([[comp.flux[0, 0], 0.0, 0.0, 0.0]]) s3_components = filter_skycomponents_by_flux(s3_components, 0.0, 10.0) from rascil.processing_components.image import show_image import matplotlib.pyplot as plt plt.clf() show_image(vpbeam, components=s3_components) plt.show(block=False) vpcomp = apply_voltage_pattern_to_skycomponent(s3_components, vpbeam) bvis.data['vis'][...] = 0.0 + 0.0j bvis = dft_skycomponent_visibility(bvis, vpcomp) rec_comp = idft_visibility_skycomponent(bvis, vpcomp)[0] stokes_comp = list() for comp in rec_comp: stokes_comp.append( convert_pol_frame(comp.flux[0], PolarisationFrame("linear"), PolarisationFrame("stokesIQUV"))) stokesI = numpy.abs( numpy.array([comp_flux[0] for comp_flux in stokes_comp]).real) stokesQ = numpy.abs( numpy.array([comp_flux[1] for comp_flux in stokes_comp]).real) stokesU = numpy.abs( numpy.array([comp_flux[2] for comp_flux in stokes_comp]).real) stokesV = numpy.abs( numpy.array([comp_flux[3] for comp_flux in stokes_comp]).real) plt.clf() plt.loglog(stokesI, stokesQ, '.', label='Q') plt.loglog(stokesI, stokesU, '.', label='U') plt.loglog(stokesI, stokesV, '.', label='V') plt.xlabel("Stokes Flux I (Jy)") plt.ylabel("Flux (Jy)") plt.legend() plt.savefig('%s/test_primary_beams_pol_rsexecute_stokes_errors.png' % self.dir) plt.show(block=False) split_times = False if split_times: bvis_list = list() for rows in vis_timeslice_iter(bvis, vis_slices=8): bvis_list.append(create_visibility_from_rows(bvis, rows)) else: bvis_list = [bvis] bvis_list = rsexecute.scatter(bvis_list) model_list = \ [rsexecute.execute(create_image_from_visibility, nout=1)(bv, cellsize=cellsize, npixel=4096, phasecentre=self.phasecentre, override_cellsize=False, polarisation_frame=PolarisationFrame("stokesIQUV")) for bv in bvis_list] model_list = rsexecute.persist(model_list) bvis_list = weight_list_rsexecute_workflow(bvis_list, model_list) continuum_imaging_list = \ continuum_imaging_list_rsexecute_workflow(bvis_list, model_list, context='2d', algorithm='hogbom', facets=1, niter=1000, fractional_threshold=0.1, threshold=1e-4, nmajor=5, gain=0.1, deconvolve_facets=4, deconvolve_overlap=32, deconvolve_taper='tukey', psf_support=64, restore_facets=4, psfwidth=1.0) clean, residual, restored = rsexecute.compute(continuum_imaging_list, sync=True) centre = 0 if self.persist: export_image_to_fits( clean[centre], '%s/test_primary_beams_pol_rsexecute_clean.fits' % self.dir) export_image_to_fits( residual[centre][0], '%s/test_primary_beams_pol_rsexecute_residual.fits' % self.dir) export_image_to_fits( restored[centre], '%s/test_primary_beams_pol_rsexecute_restored.fits' % self.dir) plt.clf() show_image(restored[centre]) plt.show(block=False) qa = qa_image(restored[centre]) assert numpy.abs(qa.data['max'] - 0.9953017707113947) < 1.0e-7, str(qa) assert numpy.abs(qa.data['min'] + 0.0036396480874570846) < 1.0e-7, str(qa)