MAOS
Multithreaded Adaptive Optics Simulator
|
Functions | |
recon_t * | setup_recon_prep (const parms_t *parms, const aper_t *aper, const powfs_t *powfs) |
void | setup_recon_prep_ga (recon_t *recon, const parms_t *parms, const aper_t *aper, const powfs_t *powfs) |
void | setup_recon_saneai (recon_t *recon, const parms_t *parms, const powfs_t *powfs) |
void | setup_recon_fit (recon_t *recon, const parms_t *parms) |
void | setup_powfs_fit (powfs_t *powfs, const recon_t *recon, const parms_t *parms) |
void | free_powfs_fit (powfs_t *powfs, const parms_t *parms) |
void | free_fit (fit_t *fit, int nfit) |
void | setup_recon_tomo_prep (recon_t *recon, const parms_t *parms) |
void | setup_recon_tomo_matrix (recon_t *recon, const parms_t *parms) |
void | setup_recon_update_cn2 (recon_t *recon, const parms_t *parms) |
void | setup_recon_tomo (recon_t *recon, const parms_t *parms, const powfs_t *powfs) |
void | setup_recon_lsr (recon_t *recon, const parms_t *parms) |
void | setup_recon_mvm (const parms_t *parms, recon_t *recon, const powfs_t *powfs) |
void | setup_recon_control (recon_t *recon, const parms_t *parms, const powfs_t *powfs) |
void | setup_recon_post (recon_t *recon, const parms_t *parms, const aper_t *aper, const powfs_t *powfs) |
void | setup_recon_HXW (recon_t *recon, const parms_t *parms, mapcell *atm) |
void | free_recon (const parms_t *parms, recon_t *recon) |
void | free_recon_unused (const parms_t *parms, recon_t *recon) |
void | setup_recon_mvst (recon_t *recon, const parms_t *parms) |
void | setup_recon_dither_dm (recon_t *recon, const powfs_t *powfs, const parms_t *parms) |
void | tomofit (dcell **dmout, sim_t *simu, dcell *gradin) |
void * | reconstruct (sim_t *simu) |
Carry out wavefront reconstruction and DM fitting.
Create reconstruction parameters that are related to the geometry only, and will not be updated when estimated WFS measurement noise changes.
This can be used to do NCPA calibration.
void setup_recon_prep_ga | ( | recon_t * | recon, |
const parms_t * | parms, | ||
const aper_t * | aper, | ||
const powfs_t * | powfs | ||
) |
That may depend on GPU data.
Setup the matrix of the inverse of gradient measurement noise equivalent angle covariance matrix. For physical optics wfs, the NEA is computed using matched filter output. For geometric optics, the NEA is from input.
Setting fitting parameter for turbulence to WFS lenslet grid.
Call this instead of free_fit directly to free powfs fit parameters.
Prepares for tomography. ALlow it to be called multiple times for Cn2 update.
assemble tomography matrix. In CG mode, this function is not executed if tomo.assemble=0, Instead, the algorithm is contained in recon.c. When you modify anything, make sure you also do it there.
For integrated tomograhy:
\(\hat{x}=(G_{lgs}^{T}C_{lgs}^{-1}G_{lgs}+C_{x}^{-1}+G_{ngs}^{T}C_{ngs}^{-1} G_{ngs})^{-1}(G_{lgs}^{T}C_{lgs}^{-1}s_{lgs}+G_{ngs}^{T}C_{ngs}^{-1}s_{ngs}){\equiv}RL^{-1}RR s.\)
For split tomography, the terms regarding NGS are dropped:
\(\hat{x}_{lgs}=(G_{lgs}^{T}C_{lgs}^{-1}G_{lgs}+C_{x}^{-1})^{-1} G_{lgs}^{T}C_{lgs}^{-1}s_{lgs}\equiv R_L^{-1} R_R s.\)
The left hand side of linear equation (inside the inverse) is stored in RL. The right hand side of the linear equation (outside of the inverse) is stored in RR. The terms regarding the NGS are handled using low rank terms. The gradients from LGS and NGS and concatenated to \(s\).
In the LGS part, there is a global tip/tilt removal operator because LGS is insensitive to global tip/tilts.
For details see www.opticsinfobase.org/abstract.cfm?URI=josaa-19-9-1803
Update assembled tomography matrix with new L2. Called from cn2est when new profiles are available.
Sets up the tomogrpahy turbulence reconstruction structs including wavefront reconstructor and DM fitting operator
Calls setup_recon_tomo_matrix() and setup_recon_fit_matrix() to setup the tomography and DM fitting matrix.
AHST is handled in setup_ngsmod().
MVST is handled in setup_recon_mvst().
MOAO is handled in setup_recon_moao().
Setup the least square reconstruct by directly inverting GA matrix. The reconstructor is simply the pseudo inverse of GA matrix:
\[\hat{x}=(G_a^TC_g^{-1}G_a)^{-1}G_a^TC_g^{-1}\]
This is very close to RR except replacing GX with GA.
We use the tomograhy parameters for lsr, since lsr is simply "tomography" onto DM directly.
Assemble matrix to do matrix vector multiply. Split from setup_recon because GPU may be used.
Setup either the control matrix using either minimum variance reconstructor by calling setup_recon_mvr() or least square reconstructor by calling setup_recon_lsr() It can be called repeatedly to update the reconstructor when NEA updates.
The results are measurement noise dependent and may be updated during simulation.
void setup_recon_post | ( | recon_t * | recon, |
const parms_t * | parms, | ||
const aper_t * | aper, | ||
const powfs_t * | powfs | ||
) |
A few further operations that needs MVM.
Setup ray tracing operator from xloc to ploc for guide stars.
recon | |
parms | |
atm | If set, atm[].vx,vy will be used for predictive offsetting. |
Free unused object in recon struct after preparation is done.
compute the MVST split tomography NGS mode reconstructor.
2010-03-16: New Implementation.
Definition:
We have
\begin{eqnarray*} \hat{x}_{lgs}&=&A^{-1}G^{T}_{lgs}C_{lgs}^{-1}s_{lgs}\\ Uw&=&A^{-1}G_{ngs}^TC_{ngs}^{-1}\\ \hat{x}_{ngs}&=&Uw(1+G_{ngs}Uw)^{-1}(s_{ngs}-G_{ngs}\hat{x}_{lgs});\\ a_{NGS}&=&F\hat{x}_{ngs}\\ MVModes&=&F\cdot Uw\\ MVRngs&=&(1+G_{ngs}Uw)^{-1} \end{eqnarray*}
If we want to orthnormalize the NGS modes. Propagate the NGS modes \(F\cdot Uw\) to fitting directions and compute the cross-coupling matrix with weighting on \(W_0, W_1\), we have
\begin{eqnarray*} Q&=&H_A\cdot F \cdot Uw\\ M_{CC}&=&Q^T(W_0-W_1 W_1^T)Q\\ \end{eqnarray*}
Do SVD on \(MCC\) we have \(M_{CC}=u\Sigma v^T \) where \(u{\equiv}v\) because \(MCC\) is symmetric. Redefine the NGS modes and reconstructor as
\begin{eqnarray*} MVModes&=&F\cdot Uw\cdot u \Sigma^{-1/2}\\ MVRngs&=&\Sigma^{1/2}u^T (1+G_{ngs}A^{-1}G_{ngs}^T C_{ngs}^{-1})^{-1} \end{eqnarray*}
Dither using command path (DM) aberration
Calls tomo() and fit() to do the tomography and DM fit. Do error signal and split tomography. In closedloop mode, the gradients are from time step isim-1.