Ejemplo n.º 1
0
def exer5(odir):
    """ Create ripples of phase aberration in the pupil
        to simulate the Ruffles Potato Chip Telescope (RPCT) """
    # instantiate an mft object:
    ft = matrixDFT.MatrixFourierTransform()
    npup = 100
    radius = 20.0
    fp_size_reselt = 100
    pupil = utils.makedisk(npup, radius=radius)
    number_of_waves_across = (2,3,4,5,6)

    peaks = (0.0, 0.1, 0.3, 1.0, 3.0) # radians, amplitude of phase ripple
    arrayshape = (npup, npup)
    diam_pupil = radius*2 # pixels
    center = (arrayshape[0]/2, arrayshape[1]/2)
    for nwaves in number_of_waves_across:
        for peak in peaks: 
            for offset in (0, ):
                for angle in (0, ):
                    spatialwavelen = diam_pupil / nwaves
                    offset = offset * np.pi/180.0
                    khat = np.array((np.sin(angle*np.pi/180.0), np.cos(angle*np.pi/180.0))) # unit vector
                    kwavedata = np.fromfunction(utils.kwave2d, arrayshape,
                            spatialwavelen=spatialwavelen,
                            center=center,
                            offset=offset,
                            khat=khat)
                    pupilarray = pupil * np.exp(1j * peak * kwavedata)
                    imagefield = ft.perform(pupilarray, fp_size_reselt, npup)
                    image_intensity = (imagefield*imagefield.conj()).real
                    psf = image_intensity / image_intensity.max()  # peak intensity unity
                    fits.PrimaryHDU(peak*kwavedata).writeto( 
                         odir+"/ex5_pupilarrayripple_{0:d}acrossD_peak_{1:.1f}.fits".format(nwaves,peak), 
                         overwrite=True)
Ejemplo n.º 2
0
def exer4(odir):
    """ Precise control of image positioning using phase slopes? """
    # instantiate an mft object:
    ft = matrixDFT.MatrixFourierTransform()
    npup = 100
    radius = 20.0
    fp_size_reselt = 100
    pupil = utils.makedisk(npup, radius=radius)
    tilts = ((0,0),
             (0*np.pi/npup, 0), 
             (1*np.pi/npup, 0),
             (2*np.pi/npup, 0),
             (3*np.pi/npup, 0),
             (4*np.pi/npup, 0),
             (95*np.pi/npup, 0),
             (npup*np.pi/npup, 0))
    phases = utils.phasearrays(npup, tilts)
    for nt, ph in enumerate(phases):
        pupilarray = pupil * np.exp(1j * ph)
        imagefield = ft.perform(pupilarray, fp_size_reselt, npup)
        image_intensity = (imagefield*imagefield.conj()).real
        psf = image_intensity / image_intensity.max()  # peak intensity unity
        fits.PrimaryHDU(psf).writeto( 
            odir+"/ex4_tiltPi_a_{0:.3f}_b_{1:.3f}.fits".format(tilts[nt][0],tilts[nt][1]), 
            overwrite=True)
Ejemplo n.º 3
0
def find_centroid(a, thresh):
    """
    Short Summary
    -------------
    Calculate the centroid of the image

    Long Summary
    ------------
    Original domain a, Fourier domain CV
    sft square image a to CV array, no loss or oversampling - like an fft.
    Normalize peak of abs(CV) to unity
    Create 'live area' mask from abs(CV) with slight undersizing
        (allow for 1 pixel shifts to live data still)
        (splodges, or full image a la KP)
    Calculate phase slopes using CV.angle() phase array
    Calculate mean of phase slopes over mask
    Normalize phase slopes to reflect image centroid location in pixels

    XY conventions meshed to lg_model conventions:
    if you simulate a psf with pixel_offset = ( (0.2, 0.4), ) then blind
        application  centroid = utils.find_centroid()

    returns the image centroid (0.40036, 0.2000093) pixels in image space. To
    use this in lg_model, nrm_core,... you will want to calculate the new image
    center using:
    image_center = utils.centerpoint(s) + np.array((centroid[1], centroid[0])
    and everything holds together sensibly looking at DS9 images of a.

    Parameters
    ----------
    a: 2D square float
        input image array

    thresh: float
        Threshold for the absolute value of the FT(a). Normalize abs(CV = FT(a))
        for unity peak, and define the support of good CV when this is above
        threshold, then find the phase slope of the CV only over this support.

    Returns
    -------
    htilt, vtilt: float, float
        Centroid of a, as offset from array center, as calculated by the DFT's.
    """
    ft = matrixDFT.MatrixFourierTransform()
    cv = ft.perform(a, a.shape[0], a.shape[0])
    cvmod, cvpha = np.abs(cv), np.angle(cv)

    cvmod = cvmod / cvmod.max()  # normalize to unity peak

    cvmask = np.where(cvmod >= thresh)

    cvmask_edgetrim = trim(cvmask, a.shape[0])

    htilt, vtilt = findslope(cvpha, cvmask_edgetrim)

    return htilt, vtilt
Ejemplo n.º 4
0
def create_input_datafiles(rfn=None):
    """
        Pupil and image data, true (zero mean) phase map (rad) 
        No de-tilting of the thase done.
        Didactic case, no input parameters: create the pupil & monochromatic image 
        on appropriate pixel scales.  
        Pupil and image arrays of same size, ready to FT into each 
        other losslessly.  For real data you'll need to worry about adjusting
        sizes to satify this sampling relation between the two planes.
        For finite bandwidth data image simulation will loop over wavelengths.
        For coarsely sampled pixels image simulation will need to use finer
        sampling and rebin to detector pixel sizes.
        [email protected] Jan 2019
    """
    pupil0 = utils.makedisk(250, radius=50)  # D=100 pix, array 250 pix
    pupil1 = (pupil0 -
              utils.makedisk(250, radius=15, ctr=(125.0, 75.0))) * pupil0
    pupil2 = (pupil1 -
              utils.makedisk(250, radius=15, ctr=(90.0, 90.0))) * pupil0
    pupil3 = (pupil2 -
              utils.makedisk(250, radius=15, ctr=(75.0, 125.0))) * pupil0
    fits.writeto(rfn + "0__input_pup.fits", pupil0, overwrite=True)
    fits.writeto(rfn + "1__input_pup.fits", pupil1, overwrite=True)
    fits.writeto(rfn + "2__input_pup.fits", pupil2, overwrite=True)
    fits.writeto(rfn + "3__input_pup.fits", pupil3, overwrite=True)

    for pnum in (0, 1, 2, 3):
        pupil = fits.getdata(rfn + "{:1d}__input_pup.fits".format(pnum))
        offctrbump = (pupil.shape[0] / 2.0 + 30, pupil.shape[1] / 2.0 + 20.0
                      )  # off-center bump
        ctrbump = (pupil.shape[0] / 2.0, pupil.shape[1] / 2.0)  # central bump
        #quarter wave phase bump, off center
        phase = (2.0 * np.pi / 4.0) * utils.makegauss(
            pupil.shape[0], ctr=offctrbump, sigma=10.0) * pupil
        phase = de_mean(phase, pupil)  # zero mean phase - doesn't change image
        fits.writeto(rfn + "{:1d}__input_truepha.fits".format(pnum),
                     phase,
                     overwrite=True)
        mft = matrixDFT.MatrixFourierTransform()
        imagefield = mft.perform(pupil * np.exp(1j * phase), pupil.shape[0],
                                 pupil.shape[0])
        image = (imagefield * imagefield.conj()).real
        fits.writeto(rfn + "{:1d}__input_img.fits".format(pnum),
                     image / image.sum(),
                     overwrite=True)
        del mft
