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