8#ifndef __imageFilters_hpp__
9#define __imageFilters_hpp__
13#include "../math/gslInterpolator.hpp"
14#include "../math/vectorUtils.hpp"
15#include "../math/geo.hpp"
28template <
typename _arrayT,
size_t _kernW = 4>
31 typedef _arrayT arrayT;
32 typedef typename _arrayT::Scalar arithT;
33 static const int kernW = _kernW;
43 int w = kernW * _fwhm;
48 kernel.resize( w, w );
50 arithT kcen = 0.5 * ( w - 1.0 );
52 arithT sig2 = _fwhm / 2.354820045030949327;
56 for(
int i = 0; i < w; ++i )
58 for(
int j = 0; j < w; ++j )
60 r2 = pow( i - kcen, 2 ) + pow( j - kcen, 2 );
61 kernel( i, j ) = exp( -r2 / ( 2.0 * sig2 ) );
65 kernel /= kernel.sum();
70 return _kernW * _fwhm;
73 void setKernel( arithT x, arithT y, arrayT &kernelArray )
76 static_cast<void>( x );
77 static_cast<void>( y );
88template <
typename _arrayT,
size_t _kernW = 2>
91 typedef _arrayT arrayT;
92 typedef typename _arrayT::Scalar arithT;
94 static const int kernW = _kernW;
120 maxAz = fabs( maxAz );
133 arithT mx1 = kernW * ( (int)( fabs(
m_azWidth * sin( qmax ) ) + fabs(
m_radWidth * cos( qmax ) ) ) + 1 );
136 arithT mx2 = kernW * ( (int)( fabs(
m_azWidth * cos( qmax ) ) + fabs(
m_radWidth * sin( qmax ) ) ) + 1 );
138 m_maxWidth = 0.5 * std::max( mx1, mx2 );
146 void setKernel( arithT x, arithT y, arrayT &kernel )
148 arithT rad0 = sqrt( (arithT)( x * x + y * y ) );
150 arithT sinq = y / rad0;
151 arithT cosq = x / rad0;
156 q = atan2( sinq, cosq );
161 if( w > m_maxWidth * 2 )
163 std::stringstream errs;
164 errs <<
"|" << kernW <<
"|" <<
m_radWidth <<
"|" <<
m_azWidth <<
"|" << m_maxWidth <<
"|" << x <<
"|" << y
165 <<
"|" << rad0 <<
"|" << sinq <<
"|" << cosq;
166 errs <<
"|" << w <<
"|" << h <<
"|";
168 "azBoxKernel::setKernel",
169 "width bigger than 2*maxWidth. This is a bug. Details: " + errs.str() );
172 if( h > m_maxWidth * 2 )
174 std::stringstream errs;
175 errs <<
"|" << kernW <<
"|" <<
m_radWidth <<
"|" <<
m_azWidth <<
"|" << m_maxWidth <<
"|" << x <<
"|" << y
176 <<
"|" << rad0 <<
"|" << sinq <<
"|" << cosq;
177 errs <<
"|" << w <<
"|" << h <<
"|";
178 mxThrowException( err::sizeerr,
179 "azBoxKernel::setKernel",
180 "height bigger than 2*maxWidth. This is a bug. Details: " + errs.str() );
183 kernel.resize( w, h );
185 arithT xcen = 0.5 * ( w - 1.0 );
186 arithT ycen = 0.5 * (
h - 1.0 );
190 arithT sinq2, cosq2, sindq, q2, dq;
191 for(
int j = 0; j <
h; ++j )
194 for(
int i = 0; i < w; ++i )
197 rad = sqrt( pow( xP - xcen, 2 ) + pow( yP - ycen, 2 ) );
198 radP = sqrt( pow( x + xP - xcen, 2 ) + pow( y + yP - ycen, 2 ) );
206 sinq2 = ( yP - ycen ) / rad;
207 cosq2 = ( xP - xcen ) / rad;
209 sindq = sinq2 * cosq - cosq2 * sinq;
214 q2 = atan2( y + yP - ycen, x + xP - xcen );
229 arithT ksum = kernel.sum();
232 mxThrowException( err::invalidconfig,
233 "azBoxKernel::setKernel",
234 "kernel sum 0 at " + std::to_string( x ) +
"," + std::to_string( y ) );
286template <
typename imageOutT,
typename imageInT,
typename kernelT>
287void filterImage( imageOutT &fim, imageInT im, kernelT kernel,
int maxr = 0 )
289 fim.resize( im.rows(), im.cols() );
291 float xcen = 0.5 * ( im.rows() - 1 );
292 float ycen = 0.5 * ( im.cols() - 1 );
295 maxr = 0.5 * im.rows() - kernel.maxWidth();
297 int mini = 0.5 * im.rows() - maxr;
298 int maxi = 0.5 * im.rows() + maxr;
299 int minj = 0.5 * im.cols() - maxr;
300 int maxj = 0.5 * im.cols() + maxr;
302 typename kernelT::arrayT kernelArray;
304#pragma omp parallel private( kernelArray )
306 int im_i, im_j, im_p, im_q;
307 int kern_i, kern_j, kern_p, kern_q;
308 typename imageOutT::Scalar norm;
311 for(
int i = 0; i < im.rows(); ++i )
313 for(
int j = 0; j < im.cols(); ++j )
315 if( ( i >= mini && i < maxi ) && ( j >= minj && j < maxj ) )
317 kernel.setKernel( i - xcen, j - ycen, kernelArray );
318 fim( i, j ) = ( im.block( i - 0.5 * ( kernelArray.rows() - 1 ),
319 j - 0.5 * ( kernelArray.cols() - 1 ),
321 kernelArray.cols() ) *
327 kernel.setKernel( i - xcen, j - ycen, kernelArray );
329 im_i = i - 0.5 * ( kernelArray.rows() - 1 );
333 im_j = j - 0.5 * ( kernelArray.cols() - 1 );
337 im_p = im.rows() - im_i;
338 if( im_p > kernelArray.rows() )
339 im_p = kernelArray.rows();
341 im_q = im.cols() - im_j;
342 if( im_q > kernelArray.cols() )
343 im_q = kernelArray.cols();
345 kern_i = 0.5 * ( kernelArray.rows() - 1 ) - i;
349 kern_j = 0.5 * ( kernelArray.cols() - 1 ) - j;
353 kern_p = kernelArray.rows() - kern_i;
354 if( kern_p > kernelArray.rows() )
355 kern_p = kernelArray.rows();
357 kern_q = kernelArray.cols() - kern_j;
358 if( kern_q > kernelArray.cols() )
359 kern_q = kernelArray.cols();
367 norm = kernelArray.block( kern_i, kern_j, kern_p, kern_q ).sum();
370 ( im.block( im_i, im_j, kern_p, kern_q ) * kernelArray.block( kern_i, kern_j, kern_p, kern_q ) )
373 if( !std::isfinite( fim( i, j ) ) )
287void filterImage( imageOutT &fim, imageInT im, kernelT kernel,
int maxr = 0 ) {
…}
381template <
typename imageOutT,
typename imageInT,
typename kernelT>
382void medianFilterImage( imageOutT &fim, imageInT im, kernelT kernel,
int maxr = 0,
int maxrproc = 1 )
384 fim.resize( im.rows(), im.cols() );
386 float xcen = 0.5 * ( im.rows() - 1 );
387 float ycen = 0.5 * ( im.cols() - 1 );
390 maxr = 0.5 * im.rows() - kernel.maxWidth();
392 int mini = 0.5 * im.rows() - maxr;
393 int maxi = 0.5 * im.rows() + maxr;
394 int minj = 0.5 * im.cols() - maxr;
395 int maxj = 0.5 * im.cols() + maxr;
397 typename kernelT::arrayT kernelArray;
399#pragma omp parallel private( kernelArray )
401 int im_i, im_j, im_p, im_q;
402 int kern_i, kern_j, kern_p, kern_q;
403 typename imageOutT::Scalar norm;
405 std::vector<typename imageOutT::Scalar> pixels;
408 for(
int i = 0; i < im.rows(); ++i )
410 for(
int j = 0; j < im.cols(); ++j )
412 if( ( i >= mini && i < maxi ) && ( j >= minj && j < maxj ) )
414 kernel.setKernel( i - xcen, j - ycen, kernelArray );
417 for(
int cc = 0; cc < kernelArray.cols(); ++cc )
419 for(
int rr = 0; rr < kernelArray.rows(); ++rr )
421 if( kernelArray( rr, cc ) != 0 )
422 pixels.push_back( im.block( i - 0.5 * ( kernelArray.rows() - 1 ),
423 j - 0.5 * ( kernelArray.cols() - 1 ),
425 kernelArray.cols() )( rr, cc ) );
438 kernel.setKernel( i - xcen, j - ycen, kernelArray );
440 im_i = i - 0.5 * ( kernelArray.rows() - 1 );
444 im_j = j - 0.5 * ( kernelArray.cols() - 1 );
448 im_p = im.rows() - im_i;
449 if( im_p > kernelArray.rows() )
450 im_p = kernelArray.rows();
452 im_q = im.cols() - im_j;
453 if( im_q > kernelArray.cols() )
454 im_q = kernelArray.cols();
456 kern_i = 0.5 * ( kernelArray.rows() - 1 ) - i;
460 kern_j = 0.5 * ( kernelArray.cols() - 1 ) - j;
464 kern_p = kernelArray.rows() - kern_i;
465 if( kern_p > kernelArray.rows() )
466 kern_p = kernelArray.rows();
468 kern_q = kernelArray.cols() - kern_j;
469 if( kern_q > kernelArray.cols() )
470 kern_q = kernelArray.cols();
478 norm = kernelArray.block( kern_i, kern_j, kern_p, kern_q ).sum();
481 for(
int cc = 0; cc < kern_q; ++cc )
483 for(
int rr = 0; rr < kern_p; ++rr )
485 if( kernelArray.block( kern_i, kern_j, kern_p, kern_q )( rr, cc ) != 0 )
487 pixels.push_back( im.block( im_i, im_j, kern_p, kern_q )( rr, cc ) );
522template <
typename imageTout,
typename imageTin>
524 const imageTin &imIn,
526 bool rejectMinMax =
false
529 typedef typename imageTout::Scalar scalarT;
531 int buff = 0.5 * meanFullWidth;
533 int nPix = meanFullWidth * meanFullWidth;
538 for(
int jj = buff; jj < imIn.cols() - buff; ++jj )
540 for(
int ii = buff; ii < imIn.rows() - buff; ++ii )
545 for(
int ll = jj - buff; ll < jj + buff + 1; ++ll )
547 for(
int kk = ii - buff; kk < ii + buff + 1; ++kk )
550 sum += imIn( kk, ll );
551 if( imIn( kk, ll ) > max )
552 max = imIn( kk, ll );
553 if( imIn( kk, ll ) < min )
554 min = imIn( kk, ll );
557 imOut( ii, jj ) = ( sum - max - min ) / nPix;
563 for(
int jj = buff; jj < imIn.cols() - buff; ++jj )
565 for(
int ii = buff; ii < imIn.rows() - buff; ++ii )
568 for(
int ll = jj - buff; ll < jj + buff + 1; ++ll )
570 for(
int kk = ii - buff; kk < ii + buff + 1; ++kk )
572 sum += imIn( kk, ll );
575 imOut( ii, jj ) = sum / nPix;
607template <
typename imageTout,
typename imageTin>
611 typename imageTout::Scalar &pMax,
612 const imageTin &imIn,
614 bool rejectMinMax =
false
617 typedef typename imageTout::Scalar scalarT;
619 int buff = 0.5 * meanFullWidth;
621 int nPix = meanFullWidth * meanFullWidth;
623 pMax = std::numeric_limits<scalarT>::lowest();
628 for(
int jj = buff; jj < imIn.cols() - buff; ++jj )
630 for(
int ii = buff; ii < imIn.rows() - buff; ++ii )
635 for(
int ll = jj - buff; ll < jj + buff + 1; ++ll )
637 for(
int kk = ii - buff; kk < ii + buff + 1; ++kk )
640 sum += imIn( kk, ll );
641 if( imIn( kk, ll ) > max )
642 max = imIn( kk, ll );
643 if( imIn( kk, ll ) < min )
644 min = imIn( kk, ll );
647 imOut( ii, jj ) = ( sum - max - min ) / nPix;
648 if( imOut( ii, jj ) > pMax )
650 pMax = imOut( ii, jj );
659 for(
int jj = buff; jj < imIn.cols() - buff; ++jj )
661 for(
int ii = buff; ii < imIn.rows() - buff; ++ii )
664 for(
int ll = jj - buff; ll < jj + buff + 1; ++ll )
666 for(
int kk = ii - buff; kk < ii + buff + 1; ++kk )
668 sum += imIn( kk, ll );
671 imOut( ii, jj ) = sum / nPix;
672 if( imOut( ii, jj ) > pMax )
674 pMax = imOut( ii, jj );
705template <
typename imageTout,
typename imageTin>
710 typename imageTout::Scalar &pMax,
711 const imageTin &imIn,
715 typedef typename imageTout::Scalar scalarT;
717 int buff = 0.5 * medianFullWidth;
719 std::vector<scalarT> pixs( medianFullWidth * medianFullWidth );
721 pMax = std::numeric_limits<scalarT>::lowest();
723 for(
int jj = buff; jj < imIn.cols() - buff; ++jj )
725 for(
int ii = buff; ii < imIn.rows() - buff; ++ii )
728 for(
int ll = jj - buff; ll < jj + buff + 1; ++ll )
730 for(
int kk = ii - buff; kk < ii + buff + 1; ++kk )
732 pixs[n] = imIn( kk, ll );
738 if( imOut( ii, jj ) > pMax )
740 pMax = imOut( ii, jj );
767template <
typename imageTout,
typename imageTin>
770 const imageTin &imIn,
775 typename imageTout::Scalar pMax;
776 return medianSmooth( imOut, xMax, yMax, pMax, imIn, medianFullWidth );
779template <
typename eigenImT>
784 typedef typename eigenImT::Scalar realT;
786 std::vector<realT> edge( 2 * ncols );
787 for(
int rr = 0; rr < im.rows(); ++rr )
789 for(
int cc = 0; cc < ncols; ++cc )
791 edge[cc] = im( rr, cc );
794 for(
int cc = 0; cc < ncols; ++cc )
796 edge[ncols + cc] = im( rr, im.cols() - ncols + cc );
806template <
typename eigenImT>
811 typedef typename eigenImT::Scalar realT;
813 std::vector<realT> edge( 2 * nrows );
814 for(
int cc = 0; cc < im.cols(); ++cc )
816 for(
int rr = 0; rr < nrows; ++rr )
818 edge[rr] = im( rr, cc );
821 for(
int rr = 0; rr < nrows; ++rr )
823 edge[nrows + rr] = im( im.rows() - nrows + rr, cc );
835template <
typename floatT>
842template <
typename floatT>
845 bool operator()( radval<floatT> rv1, radval<floatT> rv2 )
847 return ( rv1.r < rv2.r );
851template <
typename floatT>
854 bool operator()( radval<floatT> rv1, radval<floatT> rv2 )
856 return ( rv1.v < rv2.v );
871template <
typename vecT,
typename eigenImT1,
typename eigenImT2,
typename eigenImT3>
876 const eigenImT2 &radim,
877 const eigenImT3 *mask,
880 typename eigenImT1::Scalar minr = 0 )
882 typedef typename eigenImT1::Scalar floatT;
884 int dim1 = im.rows();
885 int dim2 = im.cols();
901 std::vector<radval<floatT>> rv( nPix );
905 for(
int c = 0; c < im.cols(); ++c )
907 for(
int r = 0; r < im.rows(); ++r )
911 if( ( *mask )( r, c ) == 0 )
915 rv[i].r = radim( r, c );
916 rv[i].v = im( r, c );
921 sort( rv.begin(), rv.end(), radvalRadComp<floatT>() );
937 floatT r1 = minr + dr;
944 while( rv[i1].r < r0 )
947 while( rv[i2].r <= r1 )
953 for(
int in = i1; in < i2; ++in )
959 n = 0.5 * ( i2 - i1 );
961 std::nth_element( rv.begin() + i1, rv.begin() + i1 + n, rv.begin() + i2, radvalValComp<floatT>() );
963 med = ( rv.begin() + i1 + n )->v;
966 if( ( i2 - i1 ) % 2 == 0 )
970 ( med + ( *std::max_element( rv.begin() + i1, rv.begin() + i1 + n, radvalValComp<floatT>() ) ).v );
974 rad.push_back( .5 * ( r0 + r1 ) );
975 prof.push_back( med );
995template <
typename vecT,
typename eigenImT1,
typename eigenImT2>
1000 const eigenImT2 &mask,
1005 radim.resize( im.cols(), im.rows() );
1009 radprof( rad, prof, im, radim, &mask, mean );
1023template <
typename vecT,
typename eigenImT1>
1027 const eigenImT1 &im,
1032 radim.resize( im.cols(), im.rows() );
1050template <
typename radprofT,
typename eigenImT1,
typename eigenImT2,
typename eigenImT3>
1053 const eigenImT2 &rad,
1054 const eigenImT3 *mask,
1062 std::vector<double> med_r, med_v;
1064 radprof( med_r, med_v, im, rad, mask );
1067 radprofIm.resize( im.rows(), im.cols() );
1071 for(
int c = 0; c < im.cols(); ++c )
1073 for(
int r = 0; r < im.rows(); ++r )
1077 if( ( *mask )( r, c ) == 0 )
1079 radprofIm( r, c ) = 0;
1084 radprofIm( r, c ) = interp( ( (
double)rad( r, c ) ) );
1086 im( r, c ) -= radprofIm( r, c );
1100template <
typename radprofT,
typename eigenImT>
1104 bool subtract =
false,
1109 rad.resize( im.rows(), im.cols() );
1127template <
typename eigenImT,
typename eigenImT1,
typename eigenImT2,
typename eigenImT3>
1129 const eigenImT1 &im,
1130 const eigenImT2 &rad,
1131 const eigenImT3 &mask,
1132 typename eigenImT::Scalar minRad,
1133 typename eigenImT::Scalar maxRad,
1138 typedef typename eigenImT::Scalar floatT;
1140 int dim1 = im.cols();
1141 int dim2 = im.rows();
1143 floatT mr = rad.maxCoeff();
1146 std::vector<radval<floatT>> rv( dim1 * dim2 );
1148 for(
int i = 0; i < rv.size(); ++i )
1150 if( mask( i ) == 0 )
1157 sort( rv.begin(), rv.end(), radvalRadComp<floatT>() );
1167 std::vector<double> std_r, std_v;
1171 while( rv[i1].r < r0 && i1 < rv.size() )
1175 if( i1 == rv.size() )
1181 while( rv[i2].r <= r1 && i2 < rv.size() )
1183 n = 0.5 * ( i2 - i1 );
1185 std::vector<double> vals;
1187 for(
int i = i1; i < i2; ++i )
1189 vals.push_back( rv[i].v );
1192 std_r.push_back( .5 * ( r0 + r1 ) );
1201 stdIm.resize( dim1, dim2 );
1204 for(
int i = 0; i < dim1; ++i )
1206 for(
int j = 0; j < dim2; ++j )
1208 if( rad( i, j ) < minRad || rad( i, j ) > maxRad )
1214 stdIm( i, j ) = interp( ( (
double)rad( i, j ) ) );
1216 stdIm( i, j ) = im( i, j ) / stdIm( i, j );
1233template <
typename eigenImT,
typename eigenImT1,
typename eigenImT2>
1235 const eigenImT1 &im,
1236 const eigenImT2 &mask,
1237 typename eigenImT::Scalar minRad,
1238 typename eigenImT::Scalar maxRad,
1243 int dim1 = im.cols();
1244 int dim2 = im.rows();
1247 rad.resize( dim1, dim2 );
1250 stddevImage( stdIm, im, rad, mask, minRad, maxRad, divide );
1262template <
typename eigenCubeT,
typename eigenCubeT1,
typename eigenCubeT2,
typename radT1,
typename radT2>
1265 const eigenCubeT1 &imc,
1266 const eigenCubeT2 &maskCube,
1273 int dim1 = imc.cols();
1274 int dim2 = imc.rows();
1276 typename eigenCubeT::Scalar minRadF = minRad;
1277 typename eigenCubeT::Scalar maxRadF = maxRad;
1280 rad.resize( dim1, dim2 );
1284 stdImc.resize( imc.rows(), imc.cols(), imc.planes() );
1287 for(
int i = 0; i < imc.planes(); ++i )
1291 im = imc.image( i );
1292 mask = maskCube.image( i );
1294 stddevImage( stdIm, im, rad, mask, minRadF, maxRadF, divide );
1296 stdImc.image( i ) = stdIm;
mxException for a size error
Class to manage interpolation using the GSL interpolation library.
constexpr units::realT h()
Planck Constant.
Eigen::Array< scalarT, -1, -1 > eigenImage
Definition of the eigenImage type, which is an alias for Eigen::Array.
constexpr floatT six_fifths()
Return 6/5 in the specified precision.
realT dtor(realT q)
Convert from degrees to radians.
int meanSmooth(imageTout &imOut, const imageTin &imIn, int meanFullWidth, bool rejectMinMax=false)
Smooth an image using the mean in a rectangular box, optionally rejecting the highest and lowest valu...
int medianSmooth(imageTout &imOut, int &xMax, int &yMax, typename imageTout::Scalar &pMax, const imageTin &imIn, int medianFullWidth)
void filterImage(imageOutT &fim, imageInT im, kernelT kernel, int maxr=0)
Filter an image with a mean kernel.
void radiusImage(eigenT &m, typename eigenT::Scalar xc, typename eigenT::Scalar yc, typename eigenT::Scalar scale=1)
Fills in the cells of an Eigen 2D Array with their radius from the center.
void radprof(vecT &rad, vecT &prof, const eigenImT1 &im, const eigenImT2 &radim, const eigenImT3 *mask, bool mean=false, typename eigenImT1::Scalar minr=0)
Calculate the the radial profile.
void radprofim(radprofT &radprofIm, eigenImT1 &im, const eigenImT2 &rad, const eigenImT3 *mask, bool subtract, bool mean=false)
Form a radial profile image, and optionally subtract it from the input.
void stddevImage(eigenImT &stdIm, const eigenImT1 &im, const eigenImT2 &rad, const eigenImT3 &mask, typename eigenImT::Scalar minRad, typename eigenImT::Scalar maxRad, bool divide)
Form a standard deviation image, and optionally divide the input by it to form a S/N map.
void stddevImageCube(eigenCubeT &stdImc, const eigenCubeT1 &imc, const eigenCubeT2 &maskCube, radT1 minRad, radT2 maxRad, bool divide=false)
vectorT::value_type vectorMedianInPlace(vectorT &vec)
Calculate median of a vector in-place, altering the vector.
valueT vectorVariance(const valueT *vec, size_t sz, valueT mean)
Calculate the variance of a vector relative to a supplied mean value.
vectorT::value_type vectorMedian(const vectorT &vec, vectorT *work=0)
Calculate median of a vector, leaving the vector unaltered.
void colEdgeMedSubtract(eigenImT &im, int nrows)
void rowEdgeMedSubtract(eigenImT &im, int ncols)
Declares and defines functions to work with image masks.
Azimuthally variable boxcare kernel.
void setMaxWidth()
Sets the max width based on the configured az and rad widths.
azBoxKernel(arithT radWidth, arithT azWidth, arithT maxAz)
arithT m_radWidth
the half-width of the averaging box, in the radial direction, in pixels.
arithT m_azWidth
the half-width of the averaging box, in the azimuthal direction, in pixels.
azBoxKernel(arithT radWidth, arithT azWidth)
Symetric Gaussian smoothing kernel.