13#pragma GCC system_header
16#include "../math/vectorUtils.hpp"
28template <
typename dataT>
32 typedef bool is_eigenCube;
35 typedef typename Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>::Index Index;
37 typedef Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>> imageRef;
39 typedef Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic> imageT;
46 dataT *m_data{
nullptr };
82 void resize(
int r,
int c,
int p );
84 void resize(
int r,
int c );
88 const dataT *data()
const;
90 const Index rows()
const;
92 const Index cols()
const;
94 const Index planes()
const;
97 Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>>
cube();
100 Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>>
image( Index n );
103 const Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>>
107 Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>,
109 Eigen::Stride<Eigen::Dynamic, Eigen::Dynamic>>
113 Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>>
asVectors();
116 void Covar( Eigen::Matrix<dataT, Eigen::Dynamic, Eigen::Dynamic> &cv );
122 template <
typename eigenT>
129 template <
typename eigenT>
137 template <
typename eigenT,
typename eigenCubeT>
141 double minGoodFract = 0.0
149 template <
typename eigenT>
151 std::vector<dataT> &weights
159 template <
typename eigenT,
typename eigenCubeT>
161 std::vector<dataT> &weights,
163 double minGoodFract =
171 template <
typename eigenT>
178 template <
typename eigenT>
188 template <
typename eigenT,
typename eigenCubeT>
193 double minGoodFract =
201 template <
typename eigenT>
203 std::vector<dataT> &weights,
212 template <
typename eigenT,
typename eigenCubeT>
215 std::vector<dataT> &weights,
218 double minGoodFract =
223template <
typename dataT>
232template <
typename dataT>
239 m_data =
new dataT[_rows * _cols * _planes];
243template <
typename dataT>
254template <
typename dataT>
257 if( _owner && m_data )
263template <
typename dataT>
264void eigenCube<dataT>::setZero()
266 int N = _rows * _cols * _planes;
268 for(
int i = 0; i < N; ++i )
269 m_data[i] = ( (dataT)0 );
272template <
typename dataT>
273eigenCube<dataT> &eigenCube<dataT>::operator=(
const eigenCube<dataT> &ec )
275 resize( ec.rows(), ec.cols(), ec.planes() );
277 int N = _rows * _cols * _planes;
279 for(
int i = 0; i < N; ++i )
280 m_data[i] = ec.m_data[i];
285template <
typename dataT>
286void eigenCube<dataT>::shallowCopy( eigenCube<dataT> &src,
bool takeOwner )
288 if( _owner && m_data )
296 _planes = src._planes;
299 if( takeOwner ==
true )
310template <
typename dataT>
313 if( _owner && m_data )
323template <
typename dataT>
326 if( _owner && m_data )
335 m_data =
new dataT[_rows * _cols * _planes];
339template <
typename dataT>
340void eigenCube<dataT>::resize(
int r,
int c )
345template <
typename dataT>
346dataT *eigenCube<dataT>::data()
351template <
typename dataT>
352const dataT *eigenCube<dataT>::data()
const
357template <
typename dataT>
358const typename eigenCube<dataT>::Index eigenCube<dataT>::rows()
const
363template <
typename dataT>
364const typename eigenCube<dataT>::Index eigenCube<dataT>::cols()
const
369template <
typename dataT>
370const typename eigenCube<dataT>::Index eigenCube<dataT>::planes()
const
375template <
typename dataT>
378 return Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>>( m_data, _rows * _cols, _planes );
381template <
typename dataT>
384 return Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>>( m_data + n * _rows * _cols, _rows, _cols );
387template <
typename dataT>
390 return Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>>( m_data + n * _rows * _cols, _rows, _cols );
393template <
typename dataT>
394Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>,
396 Eigen::Stride<Eigen::Dynamic, Eigen::Dynamic>>
399 return Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>,
401 Eigen::Stride<Eigen::Dynamic, Eigen::Dynamic>>(
402 m_data + j * _rows + i,
405 Eigen::Stride<Eigen::Dynamic, Eigen::Dynamic>( 0, _rows * _cols ) );
408template <
typename dataT>
411 return Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>>( m_data, _rows * _cols, _planes );
414template <
typename dataT>
417 cv = asVectors().matrix().transpose() * asVectors().matrix();
420template <
typename dataT>
421template <
typename eigenT>
424 mim.resize( _rows, _cols );
426#pragma omp parallel for schedule( static, 1 ) num_threads( Eigen::nbThreads() )
427 for( Index i = 0; i < _rows; ++i )
429 for( Index j = 0; j < _cols; ++j )
431 mim( i, j ) = pixel( i, j ).sum();
436template <
typename dataT>
437template <
typename eigenT>
440 mim.resize( _rows, _cols );
442#pragma omp parallel for schedule( static, 1 ) num_threads( Eigen::nbThreads() )
443 for( Index i = 0; i < _rows; ++i )
445 for( Index j = 0; j < _cols; ++j )
447 mim( i, j ) = pixel( i, j ).mean();
452template <
typename dataT>
453template <
typename eigenT,
typename eigenCubeT>
456 mim.resize( _rows, _cols );
460 std::vector<dataT> work;
463 for( Index i = 0; i < _rows; ++i )
465 for( Index j = 0; j < _cols; ++j )
470 for( Index k = 0; k < _planes; ++k )
472 if( ( mask.pixel( i, j ) )( k, 0 ) == 1 )
474 work.push_back( ( pixel( i, j ) )( k, 0 ) );
478 if( work.size() > minGoodFract * _planes )
484 mim( i, j ) = invalidNumber<float>();
491template <
typename dataT>
492template <
typename eigenT>
495 mim.resize( _rows, _cols );
497#pragma omp parallel num_threads( Eigen::nbThreads() )
499 std::vector<Scalar> work;
501#pragma omp for schedule( static, 10 )
502 for( Index i = 0; i < _rows; ++i )
504 for( Index j = 0; j < _cols; ++j )
506 work.resize( _planes );
508 for(
int k = 0; k < _planes; ++k )
509 work[k] = ( pixel( i, j ) )( k, 0 );
517template <
typename dataT>
518template <
typename eigenT,
typename eigenCubeT>
521 mim.resize( _rows, _cols );
523#pragma omp parallel num_threads( Eigen::nbThreads() )
525 std::vector<Scalar> work, wwork;
527#pragma omp for schedule( static, 10 )
528 for( Index i = 0; i < _rows; ++i )
530 for( Index j = 0; j < _cols; ++j )
537 for( Index k = 0; k < _planes; ++k )
539 if( ( mask.pixel( i, j ) )( k, 0 ) == 1 )
541 work.push_back( ( pixel( i, j ) )( k, 0 ) );
542 wwork.push_back( weights[k] );
545 if( work.size() > minGoodFract * _planes )
550 mim( i, j ) = invalidNumber<float>();
556template <
typename dataT>
557template <
typename eigenT>
560 mim.resize( _rows, _cols );
562#pragma omp parallel for schedule( static, 10 ) num_threads( Eigen::nbThreads() )
563 for( Index i = 0; i < _rows; ++i )
566 std::vector<Scalar> work;
567 for( Index j = 0; j < _cols; ++j )
574template <
typename dataT>
575template <
typename eigenT>
578 mim.resize( _rows, _cols );
583 std::vector<Scalar> work;
587 for( Index i = 0; i < _rows; ++i )
589 for( Index j = 0; j < _cols; ++j )
591 work.resize( _planes );
592 for(
int k = 0; k < _planes; ++k )
594 work[k] = ( pixel( i, j ) )( k, 0 );
603template <
typename dataT>
604template <
typename eigenT,
typename eigenCubeT>
607 mim.resize( _rows, _cols );
611 std::vector<Scalar> work;
614 for( Index i = 0; i < _rows; ++i )
616 for( Index j = 0; j < _cols; ++j )
621 for( Index k = 0; k < _planes; ++k )
623 if( ( mask.pixel( i, j ) )( k, 0 ) == 1 )
625 work.push_back( ( pixel( i, j ) )( k, 0 ) );
629 if( work.size() > minGoodFract * _planes )
635 mim( i, j ) = invalidNumber<float>();
642template <
typename dataT>
643template <
typename eigenT>
646 mim.resize( _rows, _cols );
648#pragma omp parallel num_threads( Eigen::nbThreads() )
650 std::vector<Scalar> work;
652#pragma omp for schedule( static, 10 )
653 for( Index i = 0; i < _rows; ++i )
655 for( Index j = 0; j < _cols; ++j )
657 work.resize( _planes );
659 for(
int k = 0; k < _planes; ++k )
660 work[k] = ( pixel( i, j ) )( k, 0 );
668template <
typename dataT>
669template <
typename eigenT,
typename eigenCubeT>
671 eigenT &mim, std::vector<dataT> &weights, eigenCubeT &mask, dataT sigma,
double minGoodFract )
673 mim.resize( _rows, _cols );
675#pragma omp parallel num_threads( Eigen::nbThreads() )
677 std::vector<Scalar> work, wwork;
679#pragma omp for schedule( static, 10 )
680 for( Index i = 0; i < _rows; ++i )
682 for( Index j = 0; j < _cols; ++j )
689 for( Index k = 0; k < _planes; ++k )
691 if( ( mask.pixel( i, j ) )( k, 0 ) == 1 )
693 work.push_back( ( pixel( i, j ) )( k, 0 ) );
694 wwork.push_back( weights[k] );
697 if( work.size() > minGoodFract * _planes )
702 mim( i, j ) = invalidNumber<float>();
An image cube with an Eigen-like API.
void sigmaMean(eigenT &mim, std::vector< dataT > &weights, eigenCubeT &mask, Scalar sigma, double minGoodFract=0.0)
Calculate the sigma clipped weighted mean image of the cube, with a mask cube.
void Covar(Eigen::Matrix< dataT, Eigen::Dynamic, Eigen::Dynamic > &cv)
Calculate the covariance matrix of the images in the cube.
void sigmaMean(eigenT &mim, Scalar sigma)
Calculate the sigma clipped mean image of the cube.
void mean(eigenT &mim, std::vector< dataT > &weights, eigenCubeT &mask, double minGoodFract=0.0)
Calculate the weighted mean image of the cube, with a mask cube.
void clear()
De-allocate and set all sizes to 0.
const Eigen::Map< Eigen::Array< dataT, Eigen::Dynamic, Eigen::Dynamic > > image(Index n) const
Returns a 2D Eigen::Eigen::Map pointed at the specified image.
void sum(eigenT &mim)
Calculate the sum image of the cube.
void sigmaMean(eigenT &mim, std::vector< dataT > &weights, Scalar sigma)
Calculate the sigma clipped weighted mean image of the cube.
Eigen::Map< Eigen::Array< dataT, Eigen::Dynamic, Eigen::Dynamic >, Eigen::Unaligned, Eigen::Stride< Eigen::Dynamic, Eigen::Dynamic > > pixel(Index i, Index j)
Returns an Eigen::Eigen::Map-ed vector of the pixels at the given coordinate.
void median(eigenT &mim)
Calculate the median image of the cube.
void mean(eigenT &mim, std::vector< dataT > &weights)
Calculate the weighted mean image of the cube.
Eigen::Map< Eigen::Array< dataT, Eigen::Dynamic, Eigen::Dynamic > > cube()
Returns a 2D Eigen::Eigen::Map pointed at the entire cube.
void mean(eigenT &mim, eigenCubeT &mask, double minGoodFract=0.0)
Calculate the mean image of the cube with a mask.
void mean(eigenT &mim)
Calculate the mean image of the cube.
eigenCube(dataT *ndata, size_t nrows, size_t ncols, size_t nplanes)
C'tor taking an existing array as an argument.
void sigmaMean(eigenT &mim, eigenCubeT &mask, Scalar sigma, double minGoodFract=0.0)
Calculate the sigma clipped mean image of the cube, with a mask cube.
eigenCube(Index nrows, Index ncols, Index nplanes)
C'tor which will allocate space..
Eigen::Map< Eigen::Array< dataT, Eigen::Dynamic, Eigen::Dynamic > > asVectors()
Return an Eigen::Eigen::Map of the cube where each image is a vector.
Eigen::Map< Eigen::Array< dataT, Eigen::Dynamic, Eigen::Dynamic > > image(Index n)
Returns a 2D Eigen::Eigen::Map pointed at the specified image.
Tools for using the eigen library for image processing.
imageT::Scalar imageMedian(const imageT &mat, const maskT *mask, std::vector< typename imageT::Scalar > *work=0)
Calculate the median of an Eigen-like array.
vectorT::value_type vectorSigmaMean(const vectorT &vec, const vectorT *weights, const sigmaT &sigma, int &maxPasses)
Calculate the sigma-clipped mean of a vector.
valueT vectorMean(const valueT *vec, size_t sz)
Calculate the mean of a vector.
Header for the image processing utilities.