mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
Loading...
Searching...
No Matches
imCenterRadon.hpp
Go to the documentation of this file.
1/** \file imCenterRadon.hpp
2 * \brief A class to find PSF centers using circular symmetry.
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 __imCenterRadon_hpp__
28#define __imCenterRadon_hpp__
29
30#include "../math/vectorUtils.hpp"
31#include "../math/fit/fitGaussian.hpp"
32#include "imageFilters.hpp"
33#include "eigenImage.hpp"
34#include "imagePeakInterp.hpp"
35
36// #define ICCS_OMP
37namespace mx
38{
39namespace improc
40{
41
42/// Find the center of the image of a point source using circular symmetry of the PSF.
43/** Conducts a grid search through possible center positions. Finds the point with the minimum
44 * total variance around concentric circles about that point. The cost grid is fit with a Gaussian.
45 *
46 * \tparam realT is a real floating point type.
47 *
48 * \ingroup image_reg
49 */
50template <typename transformT>
52{
53
54 typedef typename transformT::arithT realT;
55
56 eigenImage<realT> m_grid; ///< The cost grid.
57
58 realT m_maxPix{ 2 }; ///< The maximum range of the grid search, the grid will be +/- maxPix. Default: 5.0.
59 realT m_dPix{ 0.1 }; ///< The spacing of the points in the grid search. Default: 0.1.
60 realT m_peakTol{ 0.1 }; ///< The tolerance for the peak finding. This is multiplicative with m_dPix.
61
62 realT m_minRad{ 1.0 }; ///< The minimum radius to include in the transform. Default: 1.0.
63 realT m_maxRad{ 10 }; ///< The maximum radius to include in the transform. Default: 10.0.
64 realT m_dRad{ 0.1 }; ///< The spacing of the radii. Default: 0.1.
65
67 1.0 }; ///< The angle tolerance in pixels at the maximum radius, sets the angular step size. Default is 1.0.
68
69 realT m_guessWidth{ 0.5 }; ///< Guess in original pixels for the FWHM of the cost grid. Default is 0.5, which seems
70 ///< to work well in most cases.
71
72 /// The fitter for finding the centroid of the grid. You can use this to inspect the details of the fit if desired.
73 // mx::math::fit::fitGaussian2D<mx::gaussian2D_gen_fitter<realT>> fit;
75
76 realT m_x0; ///< Internal storage of the guess center position
77 realT m_y0; ///< Internal storage of the guess center position.
78
80 {
81 }
82
83 /// Peform the grid search and fit.
84 int center( realT &x, ///< [out] the x-coordinate of the estimated center of rotational symmetry
85 realT &y, ///< [out] the y-coordinate estimated center of rotational symmetry
86 eigenImage<realT> &im, ///< [in] the image to centroid
87 realT x0, ///< [in] an initial guess at the x-coordinate, the grid is centered about this point.
88 realT y0 ///< [in] an initial guess at the y-coordinate, the grid is centered about this point.
89 )
90 {
91 m_x0 = x0;
92 m_y0 = y0;
93
94 const int lbuff = transformT::lbuff;
95 const int width = transformT::width;
96
97 transformT trans;
98
99 m_grid.resize( 2 * ( m_maxPix / m_dPix ) + 1, 2 * ( m_maxPix / m_dPix ) + 1 );
100
101 int N = 2 * ( m_maxPix / m_dPix ) + 1;
102
103 realT startX, startY;
104
106
107 for( int i = 0; i < N; ++i )
108 {
109 startX = x0 - m_maxPix + i * m_dPix;
110 for( int j = 0; j < N; ++j )
111 {
112 startY = y0 - m_maxPix + j * m_dPix;
113
114 m_grid( i, j ) = 0;
115
116 realT dang = math::rtod( atan( 1.0 / m_maxRad ) );
117 int Nang = 360. / dang + 1;
118 dang = 360. / Nang;
119
120#ifdef MX_IMPROC_IMCENTERRADON_OMP
121#pragma omp parallel for
122#endif
123 for( int iang = 0; iang < Nang; ++iang )
124 {
125 realT ang = iang * dang;
126
128 kern.resize( width, width );
129
130 realT cosang = cos( ang * 3.14159 / 180. );
131 realT sinang = sin( ang * 3.14159 / 180. );
132 realT lineint = 0;
133 for( realT rad = m_minRad; rad < m_maxRad; rad += m_dRad )
134 {
135 realT xx = startX + rad * cosang;
136 realT yy = startY + rad * sinang;
137
138 int i0 = xx;
139 int j0 = yy;
140 realT dx = xx - i0;
141 realT dy = yy - j0;
142
143 trans( kern, dx, dy );
144 lineint += ( im.block( i0 - lbuff, j0 - lbuff, width, width ) * kern ).sum();
145 }
146#pragma omp critical
147 m_grid( i, j ) += lineint * lineint * pow( m_dRad, 2 );
148 }
149 m_grid( i, j ) /= Nang;
150 } // j
151 } // i
152
153 return fitGrid( x, y, m_grid );
154 }
155
156 /// Fit the transform grid
157 /** A separate function for testing, does not normally need to be called independently
158 */
159 int fitGrid( realT &x, ///< [out] the x-coordinate of the estimated center of rotational symmetry
160 realT &y, ///< [out] the y-coordinate estimated center of rotational symmetry
161 eigenImage<realT> grid )
162 {
164
165 ipi( x, y, grid );
166
167 x = m_x0 - m_maxPix + x * m_dPix;
168 y = m_y0 - m_maxPix + y * m_dPix;
169
170 return 0;
171
172 // Get initial guess for the fit.
173 realT mx, mn, gx, gy;
174
175 mx = grid.maxCoeff( &gx, &gy );
176 mn = grid.minCoeff();
177
178 //------------- Now fit the grid to a 2D elliptical Gaussian ------------
179 /** \todo This needs to be able to handle highly non-symmetric peaks much better.
180 */
181 m_fit.setArray( grid.data(), grid.rows(), grid.cols() );
182 m_fit.setGuess( mn, mx - mn, gx, gy, m_guessWidth / m_dPix, m_guessWidth / m_dPix, 0 );
183
184 m_fit.fit();
185
186 // The estimated enter of circular symmetry from the fit.
187 x = m_x0 - m_maxPix + m_fit.x0() * m_dPix;
188 y = m_y0 - m_maxPix + m_fit.y0() * m_dPix;
189
190 return 0;
191 }
192
193 /// Output the results to a stream
194 /** Prints a result summary to the input stream.
195 *
196 * \tparam iosT is a std::ostream-like type.
197 * \tparam comment is a comment character to start each line. Can be '\0'.
198 */
199 template <typename iosT, char comment = '#'>
200 iosT &dumpResults( iosT &ios )
201 {
202 char c[] = { comment, '\0' };
203
204 ios << c << "--------------------------------------\n";
205 ios << c << "mx::improc::imCenterRadon Results \n";
206 ios << c << "--------------------------------------\n";
207 ios << c << "Estimated x-center: " << m_x0 - m_maxPix + m_fit.x0() * m_dPix << "\n";
208 ios << c << "Estimated y-center: " << m_y0 - m_maxPix + m_fit.y0() * m_dPix << "\n";
209 ios << c << "Cost half-width: "
210 << 0.5 * math::func::sigma2fwhm<realT>( sqrt( pow( m_fit.sigma_x(), 2 ) + pow( m_fit.sigma_y(), 2 ) ) *
211 m_dPix )
212 << " pixels\n";
213 ios << c << "--------------------------------------\n";
214 ios << c << "Setup:\n";
215 ios << c << "--------------------------------------\n";
216 ios << c << "maxPix: " << m_maxPix << "\n";
217 ios << c << "dPix: " << m_dPix << "\n";
218 ios << c << "minRad: " << m_minRad << "\n";
219 ios << c << "maxRad: " << m_maxRad << "\n";
220 ios << c << "dRad: " << m_dRad << "\n";
221 ios << c << "angTolPix: " << m_angTolPix << "\n";
222 ios << c << "guessWidth: " << m_guessWidth << "\n";
223 ios << c << "--------------------------------------\n";
224 ios << c << "Fit results:\n";
225 m_fit.dumpReport( ios );
226 ios << c << "--------------------------------------\n";
227
228 return ios;
229 }
230
231 /// Output the results to std::cout
232 /** Prints a result summary to std::cout
233 *
234 * \tparam comment is a comment character to start each line. Can be '\0'.
235 */
236 template <char comment = '#'>
237 std::ostream &dumpResults()
238 {
239 return dumpResults<std::ostream, comment>( std::cout );
240 }
241};
242
243} // namespace improc
244} // namespace mx
245
246#endif //__imCenterRadon_hpp__
Tools for using the eigen library for image processing.
Eigen::Array< scalarT, -1, -1 > eigenImage
Definition of the eigenImage type, which is an alias for Eigen::Array.
realT rtod(realT q)
Convert from radians to degrees.
Definition geo.hpp:153
Image filters (smoothing, radial profiles, etc.)
A class to find the location of a peak using interpolation.
The mxlib c++ namespace.
Definition mxError.hpp:106
Transformation by cubic convolution interpolation.
Find the center of the image of a point source using circular symmetry of the PSF.
int center(realT &x, realT &y, eigenImage< realT > &im, realT x0, realT y0)
Peform the grid search and fit.
realT m_x0
Internal storage of the guess center position.
realT m_peakTol
The tolerance for the peak finding. This is multiplicative with m_dPix.
int fitGrid(realT &x, realT &y, eigenImage< realT > grid)
Fit the transform grid.
iosT & dumpResults(iosT &ios)
Output the results to a stream.
realT m_dRad
The spacing of the radii. Default: 0.1.
realT m_maxRad
The maximum radius to include in the transform. Default: 10.0.
realT m_dPix
The spacing of the points in the grid search. Default: 0.1.
mx::math::fit::fitGaussian2Dgen< realT > m_fit
The fitter for finding the centroid of the grid. You can use this to inspect the details of the fit i...
realT m_minRad
The minimum radius to include in the transform. Default: 1.0.
realT m_angTolPix
The angle tolerance in pixels at the maximum radius, sets the angular step size. Default is 1....
realT m_maxPix
The maximum range of the grid search, the grid will be +/- maxPix. Default: 5.0.
eigenImage< realT > m_grid
The cost grid.
realT m_y0
Internal storage of the guess center position.
std::ostream & dumpResults()
Output the results to std::cout.
Find the peak of an image using interpolation.