mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
Loading...
Searching...
No Matches
fitAiry.hpp
Go to the documentation of this file.
1/** \file fitAiry.hpp
2 * \author Jared R. Males
3 * \brief Tools for fitting the Airy pattern to PSF images.
4 * \ingroup fitting_files
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 fitAiry_hpp
28#define fitAiry_hpp
29
30#include "levmarInterface.hpp"
31#include "../func/airyPattern.hpp"
32
33namespace mx
34{
35namespace math
36{
37namespace fit
38{
39
40template <typename realT>
41struct array2FitAiry;
42
43template <typename _realT>
44struct airy2D_obs_fitter;
45
46template <typename _realT>
47struct airy2D_obs_fitter_ps;
48
49/** \defgroup airy_peak_fit Airy Patterns
50 * \brief Fitting Airy patterns to data.
51 *
52 * The Airy pattern is fit to data, typically diffraction-limited images of point sources.
53 *
54 * \ingroup peak_fit
55 */
56
57/// Class to manage fitting a 2D Airy pattern to data via the \ref levmarInterface
58/** In addition to the requirements on fitterT specified by \ref levmarInterface
59 * this class also requires this definition in fitterT
60 * \code
61 * static const int nparams = 5;
62 * \endcode
63 * where the number 5 is replaced by the number of parameters that fitterT expects to fit.
64 *
65 * \tparam fitterT a type meeting the above requirements.
66 *
67 * \ingroup airy_peak_fit
68 *
69 */
70template <typename fitterT>
71class fitAiry2D : public levmarInterface<fitterT>
72{
73
74 public:
75 typedef typename fitterT::realT realT;
76
77 static const int nparams = fitterT::nparams;
78
80
81 void initialize()
82 {
83 this->allocate_params( nparams );
84 this->adata = &arr;
85 }
86
87 fitAiry2D()
88 {
89 initialize();
90 }
91
93 {
94 }
95
96 /// Set the initial guess when platescale and central obscuration are fixed.
97 void setGuess( realT A0, ///< [in] the constant background
98 realT A, ///< [in] the peak scaling
99 realT x0, ///< [in] the center x-coordinate [pixels]
100 realT y0 ///< [in] the center y-coordinate [pixels]
101 )
102 {
103 static_assert( nparams == 4, "fitAiry2D: Wrong setGuess called for platescale and/or cen-obs variable." );
104
105 this->p[0] = A0;
106 this->p[1] = A;
107 this->p[2] = x0;
108 this->p[3] = y0;
109 }
110
111 /// Set the initial guess when platescale is variable, and central obscuration is fixed.
112 void setGuess( realT A0, ///< [in] the constant background
113 realT A, ///< [in] the peak scaling
114 realT x0, ///< [in] the center x-coordinate [pixels]
115 realT y0, ///< [in] the center y-coordinate [pixels]
116 realT ps ///< [in] the platescale [ (lambda/D) / pixel]
117 )
118 {
119 static_assert( nparams == 5, "fitAiry2D: Wrong setGuess called for only platescale variable." );
120
121 this->p[0] = A0;
122 this->p[1] = A;
123 this->p[2] = x0;
124 this->p[3] = y0;
125 this->p[4] = ps;
126 }
127
128 /// Set the initial guess when central-obscuration is variable.
129 void setGuess( realT A0, ///< [in] the constant background
130 realT A, ///< [in] the peak scaling
131 realT x0, ///< [in] the center x-coordinate [pixels]
132 realT y0, ///< [in] the center y-coordinate [pixels]
133 realT ps, ///< [in] the platescale [ (lambda/D) / pixel]
134 realT co ///< [in] the central obscuration [ratio]
135 )
136 {
137 static_assert( nparams == 6, "fitAiry2D: Wrong setGuess called for cen-obs not variable." );
138
139 this->p[0] = A0;
140 this->p[1] = A;
141 this->p[2] = x0;
142 this->p[3] = y0;
143 this->p[4] = ps;
144 this->p[5] = co;
145 }
146
147 void setArray( realT *data, int nx, int ny )
148 {
149 arr.data = data;
150 arr.nx = nx;
151 arr.ny = ny;
152
153 this->n = nx * ny;
154 }
155
156 void cenObs( realT co )
157 {
158 arr.cenObs = co;
159 }
160
161 void ps( realT ps )
162 {
163 arr.ps = ps;
164 }
165
166 int fit()
167 {
169 }
170
171 realT A0()
172 {
173 return this->p[0];
174 }
175
176 realT A()
177 {
178 return this->p[1];
179 }
180
181 realT x0()
182 {
183 return this->p[2];
184 }
185
186 realT y0()
187 {
188 return this->p[3];
189 }
190
191 realT ps()
192 {
193 if( nparams < 5 )
194 {
195 return arr.ps;
196 }
197 else
198 {
199 return this->p[4];
200 }
201 }
202
203 realT cenObs()
204 {
205 if( nparams < 6 )
206 {
207 return arr.cenObs;
208 }
209 else
210 {
211 return this->p[5];
212 }
213 }
214};
215
216/// Wrapper for a native array to pass to \ref levmarInterface, with Airy details.
217/** \ingroup airy_peak_fit
218 */
219template <typename realT>
221{
222 realT *data{ nullptr }; ///< Pointer to the array
223 size_t nx{ 0 }; ///< X dimension of the array
224 size_t ny{ 0 }; ///< Y dimension of the array
225
226 realT cenObs{ 0 }; ///< is the ratio of the circular central obscuration diameter to the diameter.
227 realT ps{ 0 }; ///< the platescale in \f$ (\lambda/D)/pixel \f$
228};
229
230///\ref levmarInterface fitter structure for the centrally obscured Airy pattern.
231/**
232 * Platescale and central obscuration are fixed.
233 *
234 * \ingroup airy_peak_fit
235 *
236 */
237template <typename _realT>
239{
240 typedef _realT realT;
241
242 static const int nparams = 4;
243
244 static void func( realT *p, realT *hx, int m, int n, void *adata )
245 {
247
248 size_t idx_mat, idx_dat;
249
250 idx_dat = 0;
251
252 /*
253 p[0] = floor
254 p[1] = scale
255 p[2] = x0
256 p[3] = y0
257 */
258
259 realT r;
260
261 for( int i = 0; i < arr->nx; i++ )
262 {
263 for( int j = 0; j < arr->ny; j++ )
264 {
265 idx_mat = i + j * arr->nx;
266
267 r = sqrt( pow( i - p[2], 2 ) + pow( j - p[3], 2 ) );
268
269 hx[idx_dat] = func::airyPattern( static_cast<realT>( i ),
270 static_cast<realT>( j ),
271 p[0],
272 p[1],
273 p[2],
274 p[3],
275 arr->ps,
276 arr->cenObs ) -
277 arr->data[idx_mat];
278
279 // hx[idx_dat] *= fabs(arr->data[idx_mat]);
280
281 idx_dat++;
282 }
283 }
284 }
285};
286
287///\ref levmarInterface fitter structure for the obstructed Airy pattern, including platescale.
288/** Central obscuration is fixed.
289 *
290 * \ingroup airy_peak_fit
291 *
292 */
293template <typename _realT>
295{
296 typedef _realT realT;
297
298 static const int nparams = 5;
299
300 static void func( realT *p, realT *hx, int m, int n, void *adata )
301 {
303
304 size_t idx_mat, idx_dat;
305
306 idx_dat = 0;
307
308 /*
309 p[0] = floor
310 p[1] = scale
311 p[2] = x0
312 p[3] = y0
313 p[4] = ps [(lam/D)/pix]
314 */
315
316 realT r;
317
318 for( int i = 0; i < arr->nx; i++ )
319 {
320 for( int j = 0; j < arr->ny; j++ )
321 {
322 idx_mat = i + j * arr->nx;
323
324 r = sqrt( pow( i - p[2], 2 ) + pow( j - p[3], 2 ) );
325
326 hx[idx_dat] =
328 static_cast<realT>( i ), static_cast<realT>( j ), p[0], p[1], p[2], p[3], p[4], arr->cenObs ) -
329 arr->data[idx_mat];
330
331 // hx[idx_dat] *= fabs( pow(arr->data[idx_mat],2));
332
333 idx_dat++;
334 }
335 }
336 }
337};
338
339///\ref levmarInterface fitter structure for the obstructed Airy pattern, including fitting platescale and central
340///obscuration.
341/** \ingroup airy_peak_fit
342 *
343 */
344template <typename _realT>
346{
347 typedef _realT realT;
348
349 static const int nparams = 6;
350
351 static void func( realT *p, realT *hx, int m, int n, void *adata )
352 {
354
355 size_t idx_mat, idx_dat;
356
357 idx_dat = 0;
358
359 /*
360 p[0] = floor
361 p[1] = scale
362 p[2] = x0
363 p[3] = y0
364 p[4] = ps [(lam/D)/pix]
365 p[5] = cen-obs
366 */
367
368 // realT r;
369
370 for( int i = 0; i < arr->nx; ++i )
371 {
372 for( int j = 0; j < arr->ny; ++j )
373 {
374 idx_mat = i + j * arr->nx;
375
376 // r = sqrt( pow( i-p[2],2) + pow(j-p[3],2));
377
378 hx[idx_dat] =
380 static_cast<realT>( i ), static_cast<realT>( j ), p[0], p[1], p[2], p[3], p[4], p[5] ) -
381 arr->data[idx_mat];
382
383 // hx[idx_dat] *= fabs( pow(arr->data[idx_mat],2));
384
385 idx_dat++;
386 }
387 }
388 }
389};
390
391/// Alias for the fitAiry2D type with both platescale and cen-obs fixed.
392/** \ingroup airy_peak_fit
393 */
394template <typename realT>
396
397/// Alias for the fitAiry2D type with cen-obs fixed.
398/** \ingroup airy_peak_fit
399 */
400template <typename realT>
402
403/// Alias for the fitAiry2D type with none fixed.
404/** \ingroup airy_peak_fit
405 */
406template <typename realT>
408
409} // namespace fit
410} // namespace math
411} // namespace mx
412
413#endif // fitAiry_hpp
Class to manage fitting a 2D Airy pattern to data via the levmarInterface.
Definition fitAiry.hpp:72
void setGuess(realT A0, realT A, realT x0, realT y0, realT ps)
Set the initial guess when platescale is variable, and central obscuration is fixed.
Definition fitAiry.hpp:112
void setGuess(realT A0, realT A, realT x0, realT y0)
Set the initial guess when platescale and central obscuration are fixed.
Definition fitAiry.hpp:97
void setGuess(realT A0, realT A, realT x0, realT y0, realT ps, realT co)
Set the initial guess when central-obscuration is variable.
Definition fitAiry.hpp:129
A templatized interface to the levmar package.
void allocate_params()
Allocate parameters array based on previous call to nParams.
realT * p
Parameter array. On input is the initial estimates. On output has the estimated solution.
int n
I: measurement vector dimension.
void * adata
Pointer to possibly additional data, passed uninterpreted to func & jacf.
realT airyPattern(realT x)
The classical Airy pattern.
constexpr floatT six_fifths()
Return 6/5 in the specified precision.
A c++ interface to the templatized levmar minimization routines..
The mxlib c++ namespace.
Definition mxError.hpp:106
levmarInterface fitter structure for the obstructed Airy pattern, including platescale.
Definition fitAiry.hpp:295
levmarInterface fitter structure for the centrally obscured Airy pattern.
Definition fitAiry.hpp:239
Wrapper for a native array to pass to levmarInterface, with Airy details.
Definition fitAiry.hpp:221
size_t ny
Y dimension of the array.
Definition fitAiry.hpp:224
realT * data
Pointer to the array.
Definition fitAiry.hpp:222
realT cenObs
is the ratio of the circular central obscuration diameter to the diameter.
Definition fitAiry.hpp:226
realT ps
the platescale in
Definition fitAiry.hpp:227
size_t nx
X dimension of the array.
Definition fitAiry.hpp:223