27#ifndef improc_imageTransforms_hpp
28#define improc_imageTransforms_hpp
44template <
typename _arithT>
47 typedef _arithT arithT;
49 static const size_t width = 2;
50 static const size_t lbuff = 0;
52 template <
typename arrT,
typename arithT>
53 void operator()( arrT &kern, arithT x, arithT y )
55 kern.resize( width, width );
57 kern( 0, 0 ) = ( 1. - x ) * ( 1. - y );
58 kern( 0, 1 ) = ( 1. - x ) * y;
59 kern( 1, 0 ) = x * ( 1. - y );
89template <
typename _arithT>
94 static const size_t width = 4;
95 static const size_t lbuff = 1;
123 return (
cubic + 2. ) * d * d * d - (
cubic + 3. ) * d * d + 1.;
132 template <
typename arrT,
typename arithT>
135 arithT km2x, km1x, kp1x, kp2x;
136 arithT km2y, km1y, kp1y, kp2y;
148 kern( 0, 0 ) = km2x * km2y;
149 kern( 0, 1 ) = km2x * km1y;
150 kern( 0, 2 ) = km2x * kp1y;
151 kern( 0, 3 ) = km2x * kp2y;
153 kern( 1, 0 ) = km1x * km2y;
154 kern( 1, 1 ) = km1x * km1y;
155 kern( 1, 2 ) = km1x * kp1y;
156 kern( 1, 3 ) = km1x * kp2y;
158 kern( 2, 0 ) = kp1x * km2y;
159 kern( 2, 1 ) = kp1x * km1y;
160 kern( 2, 2 ) = kp1x * kp1y;
161 kern( 2, 3 ) = kp1x * kp2y;
163 kern( 3, 0 ) = kp2x * km2y;
164 kern( 3, 1 ) = kp2x * km1y;
165 kern( 3, 2 ) = kp2x * kp1y;
166 kern( 3, 3 ) = kp2x * kp2y;
193template <
typename transformT,
typename arrT,
typename arrT2,
typename floatT>
200 typedef typename transformT::arithT arithT;
209 const int lbuff = transformT::lbuff;
210 const int width = transformT::width;
218 transim.resize( Nrows, Ncols );
221 xcen = 0.5 * ( Nrows - 1. );
222 ycen = 0.5 * ( Ncols - 1. );
224 int xulim = Nrows - width + lbuff;
225 int yulim = Ncols - width + lbuff;
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;
232 xc_x_cosq += yc_x_sinq;
233 xc_x_sinq -= yc_x_cosq;
236#pragma omp parallel private( x0, y0, i0, j0, x, y )
239 arithT i_x_cosq, i_x_sinq;
241 kern.resize( width, width );
244#pragma omp for schedule( static, 1 )
246 for(
int i = 0; i < Nrows; ++i )
248 i_x_cosq = i * cosq - xc_x_cosq;
249 i_x_sinq = -( i * sinq - xc_x_sinq );
251 for(
int j = 0; j < Ncols; ++j )
257 x0 = i_x_cosq + j * sinq;
258 y0 = i_x_sinq + j * cosq;
264 if( i0 <= lbuff || i0 >= xulim || j0 <= lbuff || j0 >= yulim )
275 transim( i, j ) = ( im.block( i0 - lbuff, j0 - lbuff, width, width ) * kern ).sum();
293template <
typename outputArrT,
typename inputArrT>
305 int outr = out.rows();
306 int outc = out.cols();
321 for(
int cc = 0; cc < outc; ++cc )
330 for(
int rr = 0; rr < outr; ++rr )
339 out( rr, cc ) = in( x, y );
348 for(
int cc = 0; cc < outc; ++cc )
352 if( y < 0 || y >= inc )
354 for(
int rr = 0; rr < outr; ++rr )
362 for(
int rr = 0; rr < outr; ++rr )
366 if( x < 0 || x >= inr )
372 out( rr, cc ) = in( x, y );
391template <
typename outputArrT,
typename inputArrT,
typename scaleArrT>
404 int outr = out.rows();
405 int outc = out.cols();
417 for(
int cc = 0; cc < outc; ++cc )
426 for(
int rr = 0; rr < outr; ++rr )
435 out( rr, cc ) = in( x, y ) * scale( rr, cc );
461template <
typename arrOutT,
typename arrInT,
typename floatT1,
typename floatT2,
typename transformT>
469 typedef typename transformT::arithT arithT;
473 const int lbuff = transformT::lbuff;
474 const int width = transformT::width;
477 if( dx == floor( dx ) && dy == floor( dy ) )
483 int xulim = Nrows - width + lbuff;
484 int yulim = Ncols - width + lbuff;
486 transim.resize( Nrows, Ncols );
494 arithT rx = 1 - ( dx - floor( dx ) );
495 arithT ry = 1 - ( dy - floor( dy ) );
498 kern.resize( width, width );
499 trans( kern, rx, ry );
504 for(
int i = 0; i < Nrows; ++i )
511 if( i0 <= lbuff || i0 >= xulim )
513 for(
int j = 0; j < Ncols; ++j )
520 for(
int j = 0; j < Ncols; ++j )
524 if( j0 <= lbuff || j0 >= yulim )
530 transim( i, j ) = ( im.block( i0 - lbuff, j0 - lbuff, width, width ) * kern ).sum();
570template <
typename arrOutT,
typename arrInT,
typename transformT>
576 typedef typename transformT::arithT arithT;
584 const int lbuff = transformT::lbuff;
585 const int width = transformT::width;
587 Nrows = transim.rows();
588 Ncols = transim.cols();
590 int xulim = im.rows() - lbuff - 1;
591 int yulim = im.cols() - lbuff - 1;
593 arithT x_scale = ( (arithT)im.rows() - 1.0 ) / ( transim.rows() - 1.0 );
594 arithT y_scale = ( (arithT)im.cols() - 1.0 ) / ( transim.cols() - 1.0 );
596 arithT xcen = 0.5 * ( (arithT)transim.rows() - 1.0 );
597 arithT ycen = 0.5 * ( (arithT)transim.cols() - 1.0 );
599 arithT xcen0 = 0.5 * ( (arithT)im.rows() - 1.0 );
600 arithT ycen0 = 0.5 * ( (arithT)im.cols() - 1.0 );
605 kern.resize( width, width );
607 for(
int j = 0; j < Ncols; ++j )
614 y0 = ycen0 + ( j - ycen ) * y_scale;
617 if( j0 < lbuff || j0 >= yulim )
619 for(
int i = 0; i < Nrows; ++i )
627 for(
int i = 0; i < Nrows; ++i )
629 x0 = xcen0 + ( i - xcen ) * x_scale;
632 if( i0 < lbuff || i0 >= xulim )
643 transim( i, j ) = ( im.block( i0 - lbuff, j0 - lbuff, width, width ) * kern ).sum();
659template <
typename arrOutT,
typename arrInT>
670template <
typename imageOutT,
typename imageInT>
673 const imageInT &imin,
677 int rebin = imin.rows() / imout.rows();
678 if( imin.cols() / imout.cols() != rebin )
684 for(
int i = 0; i < imout.rows(); ++i )
686 for(
int j = 0; j < imout.cols(); ++j )
688 imout( i, j ) = imin.block( i * rebin, j * rebin, rebin, rebin ).sum() / N;
701template <
typename imageOutT,
typename imageInT>
706 typename imageOutT::Scalar &pMax,
707 const imageInT &imin,
711 int rebin = imin.rows() / imout.rows();
712 if( imin.cols() / imout.cols() != rebin )
721 pMax = std::numeric_limits<typename imageOutT::Scalar>::lowest();
723 for(
int i = 0; i < imout.rows(); ++i )
725 for(
int j = 0; j < imout.cols(); ++j )
727 imout( i, j ) = imin.block( i * rebin, j * rebin, rebin, rebin ).sum() / N;
728 if( imout( i, j ) > pMax )
730 pMax = imout( i, j );
743template <
typename imageOutT,
typename imageInT>
757template <
typename imageOutT,
typename imageInT>
762 typename imageOutT::Scalar &pMax,
763 const imageInT &imin,
773template <
typename imageOutT,
typename imageInT>
779 int rebin = imin.rows() / imout.rows();
780 if( imin.cols() / imout.cols() != rebin )
783 int N = rebin * rebin - 2;
785 for(
int i = 0; i < imout.rows(); ++i )
787 for(
int j = 0; j < imout.cols(); ++j )
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 )
794 for(
int l = 0; l < rebin; ++l )
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 );
803 imout( i, j ) = ( sum - max - min ) / N;
816template <
typename imageOutT,
typename imageInT>
821 typename imageOutT::Scalar &pMax,
825 int rebin = imin.rows() / imout.rows();
826 if( imin.cols() / imout.cols() != rebin )
829 int N = rebin * rebin - 2;
833 pMax = std::numeric_limits<typename imageOutT::Scalar>::lowest();
835 for(
int i = 0; i < imout.rows(); ++i )
837 for(
int j = 0; j < imout.cols(); ++j )
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 )
844 for(
int l = 0; l < rebin; ++l )
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 );
853 imout( i, j ) = ( sum - max - min ) / N;
855 if( imout( i, j ) > pMax )
857 pMax = imout( i, j );
873template <
typename imageOutT,
typename imageInT>
876 typedef typename imageOutT::Scalar Scalar;
879 Scalar inputTotal = fabs( imin.sum() );
882 int closestRebin = imin.rows() / imout.rows();
884 float sample = ( (float)imin.rows() ) / closestRebin;
886 while( sample != floor( sample ) )
889 if( closestRebin == 1 )
891 sample = ( (float)imin.rows() ) / closestRebin;
896 temp.resize( imin.rows() / closestRebin, imin.cols() / closestRebin );
898 for(
int i = 0; i < temp.rows(); ++i )
900 for(
int j = 0; j < temp.cols(); ++j )
902 temp( i, j ) = imin.block( i * closestRebin, j * closestRebin, closestRebin, closestRebin ).sum();
907 if( temp.rows() == imout.rows() && temp.cols() == imout.cols() )
919 const int lbuff = transformT::lbuff;
920 const int width = transformT::width;
922 for(
int i = 0; i < imout.rows(); ++i )
924 for(
int j = 0; j < imout.cols(); ++j )
926 double x = ( (double)i / imout.rows() ) * temp.rows();
927 double y = ( (double)j / imout.cols() ) * temp.cols();
929 trans( kern, x - floor( x ), y - floor( y ) );
931 imout( i, j ) = ( temp.block( floor( x ) - lbuff, floor( y ) - lbuff, width, width ) * kern ).sum();
936 Scalar outputTotal = fabs( imout.sum() );
937 imout *= inputTotal / outputTotal;