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 "../../math/constants.hpp"
31#include "../../math/roots.hpp"
32#include "../../mxError.hpp"
33#include "../../mxException.hpp"
34
35#include "aoConstants.hpp"
36
37#include "aoAtmosphere.hpp"
38#include "aoPSDs.hpp"
39#include "aoWFS.hpp"
40
41#define FITTING_ERROR_NO 0
42#define FITTING_ERROR_ZERO 1
43#define FITTING_ERROR_X 2
44#define FITTING_ERROR_Y 3
45
46namespace mx
47{
48namespace AO
49{
50namespace analysis
51{
52
53/// Describes an analytic adaptive optics (AO) system.
54/**
55 * Templatized by the turbulence power spectral density (PSD).
56 *
57 * \tparam realT the floating point type used for all calculations
58 * \tparam inputSpecT specifies the turbulence spatial PSD type
59 * \tparam iosT is an output stream type with operator << defined (default is std::ostream)
60 *
61 * \ingroup mxAOAnalytic
62 */
63template <typename _realT, class _inputSpectT, typename iosT = std::ostream>
65{
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
90 1 }; ///< The optical gain of the WFS, 0-1. Treated as a sensitivity 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
104 0 }; ///< The optimum binning factor. If WFS modes are used, this is the mode index (0 to N-1). If not, it is
105 ///< 1 minus the pixel binning factor. 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 mxThrowException( err::sizeerr, "aoSystem::d_min", "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 mxThrowException( err::sizeerr, "aoSystem::d_min", "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 mxThrowException( err::paramnotset, "aoSystem::beta_p", "The WFS is not assigned." );
1380
1381 return m_wfsBeta->beta_p( m, n, m_D, d_opt(), atm.r_0( m_lam_sci ) );
1382}
1383
1384template <typename realT, class inputSpectT, typename iosT>
1386{
1387 if( m_wfsBeta == 0 )
1388 mxThrowException( err::paramnotset, "aoSystem::beta_r", "The WFS is not assigned." );
1389
1390 return m_wfsBeta->beta_r( m, n, m_D, d_opt(), atm.r_0( m_lam_sci ) );
1391}
1392
1393template <typename realT, class inputSpectT, typename iosT>
1395{
1396 m_opticalGain = og;
1397 m_specsChanged = true;
1398 m_dminChanged = true;
1399}
1400
1401template <typename realT, class inputSpectT, typename iosT>
1403{
1404 return m_opticalGain;
1405}
1406
1407template <typename realT, class inputSpectT, typename iosT>
1409{
1410 m_lam_wfs = nlam;
1411 m_specsChanged = true;
1412 m_dminChanged = true;
1413}
1414
1415template <typename realT, class inputSpectT, typename iosT>
1417{
1418 return m_lam_wfs;
1419}
1420
1421template <typename realT, class inputSpectT, typename iosT>
1422void aoSystem<realT, inputSpectT, iosT>::npix_wfs( const std::vector<realT> &npix )
1423{
1424 m_npix_wfs.resize( npix.size() );
1425 for( size_t n = 0; n < npix.size(); ++n )
1426 {
1427 m_npix_wfs[n] = npix[n];
1428 }
1429
1430 m_specsChanged = true;
1431 m_dminChanged = true;
1432}
1433
1434template <typename realT, class inputSpectT, typename iosT>
1436{
1437 if( m_npix_wfs.size() < idx + 1 )
1438 {
1439 mxThrowException( err::sizeerr, "aoSystem::npix_wfs", "idx larger than m_npix_wfs" );
1440 }
1441
1442 m_npix_wfs[idx] = npix;
1443
1444 m_specsChanged = true;
1445 m_dminChanged = true;
1446}
1447
1448template <typename realT, class inputSpectT, typename iosT>
1450{
1451 if( m_npix_wfs.size() < idx + 1 )
1452 {
1453 mxThrowException( err::sizeerr, "aoSystem::npix_wfs", "idx larger than m_npix_wfs" );
1454 }
1455
1456 return m_npix_wfs[idx];
1457}
1458
1459template <typename realT, class inputSpectT, typename iosT>
1461{
1462 return m_npix_wfs;
1463}
1464
1465template <typename realT, class inputSpectT, typename iosT>
1466void aoSystem<realT, inputSpectT, iosT>::ron_wfs( const std::vector<realT> &nron )
1467{
1468 m_ron_wfs.resize( nron.size() );
1469 for( size_t n = 0; n < nron.size(); ++n )
1470 {
1471 m_ron_wfs[n] = nron[n];
1472 }
1473
1474 m_specsChanged = true;
1475 m_dminChanged = true;
1476}
1477
1478template <typename realT, class inputSpectT, typename iosT>
1480{
1481 if( m_ron_wfs.size() < idx + 1 )
1482 {
1483 mxThrowException( err::sizeerr, "aoSystem::ron_wfs", "idx larger than m_ron_wfs" );
1484 }
1485
1486 m_ron_wfs[idx] = nron;
1487
1488 m_specsChanged = true;
1489 m_dminChanged = true;
1490}
1491
1492template <typename realT, class inputSpectT, typename iosT>
1494{
1495 if( m_ron_wfs.size() < idx + 1 )
1496 {
1497 mxThrowException( err::sizeerr, "aoSystem::ron_wfs", "idx larger than m_ron_wfs" );
1498 }
1499
1500 return m_ron_wfs[idx];
1501}
1502
1503template <typename realT, class inputSpectT, typename iosT>
1505{
1506 return m_ron_wfs;
1507}
1508
1509template <typename realT, class inputSpectT, typename iosT>
1510void aoSystem<realT, inputSpectT, iosT>::Fbg( const std::vector<realT> &fbg )
1511{
1512 m_Fbg.resize( fbg.size() );
1513 for( size_t n = 0; n < fbg.size(); ++n )
1514 {
1515 m_Fbg[n] = fbg[n];
1516 }
1517
1518 m_specsChanged = true;
1519 m_dminChanged = true;
1520}
1521
1522template <typename realT, class inputSpectT, typename iosT>
1524{
1525 if( m_Fbg.size() < idx + 1 )
1526 {
1527 mxThrowException( err::sizeerr, "aoSyste,::Fbg", "idx larger than m_Fbg" );
1528 }
1529
1530 m_Fbg[idx] = fbg;
1531
1532 m_specsChanged = true;
1533 m_dminChanged = true;
1534}
1535
1536template <typename realT, class inputSpectT, typename iosT>
1538{
1539 if( m_Fbg.size() < idx + 1 )
1540 {
1541 mxThrowException( err::sizeerr, "aoSyste,::Fbg", "idx larger than m_Fbg" );
1542 }
1543
1544 return m_Fbg[idx];
1545}
1546
1547template <typename realT, class inputSpectT, typename iosT>
1549{
1550 return m_Fbg;
1551}
1552
1553template <typename realT, class inputSpectT, typename iosT>
1554void aoSystem<realT, inputSpectT, iosT>::minTauWFS( const std::vector<realT> &ntau )
1555{
1556 m_minTauWFS.resize( ntau.size() );
1557 for( size_t n = 0; n < ntau.size(); ++n )
1558 {
1559 m_minTauWFS[n] = ntau[n];
1560 }
1561
1562 m_specsChanged = true;
1563 m_dminChanged = true;
1564}
1565
1566template <typename realT, class inputSpectT, typename iosT>
1568{
1569 if( m_minTauWFS.size() < idx + 1 )
1570 {
1571 mxThrowException( err::sizeerr, "aoSystem::minTauWFS", "idx larger than m_ntau_wfs" );
1572 }
1573
1574 m_minTauWFS[idx] = ntau;
1575
1576 m_specsChanged = true;
1577 m_dminChanged = true;
1578}
1579
1580template <typename realT, class inputSpectT, typename iosT>
1582{
1583 if( m_minTauWFS.size() < idx + 1 )
1584 {
1585 mxThrowException( err::sizeerr, "aoSystem::minTauWFS", "idx larger than m_ntau_wfs" );
1586 }
1587
1588 return m_minTauWFS[idx];
1589}
1590
1591template <typename realT, class inputSpectT, typename iosT>
1593{
1594 return m_minTauWFS;
1595}
1596
1597template <typename realT, class inputSpectT, typename iosT>
1599{
1600 if( bnp != m_bin_npix )
1601 {
1602 m_bin_npix = bnp;
1603 m_specsChanged = true;
1604 m_dminChanged = true;
1605 }
1606}
1607
1608template <typename realT, class inputSpectT, typename iosT>
1610{
1611 return m_bin_npix;
1612}
1613
1614template <typename realT, class inputSpectT, typename iosT>
1616{
1617 d_opt();
1618 return m_bin_opt;
1619}
1620
1621template <typename realT, class inputSpectT, typename iosT>
1623{
1624 m_tauWFS = ntau;
1625 m_specsChanged = true;
1626 m_dminChanged = true;
1627}
1628
1629template <typename realT, class inputSpectT, typename iosT>
1631{
1632 return m_tauWFS;
1633}
1634
1635template <typename realT, class inputSpectT, typename iosT>
1637{
1638 m_deltaTau = ndel;
1639 m_specsChanged = true;
1640 m_dminChanged = true;
1641}
1642
1643template <typename realT, class inputSpectT, typename iosT>
1645{
1646 return m_deltaTau;
1647}
1648
1649template <typename realT, class inputSpectT, typename iosT>
1651{
1652 m_optTau = ot;
1653 m_specsChanged = true;
1654 m_dminChanged = true;
1655}
1656
1657template <typename realT, class inputSpectT, typename iosT>
1659{
1660 return m_optTau;
1661}
1662
1663template <typename realT, class inputSpectT, typename iosT>
1665{
1666 m_lam_sci = nlam;
1667 m_specsChanged = true;
1668 m_dminChanged = true;
1669}
1670
1671template <typename realT, class inputSpectT, typename iosT>
1673{
1674 return m_lam_sci;
1675}
1676
1677template <typename realT, class inputSpectT, typename iosT>
1679{
1680 m_zeta = nz;
1681 m_secZeta = 1 / cos( m_zeta );
1682
1683 m_specsChanged = true;
1684 m_dminChanged = true;
1685}
1686
1687template <typename realT, class inputSpectT, typename iosT>
1689{
1690 return m_zeta;
1691}
1692
1693template <typename realT, class inputSpectT, typename iosT>
1695{
1696 return m_secZeta;
1697}
1698
1699template <typename realT, class inputSpectT, typename iosT>
1701{
1702 if( mnm < 0 )
1703 mnm = 0;
1704 m_fit_mn_max = mnm;
1705}
1706
1707template <typename realT, class inputSpectT, typename iosT>
1709{
1710 return m_fit_mn_max;
1711}
1712
1713template <typename realT, class inputSpectT, typename iosT>
1715{
1716 m_circularLimit = cl;
1717 m_specsChanged = true;
1718 m_dminChanged = true;
1719}
1720
1721template <typename realT, class inputSpectT, typename iosT>
1723{
1724 return m_circularLimit;
1725}
1726
1727template <typename realT, class inputSpectT, typename iosT>
1729{
1730 m_spatialFilter_ku = fabs( kx );
1731 m_specsChanged = true; // not sure if needed
1732 m_dminChanged = true; // not sure if needed
1733}
1734
1735template <typename realT, class inputSpectT, typename iosT>
1737{
1738 return m_spatialFilter_ku;
1739}
1740
1741template <typename realT, class inputSpectT, typename iosT>
1743{
1744 m_spatialFilter_kv = fabs( ky );
1745 m_specsChanged = true; // not sure if needed
1746 m_dminChanged = true; // not sure if needed
1747}
1748
1749template <typename realT, class inputSpectT, typename iosT>
1751{
1752 return m_spatialFilter_kv;
1753}
1754
1755template <typename realT, class inputSpectT, typename iosT>
1757{
1758 m_ncp_wfe = nwfe;
1759 m_specsChanged = true;
1760 m_dminChanged = true;
1761}
1762
1763template <typename realT, class inputSpectT, typename iosT>
1765{
1766 return m_ncp_wfe;
1767}
1768
1769template <typename realT, class inputSpectT, typename iosT>
1771{
1772 m_ncp_alpha = alpha;
1773 m_specsChanged = true;
1774 m_dminChanged = true;
1775}
1776
1777template <typename realT, class inputSpectT, typename iosT>
1779{
1780 return m_ncp_alpha;
1781}
1782
1783template <typename realT, class inputSpectT, typename iosT>
1785{
1786 m_F0 = nF0;
1787 m_specsChanged = true;
1788 m_dminChanged = true;
1789}
1790
1791template <typename realT, class inputSpectT, typename iosT>
1793{
1794 return m_F0;
1795}
1796
1797template <typename realT, class inputSpectT, typename iosT>
1799{
1800 m_starMag = nmag;
1801 m_specsChanged = true;
1802 m_dminChanged = true;
1803}
1804
1805template <typename realT, class inputSpectT, typename iosT>
1807{
1808 return m_starMag;
1809}
1810
1811template <typename realT, class inputSpectT, typename iosT>
1813{
1814 return m_F0 * pow( 10.0, -0.4 * mag );
1815}
1816
1817template <typename realT, class inputSpectT, typename iosT>
1819{
1820 return Fg( m_starMag );
1821}
1822
1823template <typename realT, class inputSpectT, typename iosT>
1825{
1826 if( m_specsChanged || m_dminChanged )
1827 calcStrehl();
1828
1829 return m_controlledModes;
1830}
1831
1832template <typename realT, class inputSpectT, typename iosT>
1833realT aoSystem<realT, inputSpectT, iosT>::signal2Noise2( realT &Nph, realT tau_wfs, realT d, int b )
1834{
1835 Nph = Fg() * tau_wfs;
1836
1837 double binfact = 1.0;
1838 int binidx = 0;
1839 if( m_bin_npix )
1840 {
1841 if( m_npix_wfs.size() == 1 ) // Only if "true binning", otherwise the WFS mode vectors handle it.
1842 {
1843 binfact = 1. / pow( (realT)b + 1, 2 );
1844 }
1845 else
1846 binidx = b;
1847 }
1848
1849 return pow( Nph, 2 ) / ( m_npix_wfs[binidx] * binfact * ( m_Fbg[binidx] * tau_wfs + pow( m_ron_wfs[binidx], 2 ) ) );
1850}
1851
1852template <typename realT, class inputSpectT, typename iosT>
1853realT aoSystem<realT, inputSpectT, iosT>::measurementError( realT m, realT n, realT d, int b )
1854{
1855 if( m == 0 and n == 0 )
1856 return 0;
1857
1858 realT tau_wfs;
1859
1860 if( m_optTau )
1861 tau_wfs = optimumTauWFS( m, n, d, b );
1862 else
1863 tau_wfs = m_tauWFS;
1864
1865 if( m_wfsBeta == 0 )
1866 mxThrowException( err::paramnotset, "aoSystem::measurementError", "The WFS is not assigned." );
1867
1868 realT beta_p = m_wfsBeta->beta_p( m, n, m_D, d, atm.r_0( m_lam_wfs ) ) / sqrt( m_opticalGain );
1869 realT beta_r = m_wfsBeta->beta_r( m, n, m_D, d, atm.r_0( m_lam_wfs ) ); // sqrt(m_opticalGain);
1870 realT Nph = 0;
1871 realT snr2 = signal2Noise2( Nph, tau_wfs, d, b );
1872
1873 return ( pow( beta_r, 2 ) / snr2 + pow( beta_p, 2 ) / Nph ) * pow( m_lam_wfs / m_lam_sci, 2 );
1874}
1875
1876template <typename realT, class inputSpectT, typename iosT>
1878{
1879 if( m_specsChanged || m_dminChanged )
1880 calcStrehl();
1881 return m_wfeMeasurement;
1882}
1883
1884template <typename realT, class inputSpectT, typename iosT>
1885realT aoSystem<realT, inputSpectT, iosT>::timeDelayError( realT m, realT n, realT d, int b )
1886{
1887 if( m == 0 and n == 0 )
1888 return 0;
1889
1890 realT k = sqrt( m * m + n * n ) / m_D;
1891
1892 realT tau_wfs;
1893
1894 if( m_optTau )
1895 tau_wfs = optimumTauWFS( m, n, d, b );
1896 else
1897 tau_wfs = m_tauWFS;
1898
1899 realT tau = tau_wfs + m_deltaTau;
1900
1901 ///\todo handle multiple layers
1902 return psd( atm, 0, k, m_secZeta ) / pow( m_D, 2 ) * pow( atm.lam_0() / m_lam_sci, 2 ) *
1903 sqrt( atm.X( k, m_lam_sci, m_secZeta ) ) * pow( math::two_pi<realT>() * atm.v_wind() * k, 2 ) *
1904 pow( tau, 2 );
1905}
1906
1907template <typename realT, class inputSpectT, typename iosT>
1909{
1910 if( m_specsChanged || m_dminChanged )
1911 calcStrehl();
1912 return m_wfeTimeDelay;
1913}
1914
1915template <typename realT, class inputSpectT, typename iosT>
1917{
1918 realT k = sqrt( m * m + n * n ) / m_D;
1919
1920 ///\todo handle multiple layers
1921 return psd( atm, 0, k, m_secZeta ) / pow( m_D, 2 ) * pow( atm.lam_0() / m_lam_sci, 2 );
1922}
1923
1924template <typename realT, class inputSpectT, typename iosT>
1926{
1927 if( m_specsChanged || m_dminChanged )
1928 calcStrehl();
1929 return m_wfeFitting;
1930}
1931
1932template <typename realT, class inputSpectT, typename iosT>
1934{
1935 return C4var( m, n );
1936}
1937
1938template <typename realT, class inputSpectT, typename iosT>
1940{
1941 if( m_specsChanged || m_dminChanged )
1942 calcStrehl();
1943 return m_wfeChromScintOPD;
1944}
1945
1946template <typename realT, class inputSpectT, typename iosT>
1948{
1949 return C5var( m, n );
1950
1951#if 0
1952 int mn_max = floor(0.5*m_D/d_opt());
1953
1954 realT sum = 0;
1955
1956 for(int m = -mn_max; m <= mn_max; ++m)
1957 {
1958 for(int n = -mn_max; n <= mn_max; ++n)
1959 {
1960 if(n == 0 && m == 0) continue;
1961
1962 sum += C5var(m,n);
1963 }
1964 }
1965
1966 return sum;
1967#endif
1968}
1969
1970template <typename realT, class inputSpectT, typename iosT>
1972{
1973 return C6var( m, n );
1974}
1975
1976template <typename realT, class inputSpectT, typename iosT>
1978{
1979 if( m_specsChanged || m_dminChanged )
1980 calcStrehl();
1981 return m_wfeChromIndex;
1982}
1983
1984template <typename realT, class inputSpectT, typename iosT>
1986{
1987 return C7var( m, n );
1988}
1989
1990template <typename realT, class inputSpectT, typename iosT>
1992{
1993 if( m_specsChanged || m_dminChanged )
1994 calcStrehl();
1995 return m_wfeAnisoOPD;
1996}
1997
1998template <typename realT, class inputSpectT, typename iosT>
2000{
2001 return 0;
2002
2003#if 0
2004 int mn_max = floor(0.5*m_D/d_opt());
2005
2006 realT sum = 0;
2007
2008 for(int m = -mn_max; m <= mn_max; ++m)
2009 {
2010 for(int n = -mn_max; n <= mn_max; ++n)
2011 {
2012 if(n == 0 && m == 0) continue;
2013
2014 if( m_circularLimit )
2015 {
2016 if( m*m + n*n > mn_max*mn_max ) continue;
2017 }
2018
2019 sum += C8var(m,n);
2020 }
2021 }
2022
2023 return sum;
2024#endif
2025}
2026
2027template <typename realT, class inputSpectT, typename iosT>
2029 realT m,
2030 realT n,
2031 realT dact, // here called dact due to parameter collision with root-finding
2032 int bbin )
2033{
2034 if( m_D == 0 )
2035 {
2036 mxError( "aoSystem::optimumTauWFS", MXE_PARAMNOTSET, "Diameter (D) not set." );
2037 return -1;
2038 }
2039
2040 if( m_F0 == 0 )
2041 {
2042 mxError( "aoSystem::optimumTauWFS", MXE_PARAMNOTSET, "0-mag photon flux (F0) not set." );
2043 return -1;
2044 }
2045
2046 double binfact = 1.0;
2047 int binidx = 0; // index into WFS configurations
2048 if( m_bin_npix )
2049 {
2050 if( m_npix_wfs.size() > 1 )
2051 {
2052 binidx = bbin;
2053 }
2054 else
2055 {
2056 binfact = 1.0 / pow( (realT)( bbin + 1 ), 2 );
2057 }
2058 }
2059
2060 realT k = sqrt( m * m + n * n ) / m_D;
2061
2062 realT F = Fg();
2063
2064 if( m_wfsBeta == 0 )
2065 mxThrowException( err::paramnotset, "aoSystem::beta_p", "The WFS is not assigned." );
2066
2067 realT beta_p = m_wfsBeta->beta_p( m, n, m_D, dact, atm.r_0( m_lam_wfs ) ) / sqrt( m_opticalGain );
2068 realT beta_r = m_wfsBeta->beta_r( m, n, m_D, dact, atm.r_0( m_lam_wfs ) ) / sqrt( m_opticalGain );
2069
2070 // Set up for root finding:
2071 realT a, b, c, d, e;
2072
2073 ///\todo handle multiple layers
2074 realT Atmp = 2 * pow( atm.lam_0(), 2 ) * psd( atm, 0, k, m_secZeta ) / pow( m_D, 2 ) *
2075 ( atm.X( k, m_lam_wfs, m_secZeta ) ) * pow( math::two_pi<realT>() * atm.v_wind() * k, 2 );
2076
2077 a = Atmp;
2078 b = Atmp * m_deltaTau;
2079 c = 0;
2080 d = -1 * ( pow( m_lam_wfs, 2 ) / F ) *
2081 ( ( m_npix_wfs[binidx] * binfact * m_Fbg[binidx] / F ) * pow( beta_r, 2 ) + pow( beta_p, 2 ) );
2082 e = -2 * pow( m_lam_wfs, 2 ) * ( m_npix_wfs[binidx] * binfact ) * pow( m_ron_wfs[binidx], 2 ) * pow( beta_r, 2 ) /
2083 pow( F, 2 );
2084
2085 std::vector<std::complex<realT>> x;
2086
2087 // Get the roots
2088 math::quarticRoots( x, a, b, c, d, e );
2089
2090 // Now pick the largest positive real root
2091 realT tauopt = 0.0;
2092
2093 for( int i = 0; i < 4; i++ )
2094 {
2095 if( real( x[i] ) > 0 && imag( x[i] ) == 0 && real( x[i] ) > tauopt )
2096 tauopt = real( x[i] );
2097 }
2098
2099 if( tauopt < m_minTauWFS[binidx] )
2100 tauopt = m_minTauWFS[binidx];
2101
2102 return tauopt;
2103}
2104
2105template <typename realT, class inputSpectT, typename iosT>
2107{
2108 d_opt(); // also gets m_bin_opt;
2109 return optimumTauWFS( m, n, m_d_opt, m_bin_opt );
2110}
2111
2112template <typename realT, class inputSpectT, typename iosT>
2114{
2115 if( !m_dminChanged )
2116 return m_d_opt;
2117
2118 // Validate the modes
2119 ///\todo this should be in a separate valideModes function
2120 if( m_npix_wfs.size() < 1 )
2121 {
2122 mxThrowException( err::sizeerr, "aoSystem::d_opt", "npix_wfs must have at least one entry" );
2123 }
2124
2125 if( m_d_min.size() != m_npix_wfs.size() )
2126 {
2127 mxThrowException( err::sizeerr, "aoSystem::d_opt", "d_min must be the same size as npix_wfs" );
2128 }
2129
2130 if( m_ron_wfs.size() != m_npix_wfs.size() )
2131 {
2132 mxThrowException( err::sizeerr, "aoSystem::d_opt", "ron_wfs must be the same size as npix_wfs" );
2133 }
2134
2135 if( m_Fbg.size() != m_npix_wfs.size() )
2136 {
2137 mxThrowException( err::sizeerr, "aoSystem::d_opt", "F_bg must be the same size as npix_wfs" );
2138 }
2139
2140 if( m_minTauWFS.size() != m_npix_wfs.size() )
2141 {
2142 mxThrowException( err::sizeerr, "aoSystem::d_opt", "minTauWFS must be the same size as npix_wfs" );
2143 }
2144
2145 ///\todo investigate why we tolerate m_d_min = 0 here. Is this just a config error?
2146 if( m_d_min.size() == 1 && m_d_min[0] == 0 )
2147 {
2148 m_d_opt = 1e-50;
2149 m_bin_opt = 0;
2150 m_dminChanged = false;
2151 return m_d_opt;
2152 }
2153
2154 if( !m_optd )
2155 {
2156 m_d_opt = m_d_min[0];
2157 m_bin_opt = 0;
2158 m_dminChanged = false;
2159
2160 return m_d_opt;
2161 }
2162
2163 realT best_d = 0;
2164
2165 if( m_bin_npix )
2166 {
2167 if( m_npix_wfs.size() > 1 ) // Optimize over the WFS modes
2168 {
2169 int best_idx = 0;
2170 best_d = m_d_min[0];
2171
2172 realT bestStrehlOverall = 0;
2173
2174 for( int b = 0; b < m_npix_wfs.size(); ++b )
2175 {
2176 realT d = m_d_min[b];
2177 realT s = strehl( d, b );
2178 realT bestStrehl = s;
2179
2180 // Find the actuator pitch at which AO error is less than fitting error at Nyquist
2181 while( d < m_D / 2 )
2182 {
2183 d += m_d_min[b] / m_optd_delta;
2184 s = strehl( d, b );
2185
2186 if( s < bestStrehl ) // Strehl got worse. Can't be <= otherwise small m_optd_deltas will
2187 // immediately out
2188 {
2189 d -= m_d_min[b] / m_optd_delta; // go back to previous d
2190 break;
2191 }
2192 else
2193 bestStrehl = s;
2194 }
2195
2196 // Check if this is optimum so far
2197 if( bestStrehl >= bestStrehlOverall ) // >= to favor fewer actuators
2198 {
2199 bestStrehlOverall = bestStrehl;
2200 best_idx = b;
2201 best_d = d;
2202 }
2203 }
2204
2205 m_bin_opt = best_idx;
2206 }
2207 else
2208 {
2209 realT d = m_d_min[0];
2210 m_bin_opt = 0;
2211
2212 realT s = strehl( d, m_bin_opt );
2213 realT bestStrehl = s;
2214
2215 while( d < m_D / 2 )
2216 {
2217 d += m_d_min[0] / m_optd_delta;
2218 m_bin_opt = d / m_d_min[0] - 1;
2219 if( m_bin_opt < 0 )
2220 m_bin_opt = 0;
2221
2222 s = strehl( d, m_bin_opt );
2223
2224 if( s < bestStrehl ) // Strehl got worse. Can't be <= otherwise small m_optd_deltas will immediately
2225 // out
2226 {
2227 d -= m_d_min[0] / m_optd_delta; // go back to previous d
2228 m_bin_opt = d / m_d_min[0] - 1;
2229 if( m_bin_opt < 0 )
2230 m_bin_opt = 0;
2231 break;
2232 }
2233 else
2234 bestStrehl = s;
2235 }
2236 best_d = d;
2237 }
2238 }
2239 else
2240 {
2241 realT d = m_d_min[0];
2242 m_bin_opt = 0;
2243
2244 realT s = strehl( d, m_bin_opt );
2245 realT bestStrehl = s;
2246
2247 // Find the actuator pitch which minimizes total error (AO error is less than fitting error at Nyquist)
2248 while( d < m_D / 2 )
2249 {
2250 d += m_d_min[0] / m_optd_delta;
2251
2252 s = strehl( d, m_bin_opt );
2253 if( s < bestStrehl ) // Strehl got worse. Can't be <= otherwise small m_optd_deltas will immediately out
2254 {
2255 d -= m_d_min[0] / m_optd_delta; // go back to previous d
2256 break;
2257 }
2258 else
2259 bestStrehl = s;
2260 }
2261 best_d = d;
2262 }
2263
2264 m_d_opt = best_d;
2265 m_dminChanged = false;
2266
2267 return m_d_opt;
2268}
2269
2270template <typename realT, class inputSpectT, typename iosT>
2272{
2273 if( m == 0 and n == 0 )
2274 return 0;
2275
2276 realT k = sqrt( m * m + n * n ) / m_D;
2277
2278 realT kmax = m_fit_mn_max / m_D;
2279 realT kmin = 1. / m_D;
2280
2281 realT beta;
2282 if( m_ncp_alpha != 2 )
2283 {
2284 beta = ( m_ncp_alpha - 2 ) / ( math::two_pi<realT>() ) * 1. /
2285 ( pow( kmin, -m_ncp_alpha + 2 ) - pow( kmax, -m_ncp_alpha + 2 ) );
2286 }
2287 else
2288 {
2289 beta = 1. / ( math::two_pi<realT>() * log( kmax / kmin ) );
2290 }
2291 return beta * ncpError() * pow( k, -m_ncp_alpha ) * pow( kmin, 2 );
2292}
2293
2294template <typename realT, class inputSpectT, typename iosT>
2296{
2297 return pow( m_ncp_wfe, 2 ) * pow( math::two_pi<realT>() / m_lam_sci, 2 );
2298}
2299
2300template <typename realT, class inputSpectT, typename iosT>
2302{
2303
2304 int mn_max = floor( 0.5 * m_D / d );
2305
2306 realT wfeVar = 0;
2307
2308 bool controlled;
2309 for( int m = -m_fit_mn_max; m <= m_fit_mn_max; ++m )
2310 {
2311 for( int n = -m_fit_mn_max; n <= m_fit_mn_max; ++n )
2312 {
2313 controlled = false;
2314 if( m_circularLimit )
2315 {
2316 if( m * m + n * n <= mn_max * mn_max )
2317 controlled = true;
2318 }
2319 else
2320 {
2321 if( abs( m ) <= mn_max && abs( n ) <= mn_max )
2322 controlled = true;
2323 }
2324
2325 if( controlled )
2326 {
2327 realT wfeMeasurement = measurementError( m, n, d, b );
2328 realT wfeTimeDelay = timeDelayError( m, n, d, b );
2329 realT wfeChromScintOPD = chromScintOPDError( m, n );
2330 realT wfeChromIndex = chromIndexError( m, n );
2331 realT wfeAnisoOPD = dispAnisoOPDError( m, n );
2332
2333 realT wfeFitting = fittingError( m, n );
2334
2335 if( wfeFitting < wfeMeasurement + wfeTimeDelay + wfeChromScintOPD + wfeChromIndex + wfeAnisoOPD )
2336 {
2337 wfeVar += wfeFitting;
2338 }
2339 else
2340 {
2341 wfeVar += wfeMeasurement + wfeTimeDelay + wfeChromScintOPD + wfeChromIndex + wfeAnisoOPD;
2342 }
2343 }
2344 else
2345 {
2346 wfeVar += fittingError( m, n );
2347 }
2348 }
2349 }
2350
2351 wfeVar += ncpError();
2352 return exp( -1 * wfeVar );
2353}
2354
2355template <typename realT, class inputSpectT, typename iosT>
2357{
2358 realT d = d_opt(); // have to run this to get m_bin_opt set.
2359 int b = m_bin_opt;
2360
2361 int mn_max = floor( 0.5 * m_D / d );
2362
2363 m_wfeMeasurement = 0;
2364 m_wfeTimeDelay = 0;
2365 m_wfeChromScintOPD = 0;
2366 m_wfeChromIndex = 0;
2367 m_wfeAnisoOPD = 0;
2368
2369 m_wfeFitting = 0;
2370
2371 m_wfeNCP = 0;
2372
2373 m_controlledModes.resize( 2 * m_fit_mn_max + 1, 2 * m_fit_mn_max + 1 );
2374 m_controlledModes.setZero();
2375
2376 bool controlled;
2377 for( int m = -m_fit_mn_max; m <= m_fit_mn_max; ++m )
2378 {
2379 for( int n = -m_fit_mn_max; n <= m_fit_mn_max; ++n )
2380 {
2381 if( m == 0 && n == 0 )
2382 {
2383 continue;
2384 }
2385
2386 controlled = false;
2387 if( m_circularLimit )
2388 {
2389 if( m * m + n * n <= mn_max * mn_max )
2390 controlled = true;
2391 }
2392 else
2393 {
2394 if( abs( m ) <= mn_max && abs( n ) <= mn_max )
2395 controlled = true;
2396 }
2397
2398 if( controlled )
2399 {
2400 realT wfeMeasurement = measurementError( m, n, d, b );
2401 realT wfeTimeDelay = timeDelayError( m, n, d, b );
2402 realT wfeChromScintOPD = chromScintOPDError( m, n );
2403 realT wfeChromIndex = chromIndexError( m, n );
2404 realT wfeAnisoOPD = dispAnisoOPDError( m, n );
2405
2406 realT wfeFitting = fittingError( m, n );
2407
2408 if( wfeFitting < wfeMeasurement + wfeTimeDelay + wfeChromScintOPD + wfeChromIndex + wfeAnisoOPD )
2409 {
2410 m_wfeFitting += wfeFitting;
2411 }
2412 else // it's worth controlling this mode.
2413 {
2414 m_wfeMeasurement += wfeMeasurement;
2415 m_wfeTimeDelay += wfeTimeDelay;
2416 m_wfeChromScintOPD += wfeChromScintOPD;
2417 m_wfeChromIndex += wfeChromIndex;
2418 m_wfeAnisoOPD += wfeAnisoOPD;
2419 m_controlledModes( m_fit_mn_max + m, m_fit_mn_max + n ) = 1;
2420 }
2421 }
2422 else
2423 {
2424 m_wfeFitting += fittingError( m, n );
2425 }
2426 }
2427 }
2428
2429 m_wfeNCP = ncpError();
2430
2431 m_wfeVar = m_wfeMeasurement + m_wfeTimeDelay + m_wfeFitting + m_wfeChromScintOPD + m_wfeChromIndex + m_wfeAnisoOPD +
2432 m_wfeNCP;
2433
2434 m_strehl = exp( -1 * m_wfeVar );
2435
2436 m_specsChanged = false;
2437}
2438
2439template <typename realT, class inputSpectT, typename iosT>
2441{
2442 if( m_specsChanged || m_dminChanged )
2443 calcStrehl();
2444
2445 return m_wfeVar;
2446}
2447
2448template <typename realT, class inputSpectT, typename iosT>
2450{
2451 if( m_specsChanged || m_dminChanged )
2452 calcStrehl();
2453
2454 return m_strehl;
2455}
2456
2457template <typename realT, class inputSpectT, typename iosT>
2458template <typename varFuncT>
2459realT aoSystem<realT, inputSpectT, iosT>::C_( int m, int n, bool normStrehl, varFuncT varFunc, int doFittingError )
2460{
2461 if( m == 0 && n == 0 )
2462 return 0;
2463
2464 // we run this to optimize regardless of whether we use it
2465 // if already optimized then this is a minor hit
2466 realT S = strehl();
2467
2468 if( !normStrehl )
2469 S = 1;
2470
2471 if( doFittingError != FITTING_ERROR_NO )
2472 {
2473 int mn_max = m_D / ( 2. * d_opt() );
2474
2475 if( m_controlledModes( m_fit_mn_max + m, m_fit_mn_max + n ) == false )
2476 {
2477 if( doFittingError == FITTING_ERROR_ZERO )
2478 return 0;
2479
2480 realT fe = fittingError( m, n );
2481
2482 realT k = sqrt( m * m + n * n ) / m_D;
2483
2484 if( doFittingError == FITTING_ERROR_X )
2485 fe *= ( atm.X( k, m_lam_sci, m_secZeta ) );
2486 else if( doFittingError == FITTING_ERROR_Y )
2487 fe *= ( atm.Y( k, m_lam_sci, m_secZeta ) );
2488 else
2489 {
2490 std::cerr << "Unknown doFittingError\n";
2491 exit( -1 );
2492 }
2493
2494 return fe / S;
2495 }
2496 }
2497 // get if not doing fitting error or if inside control region:
2498
2499 realT var = ( this->*varFunc )( m, n );
2500
2501 return var / S;
2502}
2503
2504template <typename realT, class inputSpectT, typename iosT>
2505template <typename imageT, typename CfuncT>
2506void aoSystem<realT, inputSpectT, iosT>::C_Map( imageT &im, CfuncT Cfunc, bool normStrehl )
2507{
2508 int dim1 = im.rows();
2509 int dim2 = im.cols();
2510
2511 int mc = 0.5 * ( dim1 - 1 );
2512 int nc = 0.5 * ( dim2 - 1 );
2513
2514 for( int i = 0; i < dim1; ++i )
2515 {
2516 int m = i - mc;
2517
2518 for( int j = 0; j < dim2; ++j )
2519 {
2520 int n = j - nc;
2521
2522 im( i, j ) = ( this->*Cfunc )( m, n, normStrehl );
2523 }
2524 }
2525}
2526
2527template <typename realT, class inputSpectT, typename iosT>
2529{
2530 realT k = sqrt( m * m + n * n ) / m_D;
2531 ///\todo handle multiple layers
2532 return psd( atm, 0, k, m_secZeta ) / pow( m_D, 2 ) * pow( atm.lam_0() / m_lam_sci, 2 ) *
2533 ( atm.X( k, m_lam_sci, m_secZeta ) );
2534}
2535
2536template <typename realT, class inputSpectT, typename iosT>
2537realT aoSystem<realT, inputSpectT, iosT>::C0( realT m, realT n, bool normStrehl )
2538{
2539
2540 return C_( m, n, normStrehl, &aoSystem<realT, inputSpectT, iosT>::C0var, FITTING_ERROR_NO );
2541}
2542
2543template <typename realT, class inputSpectT, typename iosT>
2544template <typename imageT>
2545void aoSystem<realT, inputSpectT, iosT>::C0Map( imageT &im, bool normStrehl )
2546{
2547 C_Map( im, &aoSystem<realT, inputSpectT, iosT>::C0, normStrehl );
2548}
2549
2550template <typename realT, class inputSpectT, typename iosT>
2552{
2553 realT k = sqrt( m * m + n * n ) / m_D;
2554
2555 ///\todo handle multiple layers
2556 return psd( atm, 0, k, m_secZeta ) / pow( m_D, 2 ) * pow( atm.lam_0() / m_lam_sci, 2 ) *
2557 ( atm.Y( k, m_lam_sci, m_secZeta ) );
2558}
2559
2560template <typename realT, class inputSpectT, typename iosT>
2561realT aoSystem<realT, inputSpectT, iosT>::C1( realT m, realT n, bool normStrehl )
2562{
2563 return C_( m, n, normStrehl, &aoSystem<realT, inputSpectT, iosT>::C1var, FITTING_ERROR_NO );
2564}
2565
2566template <typename realT, class inputSpectT, typename iosT>
2567template <typename imageT>
2568void aoSystem<realT, inputSpectT, iosT>::C1Map( imageT &im, bool normStrehl )
2569{
2570 C_Map( im, &aoSystem<realT, inputSpectT, iosT>::C1, normStrehl );
2571}
2572
2573template <typename realT, class inputSpectT, typename iosT>
2575{
2576 d_opt();
2577 return measurementError( m, n, d_opt(), m_bin_opt ) + timeDelayError( m, n, d_opt(), m_bin_opt );
2578}
2579
2580template <typename realT, class inputSpectT, typename iosT>
2581realT aoSystem<realT, inputSpectT, iosT>::C2( realT m, realT n, bool normStrehl )
2582{
2583 return C_( m, n, normStrehl, &aoSystem<realT, inputSpectT, iosT>::C2var, FITTING_ERROR_X );
2584}
2585
2586template <typename realT, class inputSpectT, typename iosT>
2587template <typename imageT>
2588void aoSystem<realT, inputSpectT, iosT>::C2Map( imageT &im, bool normStrehl )
2589{
2590 C_Map( im, &aoSystem<realT, inputSpectT, iosT>::C2, normStrehl );
2591}
2592
2593template <typename realT, class inputSpectT, typename iosT>
2595{
2596 return 0; // measurementError(m, n) + timeDelayError(m,n);
2597}
2598
2599template <typename realT, class inputSpectT, typename iosT>
2600realT aoSystem<realT, inputSpectT, iosT>::C3( realT m, realT n, bool normStrehl )
2601{
2602 return C_( m, n, normStrehl, &aoSystem<realT, inputSpectT, iosT>::C3var, FITTING_ERROR_ZERO ); // FITTING_ERROR_Y);
2603}
2604
2605template <typename realT, class inputSpectT, typename iosT>
2606template <typename imageT>
2607void aoSystem<realT, inputSpectT, iosT>::C3Map( imageT &im, bool normStrehl )
2608{
2609 C_Map( im, &aoSystem<realT, inputSpectT, iosT>::C3, normStrehl );
2610}
2611
2612template <typename realT, class inputSpectT, typename iosT>
2614{
2615 realT k = sqrt( m * m + n * n ) / m_D;
2616
2617 ///\todo handle multiple layers
2618 // This does not need to be divided by X b/c we haven't multiplied by it, this is is C0/X.
2619 return psd( atm, 0, k, m_secZeta ) / pow( m_D, 2 ) * pow( atm.lam_0() / m_lam_sci, 2 ) *
2620 atm.dX( k, m_lam_sci, m_lam_wfs );
2621}
2622
2623template <typename realT, class inputSpectT, typename iosT>
2624realT aoSystem<realT, inputSpectT, iosT>::C4( realT m, realT n, bool normStrehl )
2625{
2626 return C_( m, n, normStrehl, &aoSystem<realT, inputSpectT, iosT>::C4var, FITTING_ERROR_ZERO );
2627}
2628
2629template <typename realT, class inputSpectT, typename iosT>
2630template <typename imageT>
2631void aoSystem<realT, inputSpectT, iosT>::C4Map( imageT &im, bool normStrehl )
2632{
2633 C_Map( im, &aoSystem<realT, inputSpectT, iosT>::C4, normStrehl );
2634}
2635
2636template <typename realT, class inputSpectT, typename iosT>
2638{
2639 realT k = sqrt( m * m + n * n ) / m_D;
2640
2641 ///\todo handle multiple layers
2642 // This does not need to be divided by Y b/c we haven't multiplied by it, this is is C1/Y.
2643 return psd( atm, 0, k, m_secZeta ) / pow( m_D, 2 ) * pow( atm.lam_0() / m_lam_sci, 2 ) *
2644 atm.dY( k, m_lam_sci, m_lam_wfs );
2645}
2646
2647template <typename realT, class inputSpectT, typename iosT>
2648realT aoSystem<realT, inputSpectT, iosT>::C5( realT m, realT n, bool normStrehl )
2649{
2650 return C_( m, n, normStrehl, &aoSystem<realT, inputSpectT, iosT>::C5var, FITTING_ERROR_ZERO );
2651}
2652
2653template <typename realT, class inputSpectT, typename iosT>
2654template <typename imageT>
2655void aoSystem<realT, inputSpectT, iosT>::C5Map( imageT &im, bool normStrehl )
2656{
2657 C_Map( im, &aoSystem<realT, inputSpectT, iosT>::C5, normStrehl );
2658}
2659
2660template <typename realT, class inputSpectT, typename iosT>
2662{
2663 realT ni = atm.n_air( m_lam_sci );
2664 realT nw = atm.n_air( m_lam_wfs );
2665
2666 return C0var( m, n ) * pow( ( ni - nw ) / ( ni - 1 ), 2 ); // Fixes error in Eqn 29
2667}
2668
2669template <typename realT, class inputSpectT, typename iosT>
2670realT aoSystem<realT, inputSpectT, iosT>::C6( realT m, realT n, bool normStrehl )
2671{
2672 return C_( m, n, normStrehl, &aoSystem<realT, inputSpectT, iosT>::C6var, FITTING_ERROR_ZERO );
2673}
2674
2675template <typename realT, class inputSpectT, typename iosT>
2676template <typename imageT>
2677void aoSystem<realT, inputSpectT, iosT>::C6Map( imageT &im, bool normStrehl )
2678{
2679 C_Map( im, &aoSystem<realT, inputSpectT, iosT>::C6, normStrehl );
2680}
2681
2682template <typename realT, class inputSpectT, typename iosT>
2684{
2685 realT k = sqrt( m * m + n * n ) / m_D;
2686
2687 ///\todo this needs to handle multiple layers -- do we violate an assumption of L_0 isn't constant?
2688 return psd( atm, 0, k, m_secZeta ) / pow( m_D, 2 ) * pow( atm.lam_0() / m_lam_sci, 2 ) *
2689 atm.X_Z( k, m_lam_wfs, m_lam_sci, m_secZeta );
2690}
2691
2692template <typename realT, class inputSpectT, typename iosT>
2693realT aoSystem<realT, inputSpectT, iosT>::C7( realT m, realT n, bool normStrehl )
2694{
2695 return C_( m, n, normStrehl, &aoSystem<realT, inputSpectT, iosT>::C7var, FITTING_ERROR_ZERO );
2696}
2697
2698template <typename realT, class inputSpectT, typename iosT>
2699template <typename imageT>
2700void aoSystem<realT, inputSpectT, iosT>::C7Map( imageT &im, bool normStrehl )
2701{
2702 C_Map( im, &aoSystem<realT, inputSpectT, iosT>::C7, normStrehl );
2703}
2704
2705template <typename realT, class inputSpectT, typename iosT>
2707{
2708 ios << "# AO Params:\n";
2709 ios << "# D = " << D() << '\n';
2710
2711 if( m_d_min.size() > 0 )
2712 {
2713 ios << "# d_min = " << d_min( (size_t)0 );
2714 for( size_t n = 1; n < m_d_min.size(); ++n )
2715 ios << ',' << d_min( n );
2716 ios << '\n';
2717 }
2718 else
2719 ios << "# d_min = null\n";
2720
2721 ios << "# optd = " << std::boolalpha << m_optd << '\n';
2722 ios << "# d_opt_delta = " << optd_delta() << '\n';
2723 ios << "# lam_sci = " << lam_sci() << '\n';
2724 ios << "# F0 = " << F0() << '\n';
2725 ios << "# starMag = " << starMag() << '\n';
2726 ios << "# lam_sci = " << lam_sci() << '\n';
2727 ios << "# zeta = " << zeta() << '\n';
2728 ios << "# lam_wfs = " << lam_wfs() << '\n';
2729
2730 if( npix_wfs().size() > 0 )
2731 {
2732 ios << "# npix_wfs = " << npix_wfs( (size_t)0 );
2733 for( size_t n = 1; n < npix_wfs().size(); ++n )
2734 ios << ',' << npix_wfs( n );
2735 ios << '\n';
2736 }
2737 else
2738 ios << "# npix_wfs = null\n";
2739
2740 if( ron_wfs().size() > 0 )
2741 {
2742 ios << "# ron_wfs = " << ron_wfs( (size_t)0 );
2743 for( size_t n = 1; n < ron_wfs().size(); ++n )
2744 ios << ',' << ron_wfs( n );
2745 ios << '\n';
2746 }
2747 else
2748 ios << "# ron_wfs = null\n";
2749
2750 if( Fbg().size() > 0 )
2751 {
2752 ios << "# Fbg = " << Fbg( (size_t)0 );
2753 for( size_t n = 1; n < Fbg().size(); ++n )
2754 ios << ',' << Fbg( n );
2755 ios << '\n';
2756 }
2757 else
2758 ios << "# Fbg = null\n";
2759
2760 if( minTauWFS().size() > 0 )
2761 {
2762 ios << "# minTauWFS = " << minTauWFS( (size_t)0 );
2763 for( size_t n = 1; n < minTauWFS().size(); ++n )
2764 ios << ',' << minTauWFS( n );
2765 ios << '\n';
2766 }
2767 else
2768 ios << "# minTauWFS = null\n";
2769
2770 ios << "# bin_npix = " << std::boolalpha << m_bin_npix << '\n';
2771 ios << "# tauWFS = " << tauWFS() << '\n';
2772 ios << "# optTau = " << std::boolalpha << m_optTau << '\n';
2773 ios << "# deltaTau = " << deltaTau() << '\n';
2774 ios << "# fit_mn_max = " << m_fit_mn_max << '\n';
2775 ios << "# spatialFilter_ku = " << m_spatialFilter_ku << '\n';
2776 ios << "# spatialFilter_kv = " << m_spatialFilter_kv << '\n';
2777 ios << "# ncp_wfe = " << m_ncp_wfe << '\n';
2778 ios << "# ncp_alpha = " << m_ncp_alpha << '\n';
2779
2780 if( m_wfsBeta == 0 )
2781 mxThrowException( err::paramnotset, "aoSystem::beta_p", "The WFS is not assigned." );
2782 m_wfsBeta->dumpWFS( ios );
2783 psd.dumpPSD( ios );
2784 atm.dumpAtmosphere( ios );
2785
2786 dumpGitStatus( ios );
2787
2788 return ios;
2789}
2790
2791template <typename realT, class inputSpectT, typename iosT>
2793{
2794 using namespace mx::app;
2795
2796 // AO System configuration
2797 config.add( "aosys.wfs",
2798 "",
2799 "aosys.wfs",
2800 argType::Required,
2801 "aosys",
2802 "wfs",
2803 false,
2804 "string",
2805 "The WFS type: idealWFS, unmodPyWFS, asympModPyWFS, shwfs, calculatedWFS" );
2806 config.add( "aosys.wfs_beta_p",
2807 "",
2808 "aosys.wfs_beta_p",
2809 argType::Required,
2810 "aosys",
2811 "wfs_beta_p",
2812 false,
2813 "string",
2814 "The beta_p file path for calcualtedWFS" );
2815 config.add( "aosys.wfs_beta_r",
2816 "",
2817 "aosys.wfs_beta_r",
2818 argType::Required,
2819 "aosys",
2820 "wfs_beta_r",
2821 false,
2822 "string",
2823 "The beta_r file path for calcualtedWFS" );
2824 config.add( "aosys.wfs_sensitivity",
2825 "",
2826 "aosys.wfs_sensitivity",
2827 argType::Required,
2828 "aosys",
2829 "wfs_sensitivity",
2830 false,
2831 "bool",
2832 "Flag indicating that beta_p/beta_r are sensitivities (inverse) [default false]" );
2833 config.add(
2834 "aosys.D", "", "aosys.D", argType::Required, "aosys", "D", false, "real", "The telescope diameter [m]" );
2835 config.add( "aosys.d_min",
2836 "",
2837 "aosys.d_min",
2838 argType::Required,
2839 "aosys",
2840 "d_min",
2841 false,
2842 "real",
2843 "The minimum actuator spacing [m]" );
2844 config.add( "aosys.optd",
2845 "",
2846 "aosys.optd",
2847 argType::Optional,
2848 "aosys",
2849 "optd",
2850 false,
2851 "bool",
2852 "Whether or not the actuator spacing is optimized" );
2853 config.add( "aosys.optd_delta",
2854 "",
2855 "aosys.optd_delta",
2856 argType::Required,
2857 "aosys",
2858 "optd_delta",
2859 false,
2860 "bool",
2861 "The fractional change from d_min used in optimization. Set to 1 (default) for integer binnings, > 1 "
2862 "for finer sampling." );
2863 config.add( "aosys.F0",
2864 "",
2865 "aosys.F0",
2866 argType::Required,
2867 "aosys",
2868 "F0",
2869 false,
2870 "real",
2871 "Zero-mag photon flux, [photons/sec]" );
2872 config.add( "aosys.lam_wfs",
2873 "",
2874 "aosys.lam_wfs",
2875 argType::Required,
2876 "aosys",
2877 "lam_wfs",
2878 false,
2879 "real",
2880 "WFS wavelength [m]" );
2881 config.add( "aosys.npix_wfs",
2882 "",
2883 "aosys.npix_wfs",
2884 argType::Required,
2885 "aosys",
2886 "npix_wfs",
2887 false,
2888 "vector<real>",
2889 "The number of pixels in the WFS" );
2890 config.add( "aosys.ron_wfs",
2891 "",
2892 "aosys.ron_wfs",
2893 argType::Required,
2894 "aosys",
2895 "ron_wfs",
2896 false,
2897 "vector<real>",
2898 "WFS readout noise [photons/read]" );
2899 config.add( "aosys.bin_npix",
2900 "",
2901 "aosys.bin_npix",
2902 argType::Required,
2903 "aosys",
2904 "bin_npix",
2905 false,
2906 "bool",
2907 "Whether or not WFS pixels are re-binned along with actuator spacing optimization" );
2908 config.add( "aosys.Fbg",
2909 "",
2910 "aosys.Fbg",
2911 argType::Required,
2912 "aosys",
2913 "Fbg",
2914 false,
2915 "vector<real>",
2916 "Background counts, [counts/pix/sec]" );
2917 config.add( "aosys.tauWFS",
2918 "",
2919 "aosys.tauWFS",
2920 argType::Required,
2921 "aosys",
2922 "tauWFS",
2923 false,
2924 "real",
2925 "WFS integration time [s]" );
2926 config.add( "aosys.minTauWFS",
2927 "",
2928 "aosys.minTauWFS",
2929 argType::Required,
2930 "aosys",
2931 "minTauWFS",
2932 false,
2933 "vector<real>",
2934 "Minimum WFS integration time [s]" );
2935 config.add( "aosys.deltaTau",
2936 "",
2937 "aosys.deltaTau",
2938 argType::Required,
2939 "aosys",
2940 "deltaTau",
2941 false,
2942 "real",
2943 "Loop delay [s]" );
2944 config.add( "aosys.optTau",
2945 "",
2946 "aosys.optTau",
2947 argType::Optional,
2948 "aosys",
2949 "optTau",
2950 false,
2951 "bool",
2952 "Whether or not the integration time is optimized" );
2953 config.add( "aosys.lam_sci",
2954 "",
2955 "aosys.lam_sci",
2956 argType::Required,
2957 "aosys",
2958 "lam_sci",
2959 false,
2960 "real",
2961 "Science wavelength [m]" );
2962 config.add(
2963 "aosys.zeta", "", "aosys.zeta", argType::Required, "aosys", "zeta", false, "real", "Zenith distance [rad]" );
2964 config.add( "aosys.fit_mn_max",
2965 "",
2966 "aosys.fit_mn_max",
2967 argType::Required,
2968 "aosys",
2969 "fit_mn_max",
2970 false,
2971 "real",
2972 "Maximum spatial frequency index to use for analysis" );
2973 config.add( "aosys.circularLimit",
2974 "",
2975 "aosys.circularLimit",
2976 argType::Optional,
2977 "aosys",
2978 "circularLimit",
2979 false,
2980 "bool",
2981 " Flag to indicate that the spatial frequency limit is circular, not square." );
2982 config.add( "aosys.spatialFilter_ku",
2983 "",
2984 "aosys.spatialFilter_ku",
2985 argType::Required,
2986 "aosys",
2987 "spatialFilter_ku",
2988 false,
2989 "real",
2990 "Spatial filter cutoff frequency in u [m^-1]" );
2991 config.add( "aosys.spatialFilter_kv",
2992 "",
2993 "aosys.spatialFilter_kv",
2994 argType::Required,
2995 "aosys",
2996 "spatialFilter_kv",
2997 false,
2998 "real",
2999 "Spatial filter cutoff frequency in v [m^-1]" );
3000 config.add( "aosys.ncp_wfe",
3001 "",
3002 "aosys.ncp_wfe",
3003 argType::Required,
3004 "aosys",
3005 "ncp_wfe",
3006 false,
3007 "real",
3008 "NCP WFE between 1 lambda/D and fit_mn_max [rad^2]" );
3009 config.add( "aosys.ncp_alpha",
3010 "",
3011 "aosys.ncp_alpha",
3012 argType::Required,
3013 "aosys",
3014 "ncp_alpha",
3015 false,
3016 "real",
3017 "PSD index for NCP WFE" );
3018 config.add(
3019 "aosys.starMag", "", "aosys.starMag", argType::Required, "aosys", "starMag", false, "real", "Star magnitude" );
3020 config.add( "aosys.starMags",
3021 "",
3022 "aosys.starMags",
3023 argType::Required,
3024 "aosys",
3025 "starMags",
3026 false,
3027 "real vector",
3028 "A vector of star magnitudes" );
3029
3030 atm.setupConfig( config );
3031 psd.setupConfig( config );
3032}
3033
3034template <typename realT, class inputSpectT, typename iosT>
3036{
3037 // WFS
3038 if( config.isSet( "aosys.wfs" ) )
3039 {
3040 std::string wfsStr;
3041 config( wfsStr, "aosys.wfs" );
3042
3043 if( wfsStr == "ideal" )
3044 {
3045 wfsBeta<wfs<realT, iosT>>( nullptr );
3046 }
3047 else if( wfsStr == "unmodPyWFS" )
3048 {
3049 wfsBeta<pywfsUnmod<realT, iosT>>( nullptr );
3050 }
3051 else if( wfsStr == "asympModPyWFS" )
3052 {
3053 wfsBeta<pywfsModAsymptotic<realT>>( nullptr );
3054 }
3055 else if( wfsStr == "SHWFS" )
3056 {
3057 wfsBeta<shwfs<realT, iosT>>( nullptr );
3058 }
3059 else if( wfsStr == "calculatedWFS" )
3060 {
3061 wfsBeta<calculatedWFS<realT, iosT>>( nullptr );
3062
3063 calculatedWFS<realT, iosT> *cwfs = static_cast<calculatedWFS<realT, iosT> *>( m_wfsBeta );
3064 config( cwfs->m_beta_p_file, "aosys.wfs_beta_p" );
3065 config( cwfs->m_beta_r_file, "aosys.wfs_beta_r" );
3066 bool sens = cwfs->m_sensitivity;
3067 config( sens, "aosys.wfs_sensitivity" );
3068 if( config.isSet( "aosys.wfs_sensitivity" ) )
3069 {
3070 cwfs->m_sensitivity = sens;
3071 }
3072 }
3073 else
3074 {
3075 mxThrowException( err::invalidconfig, "aoSystem::loadConfig", "unknown WFS " + wfsStr + " specified" );
3076 }
3077 }
3078
3079 // We load the default value, get the config value (which will be the default if not set), and
3080 // then if it was set, call the setter function so any side effects are captured.
3081
3082 // diameter
3083 realT nD = D();
3084 config( nD, "aosys.D" );
3085 if( config.isSet( "aosys.D" ) )
3086 D( nD );
3087
3088 // d_min
3089 std::vector<realT> nd_min = d_min();
3090 config( nd_min, "aosys.d_min" );
3091 if( config.isSet( "aosys.d_min" ) )
3092 d_min( nd_min );
3093
3094 bool noptd = optd();
3095 config( noptd, "aosys.optd" );
3096 if( config.isSet( "aosys.optd" ) )
3097 optd( noptd );
3098
3099 realT noptd_delta = optd_delta();
3100 config( noptd_delta, "aosys.optd_delta" );
3101 if( config.isSet( "aosys.optd_delta" ) )
3102 optd_delta( noptd_delta );
3103
3104 realT nlam_wfs = lam_wfs();
3105 config( nlam_wfs, "aosys.lam_wfs" );
3106 if( config.isSet( "aosys.lam_wfs" ) )
3107 lam_wfs( nlam_wfs );
3108
3109 // npix_wfs
3110 std::vector<realT> nnpix_wfs = npix_wfs();
3111 config( nnpix_wfs, "aosys.npix_wfs" );
3112 if( config.isSet( "aosys.npix_wfs" ) )
3113 npix_wfs( nnpix_wfs );
3114
3115 // ron_wfs
3116 std::vector<realT> nron_wfs = ron_wfs();
3117 config( nron_wfs, "aosys.ron_wfs" );
3118 if( config.isSet( "aosys.ron_wfs" ) )
3119 ron_wfs( nron_wfs );
3120
3121 // Fbg
3122 std::vector<realT> nFbg = Fbg();
3123 config( nFbg, "aosys.Fbg" );
3124 if( config.isSet( "aosys.Fbg" ) )
3125 Fbg( nFbg );
3126
3127 // minTauWFS
3128 std::vector<realT> nminTauWFS = minTauWFS();
3129 config( nminTauWFS, "aosys.minTauWFS" );
3130 if( config.isSet( "aosys.minTauWFS" ) )
3131 minTauWFS( nminTauWFS );
3132
3133 // bin_npix
3134 bool nbin_npix = bin_npix();
3135 config( nbin_npix, "aosys.bin_npix" );
3136 if( config.isSet( "aosys.bin_npix" ) )
3137 bin_npix( nbin_npix );
3138
3139 // tauWFS
3140 realT ntauWFS = tauWFS();
3141 config( ntauWFS, "aosys.tauWFS" );
3142 if( config.isSet( "aosys.tauWFS" ) )
3143 tauWFS( ntauWFS );
3144
3145 // deltaTau
3146 realT ndeltaTau = deltaTau();
3147 config( ndeltaTau, "aosys.deltaTau" );
3148 if( config.isSet( "aosys.deltaTau" ) )
3149 deltaTau( ndeltaTau );
3150
3151 // optTau
3152 bool noptTau = optTau();
3153 config( noptTau, "aosys.optTau" );
3154 if( config.isSet( "aosys.optTau" ) )
3155 optTau( noptTau );
3156
3157 // lam_sci
3158 realT nlam_sci = lam_sci();
3159 config( nlam_sci, "aosys.lam_sci" );
3160 if( config.isSet( "aosys.lam_sci" ) )
3161 lam_sci( nlam_sci );
3162
3163 // zeta
3164 realT nzeta = zeta();
3165 config( nzeta, "aosys.zeta" );
3166 if( config.isSet( "aosys.zeta" ) )
3167 zeta( nzeta );
3168
3169 // fit_mn_max
3170 realT fmnm = fit_mn_max();
3171 config( fmnm, "aosys.fit_mn_max" );
3172 if( config.isSet( "aosys.fit_mn_max" ) )
3173 fit_mn_max( fmnm );
3174
3175 // circularLimit
3176 bool cl = circularLimit();
3177 config( cl, "aosys.circularLimit" );
3178 if( config.isSet( "aosys.circularLimit" ) )
3179 circularLimit( cl );
3180
3181 // spatialFilter_ku
3182 realT ku = spatialFilter_ku();
3183 config( ku, "aosys.spatialFilter_ku" );
3184 if( config.isSet( "aosys.spatialFilter_ku" ) )
3185 spatialFilter_ku( ku );
3186
3187 // spatialFilter_kv
3188 realT kv = spatialFilter_kv();
3189 config( kv, "aosys.spatialFilter_kv" );
3190 if( config.isSet( "aosys.spatialFilter_kv" ) )
3191 spatialFilter_kv( kv );
3192
3193 // ncp_wfe
3194 realT nwfe = ncp_wfe();
3195 config( nwfe, "aosys.ncp_wfe" );
3196 if( config.isSet( "aosys.ncp_wfe" ) )
3197 ncp_wfe( nwfe );
3198
3199 // ncp_alpha
3200 realT na = ncp_alpha();
3201 config( na, "aosys.ncp_alpha" );
3202 if( config.isSet( "aosys.ncp_alpha" ) )
3203 ncp_alpha( na );
3204
3205 // F0
3206 realT nF0 = F0();
3207 config( nF0, "aosys.F0" );
3208 if( config.isSet( "aosys.F0" ) )
3209 F0( nF0 );
3210
3211 // star_mag
3212 realT smag = starMag();
3213 config( smag, "aosys.starMag" );
3214 if( config.isSet( "aosys.starMag" ) )
3215 starMag( smag );
3216
3217 atm.loadConfig( config );
3218 psd.loadConfig( config );
3219}
3220
3221extern template class aoSystem<float, vonKarmanSpectrum<float>, std::ostream>;
3222
3223extern template class aoSystem<double, vonKarmanSpectrum<double>, std::ostream>;
3224
3225extern template class aoSystem<long double, vonKarmanSpectrum<long double>, std::ostream>;
3226
3227#ifdef HASQUAD
3228extern template class aoSystem<__float128, vonKarmanSpectrum<__float128>, std::ostream>;
3229#endif
3230
3231} // namespace analysis
3232} // namespace AO
3233} // namespace mx
3234
3235#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:65
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 m_opticalGain
The optical gain of the WFS, 0-1. Treated as a sensitivity reduction in measurement noise....
Definition aoSystem.hpp:89
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..
mxException for invalid config settings
mxException for parameters which aren't set
mxException for a size error
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.
The mxlib c++ namespace.
Definition mxError.hpp:106
iosT & dumpGitStatus(iosT &ios)
Dump the current git status of the library to a stream.
Definition mxlib.hpp:50
The calculated WFS uses sensitivities provided by FITS files.
Definition aoWFS.hpp:276
The ideal wavefront sensor sensitivity function.
Definition aoWFS.hpp:37
Class to manage a set of configurable values, and read their values from config/ini files and the com...
void add(const configTarget &tgt)
Add a configTarget.
bool isSet(const std::string &name, std::unordered_map< std::string, configTarget > &targets)
Check if a target has been set by the configuration.