Ejemplo n.º 5
0
def exer2(odir):
    """ Are you simulating samples or pixels?  Image plane sampling effect """
    # instantiate an mft object:
    ft = matrixDFT.MatrixFourierTransform()

    npup = 100
    radius = 20.0
    fp_size_reselts = 10
    fp_npixels = (4, 8, 16, 32)  
    pupil = utils.makedisk(npup, radius=radius)
    fits.PrimaryHDU(pupil).writeto(odir+"/ex2_pupil.fits", overwrite=True) # write pupil file
    for fpnpix in fp_npixels:
        imagefield = ft.perform(pupil, fp_size_reselts, fpnpix)
        imageintensity = (imagefield*imagefield.conj()).real
        psf = imageintensity / imageintensity.max()  # normalize to peak intensity unity
        zpsf = scipy.ndimage.zoom(psf, 32//fpnpix, order=0)
        fits.PrimaryHDU(zpsf).writeto(odir+"/ex2_nfppix{}_zoom{}.fits".format(fpnpix,32//fpnpix), overwrite=True)
Ejemplo n.º 6
0
def exer3(odir):
    """ What do phase slopes - tilts - in the pupil plane do? """
    # instantiate an mft object:
    ft = matrixDFT.MatrixFourierTransform()
    npup = 100
    radius = 20.0
    fp_size_reselt = 100
    pupil = utils.makedisk(npup, radius=radius)
    tilts = ((0,0),  (0.3,0), (1.0,0), (3.0,0))
    phases = utils.phasearrays(npup, tilts)
    for nt, ph in enumerate(phases):
        pupilarray = pupil * np.exp(1j * ph)
        imagefield = ft.perform(pupilarray, fp_size_reselt, npup)
        image_intensity = (imagefield*imagefield.conj()).real
        psf = image_intensity / image_intensity.max()  # peak intensity unity
        fits.PrimaryHDU(psf).writeto( 
            odir+"/ex3_tilt_a_{0:.3f}_b_{1:.3f}.fits".format(tilts[nt][0],tilts[nt][1]), 
            overwrite=True)
Ejemplo n.º 7
0
def exer1(odir):
    """ Let's get something for nothing if we can!!! """

    # instantiate an mft object:
    ft = matrixDFT.MatrixFourierTransform()

    npup = 100
    radius = 20.0
    fp_size_reselts = (100, 200, 300, 400)
    fp_npixels = (100, 200, 300, 400)  
    pupil = utils.makedisk(npup, radius=radius)
    fits.PrimaryHDU(pupil).writeto(odir+"/ex1_pupil.fits", overwrite=True) # write pupil file
    for (fpr,fpnpix) in zip(fp_npixels,fp_size_reselts):
        imagearray = np.zeros((400,400)) # create same-sized array for all 4 FOV's we calculate.
        imagefield = ft.perform(pupil, fpr, fpnpix)
        imageintensity = (imagefield*imagefield.conj()).real
        psf = imageintensity / imageintensity.max()  # normalize to peak intensity unity
        imagearray[200-fpnpix//2:200+fpnpix//2,200-fpnpix//2:200+fpnpix//2] = psf # center image in largest array size
        fits.PrimaryHDU(imagearray).writeto(odir+"/ex1_nfppix{}_fpsize{}.fits".format(fpnpix,fpr), overwrite=True)
Ejemplo n.º 8
0
    def __init__(
            self,
            pupilfn=None,
            imagefn=None,
            truephasefn=None,
            outputtag=0,
            initphase=None,
            damp=None,
            history=None,
            threshold=0.1,  # Strehl 99% ~ 1 - threshold^2
            maxiter=None,
            fileroot=None,
            verbose=True):
        """
        Assumes square arrays, equal sized arrays for pupil, image
        Iterations are pupil - to image - to pupil
        """

        # Requested attributes:
        self.damp = damp
        self.history = history
        self.threshold = threshold
        self.maxiter = maxiter
        self.fileroot = fileroot
        self.verbose = verbose
        self.pupil, self.image, self.truephase = data_input(
            pupilfn, imagefn, truephasefn, self.fileroot)
        self.tag = outputtag

        # Now for derived attributes:
        if initphase is None:
            self.initphase = self.pupil * 0.0  # initial phase guess
        else:
            self.initphase = initphase
        self.phase = self.initphase.copy()
        self.pupilpower = (self.pupil * self.pupil).sum()
        self.mft = matrixDFT.MatrixFourierTransform()
        self.npix = self.pupil.shape[0]  # size of the arrays
        # create interior of an unapodized pupil, for calculating mean tilts
        self.interior = interior(self.pupil)

        # Create initial guess at wavefront
        self.plane = PUPIL
        self.wavefront = self.pupil * np.exp(1j * self.phase)

        # List of phases to track evolution of GS
        self.iter = 0  # running count
        self.savediters = [
            -1,
        ]  #
        self.pupilintensities = [
            self.pupil,
        ]  #
        self.pupilphases = [
            self.initphase,
        ]  #
        self.imageintensities = [
            self.image,
        ]  #

        print("\nmaxiter {:d}  convergence if sigma < {:.2e}\n".format(
            self.maxiter, self.threshold))
        print("\nExample {:s} beginning\n".format(self.tag))
Ejemplo n.º 9
0
def create_input_datafiles_bumps(rfn=None):
    """
        returns: list of mgsdtasets
        Didactic case, no input parameters: create the pupil & monochromatic image 
        on appropriate pixel scales.  
        Pupil and image arrays of same size, ready to FT into each 
        other losslessly.  For real data you'll need to worry about adjusting
        sizes to satify this sampling relation between the two planes.
        For finite bandwidth data image simulation will loop over wavelengths.
        For coarsely sampled pixels image simulation will need to use finer
        sampling and rebin to detector pixel sizes.
        [email protected] Jan 2019
    """
    mft = matrixDFT.MatrixFourierTransform()

    pupilradius = 50
    pupil = utils.makedisk(250, radius=pupilradius)  # D=100 pix, array 250 pix
    pupilindex = np.where(pupil >= 1)
    pupilfn = rfn + "__input_pup.fits"
    fits.writeto(pupilfn, pupil, overwrite=True)

    mgsdatasets = []

    defocus_list = (-12.0, 12.0, -10.0, 10.0, -8.0, 8.0, -6.0, 6.0, -4.0, 4.0,
                    -2.0, 2.0)
    print(defocus_list)
    number_of_d_across_D = range(
        1, 16)  # different aberrations - dia of bump in pupil

    rippleamp = 1.0  # radians, 1/6.3 waves, about 340nm at 2.12 micron  Bump height.  bad var name

    # for a variety of bumps across the pupil:
    for nwaves in number_of_d_across_D:  # preserve var name nwaves from ripple case - bad varname here.
        pupil = fits.getdata(pupilfn)
        mgsdataset = [pupilfn]  # an mgsdataset starts with the pupil file...

        rbump = 4.0 * float(
            pupilradius) / nwaves  # sigma of gaussian bump in pupil
        #print("{:.1e}".format(rbump))
        bump = utils.makegauss(250, ctr=(145.0, 145.0),
                               sigma=rbump)  # D=100 pix, array 250 pix

        rbump = 0.5 * float(pupilradius) / nwaves  # rad of disk bump in pupil
        #print("{:.1e}".format(rbump))
        bump = utils.makedisk(250, radius=rbump,
                              ctr=(145.0, 145.0))  # D=100 pix, array 250 pix

        bump = (1.0 / np.sqrt(2)) * bump / bump[pupilindex].std(
        )  # 0.5 variance aberration, same SR hit
        ripplephase = rippleamp * bump  # bad var name

        phasefn = rfn + "bump_{0:d}acrossD_peak_{1:.1f}_pha.fits".format(
            int(nwaves), rippleamp)
        fits.PrimaryHDU(ripplephase).writeto(phasefn, overwrite=True)
        mgsdataset.append(
            phasefn)  # an mgsdataset now pupil file, phase abberration file,

        # Now create images for each defocus in the list
        #
        # First a utility array, parabola P2V unity over pupil, zero outside pupil
        prad = pupilradius  # for parabola=unity at edge of active pupil 2% < unity
        center = utils.centerpoint(pupil.shape[0])
        unityP2Vparabola = np.fromfunction(
            utils.parabola2d, pupil.shape, cx=center[0],
            cy=center[1]) / (prad * prad) * pupil
        fits.writeto(rfn + "unityP2Vparabola.fits",
                     unityP2Vparabola,
                     overwrite=True)  # sanity check - write it out

        for defoc in defocus_list:  # defoc units are waves, Peak-to-Valley

            defocusphase = unityP2Vparabola * 2 * np.pi * defoc
            aberfn = "pup__bump_defoc_{:.1f}wav.fits".format(defoc)
            fits.writeto(rfn + aberfn, defocusphase, overwrite=True)

            aberfn = rfn + "pup_bump_{0:d}acrossD_peak_{1:.1f}_defoc_{2:.1f}wav.fits".format(
                int(nwaves), rippleamp, defoc)
            imagfn = rfn + "__input_" + "img_bump_{0:d}acrossD_peak_{1:.1f}_defoc_{2:.1f}wav.fits".format(
                int(nwaves), rippleamp, defoc)

            aber = defocusphase + ripplephase
            fits.writeto(aberfn, aber, overwrite=True)

            imagefield = mft.perform(pupil * np.exp(1j * aber), pupil.shape[0],
                                     pupil.shape[0])
            image = (imagefield * imagefield.conj()).real

            fits.writeto(imagfn, image / image.sum(), overwrite=True)

            mgsdataset.append((defoc, imagfn, aberfn))

        mgsdatasets.append(mgsdataset)

    # Prepare a quicklook at signal in the pairs of defocussed images:
    side = 2.0  # inches, quicklook at curvatue signal imshow

    ndefoc = len(defocus_list) // 2
    nspatfreq = len(number_of_d_across_D)
    magnif = 2.0
    fig = plt.figure(1, (nspatfreq * magnif, ndefoc * magnif))
    grid = ImageGrid(
        fig,
        111,  # similar to subplot(111)
        nrows_ncols=(ndefoc, nspatfreq),  # creates 2x2 grid of axes
        axes_pad=0.02,  # pad between axes in inch.
    )
    i = 0
    iwav = 0
    for nwaves in number_of_d_across_D:
        ifoc = 0
        maxlist = []
        for defoc in defocus_list:  # defoc units are waves, Peak-to-Valley
            if defoc > 0:
                imagfn_pos = rfn + "__input_" + "img_bump_{0:d}acrossD_peak_{1:.1f}_defoc_{2:.1f}wav.fits".format(
                    int(nwaves), rippleamp, defoc)
                imagfn_neg = rfn + "__input_" + "img_bump_{0:d}acrossD_peak_{1:.1f}_defoc_{2:.1f}wav.fits".format(
                    int(nwaves), rippleamp, -defoc)
                datapos = fits.getdata(imagfn_pos).flatten()
                dataneg = fits.getdata(imagfn_neg).flatten()
                quicklook = (datapos - dataneg[::-1]).reshape((250, 250)) * (
                    defoc * defoc)  # normalize to equal brightness signal
                maxlist.append(quicklook.max())
                print("max {:.1e}".format(quicklook.max())
                      )  # print me out and fix the limits for imshow by hand.
                #fits.writeto(imagfn_pos.replace("img","qlk"), quicklook, overwrite=True)
                i = iwav + (ndefoc - ifoc - 1) * nspatfreq
                grid[i].imshow(
                    quicklook,
                    origin='lower',
                    cmap=plt.get_cmap(
                        'ocean'),  # RdBu, gist_rainbow, ocean, none
                    vmin=-3.0e-3,
                    vmax=3.0e-3)  # The AxesGrid object work as a list of axes.
                grid[i].set_xticks([])
                grid[i].set_yticks([])
                grid[i].text(20,
                             220,
                             "D/{:d} dia bump".format(int(nwaves + 1)),
                             color='y',
                             weight='bold')
                grid[i].text(20,
                             20,
                             "+/-{:d}w PV".format(int(defoc)),
                             color='w',
                             weight='bold')
                #print('iwav', iwav, 'ifoc', ifoc, 'i:', i)
                ifoc += 1
        iwav += 1
    strtop = "Wavefront signal from piston phase bumps of different diameters vs. defocus either side of focus"
    strbot = "Anand S. 2019, after Dean & Bowers, JOSA A 2003 (figs. 4 & 7)"
    fig.text(0.02, 0.94, strtop, fontsize=18, weight='bold')
    fig.text(0.02, 0.05, strbot, fontsize=14)
    plt.tight_layout()
    plt.savefig("DeanBowers2003_signal_vs_defocus_bump.png",
                dpi=150,
                pad_inches=1.0)
    plt.show()

    #print("Unity P-V parabola:", arraystats(unityP2Vparabola))
    """
    for dataset in mgsdatasets: 
        print("MGS data set: pupilfn, aber, (defoc/PVwaves, imagefn, defoc+aber), (repeats)")
        for thing in dataset: 
            print("\t", thing)
        print("")
    """
    return mgsdatasets
Ejemplo n.º 10
0
def create_input_datafiles(rfn=None):
    """
        returns: list of mgsdtasets
            pupilfn: name of pupil file
        Pupil and image data, true (zero mean) phase map (rad) 
        No de-tilting of the thase done.
        Didactic case, no input parameters: create the pupil & monochromatic image 
        on appropriate pixel scales.  
        Pupil and image arrays of same size, ready to FT into each 
        other losslessly.  For real data you'll need to worry about adjusting
        sizes to satify this sampling relation between the two planes.
        For finite bandwidth data image simulation will loop over wavelengths.
        For coarsely sampled pixels image simulation will need to use finer
        sampling and rebin to detector pixel sizes.
        [email protected] Jan 2019
    """
    mft = matrixDFT.MatrixFourierTransform()

    pupilradius = 50
    pupil = utils.makedisk(250, radius=pupilradius)  # D=100 pix, array 250 pix
    pupilfn = rfn + "__input_pup.fits"
    fits.writeto(pupilfn, pupil, overwrite=True)

    mgsdatasets = []

    dfoc_max = 12
    nfoci = 8  # number of defocus steps, in geo prog
    ffac = pow(10, np.log10(dfoc_max) / nfoci)
    defocus_list = []
    for i in range(nfoci):
        defocus_list.append(ffac)
        defocus_list.append(-ffac)
        ffac *= pow(10, np.log10(dfoc_max) / nfoci)
    defocus_list.reverse()
    defocus_list = (-12.0, 12.0, -10.0, 10.0, -8.0, 8.0, -6.0, 6.0, -4.0, 4.0,
                    -2.0, 2.0)
    print(defocus_list)
    number_of_waves_across_D = range(
        1, 16)  # different aberrations - number of waves across pupil
    print(number_of_waves_across_D)

    rippleamp = 1.0  # radians, 1/6.3 waves, about 340nm at 2.12 micron
    ripplepha = 30.0  # degrees, just for fun
    rippleangle = 15.0  # degrees

    # for a variety of ripples across the pupil:
    for nwaves in number_of_waves_across_D:
        pupil = fits.getdata(pupilfn)
        mgsdataset = [pupilfn]  # an mgsdataset starts with the pupil file...
        spatialwavelen = 2.0 * pupilradius / nwaves
        khat = np.array((np.sin(rippleangle * np.pi / 180.0),
                         np.cos(rippleangle * np.pi / 180.0)))  # unit vector
        kwavedata = np.fromfunction(utils.kwave2d,
                                    pupil.shape,
                                    spatialwavelen=spatialwavelen,
                                    center=utils.centerpoint(pupil.shape[0]),
                                    offset=ripplepha,
                                    khat=khat)
        ripplephase = pupil * rippleamp * kwavedata
        #imagefield = ft.perform(pupilarray, fp_size_reselt, npup)  # remove this
        #image_intensity = (imagefield*imagefield.conj()).real  # remove this
        #psf = image_intensity / image_intensity.sum()  # total intensity unity  # remove this

        phasefn = rfn + "ripple_{0:d}acrossD_peak_{1:.1f}_pha.fits".format(
            int(nwaves), rippleamp)
        fits.PrimaryHDU(ripplephase).writeto(phasefn, overwrite=True)
        mgsdataset.append(
            phasefn)  # an mgsdataset now pupil file, phase abberration file,

        # Now create images for each defocus in the list
        #
        # First a utility array, parabola P2V unity over pupil, zero outside pupil
        prad = pupilradius  # for parabola=unity at edge of active pupil 2% < unity
        center = utils.centerpoint(pupil.shape[0])
        unityP2Vparabola = np.fromfunction(
            utils.parabola2d, pupil.shape, cx=center[0],
            cy=center[1]) / (prad * prad) * pupil
        fits.writeto(rfn + "unityP2Vparabola.fits",
                     unityP2Vparabola,
                     overwrite=True)  # sanity check - write it out

        for defoc in defocus_list:  # defoc units are waves, Peak-to-Valley

            defocusphase = unityP2Vparabola * 2 * np.pi * defoc
            aberfn = "pup_defoc_{:.1f}wav.fits".format(defoc)
            fits.writeto(rfn + aberfn, defocusphase, overwrite=True)

            aberfn = rfn + "pup_ripple_{0:d}acrossD_peak_{1:.1f}_defoc_{2:.1f}wav.fits".format(
                int(nwaves), rippleamp, defoc)
            imagfn = rfn + "__input_" + "img_ripple_{0:d}acrossD_peak_{1:.1f}_defoc_{2:.1f}wav.fits".format(
                int(nwaves), rippleamp, defoc)

            aber = defocusphase + ripplephase
            #fits.writeto(aberfn, aber, overwrite=True)

            imagefield = mft.perform(pupil * np.exp(1j * aber), pupil.shape[0],
                                     pupil.shape[0])
            image = (imagefield * imagefield.conj()).real

            fits.writeto(imagfn, image / image.sum(), overwrite=True)

            mgsdataset.append((defoc, imagfn, aberfn))

        mgsdatasets.append(mgsdataset)
        """
        phase = de_mean(phase, pupil) # zero mean phase - doesn't change image
        fits.writeto(rfn+"{:1d}__input_truepha.fits".format(pnum), phase, overwrite=True)
        mft = matrixDFT.MatrixFourierTransform()
        imagefield = mft.perform(pupil * np.exp(1j*phase), pupil.shape[0], pupil.shape[0])
        image =  (imagefield*imagefield.conj()).real 
        fits.writeto(rfn+"{:1d}__input_img.fits".format(pnum), image/image.sum(), overwrite=True)
        del mft
        """

    # Prepare a quicklook at signal in the pairs of defocussed images:
    side = 2.0  # inches, quicklook at curvatue signal imshow

    ndefoc = len(defocus_list) // 2
    nspatfreq = len(number_of_waves_across_D)
    magnif = 2.0
    fig = plt.figure(1, (nspatfreq * magnif, ndefoc * magnif))
    grid = ImageGrid(
        fig,
        111,  # similar to subplot(111)
        nrows_ncols=(ndefoc, nspatfreq),  # creates 2x2 grid of axes
        axes_pad=0.02,  # pad between axes in inch.
    )
    i = 0
    iwav = 0
    for nwaves in number_of_waves_across_D:
        ifoc = 0
        maxlist = []
        for defoc in defocus_list:  # defoc units are waves, Peak-to-Valley
            if defoc > 0:
                imagfn_pos = rfn + "__input_" + "img_ripple_{0:d}acrossD_peak_{1:.1f}_defoc_{2:.1f}wav.fits".format(
                    int(nwaves), rippleamp, defoc)
                imagfn_neg = rfn + "__input_" + "img_ripple_{0:d}acrossD_peak_{1:.1f}_defoc_{2:.1f}wav.fits".format(
                    int(nwaves), rippleamp, -defoc)
                datapos = fits.getdata(imagfn_pos).flatten()
                dataneg = fits.getdata(imagfn_neg).flatten()
                quicklook = (datapos - dataneg[::-1]).reshape((250, 250)) * (
                    defoc * defoc)  # normalize to equal brightness signal
                maxlist.append(quicklook.max())
                #print("max {:.1e}".format(quicklook.max()))
                #fits.writeto(imagfn_pos.replace("img","qlk"), quicklook, overwrite=True)
                i = iwav + (ndefoc - ifoc - 1) * nspatfreq
                grid[i].imshow(
                    quicklook,
                    origin='lower',
                    cmap=plt.get_cmap(
                        'ocean'),  # RdBu, gist_rainbow, ocean, none
                    vmin=-1.3e-2,
                    vmax=1.3e-2)  # The AxesGrid object work as a list of axes.
                grid[i].set_xticks([])
                grid[i].set_yticks([])
                grid[i].text(20,
                             220,
                             "{:d} ripples ax D".format(int(nwaves)),
                             color='y',
                             weight='bold')
                grid[i].text(20,
                             20,
                             "+/-{:d}w PV".format(int(defoc)),
                             color='w',
                             weight='bold')

                ifoc += 1
        iwav += 1
    strtop = "Wavefront signal from phase ripples across the pupil vs. defocus either side of focus"
    strbot = "Anand S. 2019, illustrating Dean & Bowers, JOSA A 2003 (figs. 4 & 7)"
    fig.text(0.02, 0.94, strtop, fontsize=18, weight='bold')
    fig.text(0.02, 0.05, strbot, fontsize=14)
    plt.tight_layout()
    plt.savefig("DeanBowers2003_signal_vs_defocus_ripple.pdf",
                dpi=150,
                pad_inches=1.0)
    plt.show()

    #print("Unity P-V parabola:", arraystats(unityP2Vparabola))
    """
    for dataset in mgsdatasets: 
        print("MGS data set: pupilfn, aber, (defoc/PVwaves, imagefn, defoc+aber), (repeats)")
        for thing in dataset: 
            print("\t", thing)
        print("")
    """
    return mgsdatasets
Ejemplo n.º 11
0
def write_zernike_cube(method="poppy"):
    """
    Function to write out a FITS file of a cube 36x100x400 where each slice is 3
    side-by-side PSFs with different standard deviations (0.3, 1, and 3) and a
    pupil image of the Zernicke used to create the PSFs

    "Explores aberrated PSF morphology at 2x Nyquist image sampling"

    Parameters
    ---------
    method : str
        Defines the method used to create the Zernikes. Either "poppy" or "original".
        "poppy" calls poppy.zernike.zernike1() function and uses Noll ordering.
        "original" is Anand's original code which reads in a previously created
        "ZfunctionsUnitVar.fits" file made from running ZernikeFitter.py.

    Returns
    -------
    psf_cube : array
        Returns a 36x100x400 cube containing the PSFs and Zernikes used to
        create them
    """

    stddevs = (0.3, 1.0, 3.0)  # in radians - aberration std devs
    nab = len(stddevs)

    if method.lower() == "original":
        zerns = fits.getdata("ZernikeFitting/ZfunctionsUnitVar.fits"
                             )  # 100 x 100 arrays in cube
    elif method.lower() == "poppy":
        zerns = np.array([
            poppy.zernike.zernike1(i, npix=100, outside=0)
            for i in range(1, 37)
        ])
    else:
        raise TypeError("Not a valid method selection")

    print("nzerns, nx, ny =", zerns.shape)

    gamma = 4  # Soummer gamma of oversampling in the Fourier domain. Use integer.
    imagefov = zerns.shape[
        1] // gamma  # in ft's results - lam/D if an image plane
    npix = gamma * imagefov

    # For storing nab eg 3 PSFs per Zernike, varying strengths and the Zernike function
    psfs = np.zeros(
        (zerns.shape[0], npix, npix * (nab + 1)))  # gamma oversampling

    ft = matrixDFT.MatrixFourierTransform()

    pupil = zerns[0, :, :].copy()

    # Perfect image
    imagefield = ft.perform(pupil, imagefov, gamma * imagefov)
    imageintensity = (imagefield * imagefield.conj()).real
    perfectpeakintensity = imageintensity.max()

    fits.writeto('perfectPSF.fits',
                 imageintensity / perfectpeakintensity,
                 overwrite=True,
                 checksum=True)

    for nz in range(zerns.shape[0]):
        for iab, ab in enumerate(stddevs):
            imagefield = ft.perform(pupil * np.exp(1j * zerns[nz, :, :] * ab),
                                    imagefov, npix)
            imageintensity = (imagefield * imagefield.conj()).real
            psfs[nz, :, iab * npix:(iab + 1) *
                 npix] = imageintensity / perfectpeakintensity
            # sfs[nz, :, (iab+1)*npix:] = scipy.ndimage.interpolation.zoom(zerns[nz,:,:], 1.0/gamma,
            #                                output=None, order=0)

        displayzern = zerns[nz, :, :] - zerns[nz, :, :].min()

        # For all non-piston Z's
        if nz != 0:
            displayzern = (zerns[nz, :, :] - zerns[nz, :, :].min()) / (
                zerns[nz, :, :].max() - zerns[nz, :, :].min())

        psfs[nz, :,
             nab * npix:] = displayzern * 0.5  # looks better w lower peak

    psf_cube = psfs.astype(np.float32)

    fits.writeto('zernedPSFcube_{}.fits'.format(method.lower()),
                 psf_cube,
                 overwrite=True)

    return psf_cube
Ejemplo n.º 12
0
def analytical_model(zernike_pol, coef, cali=False):
    """

    :param zernike_pol:
    :param coef:
    :param cali: bool; True if we already have calibration coefficients to use. False if we still need to create them.
    :return:
    """

    #-# Parameters
    dataDir = os.path.join(CONFIG_INI.get('local', 'local_data_path'),
                           'active')
    telescope = CONFIG_INI.get('telescope', 'name')
    nb_seg = CONFIG_INI.getint(telescope, 'nb_subapertures')
    tel_size_m = CONFIG_INI.getfloat(telescope, 'diameter') * u.m
    real_size_seg = CONFIG_INI.getfloat(
        telescope, 'flat_to_flat'
    )  # in m, size in meters of an individual segment flatl to flat
    size_seg = CONFIG_INI.getint(
        'numerical',
        'size_seg')  # pixel size of an individual segment tip to tip
    wvln = CONFIG_INI.getint(telescope, 'lambda') * u.nm
    inner_wa = CONFIG_INI.getint(telescope, 'IWA')
    outer_wa = CONFIG_INI.getint(telescope, 'OWA')
    tel_size_px = CONFIG_INI.getint(
        'numerical', 'tel_size_px')  # pupil diameter of telescope in pixels
    im_size_pastis = CONFIG_INI.getint(
        'numerical', 'im_size_px_pastis')  # image array size in px
    sampling = CONFIG_INI.getfloat('numerical', 'sampling')  # sampling
    size_px_tel = tel_size_m / tel_size_px  # size of one pixel in pupil plane in m
    px_sq_to_rad = (size_px_tel * np.pi / tel_size_m) * u.rad
    zern_max = CONFIG_INI.getint('zernikes', 'max_zern')
    sz = CONFIG_INI.getint('numerical', 'im_size_lamD_hcipy')

    # Create Zernike mode object for easier handling
    zern_mode = util.ZernikeMode(zernike_pol)

    #-# Mean subtraction for piston
    if zernike_pol == 1:
        coef -= np.mean(coef)

    #-# Generic segment shapes

    if telescope == 'JWST':
        # Load pupil from file
        pupil = fits.getdata(
            os.path.join(dataDir, 'segmentation', 'pupil.fits'))

        # Put pupil in randomly picked, slightly larger image array
        pup_im = np.copy(pupil)  # remove if lines below this are active
        #pup_im = np.zeros([tel_size_px, tel_size_px])
        #lim = int((pup_im.shape[1] - pupil.shape[1])/2.)
        #pup_im[lim:-lim, lim:-lim] = pupil
        # test_seg = pupil[394:,197:315]    # this is just so that I can display an individual segment when the pupil is 512
        # test_seg = pupil[:203,392:631]    # ... when the pupil is 1024
        # one_seg = np.zeros_like(test_seg)
        # one_seg[:110, :] = test_seg[8:, :]    # this is the centered version of the individual segment for 512 px pupil

        # Creat a mini-segment (one individual segment from the segmented aperture)
        mini_seg_real = poppy.NgonAperture(
            name='mini', radius=real_size_seg
        )  # creating real mini segment shape with poppy
        #test = mini_seg_real.sample(wavelength=wvln, grid_size=flat_diam, return_scale=True)   # fix its sampling with wavelength
        mini_hdu = mini_seg_real.to_fits(wavelength=wvln,
                                         npix=size_seg)  # make it a fits file
        mini_seg = mini_hdu[
            0].data  # extract the image data from the fits file

    elif telescope == 'ATLAST':
        # Create mini-segment
        pupil_grid = hcipy.make_pupil_grid(dims=tel_size_px,
                                           diameter=real_size_seg)
        focal_grid = hcipy.make_focal_grid(
            pupil_grid, sampling, sz, wavelength=wvln.to(
                u.m).value)  # fov = lambda/D radius of total image
        prop = hcipy.FraunhoferPropagator(pupil_grid, focal_grid)

        mini_seg_real = hcipy.hexagonal_aperture(circum_diameter=real_size_seg,
                                                 angle=np.pi / 2)
        mini_seg_hc = hcipy.evaluate_supersampled(
            mini_seg_real, pupil_grid, 4
        )  # the supersampling number doesn't really matter in context with the other numbers
        mini_seg = mini_seg_hc.shaped  # make it a 2D array

        # Redefine size_seg if using HCIPy
        size_seg = mini_seg.shape[0]

        # Make stand-in pupil for DH array
        pupil = fits.getdata(
            os.path.join(dataDir, 'segmentation', 'pupil.fits'))
        pup_im = np.copy(pupil)

    #-# Generate a dark hole mask
    #TODO: simplify DH generation and usage
    dh_area = util.create_dark_hole(
        pup_im, inner_wa, outer_wa, sampling
    )  # this might become a problem if pupil size is not same like pastis image size. fine for now though.
    if telescope == 'ATLAST':
        dh_sz = util.zoom_cen(dh_area, sz * sampling)

    #-# Import information form segmentation script
    Projection_Matrix = fits.getdata(
        os.path.join(dataDir, 'segmentation', 'Projection_Matrix.fits'))
    vec_list = fits.getdata(
        os.path.join(dataDir, 'segmentation', 'vec_list.fits'))  # in pixels
    NR_pairs_list = fits.getdata(
        os.path.join(dataDir, 'segmentation', 'NR_pairs_list_int.fits'))

    # Figure out how many NRPs we're dealing with
    NR_pairs_nb = NR_pairs_list.shape[0]

    #-# Chose whether calibration factors to do the calibraiton with
    if cali:
        filename = 'calibration_' + zern_mode.name + '_' + zern_mode.convention + str(
            zern_mode.index)
        ck = fits.getdata(
            os.path.join(dataDir, 'calibration', filename + '.fits'))
    else:
        ck = np.ones(nb_seg)

    coef = coef * ck

    #-# Generic coefficients
    # the coefficients in front of the non redundant pairs, the A_q in eq. 13 in Leboulleux et al. 2018
    generic_coef = np.zeros(
        NR_pairs_nb
    ) * u.nm * u.nm  # setting it up with the correct units this will have

    for q in range(NR_pairs_nb):
        for i in range(nb_seg):
            for j in range(i + 1, nb_seg):
                if Projection_Matrix[i, j, 0] == q + 1:
                    generic_coef[q] += coef[i] * coef[j]

    #-# Constant sum and cosine sum - calculating eq. 13 from Leboulleux et al. 2018
    if telescope == 'JWST':
        i_line = np.linspace(-im_size_pastis / 2., im_size_pastis / 2.,
                             im_size_pastis)
        tab_i, tab_j = np.meshgrid(i_line, i_line)
        cos_u_mat = np.zeros(
            (int(im_size_pastis), int(im_size_pastis), NR_pairs_nb))
    elif telescope == 'ATLAST':
        i_line = np.linspace(-(2 * sz * sampling) / 2.,
                             (2 * sz * sampling) / 2., (2 * sz * sampling))
        tab_i, tab_j = np.meshgrid(i_line, i_line)
        cos_u_mat = np.zeros((int((2 * sz * sampling)), int(
            (2 * sz * sampling)), NR_pairs_nb))

    # Calculating the cosine terms from eq. 13.
    # The -1 with each NR_pairs_list is because the segment names are saved starting from 1, but Python starts
    # its indexing at zero, so we have to make it start at zero here too.
    for q in range(NR_pairs_nb):
        # cos(b_q <dot> u): b_q with 1 <= q <= NR_pairs_nb is the basis of NRPS, meaning the distance vectors between
        #                   two segments of one NRP. We can read these out from vec_list.
        #                   u is the position (vector) in the detector plane. Here, those are the grids tab_i and tab_j.
        # We need to calculate the dot product between all b_q and u, so in each iteration (for q), we simply add the
        # x and y component.
        cos_u_mat[:, :, q] = np.cos(
            px_sq_to_rad *
            (vec_list[NR_pairs_list[q, 0] - 1, NR_pairs_list[q, 1] - 1, 0] *
             tab_i) + px_sq_to_rad *
            (vec_list[NR_pairs_list[q, 0] - 1, NR_pairs_list[q, 1] - 1, 1] *
             tab_j)) * u.dimensionless_unscaled

    sum1 = np.sum(
        coef**2
    )  # sum of all a_{k,l} in eq. 13 - this works only for single Zernikes (l fixed), because np.sum would sum over l too, which would be wrong.
    if telescope == 'JWST':
        sum2 = np.zeros(
            (int(im_size_pastis), int(im_size_pastis))
        ) * u.nm * u.nm  # setting it up with the correct units this will have
    elif telescope == 'ATLAST':
        sum2 = np.zeros(
            (int(2 * sz * sampling), int(2 * sz * sampling))) * u.nm * u.nm

    for q in range(NR_pairs_nb):
        sum2 = sum2 + generic_coef[q] * cos_u_mat[:, :, q]

    #-# Local Zernike
    if telescope == 'JWST':
        # Generate a basis of Zernikes with the mini segment being the support
        isolated_zerns = zern.hexike_basis(nterms=zern_max,
                                           npix=size_seg,
                                           rho=None,
                                           theta=None,
                                           vertical=False,
                                           outside=0.0)

        # Calculate the Zernike that is currently being used and put it on one single subaperture, the result is Zer
        # Apply the currently used Zernike to the mini-segment.
        if zernike_pol == 1:
            Zer = np.copy(mini_seg)
        elif zernike_pol in range(2, zern_max - 2):
            Zer = np.copy(mini_seg)
            Zer = Zer * isolated_zerns[zernike_pol - 1]

        # Fourier Transform of the Zernike - the global envelope
        mf = mft.MatrixFourierTransform()
        ft_zern = mf.perform(Zer, im_size_pastis / sampling, im_size_pastis)

    elif telescope == 'ATLAST':
        isolated_zerns = hcipy.make_zernike_basis(num_modes=zern_max,
                                                  D=real_size_seg,
                                                  grid=pupil_grid,
                                                  radial_cutoff=False)
        Zer = hcipy.Wavefront(mini_seg_hc * isolated_zerns[zernike_pol - 1],
                              wavelength=wvln.to(u.m).value)

        # Fourier transform the Zernike
        ft_zern = prop(Zer)

    #-# Final image
    if telescope == 'JWST':
        # Generating the final image that will get passed on to the outer scope, I(u) in eq. 13
        intensity = np.abs(ft_zern)**2 * (sum1.value + 2. * sum2.value)
    elif telescope == 'ATLAST':
        intensity = ft_zern.intensity.shaped * (sum1.value + 2. * sum2.value)

    # PASTIS is only valid inside the dark hole, so we cut out only that part
    if telescope == 'JWST':
        tot_dh_im_size = sampling * (outer_wa + 3)
        intensity_zoom = util.zoom_cen(
            intensity, tot_dh_im_size
        )  # zoom box is (owa + 3*lambda/D) wide, in terms of lambda/D
        dh_area_zoom = util.zoom_cen(dh_area, tot_dh_im_size)

        dh_psf = dh_area_zoom * intensity_zoom

    elif telescope == 'ATLAST':
        dh_psf = dh_sz * intensity
    """
    # Create plots.
    plt.subplot(1, 3, 1)
    plt.imshow(pupil, origin='lower')
    plt.title('JWST pupil and diameter definition')
    plt.plot([46.5, 464.5], [101.5, 409.5], 'r-')   # show how the diagonal of the pupil is defined

    plt.subplot(1, 3, 2)
    plt.imshow(mini_seg, origin='lower')
    plt.title('JWST individual mini-segment')

    plt.subplot(1, 3, 3)
    plt.imshow(dh_psf, origin='lower')
    plt.title('JWST dark hole')
    plt.show()
    """

    # dh_psf is the image of the dark hole only, the pixels outside of it are zero
    # intensity is the entire final image
    return dh_psf, intensity
Ejemplo n.º 13
0
def exer6(odir):
    """ Coronagraph train, no optimization for speed.  
    2nd order BLC, didactic example, fftlike """
    # instantiate an mft object:
    ft = matrixDFT.MatrixFourierTransform()

    npup = 250 # Size of all arrays
    radius = 50.0

    # Numerical reselts in DFT setup cf telescope reselts:
    # reselts of telescope - here its 0.4 reselts per DFT output image pixel if npup=250,radius=50.
    dftpixel = 2.0 * radius / npup
    # Jinc first zero in reselts of telescope...
    firstzero_optical_reselts = 10.0
    firstzero_numericalpixels = firstzero_optical_reselts / dftpixel
    print("Jinc firstzero_numericalpixels", firstzero_numericalpixels)

    jinc = np.fromfunction(utils.Jinc, (npup,npup),
                           c=utils.centerpoint(npup),
                           scale=firstzero_numericalpixels)
    fpm_blc2ndorder = 1 - jinc*jinc
    print("Jinc fpm min = ", fpm_blc2ndorder.min(), 
          "Jinc fpm max = ", fpm_blc2ndorder.max())

    # Pupil, Pupilphase, Apodizer, FP intensity, Intensity after FPM, 
    # Lyot intensity, Lyot Stop, Post-Lyot Stop Intensity, Final image.
    #
    # Set up optical train for a typical Lyot style or phase mask coronagraph:
    Cordict = {
        "Pupil": utils.makedisk(npup, radius=radius),
        "Pupilphase": None,
        "Apodizer": None,
        "FPintensity": None,
        "FPM": fpm_blc2ndorder,
        "LyotIntensity": None,
        "LyotStop":  utils.makedisk(npup, radius=41),
        "PostLyotStopIntensity": None,
        "FinalImage": None,
        "ContrastImage": None}


    # Propagate through the coronagraph train...
    # Start with perfect incoming plane wave, no aberrations
    efield = Cordict["Pupil"]
    # Put in phase aberrations:
    if Cordict["Pupilphase"] is not None:
        efield *= np.exp(1j*Cordict["Pupilphase"])
    # Apodize the entrance pupil:
    if Cordict["Apodizer"] is not None:
        efield *= Cordict["Apodizer"]
    # PROPAGATE TO FIRST FOCAL PLANE:
    efield = ft.perform(efield, npup, npup)
    # Store FPM intensity:
    Cordict["FPintensity"] = (efield * efield.conj()).real

    # Save no-Cor efield for normalization of cor image by peak of no-FPM image
    efield_NC = efield.copy()
    # Multiply by FPM transmission function
    # Lyot style - zero in center, phase mask style: zero integral over domain
    efield *=  Cordict["FPM"]

    # PROPAGATE TO LYOT PLANE:
    efield_NC = ft.perform(efield_NC, npup, npup)
    efield = ft.perform(efield, npup, npup)
    # Save Cor Lyot intensity;
    Cordict["LyotIntensity"] = (efield * efield.conj()).real
    # Apply Lyot stop:
    if Cordict["LyotStop"] is not None: efield_NC *= Cordict["LyotStop"]
    if Cordict["LyotStop"] is not None: efield *= Cordict["LyotStop"]
    # Save Cor Lyot intensity after applying Lyot stop;
    Cordict["PostLyotStopIntensity"] = (efield * efield.conj()).real

    # PROPAGATE TO FINAL IMAGE PLANE:
    efield_NC = ft.perform(efield_NC, npup, npup)
    efield = ft.perform(efield, npup, npup)
    final_image_intensity_NC = (efield_NC * efield_NC.conj()).real
    final_image_intensity = (efield * efield.conj()).real
    Cordict["FinalImage"] = (efield * efield.conj()).real
    Cordict["ContrastImage"] = (efield * efield.conj()).real / final_image_intensity_NC.max()

    # Write our coronagraph planes:
    planenames, cube = corcube(Cordict)
    # write planemames as fits keywords
    print(odir+"/ex6_BLC_2ndOrder.fits")
    fits.PrimaryHDU(cube).writeto(odir+"/ex6_BLC_2ndOrder.fits", overwrite=True)
    fobj = fits.open(odir+"/ex6_BLC_2ndOrder.fits")
    fobj[0].header["Pupil"] = 1
    fobj[0].header["FPI"] = (2, "focal plane Intensity")
    fobj[0].header["FPM"] = (3, "focal plane mask")
    fobj[0].header["LyotIntn"] = (4, "Lyot Intensity")
    fobj[0].header["LyotStop"] = 5
    fobj[0].header["PostLyot"] = (6, "Post Lyot Stop Intensity")
    fobj[0].header["CorIm"] = (7, "Raw cor image")
    fobj[0].header["Contrast"] = (8, "Cor image in contrast units")
    fobj.writeto(odir+"/ex6_BLC_2ndOrder.fits", overwrite=True)
Ejemplo n.º 14
0
def find_centroid(a, thresh, verbose=True):
    """Return centroid of input image
 
    Parameters
    ----------
    a: input numpy array (real, square), considered 'image space'
    thresh: Threshold for the absolute value of the FT(a).
            Normalize abs(CV = FT(a)) for unity peak, and define the support 
            of "good" CV when this is above threshold, then find the phase 
            slope of the CV only over this support.
    
    Returns
    -------
    htilt, vtilt: Centroid of a, as offset from array center, in pythonese as calculated by the DFT's.

    Original domain a, Fourier domain CV
    sft square image a to CV array, no loss or oversampling - like an fft.
    Normalize peak of abs(CV) to unity
    Create 'live area' mask from abs(CV) with slight undersizing 
        (allow for 1 pixel shifts to live data still)
        (splodges, or full image a la KP)
    Calculate phase slopes using CV.angle() phase array
    Calculate mean of phase slopes over mask
    Normalize phase slopes to reflect image centroid location in pixels
    By eye, looking at smoothness of phase slope array...
    JWST NIRISS F480 F430M F380M 0.02
    JWST NIRISS F277 0.02 - untested
    GPI - untested

    XY conventions meshed to LG_Model conventions:
    if you simulate a psf with pixel_offset = ( (0.2, 0.4), ) then blind application

        centroid = utils.find_centroid()  
        
    returns the image centroid (0.40036, 0.2000093) pixels in image space. To use this
    in LG_Model, nrm_core,... you will want to calculate the new image center using

        image_center = utils.centerpoint(s) + np.array((centroid[1], centroid[0])

    and everything holds together sensibly looking at DS9 images of a.

    [email protected] 2018.02
    """

    ft = matrixDFT.MatrixFourierTransform()
    cv = ft.perform(a, a.shape[0], a.shape[0]) 
    cvmod, cvpha = np.abs(cv), np.angle(cv)
    cvmod = cvmod/cvmod.max() # normalize to unity peak
    cvmask = np.where(cvmod >= thresh)
    cvmask_edgetrim = trim(cvmask, a.shape[0])
    htilt, vtilt = findslope(cvpha, cvmask_edgetrim, verbose)

    M = np.zeros(a.shape)
    print(">>> utils.find_centroid(): M.shape {0}, a.shape {1}".format(M.shape, a.shape))
    M[cvmask_edgetrim] = 1
    print(">>> utils.find_centroid(): M.shape {0}, cvpha.shape {1}".format(M.shape, cvpha.shape))
    if 0:
        fits.PrimaryHDU(data=a).writeto("~/gitsrc/nrm_analysis/rundir/3_test_LGmethods_data/img.fits", overwrite=True)
        fits.PrimaryHDU(data=cvmod).writeto("~/gitsrc/nrm_analysis/rundir/3_test_LGmethods_data/cvmod.fits", overwrite=True)
        fits.PrimaryHDU(data=M*cvpha*180.0/np.pi).writeto("~/gitsrc/nrm_analysis/rundir/3_test_LGmethods_data/cvpha.fits", overwrite=True)
        print("CV abs max = {}".format(cvmod.max()))
        print("utils.find_centroid: {0} locations of CV array in CV mask for this data".format(len(cvmask[0])))

    return htilt, vtilt