diff --git a/package/lib/utils.py b/package/lib/utils.py index e28bbff..7d38335 100755 --- a/package/lib/utils.py +++ b/package/lib/utils.py @@ -4,6 +4,14 @@ import numpy as np def rot2D(ang): """ Return the 2D rotation matrix of given angle in degrees + ---------- + Inputs: + ang : float + Angle in degrees + ---------- + Returns: + rot_mat : numpy.ndarray + 2D matrix of shape (2,2) to perform vector rotation at angle "ang". """ alpha = np.pi * ang / 180 return np.array([[np.cos(alpha), np.sin(alpha)], [-np.sin(alpha), np.cos(alpha)]]) @@ -13,6 +21,14 @@ def princ_angle(ang): """ Return the principal angle in the 0° to 180° quadrant as PA is always defined at p/m 180°. + ---------- + Inputs: + ang : float, numpy.ndarray + Angle in degrees. Can be an array of angles. + ---------- + Returns: + princ_ang : float, numpy.ndarray + Principal angle in the 0°-180° quadrant in the same shape as input. """ if not isinstance(ang, np.ndarray): A = np.array([ang]) @@ -32,6 +48,21 @@ def PCconf(QN, UN, QN_ERR, UN_ERR): """ Compute the confidence level for 2 parameters polarisation degree and polarisation angle from the PCUBE analysis. + ---------- + Inputs: + QN : float, numpy.ndarray + Normalized Q Stokes flux. + UN : float, numpy.ndarray + Normalized U Stokes flux. + QN_ERR : float, numpy.ndarray + Normalized error on Q Stokes flux. + UN_ERR : float, numpy.ndarray + Normalized error on U Stokes flux. + ---------- + Returns: + conf : numpy.ndarray + 2D matrix of same shape as input containing the confidence on the polarization + computation between 0 and 1 for 2 parameters of interest (Q and U Stokes fluxes). """ mask = np.logical_and(QN_ERR > 0.0, UN_ERR > 0.0) conf = np.full(QN.shape, -1.0) @@ -43,6 +74,19 @@ def PCconf(QN, UN, QN_ERR, UN_ERR): def CenterConf(mask, PA, sPA): """ Compute the confidence map for the position of the center of emission. + ---------- + Inputs: + mask : bool, numpy.ndarray + Mask of the polarization vectors from which the center of emission should be drawn. + PA : float, numpy.ndarray + 2D matrix containing the computed polarization angle. + sPA : float, numpy.ndarray + 2D matrix containing the total uncertainties on the polarization angle. + ---------- + Returns: + conf : numpy.ndarray + 2D matrix of same shape as input containing the confidence on the polarization + computation between 0 and 1 for 2 parameters of interest (Q and U Stokes fluxes). """ chi2 = np.full(PA.shape, np.nan) conf = np.full(PA.shape, -1.0) @@ -64,17 +108,36 @@ def CenterConf(mask, PA, sPA): from scipy.special import gammainc conf[np.isfinite(PA)] = 1.0 - gammainc(0.5, 0.5 * chi2[np.isfinite(PA)]) - result = minimize(chisq, np.array(PA.shape) / 2.0, bounds=[(0, PA.shape[1]), (0.0, PA.shape[0])]) + c0 = np.unravel_index(np.argmax(conf), conf.shape)[::-1] + result = minimize(chisq, c0, bounds=[(0, PA.shape[1]), (0.0, PA.shape[0])]) if result.success: print("Center of emission found") else: - print("Center of emission not found") + print("Center of emission not found", result) return conf, result.x def sci_not(v, err, rnd=1, out=str): """ - Return the scientifque error notation as a string. + Return the scientific error notation as a string. + ---------- + Inputs: + v : float + Value to be transformed into scientific notation. + err : float + Error on the value to be transformed into scientific notation. + rnd : int + Number of significant numbers for the scientific notation. + Default to 1. + out : str or other + Format in which the notation should be returned. "str" means the notation + is returned as a single string, "other" means it is returned as a list of "str". + Default to str. + ---------- + Returns: + conf : numpy.ndarray + 2D matrix of same shape as input containing the confidence on the polarization + computation between 0 and 1 for 2 parameters of interest (Q and U Stokes fluxes). """ power = -int(("%E" % v)[-3:]) + 1 output = [r"({0}".format(round(v * 10**power, rnd)), round(v * 10**power, rnd)] @@ -95,6 +158,16 @@ def wcs_PA(PC21, PC22): """ Return the position angle in degrees to the North direction of a wcs from the values of coefficient of its transformation matrix. + ---------- + Inputs: + PC21 : float + Value of the WCS matric PC[1,0] + PC22 : float + Value of the WCS matric PC[1,1] + ---------- + Returns: + orient : float + Angle in degrees between the North direction and the Up direction of the WCS. """ if (abs(PC21) > abs(PC22)) and (PC21 >= 0): orient = -np.arccos(PC22) * 180.0 / np.pi diff --git a/package/test_center.py b/package/test_center.py index 54c21d4..bedfda8 100644 --- a/package/test_center.py +++ b/package/test_center.py @@ -23,15 +23,20 @@ NGC1068snr[NGC1068["POL_DEG_ERR"].data > 0.0] = ( ) NGC1068centconf, NGC1068center = CenterConf(NGC1068conf > 0.99, NGC1068["POL_ANG"].data, NGC1068["POL_ANG_ERR"].data) +NGC1068pos = WCS(NGC1068[0].header).pixel_to_world(*NGC1068center) -figngc, axngc = plt.subplots(1, 2, layout="tight", figsize=(18,9), subplot_kw=dict(projection=WCS(NGC1068[0].header))) +figngc, axngc = plt.subplots(1, 2, layout="tight", figsize=(18, 9), subplot_kw=dict(projection=WCS(NGC1068[0].header)), sharex=True, sharey=True) axngc[0].set(xlabel="RA", ylabel="DEC", title="NGC1069 intensity map with SNR and confidence contours") -imngc = axngc[0].imshow(NGC1068["I_STOKES"].data * NGC1068["I_STOKES"].header["PHOTFLAM"], norm=LogNorm(), cmap="inferno") +vmin, vmax = ( + 0.5 * np.median(NGC1068["I_STOKES"].data[NGC1068mask]) * NGC1068[0].header["PHOTFLAM"], + np.max(NGC1068["I_STOKES"].data[NGC1068mask]) * NGC1068[0].header["PHOTFLAM"], +) +imngc = axngc[0].imshow(NGC1068["I_STOKES"].data * NGC1068["I_STOKES"].header["PHOTFLAM"], norm=LogNorm(vmin, vmax), cmap="inferno") ngcsnrcont = axngc[0].contour(NGC1068snr, levelssnr, colors="b") ngcconfcont = axngc[0].contour(NGC1068conf, levelsconf, colors="r") -ngcconfcenter = axngc[0].plot(*np.unravel_index(np.argmax(NGC1068centconf), NGC1068centconf.shape)[::-1], "k+", label="Best confidence for center") -ngcconfcentcont = axngc[0].contour(NGC1068centconf, 1.-levelsconf, colors="k") +ngcconfcenter = axngc[0].plot(*NGC1068center, marker="+",color="gray", label="Best confidence for center: {0}".format(NGC1068pos.to_string('hmsdms'))) +ngcconfcentcont = axngc[0].contour(NGC1068centconf, [0.01], colors="gray") handles, labels = axngc[0].get_legend_handles_labels() labels.append("SNR contours") handles.append(Rectangle((0, 0), 1, 1, fill=False, ec=ngcsnrcont.collections[0].get_edgecolor()[0])) @@ -43,14 +48,14 @@ axngc[0].legend(handles=handles, labels=labels) axngc[1].set(xlabel="RA", ylabel="DEC", title="Location of the nucleus confidence map") ngccent = axngc[1].imshow(NGC1068centconf, vmin=0.0, cmap="inferno") -ngccentcont = axngc[1].contour(NGC1068centconf, 1.-levelsconf, colors="grey") -ngccentcenter = axngc[1].plot(*np.unravel_index(np.argmax(NGC1068centconf), NGC1068centconf.shape)[::-1], "k+", label="Best confidence for center") +ngccentcont = axngc[1].contour(NGC1068centconf, [0.01], colors="gray") +ngccentcenter = axngc[1].plot(*NGC1068center, marker="+",color="gray", label="Best confidence for center: {0}".format(NGC1068pos.to_string('hmsdms'))) handles, labels = axngc[1].get_legend_handles_labels() labels.append("CONF99 contour") handles.append(Rectangle((0, 0), 1, 1, fill=False, ec=ngccentcont.collections[0].get_edgecolor()[0])) axngc[1].legend(handles=handles, labels=labels) -figngc.savefig("NGC1068_center.pdf",dpi=150,facecolor="None") +figngc.savefig("NGC1068_center.pdf", dpi=150, facecolor="None") ################################################################################################### @@ -68,16 +73,21 @@ MRK463Esnr[MRK463E["POL_DEG_ERR"].data > 0.0] = ( ) MRK463Ecentconf, MRK463Ecenter = CenterConf(MRK463Econf > 0.99, MRK463E["POL_ANG"].data, MRK463E["POL_ANG_ERR"].data) +MRK463Epos = WCS(MRK463E[0].header).pixel_to_world(*MRK463Ecenter) -figmrk, axmrk = plt.subplots(1, 2, layout="tight", figsize=(18,9), subplot_kw=dict(projection=WCS(MRK463E[0].header))) +figmrk, axmrk = plt.subplots(1, 2, layout="tight", figsize=(18, 9), subplot_kw=dict(projection=WCS(MRK463E[0].header)), sharex=True, sharey=True) axmrk[0].set(xlabel="RA", ylabel="DEC", title="NGC1069 intensity map with SNR and confidence contours") -immrk = axmrk[0].imshow(MRK463E["I_STOKES"].data * MRK463E["I_STOKES"].header["PHOTFLAM"], norm=LogNorm(), cmap="inferno") +vmin, vmax = ( + 0.5 * np.median(MRK463E["I_STOKES"].data[MRK463Emask]) * MRK463E[0].header["PHOTFLAM"], + np.max(MRK463E["I_STOKES"].data[MRK463Emask]) * MRK463E[0].header["PHOTFLAM"], +) +immrk = axmrk[0].imshow(MRK463E["I_STOKES"].data * MRK463E["I_STOKES"].header["PHOTFLAM"], norm=LogNorm(vmin, vmax), cmap="inferno") mrksnrcont = axmrk[0].contour(MRK463Esnr, levelssnr, colors="b") mrkconfcont = axmrk[0].contour(MRK463Econf, levelsconf, colors="r") -mrkconfcenter = axmrk[0].plot(*np.unravel_index(np.argmax(MRK463Ecentconf), MRK463Ecentconf.shape)[::-1], "k+", label="Best confidence for center") -mrkconfcentcont = axmrk[0].contour(MRK463Ecentconf, 1.-levelsconf, colors="k") -handles, labels = axmrk[1].get_legend_handles_labels() +mrkconfcenter = axmrk[0].plot(*MRK463Ecenter, marker="+",color="gray", label="Best confidence for center: {0}".format(MRK463Epos.to_string('hmsdms'))) +mrkconfcentcont = axmrk[0].contour(MRK463Ecentconf, [0.01], colors="gray") +handles, labels = axmrk[0].get_legend_handles_labels() labels.append("SNR contours") handles.append(Rectangle((0, 0), 1, 1, fill=False, ec=mrksnrcont.collections[0].get_edgecolor()[0])) labels.append("CONF99 contour") @@ -88,12 +98,12 @@ axmrk[0].legend(handles=handles, labels=labels) axmrk[1].set(xlabel="RA", ylabel="DEC", title="Location of the nucleus confidence map") mrkcent = axmrk[1].imshow(MRK463Ecentconf, vmin=0.0, cmap="inferno") -mrkcentcont = axmrk[1].contour(MRK463Ecentconf, 1.-levelsconf, colors="grey") -mrkcentcenter = axmrk[1].plot(*np.unravel_index(np.argmax(MRK463Ecentconf), MRK463Ecentconf.shape)[::-1], "k+", label="Best confidence for center") +mrkcentcont = axmrk[1].contour(MRK463Ecentconf, [0.01], colors="gray") +mrkcentcenter = axmrk[1].plot(*MRK463Ecenter, marker="+",color="gray", label="Best confidence for center: {0}".format(MRK463Epos.to_string('hmsdms'))) handles, labels = axmrk[1].get_legend_handles_labels() labels.append("CONF99 contour") handles.append(Rectangle((0, 0), 1, 1, fill=False, ec=mrkcentcont.collections[0].get_edgecolor()[0])) axmrk[1].legend(handles=handles, labels=labels) -figmrk.savefig("MRK463E_center.pdf",dpi=150,facecolor="None") +figmrk.savefig("MRK463E_center.pdf", dpi=150, facecolor="None") plt.show()