Skip to content

Latest commit

 

History

History
2163 lines (1786 loc) · 89.8 KB

plot_a_hst.rst

File metadata and controls

2163 lines (1786 loc) · 89.8 KB
.. only:: html

    .. note::
        :class: sphx-glr-download-link-note

        Click :ref:`here <sphx_glr_download_examples_plot_a_hst.py>`
        to download the full example code

.. rst-class:: sphx-glr-example-title

Hubble

Aligning HST images with JHAT.

An example HST Dataset is downloaded, and then a series of alignment methods are used. For more information on the key parameters used for alignment see :ref:`params:Useful Parameters`.

import sys,os,glob
from astropy.io import fits
from astropy.table import Table
from astropy.nddata import extract_array
from astropy.coordinates import SkyCoord
from astropy import wcs
from astropy.wcs.utils import skycoord_to_pixel
from astropy import units as u
import numpy as np
import matplotlib.pyplot as plt
from astroquery.mast import Observations
from astropy.visualization import (simple_norm,LinearStretch)

import jhat
from jhat import hst_photclass,st_wcs_align

Relative Alignment

Download some Data

For this example we download 2 HST DRZ images from MAST. They're the same filter and same field, just separated in time.

obs_table = Observations.query_criteria(obs_id='hst_16264_12_wfc3_ir_f110w_iebc12')
obs_table1 = obs_table[obs_table['filters']=='F110W']

obs_table = Observations.query_criteria(obs_id='hst_16264_15_wfc3_ir_f110w_iebc15')
obs_table2 = obs_table[obs_table['filters']=='F110W']

data_products_by_obs = Observations.get_product_list(obs_table1)
data_products_by_obs = data_products_by_obs[data_products_by_obs['calib_level']==3]
data_products_by_obs = data_products_by_obs[data_products_by_obs['productSubGroupDescription']=='DRZ'][0]
Observations.download_products(data_products_by_obs,extension='fits')

data_products_by_obs = Observations.get_product_list(obs_table2)
data_products_by_obs = data_products_by_obs[data_products_by_obs['calib_level']==3]
data_products_by_obs = data_products_by_obs[data_products_by_obs['productSubGroupDescription']=='DRZ'][0]
Observations.download_products(data_products_by_obs,extension='fits')
.. rst-class:: sphx-glr-script-out

 .. code-block:: none

    Downloading URL https://mast.stsci.edu/api/v0.1/Download/file?uri=mast:HST/product/hst_16264_12_wfc3_ir_f110w_iebc12_drz.fits to ./mastDownload/HST/hst_16264_12_wfc3_ir_f110w_iebc12/hst_16264_12_wfc3_ir_f110w_iebc12_drz.fits ... [Done]
    Downloading URL https://mast.stsci.edu/api/v0.1/Download/file?uri=mast:HST/product/hst_16264_15_wfc3_ir_f110w_iebc15_drz.fits to ./mastDownload/HST/hst_16264_15_wfc3_ir_f110w_iebc15/hst_16264_15_wfc3_ir_f110w_iebc15_drz.fits ... [Done]


Table length=1
Local PathStatusMessageURL
str95str8objectobject
./mastDownload/HST/hst_16264_15_wfc3_ir_f110w_iebc15/hst_16264_15_wfc3_ir_f110w_iebc15_drz.fitsCOMPLETENoneNone


Examine the Reference Image

files = glob.glob('mastDownload/HST/*/*drz.fits')
ref_image = files[0]
ref_fits = fits.open(ref_image)
ref_data = fits.open(ref_image)['SCI',1].data
norm1 = simple_norm(ref_data,stretch='log',min_cut=-1,max_cut=15)

plt.imshow(ref_data, origin='lower',
                      norm=norm1,cmap='gray')
plt.gca().tick_params(labelcolor='none',axis='both',color='none')
plt.show()
.. image-sg:: /examples/images/sphx_glr_plot_a_hst_001.png
   :alt: plot a hst
   :srcset: /examples/images/sphx_glr_plot_a_hst_001.png
   :class: sphx-glr-single-img





Zoom in to see the offset

Here add an artificial offset to the wcs, and then we see the same star in both images at the same ra/dec location, demonstrating a large offset between the images.

star_location = SkyCoord('21:29:40.5351','+0:04:42.697',unit=(u.hourangle,u.deg))
align_image = files[1]
align_fits = fits.open(align_image)
align_fits['SCI',1].header['CRPIX1']+=2
align_fits['SCI',1].header['CRPIX2']+=2
align_fits.writeto(align_image,overwrite=True)

align_data = fits.open(align_image)['SCI',1].data
ref_y,ref_x = skycoord_to_pixel(star_location,wcs.WCS(ref_fits['SCI',1],ref_fits))
align_y,align_x = skycoord_to_pixel(star_location,wcs.WCS(align_fits['SCI',1],align_fits))

ref_cutout = extract_array(ref_data,(11,11),(ref_x,ref_y))
align_cutout = extract_array(align_data,(11,11),(align_x,align_y))
norm1 = simple_norm(ref_cutout,stretch='log',min_cut=-1,max_cut=200)
norm2 = simple_norm(align_cutout,stretch='log',min_cut=-1,max_cut=200)
fig,axes = plt.subplots(1,2)
axes[0].imshow(ref_cutout, origin='lower',
                      norm=norm1,cmap='gray')
axes[1].imshow(align_cutout, origin='lower',
                      norm=norm2,cmap='gray')
axes[0].set_title('Reference')
axes[1].set_title('To Align')
axes[0].tick_params(labelcolor='none',axis='both',color='none')
axes[1].tick_params(labelcolor='none',axis='both',color='none')

plt.show()
.. image-sg:: /examples/images/sphx_glr_plot_a_hst_002.png
   :alt: Reference, To Align
   :srcset: /examples/images/sphx_glr_plot_a_hst_002.png
   :class: sphx-glr-single-img





Create a Photometric Catalog for Relative Alignment

We choose one of the images to be the reference image, and then create a catalog that we will use to align the other image.

