mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
Loading...
Searching...
No Matches
turbLayer.hpp
Go to the documentation of this file.
1/** \file turbLayer.hpp
2 * \brief Declaration and definition of a turbulence layer.
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_turbLayer_hpp
30#define mx_AO_sim_turbLayer_hpp
31
32#include <vector>
33
34#include <Eigen/Dense>
35
36#include "../../improc/imageTransforms.hpp"
37
38namespace mx
39{
40namespace AO
41{
42namespace sim
43{
44
45template <typename _aoSystemT>
46struct turbAtmosphere;
47
48/// Simulation of a single turbulent layer
49/** \todo document this
50 * \todo add facility for changing interpolator
51 *
52 * \ingroup mxAOSim
53 */
54template <typename _aoSystemT>
56{
57 typedef _aoSystemT aoSystemT;
58 typedef typename aoSystemT::realT realT;
59 typedef Eigen::Array<realT, -1, -1> imageT;
60
61 turbAtmosphere<aoSystemT> *m_parent{ nullptr };
62
63 uint32_t m_layerNo;
64
65 uint32_t m_scrnSz;
66
67 imageT m_freq;
68 imageT m_psd;
69
70 realT m_dx{ 0 };
71 realT m_dy{ 0 };
72
73 realT m_x0;
74 realT m_y0;
75
76 int m_last_wdx;
77 int m_last_wdy;
78
79 imageT m_phase;
80 imageT m_shiftPhaseWP;
81 imageT m_shiftPhase;
82 imageT m_shiftPhaseWork;
83
84 mx::math::uniDistT<realT> uniVar; ///< Uniform deviate, used in shiftRandom.
85
86 void setLayer( turbAtmosphere<aoSystemT> *parent, int layerNo, uint32_t scrnSz );
87
88 uint32_t layerNo();
89
90 uint32_t scrnSz();
91
92 void alloc();
93
94 /// Deallocate memory necessary for phase screen generation.
95 void genDealloc();
96
97 /// Shift to a timestep.
98 /**
99 * \param [in] dt is the new timestep.
100 */
101 void shift( realT dt );
102
103 /// Seed the uniform deviation. Call this if you intend to use shiftRandom.
104 /** This only needs to be called once.
105 */
106 void initRandom();
107
108 /// Shift by a random amount using the uniform distribution
109 /** Call initRandom() once before calling this method.
110 *
111 * \param [in] nofract if true then the fractional part is ignored, and only a whole-pixel shift is executed.
112 * Default=false
113 */
114 void shiftRandom( bool nofract = false );
115};
116
117// template<typename aoSystemT>
118// turbLayer<aoSystemT>::turbLayer()
119// {
120// }
121
122template <typename aoSystemT>
123void turbLayer<aoSystemT>::setLayer( turbAtmosphere<aoSystemT> *parent, int layerNo, uint32_t scrnSz )
124{
125 m_parent = parent;
126 m_layerNo = layerNo;
127
128 m_scrnSz = scrnSz;
129
130 if( m_parent == nullptr )
131 {
132 mxThrowException(
133 err::paramnotset, "mx::AO::sim::turbLayer::setLayer", "parent is not set (m_parent is nullptr)" );
134 }
135
136 if( m_parent->aosys() == nullptr )
137 {
138 mxThrowException( err::paramnotset,
139 "mx::AO::sim::turbLayer::setLayer",
140 "parent ao system is not set (m_parent->aosys() is nullptr)" );
141 }
142
143 realT vwind = m_parent->aosys()->atm.layer_v_wind( m_layerNo );
144 realT dwind = m_parent->aosys()->atm.layer_dir( m_layerNo );
145 realT pupD = m_parent->aosys()->D();
146 uint32_t wfSz = m_parent->wfSz();
147
148 m_dx = vwind * cos( dwind ) / ( pupD / wfSz );
149 m_dy = vwind * sin( dwind ) / ( pupD / wfSz );
150
151 alloc();
152}
153
154template <typename aoSystemT>
155uint32_t turbLayer<aoSystemT>::layerNo()
156{
157 return m_layerNo;
158}
159
160template <typename aoSystemT>
161uint32_t turbLayer<aoSystemT>::scrnSz()
162{
163 return m_scrnSz;
164}
165
166template <typename aoSystemT>
167void turbLayer<aoSystemT>::alloc()
168{
169 if( m_parent == nullptr )
170 {
171 mxThrowException(
172 err::paramnotset, "mx::AO::sim::turbLayer::alloc", "parent is not set (m_parent is nullptr)" );
173 }
174
175 if( m_parent->aosys() == nullptr )
176 {
177 mxThrowException( err::paramnotset,
178 "mx::AO::sim::turbLayer::alloc",
179 "parent ao system is not set (m_parent->aosys() is nullptr)" );
180 }
181
182 m_phase.resize( m_scrnSz, m_scrnSz );
183
184 uint32_t wfSz = m_parent->wfSz();
185
186 uint32_t buffSz = m_parent->buffSz();
187
188 m_shiftPhaseWP.resize( wfSz + 2 * buffSz, wfSz + 2 * buffSz );
189
190 m_shiftPhase.resize( wfSz + 2 * buffSz, wfSz + 2 * buffSz );
191
192 m_shiftPhaseWork.resize( wfSz + 2 * buffSz, wfSz + 2 * buffSz );
193
194 initRandom();
195
196 m_x0 = 0;
197 m_y0 = 0;
198 m_last_wdx = m_scrnSz + 1;
199 m_last_wdy = m_scrnSz + 1;
200
201 // Now set up for generation
202 realT r0 = m_parent->aosys()->atm.r_0( m_parent->aosys()->lam_sci() );
203 realT L0 = m_parent->aosys()->atm.L_0( m_layerNo );
204 realT l0 = m_parent->aosys()->atm.l_0( m_layerNo );
205
206 m_psd.resize( m_scrnSz, m_scrnSz );
207
208 m_freq.resize( m_scrnSz, m_scrnSz );
209 sigproc::frequencyGrid<imageT>( m_freq, m_parent->aosys()->D() / wfSz );
210
211 realT beta = 0.0218 / pow( r0, 5. / 3. );
212 realT sqrt_alpha = 0.5 * 11. / 3.;
213
214 realT L02;
215 if( L0 > 0 )
216 L02 = 1.0 / ( L0 * L0 );
217 else
218 L02 = 0;
219
220#pragma omp parallel for
221 for( size_t jj = 0; jj < m_scrnSz; ++jj )
222 {
223 for( size_t ii = 0; ii < m_scrnSz; ++ii )
224 {
225 realT p;
226 if( m_freq( ii, jj ) == 0 && L02 == 0 )
227 {
228 p = 0;
229 }
230 else
231 {
232 p = beta / pow( pow( m_freq( ii, jj ), 2 ) + L02, sqrt_alpha );
233 if( l0 > 0 )
234 p *= exp( -1 * pow( m_freq( ii, jj ) * l0, 2 ) );
235 }
236
237 realT Ppiston = 0;
238 realT Ptiptilt = 0;
239 if( m_parent->aosys()->psd.subPiston() )
240 {
241 Ppiston =
242 pow( 2 * math::func::jinc( math::pi<realT>() * m_freq( ii, jj ) * m_parent->aosys()->D() ), 2 );
243 }
244 if( m_parent->aosys()->psd.subTipTilt() )
245 {
246 Ptiptilt =
247 pow( 4 * math::func::jincN( 2, math::pi<realT>() * m_freq( ii, jj ) * m_parent->aosys()->D() ), 2 );
248 }
249
250 m_psd( ii, jj ) = sqrt( p * ( 1.0 - Ppiston - Ptiptilt ) );
251 }
252 }
253
254 if( m_parent->shLevel() > 0 )
255 {
256 realT onoff = 1;
257 if( m_parent->outerSubHarmonics() )
258 {
259 onoff = 0; // just make all of these 0
260 }
261
262 m_psd( 0, 0 ) = 0;
263 m_psd( 0, 1 ) *= 0.5 * 0.5 * onoff;
264 m_psd( 1, 0 ) *= 0.5 * 0.5 * onoff;
265 m_psd( m_psd.rows() - 1, 0 ) *= 0.5 * 0.5 * onoff;
266 m_psd( 0, m_psd.cols() - 1 ) *= 0.5 * 0.5 * onoff;
267
268 m_psd( 1, 1 ) *= 0.75 * 0.75 * onoff;
269 m_psd( m_psd.rows() - 1, m_psd.cols() - 1 ) *= 0.75 * 0.75 * onoff;
270 m_psd( m_psd.rows() - 1, 1 ) *= 0.75 * 0.75 * onoff;
271 m_psd( 1, m_psd.cols() - 1 ) *= 0.75 * 0.75 * onoff;
272 }
273}
274
275template <typename aoSystemT>
277{
278 m_freq.resize( 0, 0 );
279 m_psd.resize( 0, 0 );
280}
281
282template <typename aoSystemT>
284{
285 int wdx, wdy;
286 realT ddx, ddy;
287
288 ddx = m_x0 + m_dx * dt;
289 ddy = m_y0 + m_dy * dt;
290
291 wdx = (int)trunc( ddx );
292 ddx -= wdx;
293
294 wdy = (int)trunc( ddy );
295 ddy -= wdy;
296
297 wdx %= m_scrnSz;
298 wdy %= m_scrnSz;
299
300 // Check for a new whole-pixel shift
301 if( wdx != m_last_wdx || wdy != m_last_wdy )
302 {
303 // Need a whole pixel shift (this also extracts the m_wfSz + m_buffSz subarray)
304 improc::imageShiftWP( m_shiftPhaseWP, m_phase, wdx, wdy );
305 }
306
307 // Do the sub-pixel shift
308 improc::imageShift( m_shiftPhase, m_shiftPhaseWP, ddx, ddy, improc::cubicConvolTransform<realT>( -0.5 ) );
309
310 m_last_wdx = wdx;
311 m_last_wdy = wdy;
312}
313
314template <typename aoSystemT>
316{
317 uniVar.seed();
318}
319
320template <typename aoSystemT>
322{
323 int wdx, wdy;
324 realT ddx, ddy;
325
326 ddx = uniVar * ( m_scrnSz );
327 ddy = uniVar * ( m_scrnSz );
328
329 wdx = (int)trunc( ddx );
330 ddx -= wdx;
331 wdy = (int)trunc( ddy );
332 ddy -= wdy;
333
334 if( nofract )
335 {
336 ddx = 0;
337 ddy = 0;
338 }
339
340 improc::imageShiftWP( m_shiftPhaseWP, m_phase, wdx, wdy );
341
342 if( ddx != 0 || ddy != 0 )
343 {
344 improc::imageShift( m_shiftPhase, m_shiftPhaseWP, ddx, ddy, improc::cubicConvolTransform<realT>( -0.5 ) );
345 }
346 else
347 {
348 m_shiftPhase = m_shiftPhaseWP;
349 }
350}
351
352} // namespace sim
353} // namespace AO
354} // namespace mx
355
356#endif // mx_AO_sim_turbLayer_hpp
mxException for parameters which aren't set
A random number type, which functions like any other arithmetic type.
Definition randomT.hpp:64
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.
void imageShift(arrOutT &transim, const arrInT &im, floatT1 dx, floatT2 dy, transformT trans)
Shift an image.
void imageShiftWP(outputArrT &out, inputArrT &in, int dx, int dy, bool wrap=true)
Shift an image by whole pixels with (optional) wrapping.
The mxlib c++ namespace.
Definition mxError.hpp:106
A turbulent atmosphere simulator.
Simulation of a single turbulent layer.
Definition turbLayer.hpp:56
void genDealloc()
Deallocate memory necessary for phase screen generation.
void initRandom()
Seed the uniform deviation. Call this if you intend to use shiftRandom.
mx::math::uniDistT< realT > uniVar
Uniform deviate, used in shiftRandom.
Definition turbLayer.hpp:84
void shift(realT dt)
Shift to a timestep.
void shiftRandom(bool nofract=false)
Shift by a random amount using the uniform distribution.
Transformation by cubic convolution interpolation.