def init_wfs_geom(p_wfs: conf.Param_wfs, r0: float, p_tel: conf.Param_tel, p_geom: conf.Param_geom, ittime: float, verbose=1): """Compute the geometry of WFSs: valid subaps, positions of the subaps, flux per subap, etc... :parameters: p_wfs: (Param_wfs) : wfs settings r0: (float) : atmos r0 @ 0.5 microns p_tel: (Param_tel) : telescope settings geom: (Param_geom) : geom settings ittime: (float) : 1/loop frequency [s] verbose: (int) : (optional) display informations if 0 """ if p_geom.pupdiam: if p_wfs.type == scons.WFSType.SH or p_wfs.type == scons.WFSType.PYRHR or p_wfs.type == scons.WFSType.PYRLR: pdiam = p_geom.pupdiam // p_wfs.nxsub if (p_geom.pupdiam % p_wfs.nxsub > 0): pdiam += 1 else: pdiam = -1 p_wfs._pdiam = pdiam init_wfs_size(p_wfs, r0, p_tel, verbose) if not p_geom.is_init: # this is the wfs with largest # of subaps # the overall geometry is deduced from it #if not p_geom.pupdiam: if p_geom.pupdiam != 0 and p_geom.pupdiam != p_wfs._pdiam * p_wfs.nxsub: print("WARNING: pupdiam set value not correct") p_geom.pupdiam = p_wfs._pdiam * p_wfs.nxsub print("pupdiam used: ", p_geom.pupdiam) if p_wfs.type == scons.WFSType.PYRHR or p_wfs.type == scons.WFSType.PYRLR: geom_init(p_geom, p_tel, padding=p_wfs._nrebin) elif (p_wfs.type == scons.WFSType.SH): geom_init(p_geom, p_tel) else: raise RuntimeError("This WFS can not be used") if (p_wfs.type == scons.WFSType.PYRHR or p_wfs.type == scons.WFSType.PYRLR): init_pyrhr_geom(p_wfs, r0, p_tel, p_geom, ittime, verbose=1) elif (p_wfs.type == scons.WFSType.SH): init_sh_geom(p_wfs, r0, p_tel, p_geom, ittime, verbose=1) else: raise RuntimeError("This WFS can not be used")
def comp_new_fstop(wfs: Sensors, n: int, p_wfs: conf.Param_wfs, fssize: float, fstop: bytes): """ Compute a new field stop for pyrhr WFS :parameters: n : (int) : WFS index wfs : (Param_wfs) : WFS parameters fssize : (float) : field stop size [arcsec] fstop : (string) : "square" or "round" (field stop shape) """ fsradius_pixels = int(fssize / p_wfs._qpixsize / 2.) if (fstop == scons.FieldStopType.ROUND): p_wfs.fstop = fstop focmask = util.dist(p_wfs._Nfft, xc=p_wfs._Nfft / 2. - 0.5, yc=p_wfs._Nfft / 2. - 0.5) < (fsradius_pixels) # fstop_area = np.pi * (p_wfs.fssize/2.)**2. #UNUSED elif (p_wfs.fstop == scons.FieldStopType.SQUARE): p_wfs.fstop = fstop y, x = np.indices((p_wfs._Nfft, p_wfs._Nfft)) x -= (p_wfs._Nfft - 1.) / 2. y -= (p_wfs._Nfft - 1.) / 2. focmask = (np.abs(x) <= (fsradius_pixels)) * \ (np.abs(y) <= (fsradius_pixels)) # fstop_area = p_wfs.fssize**2. #UNUSED else: msg = "p_wfs " + str(n) + ". fstop must be round or square" raise ValueError(msg) # pyr_focmask = np.roll(focmask,focmask.shape[0]/2,axis=0) # pyr_focmask = np.roll(pyr_focmask,focmask.shape[1]/2,axis=1) pyr_focmask = focmask * 1.0 # np.fft.fftshift(focmask*1.0) p_wfs._submask = np.fft.fftshift(pyr_focmask).astype(np.float32) p_wfs_fssize = fssize wfs.d_wfs[n].set_submask(p_wfs._submask)
def init_sh_geom(p_wfs: conf.Param_wfs, r0: float, p_tel: conf.Param_tel, p_geom: conf.Param_geom, ittime: float, verbose: bool = True): """Compute the geometry of SH WFSs: valid subaps, positions of the subaps, flux per subap, etc... :parameters: p_wfs: (Param_wfs) : wfs settings r0: (float) : atmos r0 @ 0.5 microns p_tel: (Param_tel) : telescope settings geom: (Param_geom) : geom settings ittime: (float) : 1/loop frequency [s] verbose: (bool) : (optional) display informations if True """ p_wfs.nPupils = 1 # this is the i,j index of lower left pixel of subap in _spupil istart = ((np.linspace(0, p_geom.pupdiam, p_wfs.nxsub + 1))[:-1]).astype(np.int64) # Translation in _mupil useful for raytracing istart += 2 jstart = np.copy(istart) # sorting out valid subaps fluxPerSub = np.zeros((p_wfs.nxsub, p_wfs.nxsub), dtype=np.float32) for i in range(p_wfs.nxsub): indi = istart[i] # +2-1 (yorick->python) for j in range(p_wfs.nxsub): indj = jstart[j] # +2-1 (yorick->python) fluxPerSub[i, j] = np.sum(p_geom._mpupil[indi:indi + p_wfs._pdiam, indj:indj + p_wfs._pdiam]) # fluxPerSub[i,j] = np.where(p_geom._mpupil[indi:indi+pdiam,indj:indj+pdiam] > 0)[0].size fluxPerSub = fluxPerSub / p_wfs._pdiam**2. pupvalid = (fluxPerSub >= p_wfs.fracsub) * 1 p_wfs._isvalid = pupvalid.astype(np.int32) p_wfs._nvalid = int(np.sum(p_wfs._isvalid)) p_wfs._fluxPerSub = fluxPerSub.copy() validx = np.where(p_wfs._isvalid.T)[1].astype(np.int32) validy = np.where(p_wfs._isvalid.T)[0].astype(np.int32) p_wfs._validsubsx = validx p_wfs._validsubsy = validy p_wfs._validpuppixx = validx * p_wfs._pdiam + 2 p_wfs._validpuppixy = validy * p_wfs._pdiam + 2 # this defines how we cut the phase into subaps phasemap = np.zeros((p_wfs._pdiam * p_wfs._pdiam, p_wfs._nvalid), dtype=np.int32) X = np.indices((p_geom._n, p_geom._n)) tmp = X[1] + X[0] * p_geom._n n = p_wfs._nvalid # for i in range(p_wfs._nvalid): for i in range(n): indi = istart[p_wfs._validsubsy[i]] # +2-1 (yorick->python) indj = jstart[p_wfs._validsubsx[i]] phasemap[:, i] = tmp[indi:indi + p_wfs._pdiam, indj:indj + p_wfs._pdiam].flatten() p_wfs._phasemap = phasemap p_wfs._validsubsx *= p_wfs.npix p_wfs._validsubsy *= p_wfs.npix # this is a phase shift of 1/2 pix in x and y halfxy = np.linspace(0, 2 * np.pi, p_wfs._Nfft + 1)[0:p_wfs._pdiam] / 2. halfxy = np.tile(halfxy, (p_wfs._pdiam, 1)) halfxy += halfxy.T # p_wfs._halfxy = <float*>(halfxy*0.).data #dont work: half*0 is temp # python obj if (p_wfs.npix % 2 == 1 and p_wfs._nrebin % 2 == 1): # p_wfs._halfxy = <float*>(halfxy*0.) halfxy = np.zeros((p_wfs._pdiam, p_wfs._pdiam), dtype=np.float32) p_wfs._halfxy = halfxy.astype(np.float32) else: p_wfs._halfxy = halfxy.astype(np.float32) # this defines how we create a larger fov if required if (p_wfs._Ntot != p_wfs._Nfft): indi = int((p_wfs._Ntot - p_wfs._Nfft) / 2.) # +1 -1 (yorick>python) indj = int(indi + p_wfs._Nfft) X = np.indices((p_wfs._Nfft, p_wfs._Nfft)) + 1 # TODO: +1 ?? x = X[1] y = X[0] # hrpix tmp = np.zeros((p_wfs._Ntot, p_wfs._Ntot)) tmp[indi:indj, indi:indj] = np.roll(x + (y - 1) * p_wfs._Nfft, p_wfs._Nfft // 2, axis=0) tmp[indi:indj, indi:indj] = np.roll(tmp[indi:indj, indi:indj], p_wfs._Nfft // 2, axis=1) # hrmap=roll(hrpix) tmp = np.roll(tmp, p_wfs._Ntot // 2, axis=0) tmp = np.roll(tmp, p_wfs._Ntot // 2, axis=1) tmp = np.where(tmp.flatten())[0] p_wfs._hrmap = np.copy(tmp.astype(np.int32)) p_wfs._hrmap = np.reshape(p_wfs._hrmap, (p_wfs._hrmap.shape[0], 1)) # must be set even if unused else: tmp = np.zeros((2, 2)) p_wfs._hrmap = np.copy(tmp.astype(np.int32)) # must be set even if unused if (p_wfs._nrebin * p_wfs.npix % 2 != p_wfs._Ntot % 2): # +2-1 yorick>python indi = int((p_wfs._Ntot - p_wfs._nrebin * p_wfs.npix) / 2.) + 1 else: indi = int((p_wfs._Ntot - p_wfs._nrebin * p_wfs.npix) / 2.) + 0 indj = int(indi + p_wfs._nrebin * p_wfs.npix) X = np.indices((p_wfs._nrebin * p_wfs.npix, p_wfs._nrebin * p_wfs.npix)) x = (X[1] / p_wfs._nrebin).astype(np.int64) y = (X[0] / p_wfs._nrebin).astype(np.int64) # binindices binindices = np.zeros((p_wfs._Ntot, p_wfs._Ntot)) binindices[indi:indj, indi:indj] = x + y * p_wfs.npix + 1 binmap = np.zeros((p_wfs._nrebin * p_wfs._nrebin, p_wfs.npix * p_wfs.npix)) X = np.indices((p_wfs._Ntot, p_wfs._Ntot)) tmp = X[1] + X[0] * p_wfs._Ntot if (p_wfs.gsalt <= 0): binindices = np.roll(binindices, binindices.shape[0] // 2, axis=0) binindices = np.roll(binindices, binindices.shape[1] // 2, axis=1) for i in range(p_wfs.npix * p_wfs.npix): binmap[:, i] = tmp[np.where(binindices == i + 1)] # binmap=np.reshape(binmap.flatten("F"),(binmap.shape[0],binmap.shape[1])) p_wfs._binmap = np.copy(binmap.astype(np.int32)) dr0 = p_tel.diam / r0 * \ (0.5 / p_wfs.Lambda) ** 1.2 / \ np.cos(p_geom.zenithangle * CONST.DEG2RAD) ** 0.6 fwhmseeing = p_wfs.Lambda / \ (p_tel.diam / np.sqrt(p_wfs.nxsub ** 2. + (dr0 / 1.5) ** 2.)) / 4.848 kernelfwhm = np.sqrt(fwhmseeing**2. + p_wfs.kernel**2.) tmp = util.makegaussian(p_wfs._Ntot, kernelfwhm / p_wfs._qpixsize, p_wfs._Ntot // 2, p_wfs._Ntot // 2).astype(np.float32) tmp = np.roll(tmp, tmp.shape[0] // 2, axis=0) tmp = np.roll(tmp, tmp.shape[1] // 2, axis=1) tmp[0, 0] = 1. # this insures that even with fwhm=0, the kernel is a dirac tmp = tmp / np.sum(tmp) tmp = np.fft.fft2(tmp).astype(np.complex64) / (p_wfs._Ntot * p_wfs._Ntot) p_wfs._ftkernel = np.copy(tmp).astype(np.complex64) # dealing with photometry # telSurf = np.pi / 4. * (1 - p_tel.cobs**2.) * p_tel.diam**2. # from the guide star if (p_wfs.gsalt == 0): if (p_wfs.zerop == 0): p_wfs.zerop = 1e11 # p_wfs._nphotons = p_wfs.zerop * 10 ** (-0.4 * p_wfs.gsmag) * \ # p_wfs.optthroughput * \ # (p_tel.diam / p_wfs.nxsub) ** 2. * ittime # include throughput to WFS # for unobstructed subaperture # per iteration p_wfs._nphotons = compute_nphotons(scons.WFSType.SH, ittime, p_wfs.optthroughput, p_tel.diam, nxsub=p_wfs.nxsub, zerop=p_wfs.zerop, gsmag=p_wfs.gsmag, verbose=verbose) else: # we are dealing with a LGS # p_wfs._nphotons = p_wfs.lgsreturnperwatt * \ # p_wfs.laserpower * \ # p_wfs.optthroughput * \ # (p_tel.diam / p_wfs.nxsub) ** 2. * 1e4 * \ # ittime # detected by WFS # ... for given power # include throughput to WFS # for unobstructed subaperture # per iteration p_wfs._nphotons = compute_nphotons( scons.WFSType.SH, ittime, p_wfs.optthroughput, p_tel.diam, nxsub=p_wfs.nxsub, lgsreturnperwatt=p_wfs.lgsreturnperwatt, laserpower=p_wfs.laserpower, verbose=verbose)
def init_pyrhr_geom(p_wfs: conf.Param_wfs, r0: float, p_tel: conf.Param_tel, p_geom: conf.Param_geom, ittime: float, verbose: bool = True): """Compute the geometry of PYRHR WFSs: valid subaps, positions of the subaps, flux per subap, etc... :parameters: p_wfs: (Param_wfs) : wfs settings r0: (float) : atmos r0 @ 0.5 microns p_tel: (Param_tel) : telescope settings geom: (Param_geom) : geom settings ittime: (float) : 1/loop frequency [s] verbose: (bool) : (optional) display informations if True """ p_wfs.npix = p_wfs._pdiam p_wfs._Ntot = p_geom._n # Creating field stop mask fsradius_pixels = int(p_wfs.fssize / p_wfs._qpixsize / 2.) if (p_wfs.fstop == scons.FieldStopType.ROUND): focmask = util.dist(p_wfs._Nfft, xc=p_wfs._Nfft / 2. - 0.5, yc=p_wfs._Nfft / 2. - 0.5) < (fsradius_pixels) elif (p_wfs.fstop == scons.FieldStopType.SQUARE): X = np.indices((p_wfs._Nfft, p_wfs._Nfft)) + 1 # TODO: +1 ?? x = X[1] - (p_wfs._Nfft + 1.) / 2. y = X[0] - (p_wfs._Nfft + 1.) / 2. focmask = (np.abs(x) <= (fsradius_pixels)) * \ (np.abs(y) <= (fsradius_pixels)) else: msg = "PYRHR wfs fstop must be FieldStopType.[ROUND|SQUARE]" raise ValueError(msg) pyr_focmask = focmask * 1.0 # np.fft.fftshift(focmask*1.0) p_wfs._submask = np.fft.fftshift(pyr_focmask) # Creating pyramid mask pyrsize = p_wfs._Nfft cobs = p_tel.cobs rpup = p_geom.pupdiam / 2.0 dpup = p_geom.pupdiam nrebin = p_wfs._nrebin fracsub = p_wfs.fracsub if p_wfs.pyr_pup_sep == -1: pup_sep = int(p_wfs.nxsub) else: pup_sep = p_wfs.pyr_pup_sep # number of pix between two centers two pupil images y = np.tile(np.arange(pyrsize) - pyrsize / 2, (pyrsize, 1)) x = y.T Pangle = pup_sep * nrebin # Pyramid angle in HR pixels if p_wfs.nPupils == 0: p_wfs.nPupils = 4 # Centers is a nPupils x 2 array describing the position of the quadrant centers centers = Pangle / np.sin( np.pi / p_wfs.nPupils) * np.c_[np.cos( (2 * np.arange(p_wfs.nPupils) + 1) * np.pi / p_wfs.nPupils), np.sin( (2 * np.arange(p_wfs.nPupils) + 1) * np.pi / p_wfs.nPupils)] # In case nPupils == 4, we put the centers back in the normal A-B-C-D ordering scheme, for misalignment processing. if p_wfs.nPupils == 4: centers = np.round(centers[[1, 0, 2, 3], :]).astype(np.int32) # Add in the misalignments of quadrant centers, in LR px units if p_wfs.pyr_misalignments is not None: mis = np.asarray(p_wfs.pyr_misalignments) * nrebin else: mis = np.zeros((p_wfs.nPupils, 2), dtype=np.float32) # Pyramid as minimal intersect of 4 tilting planes pyr = 2 * np.pi / pyrsize * np.min(np.asarray( [(c[0] + m[0]) * x + (c[1] + m[1]) * y for c, m in zip(centers, mis)]), axis=0) p_wfs._halfxy = np.fft.fftshift(pyr.T) # Valid pixels identification # Generate buffer with pupil at center pup = np.zeros((pyrsize, pyrsize), dtype=np.int32) pup[pyrsize // 2 - p_geom._n // 2:pyrsize // 2 + p_geom._n // 2, pyrsize // 2 - p_geom._n // 2:pyrsize // 2 + p_geom._n // 2] = \ p_geom._mpupil # Shift the mask to obtain the geometrically valid pixels for qIdx in range(p_wfs.nPupils): quadOnCenter = np.roll( pup, tuple((centers[qIdx] + mis[qIdx]).astype(np.int32)), (0, 1)) mskRebin = util.rebin(quadOnCenter.copy(), [pyrsize // nrebin, pyrsize // nrebin]) if qIdx == 0: stackedSubap = np.roll( mskRebin >= fracsub, tuple((-centers[qIdx] / nrebin).astype(np.int32)), (0, 1)) mskRebTot = (mskRebin >= fracsub).astype(np.int32) else: stackedSubap += np.roll( mskRebin >= fracsub, tuple((-centers[qIdx] / nrebin).astype(np.int32)), (0, 1)) mskRebTot += (mskRebin >= fracsub).astype(np.int32) * (qIdx + 1) validRow, validCol = [], [] # If n == 4, we need to tweak -again- the order to feed # compass XY controllers in the appropriate order if p_wfs.nPupils == 4: centers = centers[[2, 1, 3, 0], :] # mis = mis[[2, 1, 3, 0], :] # We don't need mis anymore - but if so keep the order of centers for qIdx in range(p_wfs.nPupils): tmpWh = np.where( np.roll(stackedSubap, tuple((centers[qIdx] / nrebin).astype(np.int32)), (0, 1))) validRow += [tmpWh[0].astype(np.int32)] validCol += [tmpWh[1].astype(np.int32)] nvalid = validRow[0].size validRow = np.asarray(validRow).flatten() validCol = np.asarray(validCol).flatten() stack, index = np.unique(np.c_[validRow, validCol], axis=0, return_index=True) p_wfs._nvalid = nvalid p_wfs._validsubsx = validRow[np.sort(index)] p_wfs._validsubsy = validCol[np.sort(index)] p_wfs._hrmap = mskRebTot.astype(np.int32) if (p_wfs.pyr_pos is None): pixsize = (np.pi * p_wfs._qpixsize) / (3600 * 180) # scale_fact = 2 * np.pi / npup * \ # (p_wfs.Lambda / p_tel.diam / 4.848) / pixsize * p_wfs.pyr_ampl # Proposition de Flo scale_fact = 2 * np.pi / pyrsize * \ (p_wfs.Lambda * 1e-6 / p_tel.diam) / pixsize * p_wfs.pyr_ampl cx = scale_fact * \ np.sin((np.arange(p_wfs.pyr_npts)) * 2. * np.pi / p_wfs.pyr_npts) cy = scale_fact * \ np.cos((np.arange(p_wfs.pyr_npts)) * 2. * np.pi / p_wfs.pyr_npts) # mod_npts = p_wfs.pyr_npts #UNUSED else: if (verbose): print( "Using user-defined positions [arcsec] for the pyramid modulation" ) cx = p_wfs.pyr_pos[:, 0] / p_wfs._qpixsize cy = p_wfs.pyr_pos[:, 1] / p_wfs._qpixsize # mod_npts=cx.shape[0] #UNUSED p_wfs._pyr_cx = cx.copy() p_wfs._pyr_cy = cy.copy() # telSurf = np.pi / 4. * (1 - p_tel.cobs**2.) * p_tel.diam**2. # p_wfs._nphotons = p_wfs.zerop * \ # 10. ** (-0.4 * p_wfs.gsmag) * ittime * \ # p_wfs.optthroughput * telSurf p_wfs._nphotons = compute_nphotons(scons.WFSType.PYRHR, ittime, p_wfs.optthroughput, p_tel.diam, cobs=p_tel.cobs, zerop=p_wfs.zerop, gsmag=p_wfs.gsmag, verbose=verbose) # spatial filtering by the pixel extent: # *2/2 intended. min should be 0.40 = sinc(0.5)^2. y = np.tile(np.arange(pyrsize) - pyrsize // 2, (pyrsize, 1)) x = y.T x = x * 1. / pyrsize y = y * 1. / pyrsize sincar = np.fft.fftshift(np.sinc(x) * np.sinc(y)) # sincar = np.roll(np.pi*x*np.pi*y,x.shape[1],axis=1) p_wfs._sincar = sincar.astype(np.float32) #pup = p_geom._mpupil a = pyrsize // nrebin b = p_geom._n // nrebin pupvalid = stackedSubap[a // 2 - b // 2:a // 2 + b // 2, a // 2 - b // 2:a // 2 + b // 2] p_wfs._isvalid = pupvalid.T.astype(np.int32) validsubsx = np.where(pupvalid)[0].astype(np.int32) validsubsy = np.where(pupvalid)[1].astype(np.int32) istart = np.arange(p_wfs.nxsub + 2) * p_wfs.npix jstart = np.copy(istart) # sorting out valid subaps fluxPerSub = np.zeros((p_wfs.nxsub + 2, p_wfs.nxsub + 2), dtype=np.float32) for i in range(p_wfs.nxsub + 2): indi = istart[i] for j in range(p_wfs.nxsub + 2): indj = jstart[j] fluxPerSub[i, j] = np.sum(p_geom._mpupil[indi:indi + p_wfs.npix, indj:indj + p_wfs.npix]) fluxPerSub = fluxPerSub / p_wfs.nxsub**2. p_wfs._fluxPerSub = fluxPerSub.copy() phasemap = np.zeros((p_wfs.npix * p_wfs.npix, p_wfs._nvalid), dtype=np.int32) X = np.indices((p_geom._n, p_geom._n)) # we need c-like indice tmp = X[1] + X[0] * p_geom._n pyrtmp = np.zeros((p_geom._n, p_geom._n), dtype=np.int32) for i in range(len(validsubsx)): indi = istart[validsubsy[i]] # +2-1 (yorick->python indj = jstart[validsubsx[i]] phasemap[:, i] = tmp[indi:indi + p_wfs.npix, indj:indj + p_wfs.npix].flatten("C") pyrtmp[indi:indi + p_wfs.npix, indj:indj + p_wfs.npix] = i p_wfs._phasemap = phasemap p_wfs._pyr_offsets = pyrtmp # pshift
def init_wfs_size(p_wfs: conf.Param_wfs, r0: float, p_tel: conf.Param_tel, verbose=1): """Compute all the parameters usefull for further WFS image computation (array sizes) :parameters: p_wfs: (Param_wfs) : wfs settings r0: (float) : atmos r0 @ 0.5 microns p_tel: (Param_tel) : telescope settings verbose: (int) : (optional) display informations if 0 Scheme to determine arrays sizes sh : k = 6 p = k * d/r0 n = int(2*d*v/lambda/CONST.RAD2ARCSEC)+1 N = fft_goodsize(k*n/v*lambda/r0*CONST.RAD2ARCSEC) u = k * lambda / r0 * CONST.RAD2ARCSEC / N n = v/u - int(v/u) > 0.5 ? int(v/u)+1 : int(v/u) v = n * u Nt = v * Npix pyr : Fs = field stop radius in arcsec N size of big array for FFT in pixels P pupil diameter in pixels D diameter of telescope in m Nssp : number of pyr measurement points in the pupil Rf = Fs . N . D / lambda / P ideally we choose : Fs = lambda / D . Nssp / 2 if we want good sampling of r0 (avoid aliasing of speckles) P > D / r0 . m with m = 2 or 3 to get reasonable space between pupil images : N > P.(2 + 3S) with S close to 1 N must be a power of 2 to ease computation lets assume that Nssp is a divider of P scaling factor between high res pupil images in pyramid model and actual pupil size on camera images would be P / Nssp """ r0 = r0 * (p_wfs.Lambda * 2)**(6. / 5) if (r0 != 0): if (verbose): print("r0 for WFS :", "%3.2f" % r0, " m") # seeing = CONST.RAD2ARCSEC * (p_wfs.lambda * 1.e-6) / r0 if (verbose): print("seeing for WFS : ", "%3.2f" % (CONST.RAD2ARCSEC * (p_wfs.Lambda * 1.e-6) / r0), "\"") if (p_wfs._pdiam <= 0): # this case is usualy for the wfs with max # of subaps # we look for the best compromise between pixsize and fov if (p_wfs.type == scons.WFSType.SH): subapdiam = p_tel.diam / float(p_wfs.nxsub) # diam of subap k = 6 pdiam = int(k * subapdiam / r0) # number of phase points per subap if (pdiam < 16): pdiam = 16 # Must be even to keep ssp and actuators grids aligned in the pupil if ((pdiam * p_wfs.nxsub) % 2): pdiam += 1 nrebin = int(2 * subapdiam * p_wfs.pixsize / (p_wfs.Lambda * 1.e-6) / CONST.RAD2ARCSEC) + 1 nrebin = max(2, nrebin) # first atempt on a rebin factor # since we clipped pdiam we have to be carreful in nfft computation Nfft = util.fft_goodsize( int(pdiam / subapdiam * nrebin / p_wfs.pixsize * CONST.RAD2ARCSEC * (p_wfs.Lambda * 1.e-6))) elif (p_wfs.type == scons.WFSType.PYRHR or p_wfs.type == scons.WFSType.PYRLR): # while (pdiam % p_wfs.npix != 0) pdiam+=1; k = 3 if p_wfs.type == scons.WFSType.PYRHR else 2 pdiam = int(p_tel.diam / r0 * k) while (pdiam % p_wfs.nxsub != 0): pdiam += 1 # we choose to have a multiple of p_wfs.nxsub pdiam = pdiam // p_wfs.nxsub if (pdiam < 8): pdiam = 8 # quantum pixel size else: pdiam = p_wfs._pdiam # this case is for a wfs with fixed # of phase points Nfft = util.fft_goodsize(2 * pdiam) # size of the support in fourier domain # qpixsize = pdiam * \ # (p_wfs.Lambda * 1.e-6) / subapdiam * CONST.RAD2ARCSEC / Nfft # # quantum pixel size if (p_wfs.type == scons.WFSType.SH): subapdiam = p_tel.diam / float(p_wfs.nxsub) # diam of subap # size of the support in fourier domain # qpixsize = k * (p_wfs.Lambda*1.e-6) / r0 * CONST.RAD2ARCSEC / Nfft qpixsize = (pdiam * (p_wfs.Lambda * 1.e-6) / subapdiam * CONST.RAD2ARCSEC) / Nfft # # actual rebin factor if (p_wfs.pixsize / qpixsize - int(p_wfs.pixsize / qpixsize) > 0.5): nrebin = int(p_wfs.pixsize / qpixsize) + 1 else: nrebin = int(p_wfs.pixsize / qpixsize) # actual pixel size pixsize = nrebin * qpixsize if (pixsize * p_wfs.npix > qpixsize * Nfft): Ntot = util.fft_goodsize(int(pixsize * p_wfs.npix / qpixsize) + 1) else: Ntot = Nfft if (Ntot % 2 != Nfft % 2): Ntot += 1 elif (p_wfs.type == scons.WFSType.PYRHR or p_wfs.type == scons.WFSType.PYRLR): pdiam = pdiam * p_wfs.nxsub m = 3 # fft_goodsize( m * pdiam) Nfft = int(2**np.ceil(np.log2(m * pdiam))) nrebin = pdiam // p_wfs.nxsub while (Nfft % nrebin != 0): nrebin += 1 # we choose to have a divisor of Nfft pdiam = nrebin * p_wfs.nxsub Nfft = int(2**np.ceil(np.log2(m * pdiam))) qpixsize = (pdiam * (p_wfs.Lambda * 1.e-6) / p_tel.diam * CONST.RAD2ARCSEC) / Nfft Ntot = Nfft // pdiam * p_wfs.nxsub pixsize = qpixsize * nrebin pdiam = pdiam // p_wfs.nxsub p_wfs._pdiam = pdiam p_wfs.pixsize = pixsize p_wfs._qpixsize = qpixsize p_wfs._Nfft = Nfft p_wfs._Ntot = Ntot p_wfs._nrebin = nrebin p_wfs._subapd = p_tel.diam / p_wfs.nxsub if (verbose): if (p_wfs.type == scons.WFSType.SH): print("quantum pixsize : ", "%5.4f" % qpixsize, "\"") print("simulated FoV : ", "%3.2f" % (Ntot * qpixsize), "\" x ", "%3.2f" % (Ntot * qpixsize), "\"") print("actual pixsize : ", "%5.4f" % pixsize) print("actual FoV : ", "%3.2f" % (pixsize * p_wfs.npix), "\" x ", "%3.2f" % (pixsize * p_wfs.npix), "\"") print("number of phase points : ", p_wfs._pdiam) print("size of fft support : ", Nfft) print("size of HR spot support : ", Ntot) elif (p_wfs.type == scons.WFSType.PYRHR or p_wfs.type == scons.WFSType.PYRLR): print("quantum pixsize in pyr image : ", "%5.4f" % qpixsize, "\"") print("simulated FoV : ", "%3.2f" % (Nfft * qpixsize), "\" x ", "%3.2f" % (Nfft * qpixsize), "\"") print("number of phase points : ", p_wfs._pdiam * p_wfs.nxsub) print("size of fft support : ", Nfft)
def make_lgs_prof1d(p_wfs: conf.Param_wfs, p_tel: conf.Param_tel, prof: np.ndarray, h: np.ndarray, beam: float, center=b""): """same as prep_lgs_prof but cpu only. original routine from rico :parameters: p_tel: (Param_tel) : telescope settings prof: (np.ndarray[dtype=np.float32]) : Na profile intensity, in arbitrary units h: (np.ndarray[dtype=np.float32]) : altitude, in meters. h MUST be an array with EQUALLY spaced elements. beam: (float) : size in arcsec of the laser beam center: (string) : either "image" or "fourier" depending on where the centre should be. """ p_wfs._prof1d = prof p_wfs._profcum = np.zeros(h.shape[0] + 1, dtype=np.float32) p_wfs._profcum[1:] = prof.cumsum() subapdiam = p_tel.diam / p_wfs.nxsub # diam of subap if (p_wfs.nxsub > 1): xsubs = np.linspace((subapdiam - p_tel.diam) / 2, (p_tel.diam - subapdiam) / 2, p_wfs.nxsub).astype(np.float32) else: xsubs = np.zeros(1, dtype=np.float32) ysubs = xsubs.copy().astype(np.float32) # cdef int nP=prof.shape[0] #UNUSED hG = np.sum(h * prof) / np.sum(prof) x = np.arange(p_wfs._Ntot).astype(np.float32) - p_wfs._Ntot / 2 # x expressed in pixels. (0,0) is in the fourier-center. x = x * p_wfs._qpixsize # x expressed in arcseconds # cdef float dx=x[1]-x[0] #UNUSED # cdef float dh=h[1]-h[0] #UNUSED if (p_wfs.nxsub > 1): dOffAxis = np.sqrt((xsubs[p_wfs._validsubsy] - p_wfs.lltx)**2 + (ysubs[p_wfs._validsubsx] - p_wfs.llty)**2) else: dOffAxis = np.sqrt((xsubs - p_wfs.lltx)**2 + (ysubs - p_wfs.llty)**2) profi = np.zeros((p_wfs._Ntot, p_wfs._nvalid), dtype=np.float32) subsdone = np.ones(p_wfs._nvalid, dtype=np.int32) dif2do = np.zeros(p_wfs._nvalid, dtype=np.int32) while (np.any(subsdone)): tmp = dOffAxis[np.where(subsdone)][0] inds = np.where(dOffAxis == tmp)[0] # height, translated in arcsec due to perspective effect zhc = (h - hG) * (206265. * tmp / hG**2) dzhc = zhc[1] - zhc[0] if (p_wfs._qpixsize > dzhc): avg_zhc = np.zeros(zhc.size + 1, dtype=np.float32) avg_zhc[0] = zhc[0] avg_zhc[avg_zhc.size - 1] = zhc[zhc.size - 1] avg_zhc[1:-1] = 0.5 * (zhc[1:] + zhc[:-1]) avg_x = np.zeros(x.size + 1, dtype=np.float32) avg_x[0] = x[0] avg_x[avg_x.size - 1] = x[x.size - 1] avg_x[1:-1] = 0.5 * (x[1:] + x[:-1]) for i in range(inds.size): profi[:, inds[i]] = np.diff( np.interp(avg_x, avg_zhc, p_wfs._profcum)).astype(np.float32) else: for i in range(inds.size): profi[:, inds[i]] = np.interp(x, zhc, prof) subsdone[inds] = 0 w = beam / 2.35482005 if (w == 0): # TODO what is n n = 1 g = np.zeros(n, dtype=np.float32) if (center == b"image"): g[n / 2 - 1] = 0.5 g[n / 2] = 0.5 else: g[n / 2] = 1 else: if (center == b"image"): if ((p_wfs.npix * p_wfs._nrebin) % 2 != p_wfs._Nfft % 2): g = np.exp(-(x + p_wfs._qpixsize)**2 / (2 * w**2.)) else: g = np.exp(-(x + p_wfs._qpixsize / 2)**2 / (2 * w**2.)) else: g = np.exp(-x**2 / (2 * w**2.)) p_wfs._ftbeam = np.fft.fft(g).astype(np.complex64) p_wfs._beam = g.astype(np.float32) # convolved profile in 1D. g_extended = np.tile(g, (p_wfs._nvalid, 1)).T p1d = np.fft.ifft(np.fft.fft(profi, axis=0) * np.fft.fft(g_extended, axis=0), axis=0).real.astype(np.float32) p1d = p1d * p1d.shape[0] p1d = np.roll(p1d, int(p_wfs._Ntot / 2. + 0.5), axis=0) p1d = np.abs(p1d) im = np.zeros((p1d.shape[1], p1d.shape[0], p1d.shape[0]), dtype=np.float32) for i in range(p1d.shape[0]): im[:, i, :] = g[i] * p1d.T if (ysubs.size > 1): azimuth = np.arctan2(ysubs[p_wfs._validsubsy] - p_wfs.llty, xsubs[p_wfs._validsubsx] - p_wfs.lltx) else: azimuth = np.arctan2(ysubs - p_wfs.llty, xsubs - p_wfs.lltx) p_wfs._azimuth = azimuth if (center == b"image"): xcent = p_wfs._Ntot / 2. + 0.5 ycent = xcent else: xcent = p_wfs._Ntot / 2. + 1 ycent = xcent if (ysubs.size > 0): # TODO rotate im = util.rotate3d(im, azimuth * 180 / np.pi, xcent, ycent) max_im = np.max(im, axis=(1, 2)) im = (im.T / max_im).T else: im = util.rotate(im, azimuth * 180 / np.pi, xcent, ycent) im = im / np.max(im) p_wfs._lgskern = im.T
def prep_lgs_prof(p_wfs: conf.Param_wfs, nsensors: int, p_tel: conf.Param_tel, sensors: Sensors, center=b"", imat=0): """The function returns an image array(double,n,n) of a laser beacon elongated by perpective effect. It is obtaind by convolution of a gaussian of width "lgsWidth" arcseconds, with the line of the sodium profile "prof". The altitude of the profile is the array "h". :parameters: p_wfs: (Param_wfs) : WFS settings nsensors: (int) : wfs index p_tel: (Param_tel) : telescope settings Sensors: (Sensors) : WFS object center: (string) : either "image" or "fourier" depending on where the centre should be. Computation of LGS spot from the sodium profile: Everything is done here in 1D, because the Na profile is the result of the convolution of a function P(x,y) = profile(x) . dirac(y) by a gaussian function, for which variables x and y can be split : exp(-(x^2+y^2)/2.s^2) = exp(-x^2/2.s^2) * exp(-y^2/2.s^2) The convolution is (symbol $ denotes integral) C(X,Y) = $$ exp(-x^2/2.s^2) * exp(-y^2/2.s^2) * profile(x-X) * dirac(y-Y) dx dy First one performs the integration along y C(X,Y) = exp(-Y^2/2.s^2) $ exp(-x^2/2.s^2) * profile(x-X) dx which shows that the profile can be computed by - convolving the 1-D profile - multiplying it in the 2nd dimension by a gaussian function If one has to undersample the inital profile, then some structures may be "lost". In this case, it's better to try to "save" those structures by re-sampling the integral of the profile, and then derivating it afterwards. Now, if the initial profile is a coarse one, and that one has to oversample it, then a simple re-sampling of the profile is adequate. """ if (p_wfs.proftype is None or p_wfs.proftype == b""): p_wfs.set_proftype(scons.ProfType.GAUSS1) profilename = scons.ProfType.FILES[p_wfs.proftype] profile_path = shesha_savepath + profilename print("reading Na profile from", profile_path) prof = np.load(profile_path) make_lgs_prof1d(p_wfs, p_tel, np.mean(prof[1:, :], axis=0), prof[0, :], p_wfs.beamsize, center="image") p_wfs.set_altna(prof[0, :].astype(np.float32)) p_wfs.set_profna(np.mean(prof[1:, :], axis=0).astype(np.float32)) p_wfs._prof1d = p_wfs._profna p_wfs._profcum = np.zeros(p_wfs._profna.size + 1, dtype=np.float32) p_wfs._profcum[1:] = p_wfs._profna.cumsum() subapdiam = p_tel.diam / p_wfs.nxsub # diam of subap if (p_wfs.nxsub > 1): xsubs = np.linspace((subapdiam - p_tel.diam) / 2, (p_tel.diam - subapdiam) / 2, p_wfs.nxsub).astype(np.float32) else: xsubs = np.zeros(1, dtype=np.float32) ysubs = xsubs.copy().astype(np.float32) # center of gravity of the profile hG = np.sum(p_wfs._altna * p_wfs._profna) / np.sum(p_wfs._profna) x = np.arange(p_wfs._Ntot).astype(np.float32) - p_wfs._Ntot / 2 # x expressed in pixels. (0,0) is in the fourier-center x = x * p_wfs._qpixsize # x expressed in arcseconds # cdef float dx=x[1]-x[0] #UNUSED dh = p_wfs._altna[1] - p_wfs._altna[0] if (p_wfs.nxsub > 1): dOffAxis = np.sqrt((xsubs[p_wfs._validsubsx] - p_wfs.lltx)**2 + (ysubs[p_wfs._validsubsy] - p_wfs.llty)**2) else: dOffAxis = np.sqrt((xsubs - p_wfs.lltx)**2 + (ysubs - p_wfs.llty)**2) if (imat > 0): dOffAxis *= 0. w = p_wfs.beamsize / 2.35482005 # TODO: FIXME if (w == 0): # TODO what is n n = 1 g = np.zeros(n, dtype=np.float32) if (center == b"image"): g[n / 2 - 1] = 0.5 g[n / 2] = 0.5 else: g[n / 2] = 1 else: if (center == b"image"): g = np.exp(-(x + p_wfs._qpixsize / 2)**2 / (2 * w**2.)) else: g = np.exp(-x**2 / (2 * w**2.)) p_wfs._ftbeam = np.fft.fft(g, axis=0).astype(np.complex64) p_wfs._beam = g # convolved profile in 1D. if (xsubs.size > 1): azimuth = np.arctan2(ysubs[p_wfs._validsubsy] - p_wfs.llty, xsubs[p_wfs._validsubsx] - p_wfs.lltx) else: azimuth = np.arctan2(ysubs - p_wfs.llty, xsubs - p_wfs.lltx) p_wfs._azimuth = azimuth sensors.init_lgs(nsensors, p_wfs._prof1d.size, hG, p_wfs._altna[0], dh, p_wfs._qpixsize, dOffAxis, p_wfs._prof1d, p_wfs._profcum, p_wfs._beam, p_wfs._ftbeam, p_wfs._azimuth)
def init_pyrhr_geom(p_wfs: conf.Param_wfs, r0: float, p_tel: conf.Param_tel, p_geom: conf.Param_geom, ittime: float, verbose=1): """Compute the geometry of PYRHR WFSs: valid subaps, positions of the subaps, flux per subap, etc... :parameters: p_wfs: (Param_wfs) : wfs settings r0: (float) : atmos r0 @ 0.5 microns p_tel: (Param_tel) : telescope settings geom: (Param_geom) : geom settings ittime: (float) : 1/loop frequency [s] verbose: (int) : (optional) display informations if 0 """ p_wfs.npix = p_wfs._pdiam p_wfs._Ntot = p_geom._n # Creating field stop mask fsradius_pixels = int(p_wfs.fssize / p_wfs._qpixsize / 2.) if (p_wfs.fstop == scons.FieldStopType.ROUND): focmask = util.dist(p_wfs._Nfft, xc=p_wfs._Nfft / 2. + 0.5, yc=p_wfs._Nfft / 2. + 0.5) < (fsradius_pixels) elif (p_wfs.fstop == scons.FieldStopType.SQUARE): X = np.indices((p_wfs._Nfft, p_wfs._Nfft)) + 1 # TODO: +1 ?? x = X[1] - (p_wfs._Nfft + 1.) / 2. y = X[0] - (p_wfs._Nfft + 1.) / 2. focmask = (np.abs(x) <= (fsradius_pixels)) * \ (np.abs(y) <= (fsradius_pixels)) else: msg = "PYRHR wfs fstop must be FieldStopType.[ROUND|SQUARE]" raise ValueError(msg) pyr_focmask = focmask * 1.0 # np.fft.fftshift(focmask*1.0) p_wfs._submask = np.fft.fftshift(pyr_focmask) # Creating pyramid mask pyrsize = p_wfs._Nfft cobs = p_tel.cobs rpup = p_geom.pupdiam / 2.0 dpup = p_geom.pupdiam nrebin = p_wfs._nrebin fracsub = p_wfs.fracsub if p_wfs.pyr_pup_sep == -1: pup_sep = int(p_wfs.nxsub) else: pup_sep = p_wfs.pyr_pup_sep # number of pix betw two centers two pupil images y = np.tile(np.arange(pyrsize) - pyrsize / 2, (pyrsize, 1)) x = y.T Pangle = pup_sep * nrebin # Pyramid angle in HR pixels centers = Pangle * np.array([[-1, 1], [1, 1], [-1, -1], [1, -1]], dtype=np.float32) if p_wfs.pyr_misalignments is not None: mis = np.asarray(p_wfs.pyr_misalignments) * nrebin else: mis = np.zeros((4, 2), dtype=np.float32) # Pyrarmid as minimal intersect of 4 tilting planes pyr = 2 * np.pi / pyrsize * np.min(np.asarray( [(c[0] + m[0]) * x + (c[1] + m[1]) * y for c, m in zip(centers, mis)]), axis=0) p_wfs._halfxy = np.fft.fftshift(pyr.T) # Valid pixels identification # Generate buffer with pupil at center pup = np.zeros((pyrsize, pyrsize), dtype=np.int32) pup[pyrsize // 2 - p_geom._n // 2:pyrsize // 2 + p_geom._n // 2, pyrsize // 2 - p_geom._n // 2:pyrsize // 2 + p_geom._n // 2] = p_geom._mpupil for qIdx in range(4): quadOnCenter = np.roll( pup, tuple((centers[qIdx] + mis[qIdx]).astype(np.int32)), (0, 1)) mskRebin = util.rebin(quadOnCenter.copy(), [pyrsize // nrebin, pyrsize // nrebin]) if qIdx == 0: stackedSubap = np.roll( mskRebin >= fracsub, tuple((-centers[qIdx] / nrebin).astype(np.int32)), (0, 1)) mskRebTot = (mskRebin >= fracsub).astype(np.int32) else: stackedSubap += np.roll( mskRebin >= fracsub, tuple((-centers[qIdx] / nrebin).astype(np.int32)), (0, 1)) mskRebTot += (mskRebin >= fracsub).astype(np.int32) * (qIdx + 1) validRow = [] validCol = [] for qIdx in [ 2, 1, 3, 0 ]: # Do not change order !! Corresponds to what happens A-B-C-D x-y // row-col tmpWh = np.where( np.roll(stackedSubap, tuple((centers[qIdx] / nrebin).astype(np.int32)), (0, 1))) validRow += [tmpWh[0].astype(np.int32)] validCol += [tmpWh[1].astype(np.int32)] nvalid = validRow[0].size validRow = np.asarray(validRow).flatten() validCol = np.asarray(validCol).flatten() p_wfs._nvalid = nvalid p_wfs._validsubsx = validCol p_wfs._validsubsy = validRow p_wfs._hrmap = mskRebTot.astype(np.int32) if (p_wfs.pyr_pos == None): pixsize = (np.pi * p_wfs._qpixsize) / (3600 * 180) # scale_fact = 2 * np.pi / npup * \ # (p_wfs.Lambda / p_tel.diam / 4.848) / pixsize * p_wfs.pyr_ampl # Proposition de Flo scale_fact = 2 * np.pi / pyrsize * \ (p_wfs.Lambda * 1e-6 / p_tel.diam) / pixsize * p_wfs.pyr_ampl cx = scale_fact * \ np.sin((np.arange(p_wfs.pyr_npts)) * 2. * np.pi / p_wfs.pyr_npts) cy = scale_fact * \ np.cos((np.arange(p_wfs.pyr_npts)) * 2. * np.pi / p_wfs.pyr_npts) # mod_npts = p_wfs.pyr_npts #UNUSED else: if (verbose): print("Using user-defined positions for the pyramid modulation") cx = p_wfs.pyr_pos[:, 0] / p_wfs._qpixsize cy = p_wfs.pyr_pos[:, 1] / p_wfs._qpixsize # mod_npts=cx.shape[0] #UNUSED p_wfs._pyr_cx = cx.copy() p_wfs._pyr_cy = cy.copy() telSurf = np.pi / 4. * (1 - p_tel.cobs**2.) * p_tel.diam**2. p_wfs._nphotons = p_wfs.zerop * \ 10. ** (-0.4 * p_wfs.gsmag) * ittime * \ p_wfs.optthroughput * telSurf # spatial filtering by the pixel extent: # *2/2 intended. min should be 0.40 = sinc(0.5)^2. y = np.tile(np.arange(pyrsize) - pyrsize // 2, (pyrsize, 1)) x = y.T x = x * 1. / pyrsize y = y * 1. / pyrsize sincar = np.fft.fftshift(np.sinc(x) * np.sinc(y)) # sincar = np.roll(np.pi*x*np.pi*y,x.shape[1],axis=1) p_wfs._sincar = sincar.astype(np.float32) #pup = p_geom._mpupil pupreb = util.rebin(pup * 1., [pyrsize // nrebin, pyrsize // nrebin]) a = pyrsize // nrebin b = p_geom._n // nrebin pupreb = pupreb[a // 2 - b // 2:a // 2 + b // 2, a // 2 - b // 2:a // 2 + b // 2] wsubok = np.where(pupreb >= p_wfs.fracsub) pupvalid = pupreb * 0. pupvalid[wsubok] = 1 p_wfs._isvalid = pupvalid.T.astype(np.int32) validsubsx = np.where(pupvalid)[0].astype(np.int32) validsubsy = np.where(pupvalid)[1].astype(np.int32) istart = np.arange(p_wfs.nxsub + 2) * p_wfs.npix jstart = np.copy(istart) p_wfs._istart = istart.astype(np.int32) p_wfs._jstart = jstart.astype(np.int32) # sorting out valid subaps fluxPerSub = np.zeros((p_wfs.nxsub + 2, p_wfs.nxsub + 2), dtype=np.float32) for i in range(p_wfs.nxsub + 2): indi = istart[i] for j in range(p_wfs.nxsub + 2): indj = jstart[j] fluxPerSub[i, j] = np.sum(p_geom._mpupil[indi:indi + p_wfs.npix, indj:indj + p_wfs.npix]) fluxPerSub = fluxPerSub / p_wfs.nxsub**2. p_wfs._fluxPerSub = fluxPerSub.copy() phasemap = np.zeros((p_wfs.npix * p_wfs.npix, p_wfs._nvalid), dtype=np.int32) X = np.indices((p_geom._n, p_geom._n)) # we need c-like indice tmp = X[1] + X[0] * p_geom._n pyrtmp = np.zeros((p_geom._n, p_geom._n), dtype=np.int32) for i in range(len(validsubsx)): indi = istart[validsubsy[i]] # +2-1 (yorick->python indj = jstart[validsubsx[i]] phasemap[:, i] = tmp[indi:indi + p_wfs.npix, indj:indj + p_wfs.npix].flatten("C") pyrtmp[indi:indi + p_wfs.npix, indj:indj + p_wfs.npix] = i p_wfs._phasemap = phasemap p_wfs._pyr_offsets = pyrtmp # pshift