mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
Loading...
Searching...
No Matches
eigenCube.hpp
Go to the documentation of this file.
1/** \file eigenCube.hpp
2 * \brief An image cube with an Eigen API
3 *
4 * \author Jared R. Males (jaredmales@gmail.com)
5 *
6 * \ingroup image_processing_files
7 *
8 */
9
10#ifndef eigenCube_hpp
11#define eigenCube_hpp
12
13#pragma GCC system_header
14#include <Eigen/Dense>
15
16#include "../math/vectorUtils.hpp"
17#include "eigenImage.hpp"
18#include "imageUtils.hpp"
19
20namespace mx
21{
22namespace improc
23{
24
25/// An image cube with an Eigen-like API
26/** \ingroup eigen_image_processing
27 */
28template <typename dataT>
30{
31 public:
32 typedef bool is_eigenCube;
33 typedef dataT Scalar;
34
35 typedef typename Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>::Index Index;
36
37 typedef Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>> imageRef;
38
39 typedef Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic> imageT;
40
41 protected:
42 Index _rows;
43 Index _cols;
44 Index _planes;
45
46 dataT *m_data{ nullptr };
47
48 bool _owner;
49
50 public:
51 eigenCube();
52
53 /// C'tor which will allocate space..
54 /**
55 */
56 eigenCube( Index nrows, ///< [in] Number of rows in the cube
57 Index ncols, ///< [in] Number of columns the cube
58 Index nplanes ///< [in] Number of planes in the cube
59 );
60
61 /// C'tor taking an existing array as an argument.
62 /** The existing array is used as is, as in a map, and ownership is not taken. You are
63 * responsible for memory management, e.g. free-ing this array.
64 */
65 eigenCube( dataT *ndata, ///< [in] Allocated array with a cube of nrows x ncols x nplanes
66 size_t nrows, ///< [in] Number of rows in the cube
67 size_t ncols, ///< [in] Number of columns the cube
68 size_t nplanes ///< [in] Number of planes in the cube
69 );
70
71 ~eigenCube();
72
73 void setZero();
74
75 eigenCube<dataT> &operator=( const eigenCube<dataT> &ec );
76
77 void shallowCopy( eigenCube<dataT> &src, bool takeOwner = false );
78
79 /// De-allocate and set all sizes to 0
80 void clear();
81
82 void resize( int r, int c, int p );
83
84 void resize( int r, int c );
85
86 dataT *data();
87
88 const dataT *data() const;
89
90 const Index rows() const;
91
92 const Index cols() const;
93
94 const Index planes() const;
95
96 /// Returns a 2D Eigen::Eigen::Map pointed at the entire cube.
97 Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>> cube();
98
99 /// Returns a 2D Eigen::Eigen::Map pointed at the specified image.
100 Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>> image( Index n /**< [in] the image number */ );
101
102 /// Returns a 2D Eigen::Eigen::Map pointed at the specified image.
103 const Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>>
104 image( Index n /**< [in] the image number */ ) const;
105
106 /// Returns an Eigen::Eigen::Map-ed vector of the pixels at the given coordinate
107 Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>,
108 Eigen::Unaligned,
109 Eigen::Stride<Eigen::Dynamic, Eigen::Dynamic>>
110 pixel( Index i, Index j );
111
112 /// Return an Eigen::Eigen::Map of the cube where each image is a vector
113 Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>> asVectors();
114
115 /// Calculate the covariance matrix of the images in the cube
116 void Covar( Eigen::Matrix<dataT, Eigen::Dynamic, Eigen::Dynamic> &cv );
117
118 /// Calculate the sum image of the cube
119 /**
120 * \tparam eigenT an Eigen-like type.
121 */
122 template <typename eigenT>
123 void sum( eigenT &mim /**< [out] the resultant sum image. Is resized. */ );
124
125 /// Calculate the mean image of the cube
126 /**
127 * \tparam eigenT an Eigen-like type.
128 */
129 template <typename eigenT>
130 void mean( eigenT &mim /**< [out] the resultant mean image. Is resized. */ );
131
132 /// Calculate the mean image of the cube with a mask.
133 /**
134 * \tparam eigenT an Eigen-like type.
135 * \tparam eigenCubeT an eigenCube type.
136 */
137 template <typename eigenT, typename eigenCubeT>
138 void mean( eigenT &mim, ///< [out] the resultant mean image. Is resized.
139 eigenCubeT &mask, /**< [in] a mask cube. Only pixels with value 1 are included in the
140 mean calculation. */
141 double minGoodFract = 0.0 /**< [in] [optional] the minimum fraction of good pixels, if not met
142 then the pixel is NaN-ed. */
143 );
144
145 /// Calculate the weighted mean image of the cube
146 /**
147 * \tparam eigenT an Eigen-like type.
148 */
149 template <typename eigenT>
150 void mean( eigenT &mim, ///< [out] the resultant mean image. Is resized.
151 std::vector<dataT> &weights ///< [in] a vector of weights to use for calculating the mean
152 );
153
154 /// Calculate the weighted mean image of the cube, with a mask cube.
155 /**
156 * \tparam eigenT an Eigen-like type.
157 * \tparam eigenCubeT an eigenCube type.
158 */
159 template <typename eigenT, typename eigenCubeT>
160 void mean( eigenT &mim, ///< [out] the resultant mean image. Is resized.
161 std::vector<dataT> &weights, ///< [in] a vector of weights to use for calculating the mean
162 eigenCubeT &mask, ///< [in] a mask cube. Only pixels with value 1 are included in the mean calculation.
163 double minGoodFract =
164 0.0 ///< [in] [optional] the minimum fraction of good pixels, if not met then the pixel is NaN-ed.
165 );
166
167 /// Calculate the median image of the cube
168 /**
169 * \tparam eigenT an Eigen-like type.
170 */
171 template <typename eigenT>
172 void median( eigenT &mim /**< [out] the resultant median image. Is resized. */ );
173
174 /// Calculate the sigma clipped mean image of the cube
175 /**
176 * \tparam eigenT an Eigen-like type.
177 */
178 template <typename eigenT>
179 void sigmaMean( eigenT &mim, ///< [out] the resultant mean image
180 Scalar sigma ///< [in] the sigma value at which to clip.
181 );
182
183 /// Calculate the sigma clipped mean image of the cube, with a mask cube
184 /**
185 * \tparam eigenT an Eigen-like type.
186 * \tparam eigenCubeT an eigenCube type.
187 */
188 template <typename eigenT, typename eigenCubeT>
189 void
190 sigmaMean( eigenT &mim, ///< [out] the resultant mean image. Is resized.
191 eigenCubeT &mask, ///< [in] a mask cube. Only pixels with value 1 are included in the mean calculation.
192 Scalar sigma, ///< [in] the sigma value at which to clip.
193 double minGoodFract =
194 0.0 ///< [in] [optional] the minimum fraction of good pixels, if not met then the pixel is NaN-ed.
195 );
196
197 /// Calculate the sigma clipped weighted mean image of the cube
198 /**
199 * \tparam eigenT an Eigen-like type.
200 */
201 template <typename eigenT>
202 void sigmaMean( eigenT &mim, ///< [out] the resultant mean image. Is resized.
203 std::vector<dataT> &weights, ///< [in] a vector of weights to use for calculating the mean
204 Scalar sigma ///< [in] the sigma value at which to clip.
205 );
206
207 /// Calculate the sigma clipped weighted mean image of the cube, with a mask cube.
208 /**
209 * \tparam eigenT an Eigen-like type.
210 * \tparam eigenCubeT an eigenCube type.
211 */
212 template <typename eigenT, typename eigenCubeT>
213 void
214 sigmaMean( eigenT &mim, ///< [out] the resultant mean image. Is resized.
215 std::vector<dataT> &weights, ///< [in] a vector of weights to use for calculating the mean
216 eigenCubeT &mask, ///< [in] a mask cube. Only pixels with value 1 are included in the mean calculation.
217 Scalar sigma, ///< [in] the sigma value at which to clip.
218 double minGoodFract =
219 0.0 ///< [in] [optional] the minimum fraction of good pixels, if not met then the pixel is NaN-ed.
220 );
221};
222
223template <typename dataT>
225{
226 _rows = 0;
227 _cols = 0;
228 _planes = 0;
229 _owner = false;
230}
231
232template <typename dataT>
233eigenCube<dataT>::eigenCube( Index nrows, Index ncols, Index nplanes )
234{
235 _rows = nrows;
236 _cols = ncols;
237 _planes = nplanes;
238
239 m_data = new dataT[_rows * _cols * _planes];
240 _owner = true;
241}
242
243template <typename dataT>
244eigenCube<dataT>::eigenCube( dataT *ndata, size_t nrows, size_t ncols, size_t nplanes )
245{
246 _rows = nrows;
247 _cols = ncols;
248 _planes = nplanes;
249
250 m_data = ndata;
251 _owner = false;
252}
253
254template <typename dataT>
256{
257 if( _owner && m_data )
258 {
259 delete[] m_data;
260 }
261}
262
263template <typename dataT>
264void eigenCube<dataT>::setZero()
265{
266 int N = _rows * _cols * _planes;
267
268 for( int i = 0; i < N; ++i )
269 m_data[i] = ( (dataT)0 );
270}
271
272template <typename dataT>
273eigenCube<dataT> &eigenCube<dataT>::operator=( const eigenCube<dataT> &ec )
274{
275 resize( ec.rows(), ec.cols(), ec.planes() );
276
277 int N = _rows * _cols * _planes;
278
279 for( int i = 0; i < N; ++i )
280 m_data[i] = ec.m_data[i];
281
282 return *this;
283}
284
285template <typename dataT>
286void eigenCube<dataT>::shallowCopy( eigenCube<dataT> &src, bool takeOwner )
287{
288 if( _owner && m_data )
289 {
290 delete m_data;
291 _owner = false;
292 }
293
294 _rows = src._rows;
295 _cols = src._cols;
296 _planes = src._planes;
297 m_data = src.m_data;
298
299 if( takeOwner == true )
300 {
301 _owner = true;
302 src._owner = false;
303 }
304 else
305 {
306 _owner = false;
307 }
308}
309
310template <typename dataT>
312{
313 if( _owner && m_data )
314 {
315 delete[] m_data;
316 }
317
318 _rows = 0;
319 _cols = 0;
320 _planes = 0;
321}
322
323template <typename dataT>
324void eigenCube<dataT>::resize( int r, int c, int p )
325{
326 if( _owner && m_data )
327 {
328 delete[] m_data;
329 }
330
331 _rows = r;
332 _cols = c;
333 _planes = p;
334
335 m_data = new dataT[_rows * _cols * _planes];
336 _owner = true;
337}
338
339template <typename dataT>
340void eigenCube<dataT>::resize( int r, int c )
341{
342 resize( r, c, 1 );
343}
344
345template <typename dataT>
346dataT *eigenCube<dataT>::data()
347{
348 return m_data;
349}
350
351template <typename dataT>
352const dataT *eigenCube<dataT>::data() const
353{
354 return m_data;
355}
356
357template <typename dataT>
358const typename eigenCube<dataT>::Index eigenCube<dataT>::rows() const
359{
360 return _rows;
361}
362
363template <typename dataT>
364const typename eigenCube<dataT>::Index eigenCube<dataT>::cols() const
365{
366 return _cols;
367}
368
369template <typename dataT>
370const typename eigenCube<dataT>::Index eigenCube<dataT>::planes() const
371{
372 return _planes;
373}
374
375template <typename dataT>
376Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>> eigenCube<dataT>::cube()
377{
378 return Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>>( m_data, _rows * _cols, _planes );
379}
380
381template <typename dataT>
382Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>> eigenCube<dataT>::image( Index n )
383{
384 return Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>>( m_data + n * _rows * _cols, _rows, _cols );
385}
386
387template <typename dataT>
388const Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>> eigenCube<dataT>::image( Index n ) const
389{
390 return Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>>( m_data + n * _rows * _cols, _rows, _cols );
391}
392
393template <typename dataT>
394Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>,
395 Eigen::Unaligned,
396 Eigen::Stride<Eigen::Dynamic, Eigen::Dynamic>>
397eigenCube<dataT>::pixel( Index i, Index j )
398{
399 return Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>,
400 Eigen::Unaligned,
401 Eigen::Stride<Eigen::Dynamic, Eigen::Dynamic>>(
402 m_data + j * _rows + i, _planes, 1, Eigen::Stride<Eigen::Dynamic, Eigen::Dynamic>( 0, _rows * _cols ) );
403}
404
405template <typename dataT>
406Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>> eigenCube<dataT>::asVectors()
407{
408 return Eigen::Map<Eigen::Array<dataT, Eigen::Dynamic, Eigen::Dynamic>>( m_data, _rows * _cols, _planes );
409}
410
411template <typename dataT>
412void eigenCube<dataT>::Covar( Eigen::Matrix<dataT, Eigen::Dynamic, Eigen::Dynamic> &cv )
413{
414 cv = asVectors().matrix().transpose() * asVectors().matrix();
415}
416
417template <typename dataT>
418template <typename eigenT>
419void eigenCube<dataT>::sum( eigenT &mim )
420{
421 mim.resize( _rows, _cols );
422
423#pragma omp parallel for schedule( static, 1 ) num_threads( Eigen::nbThreads() )
424 for( Index i = 0; i < _rows; ++i )
425 {
426 for( Index j = 0; j < _cols; ++j )
427 {
428 mim( i, j ) = pixel( i, j ).sum();
429 }
430 }
431}
432
433template <typename dataT>
434template <typename eigenT>
435void eigenCube<dataT>::mean( eigenT &mim )
436{
437 mim.resize( _rows, _cols );
438
439#pragma omp parallel for schedule( static, 1 ) num_threads( Eigen::nbThreads() )
440 for( Index i = 0; i < _rows; ++i )
441 {
442 for( Index j = 0; j < _cols; ++j )
443 {
444 mim( i, j ) = pixel( i, j ).mean();
445 }
446 }
447}
448
449template <typename dataT>
450template <typename eigenT, typename eigenCubeT>
451void eigenCube<dataT>::mean( eigenT &mim, eigenCubeT &mask, double minGoodFract )
452{
453 mim.resize( _rows, _cols );
454
455 #pragma omp parallel
456 {
457 std::vector<dataT> work;
458
459 #pragma omp for
460 for( Index i = 0; i < _rows; ++i )
461 {
462 for( Index j = 0; j < _cols; ++j )
463 {
464
465 work.clear();
466
467 for( Index k = 0; k < _planes; ++k )
468 {
469 if( ( mask.pixel( i, j ) )( k, 0 ) == 1 )
470 {
471 work.push_back( ( pixel( i, j ) )( k, 0 ) );
472 }
473 }
474
475 if( work.size() > minGoodFract * _planes )
476 {
477 mim( i, j ) = math::vectorMean( work );
478 }
479 else
480 {
481 mim( i, j ) = invalidNumber<float>();
482 }
483 }
484 }
485 }
486}
487
488template <typename dataT>
489template <typename eigenT>
490void eigenCube<dataT>::mean( eigenT &mim, std::vector<dataT> &weights )
491{
492 mim.resize( _rows, _cols );
493
494#pragma omp parallel num_threads( Eigen::nbThreads() )
495 {
496 std::vector<Scalar> work;
497
498#pragma omp for schedule( static, 10 )
499 for( Index i = 0; i < _rows; ++i )
500 {
501 for( Index j = 0; j < _cols; ++j )
502 {
503 work.resize( _planes ); // work could be smaller after sigmaMean
504 // int ii=0;
505 for( int k = 0; k < _planes; ++k )
506 work[k] = ( pixel( i, j ) )( k, 0 );
507
508 mim( i, j ) = math::vectorMean( work, weights );
509 }
510 }
511 }
512}
513
514template <typename dataT>
515template <typename eigenT, typename eigenCubeT>
516void eigenCube<dataT>::mean( eigenT &mim, std::vector<dataT> &weights, eigenCubeT &mask, double minGoodFract )
517{
518 mim.resize( _rows, _cols );
519
520#pragma omp parallel num_threads( Eigen::nbThreads() )
521 {
522 std::vector<Scalar> work, wwork;
523
524#pragma omp for schedule( static, 10 )
525 for( Index i = 0; i < _rows; ++i )
526 {
527 for( Index j = 0; j < _cols; ++j )
528 {
529 work.clear();
530 wwork.clear();
531
532 // int ii=0;
533
534 for( Index k = 0; k < _planes; ++k )
535 {
536 if( ( mask.pixel( i, j ) )( k, 0 ) == 1 )
537 {
538 work.push_back( ( pixel( i, j ) )( k, 0 ) );
539 wwork.push_back( weights[k] );
540 }
541 }
542 if( work.size() > minGoodFract * _planes )
543 {
544 mim( i, j ) = math::vectorMean( work, wwork );
545 }
546 else
547 mim( i, j ) = invalidNumber<float>();
548 }
549 }
550 }
551}
552
553template <typename dataT>
554template <typename eigenT>
555void eigenCube<dataT>::median( eigenT &mim )
556{
557 mim.resize( _rows, _cols );
558
559#pragma omp parallel for schedule( static, 10 ) num_threads( Eigen::nbThreads() )
560 for( Index i = 0; i < _rows; ++i )
561 {
562 // Eigen::Eigen::Array<Scalar, Eigen::Eigen::Dynamic, Eigen::Eigen::Dynamic> work;
563 std::vector<Scalar> work;
564 for( Index j = 0; j < _cols; ++j )
565 {
566 mim( i, j ) = imageMedian( pixel( i, j ), &work );
567 }
568 }
569}
570
571template <typename dataT>
572template <typename eigenT>
573void eigenCube<dataT>::sigmaMean( eigenT &mim, dataT sigma )
574{
575 mim.resize( _rows, _cols );
576
577 #pragma omp parallel num_threads( Eigen::nbThreads() )
578 {
579 std::vector<Scalar> work;
580
581 #pragma omp for schedule( static, 10 )
582 for( Index i = 0; i < _rows; ++i )
583 {
584 for( Index j = 0; j < _cols; ++j )
585 {
586 work.resize( _planes ); // work could be smaller after sigmaMean
587 for( int k = 0; k < _planes; ++k )
588 work[k] = ( pixel( i, j ) )( k, 0 );
589
590 mim( i, j ) = math::vectorSigmaMean( work, sigma );
591 }
592 }
593 }
594}
595
596template <typename dataT>
597template <typename eigenT, typename eigenCubeT>
598void eigenCube<dataT>::sigmaMean( eigenT &mim, eigenCubeT &mask, dataT sigma, double minGoodFract )
599{
600 mim.resize( _rows, _cols );
601
602 #pragma omp parallel
603 {
604 std::vector<Scalar> work;
605
606 #pragma omp for
607 for( Index i = 0; i < _rows; ++i )
608 {
609 for( Index j = 0; j < _cols; ++j )
610 {
611 work.clear();
612 int ii = 0;
613
614 for( Index k = 0; k < _planes; ++k )
615 {
616 if( ( mask.pixel( i, j ) )( k, 0 ) == 1 )
617 {
618 work.push_back( ( pixel( i, j ) )( k, 0 ) );
619 }
620 }
621
622 if( work.size() > minGoodFract * _planes )
623 {
624 mim( i, j ) = math::vectorSigmaMean( work, sigma );
625 }
626 else
627 {
628 mim( i, j ) = invalidNumber<float>();
629 }
630 }
631 }
632 }
633}
634
635template <typename dataT>
636template <typename eigenT>
637void eigenCube<dataT>::sigmaMean( eigenT &mim, std::vector<dataT> &weights, dataT sigma )
638{
639 mim.resize( _rows, _cols );
640
641#pragma omp parallel num_threads( Eigen::nbThreads() )
642 {
643 std::vector<Scalar> work;
644
645#pragma omp for schedule( static, 10 )
646 for( Index i = 0; i < _rows; ++i )
647 {
648 for( Index j = 0; j < _cols; ++j )
649 {
650 work.resize( _planes ); // work could be smaller after sigmaMean
651 int ii = 0;
652 for( int k = 0; k < _planes; ++k )
653 work[k] = ( pixel( i, j ) )( k, 0 );
654
655 mim( i, j ) = math::vectorSigmaMean( work, weights, sigma );
656 }
657 }
658 }
659}
660
661template <typename dataT>
662template <typename eigenT, typename eigenCubeT>
664 eigenT &mim, std::vector<dataT> &weights, eigenCubeT &mask, dataT sigma, double minGoodFract )
665{
666 mim.resize( _rows, _cols );
667
668#pragma omp parallel num_threads( Eigen::nbThreads() )
669 {
670 std::vector<Scalar> work, wwork;
671
672#pragma omp for schedule( static, 10 )
673 for( Index i = 0; i < _rows; ++i )
674 {
675 for( Index j = 0; j < _cols; ++j )
676 {
677 work.clear();
678 wwork.clear();
679
680 int ii = 0;
681
682 for( Index k = 0; k < _planes; ++k )
683 {
684 if( ( mask.pixel( i, j ) )( k, 0 ) == 1 )
685 {
686 work.push_back( ( pixel( i, j ) )( k, 0 ) );
687 wwork.push_back( weights[k] );
688 }
689 }
690 if( work.size() > minGoodFract * _planes )
691 {
692 mim( i, j ) = math::vectorSigmaMean( work, wwork, sigma );
693 }
694 else
695 mim( i, j ) = invalidNumber<float>();
696 }
697 }
698 }
699}
700
701} // namespace improc
702
703} // namespace mx
704
705#endif // eigenCube_hpp
An image cube with an Eigen-like API.
Definition eigenCube.hpp:30
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.
The mxlib c++ namespace.
Definition mxError.hpp:106