![]() |
COMPASS
5.0.0
End-to-end AO simulation tool using GPU acceleration
|
Initialization of a Dms object. More...
Functions | |
Dms | dm_init (carmaWrap_context context, List[conf.Param_dm] p_dms, conf.Param_tel p_tel, conf.Param_geom p_geom, List[conf.Param_wfs] p_wfss=None, bool keepAllActu=False) |
Create and initialize a Dms object on the gpu. More... | |
def | dm_init_standalone (carmaWrap_context context, list p_dms, conf.Param_geom p_geom, diam=1., cobs=0., pupAngle=0., wfs_xpos=[0], wfs_ypos=[0]) |
Create and initialize a Dms object on the gpu. More... | |
def | make_pzt_dm (conf.Param_dm p_dm, conf.Param_geom p_geom, float cobs, float pupAngle, bool keepAllActu=False) |
Compute the actuators positions and the influence functions for a pzt DM. More... | |
def | init_custom_dm (conf.Param_dm p_dm, conf.Param_geom p_geom, float diam) |
Read Fits for influence pzt fonction and form. More... | |
def | make_tiptilt_dm (conf.Param_dm p_dm, int patchDiam, conf.Param_geom p_geom, float diam) |
Compute the influence functions for a tip-tilt DM. More... | |
None | make_kl_dm (conf.Param_dm p_dm, int patchDiam, conf.Param_geom p_geom, float cobs) |
Compute the influence function for a Karhunen-Loeve DM. More... | |
def | comp_dmgeom (conf.Param_dm p_dm, conf.Param_geom p_geom) |
Compute the geometry of a DM : positions of actuators and influence functions. More... | |
def | correct_dm (context, Dms dms, list p_dms, conf.Param_controller p_controller, conf.Param_geom p_geom, np.ndarray imat=None, dict dataBase={}, bool use_DB=False) |
Correct the geometry of the DMs using the imat (filter unseen actuators) More... | |
def | makePetalDm (p_dm, p_geom, pupAngleDegree) |
makePetalDm(p_dm, p_geom, pupAngleDegree) More... | |
def | make_petal_dm_core (pupImage, pupAngleDegree) |
<pupImage> : image of the pupil More... | |
def | build_petals (nbSeg, pupAngleDegree, i0, j0, npt) |
Initialization of a Dms object.
This file is part of COMPASS https://anr-compass.github.io/compass/
Copyright (C) 2011-2019 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.
def shesha.init.dm_init.build_petals | ( | nbSeg, | |
pupAngleDegree, | |||
i0, | |||
j0, | |||
npt | |||
) |
def shesha.init.dm_init.comp_dmgeom | ( | conf.Param_dm | p_dm, |
conf.Param_geom | p_geom | ||
) |
Compute the geometry of a DM : positions of actuators and influence functions.
:parameters: dm (Param_dm) : dm settings
geom (Param_geom) : geom settings
Definition at line 739 of file dm_init.py.
def shesha.init.dm_init.correct_dm | ( | context, | |
Dms | dms, | ||
list | p_dms, | ||
conf.Param_controller | p_controller, | ||
conf.Param_geom | p_geom, | ||
np.ndarray | imat = None , |
||
dict | dataBase = {} , |
||
bool | use_DB = False |
||
) |
Correct the geometry of the DMs using the imat (filter unseen actuators)
:parameters: context (carmaWrap_context): context dms (Dms) : Dms object p_dms (list of Param_dm) : dms settings p_controller (Param_controller) : controller settings p_geom (Param_geom) : geom settings imat (np.ndarray) : interaction matrix dataBase (dict): dictionary containing paths to files to load use_DB (bool): dataBase use flag
Definition at line 811 of file dm_init.py.
Dms shesha.init.dm_init.dm_init | ( | carmaWrap_context | context, |
List[conf.Param_dm] | p_dms, | ||
conf.Param_tel | p_tel, | ||
conf.Param_geom | p_geom, | ||
List[conf.Param_wfs] | p_wfss = None , |
||
bool | keepAllActu = False |
||
) |
Create and initialize a Dms object on the gpu.
:parameters: context (carmaWrap_context): context p_dms (list of Param_dms) : dms settings p_tel (Param_tel) : telescope settings p_geom (Param_geom) : geom settings p_wfss (list of Param_wfs) : wfs settings :return: Dms (Dms): Dms object
Definition at line 68 of file dm_init.py.
def shesha.init.dm_init.dm_init_standalone | ( | carmaWrap_context | context, |
list | p_dms, | ||
conf.Param_geom | p_geom, | ||
diam = 1. , |
|||
cobs = 0. , |
|||
pupAngle = 0. , |
|||
wfs_xpos = [0] , |
|||
wfs_ypos = [0] |
|||
) |
Create and initialize a Dms object on the gpu.
:parameters: p_dms (list of Param_dms) : dms settings
p_geom (Param_geom) : geom settings
diam (float) : diameter of telescope (default 1.)
cobs (float) : cobs of telescope (default 0.)
pupAngle (float) : pupil rotation angle (degrees, default 0.)
wfs_xpos (array) : guide star x position on sky (in arcsec).
wfs_ypos (array) : guide star y position on sky (in arcsec).
Definition at line 324 of file dm_init.py.
def shesha.init.dm_init.init_custom_dm | ( | conf.Param_dm | p_dm, |
conf.Param_geom | p_geom, | ||
float | diam | ||
) |
Read Fits for influence pzt fonction and form.
:parameters: p_dm (Param_dm) : dm settings
p_geom (Param_geom) : geom settings
diam (float) : tel diameter
Conversion. There are sereval coordinate systems. Some are coming from the input fits file, others from compass. Those systems differ by scales and offsets. Coord systems from the input fits file:
Variables will be named using f, c: define either fits or compass
Definition at line 544 of file dm_init.py.
None shesha.init.dm_init.make_kl_dm | ( | conf.Param_dm | p_dm, |
int | patchDiam, | ||
conf.Param_geom | p_geom, | ||
float | cobs | ||
) |
Compute the influence function for a Karhunen-Loeve DM.
:parameters: p_dm (Param_dm) : dm settings
patchDiam (int) : patchDiam for dm size
p_geom (Param_geom) : geom settings
cobs (float) : telescope cobs
Definition at line 691 of file dm_init.py.
def shesha.init.dm_init.make_petal_dm_core | ( | pupImage, | |
pupAngleDegree | |||
) |
<pupImage> : image of the pupil
La fonction renvoie des fn d'influence en forme de petale d'apres une image de la pupille, qui est supposee etre segmentee.
influ, i1, j1, smallsize, nbSeg = make_petal_dm_core(pupImage, 0.0)
Definition at line 912 of file dm_init.py.
def shesha.init.dm_init.make_pzt_dm | ( | conf.Param_dm | p_dm, |
conf.Param_geom | p_geom, | ||
float | cobs, | ||
float | pupAngle, | ||
bool | keepAllActu = False |
||
) |
Compute the actuators positions and the influence functions for a pzt DM.
:parameters: p_dm (Param_dm) : dm parameters
p_geom (Param_geom) : geometry parameters
cobs (float) : telescope central obstruction
:return: influ (np.ndarray(dims=3, dtype=np.float64)) : cube of the IF for each actuator
Definition at line 349 of file dm_init.py.
def shesha.init.dm_init.make_tiptilt_dm | ( | conf.Param_dm | p_dm, |
int | patchDiam, | ||
conf.Param_geom | p_geom, | ||
float | diam | ||
) |
Compute the influence functions for a tip-tilt DM.
:parameters: p_dm (Param_dm) : dm settings
patchDiam (int) : patchDiam for dm size
p_geom (Param_geom) : geom settings
diam (float) : telescope diameter :return: influ (np.ndarray(dims=3,dtype=np.float64)) : cube of the IF
Definition at line 657 of file dm_init.py.
def shesha.init.dm_init.makePetalDm | ( | p_dm, | |
p_geom, | |||
pupAngleDegree | |||
) |
makePetalDm(p_dm, p_geom, pupAngleDegree)
The function builds a DM, segmented in petals according to the pupil shape. The petals will be adapted to the EELT case only.
<p_geom> : compass object p_geom. The function requires the object p_geom in order to know what is the pupil mask, and what is the mpupil. <p_dm> : compass petal dm object p_dm to be created. The function will transform/modify in place the attributes of the object p_dm.
Definition at line 887 of file dm_init.py.