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
generalIntegrator.hpp
Go to the documentation of this file.
1/** \file generalIntegrator.hpp
2 * \author Jared R. Males (jaredmales@gmail.com)
3 * \brief Declares and defines a general integrator controller class for AO control.
4 * \ingroup mxAO_sim_files
5 *
6 */
7
8#ifndef __generalIntegrator_hpp__
9#define __generalIntegrator_hpp__
10
11namespace mx
12{
13namespace AO
14{
15namespace sim
16{
17
18template <typename _realT>
19struct wfMeasurement
20{
21 typedef _realT realT;
22
23 realT iterNo;
24
25 // typedef Eigen::Array< realT, Eigen::Dynamic, Eigen::Dynamic> commandT;
26 typedef std::vector<realT> commandT;
27
28 commandT measurement;
29};
30
31/// Implements a general integrator controller.
32/** \todo document the math for the G.I.
33 * \todo document the filter coefficients, etc.
34 *
35 * \tparam _realT is the floating point type for all calculations.
36 */
37template <typename _realT>
39{
40
41 public:
42 /// The real data type
43 typedef _realT realT;
44
45 /// The real data type
46 typedef std::complex<realT> complexT;
47
48 /// The wavefront data type
49 // typedef wavefront<realT> wavefrontT;
50
51 /// The command type
52 typedef wfMeasurement<realT> commandT;
53
54 /// The image type, used here as a general storage array
55 typedef Eigen::Array<realT, Eigen::Dynamic, Eigen::Dynamic> imageT;
56
57 /// Default c'tor.
59
60 protected:
61 int m_nModes{ 0 }; ///< The number of modes being filtered.
62
63 bool m_openLoop{ false }; ///< If true, then commands are not integrated. Default is false.
64
65 int m_openLoopDelay{ 0 }; ///< If > 0, then the loop is open for this time in time steps. Default is 0.
66
67 int m_closingRamp{ 0 }; ///< If > 0, then gains are famped linearly up to m_closingGains over this interval in time
68 ///< steps. Default is 0. This is relative m_openLoopDelay.
69
70 int m_closingDelay{ 0 }; ///< If > 0, then the simple integrator,with m_closingGains is used up to this timestep.
71 ///< This is relative m_openLoopDelay. Default = 0.
72
73 imageT m_closingGains; ///< Column-vector of gains used for loop closing as simple integrator
74
75 int m_lowOrders{ 0 }; ///< If > 0, then this sets the maximum mode number which is filtered. All remaining modes are
76 ///< set to 0. Default = 0.
77
78 imageT m_a;
79 int m_currA;
80
81 imageT m_b;
82 int m_currB;
83
84 imageT m_gains; ///< Column-vector of gains
85
86 imageT m_commandsIn;
87 imageT m_commandsOut;
88
89 public:
90 /// Allocate and initialize all state.
91 /**
92 * \returns 0 on success, a negative integer otherwise.
93 */
94 int initialize( int nModes /**< [in] the number of modes to be filtered */ );
95
96 /// Get the number of modes
97 /** nModes is only set by calling initialize.
98 *
99 * \returns the current value of m_nModes
100 */
101 int nModes();
102
103 /// Set the m_openLoop flag.
104 /** If m_openLoop is true, then commands are not filtered.
105 *
106 * \returns 0 on success, a negative integer on error.
107 */
108 int openLoop( bool ol /**< [in] the new value of m_openLoop */ );
109
110 /// Get the value of the m_openLoop flag.
111 /**
112 * \returns the current value of m_openLoop.
113 */
114 bool openLoop();
115
116 /// Set the open loop delay.
117 /** If m_openLoopDelay > 0, then the loop is open for this period in time-steps
118 * Requires m_closingDelay > 0.
119 *
120 * \returns 0 on success, a negative integer on error.
121 */
122 int openLoopDelay( int cr /**< [in] The new value of m_openLoopDelay */ );
123
124 /// Get the value of the m_openLoopDelay.
125 /**
126 * \returns the current value of m_openLoopDelay.
127 */
128 int openLoopDelay();
129
130 /// Set the closing ramp.
131 /** If m_closingRamp > 0, then the gains are ramped linearly up to m_closingGains over this timer interval in
132 * timesteps. Requires m_closingDelay > 0.
133 *
134 * \returns 0 on success, a negative integer on error.
135 */
136 int closingRamp( int cr /**< [in] The new value of m_closingRamp */ );
137
138 /// Get the value of the m_closingRamp.
139 /**
140 * \returns the current value of m_closingRamp.
141 */
142 int closingRamp();
143
144 /// Set the closing delay.
145 /** If m_closingDelay > 0, then the simple integragor is used until this timestep.
146 *
147 * \returns 0 on success, a negative integer on error.
148 */
149 int closingDelay( int cd /**< [in] The new value of m_closingDelay */ );
150
151 /// Get the value of the m_closingDelay.
152 /**
153 * \returns the current value of m_closingDelay.
154 */
155 int closingDelay();
156
157 /// Set the simple integrator gains to use during closing.
158 /**
159 * \returns 0 on success
160 * \returns < 0 on error
161 */
162 int closingGains( const std::vector<realT> &gains /**< [in] vector of gains.*/ );
163
164 /// Set the simple integrator gain to use for all modes during closing.
165 /**
166 * \returns 0 on success
167 * \returns < 0 on error
168 */
169 int closingGains( realT g );
170
171 /// Set m_lowOrders.
172 /** If m_lowOrders > 0, then this sets the maximum mode number which is filtered. All remaining modes are set to 0.
173 *
174 * \returns 0 on success, a negative integer on error.
175 */
176 int lowOrders( int lo /**< [in] The new value of m_lowOrders */ );
177
178 /// Get the value of the m_lowOrders.
179 /**
180 * \returns the current value of m_lowOrders.
181 */
182 int lowOrders();
183
184 /// Set the size of the IIR vector (the a coefficients)
185 /** This allocates m_a to be m_nModes X n in size.
186 *
187 * \returns 0 on success, negative number on error.
188 */
189 int setASize( int n /**< [in] the number of IIR coefficients */ );
190
191 /// Set the IIR coefficients for one mode.
192 /**
193 * \returns 0 on success, negative number on error.
194 */
195 int setA( int i, ///< [in] the mode number
196 const imageT &a ///< [in] the IIR coefficients for this mode
197 );
198
199 /// Set the size of the FIR vector (the b coefficients)
200 /** This allocates m_b to be m_nModes X n in size.
201 *
202 * \returns 0 on success, negative number on error.
203 */
204 int setBSize( int n /**< [in] the number of FIR coefficients */ );
205
206 /// Set the FIR coefficients for one mode.
207 /**
208 * \returns 0 on success, negative number on error.
209 */
210 int setB( int i, ///< [in] the mode number
211 const imageT &b ///< [in] the IIR coefficients for this mode
212 );
213
214 /// Get the gain for a single mode.
215 /**
216 * \returns the gain value if mode exists.
217 */
218 realT gain( int i /**< The mode number*/ );
219
220 /// Set the gain for a single mode.
221 /**
222 * \returns 0 on success, negative number on error.
223 */
224 int gain( int i, ///< The mode number
225 realT g ///< The new gain value
226 );
227
228 /// Set the gain for all modes to a single value
229 /**
230 * \returns 0 on success, negative number on error.
231 */
232 int gains( realT g /**< The new gain value*/ );
233
234 /// Set the gains for all modes, using a vector to specify each gain
235 /** The vector must be exactly as long as m_nModes.
236 *
237 * \returns 0 on success, negative number on error.
238 */
239 int gains( const std::vector<realT> &gains /**< [in] vector of gains.*/ );
240
241 /// Set the gains for all modes, using a file to specify each gain
242 /** The file format is a simple ASCII single column, with 1 gain per line.
243 * Must be exactly as long as m_nModes.
244 *
245 * \returns 0 on success, negative number on error.
246 */
247 int gains( const std::string &ogainf /**< [in] the name of the file, full path */ );
248
249 /// Allocate the provided command structures
250 /** Used by the calling system to allocate the commands being passed between components.
251 *
252 * \returns 0 on success, negative number on error.
253 */
254 int initMeasurements( commandT &filtAmps, ///< The structure to contain the filtered commands
255 commandT &rawAmps ///< The structure to contain the raw commands
256 );
257
258 int filterCommands( commandT &filtAmps, commandT &rawAmps, int iterNo );
259};
260
261template <typename realT>
265
266template <typename realT>
268{
269 m_nModes = nModes;
270
271 // If m_a has been sized, resize it
272 if( m_a.cols() > 0 )
273 {
274 int n = m_a.cols();
275 m_a.resize( m_nModes, n );
276 m_a.setZero();
277 m_currA = 0;
278
279 m_commandsOut.resize( m_nModes, n );
280 m_commandsOut.setZero();
281 }
282
283 // If m_b has been sized, resize it
284 if( m_b.cols() > 0 )
285 {
286 int n = m_b.cols();
287 m_b.resize( m_nModes, n );
288 m_b.setZero();
289 m_currB = 0;
290
291 m_commandsIn.resize( m_nModes, n );
292 m_commandsIn.setZero();
293 }
294
295 m_closingGains.resize( 1, m_nModes );
296
297 m_closingGains.setZero();
298
299 m_gains.resize( 1, m_nModes );
300
301 m_gains.setZero();
302
303 return 0;
304}
305
306template <typename realT>
308{
309 return m_nModes;
310}
311
312template <typename realT>
314{
315 m_openLoop = ol;
316 return 0;
317}
318
319template <typename realT>
321{
322 return m_openLoop;
323}
324
325template <typename realT>
327{
328 m_openLoopDelay = old;
329
330 return 0;
331}
332
333template <typename realT>
335{
336 return m_openLoopDelay;
337}
338
339template <typename realT>
341{
342 m_closingRamp = cr;
343
344 return 0;
345}
346
347template <typename realT>
349{
350 return m_closingRamp;
351}
352
353template <typename realT>
355{
356 m_closingDelay = cd;
357
358 return 0;
359}
360
361template <typename realT>
363{
364 return m_closingDelay;
365}
366
367template <typename realT>
368int generalIntegrator<realT>::closingGains( const std::vector<realT> &gains )
369{
370
371 if( gains.size() != (size_t)m_closingGains.cols() )
372 {
373 mxError( "generalIntegrator::closingGains", MXE_SIZEERR, "input gain vector not same size as number of modes" );
374 return -1; ///\retval -1 on vector size mismatch
375 }
376
377 for( int i = 0; i < m_closingGains.cols(); ++i )
378 {
379 m_closingGains( 0, i ) = gains[i];
380 }
381
382 return 0;
383}
384
385template <typename realT>
387{
388 for( int i = 0; i < m_closingGains.cols(); ++i )
389 {
390 m_closingGains( 0, i ) = g;
391 }
392
393 return 0;
394}
395
396template <typename realT>
398{
399 m_lowOrders = lo;
400
401 return 0;
402}
403
404template <typename realT>
406{
407 return m_lowOrders;
408}
409
410template <typename realT>
412{
413 // Resize with m_nModes if set
414 if( m_nModes > 0 )
415 {
416 m_a.resize( m_nModes, n );
417 m_commandsOut.resize( m_nModes, n );
418 }
419 else
420 {
421 m_a.resize( 1, n );
422 m_commandsOut.resize( 1, n );
423 }
424
425 m_a.setZero();
426 m_currA = 0;
427 m_commandsOut.setZero();
428
429 return 0;
430}
431
432template <typename realT>
434{
435 m_a.row( i ) = a;
436
437 return 0;
438}
439
440template <typename realT>
442{
443 // Resize with m_nModes if set
444 if( m_nModes > 0 )
445 {
446 m_b.resize( m_nModes, n );
447 m_commandsIn.resize( m_nModes, n );
448 }
449 else
450 {
451 m_b.resize( 1, n );
452 m_commandsIn.resize( 1, n );
453 }
454
455 m_b.setZero();
456 m_currB = 0;
457 m_commandsIn.setZero();
458
459 return 0;
460}
461
462template <typename realT>
464{
465 m_b.row( i ) = b;
466
467 return 0;
468}
469
470template <class realT>
472{
473 if( i < 0 || i >= m_nModes )
474 {
475 mxError( "generalIntegrator::gain", MXE_INVALIDARG, "mode index out of range" );
476 return 0; ///\retval 0 if the mode doesn't exist.
477 }
478
479 return m_gains( i );
480}
481
482template <typename realT>
484{
485 if( i < 0 || i >= m_nModes )
486 {
487 mxError( "generalIntegrator::gain", MXE_INVALIDARG, "mode index out of range" );
488 return 0; ///\retval 0 if the mode doesn't exist.
489 }
490
491 m_gains( 0, i ) = g;
492
493 return 0;
494}
495
496template <typename realT>
498{
499 for( int i = 0; i < m_gains.cols(); ++i )
500 {
501 m_gains( 0, i ) = g;
502 }
503
504 return 0;
505}
506
507template <typename realT>
508int generalIntegrator<realT>::gains( const std::vector<realT> &gains )
509{
510
511 if( gains.size() != (size_t)m_gains.cols() )
512 {
513 mxError( "generalIntegrator::gains", MXE_SIZEERR, "input gain vector not same size as number of modes" );
514 return -1; ///\retval -1 on vector size mismatch
515 }
516
517 for( int i = 0; i < m_gains.cols(); ++i )
518 {
519 m_gains( 0, i ) = gains[i];
520 }
521
522 return 0;
523}
524
525template <typename realT>
526int generalIntegrator<realT>::gains( const std::string &ogainf )
527{
528 std::ifstream fin;
529 fin.open( ogainf );
530
531 if( !fin.good() )
532 {
533 mxError( "generalIntegrator::gains", MXE_FILEOERR, "could not open gan file" );
534 return -1; /// \retval -1 if file open fails
535 }
536
537 std::string tmpstr;
538 realT g;
539 for( int i = 0; i < m_gains.cols(); ++i )
540 {
541 fin >> tmpstr;
542 g = ioutils::convertFromString<realT>( tmpstr );
543 gain( i, g );
544 }
545
546 return 0;
547}
548
549template <class realT>
551{
552 filtAmps.measurement.resize( m_nModes );
553 for( size_t n = 0; n < m_nModes; ++n )
554 filtAmps.measurement[n] = 0;
555 // filtAmps.measurement.setZero();
556
557 rawAmps.measurement.resize( m_nModes );
558 for( size_t n = 0; n < m_nModes; ++n )
559 rawAmps.measurement[n] = 0;
560 // rawAmps.measurement.resize(1, m_nModes);
561 // rawAmps.measurement.setZero();
562
563 return 0;
564}
565
566template <typename realT>
567int generalIntegrator<realT>::filterCommands( commandT &filtAmps, commandT &rawAmps, int iterNo )
568{
569 filtAmps.iterNo = rawAmps.iterNo;
570
571 if( m_openLoop || iterNo < m_openLoopDelay )
572 {
573 for( size_t n = 0; n < m_nModes; ++n )
574 filtAmps.measurement[n] = 0;
575 // filtAmps.measurement.setZero();
576 return 0;
577 }
578
579 realT aTot, bTot;
580
581 if( iterNo - m_openLoopDelay < m_closingDelay )
582 {
583 BREAD_CRUMB;
584
585 for( int i = 0; i < m_nModes; ++i )
586 {
587 // if( std::isnan( rawAmps.measurement(0,i) ) || !std::isfinite(rawAmps.measurement(0,i)))
588 // rawAmps.measurement(0,i) = 0.0;
589 if( std::isnan( rawAmps.measurement[i] ) || !std::isfinite( rawAmps.measurement[i] ) )
590 rawAmps.measurement[i] = 0.0;
591
592 // m_commandsIn(i, m_currB) = rawAmps.measurement(0,i);
593 m_commandsIn( i, m_currB ) = rawAmps.measurement[i];
594
595 aTot = m_commandsOut( i, m_currA );
596 // std::cerr << m_currA << " " << aTot << "\n";
597
598 int cA = m_currA + 1;
599 if( cA >= m_a.cols() )
600 cA = 0;
601
602 realT gf = 1.0;
603
604 if( iterNo - m_openLoopDelay < m_closingRamp )
605 gf = ( (realT)iterNo - m_openLoopDelay ) / m_closingRamp;
606
607 // m_commandsOut(i, cA) = aTot + gf*m_closingGains(i) * rawAmps.measurement(0,i);
608 m_commandsOut( i, cA ) = aTot + gf * m_closingGains( i ) * rawAmps.measurement[i];
609
610 // std::cerr << cA << " " << m_commandsOut(i, cA) << "\n";
611
612 if( i <= m_lowOrders || m_lowOrders <= 0 )
613 {
614 filtAmps.measurement[i] = m_commandsOut( i, cA );
615 }
616 else
617 {
618 filtAmps.measurement[i] = 0;
619 }
620 }
621
622 ++m_currB;
623 if( m_currB >= m_b.cols() )
624 m_currB = 0;
625
626 ++m_currA;
627 if( m_currA >= m_a.cols() )
628 m_currA = 0;
629
630 return 0;
631 }
632
633 for( int i = 0; i < m_nModes; ++i )
634 {
635 if( m_gains( i ) == 0 )
636 {
637 // filtAmps.measurement(0,i) = 0;
638 filtAmps.measurement[i] = 0;
639 continue;
640 }
641
642 // if( std::isnan( rawAmps.measurement(0,i) ) || !std::isfinite(rawAmps.measurement(0,i)))
643 // rawAmps.measurement(0,i) = 0.0;
644 if( std::isnan( rawAmps.measurement[i] ) || !std::isfinite( rawAmps.measurement[i] ) )
645 rawAmps.measurement[i] = 0.0;
646
647 aTot = 0;
648 for( int j = 0; j < m_a.cols(); ++j )
649 {
650 int k = m_currA - j;
651 if( k < 0 )
652 k += m_a.cols();
653 aTot += m_a( i, j ) * m_commandsOut( i, k );
654 }
655
656 m_commandsIn( i, m_currB ) = rawAmps.measurement[i];
657
658 bTot = 0;
659 for( int j = 0; j < m_b.cols(); ++j )
660 {
661 int k = m_currB - j;
662 if( k < 0 )
663 k += m_b.cols();
664
665 bTot += m_b( i, j ) * m_commandsIn( i, k );
666 }
667
668 int cA = m_currA + 1;
669 if( cA >= m_a.cols() )
670 cA = 0;
671
672 m_commandsOut( i, cA ) = aTot + m_gains( i ) * bTot;
673
674 filtAmps.measurement[i] = m_commandsOut( i, cA );
675 }
676
677 ++m_currB;
678 if( m_currB >= m_b.cols() )
679 m_currB = 0;
680
681 ++m_currA;
682 if( m_currA >= m_a.cols() )
683 m_currA = 0;
684
685 return 0;
686}
687
688} // namespace sim
689} // namespace AO
690} // namespace mx
691
692#endif //__generalIntegrator_hpp__
Implements a general integrator controller.
imageT m_closingGains
Column-vector of gains used for loop closing as simple integrator.
int nModes()
Get the number of modes.
int m_openLoopDelay
If > 0, then the loop is open for this time in time steps. Default is 0.
int closingDelay()
Get the value of the m_closingDelay.
int setA(int i, const imageT &a)
Set the IIR coefficients for one mode.
bool openLoop()
Get the value of the m_openLoop flag.
wfMeasurement< realT > commandT
The wavefront data type.
bool m_openLoop
If true, then commands are not integrated. Default is false.
int gains(realT g)
Set the gain for all modes to a single value.
int closingGains(const std::vector< realT > &gains)
Set the simple integrator gains to use during closing.
int openLoopDelay()
Get the value of the m_openLoopDelay.
std::complex< realT > complexT
The real data type.
int setB(int i, const imageT &b)
Set the FIR coefficients for one mode.
int initialize(int nModes)
Allocate and initialize all state.
Eigen::Array< realT, Eigen::Dynamic, Eigen::Dynamic > imageT
The image type, used here as a general storage array.
int initMeasurements(commandT &filtAmps, commandT &rawAmps)
Allocate the provided command structures.
_realT realT
The real data type.
int lowOrders()
Get the value of the m_lowOrders.
int setBSize(int n)
Set the size of the FIR vector (the b coefficients)
int setASize(int n)
Set the size of the IIR vector (the a coefficients)
realT gain(int i)
Get the gain for a single mode.
int closingRamp()
Get the value of the m_closingRamp.
imageT m_gains
Column-vector of gains.
int m_nModes
The number of modes being filtered.
constexpr units::realT k()
Boltzmann Constant.
Definition constants.hpp:69
The mxlib c++ namespace.
Definition mxError.hpp:106