COMPASS  5.4.4
End-to-end AO simulation tool using GPU acceleration
shesha.util.slopesCovariance Namespace Reference

Functions

def KLmodes (x, y, L0, filterPiston)
 Return modes <x, y> : x,y coordinates of the DM actuators. More...
 
def computePistonFilteringMatrix (n)
 Creates a sqaure matrix that removes (=filter out) piston mode. More...
 
def compute_Caa (vv1, vv2, s1, s2, Cn2h, L0, alt1, alt2, altitude)
 This function computes the covariance matrix between TWO wavefront-sensors, each pointing in some particular direction of a given Guide Star for a UNI turbulent layer at a given altitude. More...
 
def macdo_x56 (x, k=10)
 Computation of the Mc Donald function. More...
 
def asymp_macdo (x)
 Computes a term involved in the computation of the phase struct function with a finite outer scale according to the Von-Karman model. More...
 
def rodconan (r, L0, k=10)
 DOCUMENT rodconan(r,L0,k=) The phase structure function is computed from the expression Dphi(r) = k1 * L0**(5. More...
 
def DPHI (x, y, L0)
 
def computeTTFiltMatrix (nbslopes)
 Returns a matrix <p1> that filters the average TT on the vector of slopes <vec> with nbslopes = x+y slopes. More...
 
def filt_TTvec (nbslopes)
 Returns a matrix p that averages the slopes over x and y. More...
 
def Filt_Zer (modesfilt, miz)
 

Function Documentation

◆ asymp_macdo()

def shesha.util.slopesCovariance.asymp_macdo (   x)

Computes a term involved in the computation of the phase struct function with a finite outer scale according to the Von-Karman model.

The term involves the MacDonald function (modified bessel function of second kind) K_{5/6}(x), and the algorithm uses the asymptotic form for x ~ infinity.

Warnings
  • This function makes a floating point interrupt for x=0 and should not be used in this case.
  • Works only for x>0.

Definition at line 313 of file slopesCovariance.py.

Here is the caller graph for this function:

◆ compute_Caa()

def shesha.util.slopesCovariance.compute_Caa (   vv1,
  vv2,
  s1,
  s2,
  Cn2h,
  L0,
  alt1,
  alt2,
  altitude 
)

This function computes the covariance matrix between TWO wavefront-sensors, each pointing in some particular direction of a given Guide Star for a UNI turbulent layer at a given altitude.


— Warning! —

The simulated guide star —IS ALWAYS— TT sensitive i.e if the GS is set at a finite altitude, it simulates a LGS with the perfect TT sensing capability. If needed, the TT filtering should be performed afterwards on the covariance matrix by applying an averaging slopes matrix for the given WFS (see computeTTFiltMatrix for TT filtering or Filt_Zer for filtering other Zernikes modes)

<vv1>, <vv2> is the array (N,2) of the subaperture coordinates projected onto altitude layer. vv[:,0] is the x-coordinate, vv([:,1] is y. vv1 is for WFS number 1, vv2 for the other. It's up to the user to order the subapertures in the way he expects to find them in the final cov matrix. <s1>, <s2> are the sizes of the subapertures in the pupil (i.e. at ground level, i.e. the real physical size in the telescope pupil, not the one projected onto the layer) <Cn2h> strengh of the layer: it must be equal to r0**(-5./3) with r0 in metres. Initially this parameter was r0 but the Learn & Apply fit was non-linear. Now the Learn will be searching for <Cn2h> which makes it a linear fit and convergence is faster. <L0> the outer scale (in metres) <alt1>, <alt2>: altitude of the source (the Guide Star) of the WFS (metres). A value of alt==0 is interpreted as alt==infinity. <altitude>: altitude of the layer (metres). Used to recompute the real ssp size at the given altitude.


Example of use:

# Generating Arbitrary ssp coordinates with a 14x14 SH on a 8m telescope:
Nssp = 14
x = np.linspace(-1,1,Nssp)
x, y = np.meshgrid(x,x)
r=np.sqrt(x*x+y*y)<1.1
vv1 = np.array([x[r].flatten(), y[r].flatten()]).T
vv1*=4
# We shift the second WFS by 0.5m
vv2 =vv1 + 0.5
s1=0.5 # Physical size of the ssp in meters of wfs #1
s2=0.5 # Physical size of the ssp in meters of wfs #2
r0 = 0.5 # in meters
Cn2h = r0**(5/3)
L0 = 50 # in meters
alt1 = 0 # GS altitude 0 meters == infinite GS == NGS
alt2 = 0 # GS altitude 0 meters == infinite GS == NGS
altitude = 3500 # Altitude of the turbulence layer in meters
caa = compute_Caa(vv1, vv2, s1, s2, Cn2h, L0, alt1, alt2, altitude) # Computes the WFS covariance matrix
# --------------------------------------------
# Optional modes filtering:
# --------------------------------------------
nbslopes = caa.shape[0] # x +y slopes = total nb of slopes
FiltMatrix = computeTTFiltMatrix(nbslopes) # ex with TT filtering (use Filt_Zer for other modes filtering)
# Case of Covmat of NGS/LGS:
caaFilt = np.dot(filtTTMat, caa)
# Case of Covmat of LGS/LGS:
caaFilt = np.dot(np.dot(filtTTMat, caa), filtTTMat.T)

Definition at line 176 of file slopesCovariance.py.

Here is the call graph for this function:

◆ computePistonFilteringMatrix()

def shesha.util.slopesCovariance.computePistonFilteringMatrix (   n)

Creates a sqaure matrix that removes (=filter out) piston mode.

<n> : size of the matrix

P = computePistonFilteringMatrix(64)

Definition at line 96 of file slopesCovariance.py.

Here is the caller graph for this function:

◆ computeTTFiltMatrix()

def shesha.util.slopesCovariance.computeTTFiltMatrix (   nbslopes)

Returns a matrix <p1> that filters the average TT on the vector of slopes <vec> with nbslopes = x+y slopes.


Example of use:

nbslopes = 300
p1 = computeTTFiltMatrix(nbslopes)
vec_filt = np.dot(p1,vec)

Definition at line 400 of file slopesCovariance.py.

Here is the call graph for this function:

◆ DPHI()

def shesha.util.slopesCovariance.DPHI (   x,
  y,
  L0 
)
dphi = DPHI(x,y,L0) * r0**(-5./3)

Computes the phase structure function for a separation (x,y). The r0 is not taken into account : the final result of DPHI(x,y,L0) has to be scaled with r0**-5/3, with r0 expressed in meters, to get the right value.

Definition at line 370 of file slopesCovariance.py.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ filt_TTvec()

def shesha.util.slopesCovariance.filt_TTvec (   nbslopes)

Returns a matrix p that averages the slopes over x and y.


Example of use:

nbslopes = 300
p = filt_TTvec(nbslopes)
vec_filtered = np.dot(p,vec)

Definition at line 421 of file slopesCovariance.py.

Here is the caller graph for this function:

◆ Filt_Zer()

def shesha.util.slopesCovariance.Filt_Zer (   modesfilt,
  miz 
)
MatFiltZer = Filt_Zer(modesfilt, miz);

This function generates a matrix that filters all the zernikes modes up to <modesfilt> It requires a zernike interaction matrix <miz> of the given the WFS


Example of use:

modesfilt = 10
MatFiltZer = Filt_Zer(modesfilt, miz)

The filtered vector can be computed using:

vec_filt = np.dot(MatFiltZer,miz)

With <miz> computed using doit:

from hraa.tools.doit import doit
nssp = 14
obs = 0
nbzer = 100
miz = doit(nssp, obs, nbzer)

Definition at line 461 of file slopesCovariance.py.

◆ KLmodes()

def shesha.util.slopesCovariance.KLmodes (   x,
  y,
  L0,
  filterPiston 
)

Return modes <x, y> : x,y coordinates of the DM actuators.

Could come from the output of the function:

x, y, mask = generateActuCoordinates(nx, ny, Nactu)

OR can be any x,y DM coordinates.

<L0> : outer scale, in the same units as x and y. <filterPiston> : True/False, to filter the piston out or not

example of use:
L0 = 25
nx, ny = 64, 64
Nactu = 3228
x, y, mask = generateActuCoordinates(nx, ny, Nactu)
modes = KLmodes(x, y, L0, True)

Definition at line 61 of file slopesCovariance.py.

Here is the call graph for this function:

◆ macdo_x56()

def shesha.util.slopesCovariance.macdo_x56 (   x,
  k = 10 
)

Computation of the Mc Donald function.

f(x) = x**(5/6)*K_{5/6}(x) using a series for the esimation of K_{5/6}, taken from Rod Conan thesis : K_a(x)=1/2 \sum_{n=0}**\infty \frac{(-1)**n}{n!} \left(\Gamma(-n-a) (x/2)**{2n+a} + \Gamma(-n+a) (x/2)**{2n-a} \right) , with a = 5/6.

Setting x22 = (x/2)**2, setting uda = (1/2)**a, and multiplying by x**a, this becomes : x**a * Ka(x) = 0.5 $ -1**n / n! [ G(-n-a).uda x22**(n+a) + G(-n+a)/uda x22**n ] Then we use the following recurrence formulae on the following quantities : G(-(n+1)-a) = G(-n-a) / -a-n-1 G(-(n+1)+a) = G(-n+a) / a-n-1 (n+1)! = n! * (n+1) x22**(n+1) = x22**n * x22 and at each iteration on n, one will use the values already computed at step (n-1). The values of G(a) and G(-a) are hardcoded instead of being computed.

The first term of the series has also been skipped, as it vanishes with another term in the expression of Dphi.

Definition at line 270 of file slopesCovariance.py.

Here is the caller graph for this function:

◆ rodconan()

def shesha.util.slopesCovariance.rodconan (   r,
  L0,
  k = 10 
)

DOCUMENT rodconan(r,L0,k=) The phase structure function is computed from the expression Dphi(r) = k1 * L0**(5.

/3) * (k2 - (2.pi.r/L0)**5/6 K_{5/6}(2.pi.r/L0))

For small r, the expression is computed from a development of K_5/6 near 0. The value of k2 is not used, as this same value appears in the series and cancels with k2. For large r, the expression is taken from an asymptotic form.

Definition at line 339 of file slopesCovariance.py.

Here is the call graph for this function:
Here is the caller graph for this function: