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

Functions for modal basis (DM basis, KL, Btt, etc...) More...

Functions

def compute_KL2V (conf.Param_controller p_controller, Dms dms, list p_dms, conf.Param_geom p_geom, conf.Param_atmos p_atmos, conf.Param_tel p_tel)
 Compute the Karhunen-Loeve to Volt matrix (transfer matrix between the KL space and volt space for a pzt dm) More...
 
def compute_dm_basis (g_dm, conf.Param_dm p_dm, conf.Param_geom p_geom)
 Compute a the DM basis as a sparse matrix : More...
 
def compute_IFsparse (Dms g_dm, list p_dms, conf.Param_geom p_geom)
 Compute the influence functions of all DMs as a sparse matrix : More...
 
def command_on_Btt (Rtc rtc, Dms dms, list p_dms, conf.Param_geom p_geom, int nfilt)
 Compute a command matrix in Btt modal basis (see error breakdown) and set it on the sutra_rtc. More...
 
def compute_cmat_with_Btt (Rtc rtc, np.ndarray Btt, int nfilt)
 Compute a command matrix on the Btt basis and load it in the GPU. More...
 
def command_on_KL (Rtc rtc, Dms dms, conf.Param_controller p_controller, List[conf.Param_dm] p_dms, conf.Param_geom p_geom, conf.Param_atmos p_atmos, conf.Param_tel p_tel, int nfilt)
 Compute a command matrix in KL modal basis and set it on the sutra_rtc. More...
 
def compute_cmat_with_KL (Rtc rtc, np.ndarray KL2V, int nfilt)
 Compute a command matrix on the KL basis and load it in the GPU. More...
 
def compute_fourier (int nActu, float pitch, np.ndarray actu_x_pos, np.ndarray actu_y_pos, periodic='n')
 Values you are looking for are: More...
 
def compute_btt (IFpzt, IFtt, influ_petal=None, return_delta=False)
 Returns Btt to Volts and Volts to Btt matrices. More...
 

Detailed Description

Functions for modal basis (DM basis, KL, Btt, etc...)

Author
COMPASS Team https://github.com/ANR-COMPASS
Version
5.4.4
Date
2022/01/24

This file is part of COMPASS https://anr-compass.github.io/compass/

Copyright (C) 2011-2023 COMPASS Team https://github.com/ANR-COMPASS All rights reserved. Distributed under GNU - LGPL

COMPASS is free software: you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3 of the License, or any later version.

COMPASS: End-to-end AO simulation tool using GPU acceleration The COMPASS platform was designed to meet the need of high-performance for the simulation of AO systems.

The final product includes a software package for simulating all the critical subcomponents of AO, particularly in the context of the ELT and a real-time core based on several control approaches, with performances consistent with its integration into an instrument. Taking advantage of the specific hardware architecture of the GPU, the COMPASS tool allows to achieve adequate execution speeds to conduct large simulation campaigns called to the ELT.

The COMPASS platform can be used to carry a wide variety of simulations to both testspecific components of AO of the E-ELT (such as wavefront analysis device with a pyramid or elongated Laser star), and various systems configurations such as multi-conjugate AO.

COMPASS is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.

You should have received a copy of the GNU Lesser General Public License along with COMPASS. If not, see https://www.gnu.org/licenses/lgpl-3.0.txt.

Function Documentation

◆ command_on_Btt()

def shesha.ao.basis.command_on_Btt ( Rtc  rtc,
Dms  dms,
list  p_dms,
conf.Param_geom  p_geom,
int  nfilt 
)

Compute a command matrix in Btt modal basis (see error breakdown) and set it on the sutra_rtc.

It computes by itself the volts to Btt matrix.

Parameters
rtc(Rtc) : rtc object
dms(Dms): dms object
p_dms(list of Param_dm): dms settings
p_geom(Param_geom): geometry settings
nfilt(int): number of modes to filter

Definition at line 205 of file basis.py.

Here is the call graph for this function:

◆ command_on_KL()

def shesha.ao.basis.command_on_KL ( Rtc  rtc,
Dms  dms,
conf.Param_controller  p_controller,
List[conf.Param_dm]  p_dms,
conf.Param_geom  p_geom,
conf.Param_atmos  p_atmos,
conf.Param_tel  p_tel,
int  nfilt 
)

Compute a command matrix in KL modal basis and set it on the sutra_rtc.

It computes by itself the volts to KL matrix.

Parameters
rtc(Rtc) : rtc object
dms(Dms): dms object
p_dms(list of Param_dm): dms settings
p_geom(Param_geom): geometry settings
p_atmos(Param_atmos) : atmos parameters
p_tel(Param_tel) : telescope parameters
nfilt(int): number of modes to filter

