mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
Loading...
Searching...
No Matches
airyPattern.hpp
Go to the documentation of this file.
1/** \file airyPattern.hpp
2 * \author Jared R. Males
3 * \brief Utilities related to the Airy pattern point spread function.
4 * \ingroup gen_math_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 math_func_airyPattern_hpp
28#define math_func_airyPattern_hpp
29
30#include <gsl/gsl_integration.h>
31#include <gsl/gsl_errno.h>
32
33#include "../constants.hpp"
34#include "bessel.hpp"
35#include "jinc.hpp"
36
37namespace mx
38{
39
40namespace math
41{
42
43namespace func
44{
45
46/// The classical Airy pattern
47/** Returns the intensity distribution of the Airy pattern at a given \f$ \lambda/D \f$
48 *
49 * References: \cite born_and_wolf, \cite mahajan_1986, https://en.wikipedia.org/wiki/Airy_disk.
50 *
51 * \tparam realT is the floating point type used for arithmetic.
52 *
53 * \ingroup psfs
54 * \ingroup gen_math_airy_pattern
55 */
56template <typename realT>
57realT airyPattern( realT x /**< [in] is the separation in units of \f$ \lambda/D \f$. */ )
58{
59 return pow( 2 * jinc( pi<realT>() * x ), 2 );
60}
61
62/// The centrally obscured Airy pattern
63/** Returns the intensity distribution of the centrally obscured Airy pattern at a given \f$ \lambda/D \f$
64 *
65 * References: \cite born_and_wolf, \cite mahajan_1986, https://en.wikipedia.org/wiki/Airy_disk.
66 *
67 * \tparam realT is the floating point type used for arithmetic.
68 *
69 * \ingroup psfs
70 * \ingroup gen_math_airy_pattern
71 */
72template <typename realT>
73realT airyPattern( realT x, ///< [in] is the separation in units of \f$ \lambda/D \f$.
74 realT eps ///< [in] is the ratio of the circular central obscuration diameter to the diameter.
75)
76{
77 return ( 1. / pow( 1. - eps * eps, 2 ) ) *
78 pow( 2. * jinc( pi<realT>() * x ) - eps * eps * 2. * jinc( eps * pi<realT>() * x ), 2 );
79}
80
81/// The general centrally obscured Airy pattern, with arbitrary center and platescale.
82/** Returns the intensity distribution of the centrally obscured Airy pattern at a given \f$ \lambda/D \f$. This
83 * version allows for an arbitrary center, scaling, and platescale.
84 *
85 * References: \cite born_and_wolf, \cite mahajan_1986, https://en.wikipedia.org/wiki/Airy_disk.
86 *
87 * \tparam realT is the floating point type used for arithmetic.
88 *
89 * \ingroup psfs
90 * \ingroup gen_math_airy_pattern
91 */
92template <typename realT>
93realT airyPattern( realT x, ///< [in] is the x-coordinate in units of pixels
94 realT y, ///< [in] is the y-coordinate in units of pixels
95 realT A0, ///< [in] constant value added to the Airy pattern
96 realT A, ///< [in] peak scale of the Airy pattern.
97 realT x0, ///< [in] is the x-center in units of pixels
98 realT y0, ///< [in] is the y-center in units of pixels
99 realT ps, ///< [in] the platescale in \f$ (\lambda/D)/pixel \f$
100 realT eps ///< [in] is the ratio of the circular central obscuration diameter to the diameter.
101)
102{
103 realT r = sqrt( pow( x - x0, 2 ) + pow( y - y0, 2 ) ) * ps;
104
105 return A0 + A * airyPattern( r, eps );
106}
107
108/// Fill in an array with the 2D arbitrarily-centered classical Airy pattern
109/**
110 * At each pixel (x,y) of the array this computes:
111 *
112 * \f$ f(x,y) = A_0 + A\times Airy(\sqrt(x^2+y^2) \f$
113 *
114 * \tparam realT is the type to use for arithmetic
115 *
116 * \ingroup gen_math_airy_pattern
117 */
118template <typename realT>
119void airyPattern2D( realT *arr, ///< [out] is the allocated array to fill in
120 size_t nx, ///< [in] is the size of the x dimension of the array (rows)
121 size_t ny, ///< [in] is the size of the y dimension of the array (columns)
122 const realT A0, ///< [in] is the constant to add to the Gaussian
123 const realT A, ///< [in] is the scaling factor (peak height = A-A0)
124 const realT x0, ///< [in] is the x-coordinate of the center
125 const realT y0, ///< [in] is the y-coordinate of the center
126 realT ps ///< [in] the platescale in \f$ (\lambda/D)/pixel \f$
127)
128{
129 size_t idx;
130
131 for( size_t j = 0; j < ny; ++j )
132 {
133 for( size_t i = 0; i < nx; ++i )
134 {
135 idx = i + j * nx;
136
137 realT rad = sqrt( pow( i - x0, 2 ) + pow( j - y0, 2 ) ) * ps;
138
139 arr[idx] = A0 + A * airyPattern( rad );
140 }
141 }
142}
143
144/// Seeing Halo Profile
145/** A Moffat profile due to Roddier (1981)\cite roddier_1981, see also Racine et al (1999)\cite racine_1999, which
146 * can be used to describe the seeing limited point-spread function of a telescope imaging through turbulence, or the
147 * halo in partially corrected imaging.
148 *
149 * \returns the value of the profile at x, in units of fractional flux per unit area.
150 *
151 * \ingroup psfs
152 * \ingroup gen_math_airy_pattern
153 */
154template <typename realT>
156 realT x, ///< [in] the separation in the same units as fwhm.
157 realT fwhm ///< [in] the fwhm in arbitrary units. Note that this defines the area units of the density.
158)
159{
160 return ( 0.488 / ( fwhm * fwhm ) ) * pow( 1. + ( 11. / 6. ) * pow( x / fwhm, 2 ), -11. / 6. );
161}
162
163/// Calculate the fraction of enclosed power at a given radius for the unobscured Airy Pattern
164/** See Mahajan (1986) \cite mahajan_1986 and https://en.wikipedia.org/wiki/Airy_disk.
165 *
166 * \returns the fraction of enclosed power at radius x
167 *
168 * \ingroup psfs
169 * \ingroup gen_math_airy_pattern
170 */
171template <typename realT>
172realT airyPatternEnclosed( realT x /**< [in] the radius */ )
173{
174 realT x1 = x * pi<realT>();
175
176 realT b0 = bessel_j<realT>( 0, x1 );
177 b0 = b0 * b0;
178
179 realT b1 = bessel_j<realT>( 1, x1 );
180 b1 = b1 * b1;
181
182 realT encp = static_cast<realT>( 1 ) - b0 - b1;
183
184 return encp;
185}
186
187template <typename realT>
188realT apeInt( realT x, void *params )
189{
190 realT eps = *static_cast<realT *>( params );
191
192 return bessel_j<realT>( 1, x ) * bessel_j<realT>( 1, eps * x ) / x;
193}
194
195/// Calculate the fraction of enclosed power at a given radius for the centrally obscured Airy Pattern
196/** See Mahajan (1986) \cite mahajan_1986 and https://en.wikipedia.org/wiki/Airy_disk.
197 * If eps = 0, this calls the unobscured version.
198 *
199 * \returns the fraction of enclosed power at radius x
200 *
201 * \ingroup psfs
202 * \ingroup gen_math_airy_pattern
203 */
204template <typename realT>
205realT airyPatternEnclosed( realT x, ///< [in] the radius
206 realT eps ///< [in] the central obscuration fraction
207)
208{
209 if( eps == 0 )
210 return airyPatternEnclosed( x );
211
212 gsl_function func;
213 func.function = apeInt<realT>;
214 func.params = &eps;
215
216 realT jint;
217 realT abserr;
218 size_t neval;
219
220 realT x1 = x * pi<realT>();
221
223 gsl_integration_qng( &func, 0, x1, 1e-7, 1e-7, &jint, &abserr, &neval );
224
225 realT b0 = bessel_j<realT>( 0, x1 );
226 b0 = b0 * b0;
227 realT b1 = bessel_j<realT>( 1, x1 );
228 b1 = b1 * b1;
229
230 realT b0e = bessel_j<realT>( 0, eps * x1 );
231 b0e = b0e * b0e;
232 realT b1e = bessel_j<realT>( 1, eps * x1 );
233 b1e = b1e * b1e;
234
235 realT eps2 = pow( eps, 2 );
236
237 realT encp = static_cast<realT>( 1 ) - b0 - b1 + eps2 * ( static_cast<realT>( 1 ) - b0e - b1e );
238
239 encp = encp - 4 * eps * jint;
240 encp = encp / ( static_cast<realT>( 1 ) - eps2 );
241
242 return encp;
243}
244
245} // namespace func
246} // namespace math
247} // namespace mx
248
249#endif // math_func_airyPattern_hpp
Declares and defines Bessel functions of the first kind.
realT airyPattern(realT x)
The classical Airy pattern.
realT airyPatternEnclosed(realT x)
Calculate the fraction of enclosed power at a given radius for the unobscured Airy Pattern.
realT seeingHalo(realT x, realT fwhm)
Seeing Halo Profile.
void airyPattern2D(realT *arr, size_t nx, size_t ny, const realT A0, const realT A, const realT x0, const realT y0, realT ps)
Fill in an array with the 2D arbitrarily-centered classical Airy pattern.
T jinc(const T &x)
The Jinc function.
Definition jinc.hpp:61
constexpr floatT six_fifths()
Return 6/5 in the specified precision.
Declares and defines the Jinc and Jinc2 functions.
The mxlib c++ namespace.
Definition mxError.hpp:106