ERIS Pipeline Reference Manual 1.9.2
sc_contsub.c
Go to the documentation of this file.
1/*
2 * This file is part of the SKYCORR software package.
3 * Copyright (C) 2009-2013 European Southern Observatory
4 *
5 * This programme is free software; you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation; either version 2 of the License, or
8 * (at your option) any later version.
9 *
10 * This programme is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this programme. If not, see <http://www.gnu.org/licenses/>.
17 */
18
19
39/*****************************************************************************
40 * INCLUDES *
41 ****************************************************************************/
42
43#include <sc_contsub.h>
44
45
46/*****************************************************************************
47 * CODE *
48 ****************************************************************************/
49
50cpl_error_code sc_contsub(cpl_table *spec)
51{
91 cpl_error_code err_code = CPL_ERROR_NONE;
92
93 cpl_table *tab;
94
95 double *xr, *yr, *xo, *yo; // pointers for interpolation
96
97 long nrows = cpl_table_get_nrow(spec), // number of rows in spec
98 ncont = 0, // number of continuum data points
99 i;
100
101 int npix = 11; // pixels for continuum median
102
103 // add new columns to spec
104 if (cpl_table_has_column(spec, "cflux") == 0) {
105 cpl_table_new_column(spec, "cflux", CPL_TYPE_DOUBLE);
106 }
107 if (cpl_table_has_column(spec, "lflux") == 0) {
108 cpl_table_new_column(spec, "lflux", CPL_TYPE_DOUBLE);
109 }
110
111 // allocate space for temporary variables
112 xr = (double *)calloc(nrows, sizeof(double));
113 yr = (double *)calloc(nrows, sizeof(double));
114 xo = (double *)calloc(nrows, sizeof(double));
115 yo = (double *)calloc(nrows, sizeof(double));
116
117 // count number of continuum and line data points
118 for (i = 0; i < nrows; i++) {
119 *(xo+i) = cpl_table_get_double(spec, "lambda", i, NULL);
120 *(yo+i) = 0;
121 if (cpl_table_get_int(spec, "class", i, NULL) == 0) {
122 *(xr+ncont) = cpl_table_get_double(spec, "lambda", i, NULL);
123 *(yr+ncont) = cpl_table_get_double(spec, "flux", i, NULL);
124 ncont++;
125 }
126 }
127
128 // smooth continuum points by median filtering
129 tab = cpl_table_new(ncont);
130 cpl_table_new_column(tab, "cflux", CPL_TYPE_DOUBLE);
131 for (i = 0; i < ncont; i++) {
132 cpl_table_set_double(tab, "cflux", i, *(yr+i));
133 }
134 sc_basic_filtermedian(tab, "cflux", npix);
135 for (i = 0; i < ncont; i++) {
136 *(yr+i) = cpl_table_get_double(tab, "cflux", i, NULL);
137 }
138 cpl_table_delete(tab);
139
140 // interpolate over line regions
141 err_code = sc_basic_interpollin(xo, yo, nrows, xr, yr, ncont, 0);
142
143 // write results back to spectrum
144 for (i = 0; i < nrows; i++) {
145 if (cpl_table_get_int(spec, "class", i, NULL) == 0) {
146 // keep original flux at continuum region
147 cpl_table_set_double(spec, "cflux", i,
148 cpl_table_get_double(spec, "flux", i, NULL));
149 } else {
150 // set interpolated flux at line region
151 cpl_table_set_double(spec, "cflux", i, *(yo + i));
152 }
153 // fill line flux region
154 cpl_table_set_double(spec, "lflux", i,
155 cpl_table_get_double(spec, "flux", i, NULL));
156 }
157
158 // subtract continuum flux from line region
159 cpl_table_subtract_columns(spec, "lflux", "cflux");
160
161 // cleanup
162 free(xr);
163 free(yr);
164 free(xo);
165 free(yo);
166
167 return err_code;
168}
169
170
171cpl_error_code sc_contsub_identcont(cpl_table *spec, const cpl_table *groups,
172 const cpl_table *linetab,
173 const double fluxlim,
174 const cpl_parameterlist *parlist)
175{
199 cpl_table *selgroups, *subtab = NULL;
200 const cpl_parameter *p;
201
202 int i, j, // loop variables
203 mid, // central line wavelength
204 n = 0, // # of rows in table spec
205 m = 0, // # of rows in table groups
206 varfwhm = 0; // wavelength-dependent FWHM?
207
208 double fwhm = 0, // FWHM of lines
209 ratio = 0, // flux ratio
210 minflux = 0, // minimum flux for lines from list
211 nsig = 3., // line/cont separation extension based on line width
212 lam_spec0, lam_spec1, // neighbouring wavelengths in spec
213 lam_list=0., // wavelength in line list
214 dlam = 0, // wavelength step
215 meanlam = 0, // mean wavelength
216 flux = 0, // flux
217 lpix = 0, // # of line pixels
218 meanflux = 0., // mean list line flux
219 meanlflux = 0., // mean observed line flux
220 lflux = 0.; // observed line flux
221
222 /* Copy line group table */
223 selgroups = cpl_table_duplicate(groups);
224
225 /* Number of elements in tables */
226 n = cpl_table_get_nrow(spec);
227 m = cpl_table_get_nrow(selgroups);
228
229 /* Clip groups to proper wavelength interval */
230 lam_spec0 = cpl_table_get_column_min(spec, "lambda");
231 lam_spec1 = cpl_table_get_column_max(spec, "lambda");
232 cpl_table_unselect_all(selgroups);
233 for (i = 0; i < m; i++) {
234 lam_list = cpl_table_get_double(selgroups, "lambda", i, NULL);
235 if (lam_list < lam_spec0 || lam_list > lam_spec1) {
236 cpl_table_select_row(selgroups, i);
237 }
238 }
239 cpl_table_erase_selected(selgroups);
240
241 /* Update number of elements in selgroups */
242 m = cpl_table_get_nrow(selgroups);
243
244 /* delta(lambda): required for converting FWHM from [pixel] to [micron] */
245 dlam = cpl_table_get_double(spec, "lambda", 1, NULL) -
246 cpl_table_get_double(spec, "lambda", 0, NULL);
247
248 /* Get FWHM */
249 p = cpl_parameterlist_find_const(parlist, "fwhm");
250 fwhm = cpl_parameter_get_double(p);
251
252 /* Variable FWHM (linear change with wavelength)? */
253 p = cpl_parameterlist_find_const(parlist, "varfwhm");
254 varfwhm = cpl_parameter_get_int(p);
255
256 /* Get mean wavelength */
257 p = cpl_parameterlist_find_const(parlist, "meanlam");
258 meanlam = cpl_parameter_get_double(p);
259
260 /* Convert flux in group list from integrated to peak flux */
261 for (i = 0; i < m; i++) {
262 flux = cpl_table_get_double(selgroups, "flux", i, NULL) /
263 CPL_MATH_SQRT2PI / CPL_MATH_SIG_FWHM / fwhm / dlam;
264 /* Correct flux if FWHM is wavelength dependent */
265 if (varfwhm == 1) {
266 lam_list = cpl_table_get_double(selgroups, "lambda", i, NULL);
267 flux *= meanlam / lam_list;
268 }
269 cpl_table_set_double(selgroups, "flux", i, flux);
270 }
271
272 /* Rebin selgroups spectrum */
273 subtab = cpl_table_duplicate(spec);
274 cpl_table_new_column(subtab, "flux_list", CPL_TYPE_DOUBLE);
275 sc_basic_rebin(subtab, "lambda", "flux_list", selgroups, "lambda",
276 "flux");
277
278 /* Apply scaling */
279 meanflux = cpl_table_get_column_mean(subtab, "flux_list");
280 ratio = cpl_table_get_column_mean(subtab, "lflux") / meanflux;
281 if (ratio < 0) {
282 /* Avoid negative ratios */
283 for (meanlflux = 0., i = 0; i < n; i++) {
284 if ((lflux = cpl_table_get(subtab, "lflux", i, NULL)) > 0) {
285 meanlflux += lflux;
286 }
287 }
288 ratio = meanlflux / meanflux;
289 }
290 cpl_table_multiply_scalar(selgroups, "flux", ratio);
291
292 /* Get minimum line flux limit in spectrum */
293 if (cpl_table_get_nrow(linetab) == 0) {
294 /* Solution for no identified lines */
295 minflux = fluxlim * cpl_table_get_column_mean(selgroups, "flux");
296 } else {
297 /* Fraction of median peak flux of identified lines */
298 minflux = fluxlim * cpl_table_get_column_median(linetab, "peak_flux");
299 /* Make sure that lines from the list are above the threshold */
300 if (minflux > cpl_table_get_column_max(selgroups, "flux")) {
301 minflux = fluxlim * cpl_table_get_column_mean(selgroups, "flux");
302 }
303 }
304
305 /* Add line centres from line list to spectrum */
306
307 for (i = 0, j = 0; i < n-1; ) {
308
309 /* Check end of line list */
310 if (j >= m-1) {
311 break;
312 }
313
314 /* Get two consecutive wavelengths from spectrum */
315 lam_spec0 = cpl_table_get_double(spec, "lambda", i, NULL);
316 lam_spec1 = cpl_table_get_double(spec, "lambda", i+1, NULL);
317
318 /* Get wavelength and flux in line list */
319 lam_list = cpl_table_get_double(selgroups, "lambda", j, NULL);
320 flux = cpl_table_get_double(selgroups, "flux", j, NULL);
321
322 if (lam_list < lam_spec0) {
323 /* Skip lines at low end of wavelength range in spectrum */
324 j++;
325 continue;
326 }
327
328 if (flux < minflux) {
329 /* Use only significant lines */
330 j++;
331 continue;
332 }
333
334 if (lam_list < lam_spec1) {
335 /* Current line in list has wavelength >= lam_spec0
336 but < lam_spec1 -> find closer wavelength */
337 if (lam_list-lam_spec0 < lam_spec1-lam_list) {
338 /* Closer to first wavelength */
339 mid = i;
340 } else {
341 /* Closer to second wavelength */
342 mid = i+1;
343 }
344 /* Identify peak in mask (no change of flag for isolated lines) */
345 if (cpl_table_get_int(spec, "class", mid, NULL) != 3) {
346 cpl_table_set_int(spec, "class", mid, 2);
347 }
348
349 /* Next line in list */
350 j++;
351 }
352
353 if (lam_list > lam_spec0) {
354 /* Move to next wavelength in spectrum */
355 i++;
356 }
357
358 }
359
360 /* Loop over all pixels */
361
362 for (i = 0; i < n-1; i++) {
363
364 if (cpl_table_get_int(spec, "class", i, NULL) >= 2) {
365
366 /* Peak found */
367
368 /* Extend line regions in spectrum around peaks
369 plus/minus nsig sigma [pixel] */
370 lpix = nsig * fwhm / CPL_MATH_FWHM_SIG;
371
372 /* Get peak flux and wavelength */
373 flux = cpl_table_get_double(spec, "flux", i, NULL);
374 lam_list = cpl_table_get_double(spec, "lambda", i, NULL);
375
376 /* Correct peak flux if FWHM is wavelength dependent */
377 if (varfwhm == 1) {
378 flux *= lam_list / meanlam;
379 }
380
381 /* Strong lines are broader at the base */
382 if (flux > 10 * fluxlim) {
383 lpix *= 2;
384 }
385
386 /* Correct line regions if FWHM is wavelength dependent */
387 if (varfwhm == 1) {
388 lpix *= lam_list / meanlam;
389 }
390
391 /* Surround peak with line pixels */
392
393 for (j = 1; j < (int) floor(lpix + 0.5); j++) {
394 /* Lower wavelength end */
395 if (i-j >= 0) {
396 if (cpl_table_get_int(spec, "class", i-j, NULL) < 1) {
397 /* Change continuum pixels only */
398 cpl_table_set_int(spec, "class", i-j, 1);
399 }
400 }
401 /* Higher wavelength end */
402 if (i+j < n) {
403 if (cpl_table_get_int(spec, "class", i+j, NULL) < 1) {
404 /* Change continuum pixels only */
405 cpl_table_set_int(spec, "class", i+j, 1);
406 }
407 }
408 }
409
410 }
411
412 }
413
414 /* Cleanup */
415 cpl_table_delete(subtab);
416 cpl_table_delete(selgroups);
417
418 return cpl_error_get_code();
419}
420
421
422int sc_contsub_check(cpl_table *spec)
423{
441 cpl_table *subtab;
442 int nall = 0, ncont = 0, flag = 0;
443 double minlam = 0., maxlam = 0., range = 0., crange = 0.;
444
445 /* Parameters for checking */
446 int count = 5; // continuum pixel to be evaluated (at both margins)
447 double minfpix = 0.2; // minimum fraction of continuum pixels
448 double minfrange = 0.9; // minimum fraction of continuum wavelength range
449
450 /* Get total number of pixels */
451 nall = cpl_table_get_nrow(spec);
452 if (nall == 0) {
453 return 1;
454 }
455
456 /* Extract subtable of continuum pixels */
457 cpl_table_unselect_all(spec);
458 cpl_table_or_selected_int(spec, "class", CPL_EQUAL_TO, 0);
459 subtab = cpl_table_extract_selected(spec);
460 cpl_table_select_all(spec);
461
462 /* Check fraction of continuum pixels */
463 ncont = cpl_table_get_nrow(subtab);
464 if ((double) ncont / (double) nall < minfpix) {
465 flag = 1;
466 }
467
468 /* Get full wavelength range covered by the spectrum */
469 minlam = cpl_table_get(spec, "lambda", 0, NULL);
470 maxlam = cpl_table_get(spec, "lambda", nall-1, NULL);
471 range = maxlam - minlam;
472
473 /* Check wavelength range covered by continuum pixels */
474 if (ncont >= 2 * count) {
475 minlam = cpl_table_get(subtab, "lambda", count-1, NULL);
476 maxlam = cpl_table_get(subtab, "lambda", ncont-count, NULL);
477 crange = maxlam - minlam;
478 if (crange / range < minfrange) {
479 flag = 1;
480 }
481 } else {
482 flag = 1;
483 }
484
485 /* Cleanup */
486 cpl_table_delete(subtab);
487
488 cpl_error_code err = cpl_error_get_code();
489 if (err != CPL_ERROR_NONE) {
490 cpl_errorstate_dump(0, CPL_FALSE, NULL);
491 }
492
493 return flag;
494}
495
496
497cpl_error_code sc_contsub_skipthermir(cpl_table *spec)
498{
514 int n = 0, i = 0, mask = 0;
515 double lam = 0.;
516
517 /* Get number of pixels in spectrum */
518 n = cpl_table_get_nrow(spec);
519
520 /* Beyond SC_THERMIRLIM set all mask flags to 0 */
521 for (i = n - 1; i >= 0; i--) {
522 lam = cpl_table_get(spec, "lambda", i, NULL);
523 mask = cpl_table_get(spec, "class", i, NULL);
524 if (lam > SC_THERMIRLIM && mask > 0) {
525 cpl_table_set(spec, "class", i, 0);
526 }
527 }
528
529 return cpl_error_get_code();
530}
531
532
cpl_error_code sc_basic_rebin(cpl_table *outspec, const char *outlam, const char *outflux, const cpl_table *inspec, const char *inlam, const char *influx)
Definition: sc_basic.c:296
#define SC_THERMIRLIM
Definition: sc_basic.h:102
cpl_error_code sc_basic_filtermedian(cpl_table *spec, const char *colname, const int npix)
Definition: sc_basic.c:745
cpl_error_code sc_basic_interpollin(const double *x_out, double *y_out, const long n_out, const double *x_ref, const double *y_ref, const long n_ref, const int extrapolate)
Definition: sc_basic.c:2369
int sc_contsub_check(cpl_table *spec)
Definition: sc_contsub.c:422
cpl_error_code sc_contsub_skipthermir(cpl_table *spec)
Definition: sc_contsub.c:497
cpl_error_code sc_contsub_identcont(cpl_table *spec, const cpl_table *groups, const cpl_table *linetab, const double fluxlim, const cpl_parameterlist *parlist)
Definition: sc_contsub.c:171
cpl_error_code sc_contsub(cpl_table *spec)
Definition: sc_contsub.c:50