hst_phot = hst_photclass(psf_fwhm=1.8,aperture_radius=5)
hst_phot.run_phot(imagename=ref_image,photfilename='auto',overwrite=True)
ref_catname = ref_image.replace('.fits','.phot.txt') # the default
refcat = Table.read(ref_catname,format='ascii')
print(refcat)
.. rst-class:: sphx-glr-script-out

 .. code-block:: none

    0 mastDownload/HST/hst_16264_15_wfc3_ir_f110w_iebc15/hst_16264_15_wfc3_ir_f110w_iebc15_drz.phot.txt
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/astropy/stats/sigma_clipping.py:411: AstropyUserWarning: Input data contains invalid values (NaNs or infs), which were automatically clipped.
      warnings.warn('Input data contains invalid values (NaNs or '
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/astropy/stats/sigma_clipping.py:411: AstropyUserWarning: Input data contains invalid values (NaNs or infs), which were automatically clipped.
      warnings.warn('Input data contains invalid values (NaNs or '
    /Users/jpierel/CodeBase/jhat/jhat/simple_jwst_phot.py:1942: RuntimeWarning: invalid value encountered in log10
      phot['mag'] = -2.5*np.log10(phot['aper_sum_bkgsub'])+ee_corr+zp
    /Users/jpierel/CodeBase/jhat/jhat/simple_jwst_phot.py:1945: RuntimeWarning: invalid value encountered in log10
      phot['magerr'] = 2.5 * np.log10(1.0 + (fluxerr/phot['aper_sum_bkgsub']))
    aper_sum_5.0px annulus_median_5.0px aper_bkg_5.0px ...   x_idl     y_idl
    -------------- -------------------- -------------- ... --------- ----------
         98.701908             1.222414      96.008154 ... 10.282901 -58.751291
        112.858933             1.215701      95.480899 ...  0.099716 -56.988332
        101.387021             1.221539      95.939432 ... 11.171638 -55.244084
        103.395491             1.223192       96.06924 ... 12.050499 -54.702786
         98.175518             1.216526      95.545731 ...   18.3728 -52.885932
        105.185721              1.21896      95.736903 ... 20.125268 -52.409859
         98.467777             1.222879      96.044673 ...  8.878883 -51.562419
        101.178176             1.223382      96.084188 ...  9.447841 -51.376254
        100.858931             1.222232      95.993843 ...  9.589212 -51.288441
         97.386204             1.211678      95.164984 ... -1.323612 -51.067408
               ...                  ...            ... ...       ...        ...
        106.700574             1.228211      96.463442 ... 38.894969  98.502383
        100.028718             1.227917      96.440353 ... 38.374139   99.52547
         97.812416             1.220746      95.877198 ... 26.977886  99.641609
         97.602401             1.223992      96.132097 ... 37.761718 100.233421
        105.374585             1.219728      95.797236 ... 33.155942 100.494187
         97.821326             1.222629      96.025086 ... 36.519167 101.346028
        102.604081             1.224514      96.173087 ... 22.996227 101.696774
         97.176703             1.223466      96.090817 ... 29.490227 101.683135
        106.491719             1.221887      95.966789 ... 28.696314 105.056795
        136.421682             1.214225      95.364994 ... 36.206768 105.813417
        105.467064             1.229289      96.548109 ... 30.779631 106.393258
    Length = 769 rows




Align the second image

The plots outputted here show the various steps used by jhat to determine the true matching sources in the image, and the subsequent correction needed for optimal alignment.

wcs_align = st_wcs_align()


wcs_align.run_all(align_image,
              telescope='hst',
              outsubdir='mastDownload',
          refcat_racol='ra',
          refcat_deccol='dec',
          refcat_magcol='mag',
          refcat_magerrcol='dmag',
          overwrite=True,
          d2d_max=.5,
          showplots=2,
          refcatname=ref_catname,
          histocut_order='dxdy',
              sharpness_lim=(0.3,0.9),
              roundness1_lim=(-0.7, 0.7),
              SNR_min= 3,
              dmag_max=1.0,
              objmag_lim =(14,24))
.. rst-class:: sphx-glr-horizontal


    *

      .. image-sg:: /examples/images/sphx_glr_plot_a_hst_003.png
         :alt: Initial cut: d2d_max=0.5, dmag_max=None, Nbright=None, delta_mag_lim=(None, None)
         :srcset: /examples/images/sphx_glr_plot_a_hst_003.png
         :class: sphx-glr-multi-img

    *

      .. image-sg:: /examples/images/sphx_glr_plot_a_hst_004.png
         :alt: dx, dx, dx, slope:-9.7656250000017e-05, 3-sigma cut: 198 out of 222 left mean = 2.108 px, stdev = 0.055 px
         :srcset: /examples/images/sphx_glr_plot_a_hst_004.png
         :class: sphx-glr-multi-img

    *

      .. image-sg:: /examples/images/sphx_glr_plot_a_hst_005.png
         :alt: dy, dy, dy, slope:0.00014648437499998213, 3-sigma cut: 189 out of 198 left mean = 1.941 px, stdev = 0.058 px
         :srcset: /examples/images/sphx_glr_plot_a_hst_005.png
         :class: sphx-glr-multi-img


.. rst-class:: sphx-glr-script-out

 .. code-block:: none

    Warning: Setting aperture radius to twice the psf_fwhm (4.000000)
    0 ./mastDownload/hst_16264_12_wfc3_ir_f110w_iebc12.phot.txt
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/astropy/stats/sigma_clipping.py:411: AstropyUserWarning: Input data contains invalid values (NaNs or infs), which were automatically clipped.
      warnings.warn('Input data contains invalid values (NaNs or '
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/astropy/stats/sigma_clipping.py:411: AstropyUserWarning: Input data contains invalid values (NaNs or infs), which were automatically clipped.
      warnings.warn('Input data contains invalid values (NaNs or '
    /Users/jpierel/CodeBase/jhat/jhat/simple_jwst_phot.py:1942: RuntimeWarning: invalid value encountered in log10
      phot['mag'] = -2.5*np.log10(phot['aper_sum_bkgsub'])+ee_corr+zp
    /Users/jpierel/CodeBase/jhat/jhat/simple_jwst_phot.py:1945: RuntimeWarning: invalid value encountered in log10
      phot['magerr'] = 2.5 * np.log10(1.0 + (fluxerr/phot['aper_sum_bkgsub']))
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.instrument.name the following error occurred:
    'WFC3' is not one of ['NIRCAM', 'NIRSPEC', 'MIRI', 'TFI', 'FGS', 'NIRISS', 'ANY', 'N/A']

    Failed validating 'enum' in schema:
        OrderedDict([('title', 'Instrument used to acquire the data'),
                     ('type', 'string'),
                     ('enum',
                      ['NIRCAM',
                       'NIRSPEC',
                       'MIRI',
                       'TFI',
                       'FGS',
                       'NIRISS',
                       'ANY',
                       'N/A']),
                     ('fits_keyword', 'INSTRUME'),
                     ('blend_table', True)])

    On instance:
        'WFC3'
      warnings.warn(errmsg, ValidationWarning)
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.instrument.detector the following error occurred:
    'IR' is not one of ['NRCA1', 'NRCA2', 'NRCA3', 'NRCA4', 'NRCALONG', 'NRCB1', 'NRCB2', 'NRCB3', 'NRCB4', 'NRCBLONG', 'NRS1', 'NRS2', 'ANY', 'MIRIMAGE', 'MIRIFULONG', 'MIRIFUSHORT', 'NIS', 'GUIDER1', 'GUIDER2', 'N/A', 'MULTIPLE']

    Failed validating 'enum' in schema:
        OrderedDict([('title', 'Name of detector used to acquire the data'),
                     ('type', 'string'),
                     ('enum',
                      ['NRCA1',
                       'NRCA2',
                       'NRCA3',
                       'NRCA4',
                       'NRCALONG',
                       'NRCB1',
                       'NRCB2',
                       'NRCB3',
                       'NRCB4',
                       'NRCBLONG',
                       'NRS1',
                       'NRS2',
                       'ANY',
                       'MIRIMAGE',
                       'MIRIFULONG',
                       'MIRIFUSHORT',
                       'NIS',
                       'GUIDER1',
                       'GUIDER2',
                       'N/A',
                       'MULTIPLE']),
                     ('fits_keyword', 'DETECTOR'),
                     ('blend_table', True),
                     ('blend_rule', 'multi')])

    On instance:
        'IR'
      warnings.warn(errmsg, ValidationWarning)
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.subarray.name the following error occurred:
    False is not of type 'string'

    Failed validating 'type' in schema:
        OrderedDict([('title', 'Subarray used'),
                     ('type', 'string'),
                     ('anyOf',
                      [{'enum': ['8X8',
                                 '32X32',
                                 '128X128',
                                 '2048X64',
                                 'SUB128CNTR',
                                 'SUB128DIAG',
                                 'SUB128LL',
                                 'SUB32CNTR',
                                 'SUB32DIAG',
                                 'SUB32LL',
                                 'SUB8CNTR',
                                 'SUB8DIAG',
                                 'SUB8LL',
                                 'SUBIDSTRIPCENTER',
                                 'SUBIDSTRIPLL',
                                 'SUBTUNE32CENTERG1',
                                 'SUBTUNE32CENTERG2',
                                 'SUBTUNE32LLG1',
                                 'SUBTUNE32LLG2']},
                       {'enum': ['BRIGHTSKY',
                                 'MASK1065',
                                 'MASK1140',
                                 'MASK1550',
                                 'MASKLYOT',
                                 'SLITLESSPRISM',
                                 'SUB128',
                                 'SUB256',
                                 'SUB64',
                                 'SUBPRISM']},
                       {'enum': ['FULLP',
                                 'MASK210R',
                                 'MASK335R',
                                 'MASK430R',
                                 'MASKLWB',
                                 'MASKSWB',
                                 'SUB160',
                                 'SUB160P',
                                 'SUB320',
                                 'SUB320A335R',
                                 'SUB320A430R',
                                 'SUB320ALWB',
                                 'SUB320B335R',
                                 'SUB320B430R',
                      ...
      warnings.warn(errmsg, ValidationWarning)
    *** Note: close plot to continue!
        slope  intercept     maxval  index  d_bestguess  fwhm  multimax
    -0.000098   0.071191 183.673678     28     2.022576   0.8     False
    Keeping 222 out of 222, skippin 0 because of null values in columns d_rot_tmp
    median: 2.110833
    75.000000 percentile cut: max residual for cut: 0.076543
    median: 2.112779
    i:00 mean:2.112779(0.002520) stdev:0.032375(0.001777) X2norm:1.00 Nchanged:0 Ngood:166 Nclip:56

    mean: 2.114075
    i:01 mean:2.114075(0.002887) stdev:0.038415(0.002036) X2norm:1.00 Nchanged:12 Ngood:178 Nclip:44

    mean: 2.113588
    i:02 mean:2.113588(0.003081) stdev:0.041561(0.002172) X2norm:1.00 Nchanged:5 Ngood:183 Nclip:39

    mean: 2.113580
    i:03 mean:2.113580(0.003189) stdev:0.043262(0.002249) X2norm:1.00 Nchanged:2 Ngood:185 Nclip:37

    mean: 2.110885
    i:04 mean:2.110885(0.003397) stdev:0.046577(0.002396) X2norm:1.00 Nchanged:4 Ngood:189 Nclip:33

    mean: 2.111600
    i:05 mean:2.111600(0.003454) stdev:0.047488(0.002436) X2norm:1.00 Nchanged:1 Ngood:190 Nclip:32

    mean: 2.112334
    i:06 mean:2.112334(0.003631) stdev:0.050314(0.002561) X2norm:1.00 Nchanged:3 Ngood:193 Nclip:29

    mean: 2.110830
    i:07 mean:2.110830(0.003748) stdev:0.052203(0.002643) X2norm:1.00 Nchanged:2 Ngood:195 Nclip:27

    mean: 2.109305
    i:08 mean:2.109305(0.003863) stdev:0.054084(0.002725) X2norm:1.00 Nchanged:2 Ngood:197 Nclip:25

    mean: 2.108492
    i:09 mean:2.108492(0.003929) stdev:0.055147(0.002771) X2norm:1.00 Nchanged:1 Ngood:198 Nclip:24

       slope  intercept    maxval  index  d_bestguess  fwhm  multimax
    0.000146   -0.10686 185.08144      5     1.964436   0.8     False
    Keeping 198 out of 198, skippin 0 because of null values in columns d_rot_tmp
    median: 1.946136
    75.000000 percentile cut: max residual for cut: 0.064640
    median: 1.947841
    i:00 mean:1.947841(0.002509) stdev:0.030420(0.001768) X2norm:1.00 Nchanged:0 Ngood:148 Nclip:50

    mean: 1.939984
    i:01 mean:1.939984(0.002941) stdev:0.037665(0.002073) X2norm:1.00 Nchanged:17 Ngood:165 Nclip:33

    mean: 1.938798
    i:02 mean:1.938798(0.003106) stdev:0.040264(0.002190) X2norm:1.00 Nchanged:4 Ngood:169 Nclip:29

    mean: 1.937518
    i:03 mean:1.937518(0.003415) stdev:0.045050(0.002408) X2norm:1.00 Nchanged:6 Ngood:175 Nclip:23

    mean: 1.938289
    i:04 mean:1.938289(0.003692) stdev:0.049402(0.002604) X2norm:1.00 Nchanged:5 Ngood:180 Nclip:18

    mean: 1.939838
    i:05 mean:1.939838(0.004048) stdev:0.055060(0.002855) X2norm:1.00 Nchanged:6 Ngood:186 Nclip:12

    mean: 1.939851
    i:06 mean:1.939851(0.004173) stdev:0.057063(0.002943) X2norm:1.00 Nchanged:2 Ngood:188 Nclip:10

    mean: 1.940734
    i:07 mean:1.940734(0.004244) stdev:0.058191(0.002993) X2norm:1.00 Nchanged:1 Ngood:189 Nclip:9

    mean: 1.940734
    i:08 mean:1.940734(0.004244) stdev:0.058191(0.002993) X2norm:1.00 Nchanged:0 Ngood:189 Nclip:9
    *** Note: close plots to continue!
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.instrument.name the following error occurred:
    'WFC3' is not one of ['NIRCAM', 'NIRSPEC', 'MIRI', 'TFI', 'FGS', 'NIRISS', 'ANY', 'N/A']

    Failed validating 'enum' in schema:
        OrderedDict([('title', 'Instrument used to acquire the data'),
                     ('type', 'string'),
                     ('enum',
                      ['NIRCAM',
                       'NIRSPEC',
                       'MIRI',
                       'TFI',
                       'FGS',
                       'NIRISS',
                       'ANY',
                       'N/A']),
                     ('fits_keyword', 'INSTRUME'),
                     ('blend_table', True)])

    On instance:
        'WFC3'
      warnings.warn(errmsg, ValidationWarning)
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.instrument.detector the following error occurred:
    'IR' is not one of ['NRCA1', 'NRCA2', 'NRCA3', 'NRCA4', 'NRCALONG', 'NRCB1', 'NRCB2', 'NRCB3', 'NRCB4', 'NRCBLONG', 'NRS1', 'NRS2', 'ANY', 'MIRIMAGE', 'MIRIFULONG', 'MIRIFUSHORT', 'NIS', 'GUIDER1', 'GUIDER2', 'N/A', 'MULTIPLE']

    Failed validating 'enum' in schema:
        OrderedDict([('title', 'Name of detector used to acquire the data'),
                     ('type', 'string'),
                     ('enum',
                      ['NRCA1',
                       'NRCA2',
                       'NRCA3',
                       'NRCA4',
                       'NRCALONG',
                       'NRCB1',
                       'NRCB2',
                       'NRCB3',
                       'NRCB4',
                       'NRCBLONG',
                       'NRS1',
                       'NRS2',
                       'ANY',
                       'MIRIMAGE',
                       'MIRIFULONG',
                       'MIRIFUSHORT',
                       'NIS',
                       'GUIDER1',
                       'GUIDER2',
                       'N/A',
                       'MULTIPLE']),
                     ('fits_keyword', 'DETECTOR'),
                     ('blend_table', True),
                     ('blend_rule', 'multi')])

    On instance:
        'IR'
      warnings.warn(errmsg, ValidationWarning)
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.subarray.name the following error occurred:
    False is not of type 'string'

    Failed validating 'type' in schema:
        OrderedDict([('title', 'Subarray used'),
                     ('type', 'string'),
                     ('anyOf',
                      [{'enum': ['8X8',
                                 '32X32',
                                 '128X128',
                                 '2048X64',
                                 'SUB128CNTR',
                                 'SUB128DIAG',
                                 'SUB128LL',
                                 'SUB32CNTR',
                                 'SUB32DIAG',
                                 'SUB32LL',
                                 'SUB8CNTR',
                                 'SUB8DIAG',
                                 'SUB8LL',
                                 'SUBIDSTRIPCENTER',
                                 'SUBIDSTRIPLL',
                                 'SUBTUNE32CENTERG1',
                                 'SUBTUNE32CENTERG2',
                                 'SUBTUNE32LLG1',
                                 'SUBTUNE32LLG2']},
                       {'enum': ['BRIGHTSKY',
                                 'MASK1065',
                                 'MASK1140',
                                 'MASK1550',
                                 'MASKLYOT',
                                 'SLITLESSPRISM',
                                 'SUB128',
                                 'SUB256',
                                 'SUB64',
                                 'SUBPRISM']},
                       {'enum': ['FULLP',
                                 'MASK210R',
                                 'MASK335R',
                                 'MASK430R',
                                 'MASKLWB',
                                 'MASKSWB',
                                 'SUB160',
                                 'SUB160P',
                                 'SUB320',
                                 'SUB320A335R',
                                 'SUB320A430R',
                                 'SUB320ALWB',
                                 'SUB320B335R',
                                 'SUB320B430R',
                      ...
      warnings.warn(errmsg, ValidationWarning)
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/jwst/datamodels/util.py:234: NoTypeWarning: model_type not found. Opening mastDownload/HST/hst_16264_12_wfc3_ir_f110w_iebc12/hst_16264_12_wfc3_ir_f110w_iebc12_drz.fits as a ImageModel
      warnings.warn(f"model_type not found. Opening {file_name} as a {class_name}",
    *** Note: close plots to continue!

    0



Check the Output

The reference image has not changed, but let's read in the newly aligned image and compare with the original. subsequent correction needed for optimal alignment.

aligned_image = os.path.join('mastDownload',os.path.basename(align_image).replace('drz.fits','jhat.fits'))
aligned_fits = fits.open(aligned_image)
aligned_data = fits.open(aligned_image)['SCI',1].data
aligned_y,aligned_x = skycoord_to_pixel(star_location,wcs.WCS(aligned_fits['SCI',1],aligned_fits))
aligned_cutout = extract_array(aligned_data,(11,11),(aligned_x,aligned_y))

norm3 = simple_norm(aligned_cutout,stretch='log',min_cut=-1,max_cut=200)
fig,axes = plt.subplots(1,3)
axes[0].imshow(ref_cutout, origin='lower',
                      norm=norm1,cmap='gray')
axes[1].imshow(align_cutout, origin='lower',
                      norm=norm2,cmap='gray')
axes[2].imshow(aligned_cutout, origin='lower',
                      norm=norm3,cmap='gray')
axes[0].set_title('Reference')
axes[1].set_title('To Align')
axes[2].set_title('Aligned')
for i in range(3):
    axes[i].tick_params(labelcolor='none',axis='both',color='none')


plt.show()
.. image-sg:: /examples/images/sphx_glr_plot_a_hst_006.png
   :alt: Reference, To Align, Aligned
   :srcset: /examples/images/sphx_glr_plot_a_hst_006.png
   :class: sphx-glr-single-img





Align to Gaia

You can also align each image to the Gaia DR3 catalog, or you could replace the catalog created in step one with your own catalog of the field.

wcs_align.run_all(align_image,
              telescope='hst',
              outsubdir='mastDownload',
          overwrite=True,
          d2d_max=.5,
          showplots=0,
          refcatname='Gaia',
          histocut_order='dxdy',
              sharpness_lim=(0.3,0.9),
              roundness1_lim=(-0.7, 0.7),
              SNR_min= 3,
              dmag_max=1.0,
              objmag_lim =(14,24))

aligned_image = os.path.join('mastDownload',os.path.basename(align_image).replace('drz.fits','jhat.fits'))
aligned_fits = fits.open(aligned_image)
aligned_data = fits.open(aligned_image)['SCI',1].data
aligned_y,aligned_x = skycoord_to_pixel(star_location,wcs.WCS(aligned_fits['SCI',1],aligned_fits))
aligned_cutout = extract_array(aligned_data,(11,11),(aligned_x,aligned_y))

norm3 = simple_norm(aligned_cutout,stretch='log',min_cut=-1,max_cut=200)
fig,axes = plt.subplots(1,2)
axes[0].imshow(align_cutout, origin='lower',
                      norm=norm2,cmap='gray')
axes[1].imshow(aligned_cutout, origin='lower',
                      norm=norm3,cmap='gray')
axes[0].set_title('To Align')
axes[1].set_title('Aligned')
for i in range(2):
    axes[i].tick_params(labelcolor='none',axis='both',color='none')


plt.show()
.. image-sg:: /examples/images/sphx_glr_plot_a_hst_007.png
   :alt: To Align, Aligned
   :srcset: /examples/images/sphx_glr_plot_a_hst_007.png
   :class: sphx-glr-single-img


.. rst-class:: sphx-glr-script-out

 .. code-block:: none

    Warning: Setting aperture radius to twice the psf_fwhm (4.000000)
    0 ./mastDownload/hst_16264_12_wfc3_ir_f110w_iebc12.phot.txt
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/astropy/stats/sigma_clipping.py:411: AstropyUserWarning: Input data contains invalid values (NaNs or infs), which were automatically clipped.
      warnings.warn('Input data contains invalid values (NaNs or '
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/astropy/stats/sigma_clipping.py:411: AstropyUserWarning: Input data contains invalid values (NaNs or infs), which were automatically clipped.
      warnings.warn('Input data contains invalid values (NaNs or '
    /Users/jpierel/CodeBase/jhat/jhat/simple_jwst_phot.py:1942: RuntimeWarning: invalid value encountered in log10
      phot['mag'] = -2.5*np.log10(phot['aper_sum_bkgsub'])+ee_corr+zp
    /Users/jpierel/CodeBase/jhat/jhat/simple_jwst_phot.py:1945: RuntimeWarning: invalid value encountered in log10
      phot['magerr'] = 2.5 * np.log10(1.0 + (fluxerr/phot['aper_sum_bkgsub']))
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.instrument.name the following error occurred:
    'WFC3' is not one of ['NIRCAM', 'NIRSPEC', 'MIRI', 'TFI', 'FGS', 'NIRISS', 'ANY', 'N/A']

    Failed validating 'enum' in schema:
        OrderedDict([('title', 'Instrument used to acquire the data'),
                     ('type', 'string'),
                     ('enum',
                      ['NIRCAM',
                       'NIRSPEC',
                       'MIRI',
                       'TFI',
                       'FGS',
                       'NIRISS',
                       'ANY',
                       'N/A']),
                     ('fits_keyword', 'INSTRUME'),
                     ('blend_table', True)])

    On instance:
        'WFC3'
      warnings.warn(errmsg, ValidationWarning)
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.instrument.detector the following error occurred:
    'IR' is not one of ['NRCA1', 'NRCA2', 'NRCA3', 'NRCA4', 'NRCALONG', 'NRCB1', 'NRCB2', 'NRCB3', 'NRCB4', 'NRCBLONG', 'NRS1', 'NRS2', 'ANY', 'MIRIMAGE', 'MIRIFULONG', 'MIRIFUSHORT', 'NIS', 'GUIDER1', 'GUIDER2', 'N/A', 'MULTIPLE']

    Failed validating 'enum' in schema:
        OrderedDict([('title', 'Name of detector used to acquire the data'),
                     ('type', 'string'),
                     ('enum',
                      ['NRCA1',
                       'NRCA2',
                       'NRCA3',
                       'NRCA4',
                       'NRCALONG',
                       'NRCB1',
                       'NRCB2',
                       'NRCB3',
                       'NRCB4',
                       'NRCBLONG',
                       'NRS1',
                       'NRS2',
                       'ANY',
                       'MIRIMAGE',
                       'MIRIFULONG',
                       'MIRIFUSHORT',
                       'NIS',
                       'GUIDER1',
                       'GUIDER2',
                       'N/A',
                       'MULTIPLE']),
                     ('fits_keyword', 'DETECTOR'),
                     ('blend_table', True),
                     ('blend_rule', 'multi')])

    On instance:
        'IR'
      warnings.warn(errmsg, ValidationWarning)
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.subarray.name the following error occurred:
    False is not of type 'string'

    Failed validating 'type' in schema:
        OrderedDict([('title', 'Subarray used'),
                     ('type', 'string'),
                     ('anyOf',
                      [{'enum': ['8X8',
                                 '32X32',
                                 '128X128',
                                 '2048X64',
                                 'SUB128CNTR',
                                 'SUB128DIAG',
                                 'SUB128LL',
                                 'SUB32CNTR',
                                 'SUB32DIAG',
                                 'SUB32LL',
                                 'SUB8CNTR',
                                 'SUB8DIAG',
                                 'SUB8LL',
                                 'SUBIDSTRIPCENTER',
                                 'SUBIDSTRIPLL',
                                 'SUBTUNE32CENTERG1',
                                 'SUBTUNE32CENTERG2',
                                 'SUBTUNE32LLG1',
                                 'SUBTUNE32LLG2']},
                       {'enum': ['BRIGHTSKY',
                                 'MASK1065',
                                 'MASK1140',
                                 'MASK1550',
                                 'MASKLYOT',
                                 'SLITLESSPRISM',
                                 'SUB128',
                                 'SUB256',
                                 'SUB64',
                                 'SUBPRISM']},
                       {'enum': ['FULLP',
                                 'MASK210R',
                                 'MASK335R',
                                 'MASK430R',
                                 'MASKLWB',
                                 'MASKSWB',
                                 'SUB160',
                                 'SUB160P',
                                 'SUB320',
                                 'SUB320A335R',
                                 'SUB320A430R',
                                 'SUB320ALWB',
                                 'SUB320B335R',
                                 'SUB320B430R',
                      ...
      warnings.warn(errmsg, ValidationWarning)
    INFO: Query finished. [astroquery.utils.tap.core]
    Number of stars: 81
    ### NO propoer motion correction!!!
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/pandas/core/arraylike.py:402: RuntimeWarning: invalid value encountered in sqrt
      result = getattr(ufunc, method)(*inputs, **kwargs)
    Number of stars after removing nan's: 81
       slope  intercept   maxval  index  d_bestguess  fwhm  multimax
    0.000537  -0.391553 6.837554      6     1.974222   0.8     False
    Keeping 9 out of 9, skippin 0 because of null values in columns d_rot_tmp
    median: 1.981504
    75.000000 percentile cut: max residual for cut: 0.180525
    median: 2.017531
    i:00 mean:2.017531(0.051121) stdev:0.114310(0.032998) X2norm:0.91 Nchanged:0 Ngood:6 Nclip:3

    mean: 2.001886
    i:01 mean:2.001886(0.055165) stdev:0.145953(0.036488) X2norm:1.00 Nchanged:2 Ngood:8 Nclip:1

    mean: 2.001886
    i:02 mean:2.001886(0.055165) stdev:0.145953(0.036488) X2norm:1.00 Nchanged:0 Ngood:8 Nclip:1
       slope  intercept  maxval  index  d_bestguess  fwhm  multimax
    0.000146   -0.10686   5.236      6     2.360891   1.0     False
    Keeping 8 out of 8, skippin 0 because of null values in columns d_rot_tmp
    median: 2.342233
    75.000000 percentile cut: max residual for cut: 0.358199
    median: 2.342233
    i:00 mean:2.342233(0.100460) stdev:0.224635(0.064847) X2norm:0.91 Nchanged:0 Ngood:6 Nclip:2

    mean: 2.331019
    i:01 mean:2.331019(0.114384) stdev:0.302631(0.075658) X2norm:1.00 Nchanged:2 Ngood:8 Nclip:0

    mean: 2.331019
    i:02 mean:2.331019(0.114384) stdev:0.302631(0.075658) X2norm:1.00 Nchanged:0 Ngood:8 Nclip:0
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.instrument.name the following error occurred:
    'WFC3' is not one of ['NIRCAM', 'NIRSPEC', 'MIRI', 'TFI', 'FGS', 'NIRISS', 'ANY', 'N/A']

    Failed validating 'enum' in schema:
        OrderedDict([('title', 'Instrument used to acquire the data'),
                     ('type', 'string'),
                     ('enum',
                      ['NIRCAM',
                       'NIRSPEC',
                       'MIRI',
                       'TFI',
                       'FGS',
                       'NIRISS',
                       'ANY',
                       'N/A']),
                     ('fits_keyword', 'INSTRUME'),
                     ('blend_table', True)])

    On instance:
        'WFC3'
      warnings.warn(errmsg, ValidationWarning)
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.instrument.detector the following error occurred:
    'IR' is not one of ['NRCA1', 'NRCA2', 'NRCA3', 'NRCA4', 'NRCALONG', 'NRCB1', 'NRCB2', 'NRCB3', 'NRCB4', 'NRCBLONG', 'NRS1', 'NRS2', 'ANY', 'MIRIMAGE', 'MIRIFULONG', 'MIRIFUSHORT', 'NIS', 'GUIDER1', 'GUIDER2', 'N/A', 'MULTIPLE']

    Failed validating 'enum' in schema:
        OrderedDict([('title', 'Name of detector used to acquire the data'),
                     ('type', 'string'),
                     ('enum',
                      ['NRCA1',
                       'NRCA2',
                       'NRCA3',
                       'NRCA4',
                       'NRCALONG',
                       'NRCB1',
                       'NRCB2',
                       'NRCB3',
                       'NRCB4',
                       'NRCBLONG',
                       'NRS1',
                       'NRS2',
                       'ANY',
                       'MIRIMAGE',
                       'MIRIFULONG',
                       'MIRIFUSHORT',
                       'NIS',
                       'GUIDER1',
                       'GUIDER2',
                       'N/A',
                       'MULTIPLE']),
                     ('fits_keyword', 'DETECTOR'),
                     ('blend_table', True),
                     ('blend_rule', 'multi')])

    On instance:
        'IR'
      warnings.warn(errmsg, ValidationWarning)
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.subarray.name the following error occurred:
    False is not of type 'string'

    Failed validating 'type' in schema:
        OrderedDict([('title', 'Subarray used'),
                     ('type', 'string'),
                     ('anyOf',
                      [{'enum': ['8X8',
                                 '32X32',
                                 '128X128',
                                 '2048X64',
                                 'SUB128CNTR',
                                 'SUB128DIAG',
                                 'SUB128LL',
                                 'SUB32CNTR',
                                 'SUB32DIAG',
                                 'SUB32LL',
                                 'SUB8CNTR',
                                 'SUB8DIAG',
                                 'SUB8LL',
                                 'SUBIDSTRIPCENTER',
                                 'SUBIDSTRIPLL',
                                 'SUBTUNE32CENTERG1',
                                 'SUBTUNE32CENTERG2',
                                 'SUBTUNE32LLG1',
                                 'SUBTUNE32LLG2']},
                       {'enum': ['BRIGHTSKY',
                                 'MASK1065',
                                 'MASK1140',
                                 'MASK1550',
                                 'MASKLYOT',
                                 'SLITLESSPRISM',
                                 'SUB128',
                                 'SUB256',
                                 'SUB64',
                                 'SUBPRISM']},
                       {'enum': ['FULLP',
                                 'MASK210R',
                                 'MASK335R',
                                 'MASK430R',
                                 'MASKLWB',
                                 'MASKSWB',
                                 'SUB160',
                                 'SUB160P',
                                 'SUB320',
                                 'SUB320A335R',
                                 'SUB320A430R',
                                 'SUB320ALWB',
                                 'SUB320B335R',
                                 'SUB320B430R',
                      ...
      warnings.warn(errmsg, ValidationWarning)
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/jwst/datamodels/util.py:234: NoTypeWarning: model_type not found. Opening mastDownload/HST/hst_16264_12_wfc3_ir_f110w_iebc12/hst_16264_12_wfc3_ir_f110w_iebc12_drz.fits as a ImageModel
      warnings.warn(f"model_type not found. Opening {file_name} as a {class_name}",




Large Offsets

Sometimes the initial images are so poorly aligned, that the code fails. Here we read in the same image as in the first example, and add an additional 3 pixel offset in the wcs.

files = glob.glob('mastDownload/HST/*/*drz.fits')
align_image = files[1]
align_fits = fits.open(align_image)
align_fits['SCI',1].header['CRPIX1']+=3
align_fits['SCI',1].header['CRPIX2']+=3
align_fits.writeto(align_image,overwrite=True)

align_data = fits.open(align_image)['SCI',1].data
ref_y,ref_x = skycoord_to_pixel(star_location,wcs.WCS(ref_fits['SCI',1],ref_fits))
align_y,align_x = skycoord_to_pixel(star_location,wcs.WCS(align_fits['SCI',1],align_fits))

ref_cutout = extract_array(ref_data,(11,11),(ref_x,ref_y))
align_cutout = extract_array(align_data,(11,11),(align_x,align_y))
norm1 = simple_norm(ref_cutout,stretch='log',min_cut=-1,max_cut=200)
norm2 = simple_norm(align_cutout,stretch='log',min_cut=-1,max_cut=200)
fig,axes = plt.subplots(1,2)
axes[0].imshow(ref_cutout, origin='lower',
                      norm=norm1,cmap='gray')
axes[1].imshow(align_cutout, origin='lower',
                      norm=norm2,cmap='gray')
axes[0].set_title('Reference')
axes[1].set_title('To Align')
axes[0].tick_params(labelcolor='none',axis='both',color='none')
axes[1].tick_params(labelcolor='none',axis='both',color='none')

plt.show()

wcs_align = st_wcs_align()

try:
    wcs_align.run_all(align_image,
              telescope='hst',
              outsubdir='mastDownload',
          refcat_racol='ra',
          refcat_deccol='dec',
          refcat_magcol='mag',
          refcat_magerrcol='dmag',
          overwrite=True,
          d2d_max=.5,
          showplots=2,
          refcatname=ref_catname,
          histocut_order='dxdy',
              sharpness_lim=(0.3,0.9),
              roundness1_lim=(-0.7, 0.7),
              SNR_min= 3,
              dmag_max=1.0,
              objmag_lim =(14,24))

except:
    print('Failed for not enough matches!')
.. rst-class:: sphx-glr-horizontal


    *

      .. image-sg:: /examples/images/sphx_glr_plot_a_hst_008.png
         :alt: Reference, To Align
         :srcset: /examples/images/sphx_glr_plot_a_hst_008.png
         :class: sphx-glr-multi-img

    *

      .. image-sg:: /examples/images/sphx_glr_plot_a_hst_009.png
         :alt: Initial cut: d2d_max=0.5, dmag_max=None, Nbright=None, delta_mag_lim=(None, None)
         :srcset: /examples/images/sphx_glr_plot_a_hst_009.png
         :class: sphx-glr-multi-img

    *

      .. image-sg:: /examples/images/sphx_glr_plot_a_hst_010.png
         :alt: dx, dx, dx, slope:0.0020507812499999754, 3-sigma cut: 9 out of 13 left mean = 1.129 px, stdev = 0.094 px
         :srcset: /examples/images/sphx_glr_plot_a_hst_010.png
         :class: sphx-glr-multi-img

    *

      .. image-sg:: /examples/images/sphx_glr_plot_a_hst_011.png
         :alt: dy, dy, dy, slope:-0.0011230468750000134, 3-sigma cut: 3 out of 4 left mean = 1.080 px, stdev = 0.078 px
         :srcset: /examples/images/sphx_glr_plot_a_hst_011.png
         :class: sphx-glr-multi-img


.. rst-class:: sphx-glr-script-out

 .. code-block:: none

    Warning: Setting aperture radius to twice the psf_fwhm (4.000000)
    0 ./mastDownload/hst_16264_12_wfc3_ir_f110w_iebc12.phot.txt
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/astropy/stats/sigma_clipping.py:411: AstropyUserWarning: Input data contains invalid values (NaNs or infs), which were automatically clipped.
      warnings.warn('Input data contains invalid values (NaNs or '
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/astropy/stats/sigma_clipping.py:411: AstropyUserWarning: Input data contains invalid values (NaNs or infs), which were automatically clipped.
      warnings.warn('Input data contains invalid values (NaNs or '
    /Users/jpierel/CodeBase/jhat/jhat/simple_jwst_phot.py:1942: RuntimeWarning: invalid value encountered in log10
      phot['mag'] = -2.5*np.log10(phot['aper_sum_bkgsub'])+ee_corr+zp
    /Users/jpierel/CodeBase/jhat/jhat/simple_jwst_phot.py:1945: RuntimeWarning: invalid value encountered in log10
      phot['magerr'] = 2.5 * np.log10(1.0 + (fluxerr/phot['aper_sum_bkgsub']))
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.instrument.name the following error occurred:
    'WFC3' is not one of ['NIRCAM', 'NIRSPEC', 'MIRI', 'TFI', 'FGS', 'NIRISS', 'ANY', 'N/A']

    Failed validating 'enum' in schema:
        OrderedDict([('title', 'Instrument used to acquire the data'),
                     ('type', 'string'),
                     ('enum',
                      ['NIRCAM',
                       'NIRSPEC',
                       'MIRI',
                       'TFI',
                       'FGS',
                       'NIRISS',
                       'ANY',
                       'N/A']),
                     ('fits_keyword', 'INSTRUME'),
                     ('blend_table', True)])

    On instance:
        'WFC3'
      warnings.warn(errmsg, ValidationWarning)
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.instrument.detector the following error occurred:
    'IR' is not one of ['NRCA1', 'NRCA2', 'NRCA3', 'NRCA4', 'NRCALONG', 'NRCB1', 'NRCB2', 'NRCB3', 'NRCB4', 'NRCBLONG', 'NRS1', 'NRS2', 'ANY', 'MIRIMAGE', 'MIRIFULONG', 'MIRIFUSHORT', 'NIS', 'GUIDER1', 'GUIDER2', 'N/A', 'MULTIPLE']

    Failed validating 'enum' in schema:
        OrderedDict([('title', 'Name of detector used to acquire the data'),
                     ('type', 'string'),
                     ('enum',
                      ['NRCA1',
                       'NRCA2',
                       'NRCA3',
                       'NRCA4',
                       'NRCALONG',
                       'NRCB1',
                       'NRCB2',
                       'NRCB3',
                       'NRCB4',
                       'NRCBLONG',
                       'NRS1',
                       'NRS2',
                       'ANY',
                       'MIRIMAGE',
                       'MIRIFULONG',
                       'MIRIFUSHORT',
                       'NIS',
                       'GUIDER1',
                       'GUIDER2',
                       'N/A',
                       'MULTIPLE']),
                     ('fits_keyword', 'DETECTOR'),
                     ('blend_table', True),
                     ('blend_rule', 'multi')])

    On instance:
        'IR'
      warnings.warn(errmsg, ValidationWarning)
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.subarray.name the following error occurred:
    False is not of type 'string'

    Failed validating 'type' in schema:
        OrderedDict([('title', 'Subarray used'),
                     ('type', 'string'),
                     ('anyOf',
                      [{'enum': ['8X8',
                                 '32X32',
                                 '128X128',
                                 '2048X64',
                                 'SUB128CNTR',
                                 'SUB128DIAG',
                                 'SUB128LL',
                                 'SUB32CNTR',
                                 'SUB32DIAG',
                                 'SUB32LL',
                                 'SUB8CNTR',
                                 'SUB8DIAG',
                                 'SUB8LL',
                                 'SUBIDSTRIPCENTER',
                                 'SUBIDSTRIPLL',
                                 'SUBTUNE32CENTERG1',
                                 'SUBTUNE32CENTERG2',
                                 'SUBTUNE32LLG1',
                                 'SUBTUNE32LLG2']},
                       {'enum': ['BRIGHTSKY',
                                 'MASK1065',
                                 'MASK1140',
                                 'MASK1550',
                                 'MASKLYOT',
                                 'SLITLESSPRISM',
                                 'SUB128',
                                 'SUB256',
                                 'SUB64',
                                 'SUBPRISM']},
                       {'enum': ['FULLP',
                                 'MASK210R',
                                 'MASK335R',
                                 'MASK430R',
                                 'MASKLWB',
                                 'MASKSWB',
                                 'SUB160',
                                 'SUB160P',
                                 'SUB320',
                                 'SUB320A335R',
                                 'SUB320A430R',
                                 'SUB320ALWB',
                                 'SUB320B335R',
                                 'SUB320B430R',
                      ...
      warnings.warn(errmsg, ValidationWarning)
    *** Note: close plot to continue!
       slope  intercept   maxval  index  d_bestguess  fwhm  multimax
    0.002051   -1.49502 8.877268     34     1.149103   0.8     False
    Keeping 13 out of 13, skippin 0 because of null values in columns d_rot_tmp
    median: 1.142187
    75.000000 percentile cut: max residual for cut: 0.473790
    median: 1.142187
    i:00 mean:1.142187(0.034755) stdev:0.098302(0.023170) X2norm:0.94 Nchanged:0 Ngood:9 Nclip:4

    mean: 1.128891
    i:01 mean:1.128891(0.033318) stdev:0.094236(0.022212) X2norm:1.00 Nchanged:0 Ngood:9 Nclip:4
        slope  intercept   maxval  index  d_bestguess  fwhm  multimax
    -0.001123   0.819263 3.661515     19     1.083046   0.8     False
    Keeping 4 out of 4, skippin 0 because of null values in columns d_rot_tmp
    median: 1.045332
    75.000000 percentile cut: max residual for cut: 0.138967
    median: 1.095192
    i:00 mean:1.095192(0.064008) stdev:0.090521(0.036955) X2norm:0.79 Nchanged:0 Ngood:3 Nclip:1

    mean: 1.079986
    i:01 mean:1.079986(0.055176) stdev:0.078030(0.031856) X2norm:1.00 Nchanged:0 Ngood:3 Nclip:1
    *** Note: close plots to continue!
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.instrument.name the following error occurred:
    'WFC3' is not one of ['NIRCAM', 'NIRSPEC', 'MIRI', 'TFI', 'FGS', 'NIRISS', 'ANY', 'N/A']

    Failed validating 'enum' in schema:
        OrderedDict([('title', 'Instrument used to acquire the data'),
                     ('type', 'string'),
                     ('enum',
                      ['NIRCAM',
                       'NIRSPEC',
                       'MIRI',
                       'TFI',
                       'FGS',
                       'NIRISS',
                       'ANY',
                       'N/A']),
                     ('fits_keyword', 'INSTRUME'),
                     ('blend_table', True)])

    On instance:
        'WFC3'
      warnings.warn(errmsg, ValidationWarning)
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.instrument.detector the following error occurred:
    'IR' is not one of ['NRCA1', 'NRCA2', 'NRCA3', 'NRCA4', 'NRCALONG', 'NRCB1', 'NRCB2', 'NRCB3', 'NRCB4', 'NRCBLONG', 'NRS1', 'NRS2', 'ANY', 'MIRIMAGE', 'MIRIFULONG', 'MIRIFUSHORT', 'NIS', 'GUIDER1', 'GUIDER2', 'N/A', 'MULTIPLE']

    Failed validating 'enum' in schema:
        OrderedDict([('title', 'Name of detector used to acquire the data'),
                     ('type', 'string'),
                     ('enum',
                      ['NRCA1',
                       'NRCA2',
                       'NRCA3',
                       'NRCA4',
                       'NRCALONG',
                       'NRCB1',
                       'NRCB2',
                       'NRCB3',
                       'NRCB4',
                       'NRCBLONG',
                       'NRS1',
                       'NRS2',
                       'ANY',
                       'MIRIMAGE',
                       'MIRIFULONG',
                       'MIRIFUSHORT',
                       'NIS',
                       'GUIDER1',
                       'GUIDER2',
                       'N/A',
                       'MULTIPLE']),
                     ('fits_keyword', 'DETECTOR'),
                     ('blend_table', True),
                     ('blend_rule', 'multi')])

    On instance:
        'IR'
      warnings.warn(errmsg, ValidationWarning)
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.subarray.name the following error occurred:
    False is not of type 'string'

    Failed validating 'type' in schema:
        OrderedDict([('title', 'Subarray used'),
                     ('type', 'string'),
                     ('anyOf',
                      [{'enum': ['8X8',
                                 '32X32',
                                 '128X128',
                                 '2048X64',
                                 'SUB128CNTR',
                                 'SUB128DIAG',
                                 'SUB128LL',
                                 'SUB32CNTR',
                                 'SUB32DIAG',
                                 'SUB32LL',
                                 'SUB8CNTR',
                                 'SUB8DIAG',
                                 'SUB8LL',
                                 'SUBIDSTRIPCENTER',
                                 'SUBIDSTRIPLL',
                                 'SUBTUNE32CENTERG1',
                                 'SUBTUNE32CENTERG2',
                                 'SUBTUNE32LLG1',
                                 'SUBTUNE32LLG2']},
                       {'enum': ['BRIGHTSKY',
                                 'MASK1065',
                                 'MASK1140',
                                 'MASK1550',
                                 'MASKLYOT',
                                 'SLITLESSPRISM',
                                 'SUB128',
                                 'SUB256',
                                 'SUB64',
                                 'SUBPRISM']},
                       {'enum': ['FULLP',
                                 'MASK210R',
                                 'MASK335R',
                                 'MASK430R',
                                 'MASKLWB',
                                 'MASKSWB',
                                 'SUB160',
                                 'SUB160P',
                                 'SUB320',
                                 'SUB320A335R',
                                 'SUB320A430R',
                                 'SUB320ALWB',
                                 'SUB320B335R',
                                 'SUB320B430R',
                      ...
      warnings.warn(errmsg, ValidationWarning)
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/jwst/datamodels/util.py:234: NoTypeWarning: model_type not found. Opening mastDownload/HST/hst_16264_12_wfc3_ir_f110w_iebc12/hst_16264_12_wfc3_ir_f110w_iebc12_drz.fits as a ImageModel
      warnings.warn(f"model_type not found. Opening {file_name} as a {class_name}",
    Failed for not enough matches!




This is what a failure looks like (compare to the plots above). There are now a couple of options here. You can increase the d2d_max parameter, which increases the allowed distance between sources being matched in the reference and target images:

wcs_align = st_wcs_align()


wcs_align.run_all(align_image,
              telescope='hst',
              outsubdir='mastDownload',
          refcat_racol='ra',
          refcat_deccol='dec',
          refcat_magcol='mag',
          refcat_magerrcol='dmag',
          overwrite=True,
          d2d_max=1,
          showplots=2,
          refcatname=ref_catname,
          histocut_order='dxdy',
              sharpness_lim=(0.3,0.9),
              roundness1_lim=(-0.7, 0.7),
              SNR_min= 3,
              dmag_max=1.0,
              objmag_lim =(14,24))

aligned_image = os.path.join('mastDownload',os.path.basename(align_image).replace('drz.fits','jhat.fits'))
aligned_fits = fits.open(aligned_image)
aligned_data = fits.open(aligned_image)['SCI',1].data
aligned_y,aligned_x = skycoord_to_pixel(star_location,wcs.WCS(aligned_fits['SCI',1],aligned_fits))
aligned_cutout = extract_array(aligned_data,(11,11),(aligned_x,aligned_y))

norm3 = simple_norm(aligned_cutout,stretch='log',min_cut=-1,max_cut=200)
fig,axes = plt.subplots(1,3)
axes[0].imshow(ref_cutout, origin='lower',
                      norm=norm1,cmap='gray')
axes[1].imshow(align_cutout, origin='lower',
                      norm=norm2,cmap='gray')
axes[2].imshow(aligned_cutout, origin='lower',
                      norm=norm3,cmap='gray')
axes[0].set_title('Reference')
axes[1].set_title('To Align')
axes[2].set_title('Aligned')
for i in range(3):
    axes[i].tick_params(labelcolor='none',axis='both',color='none')


plt.show()
.. rst-class:: sphx-glr-horizontal


    *

      .. image-sg:: /examples/images/sphx_glr_plot_a_hst_012.png
         :alt: Initial cut: d2d_max=1, dmag_max=None, Nbright=None, delta_mag_lim=(None, None)
         :srcset: /examples/images/sphx_glr_plot_a_hst_012.png
         :class: sphx-glr-multi-img

    *

      .. image-sg:: /examples/images/sphx_glr_plot_a_hst_013.png
         :alt: dx, dx, dx, slope:-0.0002929687500000163, 3-sigma cut: 163 out of 177 left mean = 5.103 px, stdev = 0.076 px
         :srcset: /examples/images/sphx_glr_plot_a_hst_013.png
         :class: sphx-glr-multi-img

    *

      .. image-sg:: /examples/images/sphx_glr_plot_a_hst_014.png
         :alt: dy, dy, dy, slope:9.76562499999823e-05, 3-sigma cut: 152 out of 162 left mean = 4.944 px, stdev = 0.053 px
         :srcset: /examples/images/sphx_glr_plot_a_hst_014.png
         :class: sphx-glr-multi-img

    *

      .. image-sg:: /examples/images/sphx_glr_plot_a_hst_015.png
         :alt: Reference, To Align, Aligned
         :srcset: /examples/images/sphx_glr_plot_a_hst_015.png
         :class: sphx-glr-multi-img


.. rst-class:: sphx-glr-script-out

 .. code-block:: none

    Warning: Setting aperture radius to twice the psf_fwhm (4.000000)
    0 ./mastDownload/hst_16264_12_wfc3_ir_f110w_iebc12.phot.txt
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/astropy/stats/sigma_clipping.py:411: AstropyUserWarning: Input data contains invalid values (NaNs or infs), which were automatically clipped.
      warnings.warn('Input data contains invalid values (NaNs or '
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/astropy/stats/sigma_clipping.py:411: AstropyUserWarning: Input data contains invalid values (NaNs or infs), which were automatically clipped.
      warnings.warn('Input data contains invalid values (NaNs or '
    /Users/jpierel/CodeBase/jhat/jhat/simple_jwst_phot.py:1942: RuntimeWarning: invalid value encountered in log10
      phot['mag'] = -2.5*np.log10(phot['aper_sum_bkgsub'])+ee_corr+zp
    /Users/jpierel/CodeBase/jhat/jhat/simple_jwst_phot.py:1945: RuntimeWarning: invalid value encountered in log10
      phot['magerr'] = 2.5 * np.log10(1.0 + (fluxerr/phot['aper_sum_bkgsub']))
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.instrument.name the following error occurred:
    'WFC3' is not one of ['NIRCAM', 'NIRSPEC', 'MIRI', 'TFI', 'FGS', 'NIRISS', 'ANY', 'N/A']

    Failed validating 'enum' in schema:
        OrderedDict([('title', 'Instrument used to acquire the data'),
                     ('type', 'string'),
                     ('enum',
                      ['NIRCAM',
                       'NIRSPEC',
                       'MIRI',
                       'TFI',
                       'FGS',
                       'NIRISS',
                       'ANY',
                       'N/A']),
                     ('fits_keyword', 'INSTRUME'),
                     ('blend_table', True)])

    On instance:
        'WFC3'
      warnings.warn(errmsg, ValidationWarning)
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.instrument.detector the following error occurred:
    'IR' is not one of ['NRCA1', 'NRCA2', 'NRCA3', 'NRCA4', 'NRCALONG', 'NRCB1', 'NRCB2', 'NRCB3', 'NRCB4', 'NRCBLONG', 'NRS1', 'NRS2', 'ANY', 'MIRIMAGE', 'MIRIFULONG', 'MIRIFUSHORT', 'NIS', 'GUIDER1', 'GUIDER2', 'N/A', 'MULTIPLE']

    Failed validating 'enum' in schema:
        OrderedDict([('title', 'Name of detector used to acquire the data'),
                     ('type', 'string'),
                     ('enum',
                      ['NRCA1',
                       'NRCA2',
                       'NRCA3',
                       'NRCA4',
                       'NRCALONG',
                       'NRCB1',
                       'NRCB2',
                       'NRCB3',
                       'NRCB4',
                       'NRCBLONG',
                       'NRS1',
                       'NRS2',
                       'ANY',
                       'MIRIMAGE',
                       'MIRIFULONG',
                       'MIRIFUSHORT',
                       'NIS',
                       'GUIDER1',
                       'GUIDER2',
                       'N/A',
                       'MULTIPLE']),
                     ('fits_keyword', 'DETECTOR'),
                     ('blend_table', True),
                     ('blend_rule', 'multi')])

    On instance:
        'IR'
      warnings.warn(errmsg, ValidationWarning)
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.subarray.name the following error occurred:
    False is not of type 'string'

    Failed validating 'type' in schema:
        OrderedDict([('title', 'Subarray used'),
                     ('type', 'string'),
                     ('anyOf',
                      [{'enum': ['8X8',
                                 '32X32',
                                 '128X128',
                                 '2048X64',
                                 'SUB128CNTR',
                                 'SUB128DIAG',
                                 'SUB128LL',
                                 'SUB32CNTR',
                                 'SUB32DIAG',
                                 'SUB32LL',
                                 'SUB8CNTR',
                                 'SUB8DIAG',
                                 'SUB8LL',
                                 'SUBIDSTRIPCENTER',
                                 'SUBIDSTRIPLL',
                                 'SUBTUNE32CENTERG1',
                                 'SUBTUNE32CENTERG2',
                                 'SUBTUNE32LLG1',
                                 'SUBTUNE32LLG2']},
                       {'enum': ['BRIGHTSKY',
                                 'MASK1065',
                                 'MASK1140',
                                 'MASK1550',
                                 'MASKLYOT',
                                 'SLITLESSPRISM',
                                 'SUB128',
                                 'SUB256',
                                 'SUB64',
                                 'SUBPRISM']},
                       {'enum': ['FULLP',
                                 'MASK210R',
                                 'MASK335R',
                                 'MASK430R',
                                 'MASKLWB',
                                 'MASKSWB',
                                 'SUB160',
                                 'SUB160P',
                                 'SUB320',
                                 'SUB320A335R',
                                 'SUB320A430R',
                                 'SUB320ALWB',
                                 'SUB320B335R',
                                 'SUB320B430R',
                      ...
      warnings.warn(errmsg, ValidationWarning)
    *** Note: close plot to continue!
        slope  intercept     maxval  index  d_bestguess  fwhm  multimax
    -0.000293   0.213574 151.767497     58     5.058877   0.8     False
    Keeping 177 out of 177, skippin 0 because of null values in columns d_rot_tmp
    median: 5.100640
    75.000000 percentile cut: max residual for cut: 0.091924
    median: 5.091393
    i:00 mean:5.091393(0.004565) stdev:0.052249(0.003216) X2norm:1.00 Nchanged:0 Ngood:132 Nclip:45

    mean: 5.101872
    i:01 mean:5.101872(0.005427) stdev:0.067565(0.003825) X2norm:1.00 Nchanged:24 Ngood:156 Nclip:21

    mean: 5.103992
    i:02 mean:5.103992(0.005907) stdev:0.074955(0.004164) X2norm:1.00 Nchanged:6 Ngood:162 Nclip:15

    mean: 5.102721
    i:03 mean:5.102721(0.006008) stdev:0.076467(0.004235) X2norm:1.00 Nchanged:1 Ngood:163 Nclip:14

    mean: 5.102721
    i:04 mean:5.102721(0.006008) stdev:0.076467(0.004235) X2norm:1.00 Nchanged:0 Ngood:163 Nclip:14
       slope  intercept     maxval  index  d_bestguess  fwhm  multimax
    0.000098   -0.07124 153.147275     54     4.948694   0.8     False
    Keeping 162 out of 162, skippin 0 because of null values in columns d_rot_tmp
    median: 4.942136
    75.000000 percentile cut: max residual for cut: 0.062217
    median: 4.946718
    i:00 mean:4.946718(0.002966) stdev:0.032492(0.002089) X2norm:1.00 Nchanged:0 Ngood:121 Nclip:41

    mean: 4.942254
    i:01 mean:4.942254(0.003583) stdev:0.042396(0.002525) X2norm:1.00 Nchanged:20 Ngood:141 Nclip:21

    mean: 4.941777
    i:02 mean:4.941777(0.003869) stdev:0.046592(0.002727) X2norm:1.00 Nchanged:5 Ngood:146 Nclip:16

    mean: 4.942611
    i:03 mean:4.942611(0.004236) stdev:0.051878(0.002985) X2norm:1.00 Nchanged:5 Ngood:151 Nclip:11

    mean: 4.943574
    i:04 mean:4.943574(0.004317) stdev:0.053053(0.003043) X2norm:1.00 Nchanged:1 Ngood:152 Nclip:10

    mean: 4.943574
    i:05 mean:4.943574(0.004317) stdev:0.053053(0.003043) X2norm:1.00 Nchanged:0 Ngood:152 Nclip:10
    *** Note: close plots to continue!
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.instrument.name the following error occurred:
    'WFC3' is not one of ['NIRCAM', 'NIRSPEC', 'MIRI', 'TFI', 'FGS', 'NIRISS', 'ANY', 'N/A']

    Failed validating 'enum' in schema:
        OrderedDict([('title', 'Instrument used to acquire the data'),
                     ('type', 'string'),
                     ('enum',
                      ['NIRCAM',
                       'NIRSPEC',
                       'MIRI',
                       'TFI',
                       'FGS',
                       'NIRISS',
                       'ANY',
                       'N/A']),
                     ('fits_keyword', 'INSTRUME'),
                     ('blend_table', True)])

    On instance:
        'WFC3'
      warnings.warn(errmsg, ValidationWarning)
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.instrument.detector the following error occurred:
    'IR' is not one of ['NRCA1', 'NRCA2', 'NRCA3', 'NRCA4', 'NRCALONG', 'NRCB1', 'NRCB2', 'NRCB3', 'NRCB4', 'NRCBLONG', 'NRS1', 'NRS2', 'ANY', 'MIRIMAGE', 'MIRIFULONG', 'MIRIFUSHORT', 'NIS', 'GUIDER1', 'GUIDER2', 'N/A', 'MULTIPLE']

    Failed validating 'enum' in schema:
        OrderedDict([('title', 'Name of detector used to acquire the data'),
                     ('type', 'string'),
                     ('enum',
                      ['NRCA1',
                       'NRCA2',
                       'NRCA3',
                       'NRCA4',
                       'NRCALONG',
                       'NRCB1',
                       'NRCB2',
                       'NRCB3',
                       'NRCB4',
                       'NRCBLONG',
                       'NRS1',
                       'NRS2',
                       'ANY',
                       'MIRIMAGE',
                       'MIRIFULONG',
                       'MIRIFUSHORT',
                       'NIS',
                       'GUIDER1',
                       'GUIDER2',
                       'N/A',
                       'MULTIPLE']),
                     ('fits_keyword', 'DETECTOR'),
                     ('blend_table', True),
                     ('blend_rule', 'multi')])

    On instance:
        'IR'
      warnings.warn(errmsg, ValidationWarning)
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.subarray.name the following error occurred:
    False is not of type 'string'

    Failed validating 'type' in schema:
        OrderedDict([('title', 'Subarray used'),
                     ('type', 'string'),
                     ('anyOf',
                      [{'enum': ['8X8',
                                 '32X32',
                                 '128X128',
                                 '2048X64',
                                 'SUB128CNTR',
                                 'SUB128DIAG',
                                 'SUB128LL',
                                 'SUB32CNTR',
                                 'SUB32DIAG',
                                 'SUB32LL',
                                 'SUB8CNTR',
                                 'SUB8DIAG',
                                 'SUB8LL',
                                 'SUBIDSTRIPCENTER',
                                 'SUBIDSTRIPLL',
                                 'SUBTUNE32CENTERG1',
                                 'SUBTUNE32CENTERG2',
                                 'SUBTUNE32LLG1',
                                 'SUBTUNE32LLG2']},
                       {'enum': ['BRIGHTSKY',
                                 'MASK1065',
                                 'MASK1140',
                                 'MASK1550',
                                 'MASKLYOT',
                                 'SLITLESSPRISM',
                                 'SUB128',
                                 'SUB256',
                                 'SUB64',
                                 'SUBPRISM']},
                       {'enum': ['FULLP',
                                 'MASK210R',
                                 'MASK335R',
                                 'MASK430R',
                                 'MASKLWB',
                                 'MASKSWB',
                                 'SUB160',
                                 'SUB160P',
                                 'SUB320',
                                 'SUB320A335R',
                                 'SUB320A430R',
                                 'SUB320ALWB',
                                 'SUB320B335R',
                                 'SUB320B430R',
                      ...
      warnings.warn(errmsg, ValidationWarning)
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/jwst/datamodels/util.py:234: NoTypeWarning: model_type not found. Opening mastDownload/HST/hst_16264_12_wfc3_ir_f110w_iebc12/hst_16264_12_wfc3_ir_f110w_iebc12_drz.fits as a ImageModel
      warnings.warn(f"model_type not found. Opening {file_name} as a {class_name}",
    *** Note: close plots to continue!




Or you can apply a rough guess for the offset, and then use a smaller d2d_max for matching:

wcs_align = st_wcs_align()


wcs_align.run_all(align_image,
              telescope='hst',
              outsubdir='mastDownload',
          refcat_racol='ra',
          refcat_deccol='dec',
          refcat_magcol='mag',
          refcat_magerrcol='dmag',
          overwrite=True,
          d2d_max=.25,
          xshift=5,
          yshift=5,
          showplots=2,
          refcatname=ref_catname,
          histocut_order='dxdy',
              sharpness_lim=(0.3,0.9),
              roundness1_lim=(-0.7, 0.7),
              SNR_min= 3,
              dmag_max=1.0,
              objmag_lim =(14,24))

aligned_image = os.path.join('mastDownload',os.path.basename(align_image).replace('drz.fits','jhat.fits'))
aligned_fits = fits.open(aligned_image)
aligned_data = fits.open(aligned_image)['SCI',1].data
aligned_y,aligned_x = skycoord_to_pixel(star_location,wcs.WCS(aligned_fits['SCI',1],aligned_fits))
aligned_cutout = extract_array(aligned_data,(11,11),(aligned_x,aligned_y))

norm3 = simple_norm(aligned_cutout,stretch='log',min_cut=-1,max_cut=200)
fig,axes = plt.subplots(1,3)
axes[0].imshow(ref_cutout, origin='lower',
                      norm=norm1,cmap='gray')
axes[1].imshow(align_cutout, origin='lower',
                      norm=norm2,cmap='gray')
axes[2].imshow(aligned_cutout, origin='lower',
                      norm=norm3,cmap='gray')
axes[0].set_title('Reference')
axes[1].set_title('To Align')
axes[2].set_title('Aligned')
for i in range(3):
    axes[i].tick_params(labelcolor='none',axis='both',color='none')


plt.show()
.. rst-class:: sphx-glr-horizontal


    *

      .. image-sg:: /examples/images/sphx_glr_plot_a_hst_016.png
         :alt: Initial cut: d2d_max=0.25, dmag_max=None, Nbright=None, delta_mag_lim=(None, None)
         :srcset: /examples/images/sphx_glr_plot_a_hst_016.png
         :class: sphx-glr-multi-img

    *

      .. image-sg:: /examples/images/sphx_glr_plot_a_hst_017.png
         :alt: dx, dx, dx, slope:-4.8828125000017174e-05, 3-sigma cut: 203 out of 228 left mean = 5.109 px, stdev = 0.055 px
         :srcset: /examples/images/sphx_glr_plot_a_hst_017.png
         :class: sphx-glr-multi-img

    *

      .. image-sg:: /examples/images/sphx_glr_plot_a_hst_018.png
         :alt: dy, dy, dy, slope:0.00014648437499998213, 3-sigma cut: 191 out of 203 left mean = 4.941 px, stdev = 0.058 px
         :srcset: /examples/images/sphx_glr_plot_a_hst_018.png
         :class: sphx-glr-multi-img

    *

      .. image-sg:: /examples/images/sphx_glr_plot_a_hst_019.png
         :alt: Reference, To Align, Aligned
         :srcset: /examples/images/sphx_glr_plot_a_hst_019.png
         :class: sphx-glr-multi-img


.. rst-class:: sphx-glr-script-out

 .. code-block:: none

    Warning: Setting aperture radius to twice the psf_fwhm (4.000000)
    0 ./mastDownload/hst_16264_12_wfc3_ir_f110w_iebc12.phot.txt
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/astropy/stats/sigma_clipping.py:411: AstropyUserWarning: Input data contains invalid values (NaNs or infs), which were automatically clipped.
      warnings.warn('Input data contains invalid values (NaNs or '
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/astropy/stats/sigma_clipping.py:411: AstropyUserWarning: Input data contains invalid values (NaNs or infs), which were automatically clipped.
      warnings.warn('Input data contains invalid values (NaNs or '
    /Users/jpierel/CodeBase/jhat/jhat/simple_jwst_phot.py:1942: RuntimeWarning: invalid value encountered in log10
      phot['mag'] = -2.5*np.log10(phot['aper_sum_bkgsub'])+ee_corr+zp
    /Users/jpierel/CodeBase/jhat/jhat/simple_jwst_phot.py:1945: RuntimeWarning: invalid value encountered in log10
      phot['magerr'] = 2.5 * np.log10(1.0 + (fluxerr/phot['aper_sum_bkgsub']))
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.instrument.name the following error occurred:
    'WFC3' is not one of ['NIRCAM', 'NIRSPEC', 'MIRI', 'TFI', 'FGS', 'NIRISS', 'ANY', 'N/A']

    Failed validating 'enum' in schema:
        OrderedDict([('title', 'Instrument used to acquire the data'),
                     ('type', 'string'),
                     ('enum',
                      ['NIRCAM',
                       'NIRSPEC',
                       'MIRI',
                       'TFI',
                       'FGS',
                       'NIRISS',
                       'ANY',
                       'N/A']),
                     ('fits_keyword', 'INSTRUME'),
                     ('blend_table', True)])

    On instance:
        'WFC3'
      warnings.warn(errmsg, ValidationWarning)
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.instrument.detector the following error occurred:
    'IR' is not one of ['NRCA1', 'NRCA2', 'NRCA3', 'NRCA4', 'NRCALONG', 'NRCB1', 'NRCB2', 'NRCB3', 'NRCB4', 'NRCBLONG', 'NRS1', 'NRS2', 'ANY', 'MIRIMAGE', 'MIRIFULONG', 'MIRIFUSHORT', 'NIS', 'GUIDER1', 'GUIDER2', 'N/A', 'MULTIPLE']

    Failed validating 'enum' in schema:
        OrderedDict([('title', 'Name of detector used to acquire the data'),
                     ('type', 'string'),
                     ('enum',
                      ['NRCA1',
                       'NRCA2',
                       'NRCA3',
                       'NRCA4',
                       'NRCALONG',
                       'NRCB1',
                       'NRCB2',
                       'NRCB3',
                       'NRCB4',
                       'NRCBLONG',
                       'NRS1',
                       'NRS2',
                       'ANY',
                       'MIRIMAGE',
                       'MIRIFULONG',
                       'MIRIFUSHORT',
                       'NIS',
                       'GUIDER1',
                       'GUIDER2',
                       'N/A',
                       'MULTIPLE']),
                     ('fits_keyword', 'DETECTOR'),
                     ('blend_table', True),
                     ('blend_rule', 'multi')])

    On instance:
        'IR'
      warnings.warn(errmsg, ValidationWarning)
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.subarray.name the following error occurred:
    False is not of type 'string'

    Failed validating 'type' in schema:
        OrderedDict([('title', 'Subarray used'),
                     ('type', 'string'),
                     ('anyOf',
                      [{'enum': ['8X8',
                                 '32X32',
                                 '128X128',
                                 '2048X64',
                                 'SUB128CNTR',
                                 'SUB128DIAG',
                                 'SUB128LL',
                                 'SUB32CNTR',
                                 'SUB32DIAG',
                                 'SUB32LL',
                                 'SUB8CNTR',
                                 'SUB8DIAG',
                                 'SUB8LL',
                                 'SUBIDSTRIPCENTER',
                                 'SUBIDSTRIPLL',
                                 'SUBTUNE32CENTERG1',
                                 'SUBTUNE32CENTERG2',
                                 'SUBTUNE32LLG1',
                                 'SUBTUNE32LLG2']},
                       {'enum': ['BRIGHTSKY',
                                 'MASK1065',
                                 'MASK1140',
                                 'MASK1550',
                                 'MASKLYOT',
                                 'SLITLESSPRISM',
                                 'SUB128',
                                 'SUB256',
                                 'SUB64',
                                 'SUBPRISM']},
                       {'enum': ['FULLP',
                                 'MASK210R',
                                 'MASK335R',
                                 'MASK430R',
                                 'MASKLWB',
                                 'MASKSWB',
                                 'SUB160',
                                 'SUB160P',
                                 'SUB320',
                                 'SUB320A335R',
                                 'SUB320A430R',
                                 'SUB320ALWB',
                                 'SUB320B335R',
                                 'SUB320B430R',
                      ...
      warnings.warn(errmsg, ValidationWarning)
    *** Note: close plot to continue!
        slope  intercept     maxval  index  d_bestguess  fwhm  multimax
    -0.000049   0.035596 210.533193     11     5.118499   0.8     False
    Keeping 228 out of 228, skippin 0 because of null values in columns d_rot_tmp
    median: 5.108273
    75.000000 percentile cut: max residual for cut: 0.073336
    median: 5.111246
    i:00 mean:5.111246(0.002647) stdev:0.034510(0.001866) X2norm:1.00 Nchanged:0 Ngood:171 Nclip:57

    mean: 5.113886
    i:01 mean:5.113886(0.003071) stdev:0.041878(0.002165) X2norm:1.00 Nchanged:16 Ngood:187 Nclip:41

    mean: 5.113807
    i:02 mean:5.113807(0.003391) stdev:0.047237(0.002392) X2norm:1.00 Nchanged:8 Ngood:195 Nclip:33

    mean: 5.113819
    i:03 mean:5.113819(0.003501) stdev:0.049019(0.002470) X2norm:1.00 Nchanged:2 Ngood:197 Nclip:31

    mean: 5.113098
    i:04 mean:5.113098(0.003558) stdev:0.049936(0.002509) X2norm:1.00 Nchanged:1 Ngood:198 Nclip:30

    mean: 5.111626
    i:05 mean:5.111626(0.003673) stdev:0.051807(0.002590) X2norm:1.00 Nchanged:2 Ngood:200 Nclip:28

    mean: 5.110125
    i:06 mean:5.110125(0.003788) stdev:0.053702(0.002672) X2norm:1.00 Nchanged:2 Ngood:202 Nclip:26

    mean: 5.109359
    i:07 mean:5.109359(0.003846) stdev:0.054669(0.002713) X2norm:1.00 Nchanged:1 Ngood:203 Nclip:25

    mean: 5.109359
    i:08 mean:5.109359(0.003846) stdev:0.054669(0.002713) X2norm:1.00 Nchanged:0 Ngood:203 Nclip:25
       slope  intercept     maxval  index  d_bestguess  fwhm  multimax
    0.000146   -0.10686 187.297193      5     4.964436   0.8     False
    Keeping 203 out of 203, skippin 0 because of null values in columns d_rot_tmp
    median: 4.946334
    75.000000 percentile cut: max residual for cut: 0.066828
    median: 4.947336
    i:00 mean:4.947336(0.002530) stdev:0.031086(0.001783) X2norm:1.00 Nchanged:0 Ngood:152 Nclip:51

    mean: 4.940657
    i:01 mean:4.940657(0.002958) stdev:0.038228(0.002086) X2norm:1.00 Nchanged:16 Ngood:168 Nclip:35

    mean: 4.938915
    i:02 mean:4.938915(0.003185) stdev:0.041767(0.002245) X2norm:1.00 Nchanged:5 Ngood:173 Nclip:30

    mean: 4.937633
    i:03 mean:4.937633(0.003487) stdev:0.046521(0.002459) X2norm:1.00 Nchanged:6 Ngood:179 Nclip:24

    mean: 4.939127
    i:04 mean:4.939127(0.003714) stdev:0.050106(0.002619) X2norm:1.00 Nchanged:4 Ngood:183 Nclip:20

    mean: 4.939913
    i:05 mean:4.939913(0.004011) stdev:0.054843(0.002828) X2norm:1.00 Nchanged:5 Ngood:188 Nclip:15

    mean: 4.939926
    i:06 mean:4.939926(0.004134) stdev:0.056835(0.002916) X2norm:1.00 Nchanged:2 Ngood:190 Nclip:13

    mean: 4.940799
    i:07 mean:4.940799(0.004205) stdev:0.057955(0.002965) X2norm:1.00 Nchanged:1 Ngood:191 Nclip:12

    mean: 4.940799
    i:08 mean:4.940799(0.004205) stdev:0.057955(0.002965) X2norm:1.00 Nchanged:0 Ngood:191 Nclip:12
    *** Note: close plots to continue!
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.instrument.name the following error occurred:
    'WFC3' is not one of ['NIRCAM', 'NIRSPEC', 'MIRI', 'TFI', 'FGS', 'NIRISS', 'ANY', 'N/A']

    Failed validating 'enum' in schema:
        OrderedDict([('title', 'Instrument used to acquire the data'),
                     ('type', 'string'),
                     ('enum',
                      ['NIRCAM',
                       'NIRSPEC',
                       'MIRI',
                       'TFI',
                       'FGS',
                       'NIRISS',
                       'ANY',
                       'N/A']),
                     ('fits_keyword', 'INSTRUME'),
                     ('blend_table', True)])

    On instance:
        'WFC3'
      warnings.warn(errmsg, ValidationWarning)
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.instrument.detector the following error occurred:
    'IR' is not one of ['NRCA1', 'NRCA2', 'NRCA3', 'NRCA4', 'NRCALONG', 'NRCB1', 'NRCB2', 'NRCB3', 'NRCB4', 'NRCBLONG', 'NRS1', 'NRS2', 'ANY', 'MIRIMAGE', 'MIRIFULONG', 'MIRIFUSHORT', 'NIS', 'GUIDER1', 'GUIDER2', 'N/A', 'MULTIPLE']

    Failed validating 'enum' in schema:
        OrderedDict([('title', 'Name of detector used to acquire the data'),
                     ('type', 'string'),
                     ('enum',
                      ['NRCA1',
                       'NRCA2',
                       'NRCA3',
                       'NRCA4',
                       'NRCALONG',
                       'NRCB1',
                       'NRCB2',
                       'NRCB3',
                       'NRCB4',
                       'NRCBLONG',
                       'NRS1',
                       'NRS2',
                       'ANY',
                       'MIRIMAGE',
                       'MIRIFULONG',
                       'MIRIFUSHORT',
                       'NIS',
                       'GUIDER1',
                       'GUIDER2',
                       'N/A',
                       'MULTIPLE']),
                     ('fits_keyword', 'DETECTOR'),
                     ('blend_table', True),
                     ('blend_rule', 'multi')])

    On instance:
        'IR'
      warnings.warn(errmsg, ValidationWarning)
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/stdatamodels/validate.py:38: ValidationWarning: While validating meta.subarray.name the following error occurred:
    False is not of type 'string'

    Failed validating 'type' in schema:
        OrderedDict([('title', 'Subarray used'),
                     ('type', 'string'),
                     ('anyOf',
                      [{'enum': ['8X8',
                                 '32X32',
                                 '128X128',
                                 '2048X64',
                                 'SUB128CNTR',
                                 'SUB128DIAG',
                                 'SUB128LL',
                                 'SUB32CNTR',
                                 'SUB32DIAG',
                                 'SUB32LL',
                                 'SUB8CNTR',
                                 'SUB8DIAG',
                                 'SUB8LL',
                                 'SUBIDSTRIPCENTER',
                                 'SUBIDSTRIPLL',
                                 'SUBTUNE32CENTERG1',
                                 'SUBTUNE32CENTERG2',
                                 'SUBTUNE32LLG1',
                                 'SUBTUNE32LLG2']},
                       {'enum': ['BRIGHTSKY',
                                 'MASK1065',
                                 'MASK1140',
                                 'MASK1550',
                                 'MASKLYOT',
                                 'SLITLESSPRISM',
                                 'SUB128',
                                 'SUB256',
                                 'SUB64',
                                 'SUBPRISM']},
                       {'enum': ['FULLP',
                                 'MASK210R',
                                 'MASK335R',
                                 'MASK430R',
                                 'MASKLWB',
                                 'MASKSWB',
                                 'SUB160',
                                 'SUB160P',
                                 'SUB320',
                                 'SUB320A335R',
                                 'SUB320A430R',
                                 'SUB320ALWB',
                                 'SUB320B335R',
                                 'SUB320B430R',
                      ...
      warnings.warn(errmsg, ValidationWarning)
    /Users/jpierel/miniconda3/envs/tweakreg/lib/python3.10/site-packages/jwst/datamodels/util.py:234: NoTypeWarning: model_type not found. Opening mastDownload/HST/hst_16264_12_wfc3_ir_f110w_iebc12/hst_16264_12_wfc3_ir_f110w_iebc12_drz.fits as a ImageModel
      warnings.warn(f"model_type not found. Opening {file_name} as a {class_name}",
    *** Note: close plots to continue!





.. rst-class:: sphx-glr-timing

   **Total running time of the script:** ( 0 minutes  57.951 seconds)


.. only:: html

  .. container:: sphx-glr-footer sphx-glr-footer-example


    .. container:: sphx-glr-download sphx-glr-download-python

      :download:`Download Python source code: plot_a_hst.py <plot_a_hst.py>`

    .. container:: sphx-glr-download sphx-glr-download-jupyter

      :download:`Download Jupyter notebook: plot_a_hst.ipynb <plot_a_hst.ipynb>`


.. only:: html

 .. rst-class:: sphx-glr-signature

    `Gallery generated by Sphinx-Gallery <https://sphinx-gallery.github.io>`_