mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
Loading...
Searching...
No Matches
imageXCorrFFT.hpp
Go to the documentation of this file.
1/** \file imageXCorrFFT.hpp
2 * \brief A class to register images using the Fourier cross correlation with a peak fit.
3 * \ingroup image_processing_files
4 * \author Jared R. Males (jaredmales@gmail.com)
5 *
6 */
7
8//***********************************************************************//
9// Copyright 2020 Jared R. Males (jaredmales@gmail.com)
10//
11// This file is part of mxlib.
12//
13// mxlib is free software: you can redistribute it and/or modify
14// it under the terms of the GNU General Public License as published by
15// the Free Software Foundation, either version 3 of the License, or
16// (at your option) any later version.
17//
18// mxlib is distributed in the hope that it will be useful,
19// but WITHOUT ANY WARRANTY; without even the implied warranty of
20// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21// GNU General Public License for more details.
22//
23// You should have received a copy of the GNU General Public License
24// along with mxlib. If not, see <http://www.gnu.org/licenses/>.
25//***********************************************************************//
26
27#ifndef imageXCorrFFT_hpp
28#define imageXCorrFFT_hpp
29
30#include "../mxError.hpp"
31#include "../mxException.hpp"
32#include "../math/ft/fftT.hpp"
33#include "../math/ft/mftT.hpp"
34#include "../math/ft/ftUtils.hpp"
35#include "../math/ft/fftwEnvironment.hpp"
36
37#include "../math/fit/fitGaussian.hpp"
38
39#include "eigenImage.hpp"
40#include "imageUtils.hpp"
41#include "imageTransforms.hpp"
42#include "imageXCorr.hpp"
43
44namespace mx
45{
46namespace improc
47{
48
49/// Find the optimum shift to align two images using the FFT cross correlation.
50/** The reference image must be the same size as the target image. Both the reference image and the
51 * target image being analyzed can be masked, normalized by mean subtraction and variance division,
52 * and windowed. Seperate masks and windows for the reference and target image are used.
53 *
54 * The reference mask, reference window, and normalization flag must be set before setting the reference
55 * (unless using the defaults). If the reference mask, reference window, or normalization flag are changed
56 * then the reference must be set again. After configuring the reference, then operator() can be called
57 * repeatedly to determine the shifts for a sequence of imaages. The mask and window can be changed between
58 * calls without requiring resetting the reference. No new heap allocations take place on these calls
59 * to operator(), and the reference image is not re-normalized on each call.
60 *
61 * The shift is reported in pixels such that if the mxlib imageShift function is used
62 * to shift the input image by the negative of the shifts, it will align with the
63 * reference at the center of the array.
64 *
65 * Several peak finding methods are provided:
66 * - xcorrPeakMethod::centerOfLight uses center of light
67 * - xcorrPeakMethod::gaussFit uses Gaussian centroiding
68 * - xcorrPeakMethod::interpPeak uses interpolation to find the peak to a given tolerance.
69 * - xcorrPeakMethod::none will just return the coordinate of the max pixel
70 *
71 * This class uses the real-to-complex transforms of FFTW by default. This should offer
72 * improved speed and memory use. For testing the complex-to-complex transforms can be
73 * used by defining `XCFFT_C2C` before including this file.
74 *
75 * \tparam _realImageT is the Eigen-like array type used for image processing. Must have
76 * a real floating point Scalar type (float, double, etc.). See typedefs.
77 *
78 * \ingroup image_reg
79 */
80template <class _realImageT>
82{
83 public:
84 typedef _realImageT realImageT; ///< the Eigen-like array type used for image processing
85
86 typedef typename _realImageT::Scalar realT; ///< the scalar type of the image type
87
88 typedef std::complex<realT> complexT; ///< Complex floating point type.
89
90// defaut is the real-to-complex. Only use C2C if defined:
91#ifdef XCFFT_C2C
92 typedef complexT fftOutT;
93 typedef complexT fftInT;
94#else
95 typedef realT fftOutT;
96 typedef complexT fftInT;
97#endif
98
99 typedef eigenImage<complexT> complexArrayT; ///< Complex eigen array type with Scalar==complexT
100
101 protected:
102 /** \name Configuration Member Data
103 * @{
104 */
105 int m_rows{ 0 }; ///< The number of rows in the images
106
107 int m_cols{ 0 }; ///< The number of columns in the images
108
109 realT m_padFactor{ 1 }; ///< The padding factor for the CC
110
111 realT m_padFactorR{ 1 }; /**< The actual padding factor for rows of the CC
112 set based on the resultant array sizes.*/
113
114 realT m_padFactorC{ 1 }; /**< The actual padding factor for cols of the CC
115 set based on the resultant array sizes.*/
116
117 int m_rowsPadded{ 0 }; ///< The number of rows in the padded CC
118
119 int m_colsPadded{ 0 }; ///< The number of columns in the padded CC
120
121 realImageT m_refMaskIm; /**< Mask image to use for the reference, may be needed for proper normalization
122 even if refIm has 0 mask applied.*/
123
124 bool m_haveRefMask{ false }; ///< Flag indicating that a referece mask has been provided.
125
126 realImageT m_maskIm; /**< Mask image to use, may be needed for proper normalization even if
127 refIm has 0 mask applied.*/
128
129 bool m_haveMask{ false }; ///< Flag indicating that a mask has been provided.
130
131 realImageT m_refWinIm; ///< Window image to use for the reference
132
133 bool m_haveRefWindow{ false }; ///< Flag indicating that a window has been provided for the reference
134
135 realImageT m_winIm; ///< Window image for the target image.
136
137 bool m_haveWindow{ false }; ///< Flag indicating that a window has been provided.
138
139 bool m_normalize{ true }; /**< Flag specifying whether images should be normalized prior to performing
140 the cross-correlation*/
141
142 realImageT m_refIm; ///< The normalized reference image.
143
144 bool m_refValid{ false }; /**< Flag indicating whether the reference is valid. It can be invalidated
145 if a mask, window, or normalize is set after the reference is set.*/
146
147 int m_maxLag{ 5 }; ///< The maximum lag to consider in the initial cross-correlation. Default is 5.
148
149 realT m_tol{ 0.1 }; /**< For xcorrPeakMethod::interpPeak, the tolerance of the
150 interpolated-magnified image, in pixels.*/
151
152 realT m_magSize{ 0 }; /**< Magnified size of the ccIm when using interp.
153 Set as function of m_tol and m_maxLag.*/
154
155 realT m_mag{ 1 }; /**< The magnification, in principle is 1/m_tol but will vary
156 slightly based on integer array sizes*/
157
158 /// The peak finding method to use
159 xcorrPeakMethod m_peakMethod{ xcorrPeakMethod::interpPeak };
160
161 ///@}
162
163 /** \name Working Memory Member Data
164 * @{
165 */
166
167 realImageT m_refACIm; ///< The auto-correlation image of the reference.
168
169 realT m_refX0{ 0 }; /**< The x-shift of the reference image to itself using the selected
170 algorithm, used as coordinate origin*/
171
172 realT m_refY0{ 0 }; /**< The x-shift of the reference image to itself using the selected
173 algorithm, used as coordinate origin*/
174
175 realImageT m_normIm; ///< The normalized image.
176
177 realImageT m_ccIm; ///< The cross-correlation image
178
179 realImageT m_magIm; ///< The magnified image, used if m_peakMethod == xcorrPeakMethod::interpPeak
180
181 complexArrayT m_ftIm0; ///< Working memory for the FT of the reference image.
182
183 public:
184 complexArrayT m_ftWork; ///< Working memory for the FFT.
185
186#ifdef XCFFT_C2C
187 complexArrayT m_ftWorkIn; ///< Working memory for the FFT input.
188#endif
189
190 complexArrayT m_ftWorkPadded; ///< Working memory for the padded inverse FFT input.
191
192 complexArrayT m_ftWorkPaddedOut; ///< Working memory for the FFT and MFT output.
193
194 complexArrayT m_mtWorkPadded; ///< Working memory for the MFT for xcorrPeakMethod::mftOversamp
195
196 math::ft::fftT<fftOutT, fftInT, 2> m_fft_fwd; ///< FFT object for the forward transform.
197
198 math::ft::fftT<fftInT, fftOutT, 2> m_fft_bwd; ///< FFT object for the backward transform.
199
200 math::ft::mftT<complexT, complexT, 2> m_mft_bwd; /**< MFT object used to oversample the
201 peak for xcorrPeakMethod::mftOversamp*/
202
203 ///@}
204
205 /** \name Peak Finding Tools
206 * @{
207 */
209
210 ///@}
211
212 public:
213 /** \name Construction and Destruction
214 * @{
215 */
216
217 /// Default c'tor
219
220 ///@}
221
222 /** \name Configuration Member Access
223 * @{
224 */
225
226 /// Get the number of rows in the input images
227 /** This can only be set by a call to resize() or refIm().
228 *
229 * \returns the current value of m_rows
230 */
231 int rows();
232
233 /// Get the number of columns in the input images
234 /** This can only be set by a call to resize() or refIm().
235 *
236 * \returns the current value of m_cols
237 */
238 int cols();
239
240 /// Set the padding factor
241 /** This can also be set by a call to resize() or refIm().
242 * The padding factor must be greater than or equal to 1.
243 *
244 * If m_rows and m_cols are both greater than 0 this calls resize(int, int, realT).
245 *
246 */
247 int padFactor( realT os /**< [in] the new padding factor*/ );
248
249 /// Get the padding factor
250 /**
251 * \returns the current value of m_padFactor
252 */
254
255 /// Get the number of rows in the padded cross-correlation
256 /** This can only be set by a call to resize() or refIm().
257 *
258 * \returns the current value of m_rowsPadded
259 */
260 int rowsPadded();
261
262 /// Get the number of rows in the padded cross-correlation
263 /** This can only be set by a call to resize() or refIm().
264 *
265 * \returns the current value of m_colsPadded
266 */
267 int colsPadded();
268
269 /// Set the size of the cross-correlation problem based on input image size.
270 /** This resizes all working memory and conducts fftw planning.
271 * If there are no changes to size or padding, this returns immediately
272 * making no changes.
273 *
274 * If there are changes to size or padding the reference is invalidated.
275 *
276 * If nrows or ncols is 0, this resizes all working memory to 0.
277 *
278 * \returns 0 on success
279 * \returns -1 on error
280 */
281 int resize( int nrows, /**< [in] the number of rows in the images to register */
282 int ncols, /**< [in] the number of columns in the images to register */
283 realT padFactor /**< [in] the padding factor. Must be `>=` 1. */ );
284
285 /// Set the size of the cross-correlation problem based on input image size.
286 /** Calls resize(int, int, realT) with the current value of m_padFactor.
287 *
288 * \overload
289 *
290 * \returns 0 on success
291 * \returns -1 on error
292 */
293 int resize( int nrows, /**< [in] the number of rows in the images to register */
294 int ncols /**< [in] the number of columns in the images to register */ );
295
296 /// Set the reference mask image
297 /** This invalidates the reference
298 * \returns 0 on success
299 * \returns -1 on error
300 */
301 int refMaskIm( const realImageT &mask /**< [in] the new reference mask image */ );
302
303 /// Get a reference to the reference mask image.
304 /**
305 * \returns a const reference to the reference mask image.
306 */
307 const realImageT &refMaskIm();
308
309 /// Get the value of the haveRefMask flag
310 /**
311 * \returns the current value of m_haveRefMask
312 */
313 bool haveRefMask();
314
315 /// Set the mask image
316 /** This does not invalidate the reference
317 *
318 * \returns 0 on success
319 * \returns -1 on error
320 */
321 int maskIm( const realImageT &mask /**< [in] the new mask image */ );
322
323 /// Get a reference to the mask image.
324 /**
325 * \returns a const referance to the mask image.
326 */
327 const realImageT &maskIm();
328
329 /// Get the value of the haveMask flag
330 /**
331 * \returns the current value of m_haveMask
332 */
333 bool haveMask();
334
335 /// Set the reference window image
336 /** This invalidates the reference.
337 *
338 * \returns 0 on success
339 * \returns -1 on error
340 */
341 int refWinIm( const realImageT &win /**< [in] the new reference window image */ );
342
343 /// Get a reference to the reference window image.
344 /**
345 * \returns a const referance to the reference window image.
346 */
347 const realImageT &refWinIm();
348
349 /// Get the value of the haveRefWindow flag
350 /**
351 * \returns the current value of m_haveRefWindow
352 */
353 bool haveRefWindow();
354
355 /// Set the window image
356 /** This does not invalidate the reference.
357 *
358 * \returns 0 on success
359 * \returns -1 on error
360 */
361 int winIm( const realImageT &win /**< [in] the new window image */ );
362
363 /// Get a reference to the window image.
364 /**
365 * \returns a const referance to the window image.
366 */
367 const realImageT &winIm();
368
369 /// Get the value of the haveWindow flag
370 /**
371 * \returns the current value of m_haveWindow
372 */
373 bool haveWindow();
374
375 /// Set the normalize flag
376 /** Changing this invalidates the reference.
377 */
378 void normalize( bool no /**< [in] the new normalize flag value */ );
379
380 /// Get the current value of the normalize flag
381 /**
382 * \returns the current value of m_normalize
383 */
384 bool normalize();
385
386 /// Set the reference image
387 /** Performs the following steps
388 * - Applies the mask if supplied.
389 * - If desired, normalizes the reference image by mean subtraction and variance division.
390 * - Applies the window if supplied.
391 * - Stores the result as m_refIm
392 * - Resizes working memory with a call to resize()
393 * - Shifts to FFT order and Fourier transforms and conjugates the result.
394 * - Calls operator() to measure the reference shifts
395 * - Stores the resultant auto-correlation in m_refACIm
396 *
397 * \returns 0 on success
398 * \returns -1 on error
399 */
400 int refIm( const realImageT &im0, /**< [in] The new reference image */
401 realT padFactor /**< [in] [optional] the padding factor to use */ );
402
403 /// Set the reference image
404 /** Calls refIm(const realImageT &, realT) with the current value of m_padFactor.
405 *
406 * \overload
407 *
408 * \returns 0 on success
409 * \returns -1 on error
410 */
411 int refIm( const realImageT &im0 /**< [in] The new reference image */ );
412
413 /// Get a reference to the reference image.
414 /**
415 * \returns a const referent to m_refIm.
416 */
417 const realImageT &refIm();
418
419 /// Get the value of the reference valid flag
420 /**
421 * \returns the value of m_refValid
422 */
423 bool refValid();
424
425 /// Set the maximum lag
426 void maxLag( int ml /**< [in] the new maximum lag */ );
427
428 /// Get the current maximum lag
429 /**
430 * \returns the current value of m_maxLag
431 */
432 int maxLag();
433
434 /// Set the tolerance of the interpolated-magnified image, in pixels.
435 void tol( realT nt /**< [in] The new value of the interpolation tolerance. */ );
436
437 /// Get the tolerance of the interpolated-magnified image, in pixels.
438 /**
439 * \returns the current value of m_tol.
440 */
441 realT tol();
442
443 void peakMethod( xcorrPeakMethod xpm );
444
445 xcorrPeakMethod peakMethod();
446
447 ///@}
448
449 /** \name Working Memory Member Access
450 * @{
451 */
452
453 /// Get a reference to the reference auto-correlation image.
454 /**
455 * \returns a const referent to m_refACIm.
456 */
457 const realImageT &refACIm();
458
459 /// Get the x-shift of the reference image to itself.
460 /**
461 * \returns m_refX0
462 */
463 realT refX0();
464
465 /// Get the y-shift of the reference image to itself.
466 /**
467 * \returns m_refY0
468 */
469 realT refY0();
470
471 /// Get a reference to the normalized image.
472 /**
473 * \returns a const referent to m_normIm.
474 */
475 const realImageT &normIm();
476
477 /// Get a reference to the cross correlation image.
478 /**
479 * \returns a const referent to m_ccIm.
480 */
481 const realImageT &ccIm();
482
483 /// Get a reference to the magnified image.
484 /**
485 * \returns a const referent to m_magIm.
486 */
487 const realImageT &magIm();
488
489 ///@}
490
491 /** \name Peak Finding
492 * @{
493 */
494 protected:
495 void findPeak( realT &xShift, /**< [out] the x shift of im w.r.t. im0, in pixels */
496 realT &yShift /**< [out] the y shift of im w.r.t. im0, in pixels */ );
497
498 ///@}
499
500 /** \name Cross Correlation Operator
501 * @{
502 */
503
504 public:
505 /// Conduct the cross correlation to a specified tolerance
506 /**
507 * \returns 0 on success
508 * \returns -1 on error
509 */
510 template <class imT>
511 int operator()( realT &xShift, ///< [out] the x shift of im w.r.t. im0, in pixels
512 realT &yShift, ///< [out] the y shift of im w.r.t. im0, in pixels
513 const imT &im /**< [in] the image to cross-correlate with the reference*/ );
514
515 /// Conduct the cross correlation to a specified tolerance
516 /**
517 * \overload
518 *
519 * \returns 0 on success
520 * \returns -1 on error
521 */
522 template <class im0T, class imT>
523 int operator()( realT &xShift, ///< [out] the x shift of im w.r.t. im0, in pixels
524 realT &yShift, ///< [out] the y shift of im w.r.t. im0, in pixels
525 im0T &im0, ///< [in] a new reference image
526 imT &im ///< [in] the image to cross-correlate with the reference
527 );
528
529 ///@}
530};
531
532template <class realImageT>
534{
535 maxLag( m_maxLag ); // this sets up m_magSize.
536}
537
538template <class realImageT>
540{
541 return m_rows;
542}
543
544template <class realImageT>
546{
547 return m_cols;
548}
549
550template <class realImageT>
552{
553 if( os < 1 )
554 {
555 mxError( "imageXCorrFFT::padFactor", MXE_INVALIDARG, "padding factor can't be less than 1" );
556 return MXE_INVALIDARG;
557 }
558
559 if( m_rows > 0 && m_cols > 0 )
560 {
561 return resize( m_rows, m_cols, os );
562 }
563 else
564 {
565 m_padFactor = os;
566 return 0;
567 }
568}
569
570template <class realImageT>
575
576template <class realImageT>
578{
579 return m_rowsPadded;
580}
581
582template <class realImageT>
584{
585 return m_colsPadded;
586}
587
588template <class realImageT>
589int imageXCorrFFT<realImageT>::resize( int nrows, int ncols, realT padFactor )
590{
591 // No-op if no change
592 if( m_rows == nrows && m_cols == ncols && m_padFactor == padFactor )
593 {
594 return 0;
595 }
596
597 if( padFactor < 1 )
598 {
599 mxError( "imageXCorrFFT::resize", MXE_INVALIDARG, "padding factor can't be less than 1" );
600 return MXE_INVALIDARG;
601 }
602
603 m_rows = nrows;
604
605 m_cols = ncols;
606
607 m_padFactor = padFactor;
608
609 /* Any change here invalidates the reference.
610 * This is trivially true for rows or cols. Also must be
611 * handled for padding which will change the reference shifts.
612 */
613 m_refValid = false;
614
615 if( m_rows == 0 || m_cols == 0 )
616 {
617 m_ftIm0.resize( 0, 0 );
618 m_ftWork.resize( 0, 0 );
619#ifdef XCFFT_C2C
620 m_ftWorkIn.resize( 0, 0 );
621#endif
622 m_rowsPadded = 0;
623 m_colsPadded = 0;
624 m_ccIm.resize( 0, 0 );
625 m_ftWorkPadded.resize( 0, 0 );
626 m_mtWorkPadded.resize( 0, 0 );
627 m_ftWorkPaddedOut.resize( 0, 0 );
628
629 return 0;
630 }
631
632#ifdef XCFFT_C2C
633
634 m_ftIm0.resize( m_rows, m_cols );
635
636 m_ftWork.resize( m_rows, m_cols );
637
638 m_ftWorkIn.resize( m_rows, m_cols );
639
640#else
641
642 m_ftIm0.resize( (int)( 0.5 * m_rows ) + 1, m_cols );
643
644 m_ftWork.resize( (int)( 0.5 * m_rows ) + 1, m_cols );
645
646#endif
647
648 // fftw is row-major, eigen defaults to column-major
649 m_fft_fwd.plan( m_cols, m_rows, math::ft::dir::forward, false );
650
651 m_rowsPadded = ceil( m_rows * m_padFactor );
652 m_colsPadded = ceil( m_cols * m_padFactor );
653
654 m_padFactorR = static_cast<realT>( m_rowsPadded ) / static_cast<realT>( m_rows );
655 m_padFactorC = static_cast<realT>( m_colsPadded ) / static_cast<realT>( m_cols );
656
657 m_ccIm.resize( m_rowsPadded, m_colsPadded );
658
659#ifdef XCFFT_C2C
660
661 m_ftWorkPadded.resize( m_rowsPadded, m_colsPadded );
662 m_ftWorkPaddedOut.resize( m_rowsPadded, m_colsPadded );
663
664#else
665
666 m_ftWorkPadded.resize( (int)( 0.5 * m_rowsPadded ) + 1, m_colsPadded );
667
668 if( m_peakMethod == xcorrPeakMethod::mftOversamp )
669 {
670 m_ftWorkPaddedOut.resize( m_rowsPadded, m_colsPadded );
671 }
672 else
673 {
674 m_ftWorkPaddedOut.resize( 0, 0 );
675 }
676
677#endif
678
679 if( m_peakMethod == xcorrPeakMethod::mftOversamp )
680 {
681 m_mtWorkPadded.resize( m_rowsPadded, m_colsPadded );
682 }
683 else
684 {
685 m_mtWorkPadded.resize( 0, 0 );
686 }
687
688 m_fft_bwd.plan( m_colsPadded, m_rowsPadded, math::ft::dir::backward, false );
689
690 return 0;
691} // resize(int, int, realT)
692
693template <class realImageT>
694int imageXCorrFFT<realImageT>::resize( int nrows, int ncols )
695{
696 return resize( nrows, ncols, m_padFactor );
697}
698
699template <class realImageT>
701{
702 m_refMaskIm = mask;
703 m_haveRefMask = true;
704
705 m_refValid = false;
706
707 return 0;
708}
709
710template <class realImageT>
712{
713 return m_refMaskIm;
714}
715
716template <class realImageT>
718{
719 return m_haveRefMask;
720}
721
722template <class realImageT>
724{
725 m_maskIm = mask;
726 m_haveMask = true;
727
728 return 0;
729}
730
731template <class realImageT>
733{
734 return m_maskIm;
735}
736
737template <class realImageT>
739{
740 return m_haveMask;
741}
742
743template <class realImageT>
745{
746 m_refWinIm = win;
747 m_haveRefWindow = true;
748
749 m_refValid = false;
750
751 return 0;
752}
753
754template <class realImageT>
756{
757 return m_refWinIm;
758}
759
760template <class realImageT>
762{
763 return m_haveRefWindow;
764}
765
766template <class realImageT>
768{
769 m_winIm = win;
770 m_haveWindow = true;
771
772 return 0;
773}
774
775template <class realImageT>
777{
778 return m_winIm;
779}
780
781template <class realImageT>
783{
784 return m_haveWindow;
785}
786
787template <class realImageT>
789{
790 if( no != m_normalize )
791 {
792 m_refValid = false;
793 m_normalize = no;
794 }
795}
796
797template <class realImageT>
799{
800 return m_normalize;
801}
802
803template <class realImageT>
805{
806 realImageT im0;
807
808 // Mask if needed
809 if( m_haveRefMask )
810 {
811 if( im.rows() != m_refMaskIm.rows() && im.cols() != m_refMaskIm.cols() )
812 {
813 mxError( "imageXCorrFFT::refIm", MXE_SIZEERR, "reference and reference mask are not the same size" );
814 return MXE_SIZEERR;
815 }
816
817 // Now normalize
818 if( m_normalize )
819 {
820 realT m = imageMean<realT>( im, m_refMaskIm );
821 realT v = imageVariance<realT>( im, m, m_refMaskIm );
822 im0 = ( ( im - m ) * m_refMaskIm ) / sqrt( v );
823 }
824 else
825 {
826 im0 = im * m_refMaskIm;
827 }
828 }
829 else
830 {
831 // Normalize
832 if( m_normalize )
833 {
834 realT m = imageMean<realT>( im );
835 realT v = imageVariance<realT>( im, m );
836 im0 = ( im - m ) / sqrt( v );
837 }
838 else
839 {
840 im0 = im;
841 }
842 }
843
844 // Window if needed
845 if( m_haveRefWindow )
846 {
847 if( im0.rows() != m_refWinIm.rows() && im.cols() != m_refWinIm.cols() )
848 {
849 mxError( "imageXCorrFFT::refIm", MXE_SIZEERR, "reference and reference window are not the same size" );
850 return MXE_SIZEERR;
851 }
852
853 im0 *= m_refWinIm;
854 }
855
856 // We save refIm as the un-shifted version
857 m_refIm = im0;
858
859 /* Setup the FFTW space
860 * Note: we don't do this earlier in case we ever implement something (e.g. a filter)
861 * that changes size before this point
862 */
863 resize( im0.rows(), im0.cols(), padFactor );
864
865 // Now shift so center pixel is 0,0
866 imageShiftWP( im0, m_refIm, 0.5 * m_rows + 1, 0.5 * m_cols + 1 );
867
868// Then FT
869#ifdef XCFFT_C2C
870
871 m_ftWork = im0.template cast<complexT>();
872
873 m_fft_fwd( m_ftIm0.data(), m_ftWork.data() );
874#else
875 m_fft_fwd( m_ftIm0.data(), im0.data() );
876#endif
877
878 // Conjugate and normalize for FFTW scaling.
879 for( int c = 0; c < m_ftIm0.cols(); ++c )
880 {
881 for( int r = 0; r < m_ftIm0.rows(); ++r )
882 {
883 complexT val = std::conj( m_ftIm0( r, c ) );
884 val /= ( m_rows * m_cols );
885 m_ftIm0( r, c ) = val;
886 }
887 }
888
889 // And finally find the reference shift
890 m_refX0 = 0;
891 m_refY0 = 0;
892
893 m_refValid = true;
894
895 operator()( m_refX0, m_refY0, m_refIm );
896
897 m_refACIm = m_ccIm;
898
899 return 0;
900}
901
902template <class realImageT>
904{
905 return refIm( im, m_padFactor );
906}
907
908template <class realImageT>
910{
911 return m_refIm;
912}
913
914template <class realImageT>
916{
917 return m_refValid;
918}
919
920template <class realImageT>
922{
923 m_maxLag = ml;
924 tol( m_tol );
925}
926
927template <class realImageT>
929{
930 return m_maxLag;
931}
932
933template <class realImageT>
938
939template <class realImageT>
941{
942 m_tol = nt;
943
944 m_magSize = ceil( ( ( 2. * m_maxLag + 1 ) - 1.0 ) / nt ) + 1;
945
946 m_mag = ( m_magSize - 1.0 ) / ( ( 2. * m_maxLag + 1 ) - 1.0 );
947
948 if( m_refIm.rows() != m_rows || m_refIm.cols() != m_cols || m_rows == 0 || m_cols == 0 )
949 return;
950
951 // Find the reference shift
952 m_refX0 = 0;
953 m_refY0 = 0;
954 operator()( m_refX0, m_refY0, m_refIm );
955}
956
957template <class realImageT>
959{
960 m_peakMethod = xpm;
961
962 if( xpm == xcorrPeakMethod::mftOversamp )
963 {
964 m_refValid = false;
965 }
966
967 if( m_refIm.rows() != m_rows || m_refIm.cols() != m_cols || m_rows == 0 || m_cols == 0 )
968 {
969 return;
970 }
971
972 // Find the reference shift
973 m_refX0 = 0;
974 m_refY0 = 0;
975 operator()( m_refX0, m_refY0, m_refIm );
976}
977
978template <class realImageT>
979xcorrPeakMethod imageXCorrFFT<realImageT>::peakMethod()
980{
981 return m_peakMethod;
982}
983
984template <class realImageT>
986{
987 return m_refACIm;
988}
989
990template <class realImageT>
995
996template <class realImageT>
1001
1002template <class realImageT>
1004{
1005 return m_normIm;
1006}
1007
1008template <class realImageT>
1010{
1011 return m_ccIm;
1012}
1013
1014template <class realImageT>
1016{
1017 return m_magIm;
1018}
1019
1020template <class realImageT>
1022{
1023 int xLag0, yLag0;
1024
1025 /*realImageT sm = m_ccIm;
1026 filterImage(sm, m_ccIm, gaussKernel<eigenImage<realT>, 2>(3));
1027 m_ccIm = sm;*/
1028
1029 // m_ccIm = m_ccIm.pow(3);
1030 realT pk = m_ccIm.maxCoeff( &xLag0, &yLag0 );
1031 realT mn = m_ccIm.minCoeff();
1032
1033 if( xLag0 - m_maxLag < 0 )
1034 {
1035 m_maxLag = xLag0;
1036 }
1037 if( xLag0 + 2 * m_maxLag + 1 >= m_ccIm.rows() )
1038 {
1039 m_maxLag = ( m_ccIm.rows() - 1 - xLag0 ) / 2;
1040 }
1041 if( yLag0 - m_maxLag < 0 )
1042 {
1043 m_maxLag = yLag0;
1044 }
1045 if( yLag0 + 2 * m_maxLag + 1 >= m_ccIm.cols() )
1046 {
1047 m_maxLag = ( m_ccIm.cols() - 1 - yLag0 ) / 2;
1048 }
1049
1050 realT x0 = xLag0 - m_maxLag;
1051 realT y0 = yLag0 - m_maxLag;
1052
1053 if( m_peakMethod == xcorrPeakMethod::gaussFit )
1054 {
1055 m_magIm = m_ccIm.block( x0, y0, 2 * m_maxLag + 1, 2 * m_maxLag + 1 );
1056
1057 m_fitter.setArray( m_magIm.data(), m_magIm.rows(), m_magIm.cols() );
1058 m_magIm.maxCoeff( &xLag0, &yLag0 );
1059
1060 m_fitter.setGuess( mn, pk - mn, m_maxLag, m_maxLag, 3, 3, 0 );
1061 m_fitter.fit();
1062
1063 x0 = ( m_fitter.x0() + x0 );
1064 y0 = ( m_fitter.y0() + y0 );
1065 }
1066 else if( m_peakMethod == xcorrPeakMethod::interpPeak )
1067 {
1068 int x, y;
1069 m_magIm.resize( m_magSize, m_magSize );
1070
1071 imageMagnify( m_magIm,
1072 m_ccIm.block( x0, y0, 2 * m_maxLag + 1, 2 * m_maxLag + 1 ),
1074
1075 m_magIm.maxCoeff( &x, &y );
1076 x0 = ( x * ( 1.0 / m_mag ) + x0 );
1077 y0 = ( y * ( 1.0 / m_mag ) + y0 );
1078 }
1079 else if( m_peakMethod == xcorrPeakMethod::centerOfLight )
1080 {
1081 realT x, y;
1082
1083 m_magIm = m_ccIm.block( x0, y0, 2 * m_maxLag + 1, 2 * m_maxLag + 1 );
1084 m_magIm -= m_magIm.minCoeff(); // Must sum to > 0.
1085
1086 imageCenterOfLight( x, y, m_magIm );
1087
1088 x0 = ( x + x0 );
1089 y0 = ( y + y0 );
1090 }
1091 else if( m_peakMethod == xcorrPeakMethod::mftOversamp )
1092 {
1093 int x, y;
1094 m_ccIm.maxCoeff( &x, &y );
1095
1096 if( m_refX0 == 0 && m_refY0 == 0 ) // reference finding
1097 {
1098 x0 = x;
1099 y0 = y;
1100 }
1101 else
1102 {
1103 x = x / m_padFactorR - m_refX0;
1104 y = y / m_padFactorC - m_refY0;
1105
1106 realT mag = 1.0 / m_tol; // m_mag is the effective mag for interpPeak...
1107
1108 realT tmagR = mag * m_padFactorR; // the total mag
1109 realT tmagC = mag * m_padFactorC; // the total mag
1110
1111 realT shiftx0 = -1 * ( m_refX0 * m_padFactorR ) * ( mag - 1 ) - x * tmagR;
1112 realT shifty0 = -1 * ( m_refY0 * m_padFactorC ) * ( mag - 1 ) - y * tmagC;
1113
1114 m_mft_bwd.plan( m_rowsPadded, m_colsPadded, math::ft::dir::backward, shiftx0, shifty0, mag );
1115
1116#ifdef XCFFT_C2C
1117
1118 m_mft_bwd( m_ftWorkPaddedOut, m_ftWorkPadded );
1119
1120#else
1121
1122 m_mtWorkPadded.resize( m_rowsPadded, m_colsPadded );
1123 m_mtWorkPadded.setZero();
1124 math::ft::augmentR2CFFTOutput( m_mtWorkPadded, m_ftWork );
1125
1126 m_mft_bwd( m_ftWorkPaddedOut, m_mtWorkPadded );
1127
1128#endif
1129
1130 m_magIm = m_ftWorkPaddedOut.real();
1131
1132 int xf, yf;
1133 m_magIm.maxCoeff( &xf, &yf );
1134
1135 xShift = x + ( xf - m_refX0 * m_padFactorR ) / tmagR;
1136 yShift = y + ( yf - m_refY0 * m_padFactorC ) / tmagC;
1137
1138 return;
1139 }
1140 }
1141 else if( m_peakMethod == xcorrPeakMethod::none )
1142 {
1143 int x, y;
1144 m_ccIm.maxCoeff( &x, &y );
1145 x0 = x;
1146 y0 = y;
1147 }
1148 else
1149 {
1150 mxThrowException( mx::err::invalidconfig, "imageXCorrFFT::operator()", "unknown peak finding method" );
1151 }
1152
1153 //--> unpad here, scaling the shifts
1154 xShift = x0 / m_padFactorR - m_refX0;
1155 yShift = y0 / m_padFactorC - m_refY0;
1156}
1157
1158template <class realImageT>
1159template <class imT>
1160int imageXCorrFFT<realImageT>::operator()( realT &xShift, realT &yShift, const imT &im )
1161{
1162 if( !m_refValid )
1163 {
1164 mxThrowException( err::invalidconfig, "imageXCorrFFT", "reference image is not valid" );
1165 }
1166
1167 if( im.rows() != m_rows )
1168 {
1169 mxThrowException( err::sizeerr, "imageXCorrFFT", "image must be same size as reference (rows)" );
1170 }
1171
1172 if( im.cols() != m_cols )
1173 {
1174 mxThrowException( err::sizeerr, "imageXCorrFFT", "image must be same size as reference (rows)" );
1175 }
1176
1177 // Mask and normalize as needed
1178 if( m_haveMask )
1179 {
1180 if( m_normalize )
1181 {
1182 realT m = imageMean<realT>( im, m_maskIm );
1183 realT v = imageVariance<realT>( im, m, m_maskIm );
1184
1185 m_normIm = ( im - m ) * m_maskIm / sqrt( v );
1186 }
1187 else
1188 {
1189 m_normIm = im * m_maskIm;
1190 }
1191 }
1192 else
1193 {
1194
1195 if( m_normalize )
1196 {
1197 realT m = imageMean<realT>( m_normIm );
1198 realT v = imageVariance<realT>( m_normIm, m );
1199 m_normIm = ( im - m ) / sqrt( v );
1200 }
1201 else
1202 {
1203 m_normIm = im;
1204 }
1205 }
1206
1207 if( m_haveWindow )
1208 {
1209 m_normIm *= m_winIm;
1210 }
1211
1212#ifdef XCFFT_C2C
1213 m_ftWorkIn = m_normIm.template cast<complexT>();
1214
1215 m_fft_fwd( m_ftWork.data(), m_ftWorkIn.data() );
1216#else
1217 m_fft_fwd( m_ftWork.data(), m_normIm.data() );
1218#endif
1219
1220 // m_ftIm0 is the conjugated, fftw-normalized, reference image in the Fourier domain
1221 // So this is the FT of the cross-correlation:
1222 m_ftWork *= m_ftIm0;
1223
1224 m_ftWorkPadded.setZero();
1225
1226#ifdef XCFFT_C2C
1227
1228 math::ft::padC2CFFTOutput( m_ftWorkPadded, m_ftWork );
1229
1230 m_fft_bwd( m_ftWorkPaddedOut.data(), m_ftWorkPadded.data() );
1231
1232 m_ccIm = m_ftWorkPaddedOut.real();
1233
1234#else
1235
1236 math::ft::padR2CFFTOutput( m_ftWorkPadded, m_ftWork );
1237
1238 m_fft_bwd( m_ccIm, m_ftWorkPadded );
1239
1240#endif
1241
1242 findPeak( xShift, yShift );
1243
1244 return 0;
1245}
1246
1247template <class realImageT>
1248template <class im0T, class imT>
1249int imageXCorrFFT<realImageT>::operator()( realT &xShift, realT &yShift, im0T &im0, imT &im )
1250{
1251 setReference( im0 );
1252 return operator()( xShift, yShift, im );
1253}
1254
1255} // namespace improc
1256} // namespace mx
1257
1258#endif // imageXCorrFFT_hpp
mxException for invalid config settings
mxException for a size error
Find the optimum shift to align two images using the FFT cross correlation.
realT padFactor()
Get the padding factor.
int m_cols
The number of columns in the images.
realImageT m_refACIm
The auto-correlation image of the reference.
int m_rowsPadded
The number of rows in the padded CC.
int colsPadded()
Get the number of rows in the padded cross-correlation.
eigenImage< complexT > complexArrayT
Complex eigen array type with Scalar==complexT.
int m_maxLag
The maximum lag to consider in the initial cross-correlation. Default is 5.
complexArrayT m_mtWorkPadded
Working memory for the MFT for xcorrPeakMethod::mftOversamp.
bool haveMask()
Get the value of the haveMask flag.
bool haveWindow()
Get the value of the haveWindow flag.
int m_rows
The number of rows in the images.
math::ft::fftT< fftOutT, fftInT, 2 > m_fft_fwd
FFT object for the forward transform.
const realImageT & maskIm()
Get a reference to the mask image.
const realImageT & refIm()
Get a reference to the reference image.
realImageT m_refIm
The normalized reference image.
bool haveRefMask()
Get the value of the haveRefMask flag.
int resize(int nrows, int ncols, realT padFactor)
Set the size of the cross-correlation problem based on input image size.
int maxLag()
Get the current maximum lag.
realImageT m_magIm
The magnified image, used if m_peakMethod == xcorrPeakMethod::interpPeak.
std::complex< realT > complexT
Complex floating point type.
realT m_padFactor
The padding factor for the CC.
bool m_haveMask
Flag indicating that a mask has been provided.
complexArrayT m_ftWorkPadded
Working memory for the padded inverse FFT input.
complexArrayT m_ftIm0
Working memory for the FT of the reference image.
complexArrayT m_ftWork
Working memory for the FFT.
const realImageT & refWinIm()
Get a reference to the reference window image.
bool haveRefWindow()
Get the value of the haveRefWindow flag.
complexArrayT m_ftWorkIn
Working memory for the FFT input.
math::ft::mftT< complexT, complexT, 2 > m_mft_bwd
const realImageT & refACIm()
Get a reference to the reference auto-correlation image.
const realImageT & ccIm()
Get a reference to the cross correlation image.
complexArrayT m_ftWorkPaddedOut
Working memory for the FFT and MFT output.
int m_colsPadded
The number of columns in the padded CC.
int rowsPadded()
Get the number of rows in the padded cross-correlation.
realT refX0()
Get the x-shift of the reference image to itself.
bool refValid()
Get the value of the reference valid flag.
math::ft::fftT< fftInT, fftOutT, 2 > m_fft_bwd
FFT object for the backward transform.
bool normalize()
Get the current value of the normalize flag.
realImageT m_refWinIm
Window image to use for the reference.
realT tol()
Get the tolerance of the interpolated-magnified image, in pixels.
int operator()(realT &xShift, realT &yShift, const imT &im)
Conduct the cross correlation to a specified tolerance.
const realImageT & refMaskIm()
Get a reference to the reference mask image.
bool m_haveRefWindow
Flag indicating that a window has been provided for the reference.
const realImageT & magIm()
Get a reference to the magnified image.
realT refY0()
Get the y-shift of the reference image to itself.
_realImageT::Scalar realT
the scalar type of the image type
bool m_haveRefMask
Flag indicating that a referece mask has been provided.
realImageT m_ccIm
The cross-correlation image.
realImageT m_winIm
Window image for the target image.
realImageT m_normIm
The normalized image.
int rows()
Get the number of rows in the input images.
bool m_haveWindow
Flag indicating that a window has been provided.
const realImageT & winIm()
Get a reference to the window image.
const realImageT & normIm()
Get a reference to the normalized image.
xcorrPeakMethod m_peakMethod
The peak finding method to use.
void findPeak(realT &xShift, realT &yShift)
int cols()
Get the number of columns in the input images.
_realImageT realImageT
the Eigen-like array type used for image processing
Class to manage fitting a 2D Gaussian to data via the levmarInterface.
Tools for using the eigen library for image processing.
Eigen::Array< scalarT, -1, -1 > eigenImage
Definition of the eigenImage type, which is an alias for Eigen::Array.
void augmentR2CFFTOutput(eigenArrOutT &aout, eigenArrInT &ain)
Definition ftUtils.hpp:105
void padR2CFFTOutput(eigenArrOutT &aout, eigenArrInT &ain)
Add padding to the output of a real-to-complex (R2C) FFT.
Definition ftUtils.hpp:50
void padC2CFFTOutput(eigenArrOutT &aout, eigenArrInT &ain)
Add padding to the output of a complex-to-complex (C2C) FFT.
Definition ftUtils.hpp:74
@ backward
Specifies the backward transform.
@ forward
Specifies the forward transform.
xcorrPeakMethod
Methods for finding the cross-correlation peak.
void imageMagnify(arrOutT &transim, const arrInT &im, transformT trans)
Magnify an image.
void imageShiftWP(outputArrT &out, inputArrT &in, int dx, int dy, bool wrap=true)
Shift an image by whole pixels with (optional) wrapping.
int imageCenterOfLight(typename imageT::Scalar &x, typename imageT::Scalar &y, const imageT &im)
Calculate the center of light of an image.
Image interpolation and transformation.
Header for the image processing utilities.
CRTP base class to register images.
The mxlib c++ namespace.
Definition mxError.hpp:106
Transformation by cubic convolution interpolation.