Kstars

sep.h
1/*
2 This file is part of SEP
3
4 SPDX-FileCopyrightText: 1993-2011 Emmanuel Bertin -- IAP /CNRS/UPMC
5 SPDX-FileCopyrightText: 2014 SEP developers
6
7 SPDX-License-Identifier: LGPL-3.0-or-later
8*/
9
10/* datatype codes */
11#define SEP_TBYTE 11 /* 8-bit unsigned byte */
12#define SEP_TINT 31 /* native int type */
13#define SEP_TFLOAT 42
14#define SEP_TDOUBLE 82
15
16/* object & aperture flags */
17#define SEP_OBJ_MERGED 0x0001 /* object is result of deblending */
18#define SEP_OBJ_TRUNC 0x0002 /* object truncated at image boundary */
19#define SEP_OBJ_DOVERFLOW 0x0004 /* not currently used, but could be */
20#define SEP_OBJ_SINGU 0x0008 /* x,y fully correlated */
21#define SEP_APER_TRUNC 0x0010
22#define SEP_APER_HASMASKED 0x0020
23#define SEP_APER_ALLMASKED 0x0040
24#define SEP_APER_NONPOSITIVE 0x0080
25
26/* noise_type values in sep_image */
27#define SEP_NOISE_NONE 0
28#define SEP_NOISE_STDDEV 1
29#define SEP_NOISE_VAR 2
30
31/* input flags for aperture photometry */
32#define SEP_MASK_IGNORE 0x0004
33
34/* threshold interpretation for sep_extract */
35#define SEP_THRESH_REL 0 /* in units of standard deviations (sigma) */
36#define SEP_THRESH_ABS 1 /* absolute data values */
37
38/* filter types for sep_extract */
39#define SEP_FILTER_CONV 0
40#define SEP_FILTER_MATCHED 1
41
42/* structs ------------------------------------------------------------------*/
43
44/* sep_image
45 *
46 * Represents an image, including data, noise and mask arrays, and
47 * gain.
48 */
49typedef struct {
50 void *data; /* data array */
51 void *noise; /* noise array (can be NULL) */
52 void *mask; /* mask array (can be NULL) */
53 int dtype; /* element type of image */
54 int ndtype; /* element type of noise */
55 int mdtype; /* element type of mask */
56 int w; /* array width */
57 int h; /* array height */
58 double noiseval; /* scalar noise value; used only if noise == NULL */
59 short noise_type; /* interpretation of noise value */
60 double gain; /* (poisson counts / data unit) */
61 double maskthresh; /* pixel considered masked if mask > maskthresh */
62} sep_image;
63
64/* sep_bkg
65 *
66 * The result of sep_background() -- represents a smooth image background
67 * and its noise with splines.
68 */
69typedef struct {
70 int w, h; /* original image width, height */
71 int bw, bh; /* single tile width, height */
72 int nx, ny; /* number of tiles in x, y */
73 int n; /* nx*ny */
74 float global; /* global mean */
75 float globalrms; /* global sigma */
76 float *back; /* node data for interpolation */
77 float *dback;
78 float *sigma;
79 float *dsigma;
80} sep_bkg;
81
82/* sep_catalog
83 *
84 * The result of sep_extract(). This is a struct of arrays. Each array has
85 * one entry per detected object.
86 */
87typedef struct {
88 int nobj; /* number of objects (length of all arrays) */
89 float *thresh; /* threshold (ADU) */
90 int *npix; /* # pixels extracted (size of pix array) */
91 int *tnpix; /* # pixels above thresh (unconvolved) */
92 int *xmin, *xmax;
93 int *ymin, *ymax;
94 double *x, *y; /* barycenter (first moments) */
95 double *x2, *y2, *xy; /* second moments */
96 double *errx2, *erry2, *errxy; /* second moment errors */
97 float *a, *b, *theta; /* ellipse parameters */
98 float *cxx, *cyy, *cxy; /* ellipse parameters (alternative) */
99 float *cflux; /* total flux of pixels (convolved im) */
100 float *flux; /* total flux of pixels (unconvolved) */
101 float *cpeak; /* peak intensity (ADU) (convolved) */
102 float *peak; /* peak intensity (ADU) (unconvolved) */
103 int *xcpeak, *ycpeak; /* x, y coords of peak (convolved) pixel */
104 int *xpeak, *ypeak; /* x, y coords of peak (unconvolved) pixel */
105 short *flag; /* extraction flags */
106 int **pix; /* array giving indicies of object's pixels in */
107 /* image (linearly indexed). Length is `npix`. */
108 /* (pointer to within the `objectspix` buffer) */
109 int *objectspix; /* buffer holding pixel indicies for all objects */
110} sep_catalog;
111
112#ifdef __cplusplus
113extern "C" {
114#endif
115
116/*--------------------- global background estimation ------------------------*/
117
118/* sep_background()
119 *
120 * Create representation of spatially varying image background and variance.
121 *
122 * Note that the returned pointer must eventually be freed by calling
123 * `sep_bkg_free()`.
124 *
125 * In addition to the image mask (if present), pixels <= -1e30 and NaN
126 * are ignored.
127 *
128 * Source Extractor defaults:
129 *
130 * - bw, bh = (64, 64)
131 * - fw, fh = (3, 3)
132 * - fthresh = 0.0
133 */
134int sep_background(sep_image *image,
135 int bw, int bh, /* size of a single background tile */
136 int fw, int fh, /* filter size in tiles */
137 double fthresh, /* filter threshold */
138 sep_bkg **bkg); /* OUTPUT */
139
140
141/* sep_bkg_global[rms]()
142 *
143 * Get the estimate of the global background "median" or standard deviation.
144 */
145float sep_bkg_global(sep_bkg *bkg);
146float sep_bkg_globalrms(sep_bkg *bkg);
147
148
149/* sep_bkg_pix()
150 *
151 * Return background at (x, y).
152 * Unlike other routines, this uses simple linear interpolation.
153 */
154float sep_bkg_pix(sep_bkg *bkg, int x, int y);
155
156
157/* sep_bkg_[sub,rms]line()
158 *
159 * Evaluate the background or RMS at line `y`.
160 * Uses bicubic spline interpolation between background map vertices.
161 * The second function subtracts the background from the input array.
162 * Line must be an array with same width as original image.
163 */
164int sep_bkg_line(sep_bkg *bkg, int y, void *line, int dtype);
165int sep_bkg_subline(sep_bkg *bkg, int y, void *line, int dtype);
166int sep_bkg_rmsline(sep_bkg *bkg, int y, void *line, int dtype);
167
168
169/* sep_bkg_[sub,rms]array()
170 *
171 * Evaluate the background or RMS for entire image.
172 * Uses bicubic spline interpolation between background map vertices.
173 * The second function subtracts the background from the input array.
174 * `arr` must be an array of the same size as original image.
175 */
176int sep_bkg_array(sep_bkg *bkg, void *arr, int dtype);
177int sep_bkg_subarray(sep_bkg *bkg, void *arr, int dtype);
178int sep_bkg_rmsarray(sep_bkg *bkg, void *arr, int dtype);
179
180/* sep_bkg_free()
181 *
182 * Free memory associated with bkg.
183 */
184void sep_bkg_free(sep_bkg *bkg);
185
186/*-------------------------- source extraction ------------------------------*/
187
188/* sep_extract()
189 *
190 * Extract sources from an image. Source Extractor defaults are shown
191 * in [ ] above.
192 *
193 * Notes
194 * -----
195 * `dtype` and `ndtype` indicate the data type (float, int, double) of the
196 * image and noise arrays, respectively.
197 *
198 * If `noise` is NULL, thresh is interpreted as an absolute threshold.
199 * If `noise` is not null, thresh is interpreted as a relative threshold
200 * (the absolute threshold will be thresh*noise[i,j]).
201 *
202 */
203int sep_extract(sep_image *image,
204 float thresh, /* detection threshold [1.5] */
205 int thresh_type, /* threshold units [SEP_THRESH_REL] */
206 int minarea, /* minimum area in pixels [5] */
207 float *conv, /* convolution array (can be NULL) */
208 /* [{1 2 1 2 4 2 1 2 1}] */
209 int convw, int convh, /* w, h of convolution array [3,3] */
210 int filter_type, /* convolution (0) or matched (1) [0] */
211 int deblend_nthresh, /* deblending thresholds [32] */
212 double deblend_cont, /* min. deblending contrast [0.005] */
213 int clean_flag, /* perform cleaning? [1] */
214 double clean_param, /* clean parameter [1.0] */
215 sep_catalog **catalog); /* OUTPUT catalog */
216
217
218
219/* set and get the size of the pixel stack used in extract() */
220void sep_set_extract_pixstack(size_t val);
221size_t sep_get_extract_pixstack(void);
222
223/* free memory associated with a catalog */
224void sep_catalog_free(sep_catalog *catalog);
225
226/*-------------------------- aperture photometry ----------------------------*/
227
228
229/* Sum array values within a circular aperture.
230 *
231 * Notes
232 * -----
233 * error : Can be a scalar (default), an array, or NULL
234 * If an array, set the flag SEP_ERROR_IS_ARRAY in `inflags`.
235 * Can represent 1-sigma std. deviation (default) or variance.
236 * If variance, set the flag SEP_ERROR_IS_VARIANCE in `inflags`.
237 *
238 * gain : If 0.0, poisson noise on sum is ignored when calculating error.
239 * Otherwise, (sum / gain) is added to the variance on sum.
240 *
241 * area : Total pixel area included in sum. Includes masked pixels that were
242 * corrected. The area can differ from the exact area of a circle due
243 * to inexact subpixel sampling and intersection with array boundaries.
244 */
245int sep_sum_circle(sep_image *image,
246 double x, /* center of aperture in x */
247 double y, /* center of aperture in y */
248 double r, /* radius of aperture */
249 int subpix, /* subpixel sampling */
250 short inflags, /* input flags (see below) */
251 double *sum, /* OUTPUT: sum */
252 double *sumerr, /* OUTPUT: error on sum */
253 double *area, /* OUTPUT: area included in sum */
254 short *flag); /* OUTPUT: flags */
255
256
257int sep_sum_circann(sep_image *image,
258 double x, double y, double rin, double rout,
259 int subpix, short inflags,
260 double *sum, double *sumerr, double *area, short *flag);
261
262int sep_sum_ellipse(sep_image *image,
263 double x, double y, double a, double b, double theta,
264 double r, int subpix, short inflags,
265 double *sum, double *sumerr, double *area, short *flag);
266
267int sep_sum_ellipann(sep_image *image,
268 double x, double y, double a, double b, double theta,
269 double rin, double rout, int subpix, short inflags,
270 double *sum, double *sumerr, double *area, short *flag);
271
272/* sep_sum_circann_multi()
273 *
274 * Sum an array of circular annuli more efficiently (but with no exact mode).
275 *
276 * Notable parameters:
277 *
278 * rmax: Input radii are [rmax/n, 2*rmax/n, 3*rmax/n, ..., rmax].
279 * n: Length of input and output arrays.
280 * sum: Preallocated array of length n holding sums in annuli. sum[0]
281 * corresponds to r=[0, rmax/n], sum[n-1] to outermost annulus.
282 * sumvar: Preallocated array of length n holding variance on sums.
283 * area: Preallocated array of length n holding area summed in each annulus.
284 * maskarea: Preallocated array of length n holding masked area in each
285 annulus (if mask not NULL).
286 * flag: Output flag (non-array).
287 */
288int sep_sum_circann_multi(sep_image *im,
289 double x, double y, double rmax, int n, int subpix,
290 short inflag,
291 double *sum, double *sumvar, double *area,
292 double *maskarea, short *flag);
293
294/* sep_flux_radius()
295 *
296 * Calculate the radii enclosing the requested fraction of flux relative
297 * to radius rmax.
298 *
299 * (see previous functions for most arguments)
300 * rmax : maximum radius to analyze
301 * fluxtot : scale requested flux fractions to this. (If NULL, flux within
302 `rmax` is used.)
303 * fluxfrac : array of requested fractions.
304 * n : length of fluxfrac
305 * r : (output) result array of length n.
306 * flag : (output) scalar flag
307 */
308int sep_flux_radius(sep_image *im,
309 double x, double y, double rmax, int subpix, short inflag,
310 double *fluxtot, double *fluxfrac, int n,
311 double *r, short *flag);
312
313/* sep_kron_radius()
314 *
315 * Calculate Kron radius within an ellipse given by
316 *
317 * cxx*(x'-x)^2 + cyy*(y'-y)^2 + cxy*(x'-x)*(y'-y) < r^2
318 *
319 * The Kron radius is sum(r_i * v_i) / sum(v_i) where v_i is the value of pixel
320 * i and r_i is the "radius" of pixel i, as given by the left hand side of
321 * the above equation.
322 *
323 * Flags that might be set:
324 * SEP_APER_HASMASKED - at least one of the pixels in the ellipse is masked.
325 * SEP_APER_ALLMASKED - All pixels in the ellipse are masked. kronrad = 0.
326 * SEP_APER_NONPOSITIVE - There was a non-positive numerator or denominator.
327 * kronrad = 0.
328 */
329int sep_kron_radius(sep_image *im, double x, double y,
330 double cxx, double cyy, double cxy, double r,
331 double *kronrad, short *flag);
332
333
334/* sep_windowed()
335 *
336 * Calculate "windowed" position parameters via an iterative procedure.
337 *
338 * x, y : initial center
339 * sig : sigma of Gaussian to use for weighting. The integration
340 * radius is 4 * sig.
341 * subpix : Subpixels to use in aperture-pixel overlap.
342 * SExtractor uses 11. 0 is supported for exact overlap.
343 * xout, yout : output center.
344 * niter : number of iterations used.
345 */
346int sep_windowed(sep_image *im,
347 double x, double y, double sig, int subpix, short inflag,
348 double *xout, double *yout, int *niter, short *flag);
349
350
351/* sep_set_ellipse()
352 *
353 * Set array elements within an elliptical aperture to a given value.
354 *
355 * Ellipse: cxx*(x'-x)^2 + cyy*(y'-y)^2 + cxy*(x'-x)*(y'-y) = r^2
356 */
357void sep_set_ellipse(unsigned char *arr, int w, int h,
358 double x, double y, double cxx, double cyy, double cxy,
359 double r, unsigned char val);
360
361
362/* sep_ellipse_axes()
363 * sep_ellipse_coeffs()
364 *
365 * Convert between coefficient representation of ellipse,
366 * cxx*(x'-x)^2 + cyy*(y'-y)^2 + cxy*(x'-x)*(y'-y) = r^2,
367 * and axis representation of an ellipse. The axis representation is
368 * defined by:
369 *
370 * a = semimajor axis
371 * b = semiminor axis
372 * theta = angle in radians counter-clockwise from positive x axis
373 */
374int sep_ellipse_axes(double cxx, double cyy, double cxy,
375 double *a, double *b, double *theta);
376void sep_ellipse_coeffs(double a, double b, double theta,
377 double *cxx, double *cyy, double *cxy);
378
379/*----------------------- info & error messaging ----------------------------*/
380
381/* sep_version_string : library version (e.g., "0.2.0") */
382extern char *sep_version_string;
383
384/* sep_get_errmsg()
385 *
386 * Return a short descriptive error message that corresponds to the input
387 * error status value. The message may be up to 60 characters long, plus
388 * the terminating null character.
389 */
390void sep_get_errmsg(int status, char *errtext);
391
392
393/* sep_get_errdetail()
394 *
395 * Return a longer error message with more specifics about the problem.
396 * The message may be up to 512 characters.
397 */
398void sep_get_errdetail(char *errtext);
399
400#ifdef __cplusplus
401}
402#endif
const QList< QKeySequence > & back()
This file is part of the KDE documentation.
Documentation copyright © 1996-2024 The KDE developers.
Generated on Tue Mar 26 2024 11:19:03 by doxygen 1.10.0 written by Dimitri van Heesch, © 1997-2006

KDE's Doxygen guidelines are available online.