Circular Functions
Zernike Modes
- aotools.functions.zernike.makegammas(nzrad)[source]
 Make “Gamma” matrices which can be used to determine first derivative of Zernike matrices (Noll 1976).
- Parameters:
 nzrad – Number of Zernike radial orders to calculate Gamma matrices for
- Returns:
 Array with x, then y gamma matrices
- Return type:
 ndarray
- aotools.functions.zernike.phaseFromZernikes(zCoeffs, size, norm='noll', rot=0)[source]
 Creates an array of the sum of zernike polynomials with specified coefficeints
- Parameters:
 - Returns:
 a size x size array of summed Zernike polynomials
- Return type:
 ndarray
- aotools.functions.zernike.zernIndex(j)[source]
 Find the [n,m] list giving the radial order n and azimuthal order of the Zernike polynomial of Noll index j.
- aotools.functions.zernike.zernikeArray(J, N, norm='noll', rot=0)[source]
 Creates an array of Zernike Polynomials
- Parameters:
 maxJ (int or list) – Max Zernike polynomial to create, or list of zernikes J indices to create
N (int) – size of created arrays
norm (string, optional) – The normalisation of Zernike modes. Can be
"noll","p2v"(peak to valley), or"rms". default is"noll".rot (float) – Rotates the zernike modes by rot radians around their centre. Defaults to 0.
- Returns:
 array of Zernike Polynomials
- Return type:
 ndarray
- aotools.functions.zernike.zernikeRadialFunc(n, m, r)[source]
 Fucntion to calculate the Zernike radial function
- aotools.functions.zernike.zernike_nm(n, m, N, rot=0)[source]
 Creates the Zernike polynomial with radial index, n, and azimuthal index, m.
Karhunen Loeve Modes
Collection of routines related to Karhunen-Loeve modes.
The theory is based on the paper “Optimal bases for wave-front simulation and reconstruction on annular apertures”, Robert C. Cannon, 1996, JOSAA, 13, 4
The present implementation is based on the IDL package of R. Cannon (wavefront modelling and reconstruction). A closely similar implementation can also be find in Yorick in the YAO package.
Usage
Main routine is ‘make_kl’ to generate KL basis of dimension [dim, dim, nmax].
For Kolmogorov statistics, e.g.
kl, _, _, _ = make_kl(150, 128, ri = 0.2, stf='kolmogorov')
Warning
make_kl with von Karman stf fails. It has been implemented but KL generation failed in ‘while loop’ of gkl_fcom…
- date:
 November 2017
- aotools.functions.karhunenLoeve.gkl_azimuthal(nord, npp)[source]
 Compute the azimuthal function of the KL basis.
- aotools.functions.karhunenLoeve.gkl_basis(ri=0.25, nr=40, npp=None, nfunc=500, stf='kolstf', outerscale=None)[source]
 Wrapper to create the radial and azimuthal K-L functions.
- Parameters:
 - Returns:
 gklbasis – dictionary containing the radial and azimuthal basis + other relevant information
- Return type:
 dic
- aotools.functions.karhunenLoeve.gkl_fcom(ri, kernels, nfunc, verbose=False)[source]
 Computation of the radial eigenvectors of the KL basis.
Obtained by taking the eigenvectors from the matrix L^p. The final function corresponds to the ‘nfunc’ largest eigenvalues. See eq. 16-18 in Cannon 1996.
- Parameters:
 - Returns:
 evals (1darray) – eigenvalues
nord (int) – resulting number of azimuthal orders
npo (int)
ord (1darray)
rabas (ndarray) – radial eigenvectors of the KL basis
- aotools.functions.karhunenLoeve.gkl_kernel(ri, nr, rad, stfunc='kolmogorov', outerscale=None)[source]
 Calculation of the kernel L^p
The kernel constructed here should be simply a discretization of the continuous kernel. It needs rescaling before it is treated as a matrix for finding the eigen-values
- Parameters:
 ri (float) – radius of central obscuration radius (normalized by D/2; <1)
nr (int) – number of resolution elements
rad (1d-array) – grid of points where the kernel is evaluated
stfunc (string) – string tag of the structure function on which the kernel are computed
outerscale (float) – in unit of telescope diameter. Outer-scale for von Karman structure function.
- Returns:
 L^p – kernel L^p of dimension (nr, nr, nth), where nth is the azimuthal discretization (5 * nr)
