Quad matching in DIPOL polarimetry images#

%autoawait off
%matplotlib inline
%run 01_notebook_configuration.py

# To reduce disk usage, only a central cut of of the full field of the DIPOL 
# camera is actually saved.
# polarimetry (cut)
from IPython.display import display
from iop4lib.db import RawFit, ReducedFit
from iop4lib.utils.plotting import imshow_w_sources
rf_pol = RawFit.by_fileloc(f"OSN-T090/2023-11-06/BLLAC_R_IAR-0760.fts")
display(rf_pol)
imshow_w_sources(rf_pol.mdata)
RawFit(id=7):
- telescope: OSN-T090
- night: 2023-11-06
- filename: BLLAC_R_IAR-0760.fts
- instrument: DIPOL
- imgtype: LIGHT
- size: 1104x896
- obsmode: POLARIMETRY
- band: R
- exptime: 20.0
- flags: DOWNLOADED,CLASSIFIED,BUILT_REDUCED
../_images/b10a88ea1aca54eb557f36ec883608edbaddf3ca163d7295d0f6d8bb8f845430.png

The position of the cut is defined in the header of the raw file

rf_pol.header['XORGSUBF'], rf_pol.header['YORGSUBF']
(1496, 1000)

We can compare it with the full field image of the photometry field

# photometry filed (full field)
rf_phot = RawFit.by_fileloc(f"OSN-T090/2023-11-06/BLLac_IAR-0001R.fit")
display(rf_phot)
imshow_w_sources(rf_phot.mdata)
RawFit(id=8):
- telescope: OSN-T090
- night: 2023-11-06
- filename: BLLac_IAR-0001R.fit
- instrument: DIPOL
- imgtype: LIGHT
- size: 4144x2822
- obsmode: PHOTOMETRY
- band: R
- exptime: 300.0
- flags: DOWNLOADED,CLASSIFIED,BUILT_REDUCED
../_images/ae0d71820ef5d3a6bc04aa1aff7f8716e0fd0546930d1bb7aa9508b6dcb13c8f.png

Let’s build the photometry field first

