ERIS Pipeline Reference Manual 1.8.15
sc_specdiss.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
40/*****************************************************************************
41 * INCLUDES *
42 ****************************************************************************/
43
44#include <sc_specdiss.h>
45
46
47/*****************************************************************************
48 * DEFINES *
49 ****************************************************************************/
50
51#define EXIST F_OK /* check for file existence, unistd.h, access() */
52#define EXECUTABLE X_OK /* check for execute permissions */
53
54
55/*****************************************************************************
56 * CODE *
57 ****************************************************************************/
58
59cpl_error_code sc_specdiss_find_emissionlines(cpl_table *input_spectrum,
60 cpl_table *linetab,
61 cpl_parameterlist *parlist,
62 const cpl_table *groups)
63{
107 cpl_error_code status = CPL_ERROR_NONE;
108 char err_msg[SC_MAXLEN]; /* String containing error message */
109
110 int in_nrow=0; /* Length of input spectrum */
111 int loop=0; /* Loop variables */
112 double value=0.; /* Simple spectrum value */
113
114 double sig = 0.; /* sigma of 1st derivative */
115 double mean = 0.; /* mean of abs. val. of 1st deriv. */
116 double maxrat = 4.0; /* max. sig/mean for peaks */
117 double deriv1_lim=0.; /* 1st derivative peak criterion */
118
119 float lim_fac=1; /* threshold: lim_fac * deriv1_lim */
120 double flux_low_lim=0., flux_up_lim=0.; /* limits to search for */
121
122 /* Histogram stuff */
123 //cpl_table *hist1; /* CPL table, 2cols 'bins', 'counts'*/
124 //int n_bins=0; /* number of bins to be used */
125 //int bin_scale=30; /* # of histogram bins n_bins = #rows / bin_scale */
126 /* This approach is required to obtain a roughly */
127 /* equal sampling of the histogram */
128 /* NOTE: n_bins is an int! */
129 //char column_name[SC_MAXLEN]; /* Column name to be used for hist */
130 char col1[SC_MAXLEN], col2[SC_MAXLEN]; /* col name string */
131 char plot_title[SC_MAXLEN];
132
133 cpl_table *derivatives; /* Containing derivatives of inspec */
134 cpl_table *tmptable;
135 //cpl_parameterlist *plottags;
136 cpl_parameter *pptr;
137
138/*--------------------------------------------------------------------------*/
139/*------------------------------- INITIALISING -----------------------------*/
140/*--------------------------------------------------------------------------*/
141
142 /* checking error state */
143
144 /* Checking input spectrum */
145 if ( (cpl_table_has_column(input_spectrum, "lambda")!= 1 ) ||
146 (cpl_table_has_column(input_spectrum, "flux") != 1 ) )
147 {
148 sprintf(err_msg, "Either column \"lambda\" or \"flux\" missing in "
149 "input spectrum CPL table.\n"
150 "Cannot continue. Emergency stop.");
151 cpl_msg_debug(cpl_func, "-------------------------------------------");
152 cpl_msg_debug(cpl_func, "Incorrect input spectrum CPL table "
153 "structure:");
154 cpl_table_dump_structure(input_spectrum, NULL);
155 return cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
156 "%s", err_msg);
157 }
158 in_nrow=cpl_table_get_nrow(input_spectrum);
159 if (in_nrow == 0)
160 {
161 sprintf(err_msg, "Input spectrum has length zero.\n"
162 "Cannot continue. Emergency stop.");
163 cpl_msg_debug(cpl_func, "-------------------------------------------");
164 cpl_msg_debug(cpl_func, "Incorrect input spectrum CPL table "
165 "structure:");
166 cpl_table_dump_structure(input_spectrum, NULL);
167 return cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
168 "%s", err_msg);
169 }
170 if (cpl_error_get_code() != CPL_ERROR_NONE) {
171 return cpl_error_get_code();
172 }
173
174/*--------------------------------------------------------------------------*/
175/* ------------------------- Stats + Derivatives ---------------------------*/
176/*--------------------------------------------------------------------------*/
177
178 /* Creating histogram of spectral values */
179 //n_bins=ceil(in_nrow/bin_scale);
180 //hist1=cpl_table_new(n_bins); /* first derivative */
181
182 //cpl_table_new_column(hist1, "bins", CPL_TYPE_DOUBLE);
183 //cpl_table_new_column(hist1, "counts", CPL_TYPE_INT);
184
185 /* Creating derivatives table */
186 derivatives=cpl_table_new(in_nrow);
187 cpl_table_duplicate_column(derivatives, "lambda",
188 input_spectrum, "lambda");
189 cpl_table_duplicate_column(derivatives, "flux",
190 input_spectrum, "flux");
191 cpl_table_duplicate_column(derivatives, "weight",
192 input_spectrum, "weight");
193 cpl_table_new_column(derivatives, "first", CPL_TYPE_DOUBLE);
194 cpl_table_fill_column_window(derivatives, "first", 0, in_nrow, 0.);
195 //cpl_table_new_column(derivatives, "second", CPL_TYPE_DOUBLE);
196 if (cpl_error_get_code() != CPL_ERROR_NONE) {
197 return cpl_error_get_code();
198 }
199
200 /* Calculating first derivative: */
201 /* Defined simply as the flux difference between two consecutive pixels:*/
202 /* \Delta(flux)/pixel */
203 /* This simple approach turned out to be the most efficient method */
204 cpl_table_set_double(derivatives, "first", 0, 0.);
205 for(loop=0;loop<in_nrow-1;loop++)
206 {
207 if (cpl_table_get_double(input_spectrum, "weight", loop, NULL)
208 != 0. &&
209 cpl_table_get_double(input_spectrum, "weight", loop+1, NULL)
210 != 0. &&
211 cpl_table_get_double(input_spectrum, "lambda", loop, NULL)
212 <= SC_THERMIRLIM)
213 {
214 value=cpl_table_get_double(input_spectrum, "flux", loop+1, NULL)-
215 cpl_table_get_double(input_spectrum, "flux", loop, NULL);
216 cpl_table_set_double(derivatives, "first", loop+1, value);
217 }
218 }
219 if (cpl_error_get_code() != CPL_ERROR_NONE) {
220 return cpl_error_get_code();
221 }
222
223 /* Creating hists/stats on input spectrum/derivatives */
224 //sprintf(column_name, "%s", "first");
225 //sc_specdiss_create_hist(hist1, derivatives, column_name, n_bins);
226 //deriv1_lim=cpl_table_get_column_stdev(derivatives, "first");
227 cpl_table_unselect_all(derivatives);
228 cpl_table_or_selected_double(derivatives, "first", CPL_NOT_EQUAL_TO, 0.);
229 tmptable = cpl_table_extract_selected(derivatives);
230 cpl_table_select_all(derivatives);
231 cpl_table_abs_column(tmptable, "first");
232 sig=cpl_table_get_column_stdev(tmptable, "first");
233 mean=cpl_table_get_column_mean(tmptable, "first");
234
235 cpl_table_delete(tmptable);
236 if (sig / mean > maxrat) {
237 deriv1_lim = maxrat * mean;
238 } else {
239 deriv1_lim = sig;
240 }
241
242 /* Creating histograms */
243// plottags=cpl_parameterlist_new();
244// sc_setplottags_hist(plottags, "first derivative", "{/Symbol D}flux",
245// "counts per bin", parlist);
246// sc_plot_hist(hist1, plottags);
247// cpl_parameterlist_delete(plottags);
248
249 if (cpl_error_get_code() != CPL_ERROR_NONE) {
250 return cpl_error_get_code();
251 }
252
253/*--------------------------------------------------------------------------*/
254/* ------------------------------ Line search ------------------------------*/
255/*--------------------------------------------------------------------------*/
256
257 /* Searching for line peaks */
258 flux_low_lim=-lim_fac*deriv1_lim;
259 flux_up_lim=lim_fac*deriv1_lim;
260 sprintf(col1, "%s", "first");
261 sprintf(col2, "%s", "peaks");
262 sc_specdiss_find_localmaxima(derivatives, col2, col1, flux_low_lim,
263 flux_up_lim);
264 if (status != CPL_ERROR_NONE) {
265 //cpl_table_delete(hist1);
266 cpl_table_delete(derivatives);
267 return status;
268 }
269
270 /* Column "range" in CPL table derivatives denote the pixels regions */
271 /* belonging to lines (line pixels: value = 1) */
272 /* --> copying to input spectrum table as column "class" */
273 /* Column "peaks" in CPL table derivatives denote the line peaks */
274 /* (line peaks: value = 2) */
275
276 if (cpl_table_has_column(input_spectrum, "class") == 0) {
277 cpl_table_duplicate_column(input_spectrum, "class", derivatives,
278 "range");
279 } else {
280 cpl_table_copy_data_int(input_spectrum, "class",
281 cpl_table_get_data_int_const(derivatives, "range"));
282 }
283
284 for(loop=0;loop<in_nrow;loop++)
285 {
286 if(cpl_table_get_int(derivatives, "peaks", loop, NULL)==2)
287 {
288 cpl_table_set_int(input_spectrum, "class", loop, 2);
289 }
290 }
291
292 /* Checking for airglow lines */
293 /* Read + merging linelist and input spectrum */
294 status = sc_specdiss_merge_speclinelist(input_spectrum, groups, parlist);
295 if (status != CPL_ERROR_NONE) {
296 //cpl_table_delete(hist1);
297 cpl_table_delete(derivatives);
298 return status;
299 }
300
301 /* Identify airglow lines in spectrum */
302 status = sc_specdiss_identify_airglowlines(input_spectrum);
303 if (status != CPL_ERROR_NONE) {
304 //cpl_table_delete(hist1);
305 cpl_table_delete(derivatives);
306 return status;
307 }
308// /* Plotting spectrum with identified line peaks */
309// plottags=cpl_parameterlist_new();
310// sc_setplottags_single_spec_lines(plottags, "SINFONI",*/
311// "{/Symbol l} [micron]", "Flux", parlist);
312// sc_plot_single_spec_with_lines(input_spectrum, plottags);
313// cpl_parameterlist_delete(plottags);
314
315 /* Identifying isolated line peaks for FWHM determination with the
316 "class" column */
317 status = sc_specdiss_find_isolatedlines(linetab, input_spectrum, parlist);
318 if (status != CPL_ERROR_NONE) {
319 //cpl_table_delete(hist1);
320 cpl_table_delete(derivatives);
321 return status;
322 }
323
324 /* Plotting spectrum with identified line peaks */
325 pptr = cpl_parameterlist_find(parlist, "output_name");
326 sprintf(plot_title,"%s",cpl_parameter_get_string(pptr));
327
328// plottags=cpl_parameterlist_new();
329// sc_setplottags_single_spec_lines(plottags, plot_title,
330// "{/Symbol l} [micron]", "Flux", parlist);
331// sc_plot_single_spec_with_lines(input_spectrum, plottags);
332// cpl_parameterlist_delete(plottags);
333
334 /* Delete temporary column in input spectrum */
335 cpl_table_erase_column(input_spectrum, "linelist");
336
337 /* Cleaning */
338 //cpl_table_delete(hist1);
339 cpl_table_delete(derivatives);
340
341 return cpl_error_get_code();
342}
343
344/*--------------------------------------------------------------------------*/
345
346cpl_error_code sc_specdiss_find_isolatedlines(cpl_table *linetab,
347 cpl_table *input_spectrum,
348 cpl_parameterlist *parlist)
349{
422 int in_nrow=0; /* Number of columns / rows in input spec */
423 int loop=0; /* Loop variables */
424 int iteration=0; /* # of current iteration */
425 int peak_flag=0, peak_loc=0; /* Flag if line peak is found within line */
426 /* peak_loc=location of line peak in pix! */
427 int line_width=0, px_bef_peak=0, px_aft_peak=-1;
428 int line_dist=0, line_count=0;
429
430 int line_list_val=0;
431 double min_line_dist=0, /* Minimum distance between two lines */
432 min_line_dist_0=0;
433
434 double fwhm=0; /* initial FWHM */
435 int varfwhm=0; /* wavelength-dependent FWHM? */
436 double meanlam=0; /* mean wavelength [micron] */
437
438 int linetabsize=0;
439 int line_start_flag=0, line_px_start=0;
440 double peak_lambda=0., peak_flux=0.;
441 int found=0, index=0; /* found flag / index for writing isol lines*/
442 int dist_bef=0, dist_aft=0, symm=0; /* used for definition of criteria */
443
444 cpl_parameter *p; /* CPL parameters used to read parlist */
445
446/*--------------------------------------------------------------------------*/
447/* ------------------------------ INITIALISING -----------------------------*/
448/*--------------------------------------------------------------------------*/
449
450 /* Message of iteration */
451 p = cpl_parameterlist_find(parlist,"iteration");
452 iteration = cpl_parameter_get_int(p);
453 iteration++;
454 cpl_parameter_set_int(p,iteration);
455
456 /* Initial FWHM */
457 p = cpl_parameterlist_find(parlist, "fwhm");
458 fwhm = cpl_parameter_get_double(p);
459
460 /* Variable FWHM (linear change with wavelength)? */
461 p = cpl_parameterlist_find(parlist, "varfwhm");
462 varfwhm = cpl_parameter_get_int(p);
463
464 /* Get mean wavelength */
465 p = cpl_parameterlist_find(parlist, "meanlam");
466 meanlam = cpl_parameter_get_double(p);
467
468 /* Minimal line dist */
469 p = cpl_parameterlist_find(parlist,"min_line_dist_fac");
470 min_line_dist_0 = fwhm * cpl_parameter_get_double(p);
471
472 /* Creating histogram of spectral values */
473 in_nrow=cpl_table_get_nrow(input_spectrum);
474
475 line_dist=0;
476 for(loop=1;loop<in_nrow-1;loop++)
477 {
478 /* Initialising search over line region */
479 line_width=0;
480 px_bef_peak=0; /* # of pixels before/after the line peak */
481 px_aft_peak=-1; /* default -1, because px of peak is counted */
482 peak_loc=0;
483 line_start_flag=0;
484 line_list_val=0;
485 line_dist++;
486
487 /* Investigating region with a line (==value 1 in col "class") */
488 while(cpl_table_get_int(input_spectrum, "class", loop, NULL)>=1)
489 {
490 /* Determining first pixel of line */
491 if (line_start_flag==0)
492 {
493 line_start_flag=1;
494 line_px_start=loop;
495 }
496 line_width++;
497 found=1;
498 /* Checking location of line peak */
499 peak_flag=cpl_table_get_int(input_spectrum, "class", loop,
500 NULL);
501 if(peak_flag==2)
502 {
503 peak_loc=loop;
504 peak_lambda=cpl_table_get_double(input_spectrum, "lambda",
505 loop, NULL);
506 peak_flux=cpl_table_get_double(input_spectrum, "flux",
507 loop, NULL);
508 line_list_val=cpl_table_get_int(input_spectrum, "linelist",
509 loop, NULL);
510
511 }
512 /* Counting px before/after peak for symmetry criterion */
513 if(peak_loc==0) /* peak found if peak_loc > 0 */
514 {
515 px_bef_peak++;
516 }
517 else
518 {
519 px_aft_peak++;
520 }
521 loop++;
522 }
523
524 /* Correct peak_flux if FWHM is wavelength-dependent */
525 if (varfwhm == 1) {
526 peak_flux *= peak_lambda / meanlam;
527 }
528
529 /* Filling line property table */
530 if (found==1 && peak_flux > 0.)
531 {
532 linetabsize++;
533 cpl_table_set_size(linetab, linetabsize);
534 cpl_table_set_int(linetab, "width", linetabsize-1, line_width);
535 cpl_table_set_int(linetab, "peak_loc", linetabsize-1, peak_loc);
536 cpl_table_set_double(linetab, "peak_lam", linetabsize-1,
537 peak_lambda);
538 cpl_table_set_double(linetab, "peak_flux", linetabsize-1,
539 peak_flux);
540 cpl_table_set_int(linetab, "px_bef", linetabsize-1, px_bef_peak);
541 cpl_table_set_int(linetab, "px_aft", linetabsize-1, px_aft_peak);
542 cpl_table_set_int(linetab, "line_px_start", linetabsize-1,
543 line_px_start);
544 cpl_table_set_int(linetab, "line_px_end", linetabsize-1, loop-1);
545 if (iteration == 1) {
546 cpl_table_set_double(linetab, "fwhm", linetabsize-1, 0.);
547 } else {
548 cpl_table_set_double(linetab, "fwhm", linetabsize-1, fwhm);
549 }
550 cpl_table_set_int(linetab, "line_list", linetabsize-1,
551 line_list_val);
552 found=0;
553 line_dist=0;
554 }
555 }
556
557 /* Search for isolated lines:
558 *
559 * Applied criteria:
560 *
561 * (1) Distance to neighbouring line >=
562 * <min_line_dist_fac> * <measured FWHM> [pixel]
563 * For first line: = distance to spectrum start (in px),
564 * last line: = distance to spec end
565 * (2) # pixels left of peak = # pixels on the right hand side
566 * (+/- 1 px)
567 *
568 */
569
570 for(loop=0;loop<linetabsize;loop++)
571 {
572 peak_lambda=cpl_table_get_double(linetab, "peak_lam", loop, NULL);
573
574 /* First isolated line: Distance before=dist to spectrum start */
575 if(loop==0)
576 {
577 dist_bef=cpl_table_get_int(linetab, "line_px_start", loop, NULL);
578 cpl_table_set_int(linetab, "dist_bef", loop, dist_bef);
579 }
580 else
581 {
582 dist_bef=cpl_table_get_int(linetab, "line_px_start", loop, NULL)-
583 cpl_table_get_int(linetab, "line_px_end", loop-1, NULL);
584 cpl_table_set_int(linetab, "dist_bef", loop, dist_bef);
585 }
586
587 if(loop==linetabsize-1)
588 {
589 dist_aft=in_nrow-cpl_table_get_int(linetab, "line_px_end", loop,
590 NULL);
591 cpl_table_set_int(linetab, "dist_aft", loop, dist_aft);
592 }
593 else
594 {
595 dist_aft=cpl_table_get_int(linetab,"line_px_start", loop+1, NULL)-
596 cpl_table_get_int(linetab, "line_px_end", loop, NULL);
597 cpl_table_set_int(linetab, "dist_aft", loop, dist_aft);
598 }
599
600 symm=abs(cpl_table_get_int(linetab, "px_bef", loop, NULL)-
601 cpl_table_get_int(linetab, "px_aft", loop,
602 NULL));
603
604 /* Get min_line_dist for given wavelength depending on varfwhm flag */
605 if (varfwhm == 1) {
606 min_line_dist = min_line_dist_0 * peak_lambda / meanlam;
607 } else {
608 min_line_dist = min_line_dist_0;
609 }
610
611 /* Applying criteria */
612 if( ( dist_bef >= min_line_dist ) &&
613 ( dist_aft >= min_line_dist ) &&
614 ( symm <= 1 ) )
615 {
616 index = cpl_table_get_int(linetab, "peak_loc", loop, NULL);
617 cpl_table_set_int(input_spectrum, "class", index, 3);
618 cpl_table_set_int(linetab, "isol_flag", loop, 1);
619 line_count++;
620 }
621 else
622 {
623 cpl_table_set_int(linetab, "isol_flag", loop, 0);
624 }
625 }
626
627 return cpl_error_get_code();
628}
629
630/*--------------------------------------------------------------------------*/
631
632/* Internal routines */
633
634cpl_error_code sc_specdiss_find_localmaxima(cpl_table *input_table,
635 char *newcolname,
636 char *colname,
637 const double lower_threshold,
638 const double upper_threshold)
639{
664 int nrow=0; /* # of rows of input table */
665 int loop=0; /* Used for loopings */
666 double value=0.; /* Table value */
667 double loc_max_val=0.; /* temporary local maximum */
668 int found_flag=0;
669 int loc_max=0;
670 char col1[SC_MAXLEN], col2[SC_MAXLEN];
671
672 cpl_table *temporary; /* internal temporary table containing info not */
673 /* to be stored permanently */
674 int range_flag=0;
675
676/*--------------------------------------------------------------------------*/
677
678 nrow=cpl_table_get_nrow(input_table);
679 temporary=cpl_table_new(nrow);
680
681 /* all entries fullfilling threshold conditions */
682 cpl_table_new_column(temporary, "found", CPL_TYPE_INT);
683 cpl_table_new_column(temporary, "range", CPL_TYPE_INT);
684 cpl_table_new_column(temporary, "range_red", CPL_TYPE_INT);
685 cpl_table_new_column(temporary, "range_blue", CPL_TYPE_INT);
686
687 /* Creating new col to be added to input table & initially filled with 0*/
688 cpl_table_new_column(input_table, newcolname, CPL_TYPE_INT);
689 cpl_table_fill_column_window_int(input_table, newcolname, 0, nrow, 0);
690
691 /* Searching for values */
692 cpl_table_duplicate_column(temporary, "target_col", input_table, colname);
693/* sc_specdiss_find_valuesabove(temporary, "above_thresh",
694 "target_col", upper_threshold);*/
695 sprintf(col1, "%s", "target_col");
696 sprintf(col2, "%s", "found");
697
698 /* Searching for pixels above/below the thresh (= +1/-1, 0 otherwise) */
699 sc_specdiss_find_valuesoutside(temporary, col2, col1, lower_threshold,
700 upper_threshold);
701
702 /* Searching for line range (=range of pixels belonging to one line) */
703 /* Line range is limited by a change of the values in column "found" */
704 /* from -1 --> 0 or 1 and 1 --> 0 or -1 */
705
706 /* Applying sweep technique to avoid false line range identifications */
707 /* caused by asymmetric lines. Sweep technique = search direction first */
708 /* towards red end, towards end afterwards blue */
709
710 /* Run towards red end of spectrum */
711 for (loop=0;loop<nrow-1;loop++)
712 {
713 found_flag=cpl_table_get_int(temporary, "found", loop, NULL);
714 if ( found_flag == 1 )
715 {
716 range_flag=1;
717 }
718 if ((cpl_table_get_int(temporary, "found", loop, NULL) ==-1 &&
719 cpl_table_get_int(temporary, "found", loop+1, NULL) == 0) ||
720 (cpl_table_get_int(temporary, "found", loop, NULL) ==-1 &&
721 cpl_table_get_int(temporary, "found", loop+1, NULL) == 1))
722 {
723 range_flag=0;
724 }
725 cpl_table_set_int(temporary, "range_red", loop, range_flag);
726 }
727
728 /* Run towards blue end of spectrum */
729 for ( loop = nrow-2; loop > 0; loop-- )
730 {
731 found_flag = cpl_table_get_int(temporary, "found", loop, NULL);
732 if ( found_flag == -1 )
733 {
734 range_flag = 1;
735 }
736 if (( cpl_table_get_int(temporary, "found", loop+1, NULL) == 1 &&
737 cpl_table_get_int(temporary, "found", loop, NULL) == 0 ) ||
738 ( cpl_table_get_int(temporary, "found", loop+1, NULL) == 1 &&
739 cpl_table_get_int(temporary, "found", loop, NULL) == -1 ))
740 {
741 range_flag = 0;
742 }
743 cpl_table_set_int(temporary, "range_blue", loop, range_flag);
744 }
745 cpl_table_set_int(temporary, "range_blue", loop, 0); /* 1st row */
746
747 /* Merging the two line searches */
748 cpl_table_new_column(input_table, "range", CPL_TYPE_INT);
749
750 for (loop=0;loop<nrow;loop++)
751 {
752 if (( cpl_table_get_int(temporary, "range_red", loop, NULL) == 1 ) &&
753 ( cpl_table_get_int(temporary, "range_blue", loop, NULL) == 1 ) &&
754 ( cpl_table_get_double(input_table, "lambda", loop, NULL)
755 <= SC_THERMIRLIM ))
756 {
757 cpl_table_set_int(input_table, "range", loop, 1);
758 }
759 else
760 {
761 cpl_table_set_int(input_table, "range", loop, 0);
762 }
763 }
764
765 /* Searching for local maxima in line ranges */
766 loop=0;
767 loc_max_val=0;
768
769 while (loop < nrow)
770 {
771 range_flag = cpl_table_get_int(input_table, "range", loop, NULL);
772 if (range_flag == 0)
773 {
774 loc_max=0;
775 loc_max_val=0;
776 }
777 else
778 {
779 while ( cpl_table_get_int(input_table, "range", loop, NULL) == 1)
780 {
781 value=cpl_table_get_double(input_table, "flux", loop, NULL);
782 if (value > loc_max_val)
783 {
784 loc_max=loop;
785 loc_max_val=value;
786 }
787 loop++;
788 }
789 loop--;
790 cpl_table_set_int(input_table, newcolname, loc_max, 2);
791 }
792 loop++;
793 }
794
795 cpl_table_delete(temporary);
796 return cpl_error_get_code();
797
798}
799
800/*--------------------------------------------------------------------------*/
801
802cpl_error_code sc_specdiss_find_valuesoutside(cpl_table *input_table,
803 char *newcolname,
804 char *colname,
805 const double lowlim,
806 const double uplim)
807{
839 int nrow=0; /* # of rows of input table */
840 int loop=0; /* Used for loopings */
841 double value=0.; /* Table value */
842
843 nrow=cpl_table_get_nrow(input_table);
844
845 for(loop=0;loop<nrow;loop++)
846 {
847 cpl_table_set_int(input_table, newcolname, loop, 0); /* default=0 */
848 value=cpl_table_get(input_table, colname, loop, NULL);
849 if( (value <= lowlim) )
850 {
851 cpl_table_set_int(input_table, newcolname, loop, -1);
852 }
853 if ((value >= uplim) )
854 {
855 cpl_table_set_int(input_table, newcolname, loop, 1);
856 }
857
858 }
859
860 return CPL_ERROR_NONE;
861
862}
863
864/*--------------------------------------------------------------------------*/
865
866int sc_specdiss_count_lines(const cpl_table *intable)
867{
879 int n_lines=0;
880 int loop=0;
881 int in_nrow=0;
882
883 in_nrow=cpl_table_get_nrow(intable);
884
885 for(loop=0;loop<in_nrow-1;loop++)
886 {
887 n_lines=n_lines+cpl_table_get_int(intable, "peaks", loop, NULL);
888 }
889 return n_lines;
890}
891
892/*--------------------------------------------------------------------------*/
893
894cpl_error_code sc_specdiss_create_hist(cpl_table *hist,
895 const cpl_table *input_table,
896 char *col_name, const int n_bins)
897{
918 int tab_nrow=0; /* length of table */
919 double col_min=0., col_max=0.; /* min/max value of column */
920 int loop=0, loop2=0; /* for loopings */
921 int count=0;
922 double bin_loop=0.;
923 double data_point=0.;
924 double lowlim=0., uplim=0.;
925
926 cpl_table *bin_limits;
927
928/*--------------------------------------------------------------------------*/
929
930 tab_nrow=cpl_table_get_nrow(input_table);
931
932 col_min=cpl_table_get_column_min(input_table, col_name);
933 col_max=cpl_table_get_column_max(input_table, col_name);
934
935 /* Creating equal n_bins histogram bins between col_min and col_max */
936 bin_limits=cpl_table_new(n_bins);
937 cpl_table_new_column(bin_limits, "lowlim", CPL_TYPE_DOUBLE);
938 cpl_table_new_column(bin_limits, "uplim", CPL_TYPE_DOUBLE);
939
940 for (loop=0;loop<=n_bins-1;loop++)
941 {
942 bin_loop=col_min+loop*(col_max-col_min)/n_bins;
943 cpl_table_set_double(bin_limits, "lowlim", loop, bin_loop);
944 bin_loop=col_min+(loop+1)*(col_max-col_min)/n_bins;
945 cpl_table_set_double(bin_limits, "uplim", loop, bin_loop);
946 }
947 /* Creating hist */
948 for (loop=0;loop<=n_bins-1;loop++)
949 {
950 lowlim=cpl_table_get_double(bin_limits, "lowlim", loop, NULL);
951 uplim=cpl_table_get_double(bin_limits, "uplim", loop, NULL);
952
953 /* Searching for values in bin */
954 count=0;
955 for (loop2=0;loop2<=tab_nrow-1;loop2++)
956 {
957 data_point=cpl_table_get(input_table, col_name, loop2, NULL);
958 if( ( data_point >= lowlim ) && ( data_point < uplim ) )
959 {
960 count++;
961 }
962 }
963 cpl_table_set(hist, "bins", loop, lowlim+(uplim-lowlim)/2);
964 cpl_table_set(hist, "counts", loop, count);
965 }
966
967 /* Cleaning */
968 cpl_table_delete(bin_limits);
969 return CPL_ERROR_NONE;
970}
971
972/*--------------------------------------------------------------------------*/
973
974void sc_specdiss_init_linetab(cpl_table *linetab)
975{
1006 cpl_table_new_column(linetab, "width", CPL_TYPE_INT);
1007 cpl_table_new_column(linetab, "peak_loc", CPL_TYPE_INT);
1008 cpl_table_new_column(linetab, "peak_lam", CPL_TYPE_DOUBLE);
1009 cpl_table_new_column(linetab, "peak_flux", CPL_TYPE_DOUBLE);
1010 cpl_table_new_column(linetab, "px_bef", CPL_TYPE_INT);
1011 cpl_table_new_column(linetab, "px_aft", CPL_TYPE_INT);
1012 cpl_table_new_column(linetab, "line_px_start", CPL_TYPE_INT);
1013 cpl_table_new_column(linetab, "line_px_end", CPL_TYPE_INT);
1014 cpl_table_new_column(linetab, "dist_bef", CPL_TYPE_INT);
1015 cpl_table_new_column(linetab, "dist_aft", CPL_TYPE_INT);
1016 cpl_table_new_column(linetab, "isol_flag", CPL_TYPE_INT);
1017 cpl_table_new_column(linetab, "fwhm", CPL_TYPE_DOUBLE);
1018 cpl_table_new_column(linetab, "line_list", CPL_TYPE_INT);
1019}
1020
1021/*--------------------------------------------------------------------------*/
1022
1023cpl_error_code sc_specdiss_merge_speclinelist(cpl_table *spec,
1024 const cpl_table *groups,
1025 const cpl_parameterlist *parlist)
1026{
1047 double fwhm = 0, fwhm_crit = 0, fwhm_crit_fac = 0.1, /* FWHM criterion */
1048 fwhm_crit_0 = 0.,
1049 line_flux_sum = 0.,
1050 pxscl = 0., /* pixel scale [micron/px] */
1051 meanlam = 0.; /* mean lambda [micron] */
1052
1053 double lam_spec0 = 0., lam_spec1 = 0., /* neighbouring criterion */
1054 lam_spec2 = 0., lam_list = 0., /* wavelength in line list */
1055 lam_lim1 = 0., lam_lim2 = 0.;
1056
1057
1058 int listlength=0; /* length of original line list */
1059 int line_arr_size=0; /* length of selected line list */
1060
1061 int i=0, j=0, k=0, /* loop variables */
1062 n_loops = 10; /* maximum number of loops/iterations */
1063
1064 int peakcount1 = 2; /* # lines of linelist in pixel */
1065
1066 int varfwhm = 0; /* wavelength-dependent FWHM? */
1067
1068 cpl_array *line_arr; /* line array */
1069
1070 const cpl_parameter *p; /* CPL parameter */
1071
1072/*--------------------------------------------------------------------------*/
1073
1074 /* Initialising */
1075
1076 /* Initialising line array */
1077 listlength = cpl_table_get_nrow(groups);
1078 n_loops = cpl_table_get_nrow(spec);
1079 if (cpl_table_has_column(spec, "linelist") != 1)
1080 {
1081 cpl_table_new_column(spec, "linelist", CPL_TYPE_INT);
1082 }
1083
1084 /* Initialising new column */
1085 for (i = 0;i < n_loops; i++)
1086 {
1087 cpl_table_set_int(spec,"linelist", i, 0);
1088 }
1089
1090 /* Calculating FWHM */
1091 p = cpl_parameterlist_find_const(parlist, "fwhm");
1092 fwhm = cpl_parameter_get_double(p);
1093
1094 pxscl = (cpl_table_get_double(spec,"lambda",n_loops-1,NULL) -
1095 cpl_table_get_double(spec,"lambda",0,NULL) ) / n_loops;
1096 fwhm_crit_0 = fwhm_crit_fac * fwhm * pxscl;
1097
1098 /* Variable FWHM (linear change with wavelength)? */
1099 p = cpl_parameterlist_find_const(parlist, "varfwhm");
1100 varfwhm = cpl_parameter_get_int(p);
1101
1102 /* Get mean wavelength */
1103 p = cpl_parameterlist_find_const(parlist, "meanlam");
1104 meanlam = cpl_parameter_get_double(p);
1105
1106 /* Search & count lines */
1107 peakcount1=0;
1108 /* Counting linelist lines per spectrum pixel */
1109 for (i = 1; i < n_loops-1; i++) {
1110 const double * pdlam;
1111 cpl_size nbad;
1112
1113 /* Reset peak counters */
1114 peakcount1=0;
1115
1116 /* Counting lines from line list belonging to pixel #i */
1117 /* Get wavelength borders of pixels #i-1, #i, #i+1 */
1118 lam_spec0 = cpl_table_get_double(spec, "lambda", i-1, NULL);
1119 lam_spec1 = cpl_table_get_double(spec, "lambda", i, NULL);
1120 lam_spec2 = cpl_table_get_double(spec, "lambda", i+1, NULL);
1121
1122 /* Half distance to neighboured pixels */
1123 lam_lim1 = lam_spec1 - (lam_spec1 - lam_spec0) / 2;
1124 lam_lim2 = lam_spec1 + (lam_spec2 - lam_spec1) / 2;
1125
1126 /* Get fwhm_crit for given wavelength depending on varfwhm flag */
1127 if (varfwhm == 1) {
1128 fwhm_crit = fwhm_crit_0 * lam_spec1 / meanlam;
1129 } else {
1130 fwhm_crit = fwhm_crit_0;
1131 }
1132
1133 /* Creating list of lines within current pixel */
1134 line_arr_size = 0;
1135 line_arr = cpl_array_new(line_arr_size, CPL_TYPE_DOUBLE);
1136 line_flux_sum = 0;
1137 /* TODO probably better to use selections here */
1138 pdlam = cpl_table_get_data_double_const(groups, "lambda");
1139 nbad = cpl_table_count_invalid(groups, "lambda");
1140
1141 /* Looping over linelist */
1142 for (j = 0; j < listlength; j++)
1143 {
1144 /* get current wavelength in line list */
1145 if (nbad == 0) {
1146 lam_list = pdlam[j];
1147 }
1148 else {
1149 lam_list = cpl_table_get_double(groups, "lambda", j, NULL);
1150 }
1151
1152 if (lam_list < lam_spec0)
1153 {
1154 /* skip lines at low end of wavelength range in spec */
1155 continue;
1156 }
1157
1158 /* Counting the lines from the line list belonging to pixel #i */
1159 /* Selection criterion: distance closer to pixel #i than */
1160 /* to pixels #i-1 or #i+1 */
1161 if ( ( lam_list > lam_lim1 ) && ( lam_list <= lam_lim2 ) )
1162 {
1163 line_arr_size++;
1164 cpl_array_set_size(line_arr,line_arr_size);
1165 cpl_array_set_double(line_arr,line_arr_size-1,lam_list);
1166 line_flux_sum += cpl_table_get_double(groups,
1167 "flux",j,NULL);
1168
1169 peakcount1++;
1170 }
1171 }
1172
1173 /* Checking significance of lines */
1174 /* Check distance between lines: If 2 or more lines are too closely */
1175 /* neighboured (delta_lambda < fwhm_crit) --> assumed as being a */
1176 /* single one */
1177 if (line_arr_size > 1)
1178 {
1179 k=0;
1180 while (k < line_arr_size-1)
1181 {
1182 if ( cpl_array_get_double(line_arr,k+1,NULL) -
1183 cpl_array_get_double(line_arr,k,NULL) < fwhm_crit )
1184 {
1185 peakcount1 --;
1186 }
1187 k++;
1188 }
1189
1190 }
1191
1192 cpl_array_delete(line_arr);
1193 cpl_table_set_int(spec, "linelist", i, peakcount1);
1194 }
1195
1196 return cpl_error_get_code();
1197
1198}
1199
1200/*--------------------------------------------------------------------------*/
1201
1202cpl_error_code sc_specdiss_identify_airglowlines(cpl_table *spec)
1203{
1214 int nrow=0;
1215 int loop=0;
1216
1217 nrow = cpl_table_get_nrow(spec);
1218
1219 for (loop = 0; loop < nrow; loop++)
1220 {
1221 if( ( cpl_table_get_int(spec, "class", loop, NULL) >= 2 ) &&
1222 ( cpl_table_get_int(spec, "linelist", loop, NULL) == 5 ) )
1223 {
1224 cpl_table_set_int(spec, "linelist", loop, 3);
1225 }
1226 }
1227
1228 return cpl_error_get_code();
1229
1230}
1231
1232/*--------------------------------------------------------------------------*/
1233
1234cpl_error_code sc_specdiss_find_valuesabove(cpl_table *input_table,
1235 char *newcolname,
1236 char *colname,
1237 const double lowlim)
1238{
1263 int nrow=0; /* # of rows of input table */
1264 int loop=0; /* Used for loopings */
1265 double value=0.; /* Table value */
1266
1267 nrow=cpl_table_get_nrow(input_table);
1268 cpl_table_new_column(input_table, newcolname, CPL_TYPE_INT);
1269
1270 for(loop=0;loop<nrow;loop++)
1271 {
1272 cpl_table_set_int(input_table, newcolname, loop, 0); /* default=0 */
1273 value=cpl_table_get(input_table, colname, loop, NULL);
1274 if( (value>=lowlim) )
1275 {
1276 cpl_table_set_int(input_table, newcolname, loop, 1);
1277 }
1278 }
1279
1280 return CPL_ERROR_NONE;
1281
1282}
1283
#define SC_MAXLEN
Definition: sc_basic.h:94
#define SC_THERMIRLIM
Definition: sc_basic.h:102
cpl_error_code sc_specdiss_find_valuesoutside(cpl_table *input_table, char *newcolname, char *colname, const double lowlim, const double uplim)
Definition: sc_specdiss.c:802
void sc_specdiss_init_linetab(cpl_table *linetab)
Definition: sc_specdiss.c:974
int sc_specdiss_count_lines(const cpl_table *intable)
Definition: sc_specdiss.c:866
cpl_error_code sc_specdiss_find_emissionlines(cpl_table *input_spectrum, cpl_table *linetab, cpl_parameterlist *parlist, const cpl_table *groups)
Definition: sc_specdiss.c:59
cpl_error_code sc_specdiss_identify_airglowlines(cpl_table *spec)
Definition: sc_specdiss.c:1202
cpl_error_code sc_specdiss_merge_speclinelist(cpl_table *spec, const cpl_table *groups, const cpl_parameterlist *parlist)
Definition: sc_specdiss.c:1023
cpl_error_code sc_specdiss_create_hist(cpl_table *hist, const cpl_table *input_table, char *col_name, const int n_bins)
Definition: sc_specdiss.c:894
cpl_error_code sc_specdiss_find_valuesabove(cpl_table *input_table, char *newcolname, char *colname, const double lowlim)
Definition: sc_specdiss.c:1234
cpl_error_code sc_specdiss_find_isolatedlines(cpl_table *linetab, cpl_table *input_spectrum, cpl_parameterlist *parlist)
Definition: sc_specdiss.c:346
cpl_error_code sc_specdiss_find_localmaxima(cpl_table *input_table, char *newcolname, char *colname, const double lower_threshold, const double upper_threshold)
Definition: sc_specdiss.c:634
Header for spectral analysis library.