- Return type:
 ndarray
- aotools.functions.karhunenLoeve.gkl_radii(ri, nr)[source]
 Generate n points evenly spaced in r^2 between r_i^2 and 1
- aotools.functions.karhunenLoeve.gkl_sfi(kl_basis, i)[source]
 return the i’th function from the generalized KL basis ‘bas’. ‘bas’ must be generated by ‘gkl_basis’
- aotools.functions.karhunenLoeve.make_kl(nmax, dim, ri=0.0, nr=40, stf='kolmogorov', outerscale=None, mask=True)[source]
 Main routine to generatre a KL basis of dimension [nmax, dim, dim].
For Kolmogorov statistics, e.g.
kl, _, _, _ = make_kl(150, 128, ri = 0.2, stf='kolmogorov')
As a rule of thumbnr x npp = 50 x 250 is fine up to 500 functions60 x 300 for a thousand80 x 400 for three thousands.- Parameters:
 nmax (int) – number of KL function to generate
dim (int) – size of the KL arrays
ri (float) – radial central obscuration normalized by D/2
nr (int) – number of point on radius. npp (number of azimuthal pts) is =2 pi * nr
stf (string) – structure function tag. Default is ‘kolmogorov’
outerscale (float) – outer scale in units of telescope diameter. Releveant if von vonKarman stf. (not implemented yet)
mask (bool) – pupil masking. Default is True
- Returns:
 kl (ndarray (nmax, dim, dim)) – KL basis in cartesian coordinates
varKL (1darray) – associated variance
pupil (2darray) – pupil
polar_base (dictionary) – polar base dictionary used for the KL basis computation. as returned by ‘gkl_basis’
- aotools.functions.karhunenLoeve.pcgeom(nr, npp, ncp, ri, ncmar)[source]
 This routine builds a geom dic.
- Parameters:
 px (int) – the x, y coordinates of points in the polar arrays.
py (int) – the x, y coordinates of points in the polar arrays.
cr – the r, phi coordinates of points in the cartesian grids.
cp – the r, phi coordinates of points in the cartesian grids.
ncmar – allows the possibility that there is a margin of ncmar points in the cartesian arrays outside the region of interest.
- aotools.functions.karhunenLoeve.piston_orth(nr)[source]
 Unitary matrix used to filter out piston term. Eq. 19 in Cannon 1996.
- Parameters:
 nr (int) – number of resolution elements
- Returns:
 U – Unitary matrix
- Return type:
 2d array
- aotools.functions.karhunenLoeve.pol2car(cpgeom, pol, mask=False)[source]
 Polar to cartesian conversion.
(points not in the aperture are treated as though they were at the first or last radial polar value…)
- aotools.functions.karhunenLoeve.polang(r)[source]
 Generate an array with the same dimensions as r, but containing the azimuthal values for a polar coordinate system.
- aotools.functions.karhunenLoeve.radii(nr, npp, ri)[source]
 Use to generate a polar coordinate system.
Generate an nr x npp array with npp copies of the radial coordinate array. Radial coordinate span the range from r=ri to r=1 with successive annuli having equal areas.
ie, the area between ri and 1 is divided into nr equal rings, and the points are positioned at the half-area mark on each ring; there are no points on the border
See also
- aotools.functions.karhunenLoeve.rebin(a, newshape)[source]
 Rebin an array to a new shape. See scipy cookbook. It is intended to be similar to ‘rebin’ of IDL
- aotools.functions.karhunenLoeve.set_pctr(bas, ncp=None, ncmar=None)[source]
 call pcgeom to build the dic geom_struct with the right initializtions.
- Parameters:
 bas (dic) – gkl_basis dic built with the gkl_bas routine
- aotools.functions.karhunenLoeve.setpincs(ax, ay, px, py, ri)[source]
 determine a set of squares for interpolating from cartesian to polar coordinates, using only those points with ri <= r <= 1
- aotools.functions.karhunenLoeve.stf_kolmogorov(r)[source]
 Kolmogorov structure function with r = (D/r0)