mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
Loading...
Searching...
No Matches
turbAtmosphere.hpp
Go to the documentation of this file.
1/** \file turbAtmosphere.hpp
2 * \brief Declaration and definition of a turbulent atmosphere.
3 *
4 * \author Jared R. Males (jaredmales@gmail.com)
5 *
6 * \ingroup mxAO_sim_files
7 *
8 */
9
10//***********************************************************************//
11// Copyright 2023 Jared R. Males (jaredmales@gmail.com)
12//
13// This file is part of mxlib.
14//
15// mxlib is free software: you can redistribute it and/or modify
16// it under the terms of the GNU General Public License as published by
17// the Free Software Foundation, either version 3 of the License, or
18// (at your option) any later version.
19//
20// mxlib is distributed in the hope that it will be useful,
21// but WITHOUT ANY WARRANTY; without even the implied warranty of
22// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
23// GNU General Public License for more details.
24//
25// You should have received a copy of the GNU General Public License
26// along with mxlib. If not, see <http://www.gnu.org/licenses/>.
27//***********************************************************************//
28
29#ifndef mx_AO_sim_turbAtmosphere_hpp
30#define mx_AO_sim_turbAtmosphere_hpp
31
32#include <vector>
33#include <iostream>
34
35#include "../../sigproc/psdFilter.hpp"
36#include "../../sigproc/psdUtils.hpp"
37
38#include "../../base/changeable.hpp"
39
40#include "../../improc/milkImage.hpp"
41
42#include "../../math/constants.hpp"
43#include "../../math/func/jinc.hpp"
44#include "../../math/randomT.hpp"
45
46#include "../../ioutils/fits/fitsFile.hpp"
47#include "../../ioutils/stringUtils.hpp"
48#include "../../sys/timeUtils.hpp"
49
50#include "turbLayer.hpp"
51#include "turbSubHarmonic.hpp"
52#include "wavefront.hpp"
53
54#ifdef DEBUG
55#define BREAD_CRUMB std::cout << "DEBUG: " << __FILE__ << " " << __LINE__ << "\n";
56#else
57#define BREAD_CRUMB
58#endif
59
60namespace mx
61{
62namespace AO
63{
64namespace sim
65{
66
67/// A turbulent atmosphere simulator
68/** Generate multi-layer phase screens using PSD filtering with the FFT. Manages
69 * the shifting and combination of the phase screens (without propagation).
70 *
71 * The \ref turbSubHarmonic class provides low-frequency sub-harmonics if desired
72 *
73 * The \ref turbLayer class manages the actual interpolation.
74 *
75 * \todo add facility for switching interpolators
76 * \todo need to include wrap detection when using subharmonics (e.g. frames needs to be right)
77 *
78 * \ingroup mxAOSim
79 */
80template <typename _aoSystemT>
81struct turbAtmosphere : public base::changeable<turbAtmosphere<_aoSystemT>>
82{
83
84 public:
85 typedef _aoSystemT aoSystemT;
86 typedef typename aoSystemT::realT realT;
87 typedef Eigen::Array<realT, -1, -1> imageT;
88
89 protected:
90 /** \name Configuration Data
91 * @{
92 */
93 uint32_t m_wfSz{ 0 }; ///< Size of the wavefront in pixels.
94
95 uint32_t m_buffSz{ 0 }; ///< Buffer to apply around wavefront for interpolation.
96
97 aoSystemT *m_aosys{ nullptr }; ///< The AO system object containing all system parameters. See \ref aoSystem.
98
99 uint32_t m_shLevel{ 0 }; /**< Number of subharmonic levels. 0 means no subharmonics. Generally only
100 * 1 level is needed to obtain good Tip/Tilt results. See \ref turbSubHarmonic.
101 */
102
103 bool m_outerSubHarmonics{ true }; /**< Whether or not the outer subharmonics are included.
104 * See \ref m_turbSubHarmonic.
105 */
106 bool m_shPreCalc{ true }; /**< Whether or not to pre-calculate the subharmonic modes. This is a
107 * trade of speed vs. memory use. See \ref turbSubHarmonic.
108 */
109
110 bool m_retain{ false }; /**< Whether or not to retain working memory after screen generation. One
111 * would set to true if many screens will be calculated in monte carlo fashion.
112 * For a standard case of one set of screens being shifted by wind velocity,
113 * this should be false to minimize memory use.
114 */
115
116 bool m_forceGen{ true }; /**< Force generation of new screens if true. Note that this class does not presently
117 * do any checks of saved screens to verify they match the current configuration.
118 * Set to true with caution.
119 *
120 * \todo need check for changes in parameters to decide if saved frames are valid per layer
121 */
122
123 std::string m_dataDir; /**< Specifies a path where phase screens are stored. Leaving this unset or "" is equivalent
124 * to \ref m_forceGen` == true`.
125 */
126
127 imageT *m_pupil{ 0 }; ///< A pointer to the pupil mask.
128
129 realT m_timeStep{ 0 }; ///< Length of each iteration, in seconds.
130
131 ///@}
132
133 /** \name Internal State
134 * @{
135 */
136 std::vector<turbLayer<aoSystemT>> m_layers; ///< Vector of turbulent layers.
137
138 std::vector<turbSubHarmonic<turbAtmosphere>> m_subHarms; ///< Sub-harmonic layers
139
140 size_t m_frames{ 0 }; ///< Length of the turbulence sequence.
141
142 int m_nWf{ 0 }; ///< Number of iterations which have occurred.
143
144 math::normDistT<realT> m_normVar; ///< Normal random deviate generator. This seeded in the constructor.
145
146 public:
147 /** \name Construction and Configuration
148 * @{
149 */
150
151 /// Default c'tor
153
154 /// Setup the overall atmosphere.
155 /**
156 */
157 void setup(
158 uint32_t wfSz, ///< [in] The size of the wavefront in pixels.
159 uint32_t buffSz, ///< [in] The size of the interpolation buffer to use.
160 aoSystemT *aosys, ///< [in] Pointer to an AO System. See \ref m_aosys.
161 uint32_t shLevel ///< [in] number of subharmonic levels to use. 0 turns off subharmonics. See \ref m_shLevel.
162 );
163
164 /// Set the wavefront size
165 /**
166 * see \ref m_wfSz
167 */
168 void wfSz( uint32_t ws /**< [in] the new wavefront size */ );
169
170 /// Get the wavefront size
171 /**
172 * \returns the current value of \ref m_wfSz
173 */
174 uint32_t wfSz();
175
176 /// Set the interpolation buffer size
177 /**
178 * see \ref m_buffSz
179 */
180 void buffSz( uint32_t bs /**< [in] the new buffer size*/ );
181
182 /// Get the interpolation buffer size
183 /**
184 * \returns the current value of \ref m_buffSz
185 */
186 uint32_t buffSz();
187
188 /// Set the pointer to an AO system object
189 /**
190 * see \ref m_aosys
191 */
192 void aosys( aoSystemT *aos /**< [in] the new pointer to an AO system object*/ );
193
194 /// Get the pointer to an AO system object
195 /**
196 * \returns the current value of \ref m_aosys
197 */
198 aoSystemT *aosys();
199
200 /// Set the subharmonic level
201 /**
202 * see \ref m_shLevel
203 */
204 void shLevel( uint32_t shl /**< [in] the new subharmonic level*/ );
205
206 /// Get the subharmonic level
207 /**
208 * \returns the current value of \ref m_shLevel
209 */
210 uint32_t shLevel();
211
212 /// Set whether outer subharmonics are included
213 /**
214 * see \ref m_outerSubHarmonics
215 */
216 void outerSubHarmonics( bool osh /**< [in] the new value of flag controlling outer subharmonics */ );
217
218 /// Get whether outer subharmonics are included
219 /**
220 * \returns the current value of \ref m_outerSubHarmonics
221 */
223
224 /// Set whether subharmonic modes are pre-calculated
225 /**
226 * see \ref m_shPreCalc
227 */
228 void shPreCalc( bool shp /**< [in] the new value of flag controlling subharmonic mode precalculation */ );
229
230 /// Get whether subharmonic modes are pre-calculated
231 /**
232 * \returns the current value of \ref m_shPreCalc
233 */
234 bool shPreCalc();
235
236 /// Set whether memory for screen generation is retained
237 /**
238 * see \ref m_retain
239 */
240 void retain( bool rtn /**< [in] the new value of the flag controlling whether memory is retained*/ );
241
242 /// Get whether memory for screen generation is retained
243 /**
244 * \returns the current value of \ref m_retain
245 */
246 bool retain();
247
248 /// Set whether new screen generation is forced
249 /**
250 * see \ref m_forceGen
251 */
252 void forceGen( bool fg /**< [in] the new value of the flag controlling whether screen generation is forced*/ );
253
254 /// Get whether new screen generation is forced
255 /**
256 * \returns the current value of m_forceGen
257 */
258 bool forceGen();
259
260 /// Set the data directory for saving phase screens
261 /**
262 * see \ref m_dataDir
263 */
264 void dataDir( const std::string &dd /**< [in] the new data directory for phase screen saving*/ );
265
266 /// Get the data directory for saving phase screens
267 /**
268 * \returns the current value of m_dataDir
269 */
270 std::string dataDir();
271
272 /// Setup the layers and prepare them for phase screen generation.
273 /** The number of entries in \p scrnSz must mach the number of layers in the atmosphere
274 * of \ref m_aosys.
275 */
276 void setLayers( const std::vector<size_t> &scrnSz /**< [in] the screen sizes for each layer.*/ );
277
278 /// Setup the layers and prepare them for phase screen generation.
279 /** This sets all layer phase screens to be the same size.
280 *
281 */
282 void setLayers( const size_t scrnSz /**< [in] the screen sizes for all layers.*/ );
283
284 /// Get the number of layers
285 /**
286 * \returns the size of the layers vector \ref m_layers
287 */
288 uint32_t nLayers();
289
290 /// Get a layer
291 /**
292 * \returns a reference to on entry in \ref m_layers
293 */
295
296 /// Get the random number generator
297 /**
298 * \returns a reference to \ref m_normVar
299 */
301
302 ///@} - Construction and Configuration
303
304 uint32_t scrnLengthPixels();
305
306 uint32_t scrnLengthPixels( uint32_t n );
307
308 realT scrnLength();
309
310 realT scrnLength( uint32_t n );
311
312 uint32_t maxShift( realT dt );
313
314 uint32_t maxShift( uint32_t n, realT dt );
315
316 /** \name Screen Generation
317 * @{
318 */
319
320 /// Generate all phase screens
321 /** Loads them from disk if possible and \ref m_forceGen == false.
322 *
323 * Deallocates working memory when done, unless \ref m_retain == true.
324 */
325 void genLayers();
326
327 /// @}
328
329 int shift( improc::milkImage<realT> &phase, realT dt );
330
331 int frames( int f );
332 size_t frames();
333
334 int wfPS( realT ps ); ///< dummy function for simulatedAOSystem. Does nothing.
335
336 bool _loopClosed;
337
338 void nextWF( wavefront<realT> &wf );
339};
340
341template <typename aoSystemT>
343{
344 m_normVar.seed();
345}
346
347template <typename aoSystemT>
348void turbAtmosphere<aoSystemT>::setup( uint32_t ws, uint32_t bs, aoSystemT *aos, uint32_t shl )
349{
350 wfSz( ws );
351 buffSz( bs );
352 aosys( aos );
353 shLevel( shl );
354}
355
356template <typename aoSystemT>
358{
359 if( ws != m_wfSz )
360 {
361 m_wfSz = ws;
362 this->changed();
363 }
364}
365
366template <typename aoSystemT>
368{
369 return m_wfSz;
370}
371
372template <typename aoSystemT>
374{
375 if( bs != m_buffSz )
376 {
377 m_buffSz = bs;
378 this->changed();
379 }
380}
381
382template <typename aoSystemT>
384{
385 return m_buffSz;
386}
387
388template <typename aoSystemT>
390{
391 if( aos != m_aosys )
392 {
393 m_aosys = aos;
394 this->changed();
395 }
396}
397
398template <typename aoSystemT>
400{
401 return m_aosys;
402}
403
404template <typename aoSystemT>
406{
407 if( shl != m_shLevel )
408 {
409 m_shLevel = shl;
410 this->changed();
411 }
412}
413
414template <typename aoSystemT>
416{
417 return m_shLevel;
418}
419
420template <typename aoSystemT>
422{
423 if( osh != m_outerSubHarmonics )
424 {
425 m_outerSubHarmonics = osh;
426 this->changed();
427 }
428}
429
430template <typename aoSystemT>
432{
433 return m_outerSubHarmonics;
434}
435
436template <typename aoSystemT>
438{
439 if( shp != m_shPreCalc )
440 {
441 m_shPreCalc = shp;
442 this->changed();
443 }
444}
445
446template <typename aoSystemT>
448{
449 return m_shPreCalc;
450}
451
452template <typename aoSystemT>
454{
455 if( rtn != m_retain )
456 {
457 m_retain = rtn;
458 this->changed();
459 }
460}
461
462template <typename aoSystemT>
464{
465 return m_retain;
466}
467
468template <typename aoSystemT>
470{
471 if( fg != m_forceGen )
472 {
473 this->changed();
474 }
475
476 m_forceGen = fg;
477}
478
479template <typename aoSystemT>
481{
482 return m_forceGen;
483}
484
485template <typename aoSystemT>
486void turbAtmosphere<aoSystemT>::dataDir( const std::string &dd )
487{
488 if( dd != m_dataDir )
489 {
490 this->changed();
491 }
492
493 m_dataDir = dd;
494}
495
496template <typename aoSystemT>
498{
499 return m_dataDir;
500}
501
502template <typename aoSystemT>
503void turbAtmosphere<aoSystemT>::setLayers( const std::vector<size_t> &scrnSz )
504{
505 if( m_aosys == nullptr )
506 {
507 mxThrowException(
508 err::paramnotset, "mx::AO::sim::turbAtmosphere::setLayers", "ao system is not set (m_aosys is nullptr)" );
509 }
510
511 size_t nLayers = scrnSz.size();
512
513 if( nLayers != m_aosys->atm.n_layers() )
514 {
515 mxThrowException( err::invalidarg,
516 "mx::AO::sim::turbAtmosphere::setLayers",
517 "Size of scrnSz vector does not match atmosphere." );
518 }
519
520 if( nLayers != m_layers.size() )
521 {
522 this->changed();
523 }
524
525 m_layers.resize( nLayers );
526
527 for( size_t i = 0; i < nLayers; ++i )
528 {
529 m_layers[i].setLayer( this, i, scrnSz[i] );
530 }
531
532 if( m_shLevel > 0 )
533 {
534 if( nLayers != m_subHarms.size() )
535 {
536 this->changed();
537 }
538
539 m_subHarms.resize( nLayers );
540
541 for( size_t i = 0; i < nLayers; ++i )
542 {
543 m_subHarms[i].turbAtmo( this );
544 m_subHarms[i].preCalc( m_shPreCalc );
545 m_subHarms[i].level( m_shLevel );
546 m_subHarms[i].outerSubHarmonics( m_outerSubHarmonics );
547 m_subHarms[i].initGrid( i );
548 }
549 }
550}
551
552template <typename aoSystemT>
553void turbAtmosphere<aoSystemT>::setLayers( const size_t scrnSz )
554{
555 if( m_aosys == nullptr )
556 {
557 mxThrowException(
558 err::paramnotset, "mx::AO::sim::turbAtmosphere::setLayers", "atmosphere is not set (m_atm is nullptr)" );
559 }
560
561 size_t n = m_aosys->atm.n_layers();
562
563 setLayers( std::vector<size_t>( n, scrnSz ) );
564}
565
566template <typename aoSystemT>
568{
569 return m_layers.size();
570}
571
572template <typename aoSystemT>
574{
575 if( n >= m_layers.size() )
576 {
577 mxThrowException( err::invalidarg, "mx::AO::sim::turbAtmosphere::layer", "n too large for number of layers." );
578 }
579
580 return m_layers[n];
581}
582
583template <typename aoSystemT>
588
589template <typename aoSystemT>
591{
592 return 0;
593}
594
595template <typename aoSystemT>
596uint32_t turbAtmosphere<aoSystemT>::scrnLengthPixels( uint32_t n )
597{
598
599 return 0;
600 /* realT dx = (scrnSz-wfSz-turb.buffSz()) / cos(m_aosys->atm.layer_v_wind(n));
601
602 pow(2*pow(scrnSz-wfSz-turb.buffSz(),2), 0.5) - 1;
603 */
604}
605
606template <typename aoSystemT>
607realT turbAtmosphere<aoSystemT>::scrnLength()
608{
609 return 0;
610}
611
612template <typename aoSystemT>
613realT turbAtmosphere<aoSystemT>::scrnLength( uint32_t n )
614{
615 return 0;
616}
617
618template <typename aoSystemT>
619uint32_t turbAtmosphere<aoSystemT>::maxShift( realT dt )
620{
621 return 0;
622}
623
624template <typename aoSystemT>
625uint32_t turbAtmosphere<aoSystemT>::maxShift( uint32_t n, realT dt )
626{
627 return 0;
628}
629
630template <typename aoSystemT>
632{
633 if( m_dataDir != "" && !m_forceGen )
634 {
635 std::string fbase = m_dataDir;
636 fbase += "/";
637 fbase += "layer_";
638
640 std::string fname;
641 for( size_t i = 0; i < m_layers.size(); ++i )
642 {
643 fname = fbase + ioutils::convertToString<int>( i ) + ".fits";
644 ff.read( m_layers[i].m_phase, fname );
645 }
646
647 this->setChangePoint();
648
649 return;
650 }
651
652 sigproc::psdFilter<realT, 2> filt;
653
654 for( size_t i = 0; i < m_layers.size(); ++i )
655 {
656 if( m_layers[i].m_phase.rows() != m_layers[i].scrnSz() || m_layers[i].m_phase.cols() != m_layers[i].scrnSz() )
657 {
658 mxThrowException( err::sizeerr, "mx::AO::sim::turbAtmosphere::genLayers", "layer phase not allocated." );
659 }
660
661 if( m_layers[i].m_freq.rows() != m_layers[i].scrnSz() || m_layers[i].m_freq.cols() != m_layers[i].scrnSz() )
662 {
663 mxThrowException( err::sizeerr, "mx::AO::sim::turbAtmosphere::genLayers", "layer freq not allocated." );
664 }
665
666 if( m_layers[i].m_psd.rows() != m_layers[i].scrnSz() || m_layers[i].m_psd.cols() != m_layers[i].scrnSz() )
667 {
668 mxThrowException( err::sizeerr, "mx::AO::sim::turbAtmosphere::genLayers", "layer psd not allocated." );
669 }
670
671 for( size_t jj = 0; jj < m_layers[i].m_scrnSz; ++jj )
672 {
673 for( size_t ii = 0; ii < m_layers[i].m_scrnSz; ++ii )
674 {
675 m_layers[i].m_phase( ii, jj ) = m_normVar;
676 }
677 }
678
679 realT dkx = m_layers[i].m_freq( 1, 0 ) - m_layers[i].m_freq( 0, 0 );
680 realT dky = m_layers[i].m_freq( 0, 1 ) - m_layers[i].m_freq( 0, 0 );
681 filt.psdSqrt( m_layers[i].m_psd, dkx, dky );
682
683 filt( m_layers[i].m_phase );
684
685 if( m_shLevel > 0 )
686 {
687 m_subHarms[i].screen( m_layers[i].m_phase );
688 }
689 }
690
691 if( !m_retain )
692 {
693 for( size_t i = 0; i < m_layers.size(); ++i )
694 {
695 m_layers[i].genDealloc();
696 }
697
698 for( size_t i = 0; i < m_subHarms.size(); ++i )
699 {
700 m_subHarms[i].deinit();
701 }
702 }
703
704 if( m_dataDir != "" )
705 {
706 std::string fbase = m_dataDir;
707 fbase += "/";
708 fbase += "layer_";
709
711 std::string fname;
712 for( size_t i = 0; i < m_layers.size(); ++i )
713 {
714 fname = fbase + ioutils::convertToString<int>( i ) + ".fits";
715 ff.write( fname, m_layers[i].m_phase );
716 }
717 }
718
719 this->setChangePoint();
720}
721
722template <typename aoSystemT>
724{
725 if( this->isChanged() )
726 {
727 mxThrowException( err::invalidconfig,
728 "mx::AO::sim::turbAtmosphere::shift",
729 "configuration has changed but genLayers has not been run." );
730 }
731
732 improc::eigenMap<realT> phase( milkPhase );
733
734 if( phase.rows() != m_wfSz || phase.cols() != m_wfSz )
735 {
736 }
737
738 // Don't use OMP if no multiple layers b/c it seems to make it use multiple threads in some awful way
739 if( m_layers.size() > 1 )
740 {
741#pragma omp parallel for
742 for( size_t j = 0; j < m_layers.size(); ++j )
743 {
744 m_layers[j].shift( dt );
745 }
746 }
747 else
748 {
749 for( size_t j = 0; j < m_layers.size(); ++j )
750 {
751 m_layers[j].shift( dt );
752 }
753 }
754
755 milkPhase.setWrite();
756 phase.setZero();
757 for( size_t j = 0; j < m_layers.size(); ++j )
758 {
759 phase +=
760 sqrt( m_aosys->atm.layer_Cn2( j ) ) * m_layers[j].m_shiftPhase.block( m_buffSz, m_buffSz, m_wfSz, m_wfSz );
761 }
762 milkPhase.post();
763
764 return 0;
765}
766
767template <typename aoSystemT>
768int turbAtmosphere<aoSystemT>::frames( int f )
769{
770 m_frames = f;
771
772 return 0;
773}
774
775template <typename aoSystemT>
776size_t turbAtmosphere<aoSystemT>::frames()
777{
778 return m_frames;
779}
780
781template <typename aoSystemT>
783{
784 static_cast<void>( ps );
785 return 0;
786}
787
788template <typename aoSystemT>
790{
791 if( !m_aosys )
792 {
793 mxThrowException(
794 err::paramnotset, "mx::AO::sim::turbAtmosphere<aoSystemT>::nextWF", "the AO system pointer is not set" );
795 }
796
797 if( !m_pupil )
798 {
799 mxThrowException(
800 err::paramnotset, "mx::AO::sim::turbAtmosphere<aoSystemT>::nextWF", "the m_pupil pointer is not set" );
801 }
802
803 shift( wf.phase, m_nWf * m_timeStep );
804 ++m_nWf;
805
806 wf.phase *= *m_pupil;
807
808 realT pixVal = sqrt( m_aosys->Fg() ) * ( m_aosys->D() / m_wfSz );
809
810 wf.amplitude = pixVal * ( *m_pupil );
811}
812
813} // namespace sim
814} // namespace AO
815} // namespace mx
816
817#endif // mx_AO_sim_turbAtmosphere_hpp
@ phase
The phase or OPD.
Definition aoPSDs.hpp:52
A simple class to track member data changes.
mxException for invalid arguments
mxException for invalid config settings
mxException for parameters which aren't set
mxException for a size error
Class to manage interactions with a FITS file.
Definition fitsFile.hpp:52
int write(const dataT *im, int d1, int d2, int d3, fitsHeader *head)
Write the contents of a raw array to the FITS file.
int read(dataT *data)
Read the contents of the FITS file into an array.
Definition fitsFile.hpp:829
Class to interface with an ImageStreamIO image in shared memory.
void setWrite(bool wrflag=true)
Set the write flag.
void post()
Update the metadata and post all semaphores.
Eigen::Map< Eigen::Array< scalarT, -1, -1 > > eigenMap
Definition of the eigenMap type, which is an alias for Eigen::Map<Array>.
The mxlib c++ namespace.
Definition mxError.hpp:106
A turbulent atmosphere simulator.
uint32_t buffSz()
Get the interpolation buffer size.
aoSystemT * m_aosys
The AO system object containing all system parameters. See aoSystem.
uint32_t m_buffSz
Buffer to apply around wavefront for interpolation.
bool outerSubHarmonics()
Get whether outer subharmonics are included.
int wfPS(realT ps)
dummy function for simulatedAOSystem. Does nothing.
uint32_t nLayers()
Get the number of layers.
std::vector< turbLayer< aoSystemT > > m_layers
Vector of turbulent layers.
bool forceGen()
Get whether new screen generation is forced.
int m_nWf
Number of iterations which have occurred.
void retain(bool rtn)
Set whether memory for screen generation is retained.
std::vector< turbSubHarmonic< turbAtmosphere > > m_subHarms
Sub-harmonic layers.
void wfSz(uint32_t ws)
Set the wavefront size.
void genLayers()
Generate all phase screens.
std::string dataDir()
Get the data directory for saving phase screens.
aoSystemT * aosys()
Get the pointer to an AO system object.
uint32_t shLevel()
Get the subharmonic level.
imageT * m_pupil
A pointer to the pupil mask.
uint32_t wfSz()
Get the wavefront size.
bool retain()
Get whether memory for screen generation is retained.
realT m_timeStep
Length of each iteration, in seconds.
void setLayers(const std::vector< size_t > &scrnSz)
Setup the layers and prepare them for phase screen generation.
void forceGen(bool fg)
Set whether new screen generation is forced.
bool shPreCalc()
Get whether subharmonic modes are pre-calculated.
math::normDistT< realT > & normVar()
Get the random number generator.
void shPreCalc(bool shp)
Set whether subharmonic modes are pre-calculated.
void shLevel(uint32_t shl)
Set the subharmonic level.
void dataDir(const std::string &dd)
Set the data directory for saving phase screens.
size_t m_frames
Length of the turbulence sequence.
turbLayer< aoSystemT > & layer(uint32_t n)
Get a layer.
uint32_t m_wfSz
Size of the wavefront in pixels.
math::normDistT< realT > m_normVar
Normal random deviate generator. This seeded in the constructor.
void setup(uint32_t wfSz, uint32_t buffSz, aoSystemT *aosys, uint32_t shLevel)
Setup the overall atmosphere.
void buffSz(uint32_t bs)
Set the interpolation buffer size.
void setLayers(const size_t scrnSz)
Setup the layers and prepare them for phase screen generation.
void outerSubHarmonics(bool osh)
Set whether outer subharmonics are included.
void aosys(aoSystemT *aos)
Set the pointer to an AO system object.
Simulation of a single turbulent layer.
Definition turbLayer.hpp:56
Structure containing the phase and amplitude of a wavefront.
Definition wavefront.hpp:24
realImageT amplitude
The wavefront amplitude.
Definition wavefront.hpp:35
realImageT phase
The wavefront phase.
Definition wavefront.hpp:38
Declaration and definition of a turbulence layer.