redf_phot = ReducedFit.create(rf_phot)
redf_phot.build_file()
2024-09-24 16:03 - pid 4241 - [reducedfit.py:121] - DEBUG - DB entry of ReducedFit for OSN-T090/2023-11-06/BLLac_IAR-0001R.fit already exists, it will be used instead.
2024-09-24 16:03 - pid 4241 - [instrument.py:373] - WARNING - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: masterflat in this epoch could not be found, attemptying adjacent epochs.
2024-09-24 16:03 - pid 4241 - [instrument.py:351] - WARNING - Master Flat from epoch <Epoch 1 | OSN-T090/2023-12-13> is more than 1 week away from epoch <Epoch 3 | OSN-T090/2023-11-06>.
2024-09-24 16:03 - pid 4241 - [instrument.py:480] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: building file
2024-09-24 16:03 - pid 4241 - [dipol.py:233] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: applying masters
2024-09-24 16:03 - pid 4241 - [instrument.py:486] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: performing astrometric calibration
2024-09-24 16:03 - pid 4241 - [astrometry.py:179] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: 6 different combinations of parameters to try.
2024-09-24 16:03 - pid 4241 - [astrometry.py:182] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: attempt 1 / 6, ({'keep_n_seg': 300, 'border_margin_px': 20, 'output_logodds_threshold': 14, 'n_rms_seg': 1.5, 'bkg_filter_size': 11, 'bkg_box_size': 32, 'seg_fwhm': 1.0, 'seg_kernel_size': None, 'npixels': 32, 'allsky': False, 'd_eps': 4.0, 'dx_eps': 4.0, 'dy_eps': 2.0, 'dx_min': 150, 'dx_max': 300, 'dy_min': 0, 'dy_max': 50, 'd_min': 150, 'd_max': 250, 'bins': 400, 'hist_range': (0, 500), 'position_hint': PositionHint(ra_deg=330.6804125, dec_deg=42.27766361111111, radius_deg=0.23050000000000004), 'size_hint': SizeHint(lower_arcsec_per_pixel=0.1273, upper_arcsec_per_pixel=0.14070000000000002)}) ...
2024-09-24 16:03 - pid 4241 - [astrometry.py:275] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: len(pos_seg)=168
2024-09-24 16:03 - pid 4241 - [astrometry.py:282] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: Removing segments within 20 px from border.
2024-09-24 16:03 - pid 4241 - [astrometry.py:289] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: seg pairs -> 74 (45.7%), seg_disp_sign=[-198.42540744   -5.23306087]
2024-09-24 16:03 - pid 4241 - [astrometry.py:291] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: seg pairs best -> 59 (36.4%), seg_disp_sign_best=[-205.95556765    2.43153595]
2024-09-24 16:03 - pid 4241 - [astrometry.py:293] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: seg pairs xy -> 31, disp_sign_xy=[-209.20635217   13.86089392]
2024-09-24 16:03 - pid 4241 - [astrometry.py:295] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: seg pairs xy best -> 31 (19.1%), seg_disp_sign_xy_best=[-209.20635217   13.86089392]
2024-09-24 16:03 - pid 4241 - [astrometry.py:312] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: trying to solve astrometry with Seg Best XY Pairs (n=31) (output_logodds_threshold=14).
2024-09-24 16:03 - pid 4273 - [__init__.py:212] - INFO - loaded 240 index files
2024-09-24 16:03 - pid 4273 - [__init__.py:276] - INFO - solve 1: start
2024-09-24 16:03 - pid 4273 - [__init__.py:276] - INFO - solve 1: slice=[0, 25[ (1 / 1), index="5200/index-5200-14.fits" (1 / 5)
2024-09-24 16:03 - pid 4273 - [__init__.py:276] - INFO - solve 1: logodds=164.959, matches=26, conflicts=0, distractors=1, ra=330.682, dec=42.2766, scale=0.133698, index="5200/index-5200-14.fits"
2024-09-24 16:03 - pid 4241 - [astrometry.py:327] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: Seg Best XY Pairs (n=31) worked.
2024-09-24 16:03 - pid 4241 - [astrometry.py:328] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: bm.index_path=PosixPath('/mnt/astrometry_cache/5200/index-5200-14.fits')
2024-09-24 16:03 - pid 4241 - [astrometry.py:329] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: bm.center_ra_deg=330.6820237797853
2024-09-24 16:03 - pid 4241 - [astrometry.py:330] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: bm.center_dec_deg=42.27656696901512
2024-09-24 16:03 - pid 4241 - [astrometry.py:331] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: bm.scale_arcsec_per_pixel=0.13369775052658256
2024-09-24 16:03 - pid 4241 - [astrometry.py:332] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: bm.logodds=164.959228515625
2024-09-24 16:03 - pid 4241 - [astrometry.py:192] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: WCS built with attempt 1 / 6 ({'keep_n_seg': 300, 'border_margin_px': 20, 'output_logodds_threshold': 14, 'n_rms_seg': 1.5, 'bkg_filter_size': 11, 'bkg_box_size': 32, 'seg_fwhm': 1.0, 'seg_kernel_size': None, 'npixels': 32, 'allsky': False, 'd_eps': 4.0, 'dx_eps': 4.0, 'dy_eps': 2.0, 'dx_min': 150, 'dx_max': 300, 'dy_min': 0, 'dy_max': 50, 'd_min': 150, 'd_max': 250, 'bins': 400, 'hist_range': (0, 500), 'position_hint': PositionHint(ra_deg=330.6804125, dec_deg=42.27766361111111, radius_deg=0.23050000000000004), 'size_hint': SizeHint(lower_arcsec_per_pixel=0.1273, upper_arcsec_per_pixel=0.14070000000000002)}).
2024-09-24 16:03 - pid 4241 - [astrometry.py:204] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: building summary images.
2024-09-24 16:03 - pid 4241 - [plotting.py:476] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: plotting astrometry summary image of background substraction results
2024-09-24 16:04 - pid 4241 - [plotting.py:488] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: plotting astrometry summary image of segmentation results
2024-09-24 16:04 - pid 4241 - [plotting.py:582] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: plotting astrometry summary image of astrometry results
2024-09-24 16:04 - pid 4241 - [plotting.py:214] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: found 4 catalog sources in field: [AstroSource.objects.get(name='2200+420'), AstroSource.objects.get(name='BL Lacertae C'), AstroSource.objects.get(name='BL Lacertae test 2'), AstroSource.objects.get(name='BL Lacertae test 4')]
2024-09-24 16:04 - pid 4241 - [instrument.py:425] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: saving WCSs to FITS header.
2024-09-24 16:04 - pid 4241 - [instrument.py:497] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: astrometric calibration was successful.
2024-09-24 16:04 - pid 4241 - [instrument.py:500] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: searching for sources in field...
2024-09-24 16:04 - pid 4241 - [instrument.py:503] - DEBUG - <ReducedFit 8 | OSN-T090/2023-11-06/BLLac_IAR-0001R.fit>: found 4 sources in field.

As we can see, the astrometric calibration was done with the blind astrometry.net solver as before.

Next, if we go to the polarimetry field, we see that there are not enough stars in the image to do the blind astrometric calibration.

The calibration of this file will use the previously calibrated photometry field and compare quads of stars in both images to find the transformation between them, as we can see in the log:

redf_pol = ReducedFit.create(rf_pol)
redf_pol.build_file()
2024-09-24 16:04 - pid 4241 - [reducedfit.py:121] - DEBUG - DB entry of ReducedFit for OSN-T090/2023-11-06/BLLAC_R_IAR-0760.fts already exists, it will be used instead.
2024-09-24 16:04 - pid 4241 - [instrument.py:373] - WARNING - <ReducedFit 7 | OSN-T090/2023-11-06/BLLAC_R_IAR-0760.fts>: masterflat in this epoch could not be found, attemptying adjacent epochs.
2024-09-24 16:04 - pid 4241 - [instrument.py:480] - DEBUG - <ReducedFit 7 | OSN-T090/2023-11-06/BLLAC_R_IAR-0760.fts>: building file
2024-09-24 16:04 - pid 4241 - [dipol.py:233] - DEBUG - <ReducedFit 7 | OSN-T090/2023-11-06/BLLAC_R_IAR-0760.fts>: applying masters
2024-09-24 16:04 - pid 4241 - [instrument.py:486] - DEBUG - <ReducedFit 7 | OSN-T090/2023-11-06/BLLAC_R_IAR-0760.fts>: performing astrometric calibration
2024-09-24 16:04 - pid 4241 - [connectionpool.py:1051] - DEBUG - Starting new HTTPS connection (1): simbad.cds.unistra.fr:443
2024-09-24 16:04 - pid 4241 - [connectionpool.py:546] - DEBUG - https://simbad.cds.unistra.fr:443 "POST /simbad/sim-script HTTP/11" 200 None
2024-09-24 16:04 - pid 4241 - [dipol.py:444] - DEBUG - target_src.srctype='blazar'
2024-09-24 16:04 - pid 4241 - [dipol.py:445] - DEBUG - n_estimate=11
2024-09-24 16:04 - pid 4241 - [dipol.py:446] - DEBUG - n_estimate_centered=5
2024-09-24 16:04 - pid 4241 - [dipol.py:447] - DEBUG - redf_phot=ReducedFit.objects.get(id=8)
2024-09-24 16:04 - pid 4241 - [dipol.py:448] - DEBUG - n_expected_simbad_sources=5
2024-09-24 16:04 - pid 4241 - [dipol.py:449] - DEBUG - n_expected_calibrators=4
2024-09-24 16:04 - pid 4241 - [dipol.py:503] - DEBUG - Trying attempt: {'method': '_build_wcs_for_polarimetry_images_photo_quads', 'conds': {'redf_phot__ne': None}, 'args': {'n_seg_threshold': [1.1, 1.0, 0.9], 'npixels': [64, 32], 'min_quad_distance': [4.0, 8.0]}} for <ReducedFit 7 | OSN-T090/2023-11-06/BLLAC_R_IAR-0760.fts>.
2024-09-24 16:04 - pid 4241 - [dipol.py:515] - INFO - Attempt: {'method': '_build_wcs_for_polarimetry_images_photo_quads', 'conds': {'redf_phot__ne': None}, 'args': {'n_seg_threshold': [1.1, 1.0, 0.9], 'npixels': [64, 32], 'min_quad_distance': [4.0, 8.0]}}: _build_wcs_for_polarimetry_images_photo_quads with {'n_seg_threshold': 1.1, 'npixels': 64, 'min_quad_distance': 4.0} for <ReducedFit 7 | OSN-T090/2023-11-06/BLLAC_R_IAR-0760.fts>.
2024-09-24 16:04 - pid 4241 - [dipol.py:596] - DEBUG - Invoked with n_seg_threshold=1.1, npixels=64
2024-09-24 16:04 - pid 4241 - [dipol.py:620] - DEBUG - Using 10 sources in polarimetry field and 10 in photometry field.
2024-09-24 16:04 - pid 4241 - [dipol.py:623] - DEBUG - Building summary image for astrometry detected sources.
2024-09-24 16:04 - pid 4241 - [dipol.py:697] - DEBUG - Selected 5 quads with distance < 4.0. I will get the one with less deviation from the median linear transform.
2024-09-24 16:04 - pid 4241 - [dipol.py:700] - DEBUG - t_L=(array([100.10160072,  28.91945821]), array([100.27705537,  27.26207555]), array([99.88422651, 28.81341856]), array([99.953966  , 28.74211039]), array([100.15157897,  28.6252599 ]))
2024-09-24 16:04 - pid 4241 - [dipol.py:703] - DEBUG - Filtering out big translations (<1020px)
2024-09-24 16:04 - pid 4241 - [dipol.py:707] - DEBUG - Filtered to 5 quads with distance < 4.0 and translation < 1020px.
2024-09-24 16:04 - pid 4241 - [dipol.py:748] - DEBUG - median_t=array([100.10160072,  28.74211039])
2024-09-24 16:04 - pid 4241 - [dipol.py:754] - DEBUG - Selected the quads [20,21]
2024-09-24 16:04 - pid 4241 - [dipol.py:756] - DEBUG - t = [100.15157897  28.6252599 ]
2024-09-24 16:04 - pid 4241 - [dipol.py:760] - DEBUG - Building summary image for quad matching.
2024-09-24 16:04 - pid 4241 - [__init__.py:345] - WARNING - Less than 5 calibrated fits for DIPOL 2200+420, using all of them
2024-09-24 16:04 - pid 4241 - [dipol.py:805] - DEBUG - Using angle=181.1811506754202 for pre wcs.
2024-09-24 16:04 - pid 4241 - [dipol.py:820] - DEBUG - Building summary image for astrometry.
2024-09-24 16:04 - pid 4241 - [plotting.py:214] - DEBUG - <ReducedFit 7 | OSN-T090/2023-11-06/BLLAC_R_IAR-0760.fts>: found 3 catalog sources in field: [AstroSource.objects.get(name='2200+420'), AstroSource.objects.get(name='BL Lacertae C'), AstroSource.objects.get(name='BL Lacertae test 4')]
2024-09-24 16:04 - pid 4241 - [connectionpool.py:291] - DEBUG - Resetting dropped connection: simbad.cds.unistra.fr
2024-09-24 16:04 - pid 4241 - [connectionpool.py:546] - DEBUG - https://simbad.cds.unistra.fr:443 "POST /simbad/sim-script HTTP/11" 200 None
2024-09-24 16:04 - pid 4241 - [instrument.py:425] - DEBUG - <ReducedFit 7 | OSN-T090/2023-11-06/BLLAC_R_IAR-0760.fts>: saving WCSs to FITS header.
2024-09-24 16:04 - pid 4241 - [instrument.py:497] - DEBUG - <ReducedFit 7 | OSN-T090/2023-11-06/BLLAC_R_IAR-0760.fts>: astrometric calibration was successful.
2024-09-24 16:04 - pid 4241 - [instrument.py:500] - DEBUG - <ReducedFit 7 | OSN-T090/2023-11-06/BLLAC_R_IAR-0760.fts>: searching for sources in field...
2024-09-24 16:04 - pid 4241 - [instrument.py:503] - DEBUG - <ReducedFit 7 | OSN-T090/2023-11-06/BLLAC_R_IAR-0760.fts>: found 3 sources in field.
!ls {redf_pol.filedpropdir}

from IPython.display import Image
Image(filename=f"{redf_pol.filedpropdir}/astrometry_matched_quads.png")
BLLAC_R_IAR-0760.fts		 astrometry_matched_quads.png
astrometry_detected_sources.png  astrometry_summary.png
astrometry_info
../_images/c51edfed91ca165d681a277b8b3ef8d7a643106e8ff12ae95b971ba079604951.png

The image shows the quads of stars that were matched in both images, their distance in the hash space and the distance in the image space.

This quad matching is necessary because the number of stars in the polarimetry field is not enough for the astrometry.net solver to work, plus the scale of the images is too small for the default index files. This quad matching needs at least 4 (ideally a few more) stars in the image, and a previously calibrated photometry field. If any of these conditions are not met, other methods must be used to calibrate the image (if there are only two bright sources at the right distance, can be assumed that they are the target star, and if there are between 3 and 4-6 stars, we can try to match to known sources in the field).

The process is implemented in the DIPOL._build_wcs_for_polarimetry_images_photo_quads method, which is called by the ReducedFit.astrometric_calibration method. We could also manually follow the process:

import numpy as np
import itertools
from iop4lib.instruments import DIPOL
from iop4lib.utils.astrometry import BuildWCSResult
from photutils.aperture import CircularAperture
from pathlib import Path
import matplotlib as mplt
import matplotlib.pyplot as plt

# get the subframe of the photometry field that corresponds to this polarimetry field, (approx)
x_start = redf_pol.rawfit.header['XORGSUBF']
y_start = redf_pol.rawfit.header['YORGSUBF']

x_end = x_start + redf_pol.rawfit.header['NAXIS1']
y_end = y_start + redf_pol.rawfit.header['NAXIS2']

idx = np.s_[y_start:y_end, x_start:x_end]

photdata_subframe = redf_phot.mdata[idx]

# find 10 brightest sources in each field

sets_L = list()

for redf, data in zip([redf_pol, redf_phot], [redf_pol.mdata, photdata_subframe]):

    positions = DIPOL._estimate_positions_from_segments(redf=redf, data=data, n_seg_threshold=1.5, npixels=32, centering=None, fwhm=1.0)
    positions = positions[:10]

    sets_L.append(positions)

fig = plt.figure(figsize=(12,6), dpi=iop4conf.mplt_default_dpi)
axs = fig.subplots(nrows=1, ncols=2)

for ax, data, positions in zip(axs, [redf_pol.mdata, photdata_subframe], sets_L):
    imshow_w_sources(data, pos1=positions, ax=ax)
    candidates_aps = CircularAperture(positions[:2], r=10.0)
    candidates_aps.plot(ax, color="b")
    for i, (x,y) in enumerate(positions):
        ax.text(x, y, f"{i}", color="orange", fontdict={"size":14, "weight":"bold"})#, verticalalignment="center", horizontalalignment="center") 
    ax.plot([data.shape[1]//2], [data.shape[0]//2], '+', color='y', markersize=10)
    
axs[0].set_title("Polarimetry field")
axs[1].set_title("Photometry field")
fig.suptitle("astrometry detected_sources")
plt.show()
../_images/11bb620c70fa69343dd5f8f2012bbac64709e0198b497f0c1e2e18abd19e99f9.png
# Build the quads for each field
quads_1 = np.array(list(itertools.combinations(sets_L[0], 4)))
quads_2 = np.array(list(itertools.combinations(sets_L[1], 4)))

from iop4lib.utils.quadmatching import hash_ish, distance, order, qorder_ish, find_linear_transformation
hash_func, qorder = hash_ish, qorder_ish

hashes_1 = np.array([hash_func(quad) for quad in quads_1])
hashes_2 = np.array([hash_func(quad) for quad in quads_2])

all_indices = np.array(list(itertools.product(range(len(quads_1)),range(len(quads_2)))))
all_distances = np.array([distance(hashes_1[i], hashes_2[j]) for i,j in all_indices])

idx = np.argsort(all_distances)
all_indices = all_indices[idx]
all_distances = all_distances[idx]

indices_selected = all_indices[:6]
colors = ["r", "g", "b", "c", "m", "y"]
figs, axs = zip(*[plt.subplots(figsize=(6,4), dpi=iop4conf.mplt_default_dpi) for _ in range(2)])

for (i,j), color in list(zip(indices_selected, colors)): 
    
    tij = find_linear_transformation(qorder(quads_1[i]), qorder(quads_2[j]))[1]

    for ax, fig, data, quad, positions in zip(axs, figs, [redf_pol.mdata, photdata_subframe], [quads_1[i], quads_2[j]], sets_L):
        imshow_w_sources(data, pos1=positions, ax=ax)
        x, y = np.array(order(quad)).T
        ax.fill(x, y, edgecolor='k', fill=True, facecolor=mplt.colors.to_rgba(color, alpha=0.2))

plt.show()
../_images/2cb4e4a71a4a76e76ad49578ce644e2f4cda814eb90a316f5097e8f8c39f979c.png ../_images/dd33cbac924fe6ed398c861c44a76b01a5a12f6f608365130760a87a653c994e.png