def recipe_phangs_broad_mask( template_mask, outfile=None, list_of_masks=[], grow_xy=None, grow_v=None, return_spectral_cube=True, overwrite=False, recipe='anyscale', fraction_of_scales=0.25, ): """Task to create the PHANGS-style "broad" masks from the combination of a set of other masks. Optionally also grow the mask at the end. Parameters: ----------- template_mask : string or SpectralCube The original mask that holds the target WCS. The other masks will be reprojected onto this one. This mask is included in the final output. Keywords: --------- outfile : string Filename where the mask will be written. The mask is also returned. list_of_masks : list List of masks or spectral cubes. These will be reprojected onto the template mask and then combined via logical or to form the final mask. grow_xy : int Number of spatial dilations of the mask. grow_v : int Number of spectralwise dilations of the mask. recipe : str If 'anyscale', the broad mask include pixels that are included in any of the contributing masks. If 'somescales', this only include pixels that are in fraction_of_scales masks of the origina list. fraction_of_scales : float Set to the fraction of scales that is the minimum for inclusion in the broadmask. """ # &%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&% # Error checking and work out inputs # &%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&% # TBD error checking, dimensions, types, etc. if type(template_mask) is SpectralCube: mask = template_mask elif type(template_mask) == str: mask = SpectralCube.read(template_mask) else: logger.error("Input mask must be a SpectralCube object or a filename.") return (None) mask.allow_huge_operations = True # &%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&% # Loop over other masks # &%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&% for other_mask in list_of_masks: if type(other_mask) is SpectralCube: other_mask = other_mask elif type(template_mask) == type("hello"): other_mask = SpectralCube.read(other_mask) else: logger.error( "Input masks must be SpectralCube objects or filenames.") return (None) other_mask.allow_huge_operations = True mask = join_masks(mask, other_mask, operation='sum', order='fast_nearest_neighbor') if recipe is 'anyscale': mask_values = mask.filled_data[:].value > 0 mask = SpectralCube(mask_values * 1.0, wcs=mask.wcs, header=mask.header, meta={ 'BUNIT': ' ', 'BTYPE': 'Mask' }) mask.allow_huge_operations = True if recipe is 'somescales': mask_values = mask.filled_data[:].value > (fraction_of_scales * len(list_of_masks)) mask = SpectralCube(mask_values * 1.0, wcs=mask.wcs, header=mask.header, meta={ 'BUNIT': ' ', 'BTYPE': 'Mask' }) mask.allow_huge_operations = True # &%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&% # Grow the mask # &%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&% if (grow_xy is not None) or (grow_v is not None): mask_values = mask.filled_data[:].value if grow_xy is not None: mask_values = grow_mask(mask_values, iters_xy=grow_xy) if grow_v is not None: mask_values = grow_mask(mask_values, iters_v=grow_v) mask = SpectralCube(mask_values * 1.0, wcs=mask.wcs, header=mask.header, meta={ 'BUNIT': ' ', 'BTYPE': 'Mask' }) mask.allow_huge_operations = True # &%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&% # Write to disk and return # &%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&% # Write to disk, if desired if outfile is not None: header = mask.header header['DATAMAX'] = 1 header['DATAMIN'] = 0 header[ 'COMMENT'] = 'Produced with PHANGS-ALMA pipeline version ' + version if tableversion: header[ 'COMMENT'] = 'Galaxy properties from PHANGS sample table version ' + tableversion hdu = fits.PrimaryHDU(np.array(mask.filled_data[:], dtype=np.uint8), header=header) hdu.writeto(outfile, overwrite=overwrite) if return_spectral_cube: return (mask) else: return (mask.filled_data[:].value)
def calc_noise_in_cube(cube, masking_scheme='simple', mask=None, spatial_average_npix=None, spatial_average_nbeam=5.0, spectral_average_nchan=5, verbose=False): """ Estimate rms noise in a (continuum-subtracted) spectral cube. Parameters ---------- cube : SpectralCube object Input spectral cube (needs to be continuum-subtracted) masking_scheme : {'simple', 'user'}, optional Scheme for flagging signal in the cube. 'simple' means to flag all values above 3*rms or below -3*rms (default scheme); 'user' means to use the user-specified mask (i.e., `mask`). mask : `np.ndarray` object, optional User-specified signal mask (this parameter is ignored if `masking_scheme` is not 'user') spatial_average_npix : int, optional Size of the spatial averaging box, in terms of pixel number If not None, `spatial_average_nbeam` will be ingored. (Default: None) spatial_average_nbeam : float, optional Size of the spatial averaging box, in the unit of beam FWHM (Default: 5.0) spectral_average_nchan : int, optional Size of the spectral averaging box, in terms of channel number (Default: 5) verbose : bool, optional Whether to print the detailed processing information in terminal Default is to not print. Returns ------- rmscube : SpectralCube object Spectral cube containing the rms noise at each ppv location """ from scipy.ndimage import generic_filter from astropy.stats import mad_std if masking_scheme not in ['simple', 'user']: raise ValueError("'masking_scheme' should be specified as" "either 'simple' or 'user'") elif masking_scheme == 'user' and mask is None: raise ValueError("'masking_scheme' set to 'user', yet " "no user-specified mask found") # extract negative values (only needed if masking_scheme='simple') if masking_scheme == 'simple': if verbose: print("Extracting negative values...") negmask = cube < (0 * cube.unit) negdata = cube.with_mask(negmask).filled_data[:].value negdata = np.stack([negdata, -1 * negdata], axis=-1) else: negdata = None # find rms noise as a function of channel if verbose: print("Estimating rms noise as a function of channel...") if masking_scheme == 'user': mask_v = mask elif masking_scheme == 'simple': rms_v = mad_std(negdata, axis=(1, 2, 3), ignore_nan=True) uplim_v = (3 * rms_v * cube.unit).reshape(-1, 1, 1) lolim_v = (-3 * rms_v * cube.unit).reshape(-1, 1, 1) mask_v = (((cube - uplim_v) < (0 * cube.unit)) & ((cube - lolim_v) > (0 * cube.unit))) rms_v = cube.with_mask(mask_v).mad_std(axis=(1, 2)).quantity.value rms_v = generic_filter(rms_v, np.nanmedian, mode='constant', cval=np.nan, size=spectral_average_nchan) # find rms noise as a function of sightline if verbose: print("Estimating rms noise as a function of sightline...") if masking_scheme == 'user': mask_s = mask elif masking_scheme == 'simple': rms_s = mad_std(negdata, axis=(0, 3), ignore_nan=True) uplim_s = 3 * rms_s * cube.unit lolim_s = -3 * rms_s * cube.unit mask_s = (((cube - uplim_s) < (0 * cube.unit)) & ((cube - lolim_s) > (0 * cube.unit))) rms_s = cube.with_mask(mask_s).mad_std(axis=0).quantity.value if spatial_average_npix is None: beamFWHM_pix = (cube.beam.major.to(u.deg).value / np.abs(cube.wcs.celestial.wcs.cdelt.min())) beamFWHM_pix = np.max([beamFWHM_pix, 3.]) spatial_average_npix = int(spatial_average_nbeam * beamFWHM_pix) rms_s = generic_filter(rms_s, np.nanmedian, mode='constant', cval=np.nan, size=spatial_average_npix) # create rms noise cube from the tensor product of rms_v and rms_s if verbose: print("Creating rms noise cube (direct tensor product)...") rmscube = SpectralCube(np.einsum('i,jk', rms_v, rms_s), wcs=cube.wcs, header=cube.header.copy(strip=True)) rmscube.allow_huge_operations = cube.allow_huge_operations # correct the normalization of the rms cube if masking_scheme == 'user': mask_n = mask elif masking_scheme == 'simple': rms_n = mad_std(negdata, ignore_nan=True) uplim_n = 3 * rms_n * cube.unit lolim_n = -3 * rms_n * cube.unit mask_n = (((cube - uplim_n) < (0 * cube.unit)) & ((cube - lolim_n) > (0 * cube.unit))) rms_n = cube.with_mask(mask_n).mad_std().value rmscube /= rms_n # apply NaN mask rmscube = rmscube.with_mask(cube.mask.include()) # check unit if rmscube.unit != cube.unit: rmscube = rmscube * (cube.unit / rmscube.unit) return rmscube
def cleansplit(filename, galaxy=None, Vwindow=650 * u.km / u.s, Vgalaxy=300 * u.km / u.s, blorder=3, HanningLoops=0, maskfile=None, circleMask=True, minRadius=1.3 * u.arcmin, edgeMask=False, weightCut=0.2, spectralSetup=None, CatalogFile=None, spatialSmooth=1.0): """ Takes a raw DEGAS cube and produces individual cubes for each spectral line. Paramters --------- filename : str The file to split. Keywords -------- galaxy : Galaxy object Currently unused Vwindow : astropy.Quantity Width of the window in velocity units Vgalaxy : astropy.Quantity Line of sight velocity of the galaxy centre blorder : int Baseline order HanningLoops : int Number of times to smooth and resample the data edgeMask : bool Determine whether to apply an edgeMask weightCut : float Minimum weight value to include in the data spatialSmooth : float Factor to increase the (linear) beam size by in a convolution. spectralSetup : str String to determine how we set up the spectrum 'hcn_hcop' -- split based on HCN/HCO+ setup '13co_c18o' -- split based on 13CO/C18O setup '12co' -- don't split; assume single line """ Cube = SpectralCube.read(filename) Cube.allow_huge_operations = True if CatalogFile is None: CatalogFile = get_pkg_data_filename('./data/dense_survey.cat', package='degas') Catalog = Table.read(CatalogFile, format='ascii') # Find which galaxy in our catalog corresponds to the object we # are mapping if galaxy is None: RABound, DecBound = Cube.world_extrema match = np.zeros_like(Catalog, dtype=np.bool) for index, row in enumerate(Catalog): galcoord = SkyCoord(row['RA'], row['DEC'], unit=(u.hourangle, u.deg)) if (galcoord.ra < RABound[1] and galcoord.ra > RABound[0] and galcoord.dec < DecBound[1] and galcoord.dec > DecBound[0]): match[index] = True MatchRow = Catalog[match] galcoord = SkyCoord(MatchRow['RA'], MatchRow['DEC'], unit=(u.hourangle, u.deg)) Galaxy = MatchRow['NAME'].data[0] print("Catalog Match with " + Galaxy) V0 = MatchRow['CATVEL'].data[0] * u.km / u.s elif type(galaxy) is str: match = np.zeros_like(Catalog, dtype=np.bool) for index, row in enumerate(Catalog): if galaxy in row['NAME']: match[index] = True MatchRow = Catalog[match] galcoord = SkyCoord(MatchRow['RA'], MatchRow['DEC'], unit=(u.hourangle, u.deg)) Galaxy = MatchRow['NAME'].data[0] print("Catalog Match with " + Galaxy) V0 = MatchRow['CATVEL'].data[0] * u.km / u.s # Check spectral setups. Use the max frequencies present to # determine which spectral setup we used if not specifed. if spectralSetup is None: if (Cube.spectral_axis.max() > 105 * u.GHz and Cube.spectral_axis.max() < 113 * u.GHz): warnings.warn("assuming 13CO/C18O spectral setup") spectralSetup = '13CO_C18O' filestr = '13co_c18o' if (Cube.spectral_axis.max() > 82 * u.GHz and Cube.spectral_axis.max() < 90 * u.GHz): warnings.warn("assuming HCN/HCO+ spectral setup") spectralSetup = 'HCN_HCO+' filestr = 'hcn_hcop' if (Cube.spectral_axis.max() > 113 * u.GHz): warnings.warn("assuming 12CO spectral setup") spectralSetup = '12CO' filestr = '12co' if spectralSetup == '13CO_C18O': CEighteenO = Cube.with_spectral_unit(u.km / u.s, velocity_convention='radio', rest_value=109.78217 * u.GHz) ThirteenCO = Cube.with_spectral_unit(u.km / u.s, velocity_convention='radio', rest_value=110.20135 * u.GHz) CubeList = (CEighteenO, ThirteenCO) LineList = ('C18O', '13CO') elif spectralSetup == 'HCN_HCO+': HCN = Cube.with_spectral_unit(u.km / u.s, velocity_convention='radio', rest_value=88.631847 * u.GHz) HCOp = Cube.with_spectral_unit(u.km / u.s, velocity_convention='radio', rest_value=89.188518 * u.GHz) CubeList = (HCN, HCOp) LineList = ('HCN', 'HCOp') elif spectralSetup == '12CO': TwelveCO = Cube.with_spectral_unit(u.km / u.s, velocity_convention='radio', rest_value=115.27120180 * u.GHz) CubeList = (TwelveCO, ) LineList = ('12CO', ) for ThisCube, ThisLine in zip(CubeList, LineList): if circleMask: x0, y0, _ = ThisCube.wcs.wcs_world2pix(galcoord.ra, galcoord.dec, 0, 0) FoVFile = get_pkg_data_filename('./data/field_of_view.csv', package='degas') if FoVFile: fov_table = Table.read(FoVFile) if np.any(fov_table['Galaxy'] == Galaxy): minRadius = (fov_table[fov_table['Galaxy'] == Galaxy] ['FoV_arcsec']) * u.arcsec ThisCube = circletrim(ThisCube, filename.replace('.fits', '_wts.fits'), x0, y0, weightCut=weightCut, minRadius=minRadius) if edgeMask: ThisCube = edgetrim(ThisCube, filename.replace('.fits', '_wts.fits'), weightCut=weightCut) # Trim each cube to the specified velocity range ThisCube = ThisCube.spectral_slab(V0 - Vwindow, V0 + Vwindow) ThisCube.write(Galaxy + '_' + ThisLine + '.fits', overwrite=True) StartChan = ThisCube.closest_spectral_channel(V0 - Vgalaxy) EndChan = ThisCube.closest_spectral_channel(V0 + Vgalaxy) if maskfile is not None: maskLookup = buildMaskLookup(maskfile) shp = ThisCube.shape TmpCube = ThisCube.with_spectral_unit(u.Hz) spaxis = TmpCube.spectral_axis spaxis = spaxis.value data = ThisCube.filled_data[:].value for y in np.arange(shp[1]): for x in np.arange(shp[2]): spectrum = data[:, y, x] if np.any(np.isnan(spectrum)): continue coords = ThisCube.world[:, y, x] mask = maskLookup(coords[2].value, coords[1].value, spaxis) spectrum = robustBaseline(spectrum, blorder=blorder, baselineIndex=~mask) data[:, y, x] = spectrum ThisCube = SpectralCube(data * ThisCube.unit, ThisCube.wcs, header=ThisCube.header, meta={'BUNIT': ThisCube.header['BUNIT']}) ThisCube.write(Galaxy + '_' + ThisLine + '_rebase{0}.fits'.format(blorder), overwrite=True) else: gbtpipe.Baseline.rebaseline(Galaxy + '_' + ThisLine + '.fits', baselineRegion=[ slice(0, StartChan, 1), slice(EndChan, ThisCube.shape[0], 1) ], blorder=blorder) ThisCube = SpectralCube.read(Galaxy + '_' + ThisLine + '_rebase{0}'.format(blorder) + '.fits') ThisCube.allow_huge_operations = True # Smooth Kern = Kernel1D(array=np.array([0.5, 1.0, 0.5])) for i in range(HanningLoops): ThisCube = ThisCube.spectral_smooth(Kern) ThisCube = ThisCube[::2, :, :] # Spatial Smooth if spatialSmooth > 1.0: newBeam = Beam(major=ThisCube.beam.major * spatialSmooth, minor=ThisCube.beam.minor * spatialSmooth) ThisCube = ThisCube.convolve_to(newBeam) smoothstr = '_smooth{0}'.format(spatialSmooth) else: smoothstr = '' # apply eta_mb ## see equations 2 and 3 in GBT memo 302. ## GBT forward efficiency is 0.99 ~= 1.0. ## assumes that rest freq ~ observing freq, which should be okay for our sources. freq = ThisCube.header['RESTFRQ'] * u.Hz (eta_a, eta_mb) = calc_etamb(freq) ThisCube = ThisCube / eta_mb # Final Writeout finalFile = Galaxy + '_' + ThisLine + '_rebase{0}'.format(blorder) + smoothstr + \ '_hanning{0}.fits'.format(HanningLoops) ThisCube.write(finalFile, overwrite=True) # Get the headers right via brute force fits manipulation. finalCube = fits.open(finalFile) finalCube[0].header['ETAMB'] = eta_mb finalCube[0].header['BUNIT'] = ( 'K', 'TMB') # indicate that units are now TMB via a comment finalCube[0].writeto(finalFile, overwrite=True)
def smooth_cube( incube=None, outfile=None, angular_resolution=None, linear_resolution=None, distance=None, velocity_resolution=None, nan_treatment='interpolate', # can also be 'fill' tol=None, make_coverage_cube=False, collapse_coverage=False, coveragefile=None, coverage2dfile=None, dtype=np.float32, overwrite=True ): """ Smooth an input cube to coarser angular or spectral resolution. This lightly wraps spectral cube and some of the error checking is left to that. tol is a fraction. When the target beam is within tol of the original beam, we just copy. Optionally, also calculate a coverage footprint in which original (finite) cube coverage starts at 1.0 and the output cube shows the fraction of finite pixels. """ # &%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&% # Error checking # &%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&% # Require a valid cube or map input if type(incube) is SpectralCube: cube = incube elif type(incube) == type("hello"): cube = SpectralCube.read(incube) else: logger.error("Input must be a SpectralCube object or a filename.") # Allow huge operations. If the speed or segfaults become a huge # problem, we will adjust our strategy here. cube.allow_huge_operations = True # Check that only one target scale is set if (angular_resolution is not None) and (linear_resolution is not None): logger.error('Only one of angular_resolution or ', 'linear_resolution can be set') return(None) # Work out the target angular resolution if angular_resolution is not None: if type(angular_resolution) is str: angular_resolution = u.Quantity(angular_resolution) if linear_resolution is not None: if distance is None: logger.error('Convolution to linear resolution requires a distance.') return(None) if type(distance) is str: distance = u.Quantity(distance) if type(linear_resolution) is str: linear_resolution = u.Quantity(linear_resolution) angular_resolution = (linear_resolution / distance * u.rad).to(u.arcsec) dist_mpc_val = float(distance.to(u.pc).value) / 1e6 cube._header.append(('DIST_MPC',dist_mpc_val,'Used in convolution')) if tol is None: tol = 0.0 # &%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&% # Convolution to coarser beam # &%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&% if angular_resolution is not None: logger.info("... convolving from beam: "+str(cube.beam)) target_beam = Beam(major=angular_resolution, minor=angular_resolution, pa=0 * u.deg) logger.info("... convolving to beam: "+str(target_beam)) new_major = float(target_beam.major.to(u.arcsec).value) old_major = float(cube.beam.major.to(u.arcsec).value) delta = (new_major-old_major)/old_major logger.info("... fractional change: "+str(delta)) if make_coverage_cube: coverage = SpectralCube(np.isfinite(cube.unmasked_data[:])*1.0, wcs=cube.wcs, header=cube.header, meta={'BUNIT': ' ', 'BTYPE': 'Coverage'}) coverage = coverage.with_mask(LazyMask(np.isfinite,cube=coverage)) # Allow huge operations. If the speed or segfaults become a huge # problem, we will adjust our strategy here. coverage.allow_huge_operations = True if delta > tol: logger.info("... proceeding with convolution.") cube = cube.convolve_to(target_beam, nan_treatment=nan_treatment) if make_coverage_cube: coverage = coverage.convolve_to(target_beam, nan_treatment=nan_treatment) if np.abs(delta) < tol: logger.info("... current resolution meets tolerance.") if delta < -1.0*tol: logger.info("... resolution cannot be matched. Returning") return(None) # &%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&% # Spectral convolution # &%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&% # This is only a boxcar smooth right now and does not downsample # or update the header. if velocity_resolution is not None: if type(velocity_resolution) is str: velocity_resolution = u.Quantity(velocity_resolution) dv = scdr.channel_width(cube) nChan = (velocity_resolution / dv).to(u.dimensionless_unscaled).value if nChan > 1: cube = cube.spectral_smooth(Box1DKernel(nChan)) if make_coverage_cube: coverage = coverage.spectral_smooth(Box1DKernel(nChan)) # &%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&% # Write or return as requested # &%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&%&% if outfile is not None: # cube.write(outfile, overwrite=overwrite) hdu = fits.PrimaryHDU(np.array(cube.filled_data[:], dtype=dtype), header=cube.header) hdu.writeto(outfile, overwrite=overwrite) if make_coverage_cube: if coveragefile is not None: hdu = fits.PrimaryHDU(np.array(coverage.filled_data[:], dtype=dtype), header=coverage.header) hdu.writeto(coveragefile, overwrite=overwrite) if collapse_coverage: if coveragefile and not coverage2dfile: coverage2dfile = coveragefile.replace('.fits','2d.fits') coverage_collapser(coverage, coverage2dfile=coverage2dfile, overwrite=overwrite) # coverage.write(coveragefile, overwrite=overwrite) return(cube)
def maskmoment(img_fits, gain_fits=None, rms_fits=None, mask_fits=None, outdir='.', outname=None, snr_hi=4, snr_lo=2, minbeam=1, snr_hi_minch=1, snr_lo_minch=1, min_tot_ch=2, nguard=[0, 0], edgech=5, fwhm=None, vsm=None, vsm_type='gauss', mom1_minch=2, mom2_minch=2, altoutput=False, output_snr_cube=False, output_snr_peak=False, output_snrsm_cube=False, output_2d_mask=False, to_kelvin=True, huge_operations=True, perpixel=False): """ Produce FITS images of moment maps using a dilated masking approach. Parameters ---------- img_fits : FITS file name, required The image cube, this should be in units of K, Jy/beam, or equivalent. gain_fits : FITS file name, optional The gain cube, e.g. pb cube from CASA. This should have a value between 0 and 1, with 0 near the edges and 1 near the center of the image, and have the same dimensions as the image cube. NOTE: The gain cube is ignored if a noise cube (rms_fits) is given. rms_fits : FITS file name, optional The noise cube, providing an estimate of the rms noise at each pixel. This should have the same dimensions and units as the image cube. NOTE: If rms_fits is not given, a noise cube is generated from the image cube, after removing any gain variation using the gain cube. mask_fits : FITS file name, optional External mask cube to use. This cube should have 1's for valid pixels and 0's for excluded pixels. If this is provided then the mask generation is skipped and the program goes straight to calculating the moments. outdir : string, optional Directory to write the output files. Default: Write to the directory where img_fits resides. NOTE: Currently this directory is assumed to exist. outname : string, optional Basename for output files. For instance, outname='foo' produces files 'foo.mom0.fits.gz', etc. Default: Based on root name of img_fits. snr_hi : float, optional The high significance sigma threshold from which to begin mask dilation. Default: 4 snr_lo : float, optional The low significance sigma threshold at which to end mask dilation. Default: 2 snr_hi_minch : int, optional High significance mask is required to span at least this many channels at all pixels. Default: 1 snr_lo_minch : int, optional Low significance mask is required to span at least this many channels at all pixels. Default: 1 min_tot_ch : int, optional Each contiguous mask region must span at least this many channels (but it's not necessary that every pixel in the region span this many channels). Default: 2 minbeam : float, optional Minimum velocity-integrated area of a mask region in units of the beam size. Default: 1 nguard : tuple of two ints, optional Expand the final mask by nguard[0] pixels in the sky directions and nguard[1] channels in velocity. Currently these values must be equal if both are non-zero. If nguard[0] = 0 then no expansion is done in sky coordinates. If nguard[1] = 0 then no expansion is done in velocity. Default: [0,0] edgech : int, optional Number of channels at each end of vel axis to use for rms estimation. Default: 5 fwhm : float or :class:`~astropy.units.Quantity`, optional Spatial resolution to smooth to before generating the dilated mask. If value is not an astropy quantity, assumed to be given in arcsec. Default: No spatial smoothing is applied. vsm : float or :class:`~astropy.units.Quantity`, optional Full width of the spectral smoothing kernel (or FWHM for gaussian). If given as astropy quantity, should be given in velocity units. If not given as astropy quantity, interpreted as number of channels. Default: No spectral smoothing is applied. vsm_type : string, optional What type of spectral smoothing to employ. Currently three options: (1) 'boxcar' - 1D boxcar smoothing, vsm rounded to integer # of chans. (2) 'gauss' - 1D gaussian smoothing, vsm is the convolving gaussian FWHM. (3) 'gaussfinal' - 1D gaussian smoothing, vsm is the gaussian FWHM after convolution, assuming FWHM before convolution is 1 channel. Default: 'gauss' mom1_minch : int, optional Minimum number of unmasked channels needed to calculate moment-1. Default: 2 mom2_minch : int, optional Minimum number of unmasked channels needed to calculate moment-2. Default: 2 perpixel : boolean, optional Whether to calculate the rms per XY pixel instead of over whole image. Set to True if you know there is a sensitivity variation across the image but you don't have a gain cube - requires rms_fits and gain_fits unset. Default: False output_snr_cube : boolean, optional Output the cube in SNR units in addition to the moment maps. Default: False output_snr_peak : boolean, optional Output the peak SNR image in addition to the moment maps. Default: False output_snrsm_cube : boolean, optional Output the smoothed cube in SNR units in addition to the moment maps. Default: False output_2d_mask : boolean, optional Output the projected 2-D mask as well as the newly generated mask. The projected mask at a given pixel is valid for all channels as long as the parent mask is valid for any channel. Default: False to_kelvin : boolean, optional Output the moment maps in K units if the cube is in Jy/beam units. Default: True altoutput : boolean, optional Also output moment maps from a "direct" calculation instead of the moment method in spectral_cube. Mainly used for debugging. Default: False huge_operations : boolean, optional Allow huge operations in spectral_cube. Default: True """ np.seterr(divide='ignore', invalid='ignore') # # --- DETERMINE OUTPUT FILE NAMES # if outdir != '': pth = outdir + '/' else: img_dir = os.path.dirname(img_fits) if img_dir == '': pth = './' else: pth = img_dir + '/' if outname is not None: basename = outname else: basename = os.path.basename(img_fits) if basename.endswith('.fits.gz'): basename = os.path.splitext(os.path.splitext(basename)[0])[0] else: basename = os.path.splitext(basename)[0] print('\nOutput basename is:', basename) # # --- READ INPUT FILES, OUTPUT NOISE CUBE IF NEWLY GENERATED # image_cube = SpectralCube.read(img_fits) hd3d = image_cube.header hd2d = hd3d.copy() hd2d['WCSAXES'] = 2 for key in ['CRVAL3', 'CTYPE3', 'CRPIX3', 'CDELT3', 'CUNIT3']: del hd2d[key] has_jypbeam = all(x in (hd3d['bunit']).upper() for x in ['JY', 'B']) print('Image cube ' + img_fits + ':\n', image_cube) if rms_fits is not None: rms_cube = SpectralCube.read(rms_fits) if rms_cube.unit != image_cube.unit: raise RuntimeError('Noise cube units', rms_cube.unit, 'differ from image units', image_cube.unit) if rms_cube.shape[0] == 1: rmsarray = np.broadcast_to(rms_cube, image_cube.shape) unit = rms_cube.unit rms_cube = SpectralCube(data=rmsarray * unit, header=hd3d, wcs=wcs.WCS(hd3d)) print('Noise cube ' + rms_fits + ':\n', rms_cube) else: if gain_fits is not None: gain_cube = SpectralCube.read(gain_fits) print('Gain cube ' + gain_fits + ':\n', gain_cube) rms_cube = makenoise(image_cube, gain_cube, edge=edgech) else: rms_cube = makenoise(image_cube, edge=edgech, perpixel=perpixel) print('Noise cube:\n', rms_cube) hd3d['datamin'] = np.nanmin(rms_cube._data[0]) hd3d['datamax'] = np.nanmax(rms_cube._data[0]) fits.writeto(pth + basename + '.ecube.fits.gz', rms_cube._data.astype(np.float32), hd3d, overwrite=True) print('Wrote', pth + basename + '.ecube.fits.gz') # # --- GENERATE AND OUTPUT SNR CUBE, PEAK SNR IMAGE # if mask_fits is not None: dilatedmask = fits.getdata(mask_fits) else: if huge_operations: image_cube.allow_huge_operations = True rms_cube.allow_huge_operations = True snr_cube = image_cube / rms_cube print('SNR cube:\n', snr_cube) if output_snr_cube: hd3d['datamin'] = snr_cube.min().value hd3d['datamax'] = snr_cube.max().value hd3d['bunit'] = ' ' fits.writeto(pth + basename + '.snrcube.fits.gz', snr_cube._data.astype(np.float32), hd3d, overwrite=True) print('Wrote', pth + basename + '.snrcube.fits.gz') if output_snr_peak: snr_peak = snr_cube.max(axis=0) hd2d['datamin'] = np.nanmin(snr_peak.value) hd2d['datamax'] = np.nanmax(snr_peak.value) hd2d['bunit'] = ' ' fits.writeto(pth + basename + '.snrpk.fits.gz', snr_peak.astype(np.float32), hd2d, overwrite=True) print('Wrote', pth + basename + '.snrpk.fits.gz') # # --- GENERATE AND OUTPUT DILATED MASK # if fwhm is not None or vsm is not None: sm_snrcube = smcube(snr_cube, fwhm=fwhm, vsm=vsm, vsm_type=vsm_type, edgech=edgech) print('Smoothed SNR cube:\n', sm_snrcube) if output_snrsm_cube: hd3d['datamin'] = sm_snrcube.min().value hd3d['datamax'] = sm_snrcube.max().value hd3d['bunit'] = ' ' fits.writeto(pth + basename + '.snrsmcube.fits.gz', sm_snrcube._data.astype(np.float32), hd3d, overwrite=True) print('Wrote', pth + basename + '.snrsmcube.fits.gz') dilcube = sm_snrcube else: dilcube = snr_cube dilatedmask = dilmsk(dilcube, header=hd3d, snr_hi=snr_hi, snr_lo=snr_lo, snr_hi_minch=snr_hi_minch, snr_lo_minch=snr_lo_minch, min_tot_ch=min_tot_ch, nguard=nguard, minbeam=minbeam) hd3d['datamin'] = 0 hd3d['datamax'] = 1 hd3d['bunit'] = ' ' fits.writeto(pth + basename + '.mask.fits.gz', dilatedmask.astype(np.float32), hd3d, overwrite=True) print('Wrote', pth + basename + '.mask.fits.gz') if output_2d_mask: summed_msk = np.broadcast_to(np.sum(dilatedmask, axis=0), image_cube.shape) proj_msk = np.minimum(summed_msk, 1) fits.writeto(pth + basename + '.mask2d.fits.gz', proj_msk.astype(np.float32), hd3d, overwrite=True) print('Wrote', pth + basename + '.mask2d.fits.gz') # # --- CALCULATE FLUXES (IN ORIGINAL UNITS) # if mask_fits is None and output_2d_mask: fluxtab = findflux(image_cube, rms_cube, dilatedmask, proj_msk) else: fluxtab = findflux(image_cube, rms_cube, dilatedmask) fluxtab.write(pth + basename + '.flux.csv', delimiter=',', format='ascii.ecsv', overwrite=True) print('Wrote', pth + basename + '.flux.csv') # # --- GENERATE AND OUTPUT MOMENT MAPS # nchanimg = np.sum(dilatedmask, axis=0) print('Units of cube are', image_cube.unit) if to_kelvin and has_jypbeam: if hasattr(image_cube, 'beam'): print('Beam info:', image_cube.beam) else: print('WARNING: Beam info is missing') image_cube = image_cube.to(u.K) rms_cube = rms_cube.to(u.K) dil_mskcub = image_cube.with_mask(dilatedmask > 0) dil_mskcub_mom0 = dil_mskcub.moment(order=0).to(image_cube.unit * u.km / u.s) if hasattr(dil_mskcub_mom0, 'unit'): print('Units of mom0 map are', dil_mskcub_mom0.unit) writemom(dil_mskcub_mom0, type='mom0', filename=pth + basename, hdr=hd2d) # --- Moment 1: mean velocity must be in range of cube dil_mskcub_mom1 = dil_mskcub.moment(order=1).to(u.km / u.s) vmin = dil_mskcub.spectral_extrema[0] vmax = dil_mskcub.spectral_extrema[1] dil_mskcub_mom1[dil_mskcub_mom1 < vmin] = np.nan dil_mskcub_mom1[dil_mskcub_mom1 > vmax] = np.nan dil_mskcub_mom1[nchanimg < mom1_minch] = np.nan writemom(dil_mskcub_mom1, type='mom1', filename=pth + basename, hdr=hd2d) # --- Moment 2: require at least 2 unmasked channels at each pixel dil_mskcub_mom2 = dil_mskcub.linewidth_sigma().to(u.km / u.s) dil_mskcub_mom2[nchanimg < mom2_minch] = np.nan writemom(dil_mskcub_mom2, type='mom2', filename=pth + basename, hdr=hd2d) # # --- CALCULATE ERRORS IN MOMENTS # altmom, errmom = calc_moments(image_cube, rms_cube, dilatedmask) errmom0 = errmom[0] * abs(hd3d['CDELT3']) / 1000 * dil_mskcub_mom0.unit errmom0[nchanimg == 0] = np.nan writemom(errmom0, type='emom0', filename=pth + basename, hdr=hd2d) errmom1 = errmom[1] * u.km / u.s errmom1[dil_mskcub_mom1 < vmin] = np.nan errmom1[dil_mskcub_mom1 > vmax] = np.nan errmom1[nchanimg < mom1_minch] = np.nan writemom(errmom1, type='emom1', filename=pth + basename, hdr=hd2d) errmom2 = errmom[2] * u.km / u.s errmom2[nchanimg < mom2_minch] = np.nan writemom(errmom2, type='emom2', filename=pth + basename, hdr=hd2d) if altoutput: altmom0 = altmom[0] * abs(hd3d['CDELT3']) / 1000 * dil_mskcub_mom0.unit altmom0[nchanimg == 0] = np.nan altmom1 = altmom[1] * u.km / u.s altmom1[altmom1 < vmin] = np.nan altmom1[altmom1 > vmax] = np.nan altmom1[nchanimg < mom1_minch] = np.nan altmom2 = altmom[2] * u.km / u.s altmom2[nchanimg < mom2_minch] = np.nan writemom(altmom0, type='amom0', filename=pth + basename, hdr=hd2d) writemom(altmom1, type='amom1', filename=pth + basename, hdr=hd2d) writemom(altmom2, type='amom2', filename=pth + basename, hdr=hd2d) return
def convolve_cube(incube, newbeam, mode='datacube', res_tol=0.0, min_coverage=0.8, nan_treatment='fill', boundary='fill', fill_value=0., append_raw=False, verbose=False, suppress_error=False): """ Convolve a spectral cube or an rms noise cube to a specified beam. 'datacube' mode: This is basically a wrapper around `~spectral_cube.SpectralCube.convolve_to()`, but it treats NaN values / edge effect in a more careful way (see the description of keyword 'min_coverage' below). 'noisecube' mode: It handles rms noise cubes in a way that it correctly predicts the rms noise for the corresponding data cube convolved to the same specified beam. Parameters ---------- incube : FITS HDU object or SpectralCube object Input spectral cube newbeam : radio_beam.Beam object Target beam to convolve to mode : {'datacube', 'noisecube'}, optional Whether the input cube is a data cube or an rms noise cube. In the former case, a direct convolution is performed; in the latter case, the convolution attempts to mimic the error propagation process to the specified lower resolution. (Default: 'datacube') res_tol : float, optional Tolerance on the difference between input/output resolution By default, a convolution is performed on the input cube when its native resolution is different from (sharper than) the target resolution. Use this keyword to specify a tolerance on resolution, within which no convolution will be performed. For example, ``res_tol=0.1`` will allow a 10% tolerance. min_coverage : float or None, optional This keyword specifies a minimum beam covering fraction of valid pixels for convolution (Default: 0.8). Locations with a beam covering fraction less than this value will be overwritten to "NaN" in the convolved cube. If the user would rather use the ``preserve_nan`` mode in `astropy.convolution.convolve_fft`, set this keyword to None. nan_treatment: {'interpolate', 'fill'}, optional To be passed to `astropy.convolution.convolve_fft`. (Default: 'fill') boundary: {'fill', 'wrap'}, optional To be passed to `astropy.convolution.convolve_fft`. (Default: 'fill') fill_value : float, optional To be passed to `astropy.convolution.convolve_fft`. (Default: 0) append_raw : bool, optional Whether to append the raw convolved cube and weight cube. Default is not to append. verbose : bool, optional Whether to print the detailed processing log in terminal. Default is to not print. suppress_error : bool, optional Whether to suppress the error message when convolution is unsuccessful. Default is to not suppress. Returns ------- outcube : FITS HDU objects or SpectralCube objects Convolved spectral cubes (when append_raw=False), or a 3-tuple including a masked verson, an unmaked version, and a coverage fraction cube (when append_raw=True). The output will be the same type of objects as the input. """ if isinstance(incube, SpectralCube): cube = incube elif isinstance(incube, (fits.PrimaryHDU, fits.ImageHDU)): if 'BUNIT' in incube.header: unit = u.Unit(incube.header['BUNIT'], parse_strict='warn') if isinstance(unit, u.UnrecognizedUnit): unit = u.dimensionless_unscaled else: unit = u.dimensionless_unscaled cube = SpectralCube(data=incube.data * unit, wcs=WCS(incube.header), header=incube.header, allow_huge_operations=True).with_mask( np.isfinite(incube.data)) else: raise ValueError("`incube` needs to be either a SpectralCube object " "or a FITS HDU object") if (res_tol > 0) and (newbeam.major != newbeam.minor): raise ValueError("Cannot handle a non-zero resolution torelance " "when the target beam is not round") if min_coverage is None: # Skip coverage check and preserve NaN values. # This uses the default 'preserve_nan' scheme # implemented in 'astropy.convolution.convolve_fft' convolve_func = partial(convolve_fft, fill_value=fill_value, nan_treatment=nan_treatment, boundary=boundary, preserve_nan=True, allow_huge=True) else: # Do coverage check to determine the mask on the output convolve_func = partial(convolve_fft, fill_value=fill_value, nan_treatment=nan_treatment, boundary=boundary, allow_huge=True) convolve_func_w = partial(convolve_fft, fill_value=0., boundary='fill', allow_huge=True) tol = newbeam.major * np.array([1 - res_tol, 1 + res_tol]) if ((tol[0] < cube.beam.major < tol[1]) and (tol[0] < cube.beam.minor < tol[1])): if verbose: print("Native resolution within tolerance - " "copying original cube...") my_append_raw = False convcube = wtcube = None newcube = cube.unmasked_copy().with_mask(cube.mask.include()) else: if verbose: print("Deconvolving beam...") try: beamdiff = newbeam.deconvolve(cube.beam) except ValueError as err: if suppress_error: if verbose: print("Unsuccessful beam deconvolution: " "{}\nOld: {}\nNew: {}" "".format(err, cube.beam, newbeam)) print("Exiting...") return else: raise ValueError("Unsuccessful beam deconvolution: " "{}\nOld: {}\nNew: {}" "".format(err, cube.beam, newbeam)) if verbose: print("Convolving cube...") if mode == 'datacube': # do convolution convcube = cube.convolve_to(newbeam, convolve=convolve_func) if min_coverage is not None: my_append_raw = True wtcube = SpectralCube(cube.mask.include().astype('float'), cube.wcs, beam=cube.beam).with_mask( np.ones(cube.shape).astype('?')) wtcube.allow_huge_operations = (cube.allow_huge_operations) wtcube = wtcube.convolve_to(newbeam, convolve=convolve_func_w) # divide the raw convolved image by the weight image # to correct for filling fraction newcube = convcube / wtcube.unmasked_data[:] # mask all pixels w/ weight smaller than min_coverage threshold = min_coverage * u.dimensionless_unscaled newcube = newcube.with_mask(wtcube >= threshold) else: my_append_raw = False newcube = convcube wtcube = None elif mode == 'noisecube': # Empirically derive a noise cube at the lower resolution # Step 1: square the high resolution noise cube cubesq = cube**2 # Step 2: convolve the squared noise cube with a kernel # that is sqrt(2) times narrower than the one # used for data cube convolution (this is because # the Gaussian weight needs to be squared in # error propagation) beamdiff_small = Beam(major=beamdiff.major / np.sqrt(2), minor=beamdiff.minor / np.sqrt(2), pa=beamdiff.pa) newbeam_small = cube.beam.convolve(beamdiff_small) convcubesq = cubesq.convolve_to(newbeam_small, convolve=convolve_func) del cubesq # release memory if min_coverage is not None: my_append_raw = True wtcube = SpectralCube(cube.mask.include().astype('float'), cube.wcs, beam=cube.beam).with_mask( np.ones(cube.shape).astype('?')) wtcube.allow_huge_operations = (cube.allow_huge_operations) # divide the raw convolved cube by the weight cube # to correct for filling fraction wtcube_d = wtcube.convolve_to(newbeam_small, convolve=convolve_func_w) newcubesq = convcubesq / wtcube_d.unmasked_data[:] del convcubesq, wtcube_d # release memory # mask all pixels w/ weight smaller than min_coverage # (here I force the masking of the noise cube to be # consistent with that of the data cube) wtcube = wtcube.convolve_to(newbeam, convolve=convolve_func_w) threshold = min_coverage * u.dimensionless_unscaled newcubesq = newcubesq.with_mask(wtcube >= threshold) else: my_append_raw = False newcubesq = convcubesq wtcube = None # Step 3: find the sqrt of the convolved noise cube convcube = np.sqrt(convcubesq) del convcubesq # release memory newcube = np.sqrt(newcubesq) del newcubesq # release memory # Step 4: apply a multiplicative factor, which accounts # for the decrease in rms noise due to averaging convcube *= np.sqrt(cube.beam.sr / newbeam.sr).to('').value newcube *= np.sqrt(cube.beam.sr / newbeam.sr).to('').value else: raise ValueError("Invalid `mode` value: {}".format(mode)) if isinstance(incube, SpectralCube): if append_raw and my_append_raw: return newcube, convcube, wtcube else: return newcube elif isinstance(incube, (fits.PrimaryHDU, fits.ImageHDU)): if append_raw and my_append_raw: return newcube.hdu, convcube.hdu, wtcube.hdu else: return newcube.hdu
def common_cube(filename, target_rms=None, target_resolution=None, datadir='', outdir='', outname='common.fits', cube=None, rmsname=None, rmsmap=None, distance=None, overwrite=True, return_cube=False, huge=True): """ Convolves and adds noise to a data cube to arrive at a target linear resolution and noise level Parameters: filename : str File name of spectral data cube to process outname : str File name for output cube target_rms : astropy.Quantity Target RMS level in units of cube brightness target_resolution : astropy.Quantity Target resolution expressed as either length (linear resolution) or angle (angular resolution) datadir : str Directory for input data outdir : str Directory for output data rmsname : str Name for input RMS data cube rmsmap : np.ndarray numpy array of noise values at positions corresponding to the spectral cube. distance : astropy.Quantity Distance to target; required if a linear target resolution is requested overwrite : bool Overwrite output products if they already exist? return_cube : bool Return a SpectralCube? """ if cube is None: cube = SpectralCube.read(os.path.join(datadir, filename)) cube.allow_huge_operations = huge orig_beam_area = cube.beam.sr target_beam = None if target_resolution is not None: if target_resolution.unit.is_equivalent(u.pc): if distance is None: warnings.warn( "A distance with units must be provided to reach a target linear resolution" ) raise ValueError target_beam_size = (target_resolution / distance).to( u.dimensionless_unscaled) * u.rad elif target_resolution.is_equivalent(u.rad): target_beam_size = target_resolution else: warnings.warn( "target_resolution must be either linear distance or angle") raise ValueError target_beam = Beam(major=target_beam_size, minor=target_beam_size, pa=0 * u.deg) if target_beam > cube.beam: cube = cube.convolve_to(target_beam) else: warnings.warn('Target linear resolution unreachable.') if target_rms is not None: if target_beam is None: target_beam = cube.beam noisevals = np.random.randn(*cube.shape) null_beam = Beam(major=1e-7 * u.deg, minor=1e-7 * u.deg, pa=0 * u.deg) noisevals[~np.isfinite(cube)] = np.nan noisecube = SpectralCube(noisevals, cube.wcs, header=cube.header, beam=null_beam, meta={'BUNIT': 'K'}) noisecube = noisecube.with_mask(np.isfinite(cube)) noisecube.allow_huge_operations = huge area_ratio = (orig_beam_area / target_beam.sr).to( u.dimensionless_unscaled).value noisecube = noisecube.convolve_to(target_beam) if (rmsname is not None) and (rmsmap is None): rmsmap = fits.getdata(os.path.join(datadir, rmsname)) * u.K if rmsmap is None: warnings.warn( "One of rmsmap or rmsname must be provided for noise homogenization" ) raise FileNotFoundError rmsamplitude = rmsmap * area_ratio**0.5 if np.all(rmsamplitude > target_rms): warnings.warn("All noise values larger than target value") else: addamplitude = (target_rms**2 - rmsamplitude**2) addamplitude[addamplitude < 0] = 0.0 addamplitude.shape = addamplitude.shape noisevals = noisecube.filled_data[:] noisevals /= np.nanstd(noisevals) noisevals *= (addamplitude**0.5) cube.allow_huge_operations = huge cube += noisevals if type(outname) is str: cube.write(outdir + outname, overwrite=overwrite) if return_cube: return (cube) else: return (True)