Definition at line 264 of file basis.py.

Here is the call graph for this function:

◆ compute_btt()

def shesha.ao.basis.compute_btt (   IFpzt,
  IFtt,
  influ_petal = None,
  return_delta = False 
)

Returns Btt to Volts and Volts to Btt matrices.

Parameters
IFpzt(csr_matrix) : influence function matrix of pzt DM, sparse and arrange as (Npts in pup x nactus)
IFtt(np.ndarray(ndim=2,dtype=np.float32)) : Influence function matrix of the TT mirror arrange as (Npts in pup x 2)
influ_petal(np.ndarray) : Influence function matrix of M4 petals. Default is None, if set, the Btt produced is also orthogonal to petal modes, then only driven by petal DM
return_delta(bool, optional) : If True, returns delta instead of P. Default is False

:returns:

Parameters
Btt(np.ndarray(ndim=2,dtype=np.float32)) : Btt to Volts matrix
P(np.ndarray(ndim=2,dtype=np.float32)) : Volts to Btt matrix

Definition at line 370 of file basis.py.

Here is the caller graph for this function:

◆ compute_cmat_with_Btt()

def shesha.ao.basis.compute_cmat_with_Btt ( Rtc  rtc,
np.ndarray  Btt,
int  nfilt 
)

Compute a command matrix on the Btt basis and load it in the GPU.

Parameters
rtc(Rtc): rtc object
Btt(np.ndarray[ndim=2, dtype=np.float32]) : volts to Btt matrix
nfilt(int): number of modes to filter

Definition at line 226 of file basis.py.

Here is the caller graph for this function:

◆ compute_cmat_with_KL()

def shesha.ao.basis.compute_cmat_with_KL ( Rtc  rtc,
np.ndarray  KL2V,
int  nfilt 
)

Compute a command matrix on the KL basis and load it in the GPU.

Parameters
rtc(Rtc): rtc object
KL2V(np.ndarray[ndim=2, dtype=np.float32]) : volts to KL matrix
nfilt(int): number of modes to filter

Definition at line 281 of file basis.py.

Here is the caller graph for this function:

◆ compute_dm_basis()

def shesha.ao.basis.compute_dm_basis (   g_dm,
conf.Param_dm  p_dm,
conf.Param_geom  p_geom 
)

Compute a the DM basis as a sparse matrix :

  • push on each actuator
  • get the corresponding dm shape
  • apply pupil mask and store in a column
Parameters
g_dm(Dm) : Dm object
p_dm(Param_dm) : dm settings
p_geom(Param_geom) : geom settings

:return:

IFbasis = (csr_matrix) : DM IF basis

Definition at line 131 of file basis.py.

Here is the caller graph for this function:

◆ compute_fourier()

def shesha.ao.basis.compute_fourier ( int  nActu,
float  pitch,
np.ndarray  actu_x_pos,
np.ndarray  actu_y_pos,
  periodic = 'n' 
)

Values you are looking for are:

config.p_dms[0].nact
config.p_dms[0]._pitch
config.p_dms[0]._i1
config.p_dms[0]._j1

Definition at line 307 of file basis.py.

◆ compute_IFsparse()

def shesha.ao.basis.compute_IFsparse ( Dms  g_dm,
list  p_dms,
conf.Param_geom  p_geom 
)

Compute the influence functions of all DMs as a sparse matrix :

  • push on each actuator
  • get the corresponding dm shape
  • apply pupil mask and store in a column
Parameters
g_dm(Dms) : Dms object
p_dms(Param_dms) : dms settings
p_geom(Param_geom) : geom settings

:return:

IFbasis = (csr_matrix) : DM IF basis

Definition at line 174 of file basis.py.

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

◆ compute_KL2V()

def shesha.ao.basis.compute_KL2V ( conf.Param_controller  p_controller,
Dms  dms,
list  p_dms,
conf.Param_geom  p_geom,
conf.Param_atmos  p_atmos,
conf.Param_tel  p_tel 
)

Compute the Karhunen-Loeve to Volt matrix (transfer matrix between the KL space and volt space for a pzt dm)

Parameters
p_controller(Param_controller) : p_controller settings
dms(shesha_dms) : Dms object
p_dms(list of Param_dm) : dms settings
p_geom(Param_geom) : geometry parameters
p_atmos(Param_atmos) : atmos parameters
p_tel(Param_tel) : telescope parameters

:return:

Parameters
KL2V(np.array(np.float32,dim=2)) : KL to Volt matrix

Definition at line 72 of file basis.py.

Here is the caller graph for this function: