Image¶
-
class
ehtim.image.
Image
(image, psize, ra, dec, pa=0.0, polrep='stokes', pol_prim=None, rf=230000000000.0, pulse=<function trianglePulse2D>, source='SgrA', mjd=51544, time=0.0)[source]¶ A polarimetric image (in units of Jy/pixel).
-
pulse
¶ The function convolved with the pixel values for continuous image.
- Type
function
-
add_const_pol
(mag, angle, cmag=0, csign=1)[source]¶ Return an with constant fractional linear and circular polarization
-
add_crescent
(flux, Rp, Rn, a, b, x=0, y=0, pol=None)[source]¶ Add a crescent to an image; see Kamruddin & Dexter (2013).
- Parameters
flux (float) – the total flux contained in the crescent in Jy
Rp (float) – the larger radius in radians
Rn (float) – the smaller radius in radians
a (float) – the relative x offset of smaller disk in radians
b (float) – the relative y offset of smaller disk in radians
x (float) – the center x coordinate of the larger disk in radians
y (float) – the center y coordinate of the larger disk in radians
pol (str) – the polarization to add the flux to. None defaults to pol_prim.
- Returns
output image add_gaus
- Return type
(Image)
-
add_qu
(qimage, uimage)[source]¶ Add Stokes Q and U images. self.polrep must be ‘stokes’
- Parameters
qimage (numpy.array) – The 2D Stokes Q values in Jy/pixel array
uimage (numpy.array) – The 2D Stokes U values in Jy/pixel array
Returns:
-
add_random_pol
(mag, corr, cmag=0.0, ccorr=0.0, seed=0)[source]¶ Return an image random linear and circular polarizations with certain correlation lengths
-
add_ring_m1
(I0, I1, r0, phi, sigma, x=0, y=0, pol=None)[source]¶ Add a ring to an image with an m=1 mode
- Parameters
I0 (float) –
I1 (float) –
r0 (float) – the radius
phi (float) – angle of m1 mode
sigma (float) – the blurring size
x (float) – the center x coordinate of the larger disk in radians
y (float) – the center y coordinate of the larger disk in radians
pol (str) – the polarization to add the flux to. None defaults to pol_prim.
- Returns
output image add_gaus
- Return type
(Image)
-
add_tophat
(flux, radius, pol=None)[source]¶ Add centered tophat flux to the Stokes I image inside a given radius.
-
add_v
(vimage)[source]¶ Add Stokes V image. self.polrep must be ‘stokes’
- Parameters
vimage (numpy.array) – The 2D Stokes Q values in Jy/pixel array
-
add_zblterm
(obs, uv_min, zblval=None, new_fov=False, gauss_sz=False, gauss_sz_factor=0.75, debias=True)[source]¶ Add a large Gaussian term to account for missing flux in the zero baseline.
- Parameters
obs – an Obsdata object to determine min non-zero baseline and 0-bl flux
uv_min (float) – The cutoff in Glambada used to determine what is a 0-bl
new_fov (rad) – The size of the padded image once the Gaussian is added (if False it will be set to 3 x the gaussian fwhm)
gauss_sz (rad) – The size of the Gaussian added to add flux to the 0-bl. (if False it is computed from the min non-zero baseline)
gauss_sz_factor (float) – The fraction of the min non-zero baseline used to caluclate the Gaussian FWHM.
debias (bool) – True if you use debiased amplitudes to caluclate the 0-bl flux in Jy
- Returns
a padded image with a large Gaussian component
- Return type
(Image)
-
align_images
(im_list, pol=None, shift=True, final_fov=False, scale='lin', gamma=0.5, dynamic_range=[1000.0])[source]¶ Align all the images in im_list to the current image (self) Aligns all images by comparison of the primary pol image.
- Parameters
im_list (list) – list of images to align to the current image
shift (list) – list of manual image shifts, otherwise use the shift from maximum cross-correlation
pol (str) – which polarization image to compare. Default is self.pol_prim
final_fov (float) – fov of the comparison image (rad). If False it is the largestinput image fov
scale (str) – compare images in ‘log’,’lin’,or ‘gamma’ scale
gamma (float) – exponent for gamma scale comparison
dynamic_range (float) – dynamic range for log and gamma scale comparisons
- Returns
(im_list_shift, shifts, im0_pad)
- Return type
(tuple)
-
blur_circ
(fwhm_i, fwhm_pol=0, filttype='gauss')[source]¶ Apply a circular gaussian filter to the image, with FWHM in radians.
-
blur_gauss
(beamparams, frac=1.0, frac_pol=0)[source]¶ Blur image with a Gaussian beam w/ beamparams [fwhm_max, fwhm_min, theta] in radians.
-
blur_mf
(freqs, fwhm, fit_order=1, filttype='gauss')[source]¶ Blur image correctly across multiple frequencies WARNING: does not currently do polarization correctly!
-
property
bvec
¶ Return the B mode image vector
-
center
(pol=None)[source]¶ Center the image based on the coordinates of the centroid(). A non-integer shift is used, which wraps the image when rotating.
- Parameters
pol (str) – The polarization for which to find the image centroid
- Returns
centroid positions (x0,y0) in radians
- Return type
(np.array)
-
centroid
(pol=None)[source]¶ Compute the location of the image centroid (corresponding to the polarization pol)
- Parameters
pol (str) – The polarization for which to find the image centroid
- Returns
centroid positions (x0,y0) in radians
- Return type
(np.array)
-
property
chivec
¶ Return the fractional polarization angle for each pixel
-
circ_polfrac
()[source]¶ Return the total fractional circular polarized flux
Args:
- Returns
image fractional circular polarized flux
- Return type
(float)
-
compare_images
(im_compare, pol=None, psize=None, target_fov=None, blur_frac=0.0, beamparams=[1.0, 1.0, 1.0], metric=['nxcorr', 'nrmse', 'rssd'], blursmall=False, shift=True)[source]¶ Compare to another image by computing normalized cross correlation, normalized root mean squared error, or square root of the sum of squared differences. Returns metrics only for the primary polarization imvec!
- Parameters
im_compare (Image) – the image to compare to
pol (str) – which polarization image to compare. Default is self.pol_prim
psize (float) – pixel size of comparison image (rad). If None it is the smallest of the input image pizel sizes
target_fov (float) – fov of the comparison image (rad). If None it is twice the largest fov of the input images
beamparams (list) – the nominal Gaussian beam parameters [fovx, fovy, position angle]
blur_frac (float) – fractional beam to blur each image to before comparison
metric (list) – a list of fidelity metrics from [‘nxcorr’,’nrmse’,’rssd’]
blursmall (bool) – True to blur the unpadded image rather than the large image.
shift (int) – manual image shift, otherwise use shift from maximum cross-correlation
- Returns
[errormetric, im1_pad, im2_shift]
- Return type
(tuple)
-
contour
(contour_levels=[0.1, 0.25, 0.5, 0.75], contour_cfun=None, color='w', legend=True, show_im=True, cfun='afmhot', scale='lin', interp='gaussian', gamma=0.5, dynamic_range=1000.0, plotp=False, nvec=20, pcut=0.01, mcut=0.1, label_type='ticks', has_title=True, has_cbar=True, cbar_lims=(), cbar_unit='Jy', 'pixel', contour_im=False, power=0, beamcolor='w', export_pdf='', show=True, beamparams=None, cbar_orientation='vertical', scale_lw=1, beam_lw=1, cbar_fontsize=12, axis=None, scale_fontsize=12)[source]¶ Display the image in a contour plot.
- Parameters
contour_levels (arr) – the fractional contour levels relative to the max flux plotted
contour_cfun (pyplot colormap function) – the function used to get the RGB colors
legend (bool) – True to show a legend that says what each contour line corresponds to
cfun (str) – matplotlib.pyplot color function
scale (str) – image scaling in [‘log’,’gamma’,’lin’]
interp (str) – image interpolation ‘gauss’ or ‘lin’
gamma (float) – index for gamma scaling
dynamic_range (float) – dynamic range for log and gamma scaling
plotp (bool) – True to plot linear polarimetic image
nvec (int) – number of polarimetric vectors to plot
pcut (float) – minimum stokes P value for displaying polarimetric vectors as fraction of maximum Stokes I pixel
mcut (float) – minimum fractional polarization for plotting vectors
label_type (string) – specifies the type of axes labeling: ‘ticks’, ‘scale’, ‘none’
has_title (bool) – True if you want a title on the plot
has_cbar (bool) – True if you want a colorbar on the plot
cbar_lims (tuple) – specify the lower and upper limit of the colorbar
cbar_unit (tuple of strings) – the unit of each pixel for the colorbar: ‘Jy’, ‘m-Jy’, ‘$mu$Jy’
export_pdf (str) – path to exported PDF with plot
show (bool) – Display the plot if true
show_im (bool) – Display the image with the contour plot if True
- Returns
figure object with image
- Return type
(matplotlib.figure.Figure)
-
copy
()[source]¶ Return a copy of the Image object.
Args:
- Returns
copy of the Image.
- Return type
(Image)
-
copy_pol_images
(old_image)[source]¶ Copy polarization images from old_image over to self.
- Parameters
old_image (Image) – image object to copy from
-
display
(pol=None, cfun=False, interp='gaussian', scale='lin', gamma=0.5, dynamic_range=1000.0, plotp=False, plot_stokes=False, nvec=20, vec_cfun=None, vec_cbar_lims=(), scut=0, pcut=0.1, mcut=0.01, scale_ticks=False, log_offset=False, label_type='ticks', has_title=True, alpha=1, has_cbar=True, only_cbar=False, cbar_lims=(), cbar_unit='Jy', 'pixel', export_pdf='', pdf_pad_inches=0.0, show=True, beamparams=None, cbar_orientation='vertical', scinot=False, scale_lw=1, beam_lw=1, cbar_fontsize=12, axis=None, scale_fontsize=12, power=0, beamcolor='w', beampos='right', scalecolor='w', dpi=500)[source]¶ Display the image.
- Parameters
pol (str) – which polarization image to plot. Default is self.pol_prim pol=’spec’ will plot spectral index pol=’curv’ will plot spectral curvature
cfun (str) – matplotlib.pyplot color function. False changes with ‘pol’, but is ‘afmhot’ for most
interp (str) – image interpolation ‘gauss’ or ‘lin’
scale (str) – image scaling in [‘log’,’gamma’,’lin’]
gamma (float) – index for gamma scaling
dynamic_range (float) – dynamic range for log and gamma scaling
plotp (bool) – True to plot linear polarimetic image
plot_stokes (bool) – True to plot stokes subplots along with plotp
nvec (int) – number of polarimetric vectors to plot
vec_cfun (str) – color function for vectors colored by lin pol frac
vec_cbar_lims (tuple) – lower and upper limit of the fractional polarization colormap
scut (float) – minimum fractional stokes I value for displaying spectral index
pcut (float) – minimum fractional stokes I value for displaying polarimetric vectors
mcut (float) – minimum fractional polarization value for displaying vectors
scale_ticks (bool) – if True, scale polarization ticks by linear polarization magnitude
label_type (string) – specifies the type of axes labeling: ‘ticks’, ‘scale’, ‘none’
has_title (bool) – True if you want a title on the plot
has_cbar (bool) – True if you want a colorbar on the plot
cbar_lims (tuple) – specify the lower and upper limit of the colorbar
cbar_unit (tuple) – specifies the unit of the colorbar: e.g., (‘Jy’,’pixel’),(‘m-Jy’,’$mu$as$^2$’),[‘Tb’]
beamparams (list) – [fwhm_maj, fwhm_min, theta], set to plot beam contour
export_pdf (str) – path to exported PDF with plot
show (bool) – Display the plot if true
scinot (bool) – Display numbers/units in scientific notation
scale_lw (float) – Linewidth of the scale overlay
beam_lw (float) – Linewidth of the beam overlay
cbar_fontsize (float) – Fontsize of the text elements of the colorbar
axis (matplotlib.axes.Axes) – An axis object
scale_fontsize (float) – Fontsize of the scale label
power (float) – Passed to colorbar for division of ticks by 1e(power)
beamcolor (str) – color of the beam overlay
scalecolor (str) – color of the scale label overlay
- Returns
figure object with image
- Return type
(matplotlib.figure.Figure)
-
property
evec
¶ Return the E mode image vector
-
evpa
()[source]¶ Return the total evpa
Args:
- Returns
image average evpa (E of N) in radian
- Return type
(float)
-
property
evpavec
¶ Return the fractional polarization angle for each pixel
-
find_shift
(im_compare, pol=None, psize=None, target_fov=None, beamparams=[1.0, 1.0, 1.0], blur_frac=0.0, blursmall=False, scale='lin', gamma=0.5, dynamic_range=1000.0)[source]¶ Find image shift that maximizes normalized cross correlation with a second image im2. Finds shift only by comparison of the primary pol image.
- Parameters
im_compare (Image) – image with respect with to switch
pol (str) – which polarization image to compare. Default is self.pol_prim
psize (float) – pixel size of comparison image (rad). If None it is the smallest of the input image pizel sizes
target_fov (float) – fov of the comparison image (rad). If None it is twice the largest fov of the input images
beamparams (list) – the nominal Gaussian beam parameters [fovx, fovy, position angle]
blur_frac (float) – fractional beam to blur each image to before comparison
blursmall (bool) – True to blur the unpadded image rather than the large image.
scale (str) – compare images in ‘log’,’lin’,or ‘gamma’ scale
gamma (float) – exponent for gamma scale comparison
dynamic_range (float) – dynamic range for log and gamma scale comparisons
- Returns
(errormetric, im1_pad, im2_shift)
- Return type
(tuple)
-
fit_gauss
(units='rad')[source]¶ Determine the Gaussian parameters that short baselines would measure for the source by diagonalizing the image covariance matrix. Returns parameters only for the primary polarization!
- Parameters
units (string) – ‘rad’ returns values in radians, ‘natural’ returns FWHM in uas and PA in degrees
- Returns
a tuple (fwhm_maj, fwhm_min, theta) of the fit Gaussian parameters
- Return type
(tuple)
-
fit_gauss_empirical
(paramguess=None)[source]¶ Determine the Gaussian parameters that short baselines would measure Returns parameters only for the primary polarization!
-
fovx
()[source]¶ Return the image fov in x direction in radians.
Args:
- Returns
image fov in x direction (radian)
- Return type
(float)
-
fovy
()[source]¶ Returns the image fov in y direction in radians.
Args:
- Returns
image fov in y direction (radian)
- Return type
(float)
-
get_image_mf
(nu)[source]¶ Get image at a given frequency given the spectral information in self._mflist
-
imarr
(pol=None)[source]¶ Return the 2D image array of a given pol parameter.
- Parameters
pol (str) – I,Q,U or V for Stokes, or RR,LL,LR,RL for Circ
- Returns
2D image array of dimension (ydim, xdim)
- Return type
(numpy.array)
-
lin_polfrac
()[source]¶ Return the total fractional linear polarized flux
Args:
- Returns
image fractional linear polarized flux
- Return type
(float)
-
property
lrvec
¶ Return the imvec of LR
-
mask
(cutoff=0.05, beamparams=None, frac=0.0)[source]¶ Produce an image mask that shows all pixels above the specified cutoff frac of the max Works off the primary image
-
property
mvec
¶ Return the fractional polarization for each pixel
-
observe
(array, tint, tadv, tstart, tstop, bw, mjd=None, timetype='UTC', polrep_obs=None, elevmin=10.0, elevmax=85.0, no_elevcut_space=False, ttype='nfft', fft_pad_factor=2, fix_theta_GMST=False, sgrscat=False, add_th_noise=True, jones=False, inv_jones=False, opacitycal=True, ampcal=True, phasecal=True, frcal=True, dcal=True, rlgaincal=True, stabilize_scan_phase=False, stabilize_scan_amp=False, neggains=False, tau=0.1, taup=0.1, gain_offset=0.1, gainp=0.1, phase_std=- 1, dterm_offset=0.05, rlratio_std=0.0, rlphase_std=0.0, sigmat=None, phasesigmat=None, rlgsigmat=None, rlpsigmat=None, caltable_path=None, seed=False, verbose=True)[source]¶ Generate baselines from an array object and observe the image.
- Parameters
array (Array) – an array object containing sites with which to generate baselines
tint (float) – the scan integration time in seconds
tadv (float) – the uniform cadence between scans in seconds
tstart (float) – the start time of the observation in hours
tstop (float) – the end time of the observation in hours
bw (float) – the observing bandwidth in Hz
mjd (int) – the mjd of the observation, if set as different from the image mjd
timetype (str) – how to interpret tstart and tstop; either ‘GMST’ or ‘UTC’
polrep_obs (str) – ‘stokes’ or ‘circ’ sets the data polarimetric representation
elevmin (float) – station minimum elevation in degrees
elevmax (float) – station maximum elevation in degrees
no_elevcut_space (bool) – if True, do not apply elevation cut to orbiters
ttype (str) – “fast”, “nfft” or “dtft”
fft_pad_factor (float) – zero pad the image to fft_pad_factor * image size in the FFT
fix_theta_GMST (bool) – if True, stops earth rotation to sample fixed u,v
sgrscat (bool) – if True, the visibilites will be blurred by the Sgr A* kernel
add_th_noise (bool) – if True, baseline-dependent thermal noise is added
jones (bool) – if True, uses Jones matrix to apply mis-calibration effects otherwise uses old formalism without D-terms
inv_jones (bool) – if True, applies estimated inverse Jones matrix (not including random terms) to calibrate data
opacitycal (bool) – if False, time-dependent gaussian errors are added to opacities
ampcal (bool) – if False, time-dependent gaussian errors are added to station gains
phasecal (bool) – if False, time-dependent station-based random phases are added
frcal (bool) – if False, feed rotation angle terms are added to Jones matrix.
dcal (bool) – if False, time-dependent gaussian errors added to Jones matrix D-terms.
rlgaincal (bool) – if False, time-dependent gains are not equal for R and L pol
stabilize_scan_phase (bool) – if True, random phase errors are constant over scans
stabilize_scan_amp (bool) – if True, random amplitude errors are constant over scans
neggains (bool) – if True, force the applied gains to be <1
tau (float) – the base opacity at all sites, or a dict giving one opacity per site
taup (float) – the fractional std. dev. of the random error on the opacities
gainp (float) – the fractional std. dev. of the random error on the gains or a dict giving one std. dev. per site
gain_offset (float) – the base gain offset at all sites, or a dict giving one gain offset per site
phase_std (float) – std. dev. of LCP phase, or a dict giving one std. dev. per site a negative value samples from uniform
dterm_offset (float) – the base std. dev. of random additive error at all sites, or a dict giving one std. dev. per site
rlratio_std (float) – the fractional std. dev. of the R/L gain offset or a dict giving one std. dev. per site
rlphase_std (float) – std. dev. of R/L phase offset, or a dict giving one std. dev. per site a negative value samples from uniform
sigmat (float) – temporal std for a Gaussian Process used to generate gains. If sigmat=None then an iid gain noise is applied.
phasesigmat (float) – temporal std for a Gaussian Process used to generate phases. If phasesigmat=None then an iid gain noise is applied.
rlgsigmat (float) – temporal std deviation for a Gaussian Process used to generate R/L gain ratios. If rlgsigmat=None then an iid gain noise is applied.
rlpsigmat (float) – temporal std deviation for a Gaussian Process used to generate R/L phase diff. If rlpsigmat=None then an iid gain noise is applied.
caltable_path (string) – If not None, path and prefix for saving the applied caltable
seed (int) – seeds the random component of the noise terms. DO NOT set to 0!
verbose (bool) – print updates and warnings
- Returns
an observation object
- Return type
(Obsdata)
-
observe_same
(obs_in, ttype='nfft', fft_pad_factor=2, sgrscat=False, add_th_noise=True, jones=False, inv_jones=False, opacitycal=True, ampcal=True, phasecal=True, frcal=True, dcal=True, rlgaincal=True, stabilize_scan_phase=False, stabilize_scan_amp=False, neggains=False, taup=0.1, gain_offset=0.1, gainp=0.1, phase_std=- 1, dterm_offset=0.05, rlratio_std=0.0, rlphase_std=0.0, sigmat=None, phasesigmat=None, rlgsigmat=None, rlpsigmat=None, caltable_path=None, seed=False, verbose=True)[source]¶ Observe the image on the same baselines as an existing observation object and add noise.
- Parameters
obs_in (Obsdata) – the existing observation
ttype (str) – “fast” or “nfft” or “direct”
fft_pad_factor (float) – zero pad the image to fft_pad_factor * image size in FFT
sgrscat (bool) – if True, the visibilites will be blurred by the Sgr A* kernel
add_th_noise (bool) – if True, baseline-dependent thermal noise is added
jones (bool) – if True, uses Jones matrix to apply mis-calibration effects
inv_jones (bool) – if True, applies estimated inverse Jones matrix (not including random terms) to a priori calibrate data
opacitycal (bool) – if False, time-dependent gaussian errors are added to opacities
ampcal (bool) – if False, time-dependent gaussian errors are added to station gains
phasecal (bool) – if False, time-dependent station-based random phases are added
frcal (bool) – if False, feed rotation angle terms are added to Jones matrices.
dcal (bool) – if False, time-dependent gaussian errors added to D-terms.
rlgaincal (bool) – if False, time-dependent gains are not equal for R and L pol
stabilize_scan_phase (bool) – if True, random phase errors are constant over scans
stabilize_scan_amp (bool) – if True, random amplitude errors are constant over scans
neggains (bool) – if True, force the applied gains to be <1
taup (float) – the fractional std. dev. of the random error on the opacities
gainp (float) – the fractional std. dev. of the random error on the gains or a dict giving one std. dev. per site
gain_offset (float) – the base gain offset at all sites, or a dict giving one gain offset per site
phase_std (float) – std. dev. of LCP phase, or a dict giving one std. dev. per site a negative value samples from uniform
dterm_offset (float) – the base std. dev. of random additive error at all sites, or a dict giving one std. dev. per site
rlratio_std (float) – the fractional std. dev. of the R/L gain offset or a dict giving one std. dev. per site
rlphase_std (float) – std. dev. of R/L phase offset, or a dict giving one std. dev. per site a negative value samples from uniform
sigmat (float) – temporal std for a Gaussian Process used to generate gains. If sigmat=None then an iid gain noise is applied.
phasesigmat (float) – temporal std for a Gaussian Process used to generate phases. If phasesigmat=None then an iid gain noise is applied.
rlgsigmat (float) – temporal std deviation for a Gaussian Process used to generate R/L gain ratios. If rlgsigmat=None then an iid gain noise is applied.
rlpsigmat (float) – temporal std deviation for a Gaussian Process used to generate R/L phase diff. If rlpsigmat=None then an iid gain noise is applied.
caltable_path (string) – If not None, path and prefix for saving the applied caltable
seed (int) – seeds the random component of the noise terms. DO NOT set to 0!
verbose (bool) – print updates and warnings
- Returns
an observation object
- Return type
(Obsdata)
-
observe_same_nonoise
(obs, sgrscat=False, ttype='nfft', cache=False, fft_pad_factor=2, zero_empty_pol=True, verbose=True)[source]¶ Observe the image on the same baselines as an existing observation without noise.
- Parameters
obs (Obsdata) – the existing observation
sgrscat (bool) – if True, the visibilites will be blurred by the Sgr A* kernel
ttype (str) – “fast” or “nfft” or “direct”
cache (bool) – Use cached fft for ‘fast’ mode – deprecated, use nfft instead!
fft_pad_factor (float) – zero pad the image to fft_pad_factor * image size in FFT
zero_empty_pol (bool) – if True, returns zero vec if the polarization doesn’t exist. Otherwise return None
verbose (bool) – Boolean value controls output prints.
- Returns
an observation object with no noise
- Return type
(Obsdata)
-
observe_vex
(vex, source, t_int=0.0, tight_tadv=False, polrep_obs=None, ttype='nfft', fft_pad_factor=2, fix_theta_GMST=False, sgrscat=False, add_th_noise=True, jones=False, inv_jones=False, opacitycal=True, ampcal=True, phasecal=True, frcal=True, dcal=True, rlgaincal=True, stabilize_scan_phase=False, stabilize_scan_amp=False, neggains=False, tau=0.1, taup=0.1, gain_offset=0.1, gainp=0.1, phase_std=- 1, dterm_offset=0.05, rlratio_std=0.0, rlphase_std=0.0, sigmat=None, phasesigmat=None, rlgsigmat=None, rlpsigmat=None, caltable_path=None, seed=False, verbose=True)[source]¶ Generate baselines from a vex file and observes the image.
- Parameters
vex (Vex) – an vex object containing sites and scan information
source (str) – the source to observe
t_int (float) – if not zero, overrides the vex scan lengths
tight_tadv (float) – if True, advance right after each integration, otherwise advance after 2x the scan length
polrep_obs (str) – ‘stokes’ or ‘circ’ sets the data polarimetric representation
ttype (str) – “fast” or “nfft” or “dtft”
fft_pad_factor (float) – zero pad the image to fft_pad_factor * image size in FFT
fix_theta_GMST (bool) – if True, stops earth rotation to sample fixed u,v
sgrscat (bool) – if True, the visibilites will be blurred by the Sgr A* kernel
add_th_noise (bool) – if True, baseline-dependent thermal noise is added
jones (bool) – if True, uses Jones matrix to apply mis-calibration effects otherwise uses old formalism without D-terms
inv_jones (bool) – if True, applies estimated inverse Jones matrix (not including random terms) to calibrate data
opacitycal (bool) – if False, time-dependent gaussian errors are added to opacities
ampcal (bool) – if False, time-dependent gaussian errors are added to station gains
phasecal (bool) – if False, time-dependent station-based random phases are added
frcal (bool) – if False, feed rotation angle terms are added to Jones matrix.
dcal (bool) – if False, time-dependent gaussian errors added to Jones matrix D-terms.
rlgaincal (bool) – if False, time-dependent gains are not equal for R and L pol
stabilize_scan_phase (bool) – if True, random phase errors are constant over scans
stabilize_scan_amp (bool) – if True, random amplitude errors are constant over scans
neggains (bool) – if True, force the applied gains to be <1
tau (float) – the base opacity at all sites, or a dict giving one opacity per site
taup (float) – the fractional std. dev. of the random error on the opacities
gainp (float) – the fractional std. dev. of the random error on the gains or a dict giving one std. dev. per site
gain_offset (float) – the base gain offset at all sites, or a dict giving one gain offset per site
phase_std (float) – std. dev. of LCP phase, or a dict giving one std. dev. per site a negative value samples from uniform
dterm_offset (float) – the base std. dev. of random additive error at all sites, or a dict giving one std. dev. per site
rlratio_std (float) – the fractional std. dev. of the R/L gain offset or a dict giving one std. dev. per site
rlphase_std (float) – std. dev. of R/L phase offset, or a dict giving one std. dev. per site a negative value samples from uniform
sigmat (float) – temporal std for a Gaussian Process used to generate gains. If sigmat=None then an iid gain noise is applied.
phasesigmat (float) – temporal std for a Gaussian Process used to generate phases. If phasesigmat=None then an iid gain noise is applied.
rlgsigmat (float) – temporal std deviation for a Gaussian Process used to generate R/L gain ratios. If rlgsigmat=None then an iid gain noise is applied.
rlpsigmat (float) – temporal std deviation for a Gaussian Process used to generate R/L phase diff. If rlpsigmat=None then an iid gain noise is applied.
caltable_path (string) – If not None, path and prefix for saving the applied caltable
seed (int) – seeds the random component of the noise terms. DO NOT set to 0!
verbose (bool) – print updates and warnings
- Returns
an observation object
- Return type
(Obsdata)
-
orth_chi
()[source]¶ Rotate the EVPA 90 degrees
Args:
- Returns
image with rotated EVPA
- Return type
(Image)
-
overlay_display
(im_list, color_coding=array([[1, 0, 1], [0, 1, 0]]), export_pdf='', show=True, f=False, shift=[0, 0], final_fov=False, interp='gaussian', scale='lin', gamma=0.5, dynamic_range=[1000.0], rescale=True)[source]¶ Overlay primary polarization images of a list of images to compare structures.
- Parameters
im_list (list) – list of images to align to the current image
color_coding (numpy.array) – Color coding of each image in the composite
f (matplotlib.pyplot.figure) – Figure to overlay on top of
export_pdf (str) – path to exported PDF with plot
show (bool) – Display the plot if true
shift (list) – list of manual image shifts, otherwise use the shift from maximum cross-correlation
final_fov (float) – fov of the comparison image (rad). If False it is the largestinput image fov
scale (str) – compare images in ‘log’,’lin’,or ‘gamma’ scale
gamma (float) – exponent for gamma scale comparison
dynamic_range (float) – dynamic range for log and gamma scale comparisons
- Returns
figure object with image
- Return type
(matplotlib.figure.Figure)
-
pad
(fovx, fovy)[source]¶ Pad an image to new fov_x by fov_y in radian. :param fovx: new fov in x dimension (rad) :type fovx: float :param fovy: new fov in y dimension (rad) :type fovy: float
- Returns
padded image
- Return type
im_pad (Image)
-
property
pvec
¶ Return the polarization magnitude for each pixel
-
regrid_image
(targetfov, npix, interp='linear')[source]¶ Resample the image to new (square) dimensions.
-
resample_square
(xdim_new, ker_size=5)[source]¶ Exactly resample a square image to new dimensions using the pulse function.
-
sample_uv
(uv, polrep_obs='stokes', sgrscat=False, ttype='nfft', cache=False, fft_pad_factor=2, zero_empty_pol=True, verbose=True)[source]¶ Sample the image on the selected uv points without creating an Obsdata object.
- Parameters
uv (ndarray) – an array of uv points
polrep_obs (str) – ‘stokes’ or ‘circ’ sets the data polarimetric representation
sgrscat (bool) – if True, the visibilites will be blurred by the Sgr A* kernel
ttype (str) – “fast” or “nfft” or “direct”
cache (bool) – Use cached fft for ‘fast’ mode – deprecated, use nfft instead!
fft_pad_factor (float) – zero pad the image to fft_pad_factor * image size in FFT
zero_empty_pol (bool) – if True, returns zero vec if the polarization doesn’t exist. Otherwise return None
verbose (bool) – Boolean value controls output prints.
- Returns
a list of [I,Q,U,V] visibilities
- Return type
(list)
-
save_fits
(fname)[source]¶ Save image data to a fits file.
- Parameters
fname (str) – path to output fits file
Returns:
-
save_txt
(fname)[source]¶ Save image data to text file.
- Parameters
fname (str) – path to output text file
Returns:
-
shift_fft
(shift)[source]¶ - Shift the image by a given vector in radians.
This allows non-integer pixel shifts, via FFT.
-
sourcevec
()[source]¶ Return the source position vector in geocentric coordinates at 0h GMST.
Args:
- Returns
normal vector pointing to source in geocentric coordinates (m)
- Return type
(numpy.array)
-
switch_polrep
(polrep_out='stokes', pol_prim_out=None)[source]¶ Return a new image with the polarization representation changed :param polrep_out: the polrep of the output data :type polrep_out: str :param pol_prim_out: The default image: I,Q,U or V for Stokes, RR,LL,LR,RL for circ :type pol_prim_out: str
- Returns
new Image object with potentially different polrep
- Return type
(Image)
-
threshold
(cutoff=0.05, beamparams=None, frac=0.0, fill_val=None)[source]¶ Apply a hard threshold to the primary polarization image. Leave other polarizations untouched.
- Parameters
cutoff (float) – Mask pixels with intensities greater than cuttoff * max
beamparams (list) – either [fwhm_maj, fwhm_min, pos_ang] or a single fwhm
frac (float) – the fraction of nominal beam to blur with
fill_val (float) – masked pixels are set to this value. If fill_val==None, they are set to the min unmasked value
- Returns
output mask image
- Return type
(Image)
-
-
ehtim.image.
blur_mf
(im, freqs, kernel, fit_order=2)[source]¶ blur multifrequncy images with the same beam
-
ehtim.image.
get_specim
(imlist, reffreq, fit_order=2)[source]¶ get the spectral index/curvature from a list of images
-
ehtim.image.
load_fits
(fname, aipscc=False, pulse=<function trianglePulse2D>, polrep='stokes', pol_prim=None, zero_pol=False)[source]¶ Read in an image from a FITS file.
- Parameters
fname (str) – path to input fits file
aipscc (bool) – if True, then AIPS CC table will be loaded
pulse (function) – The function convolved with the pixel values for continuous image.
polrep (str) – polarization representation, either ‘stokes’ or ‘circ’
pol_prim (str) – The default image: I,Q,U or V for Stokes, RR,LL,LR,RL for Circular
zero_pol (bool) – If True, loads any missing polarizations as zeros
- Returns
loaded image object
- Return type
(Image)
-
ehtim.image.
load_image
(image, display=False, aipscc=False)[source]¶ Read in an image from a text, .fits, .h5, or ehtim.image.Image object
- Parameters
image (str/Image) – path to input file
display (boolean) – determine whether to display the image default
aipscc (boolean) – if True, then AIPS CC table will be loaded instead of the original brightness distribution.
- Returns
loaded image object (boolean): False if the image cannot be read
- Return type
(Image)
-
ehtim.image.
load_txt
(fname, polrep='stokes', pol_prim=None, pulse=<function trianglePulse2D>, zero_pol=True)[source]¶ Read in an image from a text file.
- Parameters
fname (str) – path to input text file
pulse (function) – The function convolved with the pixel values for continuous image.
polrep (str) – polarization representation, either ‘stokes’ or ‘circ’
pol_prim (str) – The default image: I,Q,U or V for Stokes, RR,LL,LR,RL for Circular
zero_pol (bool) – If True, loads any missing polarizations as zeros
- Returns
loaded image object
- Return type
(Image)
-
ehtim.image.
make_empty
(npix, fov, ra, dec, rf=230000000000.0, source='SgrA', polrep='stokes', pol_prim=None, pulse=<function trianglePulse2D>, mjd=51544, time=0.0)[source]¶ Make an empty square image.
- Parameters
npix (int) – the pixel size of each axis
fov (float) – the field of view of each axis in radians
ra (float) – The source Right Ascension in fractional hours
dec (float) – The source declination in fractional degrees
rf (float) – The image frequency in Hz
source (str) – The source name
polrep (str) – polarization representation, either ‘stokes’ or ‘circ’
pol_prim (str) – The default image: I,Q,U or V for Stokes, RR,LL,LR,RL for Circular
pulse (function) – The function convolved with the pixel values for continuous image.
mjd (int) – The integer MJD of the image
time (float) – The observing time of the image (UTC hours)
- Returns
an image object
- Return type
(Image)
-
ehtim.image.
make_square
(obs, npix, fov, pulse=<function trianglePulse2D>, polrep='stokes', pol_prim=None)[source]¶ Make an empty square image.
- Parameters
obs (Obsdata) – an obsdata object with the image metadata
npix (int) – the pixel size of each axis
fov (float) – the field of view of each axis in radians
pulse (function) – the function convolved with the pixel values for continuous image
polrep (str) – polarization representation, either ‘stokes’ or ‘circ’
pol_prim (str) – The default image: I,Q,U or V for Stokes, or RR,LL,LR,RL for Circular
- Returns
an image object
- Return type
(Image)