|
| 1 | +#! /usr/bin/env python3 |
| 2 | + |
| 3 | +__author__ = 'Paul Hancock' |
| 4 | +__date__ = '2022/07/08' |
| 5 | + |
| 6 | +import argparse |
| 7 | +import os |
| 8 | +import sys |
| 9 | + |
| 10 | +import AegeanTools.wcs_helpers as wcsh |
| 11 | +import numpy as np |
| 12 | +from astropy.convolution import Gaussian2DKernel |
| 13 | +from astropy.convolution.kernels import Gaussian2DKernel |
| 14 | +from astropy.io import fits |
| 15 | +from radio_beam import Beam, Beams |
| 16 | +from scipy import ndimage, optimize, signal |
| 17 | + |
| 18 | +SIGMA_TO_FWHM = np.sqrt(8*np.log(2)) |
| 19 | + |
| 20 | + |
| 21 | +def get_common_beam(fnames): |
| 22 | + """ |
| 23 | + Read headers from a list of fits files and determine the smallest |
| 24 | + common beam to which they could be convolved |
| 25 | +
|
| 26 | + parameters |
| 27 | + ---------- |
| 28 | + fnames : list([str, ...]) |
| 29 | + List of filenames |
| 30 | +
|
| 31 | + return |
| 32 | + ------ |
| 33 | + cb : `radio_beam.Beam` |
| 34 | + The smallest common beam |
| 35 | + """ |
| 36 | + headers = [fits.getheader(f) for f in fnames] |
| 37 | + beams = [Beam.from_fits_header(h) for h in headers] |
| 38 | + beams=Beams(beams=beams) |
| 39 | + cb = beams.common_beam() |
| 40 | + return cb |
| 41 | + |
| 42 | + |
| 43 | +def beam_to_kernel(beam, |
| 44 | + ra, dec, |
| 45 | + wcshelper): |
| 46 | + """ |
| 47 | + Convert a beam object into a Gaussian kernel. |
| 48 | + Since the sky->pix mapping may not be affine the sky location is required for this conversion. |
| 49 | +
|
| 50 | + parameters |
| 51 | + ---------- |
| 52 | + beam : `radio_beam.Beam` |
| 53 | + The beam object |
| 54 | + |
| 55 | + ra, dec : float, float |
| 56 | + The sky location |
| 57 | + |
| 58 | + wcshelper: `AegeanTools.wcs_helpers.wcshelper` |
| 59 | + An augmented WCS object. |
| 60 | + |
| 61 | + returns |
| 62 | + ------- |
| 63 | + kern : `astropy.convolution.Gaussian2DKernel` |
| 64 | + The kernel |
| 65 | + """ |
| 66 | + _, _, a, b, pa = wcshelper.sky2pix_ellipse((ra,dec), beam.major.value, beam.minor.value, beam.pa.value) |
| 67 | + kern = Gaussian2DKernel(a/SIGMA_TO_FWHM, y_stddev=b/SIGMA_TO_FWHM, theta=pa) |
| 68 | + return kern |
| 69 | + |
| 70 | + |
| 71 | +def find_gaussian_kernel(image_kern, target_kern): |
| 72 | + """ |
| 73 | + Find the kernel B which when convolved with image_kern will produce target_kern |
| 74 | + This is solving A*B = C for B, where A and C are known and * is the 2d convolution. |
| 75 | +
|
| 76 | + parameters |
| 77 | + ---------- |
| 78 | + image_kern, target_kern : `astropy.convolution.Gaussian2DKernel` |
| 79 | + The kernels A and C |
| 80 | + |
| 81 | + return |
| 82 | + ------ |
| 83 | + B : `astropy.convolution.Gaussian2DKernel` |
| 84 | + """ |
| 85 | + def diff(kern): |
| 86 | + return np.sum(np.abs(signal.convolve2d(kern.array, image_kern.array, mode='same')- target_kern.array)) |
| 87 | + |
| 88 | + |
| 89 | + def func(x): |
| 90 | + a,b,pa = x |
| 91 | + B = Gaussian2DKernel(a/SIGMA_TO_FWHM, y_stddev=b/SIGMA_TO_FWHM, theta=pa, |
| 92 | + x_size=image_kern.shape[0], y_size=image_kern.shape[1]) |
| 93 | + return diff(B) |
| 94 | + |
| 95 | + res = optimize.minimize(func, x0=[1,1,0]) |
| 96 | + |
| 97 | + B_fit = Gaussian2DKernel(res.x[0]/SIGMA_TO_FWHM, y_stddev=res.x[1]/SIGMA_TO_FWHM, theta=res.x[2], |
| 98 | + x_size=image_kern.shape[0], y_size=image_kern.shape[1]) |
| 99 | + return B_fit |
| 100 | + |
| 101 | + |
| 102 | +def convolve_to(psf, infile, outfile): |
| 103 | + """ |
| 104 | + Convolved an image to the target psf. |
| 105 | + Save the resulting image. |
| 106 | +
|
| 107 | + parameters |
| 108 | + ---------- |
| 109 | + psf : `radio_beam.Beam` |
| 110 | + The target psf (must be same/larger than that of the input file) |
| 111 | +
|
| 112 | + infile : str |
| 113 | + Path to read input image |
| 114 | +
|
| 115 | + outfile : str |
| 116 | + Path to save output |
| 117 | + """ |
| 118 | + hdu = fits.open(infile) |
| 119 | + |
| 120 | + # determine the current psf |
| 121 | + wcshelper = wcsh.WCSHelper.from_header(hdu[0].header) |
| 122 | + beam = Beam.from_fits_header(hdu[0].header) |
| 123 | + im_kern = beam_to_kernel(beam, hdu[0].header['CRVAL1'], hdu[0].header['CRVAL2'], wcshelper) |
| 124 | + t_kern = beam_to_kernel(psf, hdu[0].header['CRVAL1'], hdu[0].header['CRVAL2'], wcshelper) |
| 125 | + |
| 126 | + # find the kernel that will get us from im_kern -> t_kern |
| 127 | + fit_kern = find_gaussian_kernel(im_kern, t_kern) |
| 128 | + |
| 129 | + # convolve and update the image data |
| 130 | + data = signal.convolve2d(hdu[0].data, fit_kern, mode='same') |
| 131 | + hdu[0].data = data |
| 132 | + |
| 133 | + # update the header info |
| 134 | + header = hdu[0].header |
| 135 | + header['BMAJ'] = psf.major.value |
| 136 | + header['BMIN'] = psf.minor.value |
| 137 | + header['BPA'] = psf.pa.value |
| 138 | + |
| 139 | + # save the file |
| 140 | + hdu.writeto(outfile, overwrite=True) |
| 141 | + print("Wrote {0}".format(outfile)) |
| 142 | + return |
| 143 | + |
| 144 | + |
| 145 | +if __name__ == "__main__": |
| 146 | + parser = argparse.ArgumentParser() |
| 147 | + group1 = parser.add_argument_group("Convolve images") |
| 148 | + group1.add_argument("--in", dest='files', type=str, default=None, nargs='+', |
| 149 | + help="List of input images") |
| 150 | + group1.add_argument("--suffix", dest='suffix', type=str, default="convolved", |
| 151 | + help="Output suffix for filenames.") |
| 152 | + results = parser.parse_args() |
| 153 | + |
| 154 | + if results.files is None or len(results.files) <2: |
| 155 | + parser.print_help() |
| 156 | + sys.exit(0) |
| 157 | + |
| 158 | + psf = get_common_beam(results.files) |
| 159 | + for f in results.files: |
| 160 | + base, ext = os.path.splitext(f) |
| 161 | + outfile = "{0}_{1}{2}".format(base, results.suffix, ext) |
| 162 | + convolve_to(psf, f, outfile) |
| 163 | + |
0 commit comments