Loading [MathJax]/extensions/tex2jax.js
mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Modules Pages
turbSubHarmonic.hpp
1/** \file turbSubHarm.hpp
2 * \brief A class to manage low-frequency sub-harmonic phase screen generation in atmospheric turbulence.
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_turbSubHarm_hpp
30#define mx_AO_sim_turbSubHarm_hpp
31
32#include <vector>
33
34#include "turbAtmosphere.hpp"
35
36#include "../../base/changeable.hpp"
37
38#include "../../improc/eigenImage.hpp"
39#include "../../improc/eigenCube.hpp"
40
41#include "../../math/constants.hpp"
42#include "../../math/func/jinc.hpp"
43
44namespace mx
45{
46namespace AO
47{
48namespace sim
49{
50
51/// A class to manage low-frequency sub-harmonic phase screen generation in atmospheric turbulence.
52/** Implements the method of Johansson and Gavel (1994)\cite johansson_gavel_1994. This is
53 * needed to generate adequate tip/tilt variance with large outer scales.
54 *
55 * \ingroup mxAOSim
56 */
57template <typename _turbAtmosphereT>
58class turbSubHarmonic : public base::changeable<turbSubHarmonic<_turbAtmosphereT>>
59{
60
61 public:
62 typedef _turbAtmosphereT turbAtmosphereT;
63 typedef typename turbAtmosphereT::realT realT;
64
65 protected:
66 /** \name Configuration Parameters
67 * @{
68 */
69
70 turbAtmosphereT *m_turbAtmo{ nullptr }; ///< Pointer to the parent atmosphere object.
71
72 unsigned m_level{ 1 }; ///< The subharmonic level to apply.
73
74 bool m_outerSubHarmonics{ true }; ///< Whether or not to include the outer subharmonics
75
76 bool m_preCalc{ false }; ///< whether or not the modes are pre-calculated.
77
78 ///@}
79
80 /** \name Internal State
81 * @{
82 */
83
84 uint32_t m_scrnSz; ///< The wavefront screen size from the layer being simulated.
85
86 std::vector<realT> m_noise; ///< Vector of Gaussian deviates prepared for each screen generation.
87
88 std::vector<realT> m_m; ///< m-coordinate fractional spatial frequency indices of the subharmonics
89 std::vector<realT> m_n; ///< n-coordinate fractional spatial frequency indices of the subharmonics
90
91 std::vector<realT> m_sqrtPSD; // the square-root of the PSD at each point
92
93 improc::eigenCube<realT> m_modes; ///< the pre-calculated modes
94
95 ///@}
96
97 public:
98 /** \name Construction
99 * @{
100 */
102
103 ///@}
104
105 /** \name Configuration
106 * @{
107 */
108
109 /// Set the pointer to the turbulent atmosphere
110 void turbAtmo( turbAtmosphereT *atm /**< [in] the new pointer to an AO atmosphere*/ );
111
112 /// Get the pointer to the AO atmosphere
113 /**
114 * \returns the the current pointer to the AO atmosphere
115 */
116 turbAtmosphereT *turbAtmo();
117
118 /// Set the subharmonic level to apply
119 void level( uint32_t ml /**< [in] the new level*/ );
120
121 /// Get the subharmonic level
122 /**
123 * \returns the current subharmonic level
124 */
125 uint32_t level();
126
127 /// Set whether or not the outer subharmonics are included
128 void outerSubHarmonics( bool osh /**< [in] the new value of the \ref m_outerSubHarmonics flag */ );
129
130 /// Get whether or not the outer subharmonics are included
131 /**
132 * \returns the current value of the \ref m_outerSubHarmonics flag
133 */
134 bool outerSubHarmonics();
135
136 /// Set whether or not to pre-calculate the modes
137 void preCalc( bool pc /**< [in] the new value of the \ref m_preCalc flag*/ );
138
139 /// Get whether or not the modes are pre-calculated
140 /**
141 * \returns the current value of the preCalc flag
142 */
143 bool preCalc();
144
145 ///@}
146
147 /** \name Screen Generation
148 * @{
149 */
150
151 /// Allocate needed memory and initialize the subharmonic transform
152 void initGrid( uint32_t layerNo );
153
154 /// Generate a realization of the subharmonic phase screen and add it to the input screen
155 void screen( improc::eigenImage<realT> & scrn /**< [in] the input phase screen to which to add the low-frequency screen. Must be m_scrnSz X m_scrnsz*/ );
156
157 /// Deallocate memory
158 void deinit();
159
160 ///@}
161};
162
163template <typename turbAtmosphereT>
165{
166}
167
168template <typename turbAtmosphereT>
169void turbSubHarmonic<turbAtmosphereT>::turbAtmo( turbAtmosphereT *turbatm )
170{
171 if( turbatm != m_turbAtmo )
172 {
173 m_turbAtmo = turbatm;
174 this->changed();
175 }
176}
177
178template <typename turbAtmosphereT>
180{
181 return m_turbAtmo;
182}
183
184template <typename turbAtmosphereT>
186{
187 if( ml != m_level )
188 {
189 m_level = ml;
190 this->changed();
191 }
192}
193
194template <typename turbAtmosphereT>
196{
197 return m_level;
198}
199
200template <typename turbAtmosphereT>
202{
203 if( osh != m_outerSubHarmonics )
204 {
205 m_outerSubHarmonics = osh;
206 this->changed();
207 }
208}
209
210template <typename turbAtmosphereT>
212{
213 return m_outerSubHarmonics;
214}
215
216template <typename turbAtmosphereT>
218{
219 if( pc != m_preCalc )
220 {
221 m_preCalc = pc;
222 this->changed();
223 }
224}
225
226template <typename turbAtmosphereT>
228{
229 return m_preCalc;
230}
231
232template <typename turbAtmosphereT>
234{
235 int N;
236
237 if( m_turbAtmo == nullptr )
238 {
239 mxThrowException( err::paramnotset,
240 "mx::AO::sim::turbSubHarmonic::initGrid",
241 "atmosphere is not set (m_turbAtmo is nullptr)" );
242 }
243
244 if( m_turbAtmo->aosys() == nullptr )
245 {
246 mxThrowException( err::paramnotset,
247 "mx::AO::sim::turbSubHarmonic::initGrid",
248 "ao system is not set (m_turbAtmo->m_aosys is nullptr)" );
249 }
250
251 if( m_turbAtmo->nLayers() <= layerNo )
252 {
253 mxThrowException( err::invalidconfig,
254 "mx::AO::sim::turbSubHarmonic::initGrid",
255 "atmosphere is not setup (m_turbAtmo->m_layers size is <= layerNo)" );
256 }
257
258 if( m_level == 0 )
259 {
260 N = 0;
261 }
262 else
263 {
264 N = 36.0 + 32.0 * ( m_level - 1 );
265
266 if( m_outerSubHarmonics )
267 {
268 N += 20;
269 }
270 }
271
272 m_m.resize( 0 ); // resized dynamically below
273 m_n.resize( 0 );
274 m_sqrtPSD.resize( N );
275 m_noise.resize( N );
276
277 m_scrnSz = m_turbAtmo->layer( layerNo ).scrnSz();
278
279 if( N == 0 )
280 {
281
282 m_modes.clear();
283 this->setChangePoint();
284 return;
285 }
286
287 realT r0 = m_turbAtmo->aosys()->atm.r_0( m_turbAtmo->aosys()->lam_sci() );
288 realT D = m_turbAtmo->aosys()->D();
289 uint32_t wfSz = m_turbAtmo->wfSz();
290
291 realT beta = 0.0218 / pow( r0, 5. / 3. ) / pow( ( D / wfSz ) * m_scrnSz, 2 );
292
293 realT sqrt_alpha = 0.5 * 11. / 3.;
294
295 realT L0 = m_turbAtmo->aosys()->atm.L_0( layerNo );
296
297 realT L02;
298 if( L0 > 0 )
299 L02 = 1.0 / pow( L0, 2 );
300 else
301 L02 = 0;
302
303 if( m_preCalc )
304 {
305 m_modes.resize( m_scrnSz, m_scrnSz, N );
306 }
307
308 std::vector<realT> scs;
309
310 if( m_outerSubHarmonics )
311 {
312 realT sc = 0.25;
313
314 m_m = std::vector<realT>( { -1.25, -0.75, -0.25, 0.25, 0.75, 1.25, -1.25, 1.25, -1.25, 1.25,
315 -1.25, 1.25, -1.25, 1.25, -1.25, -0.75, -0.25, 0.25, 0.75, 1.25 } );
316 m_n = std::vector<realT>( { -1.25, -1.25, -1.25, -1.25, -1.25, -1.25, -0.75, -0.75, -0.25, -0.25,
317 0.25, 0.25, 0.75, 0.75, 1.25, 1.25, 1.25, 1.25, 1.25, 1.25 } );
318 scs.resize( m_m.size(), sc );
319 }
320
321 for( int nl = 1; nl <= m_level; ++nl )
322 {
323 realT sc = pow( 3.0, -nl );
324
325 for( int mp = -3; mp < 3; ++mp )
326 {
327 for( int np = -3; np < 3; ++np )
328 {
329 if( nl < m_level )
330 {
331 if( mp == -1 )
332 {
333 if( np == -1 || np == 0 )
334 continue;
335 }
336 else if( mp == 0 )
337 {
338 if( np == -1 || np == 0 )
339 continue;
340 }
341 }
342
343 m_m.push_back( sc * ( mp + 0.5 ) );
344 m_n.push_back( sc * ( np + 0.5 ) );
345 scs.push_back( sc );
346 }
347 }
348 }
349
350 int n = 0;
351 for( int n = 0; n < m_m.size(); ++n )
352 {
353 realT k = sqrt( ( pow( m_m[n], 2 ) + pow( m_n[n], 2 ) ) ) / ( ( D / wfSz ) * m_scrnSz );
354
355 realT tpsd = beta / pow( k * k + L02, sqrt_alpha );
356
357 realT Ppiston = 0;
358 realT Ptiptilt = 0;
359 if( m_turbAtmo->aosys()->psd.subPiston() )
360 {
361 Ppiston = pow( 2 * math::func::jinc( math::pi<realT>() * k * D ), 2 );
362 }
363
364 if( m_turbAtmo->aosys()->psd.subTipTilt() )
365 {
366 Ptiptilt = pow( 4 * math::func::jincN( 2, math::pi<realT>() * k * D ), 2 );
367 }
368
369 m_sqrtPSD[n] = scs[n] * sqrt( tpsd * ( 1 - Ppiston - Ptiptilt ) );
370
371 if( m_preCalc )
372 {
373 for( int cc = 0; cc < m_scrnSz; ++cc )
374 {
375 realT np = cc - 0.5 * m_scrnSz;
376 for( int rr = 0; rr < m_scrnSz; ++rr )
377 {
378 realT mp = rr - 0.5 * m_scrnSz;
379 m_modes.image( n )( rr, cc ) =
380 m_sqrtPSD[n] * cos( math::two_pi<realT>() * ( m_m[n] * mp + m_n[n] * np ) / m_scrnSz );
381 }
382 }
383 }
384 }
385
386 this->setChangePoint();
387}
388
389template <typename turbAtmosphereT>
391{
392 if( m_level == 0 )
393 {
394 return;
395 }
396
397 if( m_turbAtmo == nullptr )
398 {
399 mxThrowException(
400 err::paramnotset, "mx::AO::sim::turbSubHarmonic::screen", "atmosphere is not set (m_turbAtmo is nullptr)" );
401 }
402
403 if( this->isChanged() )
404 {
405 mxThrowException( err::invalidconfig,
406 "mx::AO::sim::turbSubHarmonic::screen",
407 "configuration has changed but not re-initialized" );
408 }
409
410 if( scrn.rows() != m_scrnSz || scrn.cols() != m_scrnSz )
411 {
412 mxThrowException( err::sizeerr, "mx::AO::sim::turbSubHarmonic::screen", "input screen is not the right size" );
413 }
414
415 // Check that we're allocated
416 if( m_preCalc )
417 {
418 if( m_modes.rows() != scrn.rows() || m_modes.cols() != scrn.cols() || m_modes.planes() != m_noise.size() )
419 {
420 mxThrowException(
421 err::sizeerr, "mx::AO::sim::turbSubHarmonic::screen", "modes cube wrong size, call initGrid()." );
422 }
423 }
424 else
425 {
426 if( m_noise.size() != m_m.size() || m_noise.size() != m_n.size() || m_sqrtPSD.size() != m_noise.size() )
427 {
428 mxThrowException(
429 err::sizeerr, "mx::AO::sim::turbSubHarmonic::screen", "vectors not allocated, call initGrid()." );
430 }
431 }
432
433 // Now fill in the noise
434 for( size_t n = 0; n < m_noise.size(); ++n )
435 {
436 m_noise[n] = 2 * m_turbAtmo->normVar();
437 }
438
439#pragma omp parallel for
440 for( int cc = 0; cc < m_scrnSz; ++cc )
441 {
442 realT np = cc - 0.5 * m_scrnSz;
443 for( int rr = 0; rr < m_scrnSz; ++rr )
444 {
445 realT mp = rr - 0.5 * m_scrnSz;
446
447 if( m_preCalc )
448 {
449 for( unsigned n = 0; n < m_noise.size(); ++n )
450 {
451 scrn( rr, cc ) += m_noise[n] * m_modes.image( n )( rr, cc );
452 }
453 }
454 else
455 {
456 for( unsigned n = 0; n < m_m.size(); ++n )
457 {
458 scrn( rr, cc ) += m_noise[n] * m_sqrtPSD[n] *
459 cos( math::two_pi<realT>() * ( m_m[n] * mp + m_n[n] * np ) / m_scrnSz );
460 }
461 }
462 }
463 }
464}
465
466template <typename turbAtmosphereT>
468{
469 m_m.clear();
470 m_n.clear();
471 m_sqrtPSD.clear();
472 m_noise.clear();
473 m_modes.resize( 0, 0, 0 );
474
475 this->changed();
476}
477
478} // namespace sim
479} // namespace AO
480} // namespace mx
481
482#endif // mx_AO_sim_turbSubHarm_hpp
A class to manage low-frequency sub-harmonic phase screen generation in atmospheric turbulence.
void deinit()
Deallocate memory.
std::vector< realT > m_noise
Vector of Gaussian deviates prepared for each screen generation.
void screen(improc::eigenImage< realT > &scrn)
Generate a realization of the subharmonic phase screen and add it to the input screen.
void initGrid(uint32_t layerNo)
Allocate needed memory and initialize the subharmonic transform.
bool m_outerSubHarmonics
Whether or not to include the outer subharmonics.
uint32_t m_scrnSz
The wavefront screen size from the layer being simulated.
bool preCalc()
Get whether or not the modes are pre-calculated.
improc::eigenCube< realT > m_modes
the pre-calculated modes
turbAtmosphereT * m_turbAtmo
Pointer to the parent atmosphere object.
std::vector< realT > m_m
m-coordinate fractional spatial frequency indices of the subharmonics
turbAtmosphereT * turbAtmo()
Get the pointer to the AO atmosphere.
uint32_t level()
Get the subharmonic level.
bool outerSubHarmonics()
Get whether or not the outer subharmonics are included.
std::vector< realT > m_n
n-coordinate fractional spatial frequency indices of the subharmonics
bool m_preCalc
whether or not the modes are pre-calculated.
unsigned m_level
The subharmonic level to apply.
A simple class to track member data changes.
mxException for invalid config settings
mxException for parameters which aren't set
mxException for a size error
An image cube with an Eigen-like API.
Definition eigenCube.hpp:30
Eigen::Array< scalarT, -1, -1 > eigenImage
Definition of the eigenImage type, which is an alias for Eigen::Array.
T2 jincN(const T1 &v, const T2 &x)
The JincN function.
Definition jinc.hpp:112
T jinc(const T &x)
The Jinc function.
Definition jinc.hpp:61
constexpr floatT six_fifths()
Return 6/5 in the specified precision.
The mxlib c++ namespace.
Definition mxError.hpp:106
Declaration and definition of a turbulent atmosphere.