|
mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
|
An image cube with an Eigen-like API.
Definition at line 29 of file eigenCube.hpp.
#include <improc/eigenCube.hpp>
Public Member Functions | |
| eigenCube (Index nrows, Index ncols, Index nplanes) | |
| C'tor which will allocate space.. | |
| eigenCube (dataT *ndata, size_t nrows, size_t ncols, size_t nplanes) | |
| C'tor taking an existing array as an argument. | |
| void | clear () |
| De-allocate and set all sizes to 0. | |
| Eigen::Map< Eigen::Array< dataT, Eigen::Dynamic, Eigen::Dynamic > > | cube () |
| Returns a 2D Eigen::Eigen::Map pointed at the entire cube. | |
| Eigen::Map< Eigen::Array< dataT, Eigen::Dynamic, Eigen::Dynamic > > | image (Index n) |
| Returns a 2D Eigen::Eigen::Map pointed at the specified image. | |
| 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. | |
| 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. | |
| Eigen::Map< Eigen::Array< dataT, Eigen::Dynamic, Eigen::Dynamic > > | asVectors () |
| Return an Eigen::Eigen::Map of the cube where each image is a vector. | |
| void | Covar (Eigen::Matrix< dataT, Eigen::Dynamic, Eigen::Dynamic > &cv) |
| Calculate the covariance matrix of the images in the cube. | |
| template<typename eigenT > | |
| void | sum (eigenT &mim) |
| Calculate the sum image of the cube. | |
| template<typename eigenT > | |
| void | mean (eigenT &mim) |
| Calculate the mean image of the cube. | |
| template<typename eigenT , typename eigenCubeT > | |
| void | mean (eigenT &mim, eigenCubeT &mask, double minGoodFract=0.0) |
| Calculate the mean image of the cube with a mask. | |
| template<typename eigenT > | |
| void | mean (eigenT &mim, std::vector< dataT > &weights) |
| Calculate the weighted mean image of the cube. | |
| template<typename eigenT , typename eigenCubeT > | |
| 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. | |
| template<typename eigenT > | |
| void | median (eigenT &mim) |
| Calculate the median image of the cube. | |
| template<typename eigenT > | |
| void | sigmaMean (eigenT &mim, Scalar sigma) |
| Calculate the sigma clipped mean image of the cube. | |
| template<typename eigenT , typename eigenCubeT > | |
| 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. | |
| template<typename eigenT > | |
| void | sigmaMean (eigenT &mim, std::vector< dataT > &weights, Scalar sigma) |
| Calculate the sigma clipped weighted mean image of the cube. | |
| template<typename eigenT , typename eigenCubeT > | |
| 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. | |
| mx::improc::eigenCube< dataT >::eigenCube | ( | Index | nrows, |
| Index | ncols, | ||
| Index | nplanes | ||
| ) |
C'tor which will allocate space..
| [in] | nrows | Number of rows in the cube |
| [in] | ncols | Number of columns the cube |
| [in] | nplanes | Number of planes in the cube |
Definition at line 233 of file eigenCube.hpp.
| mx::improc::eigenCube< dataT >::eigenCube | ( | dataT * | ndata, |
| size_t | nrows, | ||
| size_t | ncols, | ||
| size_t | nplanes | ||
| ) |
C'tor taking an existing array as an argument.
The existing array is used as is, as in a map, and ownership is not taken. You are responsible for memory management, e.g. free-ing this array.
| [in] | ndata | Allocated array with a cube of nrows x ncols x nplanes |
| [in] | nrows | Number of rows in the cube |
| [in] | ncols | Number of columns the cube |
| [in] | nplanes | Number of planes in the cube |
Definition at line 244 of file eigenCube.hpp.
| Eigen::Map< Eigen::Array< dataT, Eigen::Dynamic, Eigen::Dynamic > > mx::improc::eigenCube< dataT >::asVectors | ( | ) |
Return an Eigen::Eigen::Map of the cube where each image is a vector.
Definition at line 406 of file eigenCube.hpp.
Referenced by mx::AO::ifPInv(), and mx::AO::m2cMatrix().
| void mx::improc::eigenCube< dataT >::clear | ( | ) |
De-allocate and set all sizes to 0.
Definition at line 311 of file eigenCube.hpp.
| void mx::improc::eigenCube< dataT >::Covar | ( | Eigen::Matrix< dataT, Eigen::Dynamic, Eigen::Dynamic > & | cv | ) |
Calculate the covariance matrix of the images in the cube.
Definition at line 412 of file eigenCube.hpp.
| Eigen::Map< Eigen::Array< dataT, Eigen::Dynamic, Eigen::Dynamic > > mx::improc::eigenCube< dataT >::cube | ( | ) |
Returns a 2D Eigen::Eigen::Map pointed at the entire cube.
Definition at line 376 of file eigenCube.hpp.
| Eigen::Map< Eigen::Array< dataT, Eigen::Dynamic, Eigen::Dynamic > > mx::improc::eigenCube< dataT >::image | ( | Index | n | ) |
Returns a 2D Eigen::Eigen::Map pointed at the specified image.
| [in] | n | the image number |
Definition at line 382 of file eigenCube.hpp.
Referenced by mx::AO::influenceFunctionsGaussian(), mx::AO::analysis::fourierTemporalPSD< _realT, aosysT >::intensityPSD(), mx::wfp::lyotCoronagraph< _realT, _fpmaskFloatT >::loadFocalMask(), mx::AO::makeModfBasis(), mx::AO::makeZernikeBasis(), mx::AO::modesOnDM(), mx::AO::orthogonalizeBasis(), and SCENARIO().
| const Eigen::Map< Eigen::Array< dataT, Eigen::Dynamic, Eigen::Dynamic > > mx::improc::eigenCube< dataT >::image | ( | Index | n | ) | const |
Returns a 2D Eigen::Eigen::Map pointed at the specified image.
| [in] | n | the image number |
Definition at line 388 of file eigenCube.hpp.
| void mx::improc::eigenCube< dataT >::mean | ( | eigenT & | mim | ) |
Calculate the mean image of the cube.
| eigenT | an Eigen-like type. |
| [out] | mim | the resultant mean image. Is resized. |
Definition at line 435 of file eigenCube.hpp.
| void mx::improc::eigenCube< dataT >::mean | ( | eigenT & | mim, |
| eigenCubeT & | mask, | ||
| double | minGoodFract = 0.0 |
||
| ) |
Calculate the mean image of the cube with a mask.
| eigenT | an Eigen-like type. |
| eigenCubeT | an eigenCube type. |
| [out] | mim | the resultant mean image. Is resized. |
| [in] | mask | a mask cube. Only pixels with value 1 are included in the mean calculation. |
| [in] | minGoodFract | [optional] the minimum fraction of good pixels, if not met then the pixel is NaN-ed. |
Definition at line 451 of file eigenCube.hpp.
References mx::math::vectorMean().
| void mx::improc::eigenCube< dataT >::mean | ( | eigenT & | mim, |
| std::vector< dataT > & | weights | ||
| ) |
Calculate the weighted mean image of the cube.
| eigenT | an Eigen-like type. |
| [out] | mim | the resultant mean image. Is resized. |
| [in] | weights | a vector of weights to use for calculating the mean |
Definition at line 490 of file eigenCube.hpp.
References mx::math::vectorMean().
| void mx::improc::eigenCube< dataT >::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.
| eigenT | an Eigen-like type. |
| eigenCubeT | an eigenCube type. |
| [out] | mim | the resultant mean image. Is resized. |
| [in] | weights | a vector of weights to use for calculating the mean |
| [in] | mask | a mask cube. Only pixels with value 1 are included in the mean calculation. |
| [in] | minGoodFract | [optional] the minimum fraction of good pixels, if not met then the pixel is NaN-ed. |
Definition at line 516 of file eigenCube.hpp.
References mx::math::vectorMean().
| void mx::improc::eigenCube< dataT >::median | ( | eigenT & | mim | ) |
Calculate the median image of the cube.
| eigenT | an Eigen-like type. |
| [out] | mim | the resultant median image. Is resized. |
Definition at line 555 of file eigenCube.hpp.
References mx::improc::imageMedian().
| Eigen::Map< Eigen::Array< dataT, Eigen::Dynamic, Eigen::Dynamic >, Eigen::Unaligned, Eigen::Stride< Eigen::Dynamic, Eigen::Dynamic > > mx::improc::eigenCube< dataT >::pixel | ( | Index | i, |
| Index | j | ||
| ) |
Returns an Eigen::Eigen::Map-ed vector of the pixels at the given coordinate.
Definition at line 397 of file eigenCube.hpp.
Referenced by SCENARIO().
| void mx::improc::eigenCube< dataT >::sigmaMean | ( | eigenT & | mim, |
| eigenCubeT & | mask, | ||
| Scalar | sigma, | ||
| double | minGoodFract = 0.0 |
||
| ) |
Calculate the sigma clipped mean image of the cube, with a mask cube.
| eigenT | an Eigen-like type. |
| eigenCubeT | an eigenCube type. |
| [out] | mim | the resultant mean image. Is resized. |
| [in] | mask | a mask cube. Only pixels with value 1 are included in the mean calculation. |
| [in] | sigma | the sigma value at which to clip. |
| [in] | minGoodFract | [optional] the minimum fraction of good pixels, if not met then the pixel is NaN-ed. |
Definition at line 598 of file eigenCube.hpp.
References mx::math::vectorSigmaMean().
| void mx::improc::eigenCube< dataT >::sigmaMean | ( | eigenT & | mim, |
| Scalar | sigma | ||
| ) |
Calculate the sigma clipped mean image of the cube.
| eigenT | an Eigen-like type. |
| [out] | mim | the resultant mean image |
| [in] | sigma | the sigma value at which to clip. |
Definition at line 573 of file eigenCube.hpp.
References mx::math::vectorSigmaMean().
| void mx::improc::eigenCube< dataT >::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.
| eigenT | an Eigen-like type. |
| eigenCubeT | an eigenCube type. |
| [out] | mim | the resultant mean image. Is resized. |
| [in] | weights | a vector of weights to use for calculating the mean |
| [in] | mask | a mask cube. Only pixels with value 1 are included in the mean calculation. |
| [in] | sigma | the sigma value at which to clip. |
| [in] | minGoodFract | [optional] the minimum fraction of good pixels, if not met then the pixel is NaN-ed. |
Definition at line 663 of file eigenCube.hpp.
References mx::math::vectorSigmaMean().
| void mx::improc::eigenCube< dataT >::sigmaMean | ( | eigenT & | mim, |
| std::vector< dataT > & | weights, | ||
| Scalar | sigma | ||
| ) |
Calculate the sigma clipped weighted mean image of the cube.
| eigenT | an Eigen-like type. |
| [out] | mim | the resultant mean image. Is resized. |
| [in] | weights | a vector of weights to use for calculating the mean |
| [in] | sigma | the sigma value at which to clip. |
Definition at line 637 of file eigenCube.hpp.
References mx::math::vectorSigmaMean().
| void mx::improc::eigenCube< dataT >::sum | ( | eigenT & | mim | ) |
Calculate the sum image of the cube.
| eigenT | an Eigen-like type. |
| [out] | mim | the resultant sum image. Is resized. |
Definition at line 419 of file eigenCube.hpp.