X-shooter Pipeline Reference Manual 3.8.15
test-xsh_spectrum_detect_lines.c
Go to the documentation of this file.
1/* *
2 * This file is part of the ESO X-shooter Pipeline *
3 * Copyright (C) 2006 European Southern Observatory *
4 * *
5 * This library 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 program 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 program; if not, write to the Free Software *
17 * Foundation, 51 Franklin St, Fifth Floor, Boston, MA 02111-1307 USA *
18 * */
19
20/*
21 * $Author: amodigli $
22 * $Date: 2011-12-02 14:15:28 $
23 * $Revision: 1.13 $
24 * $Name: not supported by cvs2svn $
25 */
26
27#ifdef HAVE_CONFIG_H
28#include <config.h>
29#endif
30
31/*---------------------------------------------------------------------------*/
38/*---------------------------------------------------------------------------*/
41/*---------------------------------------------------------------------------
42 Includes
43 ----------------------------------------------------------------------------*/
44#include <math.h>
45#include <xsh_error.h>
46#include <xsh_msg.h>
47#include <xsh_eqwidth_lib.h>
48#include <xsh_utils_table.h>
49
50/*---------------------------------------------------------------------------
51 Typedefs
52 ---------------------------------------------------------------------------*/
53#define INV_DOUBLE -9999.0
54#define LIGHTSPEED1000 299792458 // Speed of light [m/s]
55
56/*---------------------------------------------------------------------------
57 Functions prototypes
58 ---------------------------------------------------------------------------*/
59
60/*---------------------------------------------------------------------------
61 Implementation
62 ---------------------------------------------------------------------------*/
63static cpl_table*
64xsh_table_unique_wave(cpl_table* tab,const char* wname, const char* fname) {
65
66 cpl_table* xtab=NULL;
67 xtab=cpl_table_duplicate(tab);
68 int norg=cpl_table_get_nrow(xtab);
69 xsh_msg("N orig %d",norg);
70 double* pwave=cpl_table_get_data_double(xtab,wname);
71 double* pflux=cpl_table_get_data_double(xtab,fname);
72 int ninv=1;
73 int nrow=0;
74 int nbad=0;
75 int j=0;
76 for(j=0;ninv>0,j<4;j++) {
77 nrow=cpl_table_get_nrow(xtab);
78 for(int i=1;i < nrow; i++) {
79 if (!(pwave[i-1] < pwave[i])) {
80 if(pflux[i-1]<= pflux[i]) {
81 cpl_table_set_invalid(xtab,wname,i-1);
82 } else {
83 cpl_table_set_invalid(xtab,wname,i);
84 }
85 nbad++;
86 }
87 }
88
89 ninv=cpl_table_count_invalid(xtab,wname);
90 xsh_msg("iter=%d nrow=%d nbad=%d ninv=%d",j,nrow,nbad,ninv);
91 if(ninv>0) {
92 cpl_table_erase_invalid(xtab);
93 } else {
94 break;
95 }
96 }
97
98
99 int nnew=cpl_table_get_nrow(xtab);
100 xsh_msg("niter=%d",j);
101 xsh_msg("N orig %d flagged %d expected %d new %d",norg,nbad,norg-nbad,nnew);
102 return xtab;
103}
104
105/* test line detection on a spectrum provided in an input FITS table */
106cpl_error_code
107test_spectrum_detect_lines(int argc, char* argv[])
108{
109
110 char fname[256];
111
112 double hwidth=0.3;//2.0;
113 double overlap=0.5;
114 int niter=10;
115 double thres=0.01;
116 double resol=0.01;
117 int smwidth=20;
118 double wdiff=0.02;
119 cpl_table* spec_tab=NULL;
120 double zem = 0.0;
121 double vel_step=300.;
122 xsh_msg_warning("argc=%d",argc);
123
124 if((size_t)argc==1){
125 return cpl_error_get_code();
126 }
127 if((size_t)argc==2){
128 sprintf(fname,argv[1]);
129 xsh_msg_warning("filename=%s",fname);
130 }
131 else if((size_t)argc>2){
132 xsh_msg_warning("reading params");
133 hwidth=atof(argv[2]);
134 overlap=atof(argv[3]);
135 niter=atoi(argv[4]);
136 thres=atof(argv[5]);
137 resol=atof(argv[6]);
138 smwidth=atoi(argv[7]);
139 wdiff=atof(argv[8]);
140
141 }
142 xsh_msg_warning("hwidth=%f overlap=%f niter=%d thres=%f resol=%f smdiff=%d wdiff=%f",hwidth,overlap,niter,thres,resol,smwidth,wdiff);
143
144 /* load input data */
145 cpl_table* org_tab=cpl_table_load(fname, 1, 0);
146 /* remove wave duplicates */
147 spec_tab=xsh_table_unique_wave(org_tab,"WAVE","FLUX");
148 //cpl_table_dump_structure(spec_tab,stdout);
149 /* convert emission in absorbtion spectrum to work with ARES routines */
150 double max=cpl_table_get_column_max(spec_tab,"FLUX");
151 check(cpl_table_add_scalar(spec_tab,"FLUX",max));
152 check(cpl_table_duplicate_column(spec_tab,"INV_FLUX",spec_tab,"FLUX"));
153 int nrow=cpl_table_get_nrow(spec_tab);
154 check(cpl_table_fill_column_window_double(spec_tab,"INV_FLUX",0,nrow,1.));
155 check(cpl_table_divide_columns(spec_tab,"INV_FLUX","FLUX"));
156 cpl_table_erase_column(spec_tab,"FLUX");
157 //cpl_table_dump_structure(spec_tab,stderr);
158 cpl_table_name_column(spec_tab,"INV_FLUX","FLUX");
159
160 cpl_table_name_column(spec_tab,"WAVE","WAVEL");
161 check(xsh_sort_table_1(spec_tab,"WAVEL",CPL_FALSE));
162 cpl_table_save(spec_tab,NULL,NULL,"spec_tab.fits",CPL_IO_DEFAULT);
163
164
165 // Chunk definition
166 double wavel_min = cpl_table_get_double(spec_tab, "WAVEL", 0, NULL);
167 double wavel_max = cpl_table_get_double(spec_tab, "WAVEL",
168 cpl_table_get_nrow(spec_tab) - 1, NULL);
169 double pos_start = wavel_min + hwidth;
170 double pos_end = wavel_max - hwidth;
171 if (wavel_min < cpl_table_get_column_min(spec_tab, "WAVEL")) {
172 pos_start = cpl_table_get_column_min(spec_tab, "WAVEL") + hwidth;
173 }
174 if (wavel_max > cpl_table_get_column_max(spec_tab, "WAVEL")) {
175 pos_end = cpl_table_get_column_max(spec_tab, "WAVEL") - hwidth;
176 }
177 xsh_msg_warning("wmin=%g wmax=%g",wavel_min,wavel_max);
178 xsh_msg_warning("pos_start=%g pos_end=%g",pos_start,pos_end);
179
180 cpl_table *spec_proc = NULL;
181 cpl_table *line = NULL;
182 cpl_size line_num = 0;
183 cpl_size line_size = 0;
184 if (pos_end < pos_start) {
185 cpl_msg_warning(cpl_func, "Chunks size is too large w.r.t. the "
186 "wavelength range of the spectrum! It should be lower "
187 "than %f nm.", 0.5 * (wavel_max - wavel_min));
188 }
189 else {
190
191 double pos_step = hwidth * 2.0 * (1.0 - overlap);
192 int chunk_num = (int)ceil((pos_end - pos_start)/pos_step);
193 cpl_table *line_chunk = NULL;
194 cpl_table *line_temp = cpl_table_new(10000);
195 int i = 0;
196 int flag = 0;
197 cpl_table_new_column(line_temp, "WAVEL", CPL_TYPE_DOUBLE);
198 cpl_table_new_column(line_temp, "PEAK", CPL_TYPE_DOUBLE);
199 cpl_table_new_column(line_temp, "CONT", CPL_TYPE_DOUBLE);
200 cpl_table_new_column(line_temp, "CONTERR", CPL_TYPE_DOUBLE);
201 cpl_table_fill_column_window_double(line_temp, "CONT", 0, 10000,
202 INV_DOUBLE);
203 cpl_table_fill_column_window_double(line_temp, "CONTERR", 0, 10000,
204 INV_DOUBLE);
205
206 // Undersample QSO spectrum
207
208 if (zem > 0.0) {
209 cpl_msg_info(cpl_func, "*** Undersampling spectrum ***");
210 double vel_step_test = 2.0 * hwidth / (wavel_min *
211 (exp(vel_step / LIGHTSPEED1000) - 1.0));
212 if (vel_step_test < 3.0) {
213 cpl_msg_warning(cpl_func, "Undersampling of the spectrum is "
214 "incompatible with the adopted chunk size!");
215 vel_step = (double)floor(LIGHTSPEED1000 *
216 log(2.0 * hwidth / (3.0 * wavel_min) + 1.0));
217 cpl_msg_info(cpl_func, "Parameter vel-step is too low: redefined "
218 "to %6.1f.", vel_step);
219 }
220 /* AMO: comment out this not to include extra function
221 cpl_ensure_code(espda_spec_rebin(spec, kappa, wavel_min, wavel_max,
222 vel_step, &spec_proc) == CPL_ERROR_NONE,
223 cpl_error_get_code());
224 */
225 }
226 else {
227 spec_proc = cpl_table_duplicate(spec_tab);
228 }
229
230 // Loop over chunks
231 cpl_msg_info(cpl_func, "*** Creating a line list ***");
232 cpl_msg_info(cpl_func, "Processing chunks...");
233 for (double pos = pos_start; pos < pos_end; pos += pos_step) {
234 i++;
235
236 cpl_table *chunk = NULL;
237 // Stellar spectra
238 if (zem == 0.0) {
239
240 // Chunk selection
241 xsh_msg_warning("pos=%g hwidth=%g",pos,hwidth);
242 cpl_ensure_code(select_local_spec(spec_tab, hwidth, pos,
243 &chunk) == CPL_ERROR_NONE,
244 cpl_error_get_code());
245 cpl_table_unselect_all(chunk);
246 //double cont_rejt = 0.995;
247 double cont_rejt = 0.95;
248 // Chunk sorting
249 cpl_propertylist *chunk_col = cpl_propertylist_new();
250 cpl_propertylist_append_bool(chunk_col, "WAVEL", FALSE);
251 cpl_table_sort(chunk, chunk_col);
252 cpl_propertylist_delete(chunk_col);
253 // Local continuum fitting
254 cpl_ensure_code(esp_fit_lcont(chunk, cont_rejt, niter)
255 == CPL_ERROR_NONE, cpl_error_get_code());
256 // Line detection
257 cpl_ensure_code(esp_det_line(chunk, thres, resol, smwidth,
258 &line_chunk) == CPL_ERROR_NONE,
259 cpl_error_get_code());
260 // Line list creation
261 cpl_table_unselect_all(line_chunk);
262 cpl_table_or_selected_double (line_chunk, "WAVEL",
263 CPL_GREATER_THAN, pos - hwidth * (1.0 - overlap));
264 cpl_table_and_selected_double(line_chunk, "WAVEL",
265 CPL_NOT_GREATER_THAN, pos + hwidth * (1.0 - overlap));
266 if (cpl_table_count_selected(line_chunk) > 0) {
267 cpl_table *line_part = cpl_table_extract_selected(line_chunk);
268 if (flag == 0 && cpl_table_get_nrow(line_part) != 0) {
269 line = cpl_table_duplicate(line_part);
270 flag = 1;
271 }
272 else {
273 cpl_table_insert(line, line_part, line_num);
274 }
275 line_num += cpl_table_get_nrow(line_part);
276 cpl_table_delete(line_part);
277 }
278
279 cpl_table_delete(line_chunk);
280 }
281
282 // QSO spectra
283 else {
284 // Chunk selection
285 cpl_ensure_code(select_local_spec(spec_proc, hwidth, pos,
286 &chunk) == CPL_ERROR_NONE,
287 cpl_error_get_code());
288 cpl_table_unselect_all(chunk);
289 //printf("%f %i\n", pos, cpl_table_get_nrow(chunk));
290 //cpl_table_dump(chunk, 0, 100, NULL);
291 // Alternative way
292 /*
293 double wavel_lo = pos - hwidth;
294 double wavel_up = pos + hwidth;
295
296 // Last chunk is expanded to reach the end of the spectrum
297 if (pos + pos_step > pos_end) {
298 wavel_up = pos_end + hwidth;
299 }
300 cpl_table_unselect_all(spec_tab);
301 cpl_table_or_selected_double (spec_tab, "WAVEL", CPL_NOT_LESS_THAN,
302 wavel_lo);
303 cpl_table_and_selected_double(spec_tab, "WAVEL", CPL_NOT_GREATER_THAN,
304 wavel_up);
305 chunk = cpl_table_extract_selected(spec_tab);
306 */
307
308 cpl_propertylist *chunk_col = cpl_propertylist_new();
309 cpl_propertylist_append_bool(chunk_col, "FLUX", FALSE);
310 cpl_table_sort(chunk, chunk_col);
311 cpl_propertylist_delete(chunk_col);
312 int med_row = floor((int)cpl_table_get_nrow(chunk) / 2.0);
313 double flux_ave = cpl_table_get_column_mean(chunk, "FLUX");
314 double flux_med =
315 (cpl_table_get_double(chunk, "FLUX", med_row, NULL) +
316 cpl_table_get_double(chunk, "FLUX", med_row + 1, NULL))
317 / 2.0;
318 if (flux_med - flux_ave > thres) {
319 cpl_table_set_double(line_temp, "WAVEL", line_num,
320 cpl_table_get_double(chunk, "WAVEL", 0, NULL));
321 cpl_table_set_double(line_temp, "PEAK", line_num,
322 cpl_table_get_double(chunk, "FLUX", 0, NULL));
323 line_num++;
324 }
325 }
326 cpl_table_delete(chunk);
327
328 // Loop printout
329 int chunk_ord = (int)ceil(log10(chunk_num));
330 // if (i % (int)round(chunk_num/pow(10.0, (double)chunk_ord - 1.0)) == 0) {
331 if ((i + 1) % (int)ceil(chunk_num/10.0) == 0 && i < chunk_num) {
332 cpl_msg_warning(cpl_func, " %-*i chunks processed; %-*i remaining... ",
333 chunk_ord, i, chunk_ord, chunk_num - i);
334 }
335 if (i == chunk_num) {
336 cpl_msg_warning(cpl_func, "%i lines found...", (int)line_num);
337 }
338 }
339 line = cpl_table_extract(line_temp, 0, line_num);
340 cpl_table_save(line,NULL,NULL,"line.fits",CPL_IO_DEFAULT);
341
342 cpl_table_delete(line_temp);
343 }
344cleanup:
345 xsh_free_table(&spec_tab);
346 xsh_free_table(&org_tab);
347 return cpl_error_get_code();
348}
349
350
351/*----------------------------------------------------------------------------*/
355/*----------------------------------------------------------------------------*/
356
357int main(int argc, char* argv[])
358{
359
360 cpl_test_init(PACKAGE_BUGREPORT, CPL_MSG_WARNING);
362
363 return cpl_test_end(0);
364}
365
int main()
Unit test of xsh_bspline_interpol.
#define check(COMMAND)
Definition: xsh_error.h:71
#define xsh_msg_warning(...)
Print an warning message.
Definition: xsh_msg.h:88
#define xsh_msg(...)
Print a message on info level.
Definition: xsh_msg.h:121
static cpl_table * xsh_table_unique_wave(cpl_table *tab, const char *wname, const char *fname)
#define LIGHTSPEED1000
cpl_error_code test_spectrum_detect_lines(int argc, char *argv[])
#define INV_DOUBLE
cpl_error_code xsh_sort_table_1(cpl_table *t, const char *column, cpl_boolean reverse)
Sort a table by one column.
void xsh_free_table(cpl_table **t)
Deallocate a table and set the pointer to NULL.
Definition: xsh_utils.c:2133
int niter
Definition: xsh_detmon_lg.c:82
cpl_error_code esp_det_line(cpl_table *spec_cont_region, double det_line_thres, double det_line_resol, int det_line_smwidth, cpl_table **line_table)
Performs a local normalization of the spectrum region.
cpl_error_code esp_fit_lcont(cpl_table *spec_region, double cont_rejt, int cont_iter)
Performs a local normalization of the spectrum region.
cpl_error_code select_local_spec(cpl_table *spec_total, double ew_space, double lambda, cpl_table **spec_region)
Select_local_spec.
#define max(a, b)