mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
Loading...
Searching...
No Matches
aoSystem.hpp
Go to the documentation of this file.
1/** \file aoSystem.hpp
2 * \author Jared R. Males (jaredmales@gmail.com)
3 * \brief Declares and defines an analytical AO system
4 * \ingroup mxAO_files
5 *
6 */
7
8//***********************************************************************//
9// Copyright 2015-2022 Jared R. Males (jaredmales@gmail.com)
10//
11// This file is part of mxlib.
12//
13// mxlib is free software: you can redistribute it and/or modify
14// it under the terms of the GNU General Public License as published by
15// the Free Software Foundation, either version 3 of the License, or
16// (at your option) any later version.
17//
18// mxlib is distributed in the hope that it will be useful,
19// but WITHOUT ANY WARRANTY; without even the implied warranty of
20// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21// GNU General Public License for more details.
22//
23// You should have received a copy of the GNU General Public License
24// along with mxlib. If not, see <http://www.gnu.org/licenses/>.
25//***********************************************************************//
26
27#ifndef aoSystem_hpp
28#define aoSystem_hpp
29
30#include "../../mxlib.hpp"
31#include "../../math/constants.hpp"
32#include "../../math/roots.hpp"
33#include "../../mxError.hpp"
34#include "../../mxException.hpp"
35
36#include "aoConstants.hpp"
37
38#include "aoAtmosphere.hpp"
39#include "aoPSDs.hpp"
40#include "aoWFS.hpp"
41
42#define FITTING_ERROR_NO 0
43#define FITTING_ERROR_ZERO 1
44#define FITTING_ERROR_X 2
45#define FITTING_ERROR_Y 3
46
47namespace mx
48{
49namespace AO
50{
51namespace analysis
52{
53
54/// Describes an analytic adaptive optics (AO) system.
55/**
56 * Templatized by the turbulence power spectral density (PSD).
57 *
58 * \tparam realT the floating point type used for all calculations
59 * \tparam inputSpecT specifies the turbulence spatial PSD type
60 * \tparam iosT is an output stream type with operator << defined (default is std::ostream)
61 *
62 * \ingroup mxAOAnalytic
63 */
64template <typename _realT, class _inputSpectT, typename iosT = std::ostream>
66{
67
68 public:
69 typedef _realT realT;
70 typedef _inputSpectT inputSpectT;
72
73 aoAtmosphereT atm;
74 inputSpectT psd;
75
76 protected:
77 realT m_D{ 0 }; ///< Telescope diameter [m]
78
79 std::vector<realT> m_d_min{ 0 }; ///< Minimum AO system actuator pitch [m]. One per WFS mode.
80
81 realT m_d_opt{ 0 }; ///< Current optimum AO system actuator pitch [m]
82 bool m_optd{ false }; ///< Flag controlling whether actuator pitch is optimized (true) or just uses m_d_min (false).
83 ///< Default: true.
84 realT m_optd_delta{ 1.0 }; ///< The fractional change from d_min used in optimization. Set to 1 for integer
85 ///< binnings, > 1 for finer sampling.
86
87 wfs<realT, iosT> *m_wfsBeta{ nullptr }; ///< The WFS beta_p class.
88 bool m_ownWfsBeta{ false }; ///< Flag indicating if the WFS beta_p pointer is owned by this instance.
89
91 1 }; ///< The optical gain of the WFS, 0-1. Treated as a sensitivity reduction in measurement noise. Default 1.
92
93 protected:
94 realT m_lam_wfs{ 0 }; ///< WFS wavelength [m]
95
96 // WFS detector parameters, which may be a function of binning
97 std::vector<realT> m_npix_wfs{ 0 }; ///< Number of WFS pixels. One per WFS mode.
98 std::vector<realT> m_ron_wfs{ 0 }; ///< WFS readout noise [electrons/pix]. One per WFS mode.
99 std::vector<realT> m_Fbg{ 0 }; ///< Background flux, [photons/sec/pixel]. One per WFS mode.
100 std::vector<realT> m_minTauWFS{ 0 }; ///< Minimum WFS exposure time [sec]. One per WFS mode.
101
102 bool m_bin_npix{ true }; ///< Flag controlling whether or not to bin WFS pixels according to the actuator spacing.
103
105 0 }; ///< The optimum binning factor. If WFS modes are used, this is the mode index (0 to N-1). If not, it is
106 ///< 1 minus the pixel binning factor. It is always 1 minus the actuator binning factor.
107
108 realT m_tauWFS{ 0 }; ///< Actual WFS exposure time [sec]
109
110 realT m_deltaTau{ 0 }; ///< Loop latency [sec]
111
112 bool m_optTau{ true }; ///< Flag controlling whether optimum integration time is calculated (true) enforcing
113 ///< m_minTauWFS, or if m_tauWFS is used (false). Default: true.
114
115 realT m_lam_sci{ 0 }; ///< Science wavelength [m]
116
117 realT m_zeta{ 0 }; ///< Zenith angle [radians]
118 realT m_secZeta{ 1 }; ///< Secant of the Zenith angle (calculated)
119
120 int m_fit_mn_max{ 100 }; ///< Maximum spatial frequency index to use for fitting error calculation.
121
122 bool m_circularLimit{ false }; ///< Flag to indicate that the spatial frequency limit is circular, not square.
123
124 realT m_spatialFilter_ku{ std::numeric_limits<realT>::max() }; ///< The spatial filter cutoff in u, [m^-1]
125 realT m_spatialFilter_kv{ std::numeric_limits<realT>::max() }; ///< The spatial filter cutoff in v, [m^-1]
126
127 realT m_ncp_wfe{ 0 }; ///< Static WFE [m rms]
128 realT m_ncp_alpha{ 2 }; ///< Power-law exponent for the NCP aberations. Default is 2.
129
130 realT m_F0{ 0 }; ///< 0 mag flux from star at WFS [photons/sec]
131
132 realT m_starMag{ 0 }; ///< The magnitude of the star.
133
134 bool m_specsChanged{ true }; ///< Flag to indicate that a specification has changed.
135 bool m_dminChanged{ true }; ///< Flag to indicate that d_min has changed.
136
137 Eigen::Array<int, -1, -1> m_controlledModes; ///< Map of which modes are under control. Set by calcStrehl.
138
139 realT m_wfeMeasurement{ 0 }; ///< Total WFE due to measurement a error [rad^2 at m_lam_sci]
140 realT m_wfeTimeDelay{ 0 }; ///< Total WFE due to time delay [rad^2 at m_lam_sci]
141 realT m_wfeFitting{ 0 }; ///< Total WFE due to fitting error [rad^2 at m_lam_sci]
142 realT m_wfeChromScintOPD{ 0 }; ///< Total WFE due to the chromaticity of scintillation OPD [rad^2 at lam_sc]
143 realT m_wfeChromIndex{ 0 }; ///< Total WFE due to the chromaticity of the index of refraction [rad^2 at lam_Sci]
144 realT m_wfeAnisoOPD{ 0 }; ///< Total WFE due to dispersive anisoplanatism OPD.
145
146 realT m_wfeNCP{ 0 }; ///< Total WFE due to NCP errors [rad^2 at m_lam_sci]
147
148 realT m_wfeVar{ 0 }; ///< The WFE variance, in meters^2. Never use this directly, instead use wfeVar().
149 realT m_strehl{ 0 }; ///< Strehl ratio, a calculated quantity. Never use this directdy, instead use strehl().
150
151 public:
152 /// Default c'tor
153 /** Calls initialize().
154 */
155 aoSystem();
156
157 /// Destructor
158 ~aoSystem();
159
160 /// Load the default parameters from Guyon, 2005 \cite guyon_2005.
161 /**
162 *
163 */
164 void loadGuyon2005();
165
166 /// Load parameters corresponding to the MagAO-X system.
167 void loadMagAOX();
168
169 /// Load parameters corresponding to the G-MagAO-X system.
170 void loadGMagAOX();
171
172 /// Set the value of the primary mirror diameter.
173 /**
174 */
175 void D( realT nD /**< [in] is the new value of m_D. */ );
176
177 /// Get the value of the primary mirror diamter
178 /**
179 * \returns the current value of m_D.
180 */
181 realT D();
182
183 /// Set the minimum subaperture sampling for each WFS mode.
184 /** This sets the vector
185 */
186 void d_min( const std::vector<realT> &nd /**< [in] is the new values of m_d_min */ );
187
188 /// Set the minimum subaperture sampling for a given WFS mode
189 /** This sets the entry in the vector, if it has the proper length.
190 */
191 void d_min( int idx, ///< [in] the index of the WFS mode
192 realT nd ///< [in] is the new value of m_d_min[idx]
193 );
194
195 /// Get the minimum subaperture sampling for a given WFS mode
196 /**
197 * \returns the current value of m_d_min[idx]
198 */
199 realT d_min( size_t idx /**< [in] the index of the WFS mode.*/ );
200
201 /// Get the minimum subaperture sampling for all WFS modes.
202 /**
203 * \returns the current value of m_d_min.
204 */
205 std::vector<realT> d_min();
206
207 /// Set whether or not the value of d is optimized or just set to m_d_min.
208 /**
209 */
210 void optd( bool od /**< [in] is the new value of m_optd */ );
211
212 /// Get the value of m_optd.
213 /**
214 * \returns the new value of m_optd_delta.
215 */
216 bool optd();
217
218 /// Set the fractional change in actuator spacing for optimization.
219 /** Sets the fraction of m_d_min by which the optimizer changes actautor spacing.
220 */
221 void optd_delta( realT odd /**< [in] is the new value of m_optd_delta */ );
222
223 /// Get the value of the fractional change in actuator spacing for optimization..
224 /**
225 * \returns the value of m_optd_delta.
226 */
227 realT optd_delta();
228
229 /// Set the WFS Beta pointer
230 /** The WFS beta parameter determines the photon noise sensitivity.
231 *
232 * \tparam wfsT is a type derived from ao::analysis::wfs
233 */
234 template <typename wfsT>
235 void wfsBeta( const wfsT &w /**< [in] an object derived from ao::analysis::wfs base class*/ );
236
237 /// Set the WFS Beta pointer
238 /** The WFS beta parameter determines the photon noise sensitivity.
239 *
240 * \tparam wfsT is a type derived from ao::analysis::wfs
241 */
242 template<typename wfsT>
243 void wfsBeta( const wfsT * w /**< [in] pointer to an object derived from ao::analysis::wfs base class. If nullptr then a new object is allocated and managed.*/);
244
245 /// Get the WFS Beta pointer
246 /**
247 * \returns a pointer to the current m_wfsBeta
248 */
250
251 /// Get the value of beta_p for a spatial frequency
252 /** beta_p is the photon noise sensitivity of the WFS.
253 *
254 * \returns beta_p as calculated by the WFS.
255 */
256 realT beta_p( realT m, ///< [in] the spatial frequency index
257 realT n ///< [in] the spatial frequency index
258 );
259
260 /// Get the value of beta_r for a spatial frequency
261 /** beta_r is the read noise sensitivity of the WFS.
262 *
263 * \returns beta_r as calculated by the WFS.
264 */
265 realT beta_r( realT m, ///< [in] the spatial frequency index
266 realT n ///< [in] the spatial frequency index
267 );
268
269 /// Set the value of the Optical Gain.
270 /**
271 */
272 void opticalGain( realT og /**< [in] the new value of m_opticalGain*/ );
273
274 /// Get the value of the optical gain
275 /** Optical gain is applied in the WFS noise calculation
276 *
277 * \returns the current value of m_opticalGain
278 */
279 realT opticalGain();
280
281 /// Set the value of the WFS wavelength.
282 /**
283 */
284 void lam_wfs( realT nlam /**< [in] is the new value of m_lam_wfs */ );
285
286 /// Get the value of the WFS wavelength.
287 /**
288 * \returns the current value of m_lam_wfs.
289 */
290 realT lam_wfs();
291
292 /// Set the number of pixels in the WFS for each WFS mode
293 /** This sets the vector.
294 */
295 void npix_wfs( const std::vector<realT> &npix /**< [in] is the new values of m_npix_wfs */ );
296
297 /// Set the number of pixels in the WFS for a given WFS mode
298 /** This sets the entry in the vector, if it has the proper length.
299 */
300 void npix_wfs( int idx, ///< [in] the index of the WFS mode
301 realT npix ///< [in] is the new value of m_npix_wfs[idx]
302 );
303
304 /// Get the number of pixels in the WFS for a given WFS mode
305 /**
306 * \returns the current value of m_npix_wfs[idx]
307 */
308 realT npix_wfs( size_t idx /**< [in] the index of the WFS mode.*/ );
309
310 /// Get the number of pixels in the WFS for all WFS modes
311 /**
312 * \returns the current value of m_npix_wfs
313 */
314 std::vector<realT> npix_wfs();
315
316 /// Set the value of the WFS readout noise for each WFS mode
317 /** This sets the vector.
318 */
319 void ron_wfs( const std::vector<realT> &nron /**< [in] is the new value of m_ron_wfs */ );
320
321 /// Set the value of the WFS readout noise for a given bining mode
322 /** This sets the entry in the vector if the vector has the proper length
323 */
324 void ron_wfs( int idx, ///< [in] the index of the WFS mode
325 realT nron ///< [in] is the new value of m_ron_wfs [idx]
326 );
327
328 /// Get the value of the WFS readout noise for a given WFS mode
329 /**
330 * \returns the current value of m_ron_wfs[idx]
331 */
332 realT ron_wfs( size_t idx /**< [in] the index of the WFS mode.*/ );
333
334 /// Get the value of the WFS readout noise for all WFS modes
335 /**
336 * \returns the current value of m_ron_wfs
337 */
338 std::vector<realT> ron_wfs();
339
340 /// Set the value of the background fluxes.
341 /** This sets the value of the background flux for each WFS mode
342 */
343 void Fbg( const std::vector<realT> &fbg /**< [in] is the new value of m_Fbg */ );
344
345 /// Set a single value of the background flux for a given WFS mode
346 /** This sets the entry in the vector if the vector has the proper length
347 */
348 void Fbg( int idx, ///< [in] the index of the WFS mode
349 realT fbg ///< [in] is the new value of m_Fbg[idx]
350 );
351
352 /// Get the value of the background flux for a given WFS mode
353 /**
354 * \returns m_Fbg[idx]
355 */
356 realT Fbg( size_t idx /**< [in] the index of the WFS mode.*/ );
357
358 /// Get the value of the background flux for all WFS modes
359 /**
360 * \returns the current value of m_Fbg
361 */
362 std::vector<realT> Fbg();
363
364 /// Set the value of the minimum WFS exposure times.
365 /** This sets the vector of the minimum WFS exposure times
366 */
367 void minTauWFS( const std::vector<realT> &ntau /**< [in] is the new value of m_minTauWFS */ );
368
369 /// Set a single value of the minimum WFS exposure time for a given WFS mode.
370 /** This sets the entry in the vector if the vector has sufficient length.
371 */
372 void minTauWFS( size_t idx, ///< [in] the index of the WFS mode
373 realT ntau ///< [in] is the new value of m_minTauWFS
374 );
375
376 /// Get the value of the minimum WFS exposure time for a given binning.
377 /**
378 * \returns the current value of m_minTauWFS[idx].
379 */
380 realT minTauWFS( size_t idx /**< [in] the index of the WFS mode.*/ );
381
382 /// Get the values of the minimum WFS exposure time.
383 /**
384 * \returns the current value of m_minTauWFS.
385 */
386 std::vector<realT> minTauWFS();
387
388 /// Set the value of the pixel binning flag
389 /**
390 */
391 void bin_npix( bool bnp /**< [in] is the new value of m_bin_npix */ );
392
393 /// Get the value of the pixel binngin flag
394 /**
395 * \returns
396 */
397 bool bin_npix();
398
399 /// Get the value of the optimum binning factor/index.
400 /** If WFS modes are used, this is the mode index (0 to N-1). If not, it is 1 minus the pixel binning factor.
401 * It is always 1 minus the actuator binning factor.
402 * This calls opt_d() to perform the optimization if needed.
403 *
404 * \returns the current value of m_bin_opt.
405 */
406 int bin_opt();
407
408 /// Set the value of the WFS exposure time.
409 /**
410 */
411 void tauWFS( realT ntau /**< [in] is the new value of m_tauWFS */ );
412
413 /// Get the value of the minimum WFS exposure time.
414 /**
415 * \returns the current value of m_tauWFS.
416 */
417 realT tauWFS();
418
419 /// Set the value of m_deltaTau.
420 /**
421 */
422 void deltaTau( realT ndel /**< [in] is the new value of m_deltaTau*/ );
423
424 /// Get the value of m_deltaTau.
425 /**
426 * \returns the current value of m_deltaTau.
427 */
428 realT deltaTau();
429
430 /// Set the value of m_optTau.
431 /**
432 */
433 void optTau( bool ot /**< [in] is the new value of m_optTau */ );
434
435 /// Get the value of m_optTau.
436 /**
437 * \returns the current value of m_optTau.
438 */
439 bool optTau();
440
441 /// Set the science wavelength.
442 /**
443 */
444 void lam_sci( realT nlam /**< [in] is the new value of m_lam_sci */ );
445
446 /// Get the science wavelength.
447 /**
448 * \returns the current value of m_lam_sci.
449 */
450 realT lam_sci();
451
452 /// Set the zenith angle, and its secant.
453 /**
454 */
455 void zeta( realT nz /**< [in] The new value of m_zeta */ );
456
457 /// Get the zenith angle
458 /**
459 * \return the current value of m_zeta
460 */
461 realT zeta();
462
463 /// Get the zecant of the zenith angle
464 /**
465 * \return the current value of m_secZeta
466 */
467 realT secZeta();
468
469 /// Set the value of m_fit_mn_max
470 /**
471 */
472 void fit_mn_max( int mnm /**< [in] is the new value of m_fit_mn_max */ );
473
474 /// Get the value of m_fit_mn_max
475 /**
476 */
477 int fit_mn_max();
478
479 /// Set the value of the circularLimit flag
480 /** If this is true, then the spatial frequency limits are treated as a circle
481 * rather than a square. This will create a circular dark hole.
482 */
483 void circularLimit( bool cl /**< [in] the new value of the circularLimit flag*/ );
484
485 /// Get the value of the circularLimit flag
486 /**
487 * \returns the current value of m_circularLimit
488 */
489 bool circularLimit();
490
491 /// Set the value of spatialFilter_ku
492 /**
493 */
494 void spatialFilter_ku( realT ku /**< [in] is the new value of spatialFilter_ku*/ );
495
496 /// Get the value of spatialFilter_ku
497 /** \returns the current value of m_spatialFilter_ku
498 */
499 realT spatialFilter_ku();
500
501 /// Set the value of spatialFilter_kv
502 /**
503 */
504 void spatialFilter_kv( realT kv /**< [in] is the new value of spatialFilter_kv*/ );
505
506 /// Get the value of spatialFilter_kv
507 /** \returns the current value of m_spatialFilter_kv
508 */
509 realT spatialFilter_kv();
510
511 /// Set the value of the non-common path WFE.
512 /**
513 */
514 void ncp_wfe( realT nwfe /**< [in] is the new value of m_ncp_wfe*/ );
515
516 /// Get the value of the non-common path WFE.
517 /**
518 * \returns the current value of m_ncp_wfe.
519 */
520 realT ncp_wfe();
521
522 /// Set the value of the non-common path WFE PSD index.
523 /**
524 */
525 void ncp_alpha( realT alpha /**< [in] is the new value of m_ncp_alpha*/ );
526
527 /// Get the value of the non-common path WFE PSD index.
528 /**
529 * \returns the current value of m_ncp_alpha.
530 */
531 realT ncp_alpha();
532
533 /// Set the value of the 0 magnitude photon rate
534 /** This is the photon rate at the WFS, photons/sec.
535 *
536 */
537 void F0( realT nF0 /**< [in] is the new value of m_F0.*/ );
538
539 /// Get the value of the 0 magnitude photon rate
540 /**
541 * \returns the current value of m_F0.
542 */
543 realT F0();
544
545 /// Set the value of the Star's magnitude
546 /**
547 */
548 void starMag( realT nmag /**< [in] is the new value of m_starMag.*/ );
549
550 /// Get the value of the Star's magnitude
551 /**
552 * \returns the current value of m_starMag
553 */
554 realT starMag();
555
556 /// The photon flux at a given star magnitude.
557 /**
558 * \returns the photon flux of the star in the bandpass assumed by F0.
559 */
560 realT Fg( realT mag /**< [in] is the magnitude of the star. */ );
561
562 /// Get the photon rate at the current Star magnitude.
563 /** Calculates \f$ F_\gamma = F_0 10^{-0.4 m} \f$ where \f$ F_0 \f$ is m_F0 and \f$ m \f$ is m_starMag.
564 *
565 * \returns the current value of the current photon rate.
566 */
567 realT Fg();
568
569 /// Access the array of controlledModes
570 /** This array indicates which modes are controlled after optimizeStrehl chooses them.
571 *
572 * \returns a const reference to m_controlledModes.
573 */
574 const Eigen::Array<int, -1, -1> &controlledModes();
575
576 /** \name Measurement Error
577 * Calculating the WFE due to WFS measurement noise.
578 * @{
579 */
580
581 /// Calculate the terms of the signal to noise ratio squared (S/N)^2 for the WFS measurement
582 /** The S/N squared is
583 \f[
584 (S/N)^2 = \frac{ F_\gamma^2 \tau_{wfs}^2 }{ F_\gamma \tau_{wfs} + n_{pix} F_{bg} \tau_{wfs} + n_{pix}
585 \sigma_{ron}^2 } \f]
586
587 *
588 * \returns the S/N squared
589 */
590 realT
591 signal2Noise2( realT &Nph, ///< [out] the number of photons
592 realT tau_wfs, ///< [in] specifies the WFS exposure time. If 0, then optimumTauWFS is used
593 realT d, ///< [in] the actuator spacing in meters, used if binning WFS pixels
594 int b ///< [in] the binning parameter. Either the WFS mode index, or the binning factor minus 1.
595 );
596
597 /// Calculate the measurement noise at a spatial frequency and specified actuator spacing
598 /** Calculates the wavefront phase variance due measurement noise at \f$ k = (m/D)\hat{u} + (n/D)\hat{v} \f$.
599 *
600 * \returns the measurement noise in rad^2 rms at the science wavelength
601 */
602 realT
603 measurementError( realT m, ///< [in] specifies the u component of the spatial frequency
604 realT n, ///< [in] specifies the v component of the spatial frequency
605 realT d, ///< [in] the actuator spacing in meters
606 int b ///< [in] the binning parameter. Either the WFS mode index, or the binning factor minus 1.
607 );
608
609 /// Calculate the total measurement error over all corrected spatial frequencies
610 /** This totals the measurement WFE at the optimum actuator spacing and binning.
611 *
612 * \overload
613 *
614 * \returns the total WFE due to measurement error.
615 */
616 realT measurementErrorTotal();
617
618 ///@}
619
620 /** \name Time Delay Error
621 * Calculating the WFE due to time delay.
622 * @{
623 */
624
625 /// Calculate the time delay at a spatial frequency at the optimum exposure time and the specified actuator spacing
626 /** Calculates the wavefront phase variance due to time delay at \f$ f = (m/D)\hat{u} + (n/D)\hat{v} \f$.
627 *
628 * \returns the measurement noise in rad^2 rms at the science wavelength
629 */
630 realT
631 timeDelayError( realT m, ///< [in] specifies the u component of the spatial frequency
632 realT n, ///< [in] specifies the v component of the spatial frequency
633 realT d, ///< [in] the actuator spacing, in meters
634 int b ///< [in] the binning parameter. Either the WFS mode index, or the binning factor minus 1.
635 );
636
637 /// Calculate the time delay error over all corrected spatial frequencies
638 /** This totals the time delay WFE at the optimum actuator spacing and binning.
639 *
640 * \overload
641 *
642 * \returns the total WFE due to time delay.
643 */
644 realT timeDelayErrorTotal();
645
646 ///@}
647
648 /** \name Fitting Error
649 * Calculating the WFE due to uncorrected spatial frequencies.
650 * @{
651 */
652
653 /// Calculate the fitting error at a specific spatial frequency.
654 /**
655 * \returns the fitting error in rad^2 rms at the science wavelength at (m,n).
656 */
657 realT fittingError( realT m, ///< [in] specifies the u component of the spatial frequency
658 realT n ///< [in] specifies the v component of the spatial frequency
659 );
660
661 /// Calculate the total fitting error over all uncorrected spatial frequencies.
662 /** This totals the fitting WFE for the optimum actuator spacing and binning.
663 *
664 * \overload
665 *
666 * \returns the total fitting error.
667 */
668 realT fittingErrorTotal();
669
670 ///@}
671
672 /** \name Chromatic Errors
673 * Calculating the WFE due to chromaticity in scintillaion, index of refraction, and dispersion.
674 * @{
675 */
676
677 /// Calculate the wavefront error due to scintillation chromaticity in the OPD at a spatial frequency.
678 /** This totals the WFE from scintillation chromaticity in the OPD.
679 *
680 * \returns the WFE in rad^2 rms at the science wavelength at (m,n).
681 */
682 realT chromScintOPDError( int m, ///< [in] specifies the u component of the spatial frequency
683 int n ///< [in] specifies the v component of the spatial frequency
684 );
685
686 /// Calculate the wavefront error due to scintillation chromaticity in the OPD over all spatial frequencies.
687 /** This totals the WFE from scintillation chromaticity in the OPD for the optimum actuator spacing.
688 *
689 * \overload
690 *
691 * \returns the WFE in rad^2 rms at the science wavelength at (m,n).
692 */
694
695 /// Calculate the wavefront error due to scintillation chromaticity in amplitude for a given spatial frequency.
696 /**
697 * \returns the WFE in rad^2 rms at the science wavelength at (m,n).
698 */
699 realT chromScintAmpError( int m, ///< [in] specifies the u component of the spatial frequency
700 int n ///< [in] specifies the v component of the spatial frequency
701 );
702
703 /// Calculate the wavefront error due to scintillation chromaticity in amplitude over all spatial frequencies.
704 /**
705 * \returns the WFE in rad^2 rms at the science wavelength at (m,n).
706 */
708
709 /// Calculate the wavefront error due to chromaticity in the index of refraction at a given spatial frequency.
710 /** This totals the WFE from chromaticity in the index of refraction.
711 *
712 * \overload
713 *
714 * \returns the WFE in rad^2 rms at the science wavelength at (m,n).
715 */
716 realT chromIndexError( int m, ///< [in] specifies the u component of the spatial frequency
717 int n ///< [in] specifies the v component of the spatial frequency
718 );
719
720 /// Calculate the wavefront error due to chromaticity in the index of refraction at a specific spatial frequency.
721 /** This totals the WFE from chromaticity in the index of refraction for the optimum actuator spacing.
722 *
723 * \overload
724 *
725 * \returns the WFE in rad^2 rms at the science wavelength at (m,n).
726 */
727 realT chromIndexErrorTotal();
728
729 /// Calculate the wavefront error due to dispersive anisoplanatism in the OPD at a given spatial frequency
730 /** This calculates the WFE from dispersive anisoplanatism in the OPD.
731 *
732 * \returns the WFE in rad^2 rms at the science wavelength at (m,n).
733 */
734 realT dispAnisoOPDError( int m, ///< [in] specifies the u component of the spatial frequency
735 int n ///< [in] specifies the v component of the spatial frequency
736 );
737
738 /// Calculate the wavefront error due to dispersive anisoplanatism in the OPD over all specific spatial frequencies.
739 /** This totals the WFE from dispersive anisoplanatism in the OPD for the optimum actuator spacing.
740 *
741 * \overload
742 *
743 * \returns the WFE in rad^2 rms at the science wavelength at (m,n).
744 */
746
747 /// Calculate the wavefront error due to dispersive anisoplanatism in the amplitude at a specific spatial frequency.
748 /**
749 * \returns the WFE in rad^2 rms at the science wavelength at (m,n).
750 */
751 realT dispAnisoAmpError();
752
753 ///@}
754
755 /** \name Optimum Parameters
756 * Functions to calculate the optimum integration time and actuator spacing.
757 *
758 * @{
759 */
760
761 /// Calculate the optimum exposure time for a given spatial frequency at a specified actuator spacing
762 /** Finds the optimum exposure time at \f$ k = (m/D)\hat{u} + (n/D)\hat{v} \f$.
763 *
764 * \todo Check inclusion of X in parameters
765 *
766 * \returns the optimum expsoure time.
767 */
768 realT
769 optimumTauWFS( realT m, ///< [in] is the spatial frequency index in u
770 realT n, ///< [in] is the spatial frequency index in v
771 realT d, ///< [in] the actuator spacing, in meters
772 int b ///< [in] the binning parameter. Either the WFS mode index, or the binning factor minus 1.
773 );
774
775 /// Calculate the optimum exposure time for a given spatial frequency at the optimum actuator spacing.
776 /** Finds the optimum exposure time at \f$ k = (m/D)\hat{u} + (n/D)\hat{v} \f$.
777 * If optd == false this uses the minimum actuator spacing.
778 *
779 * \todo Check inclusion of X in parameters
780 *
781 * \returns the optimum expsoure time.
782 */
783 realT optimumTauWFS( realT m, ///< [in] is the spatial frequency index in u
784 realT n ///< [in] is the spatial frequency index in v
785 );
786
787 /// Calculate the optimum actuator spacing.
788 /** Finds the value of m_d_opt where the fitting error is less than than the combined time delay and measurement
789 * error.
790 *
791 * \returns the current value of m_d_opt.
792 */
793 realT d_opt();
794
795 /// @}
796
797 /// Calculate the NCP variance at a spatial frequency.
798 /** Finds the NCP variance at \f$ k = (m/D)\hat{u} + (n/D)\hat{v} \f$.
799 *
800 * \returns the NCP variance at k.
801 */
802 realT ncpError( int m, ///< [in] is the spatial frequency index in u
803 int n ///< [in] is the spatial frequency index in v
804 );
805
806 /// Calculate the total NCP variance in rad^2.
807 /**
808 * \returns the total NCP variance.
809 */
810 realT ncpError();
811
812 /** \name Overall WFE and Strehl Ratio
813 * Functions to calculate the total WFE and Strehl ratio.
814 *
815 * @{
816 */
817
818 protected:
819 /// Calculate the component WFE, total WFE, and Strehl ratio.
820 /** This should only be called when something changes.
821 */
822 void calcStrehl();
823
824 public:
825 /// Get the current value of the total WFE variance.
826 /** If no changes, merely returns m_wfeVar. Calls calcStrehl if there are changes.
827 *
828 * \returns the current value of m_wfeVar;
829 */
830 realT wfeVar();
831
832 /// Get the Strehl ratio for a given actuator pitch and WFS WFS mode.
833 /**
834 * Strehl is calculated using the extended Marechal approximation:
835 *
836 \f[
837 S = e^{-\sigma_{wfe}^2}
838 \f]
839 * where \f$ \sigma_{wfe}^2 \f$ is the WFE for the spacing and binning
840 *
841 * \returns the strehl for these values
842 */
843 realT strehl( realT d, ///< [in] the actuator spacing, in meters
844 int b ///< [in] the binning parameter. Either the WFS mode index, or the binning factor minus 1.
845 );
846
847 /// Get the current Strehl ratio.
848 /** If no changes, merely returns m_strehl. Calls calcStrehl if there are changes.
849 * Strehl is calculated using the extended Marechal approximation:
850 *
851 \f[
852 S = e^{-\sigma_{wfe}^2}
853 \f]
854 * where \f$ \sigma_{wfe}^2 \f$ is the current value of m_wfeVar.
855 *
856 * \returns the current value of m_strehl.
857 */
858 realT strehl();
859
860 /// @}
861
862 /** \name Contrast
863 * Calculating contrast
864 *
865 * @{
866 */
867
868 /// Worker function for raw contrast fuctions
869 /** The logic in this calculation is the same regardless of component, except for the calculation of variance.
870 * The varFunc function pointer is used to calculate the variance, otherwise This handles the other details such
871 * as Strehl normalization, fitting error (if appropriate), and bounds checking.
872 *
873 * \note this gives the raw PSD, it must be convolved with the PSF for the true contrast.
874 *
875 * \tparam varFuncT is a function-pointer-to-member (of this class) with signature realT(*f)(realT, realT)
876 */
877 template <typename varFuncT>
878 realT C_( int m, ///< [in] is the spatial frequency index in u/
879 int n, ///< [in] is the spatial frequency index in v.
880 bool normStrehl, ///< [in] flag controls whether the contrast is normalized by Strehl ratio/
881 varFuncT varFunc, ///< [in] the variance function to use.
882 int doFittingError ///< [in] flag to describe how fitting error should be considered for this term:
883 ///< FITTING_ERROR_NO, FITTING_ERROR_ZERO, FITTING_ERROR_X, or FITTING_ERROR_Y.
884 );
885
886 /// Worker function for the contrast-map functions.
887 /** The map calculation is the same for all terms, except for the calculation of variance at each pixel.
888 * The Cfunc function pointer is used to actually get the contrast.
889 *
890 * \note this is the raw PSD, it must be convolved with the PSF for the true contrast.
891 *
892 * \tparam imageT is an Eigen-like image
893 * \tparam CfuncT is a function-pointer-to-member (of this class) with signature realT (*f)(realT, realT, bool)
894 */
895 template <typename imageT, typename CfuncT>
896 void C_Map( imageT &im, ///< [out] the map image to be filled in.
897 CfuncT Cfunc, ///< [in] the raw contrast function to use for filling in the map.
898 bool normStrehl ///< [in] flag controlling whether or not the map is normalized by Strehl
899 );
900
901 /// Calculate the residual variance due to uncorrected phase at a spatial frequency.
902 /** Used to calculate contrast \ref C0().
903 *
904 * \returns variance at (m,n).
905 */
906 realT C0var( realT m, ///< [in] is the spatial frequency index in u
907 realT n ///< [in] is the spatial frequency index in v
908 );
909
910 /// Calculate the contrast due to uncorrected phase, C0.
911 /** Contrast C0 is the uncorrected phase, with the effects of scintillation included. See Guyon (2005) \cite
912 * guyon_2005, and the updated derivation in Males \& Guyon (2018) \cite males_guyon_2018.
913 *
914 * \note this is the raw PSD, it must be convolved with the PSF for the true contrast.
915 *
916 * \returns C0.
917 */
918 realT C0( realT m, ///< [in] is the spatial frequency index in u
919 realT n, ///< [in] is the spatial frequency index in v
920 bool normStrehl = true ///< [in] flag controls whether the contrast is normalized by Strehl ratio
921 );
922
923 /// Calculate a 2D map of contrast C0
924 /**
925 * \tparam imageT is an Eigen-like image type.
926 */
927 template <typename imageT>
928 void C0Map( imageT &map, ///< [in] the map image to be filled in with contrast
929 bool normStrehl );
930
931 /// Calculate the residual variance due to uncorrected amplitude at a spatial frequency.
932 /** Used to calculate contrast \ref C1().
933 *
934 * \returns variance at (m,n).
935 */
936 realT C1var( realT m, ///< [in] is the spatial frequency index in u
937 realT n ///< [in] is the spatial frequency index in v
938 );
939
940 /// Calculate the contrast due to uncorrected amplitude, C1.
941 /** Contrast C1 is the uncorrected amplitude due to scintillation. See Guyon (2005) \cite guyon_2005, and the
942 * updated derivation in Males \& Guyon (2018) \cite males_guyon_2018.
943 *
944 * \returns C0.
945 */
946 realT C1( realT m, ///< [in] is the spatial frequency index in u
947 realT n, ///< [in] is the spatial frequency index in v
948 bool normStrehl = true ///< [in] flag controls whether the contrast is normalized by Strehl ratio
949 );
950
951 /// Calculate a 2D map of contrast C1.
952 /** The contrast is Strehl-normalized here.
953 *
954 * \note this is the raw PSD, it must be convolved with the PSF for the true contrast.
955 *
956 * \tparam imageT is an Eigen-like image type.
957 */
958 template <typename imageT>
959 void C1Map( imageT &map, ///< [in] the map image to be filled in with contrast
960 bool normStrehl );
961
962 /// Calculate the residual variance due to measurement and time delay errors in phase/OPD at a spatial frequency.
963 /** Used to calculate contrast \ref C2().
964 *
965 * \returns variance at (m,n).
966 */
967 realT C2var( realT m, ///< [in] is the spatial frequency index in u
968 realT n ///< [in] is the spatial frequency index in v
969 );
970
971 /// Calculate the contrast due to measurement and time delay errors in phase/OPD at a spatial frequency.
972 /** Contrast C2 is just the total variance due to time delay and measurement errors,
973 * divided by the Strehl ratio. See Guyon (2005) \cite guyon_2005, and the updated
974 * derivation in Males \& Guyon (2018) \cite males_guyon_2018.
975 *
976 * \returns C2.
977 */
978 realT C2( realT m, ///< [in] is the spatial frequency index in u
979 realT n, ///< [in] is the spatial frequency index in v
980 bool normStrehl = true ///< [in] flag controls whether the contrast is normalized by Strehl ratio.
981 );
982
983 /// Calculate a 2D map of contrast \ref C2().
984 /** The contrast is Strehl-normalized here.
985 *
986 * \note this is the raw PSD, it must be convolved with the PSF for the true contrast.
987 *
988 * \tparam imageT is an Eigen-like image type.
989 */
990 template <typename imageT>
991 void C2Map( imageT &map, ///< [in] the map image to be filled in with contrast
992 bool normStrehl );
993
994 /// Calculate the residual variance due to measurement and time delay errors in amplitude at a spatial frequency.
995 /** Used to calculate contrast \ref C3().
996 *
997 * \returns variance at (m,n).
998 */
999 realT C3var( realT m, ///< [in] is the spatial frequency index in u
1000 realT n ///< [in] is the spatial frequency index in v
1001 );
1002
1003 /// Calculate the contrast due to measurement and time delay errors in amplitude at a spatial frequency.
1004 /** Contrast C3 is just the total variance due to time delay and measurement errors,
1005 * divided by the Strehl ratio. See Guyon (2005) \cite guyon_2005, and the updated
1006 * derivation in Males \& Guyon (2018) \cite males_guyon_2018.
1007 *
1008 * \returns C3.
1009 */
1010 realT C3( realT m, ///< [in] is the spatial frequency index in u
1011 realT n, ///< [in] is the spatial frequency index in v
1012 bool normStrehl = true ///< [in] flag controls whether the contrast is normalized by Strehl ratio.
1013 );
1014
1015 /// Calculate a 2D map of contrast C3.
1016 /** The contrast is Strehl-normalized here.
1017 *
1018 * \note this is the raw PSD, it must be convolved with the PSF for the true contrast.
1019 *
1020 * \tparam imageT is an Eigen-like image type.
1021 */
1022 template <typename imageT>
1023 void C3Map( imageT &map, ///< [in] the map image to be filled in with contrast
1024 bool normStrehl );
1025
1026 /// Calculate the residual variance due to scintilation-OPD chromaticity.
1027 /** Used to calculate contrast \ref C4().
1028 *
1029 * \returns variance at (m,n).
1030 */
1031 realT C4var( realT m, ///< [in] is the spatial frequency index in u
1032 realT n ///< [in] is the spatial frequency index in v
1033 );
1034
1035 /// Calculate the contrast due to scintilation-OPD chromaticity.
1036 /** Contrast C4 is due to the chromaticity of scintillation, causing the ODP measurement to be slightly incorrect at
1037 * the science wavelength. See Guyon (2005) \cite guyon_2005, and the updated derivation in Males \& Guyon (2018)
1038 * \cite males_guyon_2018.
1039 *
1040 * \returns C4.
1041 */
1042 realT C4( realT m, ///< [in] is the spatial frequency index in u
1043 realT n, ///< [in] is the spatial frequency index in v
1044 bool normStrehl = true ///< [in] flag controls whether the contrast is normalized by Strehl ratio.
1045 );
1046
1047 /// Calculate a 2D map of contrast C4
1048 /** The contrast is Strehl-normalized here.
1049 *
1050 * \note this is the raw PSD, it must be convolved with the PSF for the true contrast.
1051 *
1052 * \tparam imageT is an Eigen-like image type.
1053 */
1054 template <typename imageT>
1055 void C4Map( imageT &map, ///< [in] the map image to be filled in with contrast
1056 bool normStrehl );
1057
1058 /// Calculate the residual variance due to to scintilation-amplitude chromaticity.
1059 /** Used to calculate contrast \ref C5().
1060 *
1061 * \note this is the raw PSD, it must be convolved with the PSF for the true contrast.
1062 *
1063 * \returns variance at (m,n).
1064 */
1065 realT C5var( realT m, ///< [in] is the spatial frequency index in u
1066 realT n ///< [in] is the spatial frequency index in v
1067 );
1068
1069 /// Calculate the contrast due to scintilation-amplitude chromaticity.
1070 /** Contrast C5 is due to the chromaticity of scintillation, causing the amplitude measurement to be slightly
1071 * incorrect at the science wavelength. See Guyon (2005) \cite guyon_2005, and the updated derivation in Males \&
1072 * Guyon (2018) \cite males_guyon_2018.
1073 *
1074 * \returns C4.
1075 */
1076 realT C5( realT m, ///< [in] is the spatial frequency index in u
1077 realT n, ///< [in] is the spatial frequency index in v
1078 bool normStrehl = true ///< [in] flag controls whether the contrast is normalized by Strehl ratio.
1079 );
1080
1081 /// Calculate a 2D map of contrast C5
1082 /** The contrast is Strehl-normalized here.
1083 *
1084 * \note this is the raw PSD, it must be convolved with the PSF for the true contrast.
1085 *
1086 * \tparam imageT is an Eigen-like image type.
1087 */
1088 template <typename imageT>
1089 void C5Map( imageT &map, ///< [in] the map image to be filled in with contrast
1090 bool normStrehl );
1091
1092 /// Calculate the residual variance due to chromaticity of the index of refraction of air.
1093 /** Used to calculate contrast \ref C6().
1094 *
1095 * \returns variance at (m,n).
1096 */
1097 realT C6var( realT m, ///< [in] is the spatial frequency index in u
1098 realT n ///< [in] is the spatial frequency index in v
1099 );
1100
1101 /// Calculate the contrast due to chromaticity of the index of refraction of air.
1102 /** Contrast C6 is due to the index of refraction of air being wavelength dependent, causing the ODP measurement to
1103 * be slightly incorrect at the science wavelength. See Guyon (2005) \cite guyon_2005, and the updated derivation
1104 * in Males \& Guyon (2018) \cite males_guyon_2018.
1105 *
1106 * \returns C6.
1107 */
1108 realT C6( realT m, ///< [in] is the spatial frequency index in u
1109 realT n, ///< [in] is the spatial frequency index in v
1110 bool normStrehl = true ///< flag controls whether the contrast is normalized by Strehl ratio.
1111 );
1112
1113 /// Calculate a 2D map of contrast C6
1114 /** The contrast is Strehl-normalized here.
1115 *
1116 * \note this is the raw PSD, it must be convolved with the PSF for the true contrast.
1117 *
1118 * \tparam imageT is an Eigen-like image type.
1119 */
1120 template <typename imageT>
1121 void C6Map( imageT &map, ///< [in] the map image to be filled in with contrast
1122 bool normStrehl );
1123
1124 /// Calculate the residual variance due to dispersive anisoplanatism.
1125 /** Used to calculate contrast \ref C7().
1126 *
1127 * \returns variance at (m,n).
1128 */
1129 realT C7var( realT m, ///< [in] is the spatial frequency index in u
1130 realT n ///< [in] is the spatial frequency index in v
1131 );
1132
1133 /// Calculate the contrast due to dispersive anisoplanatism.
1134 /** Contrast C7 is due to atmospheric dispersion causing light at different wavelengths to take different paths
1135 * through atmospheric turbulence. See Fitzgerald (2017, in prep) \cite fitzgerald_2017
1136 *
1137 * \returns C7.
1138 */
1139 realT C7( realT m, ///< [in] is the spatial frequency index in u
1140 realT n, ///< [in] is the spatial frequency index in v
1141 bool normStrehl = true ///< flag controls whether the contrast is normalized by Strehl ratio.
1142 );
1143
1144 /// Calculate a 2D map of contrast C7
1145 /** The contrast is Strehl-normalized here.
1146 *
1147 * \note this is the raw PSD, it must be convolved with the PSF for the true contrast.
1148 *
1149 * \tparam imageT is an Eigen-like image type.
1150 */
1151 template <typename imageT>
1152 void C7Map( imageT &map, ///< [in] the map image to be filled in with contrast
1153 bool normStrehl );
1154
1155 ///@}
1156
1157 /// Output current parameters to a stream
1158 /** Outputs a formatted list of all current parameters.
1159 *
1160 */
1161 iosT &dumpAOSystem( iosT &ios /**< [in] a std::ostream-like stream. */ );
1162
1163 /// Setup the configurator to configure this class
1164 /**
1165 * todo: "\test Loading aoAtmosphere config settings \ref tests_ao_analysis_aoAtmosphere_config "[test doc]"
1166 */
1167 void setupConfig( app::appConfigurator &config /**< [in] the app::configurator object*/ );
1168
1169 /// Load the configuration of this class from a configurator
1170 /**
1171 * \todo: "\test Loading aoAtmosphere config settings \ref tests_ao_analysis_aoAtmosphere_config "[test doc]""
1172 */
1173 void loadConfig( app::appConfigurator &config /**< [in] the app::configurator object*/ );
1174};
1175
1176template <typename realT, class inputSpectT, typename iosT>
1178{
1179 wfsBeta<wfs<realT, iosT>>( nullptr ); // allocate a default ideal WFS.
1180}
1181
1182template <typename realT, class inputSpectT, typename iosT>
1184{
1185 if( m_wfsBeta && m_ownWfsBeta )
1186 {
1187 delete m_wfsBeta;
1188 }
1189}
1190
1191template <typename realT, class inputSpectT, typename iosT>
1193{
1194 atm.loadGuyon2005();
1195
1196 F0( 1.75e9 * 0.25 * math::pi<realT>() * 64. * 0.18 ); // Converting to photons/sec
1197 lam_wfs( 0.55e-6 );
1198 lam_sci( 1.6e-6 );
1199 D( 8. );
1200 starMag( 5 );
1201
1202 // The rest of Guyon 05 is very ideal
1203 npix_wfs( std::vector<realT>( { 12868 } ) );
1204 ron_wfs( std::vector<realT>( { 0.0 } ) );
1205 Fbg( std::vector<realT>( { 0.0 } ) );
1206
1207 d_min( 8.0 / 1e3 ); // Allow super fine sampling
1208 minTauWFS( (realT)( 1. / 1e9 ) ); // Just be super fast.
1209 tauWFS( 1. / 1e9 ); // Just be super fast
1210
1211 m_specsChanged = true;
1212 m_dminChanged = true;
1213}
1214
1215template <typename realT, class inputSpectT, typename iosT>
1217{
1218 atm.loadLCO();
1219 F0( 4.2e10 );
1220 lam_wfs( 0.791e-6 );
1221 lam_sci( 0.656e-6 );
1222
1223 npix_wfs( std::vector<realT>( { 9024 } ) );
1224 ron_wfs( std::vector<realT>( { 0.57 } ) );
1225 Fbg( std::vector<realT>( { 0.22 } ) );
1226
1227 d_min( std::vector<realT>( { 6.5 / 48.0 } ) );
1228 minTauWFS( std::vector<realT>( { 1. / 3622. } ) );
1229 tauWFS( 1. / 3622. );
1230
1231 D( 6.5 );
1232
1233 m_specsChanged = true;
1234 m_dminChanged = true;
1235}
1236
1237template <typename realT, class inputSpectT, typename iosT>
1239{
1240
1241 loadMagAOX();
1242
1243 F0( 7.6e10 * 368.0 / ( 0.25 * math::pi<realT>() * 6.5 * 6.5 * ( 1. - 0.29 * 0.29 ) ) ); // Scaled up.
1244
1245 D( 25.4 );
1246
1247 m_specsChanged = true;
1248 m_dminChanged = true;
1249}
1250
1251template <typename realT, class inputSpectT, typename iosT>
1253{
1254 m_D = nD;
1255 psd.D( m_D );
1256
1257 m_specsChanged = true;
1258 m_dminChanged = true;
1259}
1260
1261template <typename realT, class inputSpectT, typename iosT>
1263{
1264 return m_D;
1265}
1266
1267template <typename realT, class inputSpectT, typename iosT>
1268void aoSystem<realT, inputSpectT, iosT>::d_min( const std::vector<realT> &nd )
1269{
1270 m_d_min.assign( nd.begin(), nd.end() );
1271 m_specsChanged = true;
1272 m_dminChanged = true;
1273}
1274
1275template <typename realT, class inputSpectT, typename iosT>
1277{
1278 if( m_d_min.size() < idx + 1 )
1279 {
1280 mxThrowException( err::sizeerr, "aoSystem::d_min", "idx larger than m_d_min" );
1281 }
1282
1283 m_d_min[idx] = nd;
1284
1285 m_specsChanged = true;
1286 m_dminChanged = true;
1287}
1288
1289template <typename realT, class inputSpectT, typename iosT>
1291{
1292 if( m_d_min.size() < idx + 1 )
1293 {
1294 mxThrowException( err::sizeerr, "aoSystem::d_min", "idx larger than m_d_min" );
1295 }
1296
1297 return m_d_min[idx];
1298}
1299
1300template <typename realT, class inputSpectT, typename iosT>
1302{
1303 return m_d_min;
1304}
1305
1306template <typename realT, class inputSpectT, typename iosT>
1308{
1309 m_optd = od;
1310 m_specsChanged = true;
1311 m_dminChanged = true;
1312}
1313
1314template <typename realT, class inputSpectT, typename iosT>
1316{
1317 return m_optd;
1318}
1319
1320template <typename realT, class inputSpectT, typename iosT>
1322{
1323 m_optd_delta = odd;
1324 m_specsChanged = true;
1325 m_dminChanged = true;
1326}
1327
1328template <typename realT, class inputSpectT, typename iosT>
1330{
1331 return m_optd_delta;
1332}
1333
1334template <typename realT, class inputSpectT, typename iosT>
1335template <typename wfsT>
1337{
1338 if( m_wfsBeta && m_ownWfsBeta )
1339 {
1340 delete m_wfsBeta;
1341 m_wfsBeta = nullptr;
1342 }
1343
1344 m_wfsBeta = (wfs<realT, iosT> *)&w;
1345 m_ownWfsBeta = false;
1346}
1347
1348template <typename realT, class inputSpectT, typename iosT>
1349template <typename wfsT>
1351{
1352 if( m_wfsBeta && m_ownWfsBeta )
1353 {
1354 delete m_wfsBeta;
1355 m_wfsBeta = nullptr;
1356 }
1357
1358 if( w )
1359 {
1360 m_wfsBeta = (wfs<realT, iosT> *)w;
1361 m_ownWfsBeta = false;
1362 }
1363 else
1364 {
1365 m_wfsBeta = new wfsT;
1366 m_ownWfsBeta = true;
1367 }
1368}
1369
1370template <typename realT, class inputSpectT, typename iosT>
1375
1376template <typename realT, class inputSpectT, typename iosT>
1378{
1379 if( m_wfsBeta == 0 )
1380 mxThrowException( err::paramnotset, "aoSystem::beta_p", "The WFS is not assigned." );
1381
1382 return m_wfsBeta->beta_p( m, n, m_D, d_opt(), atm.r_0( m_lam_sci ) );
1383}
1384
1385template <typename realT, class inputSpectT, typename iosT>
1387{
1388 if( m_wfsBeta == 0 )
1389 mxThrowException( err::paramnotset, "aoSystem::beta_r", "The WFS is not assigned." );
1390
1391 return m_wfsBeta->beta_r( m, n, m_D, d_opt(), atm.r_0( m_lam_sci ) );
1392}
1393
1394template <typename realT, class inputSpectT, typename iosT>
1396{
1397 m_opticalGain = og;
1398 m_specsChanged = true;
1399 m_dminChanged = true;
1400}
1401
1402template <typename realT, class inputSpectT, typename iosT>
1404{
1405 return m_opticalGain;
1406}
1407
1408template <typename realT, class inputSpectT, typename iosT>
1410{
1411 m_lam_wfs = nlam;
1412 m_specsChanged = true;
1413 m_dminChanged = true;
1414}
1415
1416template <typename realT, class inputSpectT, typename iosT>
1418{
1419 return m_lam_wfs;
1420}
1421
1422template <typename realT, class inputSpectT, typename iosT>
1423void aoSystem<realT, inputSpectT, iosT>::npix_wfs( const std::vector<realT> &npix )
1424{
1425 m_npix_wfs.resize( npix.size() );
1426 for( size_t n = 0; n < npix.size(); ++n )
1427 {
1428 m_npix_wfs[n] = npix[n];
1429 }
1430
1431 m_specsChanged = true;
1432 m_dminChanged = true;
1433}
1434
1435template <typename realT, class inputSpectT, typename iosT>
1437{
1438 if( m_npix_wfs.size() < idx + 1 )
1439 {
1440 mxThrowException( err::sizeerr, "aoSystem::npix_wfs", "idx larger than m_npix_wfs" );
1441 }
1442
1443 m_npix_wfs[idx] = npix;
1444
1445 m_specsChanged = true;
1446 m_dminChanged = true;
1447}
1448
1449template <typename realT, class inputSpectT, typename iosT>
1451{
1452 if( m_npix_wfs.size() < idx + 1 )
1453 {
1454 mxThrowException( err::sizeerr, "aoSystem::npix_wfs", "idx larger than m_npix_wfs" );
1455 }
1456
1457 return m_npix_wfs[idx];
1458}
1459
1460template <typename realT, class inputSpectT, typename iosT>
1462{
1463 return m_npix_wfs;
1464}
1465
1466template <typename realT, class inputSpectT, typename iosT>
1467void aoSystem<realT, inputSpectT, iosT>::ron_wfs( const std::vector<realT> &nron )
1468{
1469 m_ron_wfs.resize( nron.size() );
1470 for( size_t n = 0; n < nron.size(); ++n )
1471 {
1472 m_ron_wfs[n] = nron[n];
1473 }
1474
1475 m_specsChanged = true;
1476 m_dminChanged = true;
1477}
1478
1479template <typename realT, class inputSpectT, typename iosT>
1481{
1482 if( m_ron_wfs.size() < idx + 1 )
1483 {
1484 mxThrowException( err::sizeerr, "aoSystem::ron_wfs", "idx larger than m_ron_wfs" );
1485 }
1486
1487 m_ron_wfs[idx] = nron;
1488
1489 m_specsChanged = true;
1490 m_dminChanged = true;
1491}
1492
1493template <typename realT, class inputSpectT, typename iosT>
1495{
1496 if( m_ron_wfs.size() < idx + 1 )
1497 {
1498 mxThrowException( err::sizeerr, "aoSystem::ron_wfs", "idx larger than m_ron_wfs" );
1499 }
1500
1501 return m_ron_wfs[idx];
1502}
1503
1504template <typename realT, class inputSpectT, typename iosT>
1506{
1507 return m_ron_wfs;
1508}
1509
1510template <typename realT, class inputSpectT, typename iosT>
1511void aoSystem<realT, inputSpectT, iosT>::Fbg( const std::vector<realT> &fbg )
1512{
1513 m_Fbg.resize( fbg.size() );
1514 for( size_t n = 0; n < fbg.size(); ++n )
1515 {
1516 m_Fbg[n] = fbg[n];
1517 }
1518
1519 m_specsChanged = true;
1520 m_dminChanged = true;
1521}
1522
1523template <typename realT, class inputSpectT, typename iosT>
1525{
1526 if( m_Fbg.size() < idx + 1 )
1527 {
1528 mxThrowException( err::sizeerr, "aoSyste,::Fbg", "idx larger than m_Fbg" );
1529 }
1530
1531 m_Fbg[idx] = fbg;
1532
1533 m_specsChanged = true;
1534 m_dminChanged = true;
1535}
1536
1537template <typename realT, class inputSpectT, typename iosT>
1539{
1540 if( m_Fbg.size() < idx + 1 )
1541 {
1542 mxThrowException( err::sizeerr, "aoSyste,::Fbg", "idx larger than m_Fbg" );
1543 }
1544
1545 return m_Fbg[idx];
1546}
1547
1548template <typename realT, class inputSpectT, typename iosT>
1550{
1551 return m_Fbg;
1552}
1553
1554template <typename realT, class inputSpectT, typename iosT>
1555void aoSystem<realT, inputSpectT, iosT>::minTauWFS( const std::vector<realT> &ntau )
1556{
1557 m_minTauWFS.resize( ntau.size() );
1558 for( size_t n = 0; n < ntau.size(); ++n )
1559 {
1560 m_minTauWFS[n] = ntau[n];
1561 }
1562
1563 m_specsChanged = true;
1564 m_dminChanged = true;
1565}
1566
1567template <typename realT, class inputSpectT, typename iosT>
1569{
1570 if( m_minTauWFS.size() < idx + 1 )
1571 {
1572 mxThrowException( err::sizeerr, "aoSystem::minTauWFS", "idx larger than m_ntau_wfs" );
1573 }
1574
1575 m_minTauWFS[idx] = ntau;
1576
1577 m_specsChanged = true;
1578 m_dminChanged = true;
1579}
1580
1581template <typename realT, class inputSpectT, typename iosT>
1583{
1584 if( m_minTauWFS.size() < idx + 1 )
1585 {
1586 mxThrowException( err::sizeerr, "aoSystem::minTauWFS", "idx larger than m_ntau_wfs" );
1587 }
1588
1589 return m_minTauWFS[idx];
1590}
1591
1592template <typename realT, class inputSpectT, typename iosT>
1594{
1595 return m_minTauWFS;
1596}
1597
1598template <typename realT, class inputSpectT, typename iosT>
1600{
1601 if( bnp != m_bin_npix )
1602 {
1603 m_bin_npix = bnp;
1604 m_specsChanged = true;
1605 m_dminChanged = true;
1606 }
1607}
1608
1609template <typename realT, class inputSpectT, typename iosT>
1611{
1612 return m_bin_npix;
1613}
1614
1615template <typename realT, class inputSpectT, typename iosT>
1617{
1618 d_opt();
1619 return m_bin_opt;
1620}
1621
1622template <typename realT, class inputSpectT, typename iosT>
1624{
1625 m_tauWFS = ntau;
1626 m_specsChanged = true;
1627 m_dminChanged = true;
1628}
1629
1630template <typename realT, class inputSpectT, typename iosT>
1632{
1633 return m_tauWFS;
1634}
1635
1636template <typename realT, class inputSpectT, typename iosT>
1638{
1639 m_deltaTau = ndel;
1640 m_specsChanged = true;
1641 m_dminChanged = true;
1642}
1643
1644template <typename realT, class inputSpectT, typename iosT>
1646{
1647 return m_deltaTau;
1648}
1649
1650template <typename realT, class inputSpectT, typename iosT>
1652{
1653 m_optTau = ot;
1654 m_specsChanged = true;
1655 m_dminChanged = true;
1656}
1657
1658template <typename realT, class inputSpectT, typename iosT>
1660{
1661 return m_optTau;
1662}
1663
1664template <typename realT, class inputSpectT, typename iosT>
1666{
1667 m_lam_sci = nlam;
1668 m_specsChanged = true;
1669 m_dminChanged = true;
1670}
1671
1672template <typename realT, class inputSpectT, typename iosT>
1674{
1675 return m_lam_sci;
1676}
1677
1678template <typename realT, class inputSpectT, typename iosT>
1680{
1681 m_zeta = nz;
1682 m_secZeta = 1 / cos( m_zeta );
1683
1684 m_specsChanged = true;
1685 m_dminChanged = true;
1686}
1687
1688template <typename realT, class inputSpectT, typename iosT>
1690{
1691 return m_zeta;
1692}
1693
1694template <typename realT, class inputSpectT, typename iosT>
1696{
1697 return m_secZeta;
1698}
1699
1700template <typename realT, class inputSpectT, typename iosT>
1702{
1703 if( mnm < 0 )
1704 mnm = 0;
1705 m_fit_mn_max = mnm;
1706}
1707
1708template <typename realT, class inputSpectT, typename iosT>
1710{
1711 return m_fit_mn_max;
1712}
1713
1714template <typename realT, class inputSpectT, typename iosT>
1716{
1717 m_circularLimit = cl;
1718 m_specsChanged = true;
1719 m_dminChanged = true;
1720}
1721
1722template <typename realT, class inputSpectT, typename iosT>
1724{
1725 return m_circularLimit;
1726}
1727
1728template <typename realT, class inputSpectT, typename iosT>
1730{
1731 m_spatialFilter_ku = fabs( kx );
1732 m_specsChanged = true; // not sure if needed
1733 m_dminChanged = true; // not sure if needed
1734}
1735
1736template <typename realT, class inputSpectT, typename iosT>
1738{
1739 return m_spatialFilter_ku;
1740}
1741
1742template <typename realT, class inputSpectT, typename iosT>
1744{
1745 m_spatialFilter_kv = fabs( ky );
1746 m_specsChanged = true; // not sure if needed
1747 m_dminChanged = true; // not sure if needed
1748}
1749
1750template <typename realT, class inputSpectT, typename iosT>
1752{
1753 return m_spatialFilter_kv;
1754}
1755
1756template <typename realT, class inputSpectT, typename iosT>
1758{
1759 m_ncp_wfe = nwfe;
1760 m_specsChanged = true;
1761 m_dminChanged = true;
1762}
1763
1764template <typename realT, class inputSpectT, typename iosT>
1766{
1767 return m_ncp_wfe;
1768}
1769
1770template <typename realT, class inputSpectT, typename iosT>
1772{
1773 m_ncp_alpha = alpha;
1774 m_specsChanged = true;
1775 m_dminChanged = true;
1776}
1777
1778template <typename realT, class inputSpectT, typename iosT>
1780{
1781 return m_ncp_alpha;
1782}
1783
1784template <typename realT, class inputSpectT, typename iosT>
1786{
1787 m_F0 = nF0;
1788 m_specsChanged = true;
1789 m_dminChanged = true;
1790}
1791
1792template <typename realT, class inputSpectT, typename iosT>
1794{
1795 return m_F0;
1796}
1797
1798template <typename realT, class inputSpectT, typename iosT>
1800{
1801 m_starMag = nmag;
1802 m_specsChanged = true;
1803 m_dminChanged = true;
1804}
1805
1806template <typename realT, class inputSpectT, typename iosT>
1808{
1809 return m_starMag;
1810}
1811
1812template <typename realT, class inputSpectT, typename iosT>
1814{
1815 return m_F0 * pow( 10.0, -0.4 * mag );
1816}
1817
1818template <typename realT, class inputSpectT, typename iosT>
1820{
1821 return Fg( m_starMag );
1822}
1823
1824template <typename realT, class inputSpectT, typename iosT>
1826{
1827 if( m_specsChanged || m_dminChanged )
1828 calcStrehl();
1829
1830 return m_controlledModes;
1831}
1832
1833template <typename realT, class inputSpectT, typename iosT>
1834realT aoSystem<realT, inputSpectT, iosT>::signal2Noise2( realT &Nph, realT tau_wfs, realT d, int b )
1835{
1836 Nph = Fg() * tau_wfs;
1837
1838 double binfact = 1.0;
1839 int binidx = 0;
1840 if( m_bin_npix )
1841 {
1842 if( m_npix_wfs.size() == 1 ) // Only if "true binning", otherwise the WFS mode vectors handle it.
1843 {
1844 binfact = 1. / pow( (realT)b + 1, 2 );
1845 }
1846 else
1847 binidx = b;
1848 }
1849
1850 return pow( Nph, 2 ) / ( m_npix_wfs[binidx] * binfact * ( m_Fbg[binidx] * tau_wfs + pow( m_ron_wfs[binidx], 2 ) ) );
1851}
1852
1853template <typename realT, class inputSpectT, typename iosT>
1854realT aoSystem<realT, inputSpectT, iosT>::measurementError( realT m, realT n, realT d, int b )
1855{
1856 if( m == 0 and n == 0 )
1857 return 0;
1858
1859 realT tau_wfs;
1860
1861 if( m_optTau )
1862 tau_wfs = optimumTauWFS( m, n, d, b );
1863 else
1864 tau_wfs = m_tauWFS;
1865
1866 if( m_wfsBeta == 0 )
1867 mxThrowException( err::paramnotset, "aoSystem::measurementError", "The WFS is not assigned." );
1868
1869 realT beta_p = m_wfsBeta->beta_p( m, n, m_D, d, atm.r_0( m_lam_wfs ) ) / sqrt( m_opticalGain );
1870 realT beta_r = m_wfsBeta->beta_r( m, n, m_D, d, atm.r_0( m_lam_wfs ) ); // sqrt(m_opticalGain);
1871 realT Nph = 0;
1872 realT snr2 = signal2Noise2( Nph, tau_wfs, d, b );
1873
1874 return ( pow( beta_r, 2 ) / snr2 + pow( beta_p, 2 ) / Nph ) * pow( m_lam_wfs / m_lam_sci, 2 );
1875}
1876
1877template <typename realT, class inputSpectT, typename iosT>
1879{
1880 if( m_specsChanged || m_dminChanged )
1881 calcStrehl();
1882 return m_wfeMeasurement;
1883}
1884
1885template <typename realT, class inputSpectT, typename iosT>
1886realT aoSystem<realT, inputSpectT, iosT>::timeDelayError( realT m, realT n, realT d, int b )
1887{
1888 if( m == 0 and n == 0 )
1889 return 0;
1890
1891 realT k = sqrt( m * m + n * n ) / m_D;
1892
1893 realT tau_wfs;
1894
1895 if( m_optTau )
1896 tau_wfs = optimumTauWFS( m, n, d, b );
1897 else
1898 tau_wfs = m_tauWFS;
1899
1900 realT tau = tau_wfs + m_deltaTau;
1901
1902 ///\todo handle multiple layers
1903 return psd( atm, 0, k, m_secZeta ) / pow( m_D, 2 ) * pow( atm.lam_0() / m_lam_sci, 2 ) *
1904 sqrt( atm.X( k, m_lam_sci, m_secZeta ) ) * pow( math::two_pi<realT>() * atm.v_wind() * k, 2 ) *
1905 pow( tau, 2 );
1906}
1907
1908template <typename realT, class inputSpectT, typename iosT>
1910{
1911 if( m_specsChanged || m_dminChanged )
1912 calcStrehl();
1913 return m_wfeTimeDelay;
1914}
1915
1916template <typename realT, class inputSpectT, typename iosT>
1918{
1919 realT k = sqrt( m * m + n * n ) / m_D;
1920
1921 ///\todo handle multiple layers
1922 return psd( atm, 0, k, m_secZeta ) / pow( m_D, 2 ) * pow( atm.lam_0() / m_lam_sci, 2 );
1923}
1924
1925template <typename realT, class inputSpectT, typename iosT>
1927{
1928 if( m_specsChanged || m_dminChanged )
1929 calcStrehl();
1930 return m_wfeFitting;
1931}
1932
1933template <typename realT, class inputSpectT, typename iosT>
1935{
1936 return C4var( m, n );
1937}
1938
1939template <typename realT, class inputSpectT, typename iosT>
1941{
1942 if( m_specsChanged || m_dminChanged )
1943 calcStrehl();
1944 return m_wfeChromScintOPD;
1945}
1946
1947template <typename realT, class inputSpectT, typename iosT>
1949{
1950 return C5var( m, n );
1951
1952#if 0
1953 int mn_max = floor(0.5*m_D/d_opt());
1954
1955 realT sum = 0;
1956
1957 for(int m = -mn_max; m <= mn_max; ++m)
1958 {
1959 for(int n = -mn_max; n <= mn_max; ++n)
1960 {
1961 if(n == 0 && m == 0) continue;
1962
1963 sum += C5var(m,n);
1964 }
1965 }
1966
1967 return sum;
1968#endif
1969}
1970
1971template <typename realT, class inputSpectT, typename iosT>
1973{
1974 return C6var( m, n );
1975}
1976
1977template <typename realT, class inputSpectT, typename iosT>
1979{
1980 if( m_specsChanged || m_dminChanged )
1981 calcStrehl();
1982 return m_wfeChromIndex;
1983}
1984
1985template <typename realT, class inputSpectT, typename iosT>
1987{
1988 return C7var( m, n );
1989}
1990
1991template <typename realT, class inputSpectT, typename iosT>
1993{
1994 if( m_specsChanged || m_dminChanged )
1995 calcStrehl();
1996 return m_wfeAnisoOPD;
1997}
1998
1999template <typename realT, class inputSpectT, typename iosT>
2001{
2002 return 0;
2003
2004#if 0
2005 int mn_max = floor(0.5*m_D/d_opt());
2006
2007 realT sum = 0;
2008
2009 for(int m = -mn_max; m <= mn_max; ++m)
2010 {
2011 for(int n = -mn_max; n <= mn_max; ++n)
2012 {
2013 if(n == 0 && m == 0) continue;
2014
2015 if( m_circularLimit )
2016 {
2017 if( m*m + n*n > mn_max*mn_max ) continue;
2018 }
2019
2020 sum += C8var(m,n);
2021 }
2022 }
2023
2024 return sum;
2025#endif
2026}
2027
2028template <typename realT, class inputSpectT, typename iosT>
2030 realT m,
2031 realT n,
2032 realT dact, // here called dact due to parameter collision with root-finding
2033 int bbin )
2034{
2035 if( m_D == 0 )
2036 {
2037 mxError( "aoSystem::optimumTauWFS", MXE_PARAMNOTSET, "Diameter (D) not set." );
2038 return -1;
2039 }
2040
2041 if( m_F0 == 0 )
2042 {
2043 mxError( "aoSystem::optimumTauWFS", MXE_PARAMNOTSET, "0-mag photon flux (F0) not set." );
2044 return -1;
2045 }
2046
2047 double binfact = 1.0;
2048 int binidx = 0; // index into WFS configurations
2049 if( m_bin_npix )
2050 {
2051 if( m_npix_wfs.size() > 1 )
2052 {
2053 binidx = bbin;
2054 }
2055 else
2056 {
2057 binfact = 1.0 / pow( (realT)( bbin + 1 ), 2 );
2058 }
2059 }
2060
2061 realT k = sqrt( m * m + n * n ) / m_D;
2062
2063 realT F = Fg();
2064
2065 if( m_wfsBeta == 0 )
2066 mxThrowException( err::paramnotset, "aoSystem::beta_p", "The WFS is not assigned." );
2067
2068 realT beta_p = m_wfsBeta->beta_p( m, n, m_D, dact, atm.r_0( m_lam_wfs ) ) / sqrt( m_opticalGain );
2069 realT beta_r = m_wfsBeta->beta_r( m, n, m_D, dact, atm.r_0( m_lam_wfs ) ) / sqrt( m_opticalGain );
2070
2071 // Set up for root finding:
2072 realT a, b, c, d, e;
2073
2074 ///\todo handle multiple layers
2075 realT Atmp = 2 * pow( atm.lam_0(), 2 ) * psd( atm, 0, k, m_secZeta ) / pow( m_D, 2 ) *
2076 ( atm.X( k, m_lam_wfs, m_secZeta ) ) * pow( math::two_pi<realT>() * atm.v_wind() * k, 2 );
2077
2078 a = Atmp;
2079 b = Atmp * m_deltaTau;
2080 c = 0;
2081 d = -1 * ( pow( m_lam_wfs, 2 ) / F ) *
2082 ( ( m_npix_wfs[binidx] * binfact * m_Fbg[binidx] / F ) * pow( beta_r, 2 ) + pow( beta_p, 2 ) );
2083 e = -2 * pow( m_lam_wfs, 2 ) * ( m_npix_wfs[binidx] * binfact ) * pow( m_ron_wfs[binidx], 2 ) * pow( beta_r, 2 ) /
2084 pow( F, 2 );
2085
2086 std::vector<std::complex<realT>> x;
2087
2088 // Get the roots
2089 math::quarticRoots( x, a, b, c, d, e );
2090
2091 // Now pick the largest positive real root
2092 realT tauopt = 0.0;
2093
2094 for( int i = 0; i < 4; i++ )
2095 {
2096 if( real( x[i] ) > 0 && imag( x[i] ) == 0 && real( x[i] ) > tauopt )
2097 tauopt = real( x[i] );
2098 }
2099
2100 if( tauopt < m_minTauWFS[binidx] )
2101 tauopt = m_minTauWFS[binidx];
2102
2103 return tauopt;
2104}
2105
2106template <typename realT, class inputSpectT, typename iosT>
2108{
2109 d_opt(); // also gets m_bin_opt;
2110 return optimumTauWFS( m, n, m_d_opt, m_bin_opt );
2111}
2112
2113template <typename realT, class inputSpectT, typename iosT>
2115{
2116 if( !m_dminChanged )
2117 return m_d_opt;
2118
2119 // Validate the modes
2120 ///\todo this should be in a separate valideModes function
2121 if( m_npix_wfs.size() < 1 )
2122 {
2123 mxThrowException( err::sizeerr, "aoSystem::d_opt", "npix_wfs must have at least one entry" );
2124 }
2125
2126 if( m_d_min.size() != m_npix_wfs.size() )
2127 {
2128 mxThrowException( err::sizeerr, "aoSystem::d_opt", "d_min must be the same size as npix_wfs" );
2129 }
2130
2131 if( m_ron_wfs.size() != m_npix_wfs.size() )
2132 {
2133 mxThrowException( err::sizeerr, "aoSystem::d_opt", "ron_wfs must be the same size as npix_wfs" );
2134 }
2135
2136 if( m_Fbg.size() != m_npix_wfs.size() )
2137 {
2138 mxThrowException( err::sizeerr, "aoSystem::d_opt", "F_bg must be the same size as npix_wfs" );
2139 }
2140
2141 if( m_minTauWFS.size() != m_npix_wfs.size() )
2142 {
2143 mxThrowException( err::sizeerr, "aoSystem::d_opt", "minTauWFS must be the same size as npix_wfs" );
2144 }
2145
2146 ///\todo investigate why we tolerate m_d_min = 0 here. Is this just a config error?
2147 if( m_d_min.size() == 1 && m_d_min[0] == 0 )
2148 {
2149 m_d_opt = 1e-50;
2150 m_bin_opt = 0;
2151 m_dminChanged = false;
2152 return m_d_opt;
2153 }
2154
2155 if( !m_optd )
2156 {
2157 m_d_opt = m_d_min[0];
2158 m_bin_opt = 0;
2159 m_dminChanged = false;
2160
2161 return m_d_opt;
2162 }
2163
2164 realT best_d = 0;
2165
2166 if( m_bin_npix )
2167 {
2168 if( m_npix_wfs.size() > 1 ) // Optimize over the WFS modes
2169 {
2170 int best_idx = 0;
2171 best_d = m_d_min[0];
2172
2173 realT bestStrehlOverall = 0;
2174
2175 for( int b = 0; b < m_npix_wfs.size(); ++b )
2176 {
2177 realT d = m_d_min[b];
2178 realT s = strehl( d, b );
2179 realT bestStrehl = s;
2180
2181 // Find the actuator pitch at which AO error is less than fitting error at Nyquist
2182 while( d < m_D / 2 )
2183 {
2184 d += m_d_min[b] / m_optd_delta;
2185 s = strehl( d, b );
2186
2187 if( s < bestStrehl ) // Strehl got worse. Can't be <= otherwise small m_optd_deltas will
2188 // immediately out
2189 {
2190 d -= m_d_min[b] / m_optd_delta; // go back to previous d
2191 break;
2192 }
2193 else
2194 bestStrehl = s;
2195 }
2196
2197 // Check if this is optimum so far
2198 if( bestStrehl >= bestStrehlOverall ) // >= to favor fewer actuators
2199 {
2200 bestStrehlOverall = bestStrehl;
2201 best_idx = b;
2202 best_d = d;
2203 }
2204 }
2205
2206 m_bin_opt = best_idx;
2207 }
2208 else
2209 {
2210 realT d = m_d_min[0];
2211 m_bin_opt = 0;
2212
2213 realT s = strehl( d, m_bin_opt );
2214 realT bestStrehl = s;
2215
2216 while( d < m_D / 2 )
2217 {
2218 d += m_d_min[0] / m_optd_delta;
2219 m_bin_opt = d / m_d_min[0] - 1;
2220 if( m_bin_opt < 0 )
2221 m_bin_opt = 0;
2222
2223 s = strehl( d, m_bin_opt );
2224
2225 if( s < bestStrehl ) // Strehl got worse. Can't be <= otherwise small m_optd_deltas will immediately
2226 // out
2227 {
2228 d -= m_d_min[0] / m_optd_delta; // go back to previous d
2229 m_bin_opt = d / m_d_min[0] - 1;
2230 if( m_bin_opt < 0 )
2231 m_bin_opt = 0;
2232 break;
2233 }
2234 else
2235 bestStrehl = s;
2236 }
2237 best_d = d;
2238 }
2239 }
2240 else
2241 {
2242 realT d = m_d_min[0];
2243 m_bin_opt = 0;
2244
2245 realT s = strehl( d, m_bin_opt );
2246 realT bestStrehl = s;
2247
2248 // Find the actuator pitch which minimizes total error (AO error is less than fitting error at Nyquist)
2249 while( d < m_D / 2 )
2250 {
2251 d += m_d_min[0] / m_optd_delta;
2252
2253 s = strehl( d, m_bin_opt );
2254 if( s < bestStrehl ) // Strehl got worse. Can't be <= otherwise small m_optd_deltas will immediately out
2255 {
2256 d -= m_d_min[0] / m_optd_delta; // go back to previous d
2257 break;
2258 }
2259 else
2260 bestStrehl = s;
2261 }
2262 best_d = d;
2263 }
2264
2265 m_d_opt = best_d;
2266 m_dminChanged = false;
2267
2268 return m_d_opt;
2269}
2270
2271template <typename realT, class inputSpectT, typename iosT>
2273{
2274 if( m == 0 and n == 0 )
2275 return 0;
2276
2277 realT k = sqrt( m * m + n * n ) / m_D;
2278
2279 realT kmax = m_fit_mn_max / m_D;
2280 realT kmin = 1. / m_D;
2281
2282 realT beta;
2283 if( m_ncp_alpha != 2 )
2284 {
2285 beta = ( m_ncp_alpha - 2 ) / ( math::two_pi<realT>() ) * 1. /
2286 ( pow( kmin, -m_ncp_alpha + 2 ) - pow( kmax, -m_ncp_alpha + 2 ) );
2287 }
2288 else
2289 {
2290 beta = 1. / ( math::two_pi<realT>() * log( kmax / kmin ) );
2291 }
2292 return beta * ncpError() * pow( k, -m_ncp_alpha ) * pow( kmin, 2 );
2293}
2294
2295template <typename realT, class inputSpectT, typename iosT>
2297{
2298 return pow( m_ncp_wfe, 2 ) * pow( math::two_pi<realT>() / m_lam_sci, 2 );
2299}
2300
2301template <typename realT, class inputSpectT, typename iosT>
2303{
2304
2305 int mn_max = floor( 0.5 * m_D / d );
2306
2307 realT wfeVar = 0;
2308
2309 bool controlled;
2310 for( int m = -m_fit_mn_max; m <= m_fit_mn_max; ++m )
2311 {
2312 for( int n = -m_fit_mn_max; n <= m_fit_mn_max; ++n )
2313 {
2314 controlled = false;
2315 if( m_circularLimit )
2316 {
2317 if( m * m + n * n <= mn_max * mn_max )
2318 controlled = true;
2319 }
2320 else
2321 {
2322 if( abs( m ) <= mn_max && abs( n ) <= mn_max )
2323 controlled = true;
2324 }
2325
2326 if( controlled )
2327 {
2328 realT wfeMeasurement = measurementError( m, n, d, b );
2329 realT wfeTimeDelay = timeDelayError( m, n, d, b );
2330 realT wfeChromScintOPD = chromScintOPDError( m, n );
2331 realT wfeChromIndex = chromIndexError( m, n );
2332 realT wfeAnisoOPD = dispAnisoOPDError( m, n );
2333
2334 realT wfeFitting = fittingError( m, n );
2335
2336 if( wfeFitting < wfeMeasurement + wfeTimeDelay + wfeChromScintOPD + wfeChromIndex + wfeAnisoOPD )
2337 {
2338 wfeVar += wfeFitting;
2339 }
2340 else
2341 {
2342 wfeVar += wfeMeasurement + wfeTimeDelay + wfeChromScintOPD + wfeChromIndex + wfeAnisoOPD;
2343 }
2344 }
2345 else
2346 {
2347 wfeVar += fittingError( m, n );
2348 }
2349 }
2350 }
2351
2352 wfeVar += ncpError();
2353 return exp( -1 * wfeVar );
2354}
2355
2356template <typename realT, class inputSpectT, typename iosT>
2358{
2359 realT d = d_opt(); // have to run this to get m_bin_opt set.
2360 int b = m_bin_opt;
2361
2362 int mn_max = floor( 0.5 * m_D / d );
2363
2364 m_wfeMeasurement = 0;
2365 m_wfeTimeDelay = 0;
2366 m_wfeChromScintOPD = 0;
2367 m_wfeChromIndex = 0;
2368 m_wfeAnisoOPD = 0;
2369
2370 m_wfeFitting = 0;
2371
2372 m_wfeNCP = 0;
2373
2374 m_controlledModes.resize( 2 * m_fit_mn_max + 1, 2 * m_fit_mn_max + 1 );
2375 m_controlledModes.setZero();
2376
2377 bool controlled;
2378 for( int m = -m_fit_mn_max; m <= m_fit_mn_max; ++m )
2379 {
2380 for( int n = -m_fit_mn_max; n <= m_fit_mn_max; ++n )
2381 {
2382 if( m == 0 && n == 0 )
2383 {
2384 continue;
2385 }
2386
2387 controlled = false;
2388 if( m_circularLimit )
2389 {
2390 if( m * m + n * n <= mn_max * mn_max )
2391 controlled = true;
2392 }
2393 else
2394 {
2395 if( abs( m ) <= mn_max && abs( n ) <= mn_max )
2396 controlled = true;
2397 }
2398
2399 if( controlled )
2400 {
2401 realT wfeMeasurement = measurementError( m, n, d, b );
2402 realT wfeTimeDelay = timeDelayError( m, n, d, b );
2403 realT wfeChromScintOPD = chromScintOPDError( m, n );
2404 realT wfeChromIndex = chromIndexError( m, n );
2405 realT wfeAnisoOPD = dispAnisoOPDError( m, n );
2406
2407 realT wfeFitting = fittingError( m, n );
2408
2409 if( wfeFitting < wfeMeasurement + wfeTimeDelay + wfeChromScintOPD + wfeChromIndex + wfeAnisoOPD )
2410 {
2411 m_wfeFitting += wfeFitting;
2412 }
2413 else // it's worth controlling this mode.
2414 {
2415 m_wfeMeasurement += wfeMeasurement;
2416 m_wfeTimeDelay += wfeTimeDelay;
2417 m_wfeChromScintOPD += wfeChromScintOPD;
2418 m_wfeChromIndex += wfeChromIndex;
2419 m_wfeAnisoOPD += wfeAnisoOPD;
2420 m_controlledModes( m_fit_mn_max + m, m_fit_mn_max + n ) = 1;
2421 }
2422 }
2423 else
2424 {
2425 m_wfeFitting += fittingError( m, n );
2426 }
2427 }
2428 }
2429
2430 m_wfeNCP = ncpError();
2431
2432 m_wfeVar = m_wfeMeasurement + m_wfeTimeDelay + m_wfeFitting + m_wfeChromScintOPD + m_wfeChromIndex + m_wfeAnisoOPD +
2433 m_wfeNCP;
2434
2435 m_strehl = exp( -1 * m_wfeVar );
2436
2437 m_specsChanged = false;
2438}
2439
2440template <typename realT, class inputSpectT, typename iosT>
2442{
2443 if( m_specsChanged || m_dminChanged )
2444 calcStrehl();
2445
2446 return m_wfeVar;
2447}
2448
2449template <typename realT, class inputSpectT, typename iosT>
2451{
2452 if( m_specsChanged || m_dminChanged )
2453 calcStrehl();
2454
2455 return m_strehl;
2456}
2457
2458template <typename realT, class inputSpectT, typename iosT>
2459template <typename varFuncT>
2460realT aoSystem<realT, inputSpectT, iosT>::C_( int m, int n, bool normStrehl, varFuncT varFunc, int doFittingError )
2461{
2462 if( m == 0 && n == 0 )
2463 return 0;
2464
2465 // we run this to optimize regardless of whether we use it
2466 // if already optimized then this is a minor hit
2467 realT S = strehl();
2468
2469 if( !normStrehl )
2470 S = 1;
2471
2472 if( doFittingError != FITTING_ERROR_NO )
2473 {
2474 int mn_max = m_D / ( 2. * d_opt() );
2475
2476 if( m_controlledModes( m_fit_mn_max + m, m_fit_mn_max + n ) == false )
2477 {
2478 if( doFittingError == FITTING_ERROR_ZERO )
2479 return 0;
2480
2481 realT fe = fittingError( m, n );
2482
2483 realT k = sqrt( m * m + n * n ) / m_D;
2484
2485 if( doFittingError == FITTING_ERROR_X )
2486 fe *= ( atm.X( k, m_lam_sci, m_secZeta ) );
2487 else if( doFittingError == FITTING_ERROR_Y )
2488 fe *= ( atm.Y( k, m_lam_sci, m_secZeta ) );
2489 else
2490 {
2491 std::cerr << "Unknown doFittingError\n";
2492 exit( -1 );
2493 }
2494
2495 return fe / S;
2496 }
2497 }
2498 // get if not doing fitting error or if inside control region:
2499
2500 realT var = ( this->*varFunc )( m, n );
2501
2502 return var / S;
2503}
2504
2505template <typename realT, class inputSpectT, typename iosT>
2506template <typename imageT, typename CfuncT>
2507void aoSystem<realT, inputSpectT, iosT>::C_Map( imageT &im, CfuncT Cfunc, bool normStrehl )
2508{
2509 int dim1 = im.rows();
2510 int dim2 = im.cols();
2511
2512 int mc = 0.5 * ( dim1 - 1 );
2513 int nc = 0.5 * ( dim2 - 1 );
2514
2515 for( int i = 0; i < dim1; ++i )
2516 {
2517 int m = i - mc;
2518
2519 for( int j = 0; j < dim2; ++j )
2520 {
2521 int n = j - nc;
2522
2523 im( i, j ) = ( this->*Cfunc )( m, n, normStrehl );
2524 }
2525 }
2526}
2527
2528template <typename realT, class inputSpectT, typename iosT>
2530{
2531 realT k = sqrt( m * m + n * n ) / m_D;
2532 ///\todo handle multiple layers
2533 return psd( atm, 0, k, m_secZeta ) / pow( m_D, 2 ) * pow( atm.lam_0() / m_lam_sci, 2 ) *
2534 ( atm.X( k, m_lam_sci, m_secZeta ) );
2535}
2536
2537template <typename realT, class inputSpectT, typename iosT>
2538realT aoSystem<realT, inputSpectT, iosT>::C0( realT m, realT n, bool normStrehl )
2539{
2540
2541 return C_( m, n, normStrehl, &aoSystem<realT, inputSpectT, iosT>::C0var, FITTING_ERROR_NO );
2542}
2543
2544template <typename realT, class inputSpectT, typename iosT>
2545template <typename imageT>
2546void aoSystem<realT, inputSpectT, iosT>::C0Map( imageT &im, bool normStrehl )
2547{
2548 C_Map( im, &aoSystem<realT, inputSpectT, iosT>::C0, normStrehl );
2549}
2550
2551template <typename realT, class inputSpectT, typename iosT>
2553{
2554 realT k = sqrt( m * m + n * n ) / m_D;
2555
2556 ///\todo handle multiple layers
2557 return psd( atm, 0, k, m_secZeta ) / pow( m_D, 2 ) * pow( atm.lam_0() / m_lam_sci, 2 ) *
2558 ( atm.Y( k, m_lam_sci, m_secZeta ) );
2559}
2560
2561template <typename realT, class inputSpectT, typename iosT>
2562realT aoSystem<realT, inputSpectT, iosT>::C1( realT m, realT n, bool normStrehl )
2563{
2564 return C_( m, n, normStrehl, &aoSystem<realT, inputSpectT, iosT>::C1var, FITTING_ERROR_NO );
2565}
2566
2567template <typename realT, class inputSpectT, typename iosT>
2568template <typename imageT>
2569void aoSystem<realT, inputSpectT, iosT>::C1Map( imageT &im, bool normStrehl )
2570{
2571 C_Map( im, &aoSystem<realT, inputSpectT, iosT>::C1, normStrehl );
2572}
2573
2574template <typename realT, class inputSpectT, typename iosT>
2576{
2577 d_opt();
2578 return measurementError( m, n, d_opt(), m_bin_opt ) + timeDelayError( m, n, d_opt(), m_bin_opt );
2579}
2580
2581template <typename realT, class inputSpectT, typename iosT>
2582realT aoSystem<realT, inputSpectT, iosT>::C2( realT m, realT n, bool normStrehl )
2583{
2584 return C_( m, n, normStrehl, &aoSystem<realT, inputSpectT, iosT>::C2var, FITTING_ERROR_X );
2585}
2586
2587template <typename realT, class inputSpectT, typename iosT>
2588template <typename imageT>
2589void aoSystem<realT, inputSpectT, iosT>::C2Map( imageT &im, bool normStrehl )
2590{
2591 C_Map( im, &aoSystem<realT, inputSpectT, iosT>::C2, normStrehl );
2592}
2593
2594template <typename realT, class inputSpectT, typename iosT>
2596{
2597 return 0; // measurementError(m, n) + timeDelayError(m,n);
2598}
2599
2600template <typename realT, class inputSpectT, typename iosT>
2601realT aoSystem<realT, inputSpectT, iosT>::C3( realT m, realT n, bool normStrehl )
2602{
2603 return C_( m, n, normStrehl, &aoSystem<realT, inputSpectT, iosT>::C3var, FITTING_ERROR_ZERO ); // FITTING_ERROR_Y);
2604}
2605
2606template <typename realT, class inputSpectT, typename iosT>
2607template <typename imageT>
2608void aoSystem<realT, inputSpectT, iosT>::C3Map( imageT &im, bool normStrehl )
2609{
2610 C_Map( im, &aoSystem<realT, inputSpectT, iosT>::C3, normStrehl );
2611}
2612
2613template <typename realT, class inputSpectT, typename iosT>
2615{
2616 realT k = sqrt( m * m + n * n ) / m_D;
2617
2618 ///\todo handle multiple layers
2619 // This does not need to be divided by X b/c we haven't multiplied by it, this is is C0/X.
2620 return psd( atm, 0, k, m_secZeta ) / pow( m_D, 2 ) * pow( atm.lam_0() / m_lam_sci, 2 ) *
2621 atm.dX( k, m_lam_sci, m_lam_wfs );
2622}
2623
2624template <typename realT, class inputSpectT, typename iosT>
2625realT aoSystem<realT, inputSpectT, iosT>::C4( realT m, realT n, bool normStrehl )
2626{
2627 return C_( m, n, normStrehl, &aoSystem<realT, inputSpectT, iosT>::C4var, FITTING_ERROR_ZERO );
2628}
2629
2630template <typename realT, class inputSpectT, typename iosT>
2631template <typename imageT>
2632void aoSystem<realT, inputSpectT, iosT>::C4Map( imageT &im, bool normStrehl )
2633{
2634 C_Map( im, &aoSystem<realT, inputSpectT, iosT>::C4, normStrehl );
2635}
2636
2637template <typename realT, class inputSpectT, typename iosT>
2639{
2640 realT k = sqrt( m * m + n * n ) / m_D;
2641
2642 ///\todo handle multiple layers
2643 // This does not need to be divided by Y b/c we haven't multiplied by it, this is is C1/Y.
2644 return psd( atm, 0, k, m_secZeta ) / pow( m_D, 2 ) * pow( atm.lam_0() / m_lam_sci, 2 ) *
2645 atm.dY( k, m_lam_sci, m_lam_wfs );
2646}
2647
2648template <typename realT, class inputSpectT, typename iosT>
2649realT aoSystem<realT, inputSpectT, iosT>::C5( realT m, realT n, bool normStrehl )
2650{
2651 return C_( m, n, normStrehl, &aoSystem<realT, inputSpectT, iosT>::C5var, FITTING_ERROR_ZERO );
2652}
2653
2654template <typename realT, class inputSpectT, typename iosT>
2655template <typename imageT>
2656void aoSystem<realT, inputSpectT, iosT>::C5Map( imageT &im, bool normStrehl )
2657{
2658 C_Map( im, &aoSystem<realT, inputSpectT, iosT>::C5, normStrehl );
2659}
2660
2661template <typename realT, class inputSpectT, typename iosT>
2663{
2664 realT ni = atm.n_air( m_lam_sci );
2665 realT nw = atm.n_air( m_lam_wfs );
2666
2667 return C0var( m, n ) * pow( ( ni - nw ) / ( ni - 1 ), 2 ); // Fixes error in Eqn 29
2668}
2669
2670template <typename realT, class inputSpectT, typename iosT>
2671realT aoSystem<realT, inputSpectT, iosT>::C6( realT m, realT n, bool normStrehl )
2672{
2673 return C_( m, n, normStrehl, &aoSystem<realT, inputSpectT, iosT>::C6var, FITTING_ERROR_ZERO );
2674}
2675
2676template <typename realT, class inputSpectT, typename iosT>
2677template <typename imageT>
2678void aoSystem<realT, inputSpectT, iosT>::C6Map( imageT &im, bool normStrehl )
2679{
2680 C_Map( im, &aoSystem<realT, inputSpectT, iosT>::C6, normStrehl );
2681}
2682
2683template <typename realT, class inputSpectT, typename iosT>
2685{
2686 realT k = sqrt( m * m + n * n ) / m_D;
2687
2688 ///\todo this needs to handle multiple layers -- do we violate an assumption of L_0 isn't constant?
2689 return psd( atm, 0, k, m_secZeta ) / pow( m_D, 2 ) * pow( atm.lam_0() / m_lam_sci, 2 ) *
2690 atm.X_Z( k, m_lam_wfs, m_lam_sci, m_secZeta );
2691}
2692
2693template <typename realT, class inputSpectT, typename iosT>
2694realT aoSystem<realT, inputSpectT, iosT>::C7( realT m, realT n, bool normStrehl )
2695{
2696 return C_( m, n, normStrehl, &aoSystem<realT, inputSpectT, iosT>::C7var, FITTING_ERROR_ZERO );
2697}
2698
2699template <typename realT, class inputSpectT, typename iosT>
2700template <typename imageT>
2701void aoSystem<realT, inputSpectT, iosT>::C7Map( imageT &im, bool normStrehl )
2702{
2703 C_Map( im, &aoSystem<realT, inputSpectT, iosT>::C7, normStrehl );
2704}
2705
2706template <typename realT, class inputSpectT, typename iosT>
2708{
2709 ios << "# AO Params:\n";
2710 ios << "# D = " << D() << '\n';
2711
2712 if( m_d_min.size() > 0 )
2713 {
2714 ios << "# d_min = " << d_min( (size_t)0 );
2715 for( size_t n = 1; n < m_d_min.size(); ++n )
2716 ios << ',' << d_min( n );
2717 ios << '\n';
2718 }
2719 else
2720 ios << "# d_min = null\n";
2721
2722 ios << "# optd = " << std::boolalpha << m_optd << '\n';
2723 ios << "# d_opt_delta = " << optd_delta() << '\n';
2724 ios << "# lam_sci = " << lam_sci() << '\n';
2725 ios << "# F0 = " << F0() << '\n';
2726 ios << "# starMag = " << starMag() << '\n';
2727 ios << "# lam_sci = " << lam_sci() << '\n';
2728 ios << "# zeta = " << zeta() << '\n';
2729 ios << "# lam_wfs = " << lam_wfs() << '\n';
2730
2731 if( npix_wfs().size() > 0 )
2732 {
2733 ios << "# npix_wfs = " << npix_wfs( (size_t)0 );
2734 for( size_t n = 1; n < npix_wfs().size(); ++n )
2735 ios << ',' << npix_wfs( n );
2736 ios << '\n';
2737 }
2738 else
2739 ios << "# npix_wfs = null\n";
2740
2741 if( ron_wfs().size() > 0 )
2742 {
2743 ios << "# ron_wfs = " << ron_wfs( (size_t)0 );
2744 for( size_t n = 1; n < ron_wfs().size(); ++n )
2745 ios << ',' << ron_wfs( n );
2746 ios << '\n';
2747 }
2748 else
2749 ios << "# ron_wfs = null\n";
2750
2751 if( Fbg().size() > 0 )
2752 {
2753 ios << "# Fbg = " << Fbg( (size_t)0 );
2754 for( size_t n = 1; n < Fbg().size(); ++n )
2755 ios << ',' << Fbg( n );
2756 ios << '\n';
2757 }
2758 else
2759 ios << "# Fbg = null\n";
2760
2761 if( minTauWFS().size() > 0 )
2762 {
2763 ios << "# minTauWFS = " << minTauWFS( (size_t)0 );
2764 for( size_t n = 1; n < minTauWFS().size(); ++n )
2765 ios << ',' << minTauWFS( n );
2766 ios << '\n';
2767 }
2768 else
2769 ios << "# minTauWFS = null\n";
2770
2771 ios << "# bin_npix = " << std::boolalpha << m_bin_npix << '\n';
2772 ios << "# tauWFS = " << tauWFS() << '\n';
2773 ios << "# optTau = " << std::boolalpha << m_optTau << '\n';
2774 ios << "# deltaTau = " << deltaTau() << '\n';
2775 ios << "# fit_mn_max = " << m_fit_mn_max << '\n';
2776 ios << "# spatialFilter_ku = " << m_spatialFilter_ku << '\n';
2777 ios << "# spatialFilter_kv = " << m_spatialFilter_kv << '\n';
2778 ios << "# ncp_wfe = " << m_ncp_wfe << '\n';
2779 ios << "# ncp_alpha = " << m_ncp_alpha << '\n';
2780
2781 if( m_wfsBeta == 0 )
2782 mxThrowException( err::paramnotset, "aoSystem::beta_p", "The WFS is not assigned." );
2783 m_wfsBeta->dumpWFS( ios );
2784 psd.dumpPSD( ios );
2785 atm.dumpAtmosphere( ios );
2786
2787 dumpGitStatus( ios );
2788
2789 return ios;
2790}
2791
2792template <typename realT, class inputSpectT, typename iosT>
2794{
2795 using namespace mx::app;
2796
2797 // AO System configuration
2798 config.add( "aosys.wfs",
2799 "",
2800 "aosys.wfs",
2801 argType::Required,
2802 "aosys",
2803 "wfs",
2804 false,
2805 "string",
2806 "The WFS type: idealWFS, unmodPyWFS, asympModPyWFS, shwfs, calculatedWFS" );
2807 config.add( "aosys.wfs_beta_p",
2808 "",
2809 "aosys.wfs_beta_p",
2810 argType::Required,
2811 "aosys",
2812 "wfs_beta_p",
2813 false,
2814 "string",
2815 "The beta_p file path for calcualtedWFS" );
2816 config.add( "aosys.wfs_beta_r",
2817 "",
2818 "aosys.wfs_beta_r",
2819 argType::Required,
2820 "aosys",
2821 "wfs_beta_r",
2822 false,
2823 "string",
2824 "The beta_r file path for calcualtedWFS" );
2825 config.add( "aosys.wfs_sensitivity",
2826 "",
2827 "aosys.wfs_sensitivity",
2828 argType::Required,
2829 "aosys",
2830 "wfs_sensitivity",
2831 false,
2832 "bool",
2833 "Flag indicating that beta_p/beta_r are sensitivities (inverse) [default false]" );
2834 config.add(
2835 "aosys.D", "", "aosys.D", argType::Required, "aosys", "D", false, "real", "The telescope diameter [m]" );
2836 config.add( "aosys.d_min",
2837 "",
2838 "aosys.d_min",
2839 argType::Required,
2840 "aosys",
2841 "d_min",
2842 false,
2843 "real",
2844 "The minimum actuator spacing [m]" );
2845 config.add( "aosys.optd",
2846 "",
2847 "aosys.optd",
2848 argType::Optional,
2849 "aosys",
2850 "optd",
2851 false,
2852 "bool",
2853 "Whether or not the actuator spacing is optimized" );
2854 config.add( "aosys.optd_delta",
2855 "",
2856 "aosys.optd_delta",
2857 argType::Required,
2858 "aosys",
2859 "optd_delta",
2860 false,
2861 "bool",
2862 "The fractional change from d_min used in optimization. Set to 1 (default) for integer binnings, > 1 "
2863 "for finer sampling." );
2864 config.add( "aosys.F0",
2865 "",
2866 "aosys.F0",
2867 argType::Required,
2868 "aosys",
2869 "F0",
2870 false,
2871 "real",
2872 "Zero-mag photon flux, [photons/sec]" );
2873 config.add( "aosys.lam_wfs",
2874 "",
2875 "aosys.lam_wfs",
2876 argType::Required,
2877 "aosys",
2878 "lam_wfs",
2879 false,
2880 "real",
2881 "WFS wavelength [m]" );
2882 config.add( "aosys.npix_wfs",
2883 "",
2884 "aosys.npix_wfs",
2885 argType::Required,
2886 "aosys",
2887 "npix_wfs",
2888 false,
2889 "vector<real>",
2890 "The number of pixels in the WFS" );
2891 config.add( "aosys.ron_wfs",
2892 "",
2893 "aosys.ron_wfs",
2894 argType::Required,
2895 "aosys",
2896 "ron_wfs",
2897 false,
2898 "vector<real>",
2899 "WFS readout noise [photons/read]" );
2900 config.add( "aosys.bin_npix",
2901 "",
2902 "aosys.bin_npix",
2903 argType::Required,
2904 "aosys",
2905 "bin_npix",
2906 false,
2907 "bool",
2908 "Whether or not WFS pixels are re-binned along with actuator spacing optimization" );
2909 config.add( "aosys.Fbg",
2910 "",
2911 "aosys.Fbg",
2912 argType::Required,
2913 "aosys",
2914 "Fbg",
2915 false,
2916 "vector<real>",
2917 "Background counts, [counts/pix/sec]" );
2918 config.add( "aosys.tauWFS",
2919 "",
2920 "aosys.tauWFS",
2921 argType::Required,
2922 "aosys",
2923 "tauWFS",
2924 false,
2925 "real",
2926 "WFS integration time [s]" );
2927 config.add( "aosys.minTauWFS",
2928 "",
2929 "aosys.minTauWFS",
2930 argType::Required,
2931 "aosys",
2932 "minTauWFS",
2933 false,
2934 "vector<real>",
2935 "Minimum WFS integration time [s]" );
2936 config.add( "aosys.deltaTau",
2937 "",
2938 "aosys.deltaTau",
2939 argType::Required,
2940 "aosys",
2941 "deltaTau",
2942 false,
2943 "real",
2944 "Loop delay [s]" );
2945 config.add( "aosys.optTau",
2946 "",
2947 "aosys.optTau",
2948 argType::Optional,
2949 "aosys",
2950 "optTau",
2951 false,
2952 "bool",
2953 "Whether or not the integration time is optimized" );
2954 config.add( "aosys.lam_sci",
2955 "",
2956 "aosys.lam_sci",
2957 argType::Required,
2958 "aosys",
2959 "lam_sci",
2960 false,
2961 "real",
2962 "Science wavelength [m]" );
2963 config.add(
2964 "aosys.zeta", "", "aosys.zeta", argType::Required, "aosys", "zeta", false, "real", "Zenith distance [rad]" );
2965 config.add( "aosys.fit_mn_max",
2966 "",
2967 "aosys.fit_mn_max",
2968 argType::Required,
2969 "aosys",
2970 "fit_mn_max",
2971 false,
2972 "real",
2973 "Maximum spatial frequency index to use for analysis" );
2974 config.add( "aosys.circularLimit",
2975 "",
2976 "aosys.circularLimit",
2977 argType::Optional,
2978 "aosys",
2979 "circularLimit",
2980 false,
2981 "bool",
2982 " Flag to indicate that the spatial frequency limit is circular, not square." );
2983 config.add( "aosys.spatialFilter_ku",
2984 "",
2985 "aosys.spatialFilter_ku",
2986 argType::Required,
2987 "aosys",
2988 "spatialFilter_ku",
2989 false,
2990 "real",
2991 "Spatial filter cutoff frequency in u [m^-1]" );
2992 config.add( "aosys.spatialFilter_kv",
2993 "",
2994 "aosys.spatialFilter_kv",
2995 argType::Required,
2996 "aosys",
2997 "spatialFilter_kv",
2998 false,
2999 "real",
3000 "Spatial filter cutoff frequency in v [m^-1]" );
3001 config.add( "aosys.ncp_wfe",
3002 "",
3003 "aosys.ncp_wfe",
3004 argType::Required,
3005 "aosys",
3006 "ncp_wfe",
3007 false,
3008 "real",
3009 "NCP WFE between 1 lambda/D and fit_mn_max [rad^2]" );
3010 config.add( "aosys.ncp_alpha",
3011 "",
3012 "aosys.ncp_alpha",
3013 argType::Required,
3014 "aosys",
3015 "ncp_alpha",
3016 false,
3017 "real",
3018 "PSD index for NCP WFE" );
3019 config.add(
3020 "aosys.starMag", "", "aosys.starMag", argType::Required, "aosys", "starMag", false, "real", "Star magnitude" );
3021 config.add( "aosys.starMags",
3022 "",
3023 "aosys.starMags",
3024 argType::Required,
3025 "aosys",
3026 "starMags",
3027 false,
3028 "real vector",
3029 "A vector of star magnitudes" );
3030
3031 atm.setupConfig( config );
3032 psd.setupConfig( config );
3033}
3034
3035template <typename realT, class inputSpectT, typename iosT>
3037{
3038 // WFS
3039 if( config.isSet( "aosys.wfs" ) )
3040 {
3041 std::string wfsStr;
3042 config( wfsStr, "aosys.wfs" );
3043
3044 if( wfsStr == "ideal" )
3045 {
3046 wfsBeta<wfs<realT, iosT>>( nullptr );
3047 }
3048 else if( wfsStr == "unmodPyWFS" )
3049 {
3050 wfsBeta<pywfsUnmod<realT, iosT>>( nullptr );
3051 }
3052 else if( wfsStr == "asympModPyWFS" )
3053 {
3054 wfsBeta<pywfsModAsymptotic<realT>>( nullptr );
3055 }
3056 else if( wfsStr == "SHWFS" )
3057 {
3058 wfsBeta<shwfs<realT, iosT>>( nullptr );
3059 }
3060 else if( wfsStr == "calculatedWFS" )
3061 {
3062 wfsBeta<calculatedWFS<realT, iosT>>( nullptr );
3063
3064 calculatedWFS<realT, iosT> *cwfs = static_cast<calculatedWFS<realT, iosT> *>( m_wfsBeta );
3065 config( cwfs->m_beta_p_file, "aosys.wfs_beta_p" );
3066 config( cwfs->m_beta_r_file, "aosys.wfs_beta_r" );
3067 bool sens = cwfs->m_sensitivity;
3068 config( sens, "aosys.wfs_sensitivity" );
3069 if( config.isSet( "aosys.wfs_sensitivity" ) )
3070 {
3071 cwfs->m_sensitivity = sens;
3072 }
3073 }
3074 else
3075 {
3076 mxThrowException( err::invalidconfig, "aoSystem::loadConfig", "unknown WFS " + wfsStr + " specified" );
3077 }
3078 }
3079
3080 // We load the default value, get the config value (which will be the default if not set), and
3081 // then if it was set, call the setter function so any side effects are captured.
3082
3083 // diameter
3084 realT nD = D();
3085 config( nD, "aosys.D" );
3086 if( config.isSet( "aosys.D" ) )
3087 D( nD );
3088
3089 // d_min
3090 std::vector<realT> nd_min = d_min();
3091 config( nd_min, "aosys.d_min" );
3092 if( config.isSet( "aosys.d_min" ) )
3093 d_min( nd_min );
3094
3095 bool noptd = optd();
3096 config( noptd, "aosys.optd" );
3097 if( config.isSet( "aosys.optd" ) )
3098 optd( noptd );
3099
3100 realT noptd_delta = optd_delta();
3101 config( noptd_delta, "aosys.optd_delta" );
3102 if( config.isSet( "aosys.optd_delta" ) )
3103 optd_delta( noptd_delta );
3104
3105 realT nlam_wfs = lam_wfs();
3106 config( nlam_wfs, "aosys.lam_wfs" );
3107 if( config.isSet( "aosys.lam_wfs" ) )
3108 lam_wfs( nlam_wfs );
3109
3110 // npix_wfs
3111 std::vector<realT> nnpix_wfs = npix_wfs();
3112 config( nnpix_wfs, "aosys.npix_wfs" );
3113 if( config.isSet( "aosys.npix_wfs" ) )
3114 npix_wfs( nnpix_wfs );
3115
3116 // ron_wfs
3117 std::vector<realT> nron_wfs = ron_wfs();
3118 config( nron_wfs, "aosys.ron_wfs" );
3119 if( config.isSet( "aosys.ron_wfs" ) )
3120 ron_wfs( nron_wfs );
3121
3122 // Fbg
3123 std::vector<realT> nFbg = Fbg();
3124 config( nFbg, "aosys.Fbg" );
3125 if( config.isSet( "aosys.Fbg" ) )
3126 Fbg( nFbg );
3127
3128 // minTauWFS
3129 std::vector<realT> nminTauWFS = minTauWFS();
3130 config( nminTauWFS, "aosys.minTauWFS" );
3131 if( config.isSet( "aosys.minTauWFS" ) )
3132 minTauWFS( nminTauWFS );
3133
3134 // bin_npix
3135 bool nbin_npix = bin_npix();
3136 config( nbin_npix, "aosys.bin_npix" );
3137 if( config.isSet( "aosys.bin_npix" ) )
3138 bin_npix( nbin_npix );
3139
3140 // tauWFS
3141 realT ntauWFS = tauWFS();
3142 config( ntauWFS, "aosys.tauWFS" );
3143 if( config.isSet( "aosys.tauWFS" ) )
3144 tauWFS( ntauWFS );
3145
3146 // deltaTau
3147 realT ndeltaTau = deltaTau();
3148 config( ndeltaTau, "aosys.deltaTau" );
3149 if( config.isSet( "aosys.deltaTau" ) )
3150 deltaTau( ndeltaTau );
3151
3152 // optTau
3153 bool noptTau = optTau();
3154 config( noptTau, "aosys.optTau" );
3155 if( config.isSet( "aosys.optTau" ) )
3156 optTau( noptTau );
3157
3158 // lam_sci
3159 realT nlam_sci = lam_sci();
3160 config( nlam_sci, "aosys.lam_sci" );
3161 if( config.isSet( "aosys.lam_sci" ) )
3162 lam_sci( nlam_sci );
3163
3164 // zeta
3165 realT nzeta = zeta();
3166 config( nzeta, "aosys.zeta" );
3167 if( config.isSet( "aosys.zeta" ) )
3168 zeta( nzeta );
3169
3170 // fit_mn_max
3171 realT fmnm = fit_mn_max();
3172 config( fmnm, "aosys.fit_mn_max" );
3173 if( config.isSet( "aosys.fit_mn_max" ) )
3174 fit_mn_max( fmnm );
3175
3176 // circularLimit
3177 bool cl = circularLimit();
3178 config( cl, "aosys.circularLimit" );
3179 if( config.isSet( "aosys.circularLimit" ) )
3180 circularLimit( cl );
3181
3182 // spatialFilter_ku
3183 realT ku = spatialFilter_ku();
3184 config( ku, "aosys.spatialFilter_ku" );
3185 if( config.isSet( "aosys.spatialFilter_ku" ) )
3186 spatialFilter_ku( ku );
3187
3188 // spatialFilter_kv
3189 realT kv = spatialFilter_kv();
3190 config( kv, "aosys.spatialFilter_kv" );
3191 if( config.isSet( "aosys.spatialFilter_kv" ) )
3192 spatialFilter_kv( kv );
3193
3194 // ncp_wfe
3195 realT nwfe = ncp_wfe();
3196 config( nwfe, "aosys.ncp_wfe" );
3197 if( config.isSet( "aosys.ncp_wfe" ) )
3198 ncp_wfe( nwfe );
3199
3200 // ncp_alpha
3201 realT na = ncp_alpha();
3202 config( na, "aosys.ncp_alpha" );
3203 if( config.isSet( "aosys.ncp_alpha" ) )
3204 ncp_alpha( na );
3205
3206 // F0
3207 realT nF0 = F0();
3208 config( nF0, "aosys.F0" );
3209 if( config.isSet( "aosys.F0" ) )
3210 F0( nF0 );
3211
3212 // star_mag
3213 realT smag = starMag();
3214 config( smag, "aosys.starMag" );
3215 if( config.isSet( "aosys.starMag" ) )
3216 starMag( smag );
3217
3218 atm.loadConfig( config );
3219 psd.loadConfig( config );
3220}
3221
3222extern template class aoSystem<float, vonKarmanSpectrum<float>, std::ostream>;
3223
3224extern template class aoSystem<double, vonKarmanSpectrum<double>, std::ostream>;
3225
3226extern template class aoSystem<long double, vonKarmanSpectrum<long double>, std::ostream>;
3227
3228#ifdef HASQUAD
3229extern template class aoSystem<__float128, vonKarmanSpectrum<__float128>, std::ostream>;
3230#endif
3231
3232} // namespace analysis
3233} // namespace AO
3234} // namespace mx
3235
3236#endif // aoSystem_hpp
Provides a class to specify atmosphere parameters.
Calculate and provide constants related to adaptive optics.
Spatial power spectra used in adaptive optics.
Definitions of various analytic wavefront sensors.
A class to specify atmosphere parameters and perform related calculations.
Describes an analytic adaptive optics (AO) system.
Definition aoSystem.hpp:66
bool circularLimit()
Get the value of the circularLimit flag.
realT m_ncp_wfe
Static WFE [m rms].
Definition aoSystem.hpp:127
void loadMagAOX()
Load parameters corresponding to the MagAO-X system.
bool bin_npix()
Get the value of the pixel binngin flag.
void loadGuyon2005()
Load the default parameters from Guyon, 2005 .
bool optd()
Get the value of m_optd.
void C1Map(imageT &map, bool normStrehl)
Calculate a 2D map of contrast C1.
realT m_wfeNCP
Total WFE due to NCP errors [rad^2 at m_lam_sci].
Definition aoSystem.hpp:146
realT F0()
Get the value of the 0 magnitude photon rate.
realT C5(realT m, realT n, bool normStrehl=true)
Calculate the contrast due to scintilation-amplitude chromaticity.
realT chromScintOPDErrorTotal()
Calculate the wavefront error due to scintillation chromaticity in the OPD over all spatial frequenci...
void C5Map(imageT &map, bool normStrehl)
Calculate a 2D map of contrast C5.
realT C1(realT m, realT n, bool normStrehl=true)
Calculate the contrast due to uncorrected amplitude, C1.
realT C0var(realT m, realT n)
Calculate the residual variance due to uncorrected phase at a spatial frequency.
realT C_(int m, int n, bool normStrehl, varFuncT varFunc, int doFittingError)
Worker function for raw contrast fuctions.
realT m_wfeVar
The WFE variance, in meters^2. Never use this directly, instead use wfeVar().
Definition aoSystem.hpp:148
realT d_opt()
Calculate the optimum actuator spacing.
std::vector< realT > d_min()
Get the minimum subaperture sampling for all WFS modes.
realT ncpError()
Calculate the total NCP variance in rad^2.
realT m_opticalGain
The optical gain of the WFS, 0-1. Treated as a sensitivity reduction in measurement noise....
Definition aoSystem.hpp:90
realT spatialFilter_kv()
Get the value of spatialFilter_kv.
bool m_specsChanged
Flag to indicate that a specification has changed.
Definition aoSystem.hpp:134
bool m_circularLimit
Flag to indicate that the spatial frequency limit is circular, not square.
Definition aoSystem.hpp:122
realT zeta()
Get the zenith angle.
realT Fg()
Get the photon rate at the current Star magnitude.
std::vector< realT > minTauWFS()
Get the values of the minimum WFS exposure time.
std::vector< realT > m_npix_wfs
Number of WFS pixels. One per WFS mode.
Definition aoSystem.hpp:97
void calcStrehl()
Calculate the component WFE, total WFE, and Strehl ratio.
realT lam_sci()
Get the science wavelength.
realT m_wfeFitting
Total WFE due to fitting error [rad^2 at m_lam_sci].
Definition aoSystem.hpp:141
realT m_wfeChromScintOPD
Total WFE due to the chromaticity of scintillation OPD [rad^2 at lam_sc].
Definition aoSystem.hpp:142
realT m_strehl
Strehl ratio, a calculated quantity. Never use this directdy, instead use strehl().
Definition aoSystem.hpp:149
realT deltaTau()
Get the value of m_deltaTau.
realT wfeVar()
Get the current value of the total WFE variance.
realT ncp_wfe()
Get the value of the non-common path WFE.
realT m_spatialFilter_ku
The spatial filter cutoff in u, [m^-1].
Definition aoSystem.hpp:124
realT tauWFS()
Get the value of the minimum WFS exposure time.
std::vector< realT > m_ron_wfs
WFS readout noise [electrons/pix]. One per WFS mode.
Definition aoSystem.hpp:98
bool m_ownWfsBeta
Flag indicating if the WFS beta_p pointer is owned by this instance.
Definition aoSystem.hpp:88
realT m_secZeta
Secant of the Zenith angle (calculated)
Definition aoSystem.hpp:118
realT m_starMag
The magnitude of the star.
Definition aoSystem.hpp:132
realT C4(realT m, realT n, bool normStrehl=true)
Calculate the contrast due to scintilation-OPD chromaticity.
realT timeDelayError(realT m, realT n, realT d, int b)
Calculate the time delay at a spatial frequency at the optimum exposure time and the specified actuat...
realT chromIndexErrorTotal()
Calculate the wavefront error due to chromaticity in the index of refraction at a specific spatial fr...
realT chromIndexError(int m, int n)
Calculate the wavefront error due to chromaticity in the index of refraction at a given spatial frequ...
realT m_d_opt
Current optimum AO system actuator pitch [m].
Definition aoSystem.hpp:81
realT C4var(realT m, realT n)
Calculate the residual variance due to scintilation-OPD chromaticity.
realT m_ncp_alpha
Power-law exponent for the NCP aberations. Default is 2.
Definition aoSystem.hpp:128
realT beta_r(realT m, realT n)
Get the value of beta_r for a spatial frequency.
realT spatialFilter_ku()
Get the value of spatialFilter_ku.
realT m_spatialFilter_kv
The spatial filter cutoff in v, [m^-1].
Definition aoSystem.hpp:125
void C7Map(imageT &map, bool normStrehl)
Calculate a 2D map of contrast C7.
realT C7var(realT m, realT n)
Calculate the residual variance due to dispersive anisoplanatism.
realT fittingErrorTotal()
Calculate the total fitting error over all uncorrected spatial frequencies.
bool m_bin_npix
Flag controlling whether or not to bin WFS pixels according to the actuator spacing.
Definition aoSystem.hpp:102
void C0Map(imageT &map, bool normStrehl)
Calculate a 2D map of contrast C0.
bool optTau()
Get the value of m_optTau.
realT m_tauWFS
Actual WFS exposure time [sec].
Definition aoSystem.hpp:108
realT beta_p(realT m, realT n)
Get the value of beta_p for a spatial frequency.
std::vector< realT > m_minTauWFS
Minimum WFS exposure time [sec]. One per WFS mode.
Definition aoSystem.hpp:100
realT C5var(realT m, realT n)
Calculate the residual variance due to to scintilation-amplitude chromaticity.
wfs< realT, iosT > * m_wfsBeta
The WFS beta_p class.
Definition aoSystem.hpp:87
realT secZeta()
Get the zecant of the zenith angle.
realT measurementError(realT m, realT n, realT d, int b)
Calculate the measurement noise at a spatial frequency and specified actuator spacing.
realT dispAnisoOPDErrorTotal()
Calculate the wavefront error due to dispersive anisoplanatism in the OPD over all specific spatial f...
realT m_lam_wfs
WFS wavelength [m].
Definition aoSystem.hpp:94
void C3Map(imageT &map, bool normStrehl)
Calculate a 2D map of contrast C3.
realT m_F0
0 mag flux from star at WFS [photons/sec]
Definition aoSystem.hpp:130
realT C2var(realT m, realT n)
Calculate the residual variance due to measurement and time delay errors in phase/OPD at a spatial fr...
bool m_dminChanged
Flag to indicate that d_min has changed.
Definition aoSystem.hpp:135
std::vector< realT > Fbg()
Get the value of the background flux for all WFS modes.
realT timeDelayErrorTotal()
Calculate the time delay error over all corrected spatial frequencies.
realT dispAnisoOPDError(int m, int n)
Calculate the wavefront error due to dispersive anisoplanatism in the OPD at a given spatial frequenc...
realT m_zeta
Zenith angle [radians].
Definition aoSystem.hpp:117
realT signal2Noise2(realT &Nph, realT tau_wfs, realT d, int b)
Calculate the terms of the signal to noise ratio squared (S/N)^2 for the WFS measurement.
realT C3var(realT m, realT n)
Calculate the residual variance due to measurement and time delay errors in amplitude at a spatial fr...
realT starMag()
Get the value of the Star's magnitude.
void loadGMagAOX()
Load parameters corresponding to the G-MagAO-X system.
realT m_D
Telescope diameter [m].
Definition aoSystem.hpp:77
realT chromScintAmpError()
Calculate the wavefront error due to scintillation chromaticity in amplitude over all spatial frequen...
std::vector< realT > ron_wfs()
Get the value of the WFS readout noise for all WFS modes.
std::vector< realT > npix_wfs()
Get the number of pixels in the WFS for all WFS modes.
realT C7(realT m, realT n, bool normStrehl=true)
Calculate the contrast due to dispersive anisoplanatism.
realT measurementErrorTotal()
Calculate the total measurement error over all corrected spatial frequencies.
void setupConfig(app::appConfigurator &config)
Setup the configurator to configure this class.
realT strehl()
Get the current Strehl ratio.
realT C0(realT m, realT n, bool normStrehl=true)
Calculate the contrast due to uncorrected phase, C0.
realT C1var(realT m, realT n)
Calculate the residual variance due to uncorrected amplitude at a spatial frequency.
realT ncp_alpha()
Get the value of the non-common path WFE PSD index.
const Eigen::Array< int, -1, -1 > & controlledModes()
Access the array of controlledModes.
realT m_wfeTimeDelay
Total WFE due to time delay [rad^2 at m_lam_sci].
Definition aoSystem.hpp:140
realT chromScintOPDError(int m, int n)
Calculate the wavefront error due to scintillation chromaticity in the OPD at a spatial frequency.
void C_Map(imageT &im, CfuncT Cfunc, bool normStrehl)
Worker function for the contrast-map functions.
int m_fit_mn_max
Maximum spatial frequency index to use for fitting error calculation.
Definition aoSystem.hpp:120
iosT & dumpAOSystem(iosT &ios)
Output current parameters to a stream.
std::vector< realT > m_Fbg
Background flux, [photons/sec/pixel]. One per WFS mode.
Definition aoSystem.hpp:99
void C4Map(imageT &map, bool normStrehl)
Calculate a 2D map of contrast C4.
realT C6var(realT m, realT n)
Calculate the residual variance due to chromaticity of the index of refraction of air.
realT C2(realT m, realT n, bool normStrehl=true)
Calculate the contrast due to measurement and time delay errors in phase/OPD at a spatial frequency.
realT m_wfeMeasurement
Total WFE due to measurement a error [rad^2 at m_lam_sci].
Definition aoSystem.hpp:139
void C2Map(imageT &map, bool normStrehl)
Calculate a 2D map of contrast C2().
realT dispAnisoAmpError()
Calculate the wavefront error due to dispersive anisoplanatism in the amplitude at a specific spatial...
Eigen::Array< int, -1, -1 > m_controlledModes
Map of which modes are under control. Set by calcStrehl.
Definition aoSystem.hpp:137
std::vector< realT > m_d_min
Minimum AO system actuator pitch [m]. One per WFS mode.
Definition aoSystem.hpp:79
void loadConfig(app::appConfigurator &config)
Load the configuration of this class from a configurator.
realT C6(realT m, realT n, bool normStrehl=true)
Calculate the contrast due to chromaticity of the index of refraction of air.
int fit_mn_max()
Get the value of m_fit_mn_max.
realT opticalGain()
Get the value of the optical gain.
int bin_opt()
Get the value of the optimum binning factor/index.
realT m_wfeChromIndex
Total WFE due to the chromaticity of the index of refraction [rad^2 at lam_Sci].
Definition aoSystem.hpp:143
realT m_lam_sci
Science wavelength [m].
Definition aoSystem.hpp:115
realT optimumTauWFS(realT m, realT n, realT d, int b)
Calculate the optimum exposure time for a given spatial frequency at a specified actuator spacing.
realT m_deltaTau
Loop latency [sec].
Definition aoSystem.hpp:110
void C6Map(imageT &map, bool normStrehl)
Calculate a 2D map of contrast C6.
realT lam_wfs()
Get the value of the WFS wavelength.
realT D()
Get the value of the primary mirror diamter.
realT m_wfeAnisoOPD
Total WFE due to dispersive anisoplanatism OPD.
Definition aoSystem.hpp:144
realT C3(realT m, realT n, bool normStrehl=true)
Calculate the contrast due to measurement and time delay errors in amplitude at a spatial frequency.
realT fittingError(realT m, realT n)
Calculate the fitting error at a specific spatial frequency.
wfs< realT, iosT > * wfsBeta()
Get the WFS Beta pointer.
realT optd_delta()
Get the value of the fractional change in actuator spacing for optimization..
mxException for invalid config settings
mxException for parameters which aren't set
mxException for a size error
#define mxError(esrc, ecode, expl)
This reports an mxlib specific error.
#define mxThrowException(extype, src, expl)
Throw an exception. This macro takes care of the file and line.
void quarticRoots(std::vector< std::complex< realT > > &x, realT a, realT b, realT c, realT d, realT e)
Find the roots of the general quartic equation.
Definition roots.hpp:126
constexpr floatT six_fifths()
Return 6/5 in the specified precision.
#define MXE_PARAMNOTSET
A parameter was not set.
The mxlib c++ namespace.
Definition mxError.hpp:40
iosT & dumpGitStatus(iosT &ios, const std::string &repoName, const std::string &branch, const std::string &sha1, const bool &modified, const std::string &section="")
Dump the git status of a repository to a stream.
Definition mxlib.hpp:52
The calculated WFS uses sensitivities provided by FITS files.
Definition aoWFS.hpp:276
The ideal wavefront sensor sensitivity function.
Definition aoWFS.hpp:37
Class to manage a set of configurable values, and read their values from config/ini files and the com...
void add(const configTarget &tgt)
Add a configTarget.
bool isSet(const std::string &name, std::unordered_map< std::string, configTarget > &targets)
Check if a target has been set by the configuration.