29#ifndef mx_AO_sim_turbAtmosphere_hpp
30#define mx_AO_sim_turbAtmosphere_hpp
35#include "../../sigproc/psdFilter.hpp"
36#include "../../sigproc/psdUtils.hpp"
38#include "../../base/changeable.hpp"
40#include "../../improc/milkImage.hpp"
42#include "../../math/constants.hpp"
43#include "../../math/func/jinc.hpp"
44#include "../../math/randomT.hpp"
46#include "../../ioutils/fits/fitsFile.hpp"
47#include "../../ioutils/stringUtils.hpp"
48#include "../../sys/timeUtils.hpp"
51#include "turbSubHarmonic.hpp"
52#include "wavefront.hpp"
55#define BREAD_CRUMB std::cout << "DEBUG: " << __FILE__ << " " << __LINE__ << "\n";
80template <
typename _aoSystemT,
class _verboseT = mx::verbose::d>
85 typedef _aoSystemT aoSystemT;
86 typedef _verboseT verboseT;
88 typedef typename aoSystemT::realT realT;
89 typedef Eigen::Array<realT, -1, -1> imageT;
307 uint32_t scrnLengthPixels();
309 uint32_t scrnLengthPixels( uint32_t n );
313 realT scrnLength( uint32_t n );
315 uint32_t maxShift( realT dt );
317 uint32_t maxShift( uint32_t n, realT dt );
332 int shift( improc::milkImage<realT> &phase, realT dt );
344template <
typename aoSystemT,
class verboseT>
350template <
typename aoSystemT,
class verboseT>
359template <
typename aoSystemT,
class verboseT>
369template <
typename aoSystemT,
class verboseT>
375template <
typename aoSystemT,
class verboseT>
385template <
typename aoSystemT,
class verboseT>
391template <
typename aoSystemT,
class verboseT>
401template <
typename aoSystemT,
class verboseT>
407template <
typename aoSystemT,
class verboseT>
410 if( shl != m_shLevel )
417template <
typename aoSystemT,
class verboseT>
423template <
typename aoSystemT,
class verboseT>
426 if( osh != m_outerSubHarmonics )
428 m_outerSubHarmonics = osh;
433template <
typename aoSystemT,
class verboseT>
436 return m_outerSubHarmonics;
439template <
typename aoSystemT,
class verboseT>
442 if( shp != m_shPreCalc )
449template <
typename aoSystemT,
class verboseT>
455template <
typename aoSystemT,
class verboseT>
458 if( rtn != m_retain )
465template <
typename aoSystemT,
class verboseT>
471template <
typename aoSystemT,
class verboseT>
474 if( fg != m_forceGen )
482template <
typename aoSystemT,
class verboseT>
488template <
typename aoSystemT,
class verboseT>
491 if( dd != m_dataDir )
499template <
typename aoSystemT,
class verboseT>
505template <
typename aoSystemT,
class verboseT>
508 if( m_aosys ==
nullptr )
513 size_t nLayers = scrnSz.size();
515 if( nLayers != m_aosys->atm.n_layers() )
519 "Size of scrnSz vector does not match atmosphere." );
522 if( nLayers != m_layers.size() )
527 m_layers.resize( nLayers );
529 for(
size_t i = 0; i < nLayers; ++i )
531 m_layers[i].setLayer(
this, i, scrnSz[i] );
536 if( nLayers != m_subHarms.size() )
541 m_subHarms.resize( nLayers );
543 for(
size_t i = 0; i < nLayers; ++i )
545 m_subHarms[i].turbAtmo(
this );
546 m_subHarms[i].preCalc( m_shPreCalc );
547 m_subHarms[i].level( m_shLevel );
548 m_subHarms[i].outerSubHarmonics( m_outerSubHarmonics );
549 m_subHarms[i].initGrid( i );
554template <
typename aoSystemT,
class verboseT>
557 if( m_aosys ==
nullptr )
561 "atmosphere is not set (m_atm is nullptr)" );
564 size_t n = m_aosys->atm.n_layers();
566 setLayers( std::vector<size_t>( n, scrnSz ) );
569template <
typename aoSystemT,
class verboseT>
572 return m_layers.size();
575template <
typename aoSystemT,
class verboseT>
578 if( n >= m_layers.size() )
586template <
typename aoSystemT,
class verboseT>
592template <
typename aoSystemT,
class verboseT>
598template <
typename aoSystemT,
class verboseT>
599uint32_t turbAtmosphere<aoSystemT, verboseT>::scrnLengthPixels( uint32_t n )
609template <
typename aoSystemT,
class verboseT>
610turbAtmosphere<aoSystemT, verboseT>::realT turbAtmosphere<aoSystemT, verboseT>::scrnLength()
615template <
typename aoSystemT,
class verboseT>
616turbAtmosphere<aoSystemT, verboseT>::realT turbAtmosphere<aoSystemT, verboseT>::scrnLength( uint32_t n )
621template <
typename aoSystemT,
class verboseT>
622uint32_t turbAtmosphere<aoSystemT, verboseT>::maxShift( realT dt )
627template <
typename aoSystemT,
class verboseT>
628uint32_t turbAtmosphere<aoSystemT, verboseT>::maxShift( uint32_t n, realT dt )
633template <
typename aoSystemT,
class verboseT>
636 if( m_dataDir !=
"" && !m_forceGen )
638 std::string fbase = m_dataDir;
644 for(
size_t i = 0; i < m_layers.size(); ++i )
646 fname = fbase + ioutils::convertToString<int>( i ) +
".fits";
647 ff.
read( m_layers[i].m_phase, fname );
650 this->setChangePoint();
655 sigproc::psdFilter<realT, 2> filt;
657 for(
size_t i = 0; i < m_layers.size(); ++i )
659 if( m_layers[i].m_phase.rows() != m_layers[i].scrnSz() || m_layers[i].m_phase.cols() != m_layers[i].scrnSz() )
664 if( m_layers[i].m_freq.rows() != m_layers[i].scrnSz() || m_layers[i].m_freq.cols() != m_layers[i].scrnSz() )
669 if( m_layers[i].m_psd.rows() != m_layers[i].scrnSz() || m_layers[i].m_psd.cols() != m_layers[i].scrnSz() )
674 for(
size_t jj = 0; jj < m_layers[i].m_scrnSz; ++jj )
676 for(
size_t ii = 0; ii < m_layers[i].m_scrnSz; ++ii )
678 m_layers[i].m_phase( ii, jj ) = m_normVar;
682 realT dkx = m_layers[i].m_freq( 1, 0 ) - m_layers[i].m_freq( 0, 0 );
683 realT dky = m_layers[i].m_freq( 0, 1 ) - m_layers[i].m_freq( 0, 0 );
684 filt.psdSqrt( m_layers[i].m_psd, dkx, dky );
686 filt( m_layers[i].m_phase );
690 m_subHarms[i].screen( m_layers[i].m_phase );
696 for(
size_t i = 0; i < m_layers.size(); ++i )
698 m_layers[i].genDealloc();
701 for(
size_t i = 0; i < m_subHarms.size(); ++i )
703 m_subHarms[i].deinit();
707 if( m_dataDir !=
"" )
709 std::string fbase = m_dataDir;
715 for(
size_t i = 0; i < m_layers.size(); ++i )
717 fname = fbase + ioutils::convertToString<int>( i ) +
".fits";
718 ff.
write( fname, m_layers[i].m_phase );
722 this->setChangePoint();
725template <
typename aoSystemT,
class verboseT>
728 if( this->isChanged() )
735 if( phase.rows() != m_wfSz || phase.cols() != m_wfSz )
740 if( m_layers.size() > 1 )
743 #pragma omp parallel for
744 for(
size_t j = 0; j < m_layers.size(); ++j )
746 m_layers[j].shift( dt );
751 for(
size_t j = 0; j < m_layers.size(); ++j )
753 m_layers[j].shift( dt );
757 milkPhase.setWrite();
759 for(
size_t j = 0; j < m_layers.size(); ++j )
762 sqrt( m_aosys->atm.layer_Cn2( j ) ) * m_layers[j].m_shiftPhase.block( m_buffSz, m_buffSz, m_wfSz, m_wfSz );
769template <
typename aoSystemT,
class verboseT>
770int turbAtmosphere<aoSystemT, verboseT>::frames(
int f )
777template <
typename aoSystemT,
class verboseT>
778size_t turbAtmosphere<aoSystemT, verboseT>::frames()
783template <
typename aoSystemT,
class verboseT>
786 static_cast<void>( ps );
790template <
typename aoSystemT,
class verboseT>
803 shift( wf.
phase, m_nWf * m_timeStep );
806 wf.
phase *= *m_pupil;
808 realT pixVal = sqrt( m_aosys->Fg() ) * ( m_aosys->D() / m_wfSz );
A simple class to track member data changes.
Augments an exception with the source file and line.
Class to manage interactions with a FITS file.
error_t read(dataT *data)
Read the contents of the FITS file into an array.
error_t write(const dataT *im, int d1, int d2, int d3, fitsHeader< verboseT > *head)
Write the contents of a raw array to the FITS file.
Eigen::Map< Eigen::Array< scalarT, -1, -1 > > eigenMap
Definition of the eigenMap type, which is an alias for Eigen::Map<Array>.
@ sizeerr
A size was invalid or calculated incorrectly.
@ paramnotset
A parameter was not set.
@ invalidconfig
A config setting was invalid.
@ invalidarg
An argument was invalid.
A turbulent atmosphere simulator.
void setLayers(const size_t scrnSz)
Setup the layers and prepare them for phase screen generation.
realT m_timeStep
Length of each iteration, in seconds.
uint32_t nLayers()
Get the number of layers.
bool forceGen()
Get whether new screen generation is forced.
uint32_t buffSz()
Get the interpolation buffer size.
turbAtmosphere()
Default c'tor.
void shPreCalc(bool shp)
Set whether subharmonic modes are pre-calculated.
math::normDistT< realT > & normVar()
Get the random number generator.
std::vector< turbLayer< aoSystemT > > m_layers
Vector of turbulent layers.
uint32_t m_wfSz
Size of the wavefront in pixels.
bool retain()
Get whether memory for screen generation is retained.
std::string dataDir()
Get the data directory for saving phase screens.
turbLayer< aoSystemT > & layer(uint32_t n)
Get a layer.
uint32_t shLevel()
Get the subharmonic level.
void wfSz(uint32_t ws)
Set the wavefront size.
aoSystemT * m_aosys
The AO system object containing all system parameters. See aoSystem.
void outerSubHarmonics(bool osh)
Set whether outer subharmonics are included.
std::vector< turbSubHarmonic< turbAtmosphere > > m_subHarms
Sub-harmonic layers.
bool shPreCalc()
Get whether subharmonic modes are pre-calculated.
int m_nWf
Number of iterations which have occurred.
math::normDistT< realT > m_normVar
Normal random deviate generator. This seeded in the constructor.
bool outerSubHarmonics()
Get whether outer subharmonics are included.
uint32_t m_buffSz
Buffer to apply around wavefront for interpolation.
void aosys(aoSystemT *aos)
Set the pointer to an AO system object.
void shLevel(uint32_t shl)
Set the subharmonic level.
int wfPS(realT ps)
dummy function for simulatedAOSystem. Does nothing.
void setup(uint32_t wfSz, uint32_t buffSz, aoSystemT *aosys, uint32_t shLevel)
Setup the overall atmosphere.
aoSystemT * aosys()
Get the pointer to an AO system object.
void dataDir(const std::string &dd)
Set the data directory for saving phase screens.
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.
size_t m_frames
Length of the turbulence sequence.
imageT * m_pupil
A pointer to the pupil mask.
void retain(bool rtn)
Set whether memory for screen generation is retained.
void buffSz(uint32_t bs)
Set the interpolation buffer size.
void genLayers()
Generate all phase screens.
uint32_t wfSz()
Get the wavefront size.
Simulation of a single turbulent layer.
Structure containing the phase and amplitude of a wavefront.
realImageT amplitude
The wavefront amplitude.
realImageT phase
The wavefront phase.
Declaration and definition of a turbulence layer.