mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
Loading...
Searching...
No Matches
eigenImage.hpp
Go to the documentation of this file.
1/** \file eigenImage.hpp
2 * \brief Tools for using the eigen library for image processing
3 * \ingroup image_processing_files
4 * \author Jared R. Males (jaredmales@gmail.com)
5 *
6 */
7
8//***********************************************************************//
9// Copyright 2015, 2016, 2017 Jared R. Males (jaredmales@gmail.com)
10//
11// This file is part of mxlib.
12//
13// mxlib is free software: you can redistribute it and/or modify
14// it under the terms of the GNU General Public License as published by
15// the Free Software Foundation, either version 3 of the License, or
16// (at your option) any later version.
17//
18// mxlib is distributed in the hope that it will be useful,
19// but WITHOUT ANY WARRANTY; without even the implied warranty of
20// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21// GNU General Public License for more details.
22//
23// You should have received a copy of the GNU General Public License
24// along with mxlib. If not, see <http://www.gnu.org/licenses/>.
25//***********************************************************************//
26
27#ifndef improc_eigenImage_hpp
28#define improc_eigenImage_hpp
29
30#pragma GCC system_header
31#include <Eigen/Dense>
32
33#include "../math/vectorUtils.hpp"
34
35namespace mx
36{
37namespace improc
38{
39
40/// Definition of the eigenImage type, which is an alias for Eigen::Array.
41/** \ingroup eigen_image_processing
42 */
43template <typename scalarT>
44using eigenImage = Eigen::Array<scalarT, -1, -1>;
45
46/// Definition of the eigenMap type, which is an alias for Eigen::Map<Array>.
47/** \ingroup eigen_image_processing
48 */
49template <typename scalarT>
50using eigenMap = Eigen::Map<Eigen::Array<scalarT, -1, -1>>;
51
52/// Test whether a type is an eigenCube by testing whether it has a typedef of "is_eigenCube"
53/** Used for compile-time determination of type
54 * Example usage:
55 * \code
56 * bool is_eC = is_eigenCube<eigenCube<float> >; //Evaluates to true
57 * bool is_not_eC = is_eigenCube<eigenImagef>; //Evaluates to false
58 * \endcode
59 *
60 * This was taken directly from the example at http://en.wikipedia.org/wiki/Substitution_failure_is_not_an_error
61 *
62 * \ingroup eigen_image_processing
63 */
64template <typename T>
66{
67 // Types "yes" and "no" are guaranteed to have different sizes,
68 // specifically sizeof(yes) == 1 and sizeof(no) == 2.
69 typedef char yes[1];
70 typedef char no[2];
71
72 template <typename imageT>
73 static yes &test( typename imageT::is_eigenCube * );
74
75 template <typename>
76 static no &test( ... );
77
78 // If the "sizeof" of the result of calling test<T>(0) would be equal to sizeof(yes),
79 // the first overload worked and T has a nested type named "is_eigenCube".
80 static const bool value = sizeof( test<T>( 0 ) ) == sizeof( yes );
81};
82
83/// Function object to return the number of planes for any Eigen like object, whether 2D or a 3D cube.
84/** Uses SFINAE to check for 3D eigenCube.
85 *
86 * \ingroup eigen_image_processing
87 */
88template <typename arrT, bool isCube = is_eigenCube<arrT>::value>
90{
91 // If it's an eigenCube, call planes planes()
92 int operator()( const arrT &arr )
93 {
94 return arr.planes();
95 }
96};
97
98template <typename arrT>
99struct eigenArrPlanes<arrT, false>
100{
101 // If it's not an eigenCube, never call planes()
102 int operator()( const arrT &arr )
103 {
104 return 1;
105 }
106};
107
108/// Calculate the median of an Eigen-like array.
109/** Calculates the median of the entire array, allowing for some pixels to be ignored using a mask.
110 * Working memory can be retained between calls.
111 *
112 * \tparam imageT is an Eigen-like type
113 * \tparam maskT is an Eigen-like type
114 *
115 * \returns the median of the unmasked pixels of mat, using \ref vectorMedianInPlace().
116 *
117 * \ingroup eigen_image_processing
118 *
119 */
120template <typename imageT, typename maskT = imageT>
121typename imageT::Scalar
122imageMedian( const imageT &mat, ///< [in] the image to take the median of
123 const maskT *mask, ///< [in] if non-0, a 1/0 mask where 0 pixels are ignored.
124 std::vector<typename imageT::Scalar> *work =
125 0 ///< [in] [optional] working memory can be retained and re-passed. Is resized.
126)
127{
128 typename imageT::Scalar med;
129
130 bool localWork = false;
131 if( work == 0 )
132 {
133 work = new std::vector<typename imageT::Scalar>;
134 localWork = true;
135 }
136
137 int sz = mat.size();
138
139 if( mask )
140 {
141 sz = mask->sum();
142 }
143
144 work->resize( sz );
145
146 int ii = 0;
147 for( int i = 0; i < mat.rows(); ++i )
148 {
149 for( int j = 0; j < mat.cols(); ++j )
150 {
151 if( mask )
152 {
153 if( ( *mask )( i, j ) == 0 )
154 continue;
155 }
156
157 ( *work )[ii] = mat( i, j );
158 ++ii;
159 }
160 }
161
162 med = math::vectorMedianInPlace( *work );
163
164 if( localWork )
165 delete work;
166
167 return med;
168}
169
170/// Calculate the median of an Eigen-like array.
171/** Calculates the median of the entire array.
172 * Working memory can be retained between calls.
173 *
174 * \tparam imageT is an Eigen-like type
175 *
176 * \returns the median of the unmasked pixels of mat, using \ref vectorMedianInPlace().
177 *
178 * \ingroup eigen_image_processing
179 *
180 */
181template <typename imageT>
182typename imageT::Scalar imageMedian(
183 const imageT &mat, ///< [in] the image to take the median of
184 std::vector<typename imageT::Scalar> *work = 0 ///< [in] [optional] working memory can be retained and re-passed.
185)
186{
187 return imageMedian( mat, (Eigen::Array<typename imageT::Scalar, -1, -1> *)0, work );
188}
189
190} // namespace improc
191} // namespace mx
192
193#endif // improc_eigenImage_hpp
imageT::Scalar imageMedian(const imageT &mat, const maskT *mask, std::vector< typename imageT::Scalar > *work=0)
Calculate the median of an Eigen-like array.
Eigen::Array< scalarT, -1, -1 > eigenImage
Definition of the eigenImage type, which is an alias for Eigen::Array.
Eigen::Map< Eigen::Array< scalarT, -1, -1 > > eigenMap
Definition of the eigenMap type, which is an alias for Eigen::Map<Array>.
vectorT::value_type vectorMedianInPlace(vectorT &vec)
Calculate median of a vector in-place, altering the vector.
The mxlib c++ namespace.
Definition mxError.hpp:106
Function object to return the number of planes for any Eigen like object, whether 2D or a 3D cube.
Test whether a type is an eigenCube by testing whether it has a typedef of "is_eigenCube".