mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
Loading...
Searching...
No Matches
pyramidSensor.hpp
Go to the documentation of this file.
1/** \file pyramidSensor.hpp
2 * \author Jared R. Males (jaredmales@gmail.com)
3 * \brief Declaration and definition of a standard 4 quadrant pyramid WFS.
4 * \ingroup mxAO_sim_files
5 *
6 */
7
8#ifndef mx_AO_sim_pyramidSensor_hpp
9#define mx_AO_sim_pyramidSensor_hpp
10
11#include "../../mxException.hpp"
12
13#include "../../wfp/imagingUtils.hpp"
14#include "../../wfp/fraunhoferPropagator.hpp"
15#include "../../sys/timeUtils.hpp"
16
17#include "../../improc/eigenImage.hpp"
18#include "../../improc/imageMasks.hpp"
19
20#include "../../math/constants.hpp"
21#include "../../math/geo.hpp"
22
23#include "wavefront.hpp"
24
25#ifdef DEBUG
26#define BREAD_CRUMB std::cout << "DEBUG: " << __FILE__ << " " << __LINE__ << "\n";
27#else
28#define BREAD_CRUMB
29#endif
30
31namespace mx
32{
33namespace AO
34{
35namespace sim
36{
37
38template <typename _realT>
39struct wfsImageT
40{
41 typedef _realT realT;
42
43 unsigned iterNo;
44
45 /// The wavefront sensor detector image type
46 typedef Eigen::Array<realT, Eigen::Dynamic, Eigen::Dynamic> imageT;
47
48 imageT image;
49
50 imageT tipImage;
51};
52
53/// A Pyramid Sensor Simulation
54/**
55 *
56 * \tparam _realT is the real floating point type used for calculations
57 * \tparam _detectorT is the detector used to record the PyWFS image. Must conform to the mx::AO::sim detector
58 * interface specifications.
59 */
60template <typename _realT, typename _detectorT>
61class pyramidSensor
62{
63 public:
64 /// The real floating point type used for calculations
65 typedef _realT realT;
66
67 /// The complex floating point type used for calculations
68 typedef std::complex<realT> complexT;
69
70 /// The wavefront data type
71 typedef wavefront<realT> wavefrontT;
72
73 /// The wavefront complex field type
74 typedef wfp::imagingArray<std::complex<realT>, wfp::fftwAllocator<std::complex<realT>>, 0> complexFieldT;
75
76 /// The wavefront sensor detector image type
77 typedef _detectorT detectorT;
78
79 public:
80 /// Default c'tor
81 pyramidSensor();
82
83 /** \name Standard WFS Interface
84 *
85 * @{
86 */
87
88 protected:
89 /* Standard WFS Interface: */
90 uint32_t m_wfSz{ 0 }; ///< Size of the wavefront in pixels
91
92 uint32_t m_detRows{ 0 }; ///< The number of rows of the WFS detector. After forming the image the WFS detector
93 ///< plane is re-binned to this.
94
95 uint32_t m_detCols{ 0 }; ///< The number of columns of the WFS detector. After forming the image the WFS detector
96 ///< plane is re-binned to this.
97
98 realT m_lambda{ 0 }; ///< Central wavelength, in meters
99
100 /// \todo when the filter should be set with astrospectrum, and should re-calculate the central wavelength.
101 /// \todo need to verify that the wavefront propagation is appropriately chromatic
102 std::vector<realT> m_wavelengths; ///< Vector of wavelengths in the WFS bandpass
103 std::vector<realT> _wavelengthWeights; ///< The relative weights of the wavelengths
104
105 int m_iTime{ 1 }; ///< Integration time in loop steps
106
107 int m_roTime{ 1 }; ///< Readout time in loop steps
108
109 realT m_simStep{ 0.001 }; ///< The simulation stepsize in seconds.
110
111 public:
112 /// Get the wavefront size in pixels
113 /**
114 * \returns the wavefront size in pixels
115 */
116 int wfSz();
117
118 /// Set the wavefront size in pixels.
119 /**
120 */
121 void wfSz(const uint32_t & sz /**< the new size*/);
122 void wfSz(const uint32_t & sz /**< the new size*/);
123
124 /// Get the detector rows in pixels
125 /**
126 * \returns m_detRows
127 */
128 uint32_t detRows();
129
130 /// Get the detector columns in pixels
131 /**
132 * \returns m_detCols
133 */
134 uint32_t detCols();
135
136 /// Set the detector columns in pixels.
137 /**
138 */
139 void detSize( const uint32_t &nrows, ///< The number of rows
140 const uint32_t &ncols ///< The number of columns
141 );
142
143 /// Get the PyWFS central wavelength
144 /**
145 * \returns the central wavelength in meters
146 */
147 realT lambda();
148
149 /// Set the PyWFS central wavelength
150 /**
151 */
152 void lambda( const realT &l /**< The central wavelength, in meters*/ );
153
154 /// Get the PyWFS integration time, in time steps
155 int iTime();
156
157 /// Set the PyWFS integration time, in time steps
158 void iTime( const uint32_t &it /**< the new integration time*/ );
159
160 /// Get the PyWFS detector readout time, in time steps
161 int roTime();
162
163 /// Set the PyWFS detector readout time, in time steps
164 void roTime( const uint32_t &rt /**< the new readout time*/ );
165
166 /// Get the simulation step-size, in seconds.
167 realT simStep();
168
169 /// Set the simulation step-size, in seconds.
170 void simStep( const realT &st /**< the new simulation step size*/ );
171
172 /// Link this WFS to an AO system simulation
173 template <typename AOSysT>
174 void linkSystem( AOSysT &AOSys /**< The AO system simulation to link to*/ );
175
176 /// Sense the wavefront aberrations
177 /** \returns true if a new wavefront measurement is ready.
178 * \returns false if still integrating.
179 */
180 bool senseWavefront( wavefrontT &pupilPlane /**< The input wavefront to be sensed*/ );
181
182 /// Sense the wavefront aberrations in calibration mode
183 /** Allows for faster calibrations.
184 */
185 bool senseWavefrontCal( wavefrontT &pupilPlane /**< The input wavefront to be sensed*/ );
186
187 public:
188 /// The WFS detector
189 detectorT detector;
190
191 /// The image on the detector, resized from m_wfsImage
192 wfsImageT<realT> detectorImage;
193
194 /// @}
195
196 /** \name Pyramid Sensor Interface
197 *
198 * @{
199 */
200protected:
201 uint32_t m_nSides {4}; ///< Number of sides in the pyramid
202 uint32_t m_nSides {4}; ///< Number of sides in the pyramid
203
204 /// The size of the pupil in wavefront pixels.
205 /** This is the maximum diameter of the pupil in wavefront pixels.
206 *
207 */
208 uint32_t m_pupilSz{ 0 };
209
210 /// The separation of the pupil images in fraction of a pupil. 0 <= m_pupilSep, default 1.
211 /// The separation of the pupil images in fraction of a pupil. 0 <= m_pupilSep, default 1.
212 /** This sets the center-to-center separation of the pupils images in the focal plane wavefront.
213 * Note that the separation in detector pixels then depends on the scaling between wavefront pixels
214 * (m_wfSz) and detector pixels (m_detRows and m_detCols).
215 *
216 *
217 *
218 * This sets the size of the region in the pre-detection image that each pupil image
219 * takes up, and therefore the size of the pre-detection image.
220 * If the pupil (as defined in the input wavefront) is 60 pixels across
221 * and m_pupilSep is set to 1.06667, then there will be a 2 pixel pad around each pupil image,
222 * and m_pupilSep is set to 1.06667, then there will be a 2 pixel pad around each pupil image,
223 * resulting in 4 pixels between each geometric pupil image.
224 *
225 * For a standard 4-sided pyramid, the pre-detection image will be
226 * 2*m_pupilSep*m_pupilSz across. For other n-sided pyramids, m_pupilSep still specifies the size of the pupil
227 * 2*m_pupilSep*m_pupilSz across. For other n-sided pyramids, m_pupilSep still specifies the size of the pupil
228 * image region, but the total image size will be a function of the resultant pupil positions.
229 *
230 * If m_pupilSep is less than 1, this will produce the "flattened pyramid", with overlap between
231 * If m_pupilSep is less than 1, this will produce the "flattened pyramid", with overlap between
232 * the pupil images. In this case, image size will also be set by pupilSz to ensure that there are enough
233 * pixels included to show all pupils.
234 */
235 realT m_pupilSep {1};
236
237 /// The angle by which to offset the pupils, in degrees. Default is 0.
238 /** If this is 0, then a 4-sided pyramid makes a square as usual. If this is set
239 * to 45 degrees, then a 4-sided pyramid makes a diamond.
240 *
241 */
242 realT m_angleOffset {0};
243 realT m_pupilSep {1};
244
245 /// The angle by which to offset the pupils, in degrees. Default is 0.
246 /** If this is 0, then a 4-sided pyramid makes a square as usual. If this is set
247 * to 45 degrees, then a 4-sided pyramid makes a diamond.
248 *
249 */
250 realT m_angleOffset {0};
251
252 /// The size of the resulting PyWFS image in wavefront pixels.
253 /** If \ref m_imageSzAuto is true, this is determined by number of sides (\ref m_nSides), the pupil size (\ref m_pupilSz), and the
254 * pupil separation (\ref m_pupilSep). For a 4 sided pyramid this will be the larger of
255 * 2*m_pupilSep*m_pupilSz and 2*m_pupilSz.
256 * 2*m_pupilSep*m_pupilSz and 2*m_pupilSz.
257 *
258 * If , then this is used regardless of the optimum size.
259 */
260 uint32_t m_imageSz {0};
261
262 bool m_imageSzAuto{ true }; ///< Flag to track if \ref m_imageSz should be set to 0.
263
264 realT m_wfPS {0}; ///< Wavefront pixel scale, in meters/pixel
265 realT m_wfPS {0}; ///< Wavefront pixel scale, in meters/pixel
266
267 realT m_D {0}; ///< Telescope diameter, in meters
268 realT m_D {0}; ///< Telescope diameter, in meters
269
270 uint32_t m_modSteps{ 20 }; ///< Number of modulation steps in one integration. Can be set explicitly, but will be
271 ///< calculated if \ref m_perStep is set.
272
273 realT m_perStep{
274 1 }; ///< The minimum number of lamba/D per step in the modulation. Smaller will result in more steps.
275
276 realT m_modRadius{ 3.0 }; ///< Radius of the modulation in pixels
277
278 public:
279 /// Get the number of pyramid sides
280 /**
281 * \returns the number of sides on the pyramid
282 */
283 int nSides();
284
285 /// Set the number of sides on the pyramid
286 /**
287 */
288 void nSides(const uint32_t & ns /**< The new number of sides on the pyramid*/);
289 void nSides(const uint32_t & ns /**< The new number of sides on the pyramid*/);
290
291 /// Get the minimum number of modulation steps
292 /**
293 * \returns m_perStep;
294 */
295 realT perStep();
296
297 /// Set the minimum number of modulation steps
298 /**
299 * \param mSt is the new number of modulation steps
300 */
301 void perStep( const realT &prStp /**< The minimum number of lamba/D per step to take*/ );
302
303 /// Get the number of modulation steps
304 /**
305 * \returns m_modSteps, which is defined by perStep.
306 */
307 int modSteps();
308
309 /// Get the radius of modulation
310 /**
311 * \returns m_modRadius;
312 */
313 realT modRadius();
314
315 /// Set the modulation radius
316 /**
317 * \param mR is the new modulation radius in lambda/D
318 */
319 void modRadius( const realT &mR /**< [in] the new value of modulation radius */ );
320
321 /// Get the wavefront pixel scale in meters per pixel
322 /**
323 * \returns the wavefront pixel scale in meters/pixel
324 */
325 realT wfPS();
326
327 /// Get the telescope diameter
328 /**
329 * \returns the telescope diameter in meters
330 */
331 realT D();
332
333 /// Set the telescope diameter
334 /**
335 * \param d is the new size in meters
336 */
337 void D( const realT &d /**< */ );
338
339 /// Get the pupil size in pixels
340 /** This is the pupil size in un-binned wavefront space
341 *
342 * \returns m_pupilSz
343 * \returns m_pupilSz
344 */
345 uint32_t pupilSz();
346
347 /// Set the pupil size in pixels.
348 /** This is the size of the pupils in un-binned wavefront space.
349 * See \ref m_pupilSz.
350 *
351 */
352 void pupilSz( const uint32_t &sz /**< the new pupil size.*/ );
353
354 /// Get the pupil separation as a fraction of pupil size
355 /// Get the pupil separation as a fraction of pupil size
356 /** This is the pupil separation in un-binned wavefront space
357 *
358 * \returns m_pupilSep
359 */
360 realT pupilSep();
361 realT pupilSep();
362
363 /// Set the pupil separation as a fraction of pupil size
364 /// Set the pupil separation as a fraction of pupil size
365 /** This is the separation of the pupils in un-binned wavefront space.
366 * See \ref m_pupilSep.
367 *
368 */
369 void pupilSep(const realT & sz /**< the new pupil separation.*/);
370
371 /// Get the angle offset
372 /** See \ref m_angleOffset
373 *
374 * \returns m_angleOffset
375 */
376 realT angleOffset();
377
378 /// Set the angle offset
379 /** See \ref m_angleOffset.
380 *
381 */
382 void angleOffset(const realT & ao /**< the new angle offset.*/);
383 void pupilSep(const realT & sz /**< the new pupil separation.*/);
384
385 /// Get the angle offset
386 /** See \ref m_angleOffset
387 *
388 * \returns m_angleOffset
389 */
390 realT angleOffset();
391
392 /// Set the angle offset
393 /** See \ref m_angleOffset.
394 *
395 */
396 void angleOffset(const realT & ao /**< the new angle offset.*/);
397
398 /// Get the image size in wavefront pixels
399 /** This is the size of the image in un-binned wavefront space
400 *
401 * \returns m_imageSz
402 */
403 uint32_t imageSz();
404
405 /// Set the image size in wavefront pixels
406 /** This is the size of the image in un-binned wavefront space
407 * Setting a non-zero value also sets m_imageSizeAuto to false.
408 * Setting 0 also sets m_imageSizeAuto to true.
409 */
410 void imageSz( const uint32_t &is );
411
412 /// Get the value of the image size auto flag
413 /** This controls whether image size is set automatically
414 *
415 * \returns m_imageSz
416 */
417 bool imageSzAuto();
418
419 /// Set the value of the image size auto flag
420 /** This controls whether image size is set automatically
421 *
422 */
423 void imageSzAuto( const bool &ia );
424 ///@}
425
426 wfp::fraunhoferPropagator<complexFieldT> m_frProp;
427
428 bool m_opdMaskMade{ false };
429 complexFieldT m_opdMask;
430
431 bool m_tiltsMade{ false };
432 std::vector<complexFieldT> m_tilts;
433
434 bool m_preAllocated{ false };
435 complexFieldT m_pupilPlaneCF;
436
437 // Pre-allocated working memory:
438
439 std::vector<complexFieldT> m_th_tiltedPlane; ///< Thread-local modulated wavefront
440
441 std::vector<complexFieldT> m_th_focalPlane; ///< Thread-local tip wavefront, used for FFT tilting
442
443 std::vector<typename wfsImageT<realT>::imageT> m_th_focalImage; ///< Thread-local tip image
444
445 std::vector<complexFieldT> m_th_sensorPlane; ///< Thread-local sensor-pupil-plane wavefront
446
447 std::vector<typename wfsImageT<realT>::imageT>
448 m_th_sensorImage; ///< Thread-local sensor-pupil-plane intensity image
449
450 int m_iTime_counter{ 0 };
451
452 int m_reading{ 0 };
453
454 int m_roTime_counter{ 0 };
455
456 std::vector<wavefrontT> _wavefronts;
457
458 int m_lastWavefront{ 0 };
459
460 /// The image formed by the WFS
461 wfsImageT<realT> m_wfsImage;
462
463 public:
464 wfsImageT<realT> wfsTipImage;
465
466 protected:
467 void makeOpdMask();
468
469 void makeTilts();
470
471 void allocThreadMem();
472
473 void preAllocate();
474
475 void doSenseWavefront();
476 void doSenseWavefront( wavefrontT & /**< */ );
477 void doSenseWavefrontNoMod( wavefrontT & /**< */ );
478
479 bool m_firstRun{ true };
480};
481
482template <typename realT, typename detectorT>
483pyramidSensor<realT, detectorT>::pyramidSensor()
484{
485 iTime( m_iTime );
486
487 m_frProp.wholePixel( 0 );
488}
489
490template <typename realT, typename detectorT>
491int pyramidSensor<realT, detectorT>::wfSz()
492{
493 return m_wfSz;
494}
495
496template <typename realT, typename detectorT>
497void pyramidSensor<realT, detectorT>::wfSz(const uint32_t & sz)
498void pyramidSensor<realT, detectorT>::wfSz(const uint32_t & sz)
499{
500 if( m_wfSz == sz )
501 {
502 return;
503 }
504
505 m_wfSz = sz;
506
507 m_tiltsMade = false;
508 m_opdMaskMade = false;
509 m_preAllocated = false;
510}
511
512template <typename realT, typename detectorT>
513uint32_t pyramidSensor<realT, detectorT>::detRows()
514uint32_t pyramidSensor<realT, detectorT>::detRows()
515{
516 return m_detRows;
517}
518
519template <typename realT, typename detectorT>
520uint32_t pyramidSensor<realT, detectorT>::detCols()
521uint32_t pyramidSensor<realT, detectorT>::detCols()
522{
523 return m_detCols;
524}
525
526template <typename realT, typename detectorT>
527void pyramidSensor<realT, detectorT>::detSize( const uint32_t &nrows, const uint32_t &ncols )
528{
529 if( m_detRows == nrows && m_detCols == ncols )
530 return;
531
532 m_detRows = nrows;
533 m_detCols = ncols;
534
535 detector.setSize(m_detRows, m_detCols);
536 detectorImage.image.resize(m_detRows, m_detCols);
537
538 m_opdMaskMade = false; //make sure size check is run on current settings
539
540 m_opdMaskMade = false; //make sure size check is run on current settings
541}
542
543template <typename realT, typename detectorT>
544realT pyramidSensor<realT, detectorT>::lambda()
545{
546 return m_lambda;
547}
548
549template <typename realT, typename detectorT>
550void pyramidSensor<realT, detectorT>::lambda( const realT &l )
551{
552 m_lambda = l;
553
554 //---------------------------------------
555 // Check if wavelength vector is filled
556 //---------------------------------------
557 if( m_wavelengths.size() == 0 )
558 {
559 m_wavelengths.resize( 1, m_lambda );
560 _wavelengthWeights.resize( 1, 1.0 );
561 }
562}
563
564template <typename realT, typename detectorT>
565int pyramidSensor<realT, detectorT>::iTime()
566{
567 return m_iTime;
568}
569
570template <typename realT, typename detectorT>
571void pyramidSensor<realT, detectorT>::iTime( const uint32_t &it )
572{
573 if( it < 1 )
574 {
575 mxThrowException( mx::err::invalidconfig, "pyramidSensor::iTime", "iTime must be >= 1" );
576 }
577
578 m_iTime = it;
579
580 _wavefronts.resize( m_iTime + 2 );
581 m_lastWavefront = -1;
582
583 detector.expTime( m_simStep * m_iTime );
584}
585
586template <typename realT, typename detectorT>
587int pyramidSensor<realT, detectorT>::roTime()
588{
589 return roTime;
590}
591
592template <typename realT, typename detectorT>
593void pyramidSensor<realT, detectorT>::roTime( const uint32_t &rt )
594{
595 if( rt < 1 )
596 {
597 mxThrowException( mx::err::invalidconfig, "pyramidSensor::roTime", "roTime must be >= 1" );
598 }
599
600 m_roTime = rt;
601}
602
603template <typename realT, typename detectorT>
604realT pyramidSensor<realT, detectorT>::simStep()
605{
606 return m_simStep;
607}
608
609template <typename realT, typename detectorT>
610void pyramidSensor<realT, detectorT>::simStep( const realT &st )
611{
612
613 m_simStep = st;
614
615 detector.expTime( m_simStep * m_iTime );
616}
617
618template <typename realT, typename detectorT>
619template <typename AOSysT>
620void pyramidSensor<realT, detectorT>::linkSystem( AOSysT &AOSys )
621{
622 AOSys.wfs.wfPS( AOSys.m_wfPS );
623 AOSys.wfs.D( AOSys.m_D );
624}
625
626template <typename realT, typename detectorT>
627bool pyramidSensor<realT, detectorT>::senseWavefront( wavefrontT &pupilPlane )
628{
629
630 ++m_lastWavefront;
631 if( m_lastWavefront >= _wavefronts.size() )
632 m_lastWavefront = 0;
633 _wavefronts[m_lastWavefront].amplitude = pupilPlane.amplitude;
634 _wavefronts[m_lastWavefront].phase = pupilPlane.phase;
635 _wavefronts[m_lastWavefront].iterNo = pupilPlane.iterNo;
636
637 // Always skip the first one for averaging to center of iTime.
638 if( m_firstRun )
639 {
640 m_firstRun = false;
641 return false;
642 }
643
644 ++m_iTime_counter;
645
646 bool rv = false;
647
648 if( m_reading )
649 {
650 ++m_roTime_counter;
651
652 if( m_roTime_counter >= m_roTime )
653 {
654 detector.exposeImage( detectorImage.image, m_wfsImage.image );
655
656 detectorImage.tipImage = wfsTipImage.image;
657 detectorImage.iterNo = m_wfsImage.iterNo;
658
659 m_roTime_counter = 0;
660 m_reading = 0;
661 rv = true;
662 }
663 }
664
665 if( m_iTime_counter >= m_iTime )
666 {
667 doSenseWavefront();
668 m_iTime_counter = 0;
669
670 m_reading = 1;
671 m_roTime_counter = 0;
672 }
673
674 return rv;
675}
676
677template <typename realT, typename detectorT>
678bool pyramidSensor<realT, detectorT>::senseWavefrontCal( wavefrontT &pupilPlane )
679{
680
681 BREAD_CRUMB;
682
683 doSenseWavefront(pupilPlane);
684
685
686 BREAD_CRUMB;
687
688 detector.exposeImage( detectorImage.image, m_wfsImage.image );
689
690 BREAD_CRUMB;
691
692 detectorImage.tipImage = wfsTipImage.image;
693
694 BREAD_CRUMB;
695
696 return true;
697}
698
699/* Pyramid Specifics */
700
701template <typename realT, typename detectorT>
702int pyramidSensor<realT, detectorT>::nSides()
703{
704 return m_nSides;
705}
706
707template <typename realT, typename detectorT>
708void pyramidSensor<realT, detectorT>::nSides(const uint32_t & ns)
709void pyramidSensor<realT, detectorT>::nSides(const uint32_t & ns)
710{
711 m_nSides = ns;
712 m_opdMaskMade = false;
713}
714
715template <typename realT, typename detectorT>
716realT pyramidSensor<realT, detectorT>::wfPS()
717{
718 return m_wfPS;
719}
720
721template <typename realT, typename detectorT>
722realT pyramidSensor<realT, detectorT>::D()
723{
724 return m_D;
725}
726
727template <typename realT, typename detectorT>
728void pyramidSensor<realT, detectorT>::D( const realT &d )
729{
730 m_D = d;
731
732 if(m_pupilSz > 0) //Avoid making inf or nan so m_wfPS remains unset. note that fast-math make detecting inf and nan hard.
733 {
734 m_wfPS = m_D/m_pupilSz;
735 }
736 else
737 {
738 m_wfPS = 0;
739 }
740
741 m_tiltsMade = false;
742 m_opdMaskMade = false;
743 m_preAllocated = false;
744 m_opdMaskMade = false;
745 m_preAllocated = false;
746}
747
748template <typename realT, typename detectorT>
749realT pyramidSensor<realT, detectorT>::perStep()
750{
751 return m_perStep;
752}
753
754template <typename realT, typename detectorT>
755void pyramidSensor<realT, detectorT>::perStep( const realT &prStp )
756{
757 m_perStep = prStp;
758
759 if( m_modRadius <= 0 )
760 {
761 m_modSteps = 0;
762 return;
763 }
764
765 realT radPerStep = m_perStep / m_modRadius;
766
767 // Get the minimum number of steps to meet m_perStep while having symmetry for the quadrants
768 m_modSteps = 1;
769 while( math::half_pi<realT>() / m_modSteps > radPerStep )
770 {
771 ++m_modSteps;
772 }
773
774 m_modSteps *= 4;
775
776 m_tiltsMade = false;
777}
778
779template <typename realT, typename detectorT>
780int pyramidSensor<realT, detectorT>::modSteps()
781{
782 return m_modSteps;
783}
784
785template <typename realT, typename detectorT>
786realT pyramidSensor<realT, detectorT>::modRadius()
787{
788 return m_modRadius;
789}
790
791template <typename realT, typename detectorT>
792void pyramidSensor<realT, detectorT>::modRadius( const realT &mR )
793{
794 m_modRadius = mR;
795 perStep( m_perStep ); // to calculate m_modSteps;
796
797 m_tiltsMade = false;
798}
799
800template <typename realT, typename detectorT>
801uint32_t pyramidSensor<realT, detectorT>::pupilSz()
802{
803 return m_pupilSz;
804}
805
806template <typename realT, typename detectorT>
807void pyramidSensor<realT, detectorT>::pupilSz( const uint32_t &sz )
808{
809 if( m_pupilSz == sz )
810 {
811 return;
812 }
813
814 m_pupilSz = sz;
815
816 if(m_pupilSz > 0) //Avoid making inf or nan so m_wfPS remains unset. note that fast-math make detecting inf and nan hard.
817 {
818 m_wfPS = m_D/m_pupilSz;
819 }
820 else
821 {
822 m_wfPS = 0;
823 }
824
825 if(m_pupilSz > 0) //Avoid making inf or nan so m_wfPS remains unset. note that fast-math make detecting inf and nan hard.
826 {
827 m_wfPS = m_D/m_pupilSz;
828 }
829 else
830 {
831 m_wfPS = 0;
832 }
833
834 m_tiltsMade = false;
835 m_opdMaskMade = false;
836 m_preAllocated = false;
837}
838
839template <typename realT, typename detectorT>
840realT pyramidSensor<realT, detectorT>::pupilSep()
841realT pyramidSensor<realT, detectorT>::pupilSep()
842{
843 return m_pupilSep;
844}
845
846template <typename realT, typename detectorT>
847void pyramidSensor<realT, detectorT>::pupilSep(const realT & sz)
848void pyramidSensor<realT, detectorT>::pupilSep(const realT & sz)
849{
850 if( m_pupilSep == sz )
851 {
852 return;
853 }
854
855 m_pupilSep = sz;
856
857 m_opdMaskMade = false;
858 m_preAllocated = false;
859}
860
861template <typename realT, typename detectorT>
862realT pyramidSensor<realT, detectorT>::angleOffset()
863{
864 return m_angleOffset;
865}
866
867template <typename realT, typename detectorT>
868void pyramidSensor<realT, detectorT>::angleOffset(const realT & ao)
869{
870 if (m_angleOffset == ao)
871 {
872 return;
873 }
874
875 m_angleOffset = ao;
876
877 m_opdMaskMade = false;
878 m_preAllocated = false;
879}
880
881template <typename realT, typename detectorT>
882realT pyramidSensor<realT, detectorT>::angleOffset()
883{
884 return m_angleOffset;
885}
886
887template <typename realT, typename detectorT>
888void pyramidSensor<realT, detectorT>::angleOffset(const realT & ao)
889{
890 if (m_angleOffset == ao)
891 {
892 return;
893 }
894
895 m_angleOffset = ao;
896
897 m_opdMaskMade = false;
898 m_preAllocated = false;
899}
900
901template <typename realT, typename detectorT>
902uint32_t pyramidSensor<realT, detectorT>::imageSz()
903{
904 if( !m_opdMaskMade )
905 {
906 makeOpdMask();
907 }
908
909 return m_imageSz;
910}
911
912template <typename realT, typename detectorT>
913void pyramidSensor<realT, detectorT>::imageSz( const uint32_t &sz )
914{
915 if( m_imageSz == sz )
916 {
917 return;
918 }
919
920 m_imageSz = sz;
921
922 if( m_imageSz == 0 )
923 {
924 m_imageSzAuto = true;
925 }
926 else
927 {
928 m_imageSzAuto = false;
929 }
930
931 m_opdMaskMade = false;
932 m_preAllocated = false;
933}
934
935template <typename realT, typename detectorT>
936bool pyramidSensor<realT, detectorT>::imageSzAuto()
937{
938 return m_imageSzAuto;
939}
940
941template <typename realT, typename detectorT>
942void pyramidSensor<realT, detectorT>::imageSzAuto( const bool &ia )
943{
944 m_imageSzAuto = ia;
945
946 m_opdMaskMade = false;
947
948 m_preAllocated = false;
949}
950
951template <typename realT, typename detectorT>
952void pyramidSensor<realT, detectorT>::makeOpdMask()
953{
954 complexFieldT opdMaskQ;
955
956 if(m_wfPS == 0)
957 {
958 mxThrowException(mx::err::invalidconfig, "pyramidSensor::makeOpdMask()", "wavefront platescale (m_wfPS) is 0. Must set pupilSz and D first.");
959 }
960
961 if(!std::isfinite(m_wfPS) || !std::isnormal(m_wfPS))
962 {
963 mxThrowException(mx::err::invalidconfig, "pyramidSensor::makeOpdMask()", "wavefront platescale (m_wfPS) is infinite. Must set pupilSz and D first.");
964 }
965
966 std::cerr << m_wfPS << " " << m_D << "\n";
967
968 if(m_D == 0)
969 {
970 mxThrowException(mx::err::invalidconfig, "pyramidSensor::makeOpdMask()", "pupil diameter is 0. Must set D > 0 first.");
971 }
972
973 //Setup the Fraunhoffer Propagator
974 m_frProp.setWavefrontSizePixels(m_wfSz);
975
976 m_opdMask.resize(m_wfSz, m_wfSz);
977 opdMaskQ.resize(m_wfSz, m_wfSz);
978
980
981 mask.resize(m_opdMask.rows(), m_opdMask.cols());
982 realT dang = mx::math::two_pi<realT>()/m_nSides;
983 if(m_wfPS == 0)
984 {
985 mxThrowException(mx::err::invalidconfig, "pyramidSensor::makeOpdMask()", "wavefront platescale (m_wfPS) is 0. Must set pupilSz and D first.");
986 }
987
988 if(!std::isfinite(m_wfPS) || !std::isnormal(m_wfPS))
989 {
990 mxThrowException(mx::err::invalidconfig, "pyramidSensor::makeOpdMask()", "wavefront platescale (m_wfPS) is infinite. Must set pupilSz and D first.");
991 }
992
993 std::cerr << m_wfPS << " " << m_D << "\n";
994
995 if(m_D == 0)
996 {
997 mxThrowException(mx::err::invalidconfig, "pyramidSensor::makeOpdMask()", "pupil diameter is 0. Must set D > 0 first.");
998 }
999
1000 //Setup the Fraunhoffer Propagator
1001 m_frProp.setWavefrontSizePixels(m_wfSz);
1002
1003 m_opdMask.resize(m_wfSz, m_wfSz);
1004 opdMaskQ.resize(m_wfSz, m_wfSz);
1005
1007
1008 mask.resize(m_opdMask.rows(), m_opdMask.cols());
1009 realT dang = mx::math::two_pi<realT>()/m_nSides;
1010
1011 realT minx = 0;
1012 realT maxx = 0;
1013 realT miny = 0;
1014 realT maxy = 0;
1015
1016 realT pupilRad = m_pupilSep*m_pupilSz /(2*sin(dang/2.0));
1017 realT minx = 0;
1018 realT maxx = 0;
1019 realT miny = 0;
1020 realT maxy = 0;
1021
1022 realT pupilRad = m_pupilSep*m_pupilSz /(2*sin(dang/2.0));
1023
1024 for(int n = 0; n < m_nSides; ++n)
1025 {
1026 realT ang = m_angleOffset*math::degreesT<realT>::radians + 0.5*dang + n * dang;
1027 for(int n = 0; n < m_nSides; ++n)
1028 {
1029 realT ang = m_angleOffset*math::degreesT<realT>::radians + 0.5*dang + n * dang;
1030
1031 realT dx = pupilRad * cos(ang);
1032 realT dx = pupilRad * cos(ang);
1033
1034 if(dx < minx) minx = dx;
1035 if(dx > maxx) maxx = dx;
1036 if(dx < minx) minx = dx;
1037 if(dx > maxx) maxx = dx;
1038
1039 realT dy = pupilRad * sin(ang);
1040 realT dy = pupilRad * sin(ang);
1041
1042 if(dy < miny) miny = dy;
1043 if(dy > maxy) maxy = dy;
1044 if(dy < miny) miny = dy;
1045 if(dy > maxy) maxy = dy;
1046
1047 opdMaskQ.set(std::complex<realT>(0, 1));
1048 wfp::tiltWavefront(opdMaskQ, dx, dy);
1049 mask.setZero();
1050 improc::maskWedge(mask, 0.5*(mask.rows()-1), 0.5*(mask.cols()-1), math::rtod(ang), 0.5*math::rtod(dang), 1);
1051 wfp::extractMaskedPixels(m_opdMask, opdMaskQ, mask);
1052 }
1053 opdMaskQ.set(std::complex<realT>(0, 1));
1054 wfp::tiltWavefront(opdMaskQ, dx, dy);
1055 mask.setZero();
1056 improc::maskWedge(mask, 0.5*(mask.rows()-1), 0.5*(mask.cols()-1), math::rtod(ang), 0.5*math::rtod(dang), 1);
1057 wfp::extractMaskedPixels(m_opdMask, opdMaskQ, mask);
1058 }
1059
1060 int xsz = 2*std::max( {fabs(maxx), fabs(minx)} ) + 2*std::max({(pupilRad/2),((realT)m_pupilSz/2)});
1061 int ysz = 2*std::max( {fabs(maxy), fabs(miny)} ) + 2*std::max({(pupilRad/2), ((realT)m_pupilSz/2)});
1062 int xsz = 2*std::max( {fabs(maxx), fabs(minx)} ) + 2*std::max({(pupilRad/2),((realT)m_pupilSz/2)});
1063 int ysz = 2*std::max( {fabs(maxy), fabs(miny)} ) + 2*std::max({(pupilRad/2), ((realT)m_pupilSz/2)});
1064
1065 if(m_imageSzAuto)
1066 {
1067 m_imageSz = std::max(xsz, ysz);
1068 }
1069
1070
1071 if(m_imageSz > m_wfSz)
1072 {
1073 std::string msg = "image size (m_imageSz = " + std::to_string(m_imageSz) + ") ";
1074 msg += "> wavefront size (m_wfSz = " + std::to_string(m_wfSz) + "). ";
1075 msg += "Decrease number of sides (m_nSides = " + std::to_string(m_nSides) + ") or increase wavefront size. ";
1076 mxThrowException(mx::err::invalidconfig, "pyramidSensor::makeOpdMask", msg);
1077 if(m_imageSzAuto)
1078 {
1079 m_imageSz = std::max(xsz, ysz);
1080 }
1081
1082
1083 if(m_imageSz > m_wfSz)
1084 {
1085 std::string msg = "image size (m_imageSz = " + std::to_string(m_imageSz) + ") ";
1086 msg += "> wavefront size (m_wfSz = " + std::to_string(m_wfSz) + "). ";
1087 msg += "Decrease number of sides (m_nSides = " + std::to_string(m_nSides) + ") or increase wavefront size. ";
1088 mxThrowException(mx::err::invalidconfig, "pyramidSensor::makeOpdMask", msg);
1089 }
1090
1091 m_wfsImage.image.resize(m_imageSz, m_imageSz);
1092
1093 if(m_detRows == 0 || m_detCols == 0)
1094 {
1095 mxThrowException(mx::err::invalidconfig, "pyramidSensor::makeOpdMask", "must set detector size");
1096 }
1097
1098 if(m_detRows > m_imageSz || m_detCols > m_imageSz)
1099 {
1100 mxThrowException(mx::err::invalidconfig, "pyramidSensor::makeOpdMask", "detector is larger than image size");
1101 }
1102
1103 if(m_detRows == 0 || m_detCols == 0)
1104 {
1105 mxThrowException(mx::err::invalidconfig, "pyramidSensor::makeOpdMask", "must set detector size");
1106 }
1107
1108 if(m_detRows > m_imageSz || m_detCols > m_imageSz)
1109 {
1110 mxThrowException(mx::err::invalidconfig, "pyramidSensor::makeOpdMask", "detector is larger than image size");
1111 }
1112
1113 m_opdMaskMade = true;
1114}
1115
1116template <typename realT, typename detectorT>
1117void pyramidSensor<realT, detectorT>::makeTilts()
1118{
1119 constexpr realT pi = math::pi<realT>();
1120
1121 if( m_modSteps == 0 )
1122 {
1123 mxThrowException( mx::err::invalidconfig,
1124 "pyramidSensor::makeTilts()",
1125 "number of modulation steps (m_modSteps) has not been set." );
1126 }
1127
1128 if(m_wfPS == 0)
1129 {
1130 mxThrowException( mx::err::invalidconfig,
1131 "pyramidSensor::makeTilts()",
1132 "wavefront platescale (m_wfPS) is 0. Must set pupilSz and D first." );
1133 }
1134
1135 if( !std::isfinite( m_wfPS ) )
1136 {
1137 mxThrowException( mx::err::invalidconfig,
1138 "pyramidSensor::makeTilts()",
1139 "wavefront platescale (m_wfPS) is infinite. Must set pupilSz and D first." );
1140 }
1141
1142 if(m_D == 0)
1143 {
1144 mxThrowException(mx::err::invalidconfig, "pyramidSensor::makeTilts()", "pupil diameter is 0. Must set D > 0 first.");
1145 }
1146
1147
1148 if(m_D == 0)
1149 {
1150 mxThrowException(mx::err::invalidconfig, "pyramidSensor::makeTilts()", "pupil diameter is 0. Must set D > 0 first.");
1151 }
1152
1153
1154 if(m_D == 0)
1155 {
1156 mxThrowException(mx::err::invalidconfig, "pyramidSensor::makeTilts()", "pupil diameter is 0. Must set D > 0 first.");
1157 }
1158
1159
1160 if(m_D == 0)
1161 {
1162 mxThrowException(mx::err::invalidconfig, "pyramidSensor::makeTilts()", "pupil diameter is 0. Must set D > 0 first.");
1163 }
1164
1165
1166 realT dang = 2 * pi / ( m_modSteps );
1167 realT dx, dy;
1168
1169 m_tilts.resize( m_modSteps );
1170
1171 std::cout << "WF Size: " << m_wfSz << "\n";
1172 std::cout << "WF PS: " << m_wfPS << "\n";
1173 std::cout << "Lambda: " << m_lambda << "\n";
1174 std::cout << "Pyr. PS: " << wfp::fftPlateScale<realT>( m_wfSz, m_wfPS, m_lambda ) * 206265. << " (mas/pix)\n";
1175 std::cout << "Mod. steps: " << m_modSteps << "\n";
1176 std::cout << "Mod rad: " << m_modRadius * ( m_lambda / m_D ) / wfp::fftPlateScale<realT>( m_wfSz, m_wfPS, m_lambda )
1177 << " (pixels)\n";
1178
1179 for( int i = 0; i < m_modSteps; ++i )
1180 {
1181 dx = m_modRadius * ( m_lambda / m_D ) / wfp::fftPlateScale<realT>( m_wfSz, m_wfPS, m_lambda ) *
1182 cos( 0.0 * dang + dang * i );
1183 dy = m_modRadius * ( m_lambda / m_D ) / wfp::fftPlateScale<realT>( m_wfSz, m_wfPS, m_lambda ) *
1184 sin( 0.0 * dang + dang * i );
1185
1186 m_tilts[i].resize( m_wfSz, m_wfSz );
1187 m_tilts[i].set( std::complex<realT>( 0, 1 ) );
1188
1189 wfp::tiltWavefront( m_tilts[i], dx, dy );
1190 }
1191
1192 m_tiltsMade = true;
1193}
1194
1195template <typename realT, typename detectorT>
1196void pyramidSensor<realT, detectorT>::allocThreadMem()
1197{
1198 if( !m_opdMaskMade )
1199 {
1200 makeOpdMask(); // Needed for m_imageSz
1201 }
1202
1203 m_pupilPlaneCF.resize( m_wfSz, m_wfSz );
1204
1205 int maxTh = omp_get_max_threads();
1206 m_th_tiltedPlane.resize( maxTh );
1207
1208 m_th_focalPlane.resize( maxTh );
1209
1210 m_th_focalImage.resize( maxTh );
1211
1212 m_th_sensorPlane.resize( maxTh );
1213
1214 m_th_sensorImage.resize( maxTh );
1215
1216 for( int nTh = 0; nTh < maxTh; ++nTh )
1217 {
1218 m_th_tiltedPlane[nTh].resize( m_wfSz, m_wfSz );
1219
1220 m_th_focalPlane[nTh].resize( m_wfSz, m_wfSz );
1221
1222 m_th_focalImage[nTh].resize( m_wfSz, m_wfSz );
1223
1224 m_th_sensorPlane[nTh].resize( m_wfSz, m_wfSz );
1225
1226 m_th_sensorImage[nTh].resize( m_imageSz, m_imageSz );
1227 }
1228
1229 m_preAllocated = true;
1230}
1231
1232template <typename realT, typename detectorT>
1233void pyramidSensor<realT, detectorT>::preAllocate()
1234{
1235 if( !m_tiltsMade )
1236 {
1237 makeTilts();
1238 }
1239
1240 if( !m_opdMaskMade )
1241 {
1242 makeOpdMask();
1243 }
1244
1245 if( !m_preAllocated )
1246 {
1247 allocThreadMem();
1248 }
1249}
1250
1251template <typename realT, typename detectorT>
1252void pyramidSensor<realT, detectorT>::doSenseWavefront()
1253{
1254 BREAD_CRUMB;
1255
1256 wavefrontT pupilPlane;
1257
1258 /* Here make average wavefront for now */
1259 int _firstWavefront = m_lastWavefront - m_iTime;
1260 if( _firstWavefront < 0 )
1261 _firstWavefront += _wavefronts.size();
1262
1263 pupilPlane.amplitude = _wavefronts[_firstWavefront].amplitude;
1264 pupilPlane.phase = _wavefronts[_firstWavefront].phase;
1265
1266 realT avgIt = _wavefronts[_firstWavefront].iterNo;
1267
1268 BREAD_CRUMB;
1269
1270 for( int i = 0; i < m_iTime; ++i )
1271 {
1272 ++_firstWavefront;
1273 if( (size_t)_firstWavefront >= _wavefronts.size() )
1274 _firstWavefront = 0;
1275
1276 pupilPlane.amplitude += _wavefronts[_firstWavefront].amplitude;
1277 pupilPlane.phase += _wavefronts[_firstWavefront].phase;
1278 avgIt += _wavefronts[_firstWavefront].iterNo;
1279 }
1280
1281 BREAD_CRUMB;
1282
1283 pupilPlane.amplitude /= ( m_iTime + 1 );
1284 pupilPlane.phase /= ( m_iTime + 1 );
1285
1286 avgIt /= ( m_iTime + 1.0 );
1287
1288 /*=====================================*/
1289 doSenseWavefront( pupilPlane );
1290
1291 m_wfsImage.iterNo = avgIt;
1292}
1293
1294template <typename realT, typename detectorT>
1295void pyramidSensor<realT, detectorT>::doSenseWavefront( wavefrontT &pupilPlane )
1296{
1297 if( m_modRadius == 0 )
1298 {
1299 return doSenseWavefrontNoMod( pupilPlane );
1300 }
1301
1302 if( !m_preAllocated )
1303 {
1304 preAllocate();
1305 }
1306
1307 BREAD_CRUMB;
1308
1309 m_wfsImage.image.resize( m_imageSz, m_imageSz );
1310 m_wfsImage.image.setZero();
1311
1312 wfsTipImage.image.resize( m_wfSz, m_wfSz );
1313 wfsTipImage.image.setZero();
1314
1315 for( size_t l = 0; l < m_wavelengths.size(); ++l )
1316 {
1317 pupilPlane.lambda = m_lambda;
1318 pupilPlane.getWavefront( m_pupilPlaneCF, m_wavelengths[l], m_wfSz );
1319
1320#pragma omp parallel
1321 {
1322 int nTh = omp_get_thread_num();
1323
1324 m_th_sensorImage[nTh].setZero();
1325 m_th_focalImage[nTh].setZero();
1326
1327 complexT *tpm_Data;
1328 complexT *ppm_Data;
1329 complexT *tim_Data;
1330 complexT *fpm_Data;
1331 complexT *opdm_Data;
1332
1333 ppm_Data = m_pupilPlaneCF.data();
1334 tpm_Data = m_th_tiltedPlane[nTh].data();
1335 fpm_Data = m_th_focalPlane[nTh].data();
1336 opdm_Data = m_opdMask.data();
1337
1338 int nelem = m_wfSz * m_wfSz;
1339
1340#pragma omp for
1341 for( int i = 0; i < m_modSteps; ++i )
1342 {
1343
1344 tim_Data = m_tilts[i].data();
1345
1346 //---------------------------------------------
1347 // Apply the modulating tip
1348 //---------------------------------------------
1349 for( int ii = 0; ii < nelem; ++ii )
1350 {
1351 tpm_Data[ii] = ppm_Data[ii] * tim_Data[ii];
1352 }
1353
1354 //---------------------------------------------
1355 // Propagate to Pyramid tip
1356 //---------------------------------------------
1357 m_frProp.propagatePupilToFocal( m_th_focalPlane[nTh], m_th_tiltedPlane[nTh], true );
1358
1359 //---------------------------------------------
1360 // Extract the tip image.
1361 //---------------------------------------------
1362 wfp::extractIntensityImageAccum(
1363 m_th_focalImage[nTh], 0, m_wfSz, 0, m_wfSz, m_th_focalPlane[nTh], 0, 0 );
1364
1365 //---------------------------------------------
1366 // Now apply the pyramid OPD
1367 //---------------------------------------------
1368 for( int ii = 0; ii < nelem; ++ii )
1369 {
1370 fpm_Data[ii] = fpm_Data[ii] * opdm_Data[ii];
1371 }
1372
1373 //---------------------------------------------
1374 // Propagate to sensor plane
1375 //---------------------------------------------
1376 m_frProp.propagateFocalToPupil( m_th_sensorPlane[nTh], m_th_focalPlane[nTh], true );
1377
1378 //---------------------------------------------
1379 // Extract the image.
1380 //---------------------------------------------
1381 wfp::extractIntensityImageAccum( m_th_sensorImage[nTh],
1382 0,
1383 m_imageSz,
1384 0,
1385 m_imageSz,
1386 m_th_sensorPlane[nTh],
1387 0.5 * m_wfSz - m_imageSz / 2,
1388 0.5 * m_wfSz - m_imageSz / 2 );
1389
1390 } // for
1391
1392 BREAD_CRUMB;
1393
1394#pragma omp critical
1395 {
1396 m_wfsImage.image += m_th_sensorImage[nTh] * _wavelengthWeights[l];
1397 wfsTipImage.image += m_th_focalImage[nTh] * _wavelengthWeights[l];
1398 }
1399 } // #pragma omp parallel
1400
1401 } // l for wavelength
1402
1403 BREAD_CRUMB;
1404
1405 m_wfsImage.image /= m_modSteps;
1406 wfsTipImage.image /= m_modSteps;
1407}
1408
1409template <typename realT, typename detectorT>
1410void pyramidSensor<realT, detectorT>::doSenseWavefrontNoMod( wavefrontT &pupilPlane )
1411{
1412 BREAD_CRUMB;
1413
1414 if( !m_opdMaskMade )
1415 {
1416 makeOpdMask();
1417 }
1418
1419 m_wfsImage.image.resize( m_imageSz, m_imageSz );
1420 m_wfsImage.image.setZero();
1421
1422 wfsTipImage.image.resize( m_wfSz, m_wfSz );
1423 wfsTipImage.image.setZero();
1424
1425 complexFieldT m_pupilPlaneCF;
1426
1427 pupilPlane.getWavefront( m_pupilPlaneCF, m_wfSz );
1428
1429 complexFieldT tiltedPlane;
1430 complexFieldT focalPlane;
1431 complexFieldT sensorPlane;
1432
1433 tiltedPlane.resize( m_wfSz, m_wfSz );
1434 focalPlane.resize( m_wfSz, m_wfSz );
1435 sensorPlane.resize( m_wfSz, m_wfSz );
1436
1437 int nelem = m_wfSz * m_wfSz;
1438
1439 complexT *tpm_Data = tiltedPlane.data();
1440 complexT *ppm_Data = m_pupilPlaneCF.data();
1441 complexT *opdm_Data = m_opdMask.data();
1442 complexT *fpm_Data = focalPlane.data();
1443
1444 BREAD_CRUMB;
1445
1446 //---------------------------------------------
1447 // Apply NO modulator tilt
1448 //---------------------------------------------
1449 for( int ii = 0; ii < nelem; ++ii )
1450 {
1451 tpm_Data[ii] = ppm_Data[ii];
1452 }
1453
1454 BREAD_CRUMB;
1455
1456 //---------------------------------------------
1457 // Propagate to Pyramid tip
1458 //---------------------------------------------
1459 m_frProp.propagatePupilToFocal( focalPlane, tiltedPlane, true );
1460
1461 BREAD_CRUMB;
1462
1463 //---------------------------------------------
1464 // Extract the tip image.
1465 //---------------------------------------------
1466 wfp::extractIntensityImageAccum( wfsTipImage.image, 0, m_wfSz, 0, m_wfSz, focalPlane, 0, 0 );
1467
1468 BREAD_CRUMB;
1469
1470 //---------------------------------------------
1471 // Now apply the pyramid OPD
1472 //---------------------------------------------
1473 for( int ii = 0; ii < nelem; ++ii )
1474 {
1475 fpm_Data[ii] = fpm_Data[ii] * opdm_Data[ii];
1476 }
1477
1478 BREAD_CRUMB;
1479
1480 //---------------------------------------------
1481 // Propagate to sensor plane
1482 //---------------------------------------------
1483 m_frProp.propagateFocalToPupil( sensorPlane, focalPlane, true );
1484
1485 BREAD_CRUMB;
1486
1487 //---------------------------------------------
1488 // Extract the image.
1489 //---------------------------------------------
1490 wfp::extractIntensityImageAccum( m_wfsImage.image,
1491 0,
1492 m_imageSz,
1493 0,
1494 m_imageSz,
1495 sensorPlane,
1496 0.5 * m_wfSz - m_imageSz / 2,
1497 0.5 * m_wfSz - m_imageSz / 2 );
1498
1499 BREAD_CRUMB;
1500}
1501
1502} // namespace sim
1503} // namespace AO
1504} // namespace mx
1505
1506#endif // mx_AO_sim_pyramidSensor_hpp
mxException for invalid config settings
Eigen::Array< scalarT, -1, -1 > eigenImage
Definition of the eigenImage type, which is an alias for Eigen::Array.
constexpr T pi()
Get the value of pi.
Definition constants.hpp:55
constexpr floatT six_fifths()
Return 6/5 in the specified precision.
realT rtod(realT q)
Convert from radians to degrees.
Definition geo.hpp:153
void maskWedge(arrayT &m, typename arrayT::Scalar xcen, typename arrayT::Scalar ycen, typename arrayT::Scalar angCen, typename arrayT::Scalar angHW, typename arrayT::Scalar val=0)
Mask a wedge in an image.
void tiltWavefront(wavefrontT &complexWavefront, typename wavefrontT::Scalar::value_type xTilt, typename wavefrontT::Scalar::value_type yTilt)
Apply a tilt to a wavefront.
The mxlib c++ namespace.
Definition mxError.hpp:106