mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
Loading...
Searching...
No Matches
imageTransforms.hpp
Go to the documentation of this file.
1/** \file imageTransforms.hpp
2 * \author Jared R. Males
3 * \brief Image interpolation and transformation
4 * \ingroup image_processing_files
5 *
6 */
7
8//***********************************************************************//
9// Copyright 2015, 2016, 2017, 2018 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 improc_imageTransforms_hpp
28#define improc_imageTransforms_hpp
29
30#include <cstddef>
31#include <cmath>
32
33#include <iostream>
34#include <limits>
35
36namespace mx
37{
38namespace improc
39{
40
41/// Transformation by bi-linear interpolation
42/** \ingroup image_transforms
43 */
44template <typename _arithT>
46{
47 typedef _arithT arithT;
48
49 static const size_t width = 2;
50 static const size_t lbuff = 0;
51
52 template <typename arrT, typename arithT>
53 void operator()( arrT &kern, arithT x, arithT y )
54 {
55 kern.resize( width, width );
56
57 kern( 0, 0 ) = ( 1. - x ) * ( 1. - y );
58 kern( 0, 1 ) = ( 1. - x ) * y;
59 kern( 1, 0 ) = x * ( 1. - y );
60 kern( 1, 1 ) = x * y;
61 }
62};
63
64/// Typedef for bilinearTransform with single precision
65/** \ingroup image_transforms
66 */
68
69/// Typedef for bilinearTransform with double precision
70/** \ingroup image_transforms
71 */
73
74/// Transformation by cubic convolution interpolation
75/** Uses the cubic convolution interpolation kernel. See <a
76 * href="https://en.wikipedia.org/wiki/Bicubic_interpolation">https://en.wikipedia.org/wiki/Bicubic_interpolation </a>.
77 *
78 * The parameter \ref cubic should be left as the default -0.5 in most cases, which gives the bicubic spline
79 * interpolator. See <a
80 * href="https://en.wikipedia.org/wiki/Cubic_Hermite_spline">https://en.wikipedia.org/wiki/Cubic_Hermite_spline</a>.
81 *
82 * \tparam _arithT is the type in which to do all calculations. Should be a floating point type.
83 *
84 * \test Scenario: Verify direction and accuracy of various image shifts \ref tests_improc_imageTransforms_imageShift
85 * "[test doc]"
86 *
87 * \ingroup image_transforms
88 */
89template <typename _arithT>
91{
92 typedef _arithT arithT; ///< The type in which all calculations are performed.
93
94 static const size_t width = 4;
95 static const size_t lbuff = 1;
96
97 arithT cubic{ -0.5 }; ///< The kernel parameter. The default value -0.5 gives the bicubic spline interpolator.
98
99 /// Default c'tor.
100 /**
101 * This will provide the bicubic spline interpolator
102 */
104 {
105 }
106
107 /// Construct setting the kernel parameter.
108 explicit cubicConvolTransform( arithT c /**< [in] [optiona] The kernel parameter. The default value -0.5 gives the bicubic spline interpolator. */)
109 {
110 cubic = c;
111 }
112
113 /// Copy c'tor
115 {
116 cubic = t.cubic;
117 }
118
119 /// Calculate the kernel value for a given residual.
121 {
122 if( d <= 1 )
123 return ( cubic + 2. ) * d * d * d - ( cubic + 3. ) * d * d + 1.;
124
125 if( d < 2 )
126 return cubic * d * d * d - 5. * cubic * d * d + 8. * cubic * d - 4. * cubic;
127
128 return 0;
129 }
130
131 ///\todo why is this arithT not just the class's arithT?
132 template <typename arrT, typename arithT>
133 void operator()( arrT &kern, arithT x, arithT y )
134 {
135 arithT km2x, km1x, kp1x, kp2x;
136 arithT km2y, km1y, kp1y, kp2y;
137
138 km2x = cubicConvolKernel( ( 1. + x ) );
139 km1x = cubicConvolKernel( x );
140 kp1x = cubicConvolKernel( 1. - x );
141 kp2x = cubicConvolKernel( 2. - x );
142
143 km2y = cubicConvolKernel( ( 1. + y ) );
144 km1y = cubicConvolKernel( y );
145 kp1y = cubicConvolKernel( 1. - y );
146 kp2y = cubicConvolKernel( 2. - y );
147
148 kern( 0, 0 ) = km2x * km2y;
149 kern( 0, 1 ) = km2x * km1y;
150 kern( 0, 2 ) = km2x * kp1y;
151 kern( 0, 3 ) = km2x * kp2y;
152
153 kern( 1, 0 ) = km1x * km2y;
154 kern( 1, 1 ) = km1x * km1y;
155 kern( 1, 2 ) = km1x * kp1y;
156 kern( 1, 3 ) = km1x * kp2y;
157
158 kern( 2, 0 ) = kp1x * km2y;
159 kern( 2, 1 ) = kp1x * km1y;
160 kern( 2, 2 ) = kp1x * kp1y;
161 kern( 2, 3 ) = kp1x * kp2y;
162
163 kern( 3, 0 ) = kp2x * km2y;
164 kern( 3, 1 ) = kp2x * km1y;
165 kern( 3, 2 ) = kp2x * kp1y;
166 kern( 3, 3 ) = kp2x * kp2y;
167 }
168};
169
170/// Typedef for cubicConvolTransform with single precision
171/** \ingroup image_transforms
172 */
174
175/// Typedef for cubicConvolTransform with double precision
176/** \ingroup image_transforms
177 */
179
180/** \ingroup image_transforms
181 * @{
182 */
183
184/// Rotate an image represented as an eigen array
185/** Uses the given transformation type to rotate an image.
186 *
187 * \tparam transformT specifies the transformation to use [will be resolved by compiler]
188 * \tparam arrT is the eigen array type of the output [will be resolved by compiler]
189 * \tparam arrT2 is the eigen array type of the input [will be resolved by compiler]
190 * \tparam floatT is a floating point type [will be resolved by compiler in most cases]
191 *
192 */
193template <typename transformT, typename arrT, typename arrT2, typename floatT>
194void imageRotate( arrT &transim, ///< [out] The rotated image. Must be pre-allocated.
195 const arrT2 &im, ///< [in] The image to be rotated.
196 floatT dq, ///< [in] the angle, in radians, by which to rotate in the c.c.w. direction
197 transformT trans ///< [in] is the transformation to use
198)
199{
200 typedef typename transformT::arithT arithT;
201 arithT cosq, sinq;
202 arithT x0, y0, x, y;
203 arithT xcen, ycen;
204
205 int Nrows, Ncols;
206
207 int i0, j0;
208
209 const int lbuff = transformT::lbuff;
210 const int width = transformT::width;
211
212 cosq = cos( dq );
213 sinq = sin( dq );
214
215 Nrows = im.rows();
216 Ncols = im.cols();
217
218 transim.resize( Nrows, Ncols );
219
220 // The geometric image center
221 xcen = 0.5 * ( Nrows - 1. );
222 ycen = 0.5 * ( Ncols - 1. );
223
224 int xulim = Nrows - width + lbuff; // - 1;
225 int yulim = Ncols - width + lbuff; // - 1;
226
227 arithT xc_x_cosq = xcen * cosq;
228 arithT xc_x_sinq = xcen * sinq;
229 arithT yc_x_cosq = ycen * cosq;
230 arithT yc_x_sinq = ycen * sinq;
231
232 xc_x_cosq += yc_x_sinq;
233 xc_x_sinq -= yc_x_cosq;
234
235#ifdef MXLIB_USE_OMP
236#pragma omp parallel private( x0, y0, i0, j0, x, y )
237#endif
238 {
239 arithT i_x_cosq, i_x_sinq;
240 arrT kern;
241 kern.resize( width, width );
242
243#ifdef MXLIB_USE_OMP
244#pragma omp for schedule( static, 1 )
245#endif
246 for( int i = 0; i < Nrows; ++i )
247 {
248 i_x_cosq = i * cosq - xc_x_cosq; // + xcen;
249 i_x_sinq = -( i * sinq - xc_x_sinq ); // + ycen;
250
251 for( int j = 0; j < Ncols; ++j )
252 {
253 // We are actually doing this rotation matrix:
254 // x0 = (i-xcen)*cosq + (j-ycen)*sinq;
255 // y0 = -(i-xcen)*sinq + (j-ycen)*cosq;
256 // This is the minimum-op representation of the above rotation matrix:
257 x0 = i_x_cosq + j * sinq;
258 y0 = i_x_sinq + j * cosq;
259
260 // Get lower left index
261 i0 = x0 + xcen;
262 j0 = y0 + ycen;
263
264 if( i0 <= lbuff || i0 >= xulim || j0 <= lbuff || j0 >= yulim )
265 {
266 transim( i, j ) = 0;
267 continue;
268 }
269
270 // Get the residual
271 x = x0 + xcen - i0;
272 y = y0 + ycen - j0;
273
274 trans( kern, x, y );
275 transim( i, j ) = ( im.block( i0 - lbuff, j0 - lbuff, width, width ) * kern ).sum();
276 } // for j
277 } // for i
278 } // #pragma omp parallel
279
280} // void imageRotate(arrT & transim, const arrT2 &im, floatT dq, transformT trans)
281
282/// Shift an image by whole pixels with (optional) wrapping.
283/** The output image can be smaller than the input image, in which case the wrapping (if enabled) still occurs for the
284 * input image, but only output images worth of pixels are actually shifted. This is useful, for instance, when
285 * propagating large turbulence phase screens where one only needs a small section at a time.
286 *
287 * \tparam outputArrT is the eigen array type of the output [will be resolved by compiler]
288 * \tparam inputArrT is the eigen array type of the input [will be resolved by compiler]
289 *
290 * \test Scenario: Verify direction and accuracy of various image shifts \ref tests_improc_imageTransforms_imageShift
291 * "[test doc]"
292 */
293template <typename outputArrT, typename inputArrT>
294void imageShiftWP( outputArrT &out, ///< [out] contains the shifted image. Must be pre-allocated, but can be smaller
295 ///< than the in array.
296 inputArrT &in, ///< [in] the image to be shifted.
297 int dx, ///< [in] the amount to shift in the x direction
298 int dy, ///< [in] the amount to shift in the y direction
299 bool wrap = true ///< [in] flag controlling whether or not to wrap around
300)
301{
302 dx %= in.rows();
303 dy %= in.cols();
304
305 int outr = out.rows();
306 int outc = out.cols();
307 int inr = in.rows();
308 int inc = in.cols();
309
310#ifdef MXLIB_USE_OMP
311#pragma omp parallel
312#endif
313 {
314 int x, y;
315
316 if( wrap )
317 {
318#ifdef MXLIB_USE_OMP
319#pragma omp for
320#endif
321 for( int cc = 0; cc < outc; ++cc )
322 {
323 y = cc - dy;
324
325 if( y < 0 )
326 y += inc;
327 else if( y >= inc )
328 y -= inc;
329
330 for( int rr = 0; rr < outr; ++rr )
331 {
332 x = rr - dx;
333
334 if( x < 0 )
335 x += inr;
336 else if( x >= inr )
337 x -= inr;
338
339 out( rr, cc ) = in( x, y );
340 }
341 }
342 }
343 else
344 {
345#ifdef MXLIB_USE_OMP
346#pragma omp for
347#endif
348 for( int cc = 0; cc < outc; ++cc )
349 {
350 y = cc - dy;
351
352 if( y < 0 || y >= inc )
353 {
354 for( int rr = 0; rr < outr; ++rr )
355 {
356 out( rr, cc ) = 0;
357 }
358
359 continue;
360 }
361
362 for( int rr = 0; rr < outr; ++rr )
363 {
364 x = rr - dx;
365
366 if( x < 0 || x >= inr )
367 {
368 out( rr, cc ) = 0;
369 continue;
370 }
371
372 out( rr, cc ) = in( x, y );
373 }
374 }
375 } // if(wrap)-else
376 }
377}
378
379/// Shift an image by whole pixels, wrapping around, with a scaling image applied to the shifted image.
380/** The output image can be smaller than the input image, in which case the wrapping still occurs for the input image,
381 * but only output images worth of pixels are actually shifted. This is useful, for instance, when propagating large
382 * turbulence phase screens where one only needs a small section at a time.
383 *
384 * The scaling is applied to the output image. The scale image must be the same size as the output image.
385 *
386 * \tparam outputArrT is the eigen array type of the output [will be resolved by compiler]
387 * \tparam inputArrT is the eigen array type of the input [will be resolved by compiler]
388 * \tparam scaleArrT is the eigen array type of the scale image [will be resolved by compiler]
389 *
390 */
391template <typename outputArrT, typename inputArrT, typename scaleArrT>
393 outputArrT
394 &out, ///< [out] contains the shifted image. Must be pre-allocated, but can be smaller than the in array.
395 inputArrT &in, ///< [in] the image to be shifted.
396 scaleArrT &scale, ///< [in] image of scale values applied per-pixel to the output (shifted) image, same size as out
397 int dx, ///< [in] the amount to shift in the x direction
398 int dy ///< [in] the amount to shift in the y direction
399)
400{
401 dx %= in.rows();
402 dy %= in.cols();
403
404 int outr = out.rows();
405 int outc = out.cols();
406 int inr = in.rows();
407 int inc = in.cols();
408#ifdef MXLIB_USE_OMP
409// #pragma omp parallel
410#endif
411 {
412 int x, y;
413
414#ifdef MXLIB_USE_OMP
415// #pragma omp for
416#endif
417 for( int cc = 0; cc < outc; ++cc )
418 {
419 y = cc - dy;
420
421 if( y < 0 )
422 y += inc;
423 else if( y >= inc )
424 y -= inc;
425
426 for( int rr = 0; rr < outr; ++rr )
427 {
428 x = rr - dx;
429
430 if( x < 0 )
431 x += inr;
432 else if( x >= inr )
433 x -= inr;
434
435 out( rr, cc ) = in( x, y ) * scale( rr, cc );
436 }
437 }
438 }
439}
440
441/// Shift an image.
442/** Uses the given transformation type to shift an image such that objects move by (\p dx,\p dy) pixels.
443 * The shift is such that an object located at the coordinate
444 * (\p -dx, \p -dy) from the center of the image will be moved to the center of the image. So to move an object
445 * located 2 pixels right (dx) and 2 pixels up (dy) from the center to be at the center, use \p dx = -2, \p dy = -2.
446 *
447 * Note that this does not treat the edges
448 * of the image, determined by the buffer width (lbuff) of the kernel and the size of shift. If you wish to
449 * treat the edges, you must pad the image by at least lbuff+abs(shift) pixels in each direction, and
450 * implement a strategy (zeros, mirror, wrap) prior to calling this function.
451 *
452 * \tparam arrOutT is the Eigen-like array type of the output [will be resolved by compiler]
453 * \tparam arrInT is the Eigen-like array type of the input [will be resolved by compiler]
454 * \tparam floatT1 is a floating point type [will be resolved by compiler]
455 * \tparam floatT2 is a floating point type [will be resolved by compiler]
456 * \tparam transformT specifies the transformation to use [will be resolved by compiler]
457 *
458 * \test Scenario: Verify direction and accuracy of various image shifts \ref tests_improc_imageTransforms_imageShift
459 * "[test doc]"
460 */
461template <typename arrOutT, typename arrInT, typename floatT1, typename floatT2, typename transformT>
462void imageShift( arrOutT &transim, ///< [out] Will contain the shifted image. Will be allocated.
463 const arrInT &im, ///< [in] the image to be shifted.
464 floatT1 dx, ///< [in] the amount to shift in the x direction
465 floatT2 dy, ///< [in] the amount to shift in the y direction
466 transformT trans ///< [in] trans is the transformation to use
467)
468{
469 typedef typename transformT::arithT arithT;
470
471 int Nrows, Ncols;
472
473 const int lbuff = transformT::lbuff;
474 const int width = transformT::width;
475
476 // If this is a whole pixel, just do that.
477 if( dx == floor( dx ) && dy == floor( dy ) )
478 return imageShiftWP( transim, im, dx, dy, false );
479
480 Nrows = im.rows();
481 Ncols = im.cols();
482
483 int xulim = Nrows - width + lbuff;
484 int yulim = Ncols - width + lbuff;
485
486 transim.resize( Nrows, Ncols );
487
488#ifdef MXLIB_USE_OMP
489#pragma omp parallel
490#endif
491 {
492 int i0, j0;
493 // (rx, ry) is fractional residual of shift
494 arithT rx = 1 - ( dx - floor( dx ) );
495 arithT ry = 1 - ( dy - floor( dy ) );
496
497 arrOutT kern;
498 kern.resize( width, width );
499 trans( kern, rx, ry );
500
501#ifdef MXLIB_USE_OMP
502#pragma omp for
503#endif
504 for( int i = 0; i < Nrows; ++i )
505 {
506 // (i,j) is position in new image
507 // (i0,j0) is integer position in old image
508
509 i0 = i - dx;
510
511 if( i0 <= lbuff || i0 >= xulim )
512 {
513 for( int j = 0; j < Ncols; ++j )
514 {
515 transim( i, j ) = 0;
516 }
517 continue;
518 }
519
520 for( int j = 0; j < Ncols; ++j )
521 {
522 j0 = j - dy;
523
524 if( j0 <= lbuff || j0 >= yulim )
525 {
526 transim( i, j ) = 0;
527 continue;
528 }
529
530 transim( i, j ) = ( im.block( i0 - lbuff, j0 - lbuff, width, width ) * kern ).sum();
531 } // for j
532 } // for i
533 } // #pragam omp
534
535} // imageShift
536
537/// Magnify an image.
538/** Uses the given transformation type to magnify the input image to the size of the output image.
539 *
540 * Here we assume that the image center is the mxlib standard:
541 * \code
542 x_center = 0.5*(im.rows()-1);
543 y_center = 0.5*(im.cols()-1);
544 * \endcode
545 * Some care is necessary to prevent magnification from shifting the image with respect to this center. The main
546 result is that the
547 * magnification factors (which can be different in x and y) are defined thus:
548 * \code
549 x_mag = (transim.rows()-1.0) / (im.rows()-1.0);
550 y_mag = (transim.cols()-1.0) / (im.cols()-1.0);
551 * \endcode
552 *
553 * Example:
554 * \code
555 im1.resize(512,512);
556 //add image to im1
557 im2.resize(1024,1024);
558 imageMagnify(im2,im1, cubicConvolTransform<double>());
559 \endcode
560 * In this exmple, the image in im1 will be magnified by `1023.0/511.0 = 2.002x` and placed in im2.
561 *
562 * This transform function does not handle edges. If treatment of edges is desired, you must pad the input
563 * image using the desired strategy before calling this function. Note that the padded-size of the input image
564 * will affect the magnification factor.
565 *
566 * \tparam arrOutT is the eigen array type of the output.
567 * \tparam arrInT is the eigen array type of the input.
568 * \tparam transformT specifies the transformation to use.
569 */
570template <typename arrOutT, typename arrInT, typename transformT>
571void imageMagnify( arrOutT &transim, ///< [out] contains the magnified image. Must be pre-allocated.
572 const arrInT &im, ///< [in] is the image to be magnified.
573 transformT trans ///< [in] is the transformation to use
574)
575{
576 typedef typename transformT::arithT arithT;
577
578 arithT x0, y0, x, y;
579
580 int Nrows, Ncols;
581
582 int i0, j0;
583
584 const int lbuff = transformT::lbuff;
585 const int width = transformT::width;
586
587 Nrows = transim.rows();
588 Ncols = transim.cols();
589
590 int xulim = im.rows() - lbuff - 1;
591 int yulim = im.cols() - lbuff - 1;
592
593 arithT x_scale = ( (arithT)im.rows() - 1.0 ) / ( transim.rows() - 1.0 ); // this is 1/x_mag
594 arithT y_scale = ( (arithT)im.cols() - 1.0 ) / ( transim.cols() - 1.0 ); // this is 1/y_mag
595
596 arithT xcen = 0.5 * ( (arithT)transim.rows() - 1.0 );
597 arithT ycen = 0.5 * ( (arithT)transim.cols() - 1.0 );
598
599 arithT xcen0 = 0.5 * ( (arithT)im.rows() - 1.0 );
600 arithT ycen0 = 0.5 * ( (arithT)im.cols() - 1.0 );
601
602 // #pragma omp parallel private(x0,y0,i0,j0,x,y) num_threads(4)
603 {
604 arrOutT kern;
605 kern.resize( width, width );
606
607 for( int j = 0; j < Ncols; ++j )
608 {
609 // (i,j) is position in new image
610 // (x0,y0) is true position in old image
611 // (i0,j0) is integer position in old image
612 // (x, y) is fractional residual of (x0-i0, y0-j0)
613
614 y0 = ycen0 + ( j - ycen ) * y_scale;
615 j0 = y0;
616
617 if( j0 < lbuff || j0 >= yulim )
618 {
619 for( int i = 0; i < Nrows; ++i )
620 {
621 transim( i, j ) = 0;
622 }
623 continue;
624 }
625
626 // #pragma omp for
627 for( int i = 0; i < Nrows; ++i )
628 {
629 x0 = xcen0 + ( i - xcen ) * x_scale;
630 i0 = x0; // just converting to int
631
632 if( i0 < lbuff || i0 >= xulim )
633 {
634 transim( i, j ) = 0;
635 continue;
636 }
637
638 // Get the residual
639 x = x0 - i0;
640 y = y0 - j0;
641
642 trans( kern, x, y );
643 transim( i, j ) = ( im.block( i0 - lbuff, j0 - lbuff, width, width ) * kern ).sum();
644 } // for j
645 } // for i
646 } // #pragma omp
647}
648
649/// Magnify an image with the cubic convolution interpolator.
650/** Uses the cubic convolution interpolator to magnify the input image to the size of the output image.
651 *
652 * This is a wrapper for imageMagnify with the transform type specified.
653 *
654 * \tparam arrOutT is the eigen array type of the output.
655 * \tparam arrInT is the eigen array type of the input.
656 *
657 * \overload
658 */
659template <typename arrOutT, typename arrInT>
660void imageMagnify( arrOutT &transim, ///< [out] contains the magnified image. Must be pre-allocated.
661 const arrInT &im ///< [in] is the image to be magnified.
662)
663{
665}
666
667/// Re-bin an image using the sum, reducing its size while conserving the total flux.
668/** Optionally this can be the mean instead of the sum filter, in which case total flux is not conserved.
669 */
670template <typename imageOutT, typename imageInT>
671int imageRebinSum( imageOutT &imout, ///< [out] the re-binned image. Must be allocated to size which is an integer
672 ///< factor smaller than imin.
673 const imageInT &imin, ///< [in] the image to rebin
674 bool mean = false ///< [in] if true the output is the mean rather than the sum.
675)
676{
677 int rebin = imin.rows() / imout.rows();
678 if( imin.cols() / imout.cols() != rebin )
679 return -1;
680
681 int N = 1;
682 if( mean )
683 N = rebin * rebin;
684 for( int i = 0; i < imout.rows(); ++i )
685 {
686 for( int j = 0; j < imout.cols(); ++j )
687 {
688 imout( i, j ) = imin.block( i * rebin, j * rebin, rebin, rebin ).sum() / N;
689 }
690 }
691
692 return 0;
693}
694
695/// Re-bin an image using the sum, reducing its size while conserving the total flux. Records the value and position of
696/// the re-binned max pixel.
697/** Optionally this can be the mean instead of the sum filter, in which case total flux is not conserved.
698 *
699 * \overload
700 */
701template <typename imageOutT, typename imageInT>
702int imageRebinSum( imageOutT &imout, ///< [out] the re-binned image. Must be allocated to size which is an integer
703 ///< factor smaller than imin.
704 int &xMax, ///< [out] the x-locatioin of the max pixel
705 int &yMax, ///< [out] the y-locatioin of the max pixel
706 typename imageOutT::Scalar &pMax, ///< [out] the value of the max pixel
707 const imageInT &imin, ///< [in] the image to rebin
708 bool mean = false ///< [in] if true the output is the mean rather than the sum.
709)
710{
711 int rebin = imin.rows() / imout.rows();
712 if( imin.cols() / imout.cols() != rebin )
713 return -1;
714
715 int N = 1;
716 if( mean )
717 N = rebin * rebin;
718
719 xMax = 0;
720 yMax = 0;
721 pMax = std::numeric_limits<typename imageOutT::Scalar>::lowest();
722
723 for( int i = 0; i < imout.rows(); ++i )
724 {
725 for( int j = 0; j < imout.cols(); ++j )
726 {
727 imout( i, j ) = imin.block( i * rebin, j * rebin, rebin, rebin ).sum() / N;
728 if( imout( i, j ) > pMax )
729 {
730 pMax = imout( i, j );
731 xMax = i;
732 yMax = j;
733 }
734 }
735 }
736
737 return 0;
738}
739
740/// Re-bin an image using the mean.
741/** This is a wrapper for imageRebinSum with `mean=true`.
742 */
743template <typename imageOutT, typename imageInT>
744int imageRebinMean( imageOutT &imout, ///< [out] the re-binned image. Must be allocated to size which is an integer
745 ///< factor smaller than imin.
746 const imageInT &imin ///< [in] the image to rebin
747)
748{
749 return imageRebinSum( imout, imin, true );
750}
751
752/// Re-bin an image using the mean. Records the value and position of the re-binned max pixel.
753/** This is a wrapper for imageRebinSum with `mean=true`.
754 *
755 * \overload
756 */
757template <typename imageOutT, typename imageInT>
758int imageRebinMean( imageOutT &imout, ///< [out] the re-binned image. Must be allocated to size which is an integer
759 ///< factor smaller than imin.
760 int &xMax, ///< [out] the x-locatioin of the max pixel
761 int &yMax, ///< [out] the y-locatioin of the max pixel
762 typename imageOutT::Scalar &pMax, ///< [out] the value of the max pixel
763 const imageInT &imin, ///< [in] the image to rebin
764 bool mean = false ///< [in] if true the output is the mean rather than the sum.
765)
766{
767 return imageRebinSum( imout, xMax, yMax, pMax, imin, true );
768}
769
770/// Re-bin an image, takes the mean with a min/max rejection.
771/** The mean is calculated after rejecting the minimuma and maximum value.
772 */
773template <typename imageOutT, typename imageInT>
774int imageRebinMeanReject( imageOutT &imout, ///< [out] the re-binned image. Must be allocated to size which is an
775 ///< integer factor smaller than imin.
776 const imageInT &imin ///< [in] the image to rebin
777)
778{
779 int rebin = imin.rows() / imout.rows();
780 if( imin.cols() / imout.cols() != rebin )
781 return -1;
782
783 int N = rebin * rebin - 2;
784
785 for( int i = 0; i < imout.rows(); ++i )
786 {
787 for( int j = 0; j < imout.cols(); ++j )
788 {
789 typename imageOutT::Scalar sum = 0;
790 typename imageOutT::Scalar max = imin( i * rebin, j * rebin );
791 typename imageOutT::Scalar min = imin( i * rebin, j * rebin );
792 for( int k = 0; k < rebin; ++k )
793 {
794 for( int l = 0; l < rebin; ++l )
795 {
796 sum += imin( i * rebin + k, j * rebin + l );
797 if( imin( i * rebin + k, j * rebin + l ) > max )
798 max = imin( i * rebin + k, j * rebin + l );
799 if( imin( i * rebin + k, j * rebin + l ) < min )
800 min = imin( i * rebin + k, j * rebin + l );
801 }
802 }
803 imout( i, j ) = ( sum - max - min ) / N; /**/
804 }
805 }
806
807 return 0;
808}
809
810/// Re-bin an image, takes the mean with a min/max rejection. Records the value and position of the re-binned max
811/// pixel.
812/** The mean is calculated after rejecting the minimuma and maximum value.
813 *
814 * \overload
815 */
816template <typename imageOutT, typename imageInT>
817int imageRebinMeanReject( imageOutT &imout, ///< [out] the re-binned image. Must be allocated to size which is an
818 ///< integer factor smaller than imin.
819 int &xMax, ///< [out] the x-locatioin of the max pixel
820 int &yMax, ///< [out] the y-locatioin of the max pixel
821 typename imageOutT::Scalar &pMax, ///< [out] the value of the max pixel
822 const imageInT &imin ///< [in] the image to rebin
823)
824{
825 int rebin = imin.rows() / imout.rows();
826 if( imin.cols() / imout.cols() != rebin )
827 return -1;
828
829 int N = rebin * rebin - 2;
830
831 xMax = 0;
832 yMax = 0;
833 pMax = std::numeric_limits<typename imageOutT::Scalar>::lowest();
834
835 for( int i = 0; i < imout.rows(); ++i )
836 {
837 for( int j = 0; j < imout.cols(); ++j )
838 {
839 typename imageOutT::Scalar sum = 0;
840 typename imageOutT::Scalar max = imin( i * rebin, j * rebin );
841 typename imageOutT::Scalar min = imin( i * rebin, j * rebin );
842 for( int k = 0; k < rebin; ++k )
843 {
844 for( int l = 0; l < rebin; ++l )
845 {
846 sum += imin( i * rebin + k, j * rebin + l );
847 if( imin( i * rebin + k, j * rebin + l ) > max )
848 max = imin( i * rebin + k, j * rebin + l );
849 if( imin( i * rebin + k, j * rebin + l ) < min )
850 min = imin( i * rebin + k, j * rebin + l );
851 }
852 }
853 imout( i, j ) = ( sum - max - min ) / N;
854
855 if( imout( i, j ) > pMax )
856 {
857 pMax = imout( i, j );
858 xMax = i;
859 yMax = j;
860 }
861 }
862 }
863
864 return 0;
865}
866
867/// Down-sample an image, reducing its size while conserving the total flux.
868/** If the old size is an integer multiple of the new size, this is just a re-bin. If not an integer multiple,
869 * the image is interpolated after performing the closest re-bin, and then re-normalized to conserve flux.
870 *
871 * \todo Allow selection of interpolator, providing a default version.
872 */
873template <typename imageOutT, typename imageInT>
874void imageDownSample( imageOutT &imout, const imageInT &imin )
875{
876 typedef typename imageOutT::Scalar Scalar;
877
878 // Record this for normalization later
879 Scalar inputTotal = fabs( imin.sum() );
880
881 // As a first step, rebin to nearest whole pixel factor which is larger than the desired output size
882 int closestRebin = imin.rows() / imout.rows(); //, imin.cols()/imout.cols() );
883
884 float sample = ( (float)imin.rows() ) / closestRebin;
885
886 while( sample != floor( sample ) )
887 {
888 --closestRebin;
889 if( closestRebin == 1 )
890 break;
891 sample = ( (float)imin.rows() ) / closestRebin;
892 }
893
894 // Eigen::Array<Scalar, Eigen::Dynamic, Eigen::Dynamic> temp;
895 imageOutT temp;
896 temp.resize( imin.rows() / closestRebin, imin.cols() / closestRebin );
897
898 for( int i = 0; i < temp.rows(); ++i )
899 {
900 for( int j = 0; j < temp.cols(); ++j )
901 {
902 temp( i, j ) = imin.block( i * closestRebin, j * closestRebin, closestRebin, closestRebin ).sum();
903 }
904 }
905
906 // If the output image is now the requested size return.
907 if( temp.rows() == imout.rows() && temp.cols() == imout.cols() )
908 {
909 imout = temp;
910 return;
911 }
912 // Otherwise, re-sample using bilinear interpolation.
913 typedef bilinearTransform<Scalar> transformT;
914
915 transformT trans;
916 // Eigen::Array<Scalar, -1,-1> kern;
917 imageOutT kern;
918
919 const int lbuff = transformT::lbuff;
920 const int width = transformT::width;
921
922 for( int i = 0; i < imout.rows(); ++i )
923 {
924 for( int j = 0; j < imout.cols(); ++j )
925 {
926 double x = ( (double)i / imout.rows() ) * temp.rows();
927 double y = ( (double)j / imout.cols() ) * temp.cols();
928
929 trans( kern, x - floor( x ), y - floor( y ) );
930
931 imout( i, j ) = ( temp.block( floor( x ) - lbuff, floor( y ) - lbuff, width, width ) * kern ).sum();
932 }
933 }
934
935 // Normalize
936 Scalar outputTotal = fabs( imout.sum() );
937 imout *= inputTotal / outputTotal;
938}
939
940///@}
941
942} // namespace improc
943} // namespace mx
944
945#endif // improc_imageTransforms_hpp
bilinearTransform< double > bilinearTransd
Typedef for bilinearTransform with double precision.
void imageRotate(arrT &transim, const arrT2 &im, floatT dq, transformT trans)
Rotate an image represented as an eigen array.
cubicConvolTransform< double > cubicConvolTransd
Typedef for cubicConvolTransform with double precision.
void imageMagnify(arrOutT &transim, const arrInT &im, transformT trans)
Magnify an image.
int imageRebinMeanReject(imageOutT &imout, const imageInT &imin)
Re-bin an image, takes the mean with a min/max rejection.
void imageShift(arrOutT &transim, const arrInT &im, floatT1 dx, floatT2 dy, transformT trans)
Shift an image.
bilinearTransform< float > bilinearTransf
Typedef for bilinearTransform with single precision.
int imageRebinSum(imageOutT &imout, const imageInT &imin, bool mean=false)
Re-bin an image using the sum, reducing its size while conserving the total flux.
int imageRebinMean(imageOutT &imout, const imageInT &imin)
Re-bin an image using the mean.
cubicConvolTransform< float > cubicConvolTransf
Typedef for cubicConvolTransform with single precision.
void imageDownSample(imageOutT &imout, const imageInT &imin)
Down-sample an image, reducing its size while conserving the total flux.
void imageShiftWPScale(outputArrT &out, inputArrT &in, scaleArrT &scale, int dx, int dy)
Shift an image by whole pixels, wrapping around, with a scaling image applied to the shifted 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
Transformation by bi-linear interpolation.
Transformation by cubic convolution interpolation.
void operator()(arrT &kern, arithT x, arithT y)
arithT cubic
The kernel parameter. The default value -0.5 gives the bicubic spline interpolator.
cubicConvolTransform(const cubicConvolTransform &t)
Copy c'tor.
arithT cubicConvolKernel(arithT d)
Calculate the kernel value for a given residual.
cubicConvolTransform(arithT c)
Construct setting the kernel parameter.
_arithT arithT
The type in which all calculations are performed.