CR2RE Pipeline Reference Manual 1.6.8
cr2res_wave.c
1
19#ifdef HAVE_CONFIG_H
20#include <config.h>
21#endif
22
23/*-----------------------------------------------------------------------------
24 Includes
25 -----------------------------------------------------------------------------*/
26
27#include <math.h>
28#include <string.h>
29
30#include <cpl.h>
31
32#include "irplib_wlxcorr.h"
33
34#include "cr2res_dfs.h"
35#include "cr2res_wave.h"
36#include "cr2res_trace.h"
37#include "cr2res_extract.h"
38#include "cr2res_io.h"
39#include "cr2res_pfits.h"
40#include "cr2res_utils.h"
41#include "cr2res_etalon.h"
42#include "cr2res_qc.h"
43#include "cr2res_pfits.h"
44
45/*-----------------------------------------------------------------------------
46 Defines
47 -----------------------------------------------------------------------------*/
48
49// The maximum difference to the neighbouring pixels to be considered a bad pix
50// For the fitting of the gaussian line centers
51#define MAX_DEVIATION_FOR_BAD_PIXEL 300
52#define SPEED_OF_LIGHT 299792.458
53
54/*-----------------------------------------------------------------------------
55 Functions prototypes
56 -----------------------------------------------------------------------------*/
57
58static cpl_vector * cr2res_wave_clean_spectrum(
59 const cpl_vector * spec_intens,
60 const cpl_polynomial * wl_poly,
61 const cpl_bivector * catalog,
62 double wl_error) ;
63static cpl_table * cr2res_wave_recompute_wl(
64 const cpl_table * spectra,
65 const cpl_table * tw) ;
66static cpl_bivector * cr2res_wave_gen_lines_spectrum(
67 const char * catalog,
68 cpl_polynomial * wavesol_init,
69 double wl_error,
70 double max_intensity,
71 int log_flag) ;
72static cpl_polynomial * cr2res_wave_line_fitting(
73 cpl_bivector * spectrum,
74 cpl_bivector * spectrum_err,
75 cpl_polynomial * wavesol_init,
76 const cpl_array * wave_error_init,
77 int order,
78 int trace_nb,
79 cpl_bivector * lines_list,
80 int degree,
81 int display,
82 cpl_vector ** sigma_fit,
83 cpl_array ** wavelength_error,
84 cpl_table ** lines_diagnostics) ;
85static cpl_polynomial * cr2res_wave_polyfit_1d(
86 cpl_matrix * px,
87 cpl_vector * py,
88 cpl_vector *sigma_py,
89 int degree,
90 const cpl_polynomial * solution_init,
91 cpl_array ** wavelength_error,
92 cpl_vector ** sigma_fit,
93 cpl_matrix ** cov) ;
94static int cr2res_poly(const double x[], const double a[], double *
95 result) ;
96static int cr2res_deriv_poly(const double x[], const double a[], double
97 * result) ;
98static int cr2res_gauss(const double x[], const double a[], double
99 *result) ;
100static int cr2res_gauss_derivative(const double x[], const double a[],
101 double result[]) ;
102static cpl_bivector * cr2res_wave_etalon_assign_fringes(
103 const cpl_vector * li,
104 const cpl_vector * li_true) ;
105static double cr2res_wave_etalon_get_x0(
106 cpl_vector * li,
107 double trueD) ;
108static double cr2res_wave_etalon_get_D(
109 cpl_vector * li) ;
110
111/*----------------------------------------------------------------------------*/
115/*----------------------------------------------------------------------------*/
116
119/*----------------------------------------------------------------------------*/
149/*----------------------------------------------------------------------------*/
151 cpl_table * tw_in,
152 cpl_table * spectra_tab,
153 const cpl_frame * catalog_frame,
154 int reduce_order,
155 int reduce_trace,
156 cr2res_wavecal_type wavecal_type,
157 int degree,
158 int xdegree, // X-disp degree for 2D
159 double wl_start,
160 double wl_end,
161 double wl_err,
162 double wl_shift,
163 int log_flag,
164 int fallback_input_wavecal_flag,
165 int keep_higher_degrees_flag,
166 int clean_spectrum,
167 int display,
168 double display_wmin,
169 double display_wmax,
170 int zp_order,
171 int grat1_order,
172 cpl_propertylist ** qcs,
173 cpl_table ** lines_diagnostics,
174 cpl_table ** extracted_out,
175 cpl_table ** trace_wave_out)
176{
177 const char * catalog_fname ;
178 cpl_bivector ** spectra ;
179 cpl_bivector ** spectra_err ;
180 cpl_polynomial ** wavesol_init ;
181 cpl_array ** wavesol_init_error ;
182 int * orders ;
183 int * traces_nb ;
184 int nb_traces ;
185 cpl_table * tw_out ;
186 cpl_table * extracted_out_loc ;
187 cpl_array * wl_array ;
188 int flag ;
189 cpl_polynomial * wave_sol_2d ;
190 cpl_table * lines_diagnostics_tmp ;
191 cpl_table * lines_diagnostics_loc ;
192 cpl_propertylist * qcs_plist ;
193 cpl_polynomial * wave_sol_1d ;
194 cpl_array * wl_err_array ;
195 int trace_id, order, i;
196 double best_xcorr, lines_resol_fwhm, lines_pos,
197 lines_resol, lines_intens, qc_wlen_central,
198 qc_wlen_disp ;
199 int out_wl_array_size, init_wl_array_size, order_idx,
200 degree_out ;
201
202 /* Check Entries */
203 if (wavecal_type != CR2RES_XCORR && wavecal_type != CR2RES_LINE1D &&
204 wavecal_type != CR2RES_LINE2D && wavecal_type != CR2RES_ETALON)
205 return -1 ;
206 if (catalog_frame == NULL && wavecal_type != CR2RES_ETALON)
207 return -1 ;
208 // Both LINE2D and ETALON only work with all orders
209 if (reduce_order > -1 &&
210 (wavecal_type == CR2RES_LINE2D || wavecal_type == CR2RES_ETALON))
211 return -1;
212
213
214 /* Initialise */
215 if (catalog_frame != NULL)
216 catalog_fname = cpl_frame_get_filename(catalog_frame) ;
217 else
218 catalog_fname = NULL ;
219 lines_diagnostics_loc = NULL ;
220
221 /* Load the spectra */
222 nb_traces = cpl_table_get_nrow(tw_in) ;
223 spectra = cpl_malloc(nb_traces * sizeof(cpl_bivector *));
224 spectra_err = cpl_malloc(nb_traces * sizeof(cpl_bivector *));
225 wavesol_init = cpl_malloc(nb_traces * sizeof(cpl_polynomial *));
226 wavesol_init_error = cpl_malloc(nb_traces * sizeof(cpl_array *));
227 orders = cpl_malloc(nb_traces * sizeof(int));
228 traces_nb = cpl_malloc(nb_traces * sizeof(int));
229
230 /* Loop over the traces spectra */
231 for (i=0 ; i<nb_traces ; i++) {
232 /* Initialise */
233 spectra[i] = spectra_err[i] = NULL ;
234 wavesol_init[i] = NULL ;
235 wavesol_init_error[i] = NULL ;
236 orders[i] = -1 ;
237 traces_nb[i] = -1 ;
238
239 /* Get Order and trace id */
240 order = cpl_table_get(tw_in, CR2RES_COL_ORDER, i, NULL) ;
241 trace_id = cpl_table_get(tw_in, CR2RES_COL_TRACENB, i, NULL) ;
242
243 /* Check if this order needs to be skipped */
244 if (reduce_order > -1 && order != reduce_order) {
245 continue ;
246 }
247
248 /* Check if this trace needs to be skipped */
249 if (reduce_trace > -1 && trace_id != reduce_trace) {
250 continue ;
251 }
252
253 cpl_msg_info(__func__, "Retrieve Order %d/Trace %d Spectrum",
254 order, trace_id) ;
255
256 /* Get the extracted spectrum */
257 if (cr2res_extract_EXTRACT1D_get_spectrum(spectra_tab, order, trace_id,
258 &(spectra[i]), &(spectra_err[i]))) {
259 cpl_msg_error(__func__, "Cannot get the extracted spectrum") ;
260 continue ;
261 }
262
263 /* Get the wavelength guess */
264 if (wl_start>0.0 && wl_end>0.0) {
265 wavesol_init[i] = cr2res_wave_estimate_compute(wl_start, wl_end) ;
266 } else {
267 if ((wavesol_init[i]=cr2res_get_trace_wave_poly(tw_in,
268 CR2RES_COL_WAVELENGTH, order, trace_id)) == NULL) {
269 cpl_msg_warning(__func__, "Cannot get the WL guess") ;
270 cpl_bivector_delete(spectra[i]);
271 cpl_bivector_delete(spectra_err[i]);
272 spectra[i] = spectra_err[i] = NULL ;
273 continue ;
274 }
275 }
276
277 /* Get the wavelength error */
278 if (wl_err>0.0) {
279 wavesol_init_error[i] = cpl_array_new(2, CPL_TYPE_DOUBLE);
280 cpl_array_set_double(wavesol_init_error[i], 0, wl_err) ;
281 cpl_array_set_double(wavesol_init_error[i], 1, wl_err) ;
282 } else {
283 if ((wavesol_init_error[i]=cpl_array_duplicate(
284 cpl_table_get_array(tw_in,
285 CR2RES_COL_WAVELENGTH_ERROR,
286 cr2res_get_trace_table_index(tw_in, order,
287 trace_id)))) == NULL) {
288 cpl_msg_warning(__func__, "Cannot get the WL ERROR guess") ;
289 cpl_bivector_delete(spectra[i]);
290 cpl_bivector_delete(spectra_err[i]);
291 cpl_polynomial_delete(wavesol_init[i]) ;
292 spectra[i] = spectra_err[i] = NULL ;
293 wavesol_init[i] = NULL ;
294 continue ;
295 }
296 }
297
298 /* Store the order / trace_nb */
299 orders[i] = order ;
300 traces_nb[i] = trace_id ;
301
302 /* Apply the shift */
303 if (fabs(wl_shift) > 1e-3) {
304 const cpl_size power = 0;
305 cpl_polynomial_set_coeff(wavesol_init[i], &power,
306 cpl_polynomial_get_coeff(wavesol_init[i], &power)+wl_shift);
307 }
308 }
309
310 /* Determine the output WL array size */
311 init_wl_array_size = -1 ;
312 for (i=0; i<nb_traces ; i++) {
313 if (wavesol_init[i] != NULL) {
314 init_wl_array_size = cpl_polynomial_get_degree(wavesol_init[i]) ;
315 init_wl_array_size++ ;
316 /* keep first not NULL */
317 if (init_wl_array_size > 0) break ;
318 }
319 }
320
321 /* Only if we use fallback or keep higher degrees is the ouput */
322 /* degree potentially higher than the requested one */
323 if (!fallback_input_wavecal_flag && !keep_higher_degrees_flag)
324 out_wl_array_size = degree+1 ;
325 else
326 if (degree+1 < init_wl_array_size)
327 out_wl_array_size = init_wl_array_size ;
328 else
329 out_wl_array_size = degree+1 ;
330 cpl_msg_debug(__func__, "Input/Output WL size (req. degree): %d/%d (%d)",
331 init_wl_array_size, out_wl_array_size, degree) ;
332 degree_out = out_wl_array_size - 1 ;
333
334 /* Output TRACE_WAVE copied from the input one */
335 tw_out = cpl_table_duplicate(tw_in) ;
336
337 /* Rename old Wavelength column and create new with right degree */
338 cpl_table_name_column(tw_out, CR2RES_COL_WAVELENGTH, "TMP_WL");
339 cpl_table_new_column_array(tw_out, CR2RES_COL_WAVELENGTH, CPL_TYPE_DOUBLE,
340 degree_out+1) ;
341 cpl_table_name_column(tw_out, CR2RES_COL_WAVELENGTH_ERROR, "TMP_WL_ERR");
342 cpl_table_new_column_array(tw_out, CR2RES_COL_WAVELENGTH_ERROR,
343 CPL_TYPE_DOUBLE, 2) ;
344
345 if (cpl_error_get_code()) {
346 cpl_msg_error(__func__, "Table error renaming old columns: %d",
347 cpl_error_get_code());
348 /* De-allocate */
349 for (i=0 ; i<nb_traces ; i++) {
350 if (spectra[i] != NULL) cpl_bivector_delete(spectra[i]) ;
351 if (spectra_err[i] != NULL)
352 cpl_bivector_delete(spectra_err[i]) ;
353 if (wavesol_init[i]!=NULL)
354 cpl_polynomial_delete(wavesol_init[i]) ;
355 if (wavesol_init_error[i]!=NULL)
356 cpl_array_delete(wavesol_init_error[i]) ;
357 }
358 cpl_free(spectra) ;
359 cpl_free(spectra_err) ;
360 cpl_free(wavesol_init) ;
361 cpl_free(wavesol_init_error) ;
362 cpl_free(orders) ;
363 cpl_free(traces_nb) ;
364 cpl_table_delete(tw_out) ;
365 return -1;
366 }
367
368 /* Copy incoming solution into output If explicitly requested */
369 if (fallback_input_wavecal_flag) {
370 for (i = 0; i < nb_traces; i++) {
371 int j;
372 /* Propagate WAVELENGTH */
373 const cpl_array * wl_array_tmp ;
374 wl_array_tmp = cpl_table_get_array(tw_out, "TMP_WL", i);
375 wl_array = cpl_array_new(degree_out+1, CPL_TYPE_DOUBLE);
376 for (j=0; j <= degree_out; j++) {
377 if ( j+1 > cpl_array_get_size(wl_array_tmp)){
378 cpl_array_set(wl_array, j, 0.0);
379 } else {
380 double coeff;
381 coeff = cpl_array_get_double(wl_array_tmp, j, &flag);
382 if (flag != 0)
383 cpl_msg_debug(__func__,"%d, %d, %s",j,flag,
384 cpl_error_get_where());
385 else cpl_array_set(wl_array, j, coeff);
386 }
387 }
388 cpl_table_set_array(tw_out, CR2RES_COL_WAVELENGTH, i, wl_array);
389 cpl_array_delete(wl_array);
390
391 /* Propagate WAVELENGTH_ERROR */
392 wl_array_tmp = cpl_table_get_array(tw_out, "TMP_WL_ERR", i);
393 cpl_table_set_array(tw_out, CR2RES_COL_WAVELENGTH_ERROR, i,
394 wl_array_tmp);
395 }
396 }
397 cpl_table_erase_column(tw_out, "TMP_WL");
398 cpl_table_erase_column(tw_out, "TMP_WL_ERR");
399 if (cpl_error_get_code()){
400 cpl_msg_error(__func__, "Table error after (not) propagating: %d",
401 cpl_error_get_code());
402 /* De-allocate */
403 for (i=0 ; i<nb_traces ; i++) {
404 if (spectra[i] != NULL) cpl_bivector_delete(spectra[i]) ;
405 if (spectra_err[i] != NULL)
406 cpl_bivector_delete(spectra_err[i]) ;
407 if (wavesol_init[i]!=NULL)
408 cpl_polynomial_delete(wavesol_init[i]) ;
409 if (wavesol_init_error[i]!=NULL)
410 cpl_array_delete(wavesol_init_error[i]) ;
411 }
412 cpl_free(spectra) ;
413 cpl_free(spectra_err) ;
414 cpl_free(wavesol_init) ;
415 cpl_free(wavesol_init_error) ;
416 cpl_free(orders) ;
417 cpl_free(traces_nb) ;
418 cpl_table_delete(tw_out) ;
419 return -1;
420 }
421
422 /* Allocate qcs_plist */
423 qcs_plist = cpl_propertylist_new() ;
424
425 /* Actual calibration */
426 if (wavecal_type == CR2RES_LINE2D) {
427 /* 2D Calibration */
428 if ((wave_sol_2d=cr2res_wave_2d(spectra, spectra_err, wavesol_init,
429 wavesol_init_error, orders, traces_nb, nb_traces,
430 catalog_fname, degree, xdegree, -1, 1, zp_order,
431 display,
432 &wl_err_array,
433 &lines_diagnostics_loc)) == NULL) {
434 cpl_msg_error(__func__, "Failed to compute 2d Wavelength solution");
435 /* De-allocate */
436 for (i=0 ; i<nb_traces ; i++) {
437 if (spectra[i] != NULL) cpl_bivector_delete(spectra[i]) ;
438 if (spectra_err[i] != NULL)
439 cpl_bivector_delete(spectra_err[i]) ;
440 if (wavesol_init[i]!=NULL)
441 cpl_polynomial_delete(wavesol_init[i]) ;
442 if (wavesol_init_error[i]!=NULL)
443 cpl_array_delete(wavesol_init_error[i]) ;
444 }
445 cpl_free(spectra) ;
446 cpl_free(spectra_err) ;
447 cpl_free(wavesol_init) ;
448 cpl_free(wavesol_init_error) ;
449 cpl_free(orders) ;
450 cpl_free(traces_nb) ;
451 cpl_propertylist_delete(qcs_plist);
452 cpl_table_delete(tw_out) ;
453 return -1 ;
454 }
455
456 /* Store the Solution in the table */
457 for (i = 0; i < nb_traces; i++) {
458 wave_sol_1d = cr2res_wave_poly_2d_to_1d(wave_sol_2d,
459 orders[i] + zp_order);
460 wl_array=cr2res_convert_poly_to_array(wave_sol_1d, degree_out+1);
461 cpl_polynomial_delete(wave_sol_1d);
462 if (wl_array != NULL) {
463 cpl_table_set_array(tw_out, CR2RES_COL_WAVELENGTH, i, wl_array);
464 cpl_array_delete(wl_array) ;
465 }
466 if (wl_err_array != NULL) {
467 cpl_table_set_array(tw_out, CR2RES_COL_WAVELENGTH_ERROR, i,
468 wl_err_array);
469 }
470 }
471 cpl_array_delete(wl_err_array) ;
472 cpl_polynomial_delete(wave_sol_2d);
473 } else if (wavecal_type == CR2RES_ETALON) {
474 /* 2D Etalon */
475 if ((wave_sol_2d = cr2res_etalon_wave_2d(
476 spectra, spectra_err, wavesol_init,
477 wavesol_init_error, orders, traces_nb, nb_traces,
478 degree, xdegree, zp_order, display,
479 &wl_err_array, &lines_diagnostics_loc)) == NULL) {
480 cpl_msg_error(__func__, "Failed to compute 2d Etalon solution");
481 /* De-allocate */
482 for (i=0 ; i<nb_traces ; i++) {
483 if (spectra[i] != NULL) cpl_bivector_delete(spectra[i]) ;
484 if (spectra_err[i] != NULL)
485 cpl_bivector_delete(spectra_err[i]) ;
486 if (wavesol_init[i]!=NULL)
487 cpl_polynomial_delete(wavesol_init[i]) ;
488 if (wavesol_init_error[i]!=NULL)
489 cpl_array_delete(wavesol_init_error[i]) ;
490 }
491 cpl_free(spectra) ;
492 cpl_free(spectra_err) ;
493 cpl_free(wavesol_init) ;
494 cpl_free(wavesol_init_error) ;
495 cpl_free(orders) ;
496 cpl_free(traces_nb) ;
497 cpl_propertylist_delete(qcs_plist);
498 cpl_table_delete(tw_out) ;
499 return -1 ;
500 }
501
502 /* Store the Solution in the table */
503 for (i = 0; i < nb_traces; i++) {
504 wave_sol_1d = cr2res_wave_poly_2d_to_1d(wave_sol_2d,
505 orders[i] + zp_order);
506 wl_array=cr2res_convert_poly_to_array(wave_sol_1d, degree_out+1);
507 cpl_polynomial_delete(wave_sol_1d);
508 if (wl_array != NULL) {
509 cpl_table_set_array(tw_out, CR2RES_COL_WAVELENGTH, i, wl_array);
510 cpl_array_delete(wl_array) ;
511 }
512 if (wl_err_array != NULL) {
513 cpl_table_set_array(tw_out, CR2RES_COL_WAVELENGTH_ERROR, i,
514 wl_err_array);
515 }
516 }
517 cpl_array_delete(wl_err_array) ;
518 cpl_polynomial_delete(wave_sol_2d);
519
520 } else {
521 /* 1D Calibration */
522 /* Loop over the traces spectra */
523
524 for (i=0 ; i<nb_traces ; i++) {
525 /* Get Order and trace id */
526 order = cpl_table_get(tw_out, CR2RES_COL_ORDER, i, NULL) ;
527 trace_id = cpl_table_get(tw_out, CR2RES_COL_TRACENB, i, NULL) ;
528
529 /* Check if this order needs to be skipped */
530 if (reduce_order > -1 && order != reduce_order) {
531 continue ;
532 }
533
534 /* Check if this trace needs to be skipped */
535 if (reduce_trace > -1 && trace_id != reduce_trace) {
536 continue ;
537 }
538
539 cpl_msg_info(__func__, "Process Order %d/Trace %d",
540 order, trace_id) ;
541
542 /* Call the Wavelength Calibration */
543 lines_diagnostics_tmp = NULL ;
544 if ((wave_sol_1d = cr2res_wave_1d(spectra[i], spectra_err[i],
545 wavesol_init[i], wavesol_init_error[i], order,
546 trace_id, wavecal_type, catalog_fname,
547 degree, clean_spectrum, log_flag,
548 keep_higher_degrees_flag,
549 display, display_wmin, display_wmax, &best_xcorr,
550 &lines_resol_fwhm, &lines_pos, &lines_resol,
551 &lines_intens,
552 &wl_err_array, &lines_diagnostics_tmp)) == NULL) {
553 cpl_msg_warning(__func__, "Cannot calibrate in Wavelength") ;
554 cpl_error_reset() ;
555 continue ;
556 }
557
558 /* Add The QC parameters */
559 if (best_xcorr > 0.0) {
560 char * qc_name = cpl_sprintf("%s-%02d-%02d",
561 CR2RES_HEADER_QC_WAVE_BESTXCORR, order, trace_id) ;
562 cpl_propertylist_append_double(qcs_plist, qc_name, best_xcorr);
563 cpl_free(qc_name) ;
564 }
565 if (lines_intens > 0.0) {
566 char * qc_name = cpl_sprintf("%s-%02d-%02d",
567 CR2RES_HEADER_QC_WAVE_LAMP_EFFIC, order, trace_id) ;
568 cpl_propertylist_append_double(qcs_plist, qc_name,
569 lines_intens);
570 cpl_free(qc_name) ;
571 }
572 if (lines_resol > 0.0) {
573 char * qc_name = cpl_sprintf("%s-%02d-%02d",
574 CR2RES_HEADER_QC_WAVE_RESOL, order, trace_id) ;
575 cpl_propertylist_append_double(qcs_plist, qc_name,
576 lines_resol);
577 cpl_free(qc_name) ;
578 }
579 if (lines_resol_fwhm > 0.0) {
580 char * qc_name = cpl_sprintf("%s-%02d-%02d",
581 CR2RES_HEADER_QC_WAVE_RESOL_FWHM, order, trace_id) ;
582 cpl_propertylist_append_double(qcs_plist, qc_name,
583 lines_resol_fwhm);
584 cpl_free(qc_name) ;
585 }
586 if (lines_pos > 0.0) {
587 char * qc_name = cpl_sprintf("%s-%02d-%02d",
588 CR2RES_HEADER_QC_WAVE_POS, order, trace_id) ;
589 cpl_propertylist_append_double(qcs_plist, qc_name,
590 lines_pos);
591 cpl_free(qc_name) ;
592 }
593 //CR2RES_HEADER_QC_WAVE_LAMP_EFFIC ?
594
595 /* Merge the lines_diagnostics */
596 if (lines_diagnostics_loc == NULL) {
597 /* First trace - Initial table */
598 lines_diagnostics_loc = lines_diagnostics_tmp ;
599 lines_diagnostics_tmp = NULL ;
600 } else if (lines_diagnostics_tmp != NULL) {
601 /* Merge with previous */
602 cpl_table_insert(lines_diagnostics_loc,
603 lines_diagnostics_tmp,
604 cpl_table_get_nrow(lines_diagnostics_loc)) ;
605 cpl_table_delete(lines_diagnostics_tmp) ;
606 }
607
608 /* Store the Solution in the table */
609 cpl_msg_debug(__func__, "Saving poly size %d, degree_out %d",
610 (int)cpl_polynomial_get_degree(wave_sol_1d), degree_out);
611 wl_array = cr2res_convert_poly_to_array(wave_sol_1d, degree_out+1) ;
612 if (wl_array != NULL) {
613 cpl_table_set_array(tw_out, CR2RES_COL_WAVELENGTH, i, wl_array);
614 cpl_array_delete(wl_array) ;
615 }
616 if (wl_err_array != NULL) {
617 cpl_table_set_array(tw_out, CR2RES_COL_WAVELENGTH_ERROR, i,
618 wl_err_array);
619 cpl_array_delete(wl_err_array) ;
620 }
621 cpl_polynomial_delete(wave_sol_1d);
622 }
623 }
624
625 /* Store the central Wavelength QC */
626 order_idx = cr2res_order_real_to_idx(grat1_order, zp_order) ;
627 qc_wlen_central = cr2res_qc_wave_central(tw_out, order_idx) ;
628 cpl_propertylist_append_double(qcs_plist,
629 CR2RES_HEADER_QC_WAVE_CENTWL, qc_wlen_central);
630 qc_wlen_disp = cr2res_qc_wave_disp(tw_out, order_idx) ;
631 cpl_propertylist_append_double(qcs_plist,
632 CR2RES_HEADER_QC_WAVE_DISPWL, qc_wlen_disp);
633
634 /* Recompute the extracted table wavelengths with the results */
635 extracted_out_loc = cr2res_wave_recompute_wl(spectra_tab, tw_out) ;
636
637 /* De-allocate */
638 for (i=0 ; i<nb_traces ; i++) {
639 if (spectra[i] != NULL) cpl_bivector_delete(spectra[i]) ;
640 if (spectra_err[i] != NULL) cpl_bivector_delete(spectra_err[i]) ;
641 if (wavesol_init[i]!=NULL) cpl_polynomial_delete(wavesol_init[i]) ;
642 if (wavesol_init_error[i]!=NULL)
643 cpl_array_delete(wavesol_init_error[i]) ;
644 }
645 cpl_free(spectra) ;
646 cpl_free(spectra_err) ;
647 cpl_free(wavesol_init) ;
648 cpl_free(wavesol_init_error) ;
649 cpl_free(orders) ;
650 cpl_free(traces_nb) ;
651
652 if (qcs != NULL) *qcs = qcs_plist ;
653 else cpl_propertylist_delete(qcs_plist) ;
654 if (extracted_out != NULL) *extracted_out = extracted_out_loc ;
655 else cpl_table_delete(extracted_out_loc) ;
656 *lines_diagnostics = lines_diagnostics_loc ;
657 *trace_wave_out = tw_out ;
658 return 0 ;
659}
660
661/*----------------------------------------------------------------------------*/
688/*----------------------------------------------------------------------------*/
689cpl_polynomial * cr2res_wave_1d(
690 cpl_bivector * spectrum,
691 cpl_bivector * spectrum_err,
692 cpl_polynomial * wavesol_init,
693 const cpl_array * wave_error_init,
694 int order,
695 int trace_nb,
696 cr2res_wavecal_type wavecal_type,
697 const char * catalog,
698 int degree,
699 int clean_spectrum,
700 int log_flag,
701 int keep_higher_degrees_flag,
702 int display,
703 double display_wmin,
704 double display_wmax,
705 double * best_xcorr,
706 double * lines_resol_fwhm,
707 double * lines_pos,
708 double * lines_resol,
709 double * lines_intens,
710 cpl_array ** wavelength_error,
711 cpl_table ** lines_diagnostics)
712{
713 cpl_bivector * spectrum_local ;
714 cpl_polynomial * solution ;
715 cpl_bivector * ref_spectrum ;
716 double wl_error_nm, pos ;
717
718 /* Check Inputs */
719 if (spectrum == NULL || spectrum_err == NULL || wavesol_init == NULL ||
720 wavelength_error == NULL || lines_diagnostics == NULL)
721 return NULL ;
722 if ((wavecal_type == CR2RES_XCORR || wavecal_type == CR2RES_LINE1D) &&
723 catalog == NULL) return NULL ;
724
725 /* Initialise */
726 solution = NULL ;
727 wl_error_nm = cpl_array_get_double(wave_error_init, 0, NULL) ;
728 *wavelength_error = NULL ;
729 *lines_diagnostics = NULL ;
730 if (lines_resol_fwhm != NULL) *lines_resol_fwhm = -1.0 ;
731 if (lines_resol != NULL) *lines_resol = -1.0 ;
732 if (lines_pos != NULL) *lines_pos = -1.0 ;
733 if (lines_intens != NULL) *lines_intens = -1.0 ;
734
735 /* Create the lines spectrum from the lines list */
736 if (catalog != NULL){
737 if ((ref_spectrum = cr2res_wave_gen_lines_spectrum(catalog,
738 wavesol_init, wl_error_nm, -1.0, log_flag)) == NULL){
739 cpl_msg_warning(__func__, "Could not read the line list");
740 ref_spectrum = NULL;
741 }
742 } else {
743 ref_spectrum = NULL ;
744 }
745
746 /* Clean the input spectrum if requested */
747 if (clean_spectrum && ref_spectrum != NULL) {
748 cpl_vector * clean_spec_vec ;
749 /* Clean the spectrum */
750 if ((clean_spec_vec = cr2res_wave_clean_spectrum(
751 cpl_bivector_get_y(spectrum), wavesol_init,
752 ref_spectrum, wl_error_nm)) == NULL) {
753 cpl_msg_warning(__func__, "Failed cleaning the spectrum lines") ;
754 cpl_error_reset() ;
755 spectrum_local = cpl_bivector_duplicate(spectrum) ;
756 } else {
757 /* Replace the input spectrum by the clean one */
758 spectrum_local = cpl_bivector_duplicate(spectrum);
759 cpl_vector_copy(cpl_bivector_get_y(spectrum_local), clean_spec_vec);
760 cpl_vector_delete(clean_spec_vec);
761 }
762 } else {
763 spectrum_local = cpl_bivector_duplicate(spectrum) ;
764 }
765
766 /* Switch on the possible methods */
767 if (wavecal_type == CR2RES_XCORR) {
768 /* Currently Hardcoded */
769 double slit_width = 2.0 ;
770 double fwhm = 2.0 ;
771 int cleaning_filter_size = 9 ;
772 solution = cr2res_wave_xcorr(spectrum_local, wavesol_init, wl_error_nm,
773 ref_spectrum, degree, keep_higher_degrees_flag,
774 cleaning_filter_size, slit_width, fwhm, display, best_xcorr,
775 wavelength_error) ;
776
777 if (cpl_error_get_code() != CPL_ERROR_NONE) {
778 cpl_error_reset();
779 } else {
780 cpl_bivector * spectrum_corrected ;
781 double * px ;
782 int i ;
783 double lamb0, disp ;
784 /* Recompute the wavelengths for the spectrum */
785 spectrum_corrected = cpl_bivector_duplicate(spectrum_local) ;
786 px = cpl_bivector_get_x_data(spectrum_corrected) ;
787 for (i=0 ; i<cpl_bivector_get_size(spectrum_corrected) ; i++)
788 px[i] = cpl_polynomial_eval_1d(wavesol_init,(double)(i+1),NULL);
789
790 /* Run the plot */
791 if (display) {
793 spectrum_corrected,
794 ref_spectrum,
795 " (INITIAL GUESS)",
796 display_wmin,
797 display_wmax) ;
798 }
799 cpl_bivector_delete(spectrum_corrected) ;
800
801 /* Recompute the wavelengths for the spectrum */
802 spectrum_corrected = cpl_bivector_duplicate(spectrum_local) ;
803 px = cpl_bivector_get_x_data(spectrum_corrected) ;
804 for (i=0 ; i<cpl_bivector_get_size(spectrum_corrected) ; i++)
805 px[i] = cpl_polynomial_eval_1d(solution, (double)(i+1), NULL) ;
806
807 /* Run the plot */
808 if (display) {
810 spectrum_corrected,
811 ref_spectrum,
812 " (COMPUTED SOLUTION)",
813 display_wmin,
814 display_wmax) ;
815 }
816
817 /* Compute QC.WAVE LAMP EFFIC */
818 if (lines_intens != NULL)
819 *lines_intens =
820 cr2res_qc_wave_lamp_effic(spectrum_corrected) ;
821
822 /* Compute QC.RESOLFWHM */
823 if (lines_resol_fwhm != NULL)
824 *lines_resol_fwhm =
825 cr2res_qc_wave_resol_fwhm(spectrum_corrected, &pos) ;
826 lamb0 = cpl_vector_get(cpl_bivector_get_x(spectrum_corrected), 0);
827 disp = (cpl_vector_get(cpl_bivector_get_x(spectrum_corrected),
828 cpl_bivector_get_size(spectrum_corrected)-1) - lamb0)
829 / cpl_bivector_get_size(spectrum_corrected) ;
830 cpl_bivector_delete(spectrum_corrected) ;
831
832 /* Compute QC.RESOL */
833 if (lines_resol != NULL)
834 *lines_resol = lamb0 / ( *lines_resol_fwhm * disp );
835 if (lines_pos != NULL) *lines_pos = pos ;
836
837// TODO : GENERATE THE LINES STATISTICS ---> *lines_diagnostics
838 }
839 } else if (wavecal_type == CR2RES_LINE1D) {
840 solution = cr2res_wave_line_fitting(spectrum_local, spectrum_err,
841 wavesol_init, wave_error_init, order, trace_nb,
842 ref_spectrum, degree, display, NULL, wavelength_error,
843 lines_diagnostics) ;
844 } else if (wavecal_type == CR2RES_ETALON) {
845
846// TODO : This should not happen any more ---> remove the if ()
847
848 solution = cr2res_wave_etalon(spectrum_local, spectrum_err,
849 wavesol_init, degree, wavelength_error);
850 }
851
852 cpl_bivector_delete(spectrum_local) ;
853 if (ref_spectrum != NULL) cpl_bivector_delete(ref_spectrum) ;
854
855 return solution ;
856}
857
858/*----------------------------------------------------------------------------*/
888/*----------------------------------------------------------------------------*/
889cpl_polynomial * cr2res_wave_2d(
890 cpl_bivector ** spectra,
891 cpl_bivector ** spectra_err,
892 cpl_polynomial ** wavesol_init,
893 cpl_array ** wavesol_init_err,
894 int * orders,
895 int * traces_nb,
896 int ninputs,
897 const char * catalog,
898 cpl_size degree,
899 cpl_size degree_y,
900 double threshold,
901 int n_iterations,
902 int zp_order,
903 int display,
904 cpl_array ** wavelength_error,
905 cpl_table ** lines_diagnostics)
906{
907 cpl_table * lines_diagnostics_loc ;
908 cpl_bivector * catalog_spec ;
909 cpl_vector * diff;
910 cpl_size old, new, i, j, k, nlines ;
911 cpl_polynomial * result ;
912 cpl_size degree_2d[2];
913 cpl_matrix * tmp_x;
914 cpl_vector * tmp_y;
915 cpl_vector * tmp_sigma;
916 cpl_vector * pos;
917 cpl_vector * pxa ;
918 cpl_vector * pxb ;
919 cpl_vector * py ;
920 cpl_vector * sigma_py ;
921 cpl_vector * heights;
922 cpl_vector * fit_errors;
923 cpl_vector * loc;
924 cpl_polynomial ** wavesol;
925 double pix_pos, lambda_cat, lambda_meas, line_width,
926 line_intens, fit_error, value;
927 int n, badpix;
928
929 /* Check Inputs */
930 if (spectra==NULL || spectra_err==NULL || wavesol_init==NULL ||
931 orders==NULL || catalog==NULL || wavelength_error == NULL ||
932 lines_diagnostics==NULL)
933 return NULL ;
934
935 /* Initialise */
936 pxa = pxb = NULL;
937 py = sigma_py = NULL ;
938 *lines_diagnostics = NULL ;
939 *wavelength_error = NULL;
940 result = NULL ;
941
942 /* Load Catalog */
943 if ((catalog_spec = cr2res_io_load_EMISSION_LINES(catalog)) == NULL) {
944 cpl_msg_error(__func__,"Failed to load the catalog");
945 return NULL ;
946 }
947 // n = cpl_bivector_get_size(catalog_spec);
948 // result = cpl_polynomial_new(2);
949 wavesol = cpl_malloc(ninputs * sizeof(cpl_vector *));
950 for (i = 0; i < ninputs; i++){
951 wavesol[i] = cpl_polynomial_duplicate(wavesol_init[i]);
952 }
953
954 for (k = 0; k < n_iterations; k++){
955
956 /* Loop on the input spectra */
957 for (i = 0; i < ninputs; i++){
958 // extract line data in 1 spectrum
959 if (cr2res_wave_extract_lines(spectra[i], spectra_err[i],
960 wavesol[i], wavesol_init_err[i], catalog_spec, -1, -1, display,
961 &tmp_x, &tmp_y, &tmp_sigma, &heights, &fit_errors) == -1) {
962 cpl_msg_warning(__func__, "Could not extract lines");
963 continue;
964 }
965
966 /* Create / Fill / Merge the lines diagnostics table */
967 if (lines_diagnostics != NULL && tmp_x != NULL) {
968 nlines = cpl_matrix_get_nrow(tmp_x) ;
969 cpl_msg_debug(__func__, "Number of lines: %"CPL_SIZE_FORMAT,nlines);
970 /* Create */
971 lines_diagnostics_loc =
973 /* Fill */
974 for (j=0 ; j<nlines ; j++) {
975 pix_pos = cpl_matrix_get(tmp_x, j, 0) ;
976 lambda_cat = cpl_vector_get(tmp_y, j) ;
977 lambda_meas = cpl_polynomial_eval_1d(wavesol_init[i], pix_pos,
978 NULL) ;
979 line_width = cpl_vector_get(tmp_sigma, j) ;
980 line_intens = cpl_vector_get(heights, j) ;
981 fit_error = cpl_vector_get(fit_errors, j) ;
982
983 cpl_table_set_int(lines_diagnostics_loc,
984 CR2RES_COL_ORDER, j, orders[i]) ;
985 cpl_table_set_int(lines_diagnostics_loc,
986 CR2RES_COL_TRACENB, j, traces_nb[i]) ;
987 cpl_table_set_double(lines_diagnostics_loc,
988 CR2RES_COL_MEASURED_LAMBDA, j, lambda_meas) ;
989 cpl_table_set_double(lines_diagnostics_loc,
990 CR2RES_COL_CATALOG_LAMBDA, j, lambda_cat);
991 cpl_table_set_double(lines_diagnostics_loc,
992 CR2RES_COL_DELTA_LAMBDA, j, lambda_cat-lambda_meas);
993 cpl_table_set_double(lines_diagnostics_loc,
994 CR2RES_COL_MEASURED_PIXEL, j, pix_pos);
995 cpl_table_set_double(lines_diagnostics_loc,
996 CR2RES_COL_LINE_WIDTH, j, line_width) ;
997 cpl_table_set_double(lines_diagnostics_loc,
998 CR2RES_COL_FIT_QUALITY, j, fit_error) ;
999 cpl_table_set_double(lines_diagnostics_loc,
1000 CR2RES_COL_INTENSITY, j, line_intens) ;
1001 }
1002
1003 /* Merge */
1004 if (*lines_diagnostics == NULL) {
1005 *lines_diagnostics = lines_diagnostics_loc ;
1006 lines_diagnostics_loc = NULL ;
1007 } else if (lines_diagnostics_loc != NULL) {
1008 /* Merge with previous */
1009 cpl_table_insert(*lines_diagnostics, lines_diagnostics_loc,
1010 cpl_table_get_nrow(*lines_diagnostics)) ;
1011 cpl_table_delete(lines_diagnostics_loc) ;
1012 }
1013 }
1014
1015 // append new data onto existing vectors/matrices
1016 if (tmp_y != NULL) new = cpl_vector_get_size(tmp_y);
1017 else new = 0;
1018
1019 if (pxa == NULL){
1020 // First order to run
1021 pxa = cpl_vector_new(new);
1022 pxb = cpl_vector_new(new);
1023 py = cpl_vector_new(new);
1024 sigma_py = cpl_vector_new(new);
1025 old = 0;
1026 } else {
1027 old = cpl_vector_get_size(py);
1028 cpl_vector_set_size(py, old + new);
1029 cpl_vector_set_size(sigma_py, old + new);
1030 cpl_vector_set_size(pxa, old + new);
1031 cpl_vector_set_size(pxb, old + new);
1032 }
1033
1034 for (j = 0; j < new; j++){
1035 cpl_vector_set(py, old + j, cpl_vector_get(tmp_y, j));
1036 cpl_vector_set(sigma_py, old + j, cpl_vector_get(tmp_sigma, j));
1037 cpl_vector_set(pxa, old + j, cpl_matrix_get(tmp_x, j, 0));
1038 cpl_vector_set(pxb, old + j, orders[i] + zp_order);
1039 }
1040 cpl_matrix_delete(tmp_x);
1041 cpl_vector_delete(tmp_y);
1042 cpl_vector_delete(tmp_sigma);
1043 cpl_vector_delete(heights);
1044 cpl_vector_delete(fit_errors);
1045 }
1046
1047
1048 if (pxa == NULL){
1049 // No orders ran successfully
1050 cpl_msg_error(__func__, "No lines could be extracted in any order");
1051 cpl_bivector_delete(catalog_spec) ;
1052 // cpl_polynomial_delete(result);
1053 return NULL;
1054 }
1055
1056 n = cpl_vector_get_size(py);
1057 for (i = 0; i < n - 1; i++){
1058 degree_2d[0] = degree ;
1059 degree_2d[1] = degree_y ;
1060 // error = cpl_polynomial_fit(result, px, NULL, py, NULL, TRUE, NULL,
1061 // degree_2d);
1062 result = cr2res_polyfit_2d(pxa, pxb, py, degree_2d);
1063
1064 // in case something went wrong during fitting
1065 if (result == NULL){
1066 cpl_msg_error(__func__, "Polynomial fit failed. Error: %s",
1067 cpl_error_get_message());
1068 // cpl_polynomial_delete(result);
1069 cpl_vector_delete(pxa);
1070 cpl_vector_delete(pxb);
1071 cpl_vector_delete(py);
1072 cpl_vector_delete(sigma_py);
1073 cpl_bivector_delete(catalog_spec) ;
1074 cpl_error_reset();
1075 return NULL;
1076 }
1077
1078 // reject outliers and do it all over again!
1079 // Determine Residuals of each point
1080 loc = cpl_vector_new(2);
1081 pos = cpl_vector_new(cpl_vector_get_size(py));
1082
1083 for (j = 0; j < cpl_vector_get_size(py); j++)
1084 {
1085 cpl_vector_set(loc, 0, cpl_vector_get(pxa, j));
1086 cpl_vector_set(loc, 1, cpl_vector_get(pxb, j));
1087 value = cpl_polynomial_eval(result, loc);
1088 cpl_vector_set(pos, j, value);
1089 }
1090 cpl_vector_delete(loc);
1091
1092 // Convert wavelength to km/s
1093 cpl_vector_subtract(pos, py);
1094 cr2res_vector_abs(pos);
1095 cpl_vector_divide(pos, py);
1096 cpl_vector_multiply_scalar(pos, SPEED_OF_LIGHT);
1097
1098 // Reject largest value if above the threshold
1099 // threshold is negative skip ahead
1100 if (threshold < 0) {
1101 cpl_vector_delete(pos);
1102 break;
1103 }
1104 if (cpl_vector_get_max(pos) > threshold){
1105 j = cpl_vector_get_maxpos(pos);
1106 cr2res_vector_erase_element(pxa, j);
1107 cr2res_vector_erase_element(pxb, j);
1108 cr2res_vector_erase_element(py, j);
1109 cr2res_vector_erase_element(sigma_py, j);
1110 cpl_vector_delete(pos);
1111 } else {
1112 cpl_vector_delete(pos);
1113 break;
1114 }
1115 }
1116
1117 if (k != n_iterations - 1){
1118 cpl_vector_delete(pxa);
1119 cpl_vector_delete(pxb);
1120 cpl_vector_delete(py);
1121 cpl_vector_delete(sigma_py);
1122 pxa = pxb = NULL;
1123 }
1124
1125 for (i = 0; i < ninputs; i++){
1126 cpl_polynomial_delete(wavesol[i]);
1127 wavesol[i] = cr2res_wave_poly_2d_to_1d(result, orders[i] + zp_order);
1128 }
1129 }
1130
1131 n = cpl_vector_get_size(py);
1132 // Calculate absolute difference between polynomial and
1133 // catalog value for each line
1134 // use px and py, so that only good lines are used
1135 diff = cpl_vector_new(n);
1136 pos = cpl_vector_new(2);
1137 for (i = 0; i < n; i++){
1138 cpl_vector_set(pos, 0, cpl_vector_get(pxa, i));
1139 cpl_vector_set(pos, 1, cpl_vector_get(pxb, i));
1140 cpl_vector_set(diff, i, fabs(
1141 cpl_polynomial_eval(result, pos)
1142 - cpl_vector_get(py, i)));
1143 }
1144
1145 // Set wavelength_error to mean and max difference
1146 if (*wavelength_error == NULL)
1147 *wavelength_error = cpl_array_new(2, CPL_TYPE_DOUBLE);
1148 cpl_array_set_double(*wavelength_error, 0,
1149 cpl_vector_get_mean(diff));
1150 cpl_array_set_double(*wavelength_error, 1,
1151 cpl_vector_get_max(diff));
1152
1153 /* Update the lines diagnostics table with the new solution*/
1154 if (lines_diagnostics != NULL) {
1155 nlines = cpl_table_get_nrow(*lines_diagnostics);
1156 /* Fill */
1157 i = 0;
1158 for (j=0 ; j<nlines ; j++) {
1159 int order, trace;
1160 pix_pos = cpl_table_get_double(*lines_diagnostics,
1161 CR2RES_COL_MEASURED_PIXEL, j, &badpix);
1162 lambda_cat = cpl_table_get_double(*lines_diagnostics,
1163 CR2RES_COL_CATALOG_LAMBDA, j, &badpix);
1164 order = cpl_table_get_int(*lines_diagnostics,
1165 CR2RES_COL_ORDER, j, &badpix);
1166 trace = cpl_table_get_int(*lines_diagnostics,
1167 CR2RES_COL_TRACENB, j, &badpix);
1168
1169 // Assure we pick the correct wave solution
1170 while(!(orders[i] == order && traces_nb[i] == trace)){
1171 i++;
1172 // reset it in case we circle around, but this should not happen
1173 if (i == ninputs) i = 0;
1174 }
1175 // Calculate the new lambda and update the table
1176 lambda_meas = cpl_polynomial_eval_1d(wavesol[i], pix_pos, NULL) ;
1177 cpl_table_set_double(*lines_diagnostics,
1178 CR2RES_COL_CATALOG_LAMBDA, j, lambda_cat);
1179 cpl_table_set_double(*lines_diagnostics,
1180 CR2RES_COL_DELTA_LAMBDA, j, lambda_cat-lambda_meas);
1181 }
1182 }
1183
1184 // Cleanup
1185 cpl_vector_delete(diff);
1186 cpl_vector_delete(pos);
1187 cpl_bivector_delete(catalog_spec) ;
1188 for (i = 0; i < ninputs; i++){
1189 cpl_polynomial_delete(wavesol[i]);
1190 }
1191 cpl_free(wavesol);
1192 cpl_vector_delete(pxa);
1193 cpl_vector_delete(pxb);
1194 cpl_vector_delete(py);
1195 cpl_vector_delete(sigma_py);
1196
1197 return result;
1198}
1199
1200/*----------------------------------------------------------------------------*/
1220/*----------------------------------------------------------------------------*/
1221cpl_polynomial * cr2res_wave_xcorr(
1222 cpl_bivector * spectrum,
1223 cpl_polynomial * wavesol_init,
1224 double wl_error,
1225 cpl_bivector * lines_list,
1226 int degree,
1227 int keep_higher_degrees_flag,
1228 int cleaning_filter_size,
1229 double slit_width,
1230 double fwhm,
1231 int display,
1232 double * best_xcorr,
1233 cpl_array ** wavelength_error)
1234{
1235 cpl_vector * wl_errors ;
1236 cpl_polynomial * sol ;
1237 cpl_polynomial * sol_guess ;
1238 cpl_vector * spec_clean ;
1239 double * pspec_clean ;
1240 cpl_bivector * lines_list_filtered ;
1241 cpl_vector * xcorrs ;
1242 double wl_min, wl_max, wl_error_nm, wl_error_pix ;
1243 int i, nsamples, degree_loc, npolys ;
1244
1245 /* Check Entries */
1246 if (spectrum == NULL || wavesol_init == NULL || lines_list == NULL
1247 || best_xcorr == NULL)
1248 return NULL ;
1249
1250 /* Compute wl boundaries */
1251 wl_min = cpl_polynomial_eval_1d(wavesol_init, 1, NULL);
1252 wl_max = cpl_polynomial_eval_1d(wavesol_init, CR2RES_DETECTOR_SIZE, NULL);
1253
1254 cpl_msg_info(__func__, "Wl Range Input : %g - %g", wl_min, wl_max) ;
1255
1256 /* Clean the spectrum from the low frequency signal if requested */
1257 if (cleaning_filter_size > 0) {
1258 cpl_vector * filtered ;
1259 cpl_msg_info(__func__, "Low Frequency removal from spectrum") ;
1260 cpl_msg_indent_more() ;
1261 /* Subtract the low frequency part */
1262 if ((filtered=cpl_vector_filter_median_create(
1263 cpl_bivector_get_y(spectrum),
1264 cleaning_filter_size))==NULL){
1265 cpl_msg_error(__func__, "Cannot filter the spectrum") ;
1266 spec_clean = cpl_vector_duplicate(cpl_bivector_get_y(spectrum)) ;
1267 } else {
1268 spec_clean = cpl_vector_duplicate(cpl_bivector_get_y(spectrum)) ;
1269 cpl_vector_subtract(spec_clean, filtered) ;
1270 cpl_vector_delete(filtered) ;
1271 }
1272 cpl_msg_indent_less() ;
1273 } else {
1274 spec_clean = cpl_vector_duplicate(cpl_bivector_get_y(spectrum)) ;
1275 }
1276
1277 /* Currently the lines list is not filtered */
1278 lines_list_filtered = cpl_bivector_duplicate(lines_list) ;
1279
1280 /* Remove Negative values */
1281 pspec_clean = cpl_vector_get_data(spec_clean) ;
1282 for (i=0 ; i<cpl_vector_get_size(spec_clean) ; i++){
1283 if (pspec_clean[i] < 0.0) pspec_clean[i] = 0 ;
1284 if (isnan(pspec_clean[i])) pspec_clean[i] = 0;
1285 }
1286
1287 /* Display */
1288 if (display) {
1289 /* Plot Catalog Spectrum */
1290 cpl_plot_bivector(
1291 "set grid;set xlabel 'Wavelength (nm)';set ylabel 'Emission';",
1292 "t 'Catalog Spectrum' w impulses", "",
1293 lines_list_filtered);
1294 /* Plot Extracted Spectrum */
1295 cpl_plot_vector(
1296 "set grid;set xlabel 'Position (Pixel)';set ylabel 'Intensity (ADU/sec)';",
1297 "t 'Cleaned Extracted spectrum' w lines", "", spec_clean);
1298 }
1299
1300 /* Prepare inputs for X-corr */
1301 degree_loc = degree ;
1302 npolys = 50000;
1303 nsamples = pow(npolys,1.0/(degree_loc+1)) +1 ;
1304 wl_error_nm = wl_error ;
1305 sol_guess = wavesol_init ;
1306 wl_error_pix = CR2RES_DETECTOR_SIZE *wl_error_nm/(wl_max-wl_min) ;
1307 wl_errors = cpl_vector_new(degree_loc+1) ;
1308 cpl_vector_fill(wl_errors, wl_error_nm) ;
1309
1310 cpl_msg_info(__func__,
1311 "XCORR: Deg:%d - Err:%g nm (%g pix) - %d samples -> %g polys",
1312 degree_loc,wl_error_nm,wl_error_pix,nsamples,pow(nsamples,
1313 degree_loc+1)) ;
1314 cpl_msg_indent_more() ;
1315 if ((sol = keep_higher_degrees_flag ?
1316 irplib_wlxcorr_best_poly_prop(spec_clean, lines_list_filtered,
1317 degree_loc, sol_guess, wl_errors, nsamples, slit_width,
1318 fwhm, best_xcorr, NULL, &xcorrs)
1319 : irplib_wlxcorr_best_poly(spec_clean, lines_list_filtered,
1320 degree_loc, sol_guess, wl_errors, nsamples, slit_width,
1321 fwhm, best_xcorr, NULL, &xcorrs)
1322 ) == NULL) {
1323 cpl_msg_error(__func__, "Cannot get the best polynomial: %d",
1324 cpl_error_get_code()) ;
1325 cpl_vector_delete(wl_errors) ;
1326 cpl_vector_delete(spec_clean) ;
1327 if (xcorrs != NULL) cpl_vector_delete(xcorrs) ;
1328 cpl_error_reset() ;
1329 cpl_msg_indent_less() ;
1330 return NULL ;
1331 }
1332 cpl_vector_delete(wl_errors) ;
1333 cpl_msg_info(__func__, "Best Cross-Correlation factor: %g",
1334 *best_xcorr) ;
1335
1336 /* Compute the strongest line distance */
1337 int x_max = -1 ;
1338 double y_max = -1.0 ;
1339 for (i=0 ; i<cpl_vector_get_size(spec_clean) ; i++) {
1340 if (pspec_clean[i] > y_max) {
1341 y_max = pspec_clean[i] ;
1342 x_max = i+1 ;
1343 }
1344 }
1345 double delta_lambda =
1346 cpl_polynomial_eval_1d(sol, x_max, NULL) -
1347 cpl_polynomial_eval_1d(sol_guess, x_max, NULL) ;
1348 double delta_pix = delta_lambda * CR2RES_DETECTOR_SIZE / (wl_max-wl_min) ;
1349 cpl_msg_debug(__func__, "Correction : %g nm / %g pixels",
1350 delta_lambda, delta_pix) ;
1351 if (cpl_msg_get_level() == CPL_MSG_DEBUG) {
1352 cpl_polynomial_dump(sol, stdout);
1353 cpl_polynomial_dump(sol_guess, stdout);
1354 }
1355
1356 /* Plot the correlation values */
1357 if (display) {
1358 cpl_plot_vector("set grid;", "t 'Correlation values' w lines",
1359 "", xcorrs) ;
1360 }
1361 if (xcorrs != NULL) cpl_vector_delete(xcorrs) ;
1362 cpl_msg_indent_less() ;
1363
1364 cpl_vector_delete(spec_clean) ;
1365 cpl_bivector_delete(lines_list_filtered) ;
1366
1367 /* Use the Line Fitting to compute the Error */
1368 if (wavelength_error != NULL){
1369 /* TODO */
1370 *wavelength_error = cpl_array_new(2, CPL_TYPE_DOUBLE);
1371 cpl_array_set_double(*wavelength_error, 0, wl_error_nm) ;
1372 cpl_array_set_double(*wavelength_error, 1, wl_error_nm) ;
1373 }
1374 return sol ;
1375}
1376
1377/*----------------------------------------------------------------------------*/
1404/*----------------------------------------------------------------------------*/
1405cpl_polynomial * cr2res_wave_etalon(
1406 cpl_bivector * spectrum,
1407 const cpl_bivector * spectrum_err,
1408 cpl_polynomial * wavesol_init,
1409 int degree,
1410 cpl_array ** wavelength_error)
1411{
1412 cpl_bivector * is_should;
1413 cpl_vector * xi;
1414 cpl_vector * li;
1415 // cpl_vector * mi;
1416 cpl_vector * li_true;
1417 cpl_matrix * px;
1418 cpl_polynomial * result;
1419 double l0, trueD;
1420 int nxi, i, npeaks;
1421
1422 if (spectrum == NULL || spectrum_err == NULL ||
1423 wavesol_init == NULL || degree < 0 || wavelength_error == NULL)
1424 return NULL;
1425
1426 // Find etalon peaks xi in spectrum
1427 // TODO: Use Spectrum Error
1428 xi = cr2res_wave_etalon_measure_fringes(
1429 cpl_bivector_get_y(spectrum));
1430 nxi=cpl_vector_get_size(xi);
1431
1432 /* apply initial solution to get wavelength li at each point xi*/
1433 li = cr2res_polynomial_eval_vector(wavesol_init, xi);
1434
1435 /* Calculate delta lambda between peaks */
1436 trueD = cr2res_wave_etalon_get_D(li);
1437 cpl_msg_debug(__func__,"trueD: %e", trueD);
1438
1439 /* Set vector with correct wavelength values */
1440 // expected number of peaks (+-2 peaks just to be sure)
1441 l0 = cr2res_wave_etalon_get_x0(li, trueD);
1442 npeaks = (cpl_vector_get(li, nxi-1) - l0) / trueD + 4;
1443 l0 -= 2 * trueD;
1444 li_true = cpl_vector_new(npeaks);
1445 for (i=0; i<npeaks; i++) {
1446 cpl_vector_set(li_true, i,l0 + (trueD*i));
1447 }
1448
1449 // For each peak find the closest expected wavelength value
1450 // initial wavelength in x
1451 // estimated wavelength in y
1452 is_should = cr2res_wave_etalon_assign_fringes(li, li_true);
1453
1454 cpl_bivector_dump(is_should, NULL);
1455
1456 // polynomial fit to points, xi, li
1457 px = cpl_matrix_wrap(nxi, 1, cpl_vector_get_data(xi));
1458 result = cr2res_wave_polyfit_1d(px, cpl_bivector_get_y(is_should), NULL,
1459 degree, wavesol_init, wavelength_error, NULL, NULL);
1460
1461 cpl_matrix_unwrap(px);
1462 cpl_vector_delete(li_true);
1463
1464 cpl_bivector_delete(is_should);
1465 cpl_vector_delete(xi);
1466 cpl_vector_delete(li);
1467
1468 return result;
1469}
1470
1471/*----------------------------------------------------------------------------*/
1480/*----------------------------------------------------------------------------*/
1482 double wmin,
1483 double wmax)
1484{
1485 cpl_polynomial * poly ;
1486 cpl_size power ;
1487 double a, b ;
1488 int nbpix = CR2RES_DETECTOR_SIZE ;
1489
1490 /* Test entries */
1491 if (wmin <= 0.0) return NULL ;
1492 if (wmax <= 0.0) return NULL ;
1493 if (wmax <= wmin) return NULL ;
1494
1495 /* Compute polynomial coeffs */
1496 a = (wmax - wmin) / (nbpix-1) ;
1497 b = wmin - a ;
1498
1499 /* Create polynomial */
1500 poly = cpl_polynomial_new(1) ;
1501 power = 0 ;
1502 cpl_polynomial_set_coeff(poly, &power, b) ;
1503 power = 1 ;
1504 cpl_polynomial_set_coeff(poly, &power, a) ;
1505
1506 return poly ;
1507}
1508
1509/*----------------------------------------------------------------------------*/
1526/*----------------------------------------------------------------------------*/
1528 const char * filename,
1529 int detector,
1530 int order)
1531{
1532 double wmin, wmax, a, b, c, b1, b3 ;
1533 cpl_array * wl ;
1534 cpl_propertylist * plist ;
1535 int wished_ext_nb ;
1536
1537 /* Check Entries */
1538 if (filename == NULL) return NULL ;
1539 if (order < 0 || detector <= 0 || detector > CR2RES_NB_DETECTORS)
1540 return NULL ;
1541
1542 /********************/
1543 /* Get b3 for order */
1544 /********************/
1545 wished_ext_nb = cr2res_io_get_ext_idx(filename, 3, 1) ;
1546 plist = cpl_propertylist_load(filename, wished_ext_nb) ;
1547 wmin = cr2res_pfits_get_wstrt(plist, order) ;
1548 wmax = cr2res_pfits_get_wend(plist, order) ;
1549 cpl_propertylist_delete(plist) ;
1550 b3 = (wmax - wmin) / (CR2RES_DETECTOR_SIZE-1) ;
1551
1552 /********************/
1553 /* Get b1 for order */
1554 /********************/
1555 wished_ext_nb = cr2res_io_get_ext_idx(filename, 1, 1) ;
1556 plist = cpl_propertylist_load(filename, wished_ext_nb) ;
1557 wmin = cr2res_pfits_get_wstrt(plist, order) ;
1558 wmax = cr2res_pfits_get_wend(plist, order) ;
1559 cpl_propertylist_delete(plist) ;
1560 b1 = (wmax - wmin) / (CR2RES_DETECTOR_SIZE-1) ;
1561
1562 if (cpl_error_get_code()) {
1563 cpl_msg_warning(__func__, "Missing WMIN/WMAX for order %d", order) ;
1564 cpl_error_reset() ;
1565 c = 0.0 ;
1566 } else {
1567 c = (b3-b1)/(4.0*(CR2RES_DETECTOR_SIZE-1)) ;
1568 }
1569
1570 /* Load the propertylist */
1571 wished_ext_nb = cr2res_io_get_ext_idx(filename, detector, 1) ;
1572 plist = cpl_propertylist_load(filename, wished_ext_nb) ;
1573
1574 /* Get the values for this order */
1575 wmin = cr2res_pfits_get_wstrt(plist, order) ;
1576 wmax = cr2res_pfits_get_wend(plist, order) ;
1577 cpl_propertylist_delete(plist) ;
1578 if (cpl_error_get_code() != CPL_ERROR_NONE) {
1579 cpl_msg_error(__func__,
1580 "Cannot get WSTRT/WEND from header for Detector %d / Order %d",
1581 detector, order) ;
1582 cpl_error_reset();
1583 return NULL ;
1584 }
1585
1586 /* Compute polynomial coefficients */
1587 b = (wmax - wmin) / (CR2RES_DETECTOR_SIZE-1) ;
1588 b -= c * CR2RES_DETECTOR_SIZE;
1589 a = wmin - b ;
1590
1591 /* Create the array */
1592 wl = cpl_array_new(3, CPL_TYPE_DOUBLE) ;
1593 cpl_array_set(wl, 0, a) ;
1594 cpl_array_set(wl, 1, b) ;
1595 cpl_array_set(wl, 2, c) ;
1596 return wl ;
1597}
1598
1599/*----------------------------------------------------------------------------*/
1607/*----------------------------------------------------------------------------*/
1609 const cpl_table * trace_wave)
1610{
1611 hdrl_image * out ;
1612 cpl_image * out_ima ;
1613 double * pout_ima ;
1614 cpl_polynomial * upper_poly ;
1615 cpl_polynomial * lower_poly ;
1616 double upper_pos, lower_pos, wavelength ;
1617 cpl_size i, j, k, nrows, nx, ny ;
1618
1619 /* Check Entries */
1620 if (trace_wave == NULL) return NULL ;
1621
1622 /* Initialise */
1623 nrows = cpl_table_get_nrow(trace_wave) ;
1624
1625 /* Create the image */
1626 out = hdrl_image_new(CR2RES_DETECTOR_SIZE, CR2RES_DETECTOR_SIZE) ;
1627 out_ima = hdrl_image_get_image(out) ;
1628 nx = cpl_image_get_size_x(out_ima) ;
1629 ny = cpl_image_get_size_y(out_ima) ;
1630 pout_ima = cpl_image_get_data_double(out_ima) ;
1631
1632 /* Loop on the traces */
1633 for (k=0 ; k<nrows ; k++) {
1634 const cpl_array * tmp_array ;
1635 cpl_polynomial * wave_poly ;
1636 /* Check if there is a Wavelength Polynomial available */
1637 tmp_array = cpl_table_get_array(trace_wave, CR2RES_COL_WAVELENGTH, k) ;
1638 wave_poly = cr2res_convert_array_to_poly(tmp_array) ;
1639 if (wave_poly != NULL) {
1640 /* Get the Upper Polynomial */
1641 tmp_array = cpl_table_get_array(trace_wave, CR2RES_COL_UPPER, k) ;
1642 upper_poly = cr2res_convert_array_to_poly(tmp_array) ;
1643
1644 /* Get the Lower Polynomial */
1645 tmp_array = cpl_table_get_array(trace_wave, CR2RES_COL_LOWER, k) ;
1646 lower_poly = cr2res_convert_array_to_poly(tmp_array) ;
1647
1648 /* Check if all Polynomials are available */
1649 if (upper_poly == NULL || lower_poly == NULL) {
1650 if (upper_poly != NULL) cpl_polynomial_delete(upper_poly) ;
1651 if (lower_poly != NULL) cpl_polynomial_delete(lower_poly) ;
1652 cpl_msg_warning(__func__, "Cannot get UPPER/LOWER information");
1653 cpl_polynomial_delete(wave_poly) ;
1654 continue ;
1655 }
1656
1657 /* Set the Pixels in the trace */
1658 for (i=0 ; i<nx ; i++) {
1659 upper_pos = cpl_polynomial_eval_1d(upper_poly, i+1, NULL) ;
1660 lower_pos = cpl_polynomial_eval_1d(lower_poly, i+1, NULL) ;
1661 wavelength = cpl_polynomial_eval_1d(wave_poly, i+1, NULL) ;
1662 for (j=0 ; j<ny ; j++) {
1663 if (j+1 >= lower_pos && j+1 <= upper_pos)
1664 pout_ima[i+j*nx] = wavelength ;
1665 }
1666 }
1667 cpl_polynomial_delete(wave_poly) ;
1668 cpl_polynomial_delete(upper_poly) ;
1669 cpl_polynomial_delete(lower_poly) ;
1670 }
1671 }
1672 return out ;
1673}
1674
1675/*----------------------------------------------------------------------------*/
1684/*----------------------------------------------------------------------------*/
1686 cpl_polynomial ** poly_1ds,
1687 int * orders,
1688 int npolys,
1689 cpl_size xdegree)
1690{
1691 cpl_polynomial * out ;
1692 cpl_polynomial * fit;
1693 cpl_matrix * samppos;
1694 cpl_vector * values;
1695 int degree;
1696 cpl_size * coef_pos;
1697 cpl_error_code error;
1698
1699 /* Check Entries */
1700 if (poly_1ds == NULL || orders == NULL) return NULL ;
1701 if (npolys <= 0) return NULL ;
1702
1703 degree = cpl_polynomial_get_degree(poly_1ds[0]);
1704 // Check that all polynomials have the same degree
1705 for (int i = 1; i < npolys; i++)
1706 {
1707 if (degree != cpl_polynomial_get_degree(poly_1ds[0]))
1708 {
1709 return NULL;
1710 }
1711 }
1712
1713 out = cpl_polynomial_new(2);
1714 fit = cpl_polynomial_new(1);
1715 samppos = cpl_matrix_new(1, npolys);
1716 values = cpl_vector_new(npolys);
1717 coef_pos = cpl_malloc(2 * sizeof(cpl_size));
1718
1719 for (cpl_size i = 0; i <= degree; i++)
1720 {
1721 for (cpl_size j = 0; j < npolys; j++)
1722 {
1723 cpl_matrix_set(samppos, 0, j, orders[j]);
1724 cpl_vector_set(values, j, cpl_polynomial_get_coeff(poly_1ds[j],&i));
1725 }
1726
1727 error = cpl_polynomial_fit(fit, samppos, NULL, values, NULL, FALSE,
1728 NULL, &xdegree);
1729
1730 if (error != CPL_ERROR_NONE) printf("%i", error);
1731
1732 for (cpl_size j = 0; j <= xdegree; j++)
1733 {
1734 coef_pos[0] = (cpl_size) i;
1735 coef_pos[1] = (cpl_size) j;
1736 error = cpl_polynomial_set_coeff(out, coef_pos,
1737 cpl_polynomial_get_coeff(fit, &j));
1738 if (error != CPL_ERROR_NONE) printf("%i", error);
1739
1740 }
1741 }
1742
1743 // for each polynomial coefficient degree (0, 1, 2 ...)
1744 // fit a polynomial of degree ? to them vs the order number
1745 // new coefficients of the 2d polynomial are then given by
1746 // the fit results
1747
1748 // in:
1749 // f^{(1)}(x) = c_0^{(1)} + c_1^{(1)} x + c_2^{(1)} x^2 + ...
1750 // f^{(2)}(x) = c_0^{(2)} + c_1^{(2)} x + c_2^{(2)} x^2 + ...
1751 // ...
1752
1753 // out:
1754 // f(x, y) = c_{00} + c_{10} x + c_{01} y + c_{11} xy + ...
1755
1756 // fit:
1757 // c_0^{(y)} = c_{00} + c_{01} y + c_{02} y^2
1758 // c_1^{(y)} = c_{10} + c_{11} y + c_{12} y^2
1759
1760 cpl_polynomial_delete(fit);
1761 cpl_matrix_delete(samppos);
1762 cpl_vector_delete(values);
1763 cpl_free(coef_pos);
1764
1765 return out ;
1766}
1767
1768/*----------------------------------------------------------------------------*/
1775/*----------------------------------------------------------------------------*/
1777 cpl_polynomial * poly_2d,
1778 int order)
1779{
1780 cpl_polynomial * out ;
1781 cpl_polynomial * tmp ;
1782 cpl_size power = 0;
1783
1784 /* Check Entries */
1785 if (poly_2d == NULL) return NULL ;
1786 if (order < 0) return NULL ;
1787
1788 // Use polynomial extract
1789 // but that requires a constant polynomial
1790 tmp = cpl_polynomial_new(1);
1791 power = 0;
1792 cpl_polynomial_set_coeff(tmp, &power, order);
1793
1794 out = cpl_polynomial_extract(poly_2d, 1, tmp);
1795
1796 cpl_polynomial_delete(tmp);
1797
1798 return out ;
1799}
1800#ifdef CR2RES_UNUSED
1801/*----------------------------------------------------------------------------*/
1807/*----------------------------------------------------------------------------*/
1808char * cr2res_wave_method_print(
1809 cr2res_wavecal_type wavecal_type)
1810{
1811 char * out_str ;
1812 if (wavecal_type == CR2RES_XCORR)
1813 out_str = cpl_sprintf("Cross-Correlation") ;
1814 else if (wavecal_type == CR2RES_LINE1D)
1815 out_str = cpl_sprintf("1d Lines Fitting") ;
1816 else if (wavecal_type == CR2RES_LINE2D)
1817 out_str = cpl_sprintf("2d Lines Fitting") ;
1818 else if (wavecal_type == CR2RES_ETALON)
1819 out_str = cpl_sprintf("Etalon") ;
1820 else if (wavecal_type == CR2RES_UNSPECIFIED)
1821 out_str = cpl_sprintf("Unspecified") ;
1822 else
1823 out_str = cpl_sprintf("Unsupported") ;
1824 return out_str ;
1825}
1826#endif
1829/*----------------------------------------------------------------------------*/
1838/*----------------------------------------------------------------------------*/
1839static cpl_vector * cr2res_wave_clean_spectrum(
1840 const cpl_vector * spec_intens,
1841 const cpl_polynomial * wl_poly,
1842 const cpl_bivector * catalog,
1843 double wl_error)
1844{
1845 cpl_vector * x;
1846 cpl_vector * wave;
1847 cpl_vector * out ;
1848 double * pout ;
1849 cpl_size ncols, nlines, i, k ;
1850 double line, line_min, line_max;
1851
1852 /* Check entries */
1853 if (spec_intens == NULL || wl_poly == NULL || catalog == NULL ||
1854 wl_error < 0.0) return NULL ;
1855
1856
1857 /* Initialise */
1858 out = cpl_vector_duplicate(spec_intens) ;
1859 pout = cpl_vector_get_data(out) ;
1860 ncols = cpl_vector_get_size(out) ;
1861 nlines = cpl_bivector_get_size(catalog);
1862
1863 // If there are no lines in the catalog
1864 // just return all NaN
1865 if (cpl_bivector_get_size(catalog) == 0){
1866 for (i=0 ; i<ncols ; i++) {
1867 pout[i] = NAN ;
1868 }
1869 return out;
1870 }
1871
1872 x = cpl_vector_new(ncols);
1873 for (i = 0; i < ncols; i++) cpl_vector_set(x, i, i);
1874 wave = cr2res_polynomial_eval_vector(wl_poly, x);
1875
1876 // Lines MUST be in wavelength order
1877 k = 0;
1878 line = cpl_vector_get(cpl_bivector_get_x_const(catalog), k);
1879 line_min = line - wl_error;
1880 line_max = line + wl_error;
1881
1882 for (i=0 ; i<ncols ; i++) {
1883 if (cpl_vector_get(wave, i) < line_min){
1884 pout[i] = NAN;
1885 } else {
1886 // Find next relevant line
1887 while (cpl_vector_get(wave, i) > line_max)
1888 {
1889 k++;
1890 if (k < nlines){
1891 line = cpl_vector_get(cpl_bivector_get_x_const(catalog), k);
1892 line_min = line - wl_error;
1893 line_max = line + wl_error;
1894 } else {
1895 break;
1896 }
1897 }
1898 if (cpl_vector_get(wave, i) > line_max)
1899 pout[i] = NAN;
1900 if (cpl_vector_get(wave, i) < line_min)
1901 pout[i] = NAN;
1902 }
1903 }
1904 if (cpl_msg_get_level() == CPL_MSG_DEBUG) {
1905 cpl_vector_save(out, "debug_cleanedspec.fits", CPL_TYPE_DOUBLE,
1906 NULL, CPL_IO_CREATE);
1907 }
1908
1909 cpl_vector_delete(x);
1910 cpl_vector_delete(wave);
1911 return out ;
1912}
1913
1914/*----------------------------------------------------------------------------*/
1921/*----------------------------------------------------------------------------*/
1922static cpl_table * cr2res_wave_recompute_wl(
1923 const cpl_table * spectra,
1924 const cpl_table * tw)
1925{
1926 cpl_table * spectra_out ;
1927 char * wave_name ;
1928 double * pwave ;
1929 int i, j, nrows;
1930
1931 /* Check entries */
1932 if (spectra == NULL || tw == NULL) return NULL ;
1933
1934 /* Initialise */
1935 nrows = cpl_table_get_nrow(tw) ;
1936
1937 /* Create output spectra */
1938 spectra_out = cpl_table_duplicate(spectra) ;
1939
1940 /* Loop on the traces */
1941 for (i=0 ; i<nrows ; i++) {
1942
1943 cpl_polynomial * wave_poly ;
1944 int order, trace_nb ;
1945 order = cpl_table_get(tw, CR2RES_COL_ORDER, i, NULL) ;
1946 trace_nb = cpl_table_get(tw, CR2RES_COL_TRACENB, i, NULL);
1947
1948 if ((wave_poly = cr2res_get_trace_wave_poly(tw,
1949 CR2RES_COL_WAVELENGTH, order, trace_nb)) != NULL) {
1950
1951 /* Get the column name */
1952 wave_name = cr2res_dfs_WAVELENGTH_colname(order, trace_nb) ;
1953
1954 /* If the column is there, update it */
1955 if (cpl_table_has_column(spectra_out, wave_name)) {
1956 if ((pwave = cpl_table_get_data_double(spectra_out,
1957 wave_name)) == NULL) {
1958 cpl_msg_error(__func__, "Cannot access the wavelength") ;
1959 cpl_free(wave_name) ;
1960 cpl_polynomial_delete(wave_poly) ;
1961 cpl_table_delete(spectra_out) ;
1962 return NULL ;
1963 }
1964
1965 /* Update the Wavelength */
1966 for (j=0 ; j<cpl_table_get_nrow(spectra_out) ; j++)
1967 pwave[j] = cpl_polynomial_eval_1d(wave_poly, j+1, NULL) ;
1968 }
1969 cpl_free(wave_name) ;
1970 cpl_polynomial_delete(wave_poly) ;
1971 }
1972 }
1973
1974 return spectra_out ;
1975}
1976
1977/*----------------------------------------------------------------------------*/
1987/*----------------------------------------------------------------------------*/
1988static cpl_bivector * cr2res_wave_gen_lines_spectrum(
1989 const char * catalog,
1990 cpl_polynomial * wavesol_init,
1991 double wl_error,
1992 double max_intensity,
1993 int log_flag)
1994{
1995 cpl_bivector * lines ;
1996 cpl_bivector * lines_sub ;
1997 double * lines_sub_wl ;
1998 double * lines_sub_intens ;
1999 double wl_min, wl_max ;
2000 int i ;
2001
2002 /* Check Entries */
2003 if (catalog == NULL || wavesol_init == NULL) return NULL ;
2004
2005 /* Load the lines */
2006 lines = cr2res_io_load_EMISSION_LINES(catalog) ;
2007
2008 /* Extract the needed spectrum */
2009 wl_min = cpl_polynomial_eval_1d(wavesol_init, 1, NULL);
2010 wl_max = cpl_polynomial_eval_1d(wavesol_init, CR2RES_DETECTOR_SIZE, NULL);
2011 if ((lines_sub = irplib_wlxcorr_cat_extract(lines, wl_min-wl_error,
2012 wl_max+wl_error)) == NULL) {
2013 cpl_bivector_delete(lines) ;
2014 cpl_msg_warning(__func__, "Cannot find lines in [%.2f %.2f]",
2015 wl_min-wl_error, wl_max+wl_error) ;
2016 cpl_error_reset();
2017 return NULL;
2018 };
2019
2020 cpl_msg_debug(__func__,
2021"Extract %"CPL_SIZE_FORMAT" catalog lines in range %g (%g-%g)nm - %g (%g+%g)nm",
2022 cpl_bivector_get_size(lines_sub),
2023 wl_min-wl_error, wl_min, wl_error,
2024 wl_max+wl_error, wl_max, wl_error) ;
2025
2026
2027 /* Zero the beginning and the end */
2028 lines_sub_wl = cpl_bivector_get_x_data(lines_sub) ;
2029 lines_sub_intens = cpl_bivector_get_y_data(lines_sub) ;
2030 for (i=0 ; i<cpl_bivector_get_size(lines_sub) ; i++) {
2031 if (wl_min > 0)
2032 if (lines_sub_wl[i] < wl_min) lines_sub_intens[i] = 0.0 ;
2033 if (wl_max > 0)
2034 if (lines_sub_wl[i] > wl_max) lines_sub_intens[i] = 0.0 ;
2035 if (lines_sub_intens[i] > max_intensity && max_intensity > 0.0)
2036 lines_sub_intens[i] = 0.0 ;
2037 if (lines_sub_intens[i] > 0.0 && log_flag)
2038 lines_sub_intens[i] = log(lines_sub_intens[i]) ;
2039 }
2040
2041 /* Free and return */
2042 cpl_bivector_delete(lines) ;
2043 return lines_sub ;
2044}
2045
2046/*----------------------------------------------------------------------------*/
2073/*----------------------------------------------------------------------------*/
2074int cr2res_wave_fit_single_line(
2075 const cpl_vector * spec,
2076 const cpl_vector * unc,
2077 cpl_size pixel_pos,
2078 cpl_size window_size,
2079 cpl_size peak_width,
2080 int display,
2081 cpl_vector ** result)
2082{
2083 // For gaussian fit of each line
2084 // gauss = A * exp((x-mu)^2/(2*sig^2)) + cont
2085
2086 cpl_matrix * x;
2087 cpl_matrix * sigma_x;
2088 cpl_vector * y;
2089 cpl_vector * sigma_y;
2090 cpl_vector * a;
2091 int * ia;
2092 double red_chisq;
2093 double value;
2094 cpl_size k, j, n, spec_size;
2095 cpl_error_code error;
2096
2097 // Check inputs
2098 if (spec == NULL) return -1;
2099 if (unc != NULL){
2100 if (cpl_vector_get_size(spec) != cpl_vector_get_size(unc)) return -1;
2101 }
2102 if (pixel_pos < 0 || pixel_pos > cpl_vector_get_size(spec)) return -1;
2103 if (window_size < 0) return -1;
2104 if (peak_width < 0) peak_width = 1;
2105
2106 x = cpl_matrix_new(window_size, 1);
2107 sigma_x = NULL;
2108 y = cpl_vector_new(window_size);
2109 sigma_y = cpl_vector_new(window_size);
2110 a = cpl_vector_new(4);
2111 ia = NULL;
2112
2113 spec_size = cpl_vector_get_size(spec);
2114 //value = 0;
2115 k = pixel_pos - window_size / 2;
2116
2117 // Extract the spectrum within the window
2118 if (k < 0 || k + window_size >= spec_size){
2119 // if the window reaches outside the spectrum don't use the line
2120 // cpl_msg_error(__func__,
2121 // "Line at pixel %lli extends past the edge of the spectrum.",
2122 // pixel_pos);
2123 // cleanup
2124 cpl_vector_delete(a);
2125 cpl_matrix_delete(x);
2126 cpl_vector_delete(y);
2127 cpl_vector_delete(sigma_y);
2128 return -1;
2129 }
2130 n = 0;
2131 for (j = 0; j < window_size; j++){
2132 double value2;
2133 k = pixel_pos - window_size / 2 + j;
2134 value = cpl_vector_get(spec, k);
2135 value2 = unc != NULL ? cpl_vector_get(unc, k) : 0;
2136 if (isnan(value) || value < 0) value = 0;
2137 if (isnan(value2) || value2 <= 0){
2138 value2 = value == 0 ? 1 : sqrt(value);
2139 }
2140 if (value2 < 1) value2 = 1;
2141 if (value != 0){
2142 cpl_matrix_set(x, n, 0, k);
2143 cpl_vector_set(y, n, value);
2144 cpl_vector_set(sigma_y, n, value2);
2145 n++;
2146 }
2147 }
2148 if (n < 5){
2149 // Need at least 5 points for the fit (4 unknowns)
2150 // cpl_msg_error(__func__, "Not enough good points for line fit");
2151 cpl_vector_delete(a);
2152 cpl_matrix_delete(x);
2153 cpl_vector_delete(y);
2154 cpl_vector_delete(sigma_y);
2155 return -1;
2156 }
2157 cpl_matrix_set_size(x, n, 1);
2158 cpl_vector_set_size(y, n);
2159 cpl_vector_set_size(sigma_y, n);
2160
2161
2162 // Filter out bad pixels
2163 // This needs at least 2 elements
2164 //if (n > 2){ n>=5 from above.
2165
2166 cpl_vector *tmp;
2167 tmp = cpl_vector_duplicate(y);
2168 // First element
2169 if (fabs(cpl_vector_get(tmp, 0) - cpl_vector_get(tmp, 1)) >
2170 MAX_DEVIATION_FOR_BAD_PIXEL)
2171 cpl_vector_set(y, 0, cpl_vector_get(tmp, 1));
2172 // Elements in between
2173 for (j = 1; j < n - 1; j++) {
2174 double diff;
2175 diff = 2 * cpl_vector_get(tmp, j) - cpl_vector_get(tmp, j - 1) -
2176 cpl_vector_get(tmp, j + 1);
2177 diff = fabs(diff);
2178 if (diff > MAX_DEVIATION_FOR_BAD_PIXEL) {
2179 value =
2180 (cpl_vector_get(tmp, j - 1) + cpl_vector_get(tmp, j + 1)) / 2.;
2181 cpl_vector_set(y, j, value);
2182 }
2183 }
2184 // Last element
2185 if (fabs(cpl_vector_get(tmp, n - 1) - cpl_vector_get(tmp, n - 2)) >
2186 MAX_DEVIATION_FOR_BAD_PIXEL)
2187 cpl_vector_set(y, n - 1, cpl_vector_get(tmp, n - 2));
2188 cpl_vector_delete(tmp);
2189
2190 //}
2191
2192 // get initial guess for gaussian fit
2193 value = pixel_pos - window_size / 2 + cpl_vector_get_maxpos(y);
2194 cpl_vector_set(a, 0, value);
2195 cpl_vector_set(a, 1, peak_width);
2196 cpl_vector_set(a, 2, cpl_vector_get_max(y) - cpl_vector_get_min(y));
2197 cpl_vector_set(a, 3, cpl_vector_get_min(y));
2198
2199 // perform the fit
2200 red_chisq = -1;
2201 error = cpl_fit_lvmq(x, sigma_x, y, sigma_y, a, ia,
2202 &cr2res_gauss, &cr2res_gauss_derivative,
2203 CPL_FIT_LVMQ_TOLERANCE, CPL_FIT_LVMQ_COUNT,
2204 CPL_FIT_LVMQ_MAXITER * 100, NULL, &red_chisq, NULL);
2205 if (error != CPL_ERROR_NONE){
2206 // cpl_msg_error(__func__, "Could not fit line at pixel %lli. %s",
2207 // pixel_pos, cpl_error_get_message_default(error));
2208 cpl_vector_delete(a);
2209 cpl_matrix_delete(x);
2210 cpl_vector_delete(y);
2211 cpl_vector_delete(sigma_y);
2212 cpl_error_reset();
2213 return -1;
2214 }
2215
2216 // Set new pixel pos based on gaussian fit
2217 *result = cpl_vector_new(5);
2218 // 0: Peak position
2219 cpl_vector_set(*result, 0, cpl_vector_get(a, 0));
2220 // 1: Peak width
2221 cpl_vector_set(*result, 1, cpl_vector_get(a, 1));
2222 // 2 : Peak height
2223 cpl_vector_set(*result, 2, cpl_vector_get(a, 2));
2224 // 3: Peak base
2225 cpl_vector_set(*result, 3, cpl_vector_get(a, 3));
2226 // 4: Peak red_chisq
2227 cpl_vector_set(*result, 4, red_chisq);
2228
2229 if (display) {
2230 /* Plot Observation and Fit */
2231 const cpl_vector ** plot;
2232 cpl_vector * fit;
2233 cpl_vector * fit_x;
2234 double dbl, res;
2235 // TODO: Currently this will create a plot window for each line
2236 // afaik, to get them all in one plot, one needs to save all the fits
2237 // and then run one big plot_vectors (maybe bivectors) command with set multiplot
2238 // alternatively create a file for each line, instead of plot window
2239 plot = cpl_malloc(3 * sizeof(cpl_vector*));
2240 fit = cpl_vector_new(n);
2241 fit_x = cpl_vector_new(n);
2242
2243 for (j = 0; j < n; j++){
2244 dbl = cpl_matrix_get(x, j, 0);
2245 cr2res_gauss(&dbl, cpl_vector_get_data_const(a), &res);
2246 cpl_vector_set(fit, j, res);
2247 // dbl = cpl_polynomial_eval_1d(wavesol_init, dbl, NULL);
2248 cpl_vector_set(fit_x, j, dbl);
2249 }
2250 plot[0] = fit_x;
2251 plot[1] = y;
2252 plot[2] = fit;
2253 cpl_plot_vectors(
2254 "set grid;set xlabel 'Position (Pixel)';set ylabel 'Intensity (ADU/sec)';",
2255 "title 'Observed' w lines", "q", plot, 3);
2256 cpl_vector_delete(fit);
2257 cpl_vector_delete(fit_x);
2258 cpl_free(plot);
2259 cpl_error_reset();
2260 }
2261
2262 // Cleanup
2263 cpl_vector_delete(a);
2264 cpl_matrix_delete(x);
2265 cpl_vector_delete(y);
2266 cpl_vector_delete(sigma_y);
2267
2268 return 0;
2269}
2270
2271/*----------------------------------------------------------------------------*/
2303/*----------------------------------------------------------------------------*/
2304int cr2res_wave_extract_lines(
2305 cpl_bivector * spectrum,
2306 cpl_bivector * spectrum_err,
2307 cpl_polynomial * wavesol_init,
2308 const cpl_array * wave_error_init,
2309 cpl_bivector * lines_list,
2310 int window_size,
2311 double peak_width,
2312 int display,
2313 cpl_matrix ** px,
2314 cpl_vector ** py,
2315 cpl_vector ** sigma_py,
2316 cpl_vector ** heights,
2317 cpl_vector ** fit_error)
2318{
2319
2320 cpl_vector * wave_vec;
2321 cpl_vector * pixel_vec;
2322 cpl_vector * width_vec;
2323 cpl_vector * flag_vec;
2324 cpl_vector * height_vec;
2325 cpl_vector * fit_error_vec;
2326 cpl_vector * tmp;
2327 cpl_vector * result;
2328 const cpl_vector * spec;
2329 const cpl_vector * unc;
2330 const double * wave;
2331 cpl_size power;
2332 cpl_size i, k, ngood, spec_size, npossible;
2333 double pixel_pos, pixel_new;
2334 double max_wl, min_wl;
2335 double peak_height;
2336 int n;
2337
2338 /* Check Entries */
2339 if (spectrum == NULL || spectrum_err == NULL || wavesol_init == NULL ||
2340 lines_list == NULL){
2341 return -1;
2342 }
2343
2344 /* set window_size using the wave_error_init, scaled by the initial guess */
2345 if (window_size < 0 && wave_error_init != NULL)
2346 if (cpl_array_get_double(wave_error_init, 1, NULL) > 0) {
2347 power = 1;
2348 window_size = ceil(cpl_array_get_double(wave_error_init, 1,
2349 NULL) / fabs(cpl_polynomial_get_coeff(wavesol_init, &power)));
2350 }
2351
2352 if (window_size < CR2RES_WAVELENGTH_MIN_FIT_PIX)
2353 window_size = CR2RES_WAVELENGTH_MIN_FIT_PIX;
2354 cpl_msg_debug(__func__, "Using window size %d pix.", window_size);
2355
2356
2357 n = cpl_bivector_get_size(lines_list);
2358 spec = cpl_bivector_get_y_const(spectrum);
2359 unc = cpl_bivector_get_y_const(spectrum_err);
2360
2361 // Prepare fit data vectors
2362 pixel_vec = cpl_vector_new(n);
2363 height_vec = cpl_vector_new(n);
2364 fit_error_vec = cpl_vector_new(n);
2365 width_vec = cpl_vector_new(n);
2366 flag_vec = cpl_vector_new(n);
2367 cpl_vector_fill(flag_vec, 1);
2368
2369 // The number of good lines, start with all
2370 ngood = 0;
2371 // The number of possible lines to fit, for debugging only
2372 npossible = 0;
2373
2374 // evaluate the initial wavelength solution for all pixels
2375 // so that we can find the closest pixel position of each line
2376 spec_size = cpl_vector_get_size(spec);
2377 wave_vec = cpl_vector_new(spec_size);
2378 for (i = 0; i< spec_size; i++){
2379 cpl_vector_set(wave_vec, i, cpl_polynomial_eval_1d(wavesol_init, i, NULL));
2380 }
2381 max_wl = cpl_vector_get_max(wave_vec);
2382 min_wl = cpl_vector_get_min(wave_vec);
2383
2384
2385 // get line data
2386 wave = cpl_bivector_get_x_data_const(lines_list);
2387 //height = cpl_bivector_get_y_data_const(lines_list);
2388 // TODO width is not provided in the catalog at the moment,
2389 // use half window size instead?
2390 if (peak_width <= 0){
2391 peak_width = 1;
2392 }
2393
2394 // for each line fit a gaussian around guessed position
2395 // and find actual pixel position
2396 for (i = 0; i < n; i++){
2397
2398 // skip lines that are outside this wavelength region
2399 if ((wave[i] < min_wl) || (wave[i] > max_wl)){
2400 cpl_vector_set(flag_vec, i, 0);
2401 continue;
2402 }
2403 // The number of possible lines to fit, for debugging
2404 npossible++;
2405
2406 // cut out a part of the spectrum around each line
2407 // assumes that the wavelength vector is ascending !!!
2408 // pixel_pos = cpl_vector_find(wave_vec, wave[i]);
2409 tmp = cpl_vector_duplicate(wave_vec);
2410 cpl_vector_subtract_scalar(tmp, wave[i]);
2411 cpl_vector_multiply(tmp, tmp);
2412 pixel_pos = cpl_vector_get_minpos(tmp);
2413 cpl_vector_delete(tmp);
2414
2415 result = NULL;
2416 if (cr2res_wave_fit_single_line(spec, unc, pixel_pos,
2417 window_size, peak_width, display, &result)){
2418 cpl_vector_set(flag_vec, i, 0);
2419 continue;
2420 }
2421
2422 // Set new pixel pos based on gaussian fit
2423 pixel_new = cpl_vector_get(result, 0);
2424 peak_height = cpl_vector_get(result, 2);
2425 cpl_vector_set(pixel_vec, i, pixel_new);
2426 // width == uncertainty of wavelength position?
2427 cpl_vector_set(width_vec, i, cpl_vector_get(result, 1));
2428 cpl_vector_set(height_vec, i, peak_height);
2429 cpl_vector_set(fit_error_vec, i, cpl_vector_get(result, 4));
2430
2431 cpl_vector_delete(result);
2432
2433 // if fit to bad set flag to 0(False)
2434 // fit is bad, when
2435 // 1) it caused an error
2436 // 2) its chi square is large
2437 // 3) the gaussian is centered outside the window
2438 // 4) the fitted height is negative
2439 // 5) Peak is smaller than the noise level (SNR > 1)
2440 // 6) sigma is too large or too small
2441
2442 if (fabs(pixel_new - pixel_pos) > window_size
2443 || peak_height < 0){
2444 // TODO: the floor is not the same as the noise
2445 // maybe one could use the std around the floor as a noise estimate?
2446 // | cpl_vector_get(a, 2) < cpl_vector_get(a, 3)
2447 // TODO: Tweak these values and make into proper parameters?
2448 // | cpl_vector_get(a, 1) < 1 // lower line width limit
2449 // | cpl_vector_get(a, 1) > 6 // upper
2450 if (fabs(pixel_new - pixel_pos) > window_size)
2451 cpl_msg_debug(__func__, "Pixel position mismatch");
2452 if (peak_height < 0)
2453 cpl_msg_debug(__func__, "Negative Peak");
2454 cpl_vector_set(flag_vec, i, 0);
2455 cpl_error_reset();
2456 continue;
2457 }
2458 ngood++;
2459 }
2460
2461 cpl_msg_debug(__func__, "Using %lli out of %lli lines", ngood, npossible);
2462
2463 if (ngood == 0)
2464 {
2465 cpl_msg_warning(__func__, "No lines");
2466 cpl_vector_delete(wave_vec);
2467 cpl_vector_delete(pixel_vec);
2468 cpl_vector_delete(width_vec);
2469 cpl_vector_delete(flag_vec);
2470 cpl_vector_delete(height_vec);
2471 cpl_vector_delete(fit_error_vec);
2472 return -1;
2473 }
2474
2475 // Set vectors/matrices for polyfit
2476 // only need space for good lines, ignoring bad ones
2477 *px = cpl_matrix_new(ngood, 1);
2478 *py = cpl_vector_new(ngood);
2479 *sigma_py = cpl_vector_new(ngood);
2480 if (heights != NULL) *heights = cpl_vector_new(ngood);
2481 if (fit_error != NULL) *fit_error = cpl_vector_new(ngood);
2482
2483 k = 0;
2484 for (i = 0; i < n; i++){
2485 // Skip bad lines
2486 if (cpl_vector_get(flag_vec, i) == 1){
2487 cpl_matrix_set(*px, k, 0, cpl_vector_get(pixel_vec, i));
2488 cpl_vector_set(*py, k, wave[i]);
2489 cpl_vector_set(*sigma_py, k, fabs(cpl_vector_get(width_vec, i)));
2490 if (heights != NULL) cpl_vector_set(*heights, k,
2491 cpl_vector_get(height_vec, i));
2492 if (fit_error != NULL) cpl_vector_set(*fit_error, k,
2493 cpl_vector_get(fit_error_vec, i));
2494 k++;
2495 }
2496 }
2497 cpl_vector_delete(wave_vec);
2498 cpl_vector_delete(pixel_vec);
2499 cpl_vector_delete(width_vec);
2500 cpl_vector_delete(flag_vec);
2501 cpl_vector_delete(height_vec);
2502 cpl_vector_delete(fit_error_vec);
2503
2504 return 0;
2505}
2506
2507/*----------------------------------------------------------------------------*/
2516/*----------------------------------------------------------------------------*/
2517cpl_vector * cr2res_wave_etalon_measure_fringes(
2518 cpl_vector * spectrum)
2519{
2520 cpl_array * peaks;
2521 cpl_vector * peak_vec;
2522 cpl_vector * spec_thresh;
2523 cpl_vector * cur_peak;
2524 cpl_vector * X_all, *X_peak;
2525 int i, k ;
2526 int smooth = 35 ; // TODO: make free parameter?
2527 // interfringe ~30 in Y, ~70 in K
2528 double thresh = 1.0 ; // TODO: derive from read-out noise
2529 int max_num_peaks = 256 ;
2530 int min_len_peak = 5 ; //TODO: tweak or make parameter?;
2531 cpl_size nx ;
2532 double x0, sigma, area, offset;
2533
2534 nx = cpl_vector_get_size(spectrum) ;
2535 spec_thresh = cr2res_threshold_spec(spectrum, smooth, thresh) ;
2536 //cpl_plot_vector("", "w lines", "", spec_thresh) ;
2537
2538 if (cpl_msg_get_level() == CPL_MSG_DEBUG) {
2539 cpl_vector_save(spec_thresh, "debug_thresh.fits", CPL_TYPE_DOUBLE,
2540 NULL, CPL_IO_CREATE);
2541 cpl_vector_save(spectrum, "debug_spectrum.fits", CPL_TYPE_DOUBLE,
2542 NULL, CPL_IO_CREATE);
2543 }
2544
2545 /*Output array, values are invalid until set.*/
2546 peaks = cpl_array_new(max_num_peaks, CPL_TYPE_DOUBLE);
2547
2548 /* X-axis to cut out from for each peak */
2549 X_all = cpl_vector_new(nx);
2550 for (i=0; i<nx; i++) cpl_vector_set(X_all, i, (double)i+1) ;
2551
2552 for (i=0; i < nx; i++){
2553 int j;
2554 j = 0;
2555 while ( cpl_vector_get(spec_thresh, i) > -1 ) {
2556 j++;
2557 i++;
2558 }
2559 if (j < min_len_peak) continue;
2560 // cpl_msg_debug(__func__, "Peak length j=%d at i=%d",j,i);
2561 cur_peak = cpl_vector_extract(spec_thresh, i-j, i, 1) ;
2562 X_peak = cpl_vector_extract(X_all, i-j, i, 1) ;
2563
2564 if (cpl_vector_fit_gaussian(X_peak, NULL, cur_peak, NULL, CPL_FIT_ALL,
2565 &x0, &sigma, &area, &offset,
2566 NULL,NULL,NULL) != CPL_ERROR_NONE ) {
2567 cpl_msg_warning(__func__, "Fit at j=%d i=%d failed",j,i);
2568 cpl_vector_delete(cur_peak);
2569 cpl_vector_delete(X_peak);
2570 cpl_error_reset();
2571 continue;
2572 }
2573
2574 //cpl_msg_debug(__func__,"Fit: %.2f, %.2f, %.2f, %.2f",
2575 // x0, sigma, area, offset);
2576 if ((k = cpl_array_count_invalid(peaks)) <1)
2577 cpl_msg_error(__func__,"Output array overflow!");
2578 //cpl_msg_debug(__func__,"k=%d, x0=%g",k,x0);
2579 cpl_array_set_double(peaks, max_num_peaks - k, x0);
2580
2581 cpl_vector_delete(cur_peak);
2582 cpl_vector_delete(X_peak);
2583 }
2584
2585 /* Copy into output array */
2586 k = max_num_peaks - cpl_array_count_invalid(peaks);
2587 peak_vec = cpl_vector_new(k) ;
2588 for (i=0; i<k; i++)
2589 cpl_vector_set(peak_vec, i, cpl_array_get(peaks, i, NULL) );
2590
2591 if (cpl_msg_get_level() == CPL_MSG_DEBUG) {
2592 cpl_vector_save(peak_vec, "debug_peakpos.fits", CPL_TYPE_DOUBLE,
2593 NULL, CPL_IO_CREATE);
2594 }
2595
2596 cpl_vector_delete(spec_thresh) ;
2597 cpl_vector_delete(X_all) ;
2598 cpl_array_delete(peaks) ;
2599 return peak_vec;
2600}
2601
2602/*----------------------------------------------------------------------------*/
2603/*------------------- LINES FITTING METHODS -----------------------------*/
2604/*----------------------------------------------------------------------------*/
2605
2606/*----------------------------------------------------------------------------*/
2630/*----------------------------------------------------------------------------*/
2631static cpl_polynomial * cr2res_wave_line_fitting(
2632 cpl_bivector * spectrum,
2633 cpl_bivector * spectrum_err,
2634 cpl_polynomial * wavesol_init,
2635 const cpl_array * wave_error_init,
2636 int order,
2637 int trace_nb,
2638 cpl_bivector * lines_list,
2639 int degree,
2640 int display,
2641 cpl_vector ** sigma_fit,
2642 cpl_array ** wavelength_error,
2643 cpl_table ** lines_diagnostics)
2644{
2645 cpl_polynomial * result;
2646 cpl_matrix * cov;
2647 cpl_matrix * px;
2648 cpl_vector * py;
2649 cpl_vector * sigma_py;
2650 cpl_vector * heights ;
2651 cpl_vector * fit_errors ;
2652 cpl_size nlines;
2653
2654 /* Check Entries */
2655 if (spectrum == NULL || spectrum_err == NULL || wavesol_init == NULL ||
2656 lines_list == NULL)
2657 return NULL;
2658
2659 // extract line data in 1 spectrum
2660 if (cr2res_wave_extract_lines(spectrum, spectrum_err, wavesol_init,
2661 wave_error_init, lines_list, -1, -1, display, &px, &py, &sigma_py,
2662 &heights, &fit_errors) != 0) {
2663 cpl_msg_error(__func__, "Cannot extract lines") ;
2664 return NULL ;
2665 }
2666
2667 // fit polynomial to data points
2668 result = cr2res_wave_polyfit_1d(px, py, sigma_py, degree, wavesol_init,
2669 wavelength_error, sigma_fit, &cov);
2670
2671 /* Create / Fill the lines diagnostics table */
2672 if (lines_diagnostics != NULL && px != NULL && result != NULL) {
2673 cpl_size j;
2674
2675 nlines = cpl_matrix_get_nrow(px) ;
2676 cpl_msg_debug(__func__, "Number of lines: %"CPL_SIZE_FORMAT, nlines);
2677 /* Create */
2678 *lines_diagnostics = cr2res_dfs_create_lines_diagnostics_table(nlines) ;
2679 /* Fill */
2680 for (j=0 ; j<nlines ; j++) {
2681
2682 double pix_pos, lambda_cat, lambda_meas, line_width,
2683 line_intens, fit_error ;
2684 pix_pos = cpl_matrix_get(px, j, 0) ;
2685 lambda_cat = cpl_vector_get(py, j) ;
2686 lambda_meas = cpl_polynomial_eval_1d(result, pix_pos, NULL) ;
2687 line_width = cpl_vector_get(sigma_py, j) ;
2688 line_intens = cpl_vector_get(heights, j) ;
2689 fit_error = cpl_vector_get(fit_errors, j) ;
2690 cpl_table_set_int(*lines_diagnostics, CR2RES_COL_ORDER, j, order) ;
2691 cpl_table_set_int(*lines_diagnostics,CR2RES_COL_TRACENB,j,trace_nb);
2692 cpl_table_set_double(*lines_diagnostics,
2693 CR2RES_COL_MEASURED_LAMBDA, j, lambda_meas) ;
2694 cpl_table_set_double(*lines_diagnostics,
2695 CR2RES_COL_CATALOG_LAMBDA, j, lambda_cat) ;
2696 cpl_table_set_double(*lines_diagnostics,
2697 CR2RES_COL_DELTA_LAMBDA, j, lambda_cat-lambda_meas) ;
2698 cpl_table_set_double(*lines_diagnostics,
2699 CR2RES_COL_MEASURED_PIXEL, j, pix_pos) ;
2700 cpl_table_set_double(*lines_diagnostics,
2701 CR2RES_COL_LINE_WIDTH, j, line_width) ;
2702 cpl_table_set_double(*lines_diagnostics,
2703 CR2RES_COL_FIT_QUALITY, j, fit_error) ;
2704 cpl_table_set_double(*lines_diagnostics,
2705 CR2RES_COL_INTENSITY, j, line_intens) ;
2706 }
2707 }
2708
2709 if (display){
2710 cpl_size n = cpl_bivector_get_size(spectrum);
2711 nlines = cpl_vector_get_size(py);
2712 const cpl_bivector ** plot = cpl_malloc(2 * sizeof(cpl_bivector*));
2713
2714 cpl_vector * wave = cpl_bivector_get_x(spectrum);
2715 for (cpl_size i = 0; i < n; i++)
2716 cpl_vector_set(wave, i, cpl_polynomial_eval_1d(result, i+1, NULL));
2717
2718 cpl_bivector * lines = cpl_bivector_new(nlines);
2719 cpl_vector * pos = cpl_bivector_get_x(lines);
2720 cpl_vector * val = cpl_bivector_get_y(lines);
2721
2722 for (cpl_size i = 0; i < nlines; i++){
2723 cpl_vector_set(pos, i, cpl_vector_get(py, i)); // Wavelength
2724 cpl_vector_set(val, i, cpl_vector_get(heights, i));
2725 }
2726 plot[0] = spectrum;
2727 plot[1] = lines;
2728
2729 const char* options[] =
2730 {"title 'Observed' w lines", "title 'Lines' w points"};
2731
2732 cpl_plot_bivectors(
2733 "set grid;set xlabel 'Position (Wavelength)';set ylabel 'Intensity (ADU/sec)';",
2734 options , "", plot, 2);
2735 cpl_bivector_delete(lines);
2736 cpl_free(plot);
2737 }
2738 cpl_vector_delete(heights);
2739 cpl_vector_delete(fit_errors);
2740 cpl_matrix_delete(px);
2741 cpl_vector_delete(py);
2742 cpl_vector_delete(sigma_py);
2743
2744 if (result != NULL) {
2745 cpl_matrix_delete(cov);
2746 }
2747 return result;
2748}
2749
2750/*----------------------------------------------------------------------------*/
2767/*----------------------------------------------------------------------------*/
2768static cpl_polynomial * cr2res_wave_polyfit_1d(
2769 cpl_matrix * px,
2770 cpl_vector * py,
2771 cpl_vector *sigma_py,
2772 int degree,
2773 const cpl_polynomial * solution_init,
2774 cpl_array ** wavelength_error,
2775 cpl_vector ** sigma_fit,
2776 cpl_matrix ** cov)
2777{
2778
2779 // For polynomial fit
2780 cpl_polynomial * result = cpl_polynomial_new(1);
2781 cpl_matrix * sigma_px = NULL;
2782 cpl_vector * pa = cpl_vector_new(degree + 1 + 1);
2783 cpl_error_code error;
2784 cpl_size i, j;
2785 int *pia;
2786
2787 // first parameter of polynomial fit is the number of degrees (the value is fixed though)
2788 // the number of parameters is then polynomial degree + 1 (constant term) + 1 (number of degrees)
2789 // i.e. pa = degrees, a0, a1, ...
2790 cpl_vector_set(pa, 0, degree);
2791 pia = cpl_malloc((degree + 1 + 1) * sizeof(int));
2792 pia[0] = 0;
2793 for (i = 1; i < degree + 1 + 1; i++){
2794 pia[i] = 1;
2795 }
2796
2797 // initial guess for polynomial fit, based on passed initial guess
2798 for (j = 0; j < degree+1; j++){
2799 cpl_vector_set(pa, j+1, cpl_polynomial_get_coeff(solution_init, &j));
2800 }
2801
2802 // polynomial fit of px, py
2803 // with px: line pixel, py: line wavelength
2804 // I would use cpl_polynomial_fit, but that does not support error estimation
2805 error = cpl_fit_lvmq(px, sigma_px, py, sigma_py, pa, pia, &cr2res_poly, &cr2res_deriv_poly,
2806 CPL_FIT_LVMQ_TOLERANCE, CPL_FIT_LVMQ_COUNT,
2807 CPL_FIT_LVMQ_MAXITER, NULL, NULL, cov);
2808
2809
2810 if (error == CPL_ERROR_NONE){
2811 // Everything is fine
2812
2813 // if there are not enough data points left for the polynomial fit
2814 // the lvmq fit will complain about a singular matrix
2815
2816 // errors of the fit are the square root of the diagonal of the covariance matrix
2817 // assuming parameters are uncorrelated, better use the whole matrix
2818 // errors on the wavelength are then given using standard error propagation
2819 // s_wl**2 = s0**2 + s1**2 * x**2 + s2**2 * x**4 + ... si**2 * x**(2*i)
2820 for (i = 0; i < degree + 1; i++){
2821 cpl_polynomial_set_coeff(result, &i, cpl_vector_get(pa, i+1));
2822 if (sigma_fit != NULL)
2823 cpl_vector_set(*sigma_fit, i, sqrt(cpl_matrix_get(*cov, i+1, i+1)));
2824 }
2825
2826 if (wavelength_error != NULL){
2827
2828 cpl_vector * diff;
2829 // Calculate absolute difference between polynomial and catalog value for each line
2830 // use px and py, so that only good lines are used
2831 diff = cpl_vector_new(cpl_vector_get_size(py));
2832 if (*wavelength_error == NULL)
2833 *wavelength_error = cpl_array_new(2, CPL_TYPE_DOUBLE);
2834
2835
2836 for (i = 0; i < cpl_vector_get_size(py); i++){
2837 cpl_vector_set(diff, i, fabs(
2838 cpl_polynomial_eval_1d(result, cpl_matrix_get(px, i, 0), NULL)
2839 - cpl_vector_get(py, i)));
2840 }
2841 // Set wavelength_error to mean and max difference
2842 cpl_array_set_double(*wavelength_error, 0, cpl_vector_get_mean(diff));
2843 cpl_array_set_double(*wavelength_error, 1, cpl_vector_get_max(diff));
2844
2845 cpl_msg_debug(__func__, "Wave Error Mean: %g", cpl_array_get(*wavelength_error, 0, NULL));
2846 cpl_msg_debug(__func__, "Wave Error Max : %g", cpl_array_get(*wavelength_error, 1, NULL));
2847
2848 cpl_vector_delete(diff);
2849 }
2850 }
2851 cpl_free(pia);
2852 cpl_vector_delete(pa);
2853
2854 if (error != CPL_ERROR_NONE){
2855 cpl_msg_info(__func__, "Can't fit Polynomial");
2856 cpl_polynomial_delete(result);
2857 cpl_error_reset();
2858 return NULL;
2859 }
2860 return result;
2861}
2862
2863// Functions for polynomial fit in _wave_catalog()
2864static int cr2res_poly(const double x[], const double a[], double * result)
2865{
2866 int j;
2867 int degree = a[0];
2868
2869 *result = 0;
2870 for (j = degree; j > 0; j--){
2871 *result = x[0] * (a[1 + j] + *result);
2872 }
2873 *result = *result + a[1];
2874
2875 return 0;
2876}
2877
2878static int cr2res_deriv_poly(const double x[], const double a[], double * result)
2879{
2880 int i, j;
2881 const double degree = a[0];
2882
2883 result[0] = 0;
2884 for (j = 0; j < degree + 1; j++){
2885 result[j+1] = 1;
2886 for (i = 0; i < j; i++){
2887 result[j+1] *= x[0];
2888 }
2889 }
2890
2891 return 0;
2892}
2893
2894// following two are shamelessly taken from cpl_vector_fit_gauss
2895
2896/*----------------------------------------------------------------------------*/
2917/*----------------------------------------------------------------------------*/
2918static int cr2res_gauss(const double x[], const double a[], double *result)
2919{
2920 const double my = a[0];
2921 const double sigma = a[1];
2922
2923 if (sigma != 0.0) {
2924
2925 const double A = a[2];
2926 const double B = a[3];
2927
2928 *result = B + A * exp(- (x[0] - my)*(x[0] - my)
2929 / (2*sigma*sigma));
2930
2931 } else {
2932
2933 /* Dirac's delta function */
2934 *result = x[0] != my ? 0.0 : DBL_MAX;
2935 }
2936
2937 return 0;
2938}
2939
2940/*----------------------------------------------------------------------------*/
2960/*----------------------------------------------------------------------------*/
2961static int cr2res_gauss_derivative(const double x[], const double a[],
2962 double result[])
2963{
2964
2965 if (a[1] != 0.0) {
2966
2967 const double my = a[0];
2968 const double sigma = a[1];
2969 const double A = a[2];
2970 /* a[3] not used */
2971
2972 /* f(x) = B + A/sqrt(2 pi s^2) exp(-(x-my)^2/2s^2)
2973 *
2974 * df/d(my) = A/sqrt(2 pi s^2) exp(-(x-my)^2/2s^2) * (x-my) / s^2
2975 * = A * fac. * (x-my) / s^2
2976 * df/ds = A/sqrt(2 pi s^2) exp(-(x-my)^2/2s^2) * ((x-my)^2/s^3 -1/s)
2977 * = A * fac. * ((x-my)^2 / s^2 - 1) / s
2978 * df/dA = 1/sqrt(2 pi s^2) exp(-(x-my)^2/2s^2)
2979 * = fac.
2980 * df/dB = 1
2981 */
2982
2983
2984 const double factor = exp( -(x[0] - my)*(x[0] - my)/(2.0*sigma*sigma) );
2985
2986 result[0] = A * factor * (x[0]-my) / (sigma*sigma);
2987 result[1] = A * factor * ((x[0]-my)*(x[0]-my) / (sigma*sigma*sigma));
2988 result[2] = factor;
2989 result[3] = 1.0;
2990
2991 } else {
2992 /* Derivative of Dirac's delta function */
2993 result[0] = 0.0;
2994 result[1] = 0.0;
2995 result[2] = 0.0;
2996 result[3] = 0.0;
2997 }
2998
2999 return 0;
3000}
3001
3002/*----------------------------------------------------------------------------*/
3003/*------------------- ETALON RELATED METHODS -----------------------------*/
3004/*----------------------------------------------------------------------------*/
3005
3006/*----------------------------------------------------------------------------*/
3012/*----------------------------------------------------------------------------*/
3013static cpl_bivector * cr2res_wave_etalon_assign_fringes(
3014 const cpl_vector * li,
3015 const cpl_vector * li_true)
3016{
3017 int i, n;
3018 cpl_bivector * is_should;
3019 cpl_vector * is;
3020 cpl_vector * should;
3021
3022 n = cpl_vector_get_size(li);
3023 is_should = cpl_bivector_new(n);
3024 is = cpl_bivector_get_x(is_should);
3025 should = cpl_bivector_get_y(is_should);
3026 for (i=0; i<n; i++) {
3027 double x;
3028 int j;
3029 x = cpl_vector_get(li, i);
3030 j = cpl_vector_find(li_true, x);
3031 cpl_vector_set(is, i, x);
3032 cpl_vector_set(should, i, cpl_vector_get(li_true, j));
3033 }
3034 return is_should;
3035}
3036
3037/*----------------------------------------------------------------------------*/
3045/*----------------------------------------------------------------------------*/
3046static double cr2res_wave_etalon_get_x0(
3047 cpl_vector * li,
3048 double trueD)
3049{
3050 cpl_vector * xs;
3051 double x0;
3052 int i, n;
3053
3054 n = cpl_vector_get_size(li);
3055 xs = cpl_vector_new(n);
3056
3057 for (i = 0; i < n; i++) {
3058 cpl_vector_set(xs, i, cpl_vector_get(li, i) - (i * trueD) );
3059 }
3060
3061 x0 = cpl_vector_get_median(xs);
3062
3063 cpl_vector_delete(xs);
3064 return x0;
3065}
3066
3067/*----------------------------------------------------------------------------*/
3073/*----------------------------------------------------------------------------*/
3074static double cr2res_wave_etalon_get_D(
3075 cpl_vector * li)
3076{
3077 int i;
3078 cpl_size nxi;
3079 double trueD=-1.0;
3080 cpl_vector * diffs;
3081
3082 nxi = cpl_vector_get_size(li);
3083 diffs = cpl_vector_new(nxi-1);
3084 for (i=1; i<nxi; i++){
3085 cpl_vector_set(diffs,i-1,
3086 cpl_vector_get(li,i) - cpl_vector_get(li,i-1) );
3087 }
3088
3089 if (cpl_msg_get_level() == CPL_MSG_DEBUG){
3090 cpl_table * tab;
3091 tab = cpl_table_new(nxi-1);
3092 cpl_table_new_column(tab, "wavediff", CPL_TYPE_DOUBLE) ;
3093 for(i=0; i<nxi-1; i++) {
3094 cpl_table_set_double(tab, "wavediff", i, cpl_vector_get(diffs, i));
3095 }
3096
3097 if ( cpl_table_save(tab, NULL, NULL, "debug_wavediffs.fits",
3098 CPL_IO_EXTEND) == CPL_ERROR_FILE_NOT_FOUND) {
3099 cpl_error_reset() ;
3100 cpl_table_save(tab, NULL, NULL, "debug_wavediffs.fits",
3101 CPL_IO_CREATE);
3102 }
3103 cpl_table_delete(tab);
3104 }
3105
3106 trueD = cpl_vector_get_median(diffs);
3107 cpl_vector_delete(diffs);
3108 return trueD;
3109}
3110
char * cr2res_dfs_WAVELENGTH_colname(int order_idx, int trace)
Get the WAVELENGTH column name for a given order/trace.
Definition: cr2res_dfs.c:396
cpl_table * cr2res_dfs_create_lines_diagnostics_table(int nrows)
Create an empty LINES DIAGNOSTICS table.
Definition: cr2res_dfs.c:553
cpl_polynomial * cr2res_etalon_wave_2d(cpl_bivector **spectra, cpl_bivector **spectra_err, cpl_polynomial **wavesol_init, cpl_array **wavesol_init_err, int *orders, int *traces_nb, int ninputs, cpl_size degree_x, cpl_size degree_y, int zp_order, int display, cpl_array **wavelength_error, cpl_table **line_diagnostics)
Create the 2d wavecal fit using etalon peaks.
int cr2res_extract_EXTRACT1D_get_spectrum(const cpl_table *tab, int order, int trace_nb, cpl_bivector **spec, cpl_bivector **spec_err)
Get a Spectrum and its error from the EXTRACT_1D table.
int cr2res_io_get_ext_idx(const char *filename, int detector, int data)
Get the wished extension number for a detector.
Definition: cr2res_io.c:644
cpl_bivector * cr2res_io_load_EMISSION_LINES(const char *filename)
Load an EMISSION_LINES bivector.
Definition: cr2res_io.c:903
double cr2res_pfits_get_wstrt(const cpl_propertylist *plist, int order_idx)
find out the Start wavelength for an order_idx (current detector)
Definition: cr2res_pfits.c:254
double cr2res_pfits_get_wend(const cpl_propertylist *plist, int order_idx)
find out the End wavelength for an order_idx (current detector)
Definition: cr2res_pfits.c:285
double cr2res_qc_wave_resol_fwhm(const cpl_bivector *spec, double *wl)
Computes the lines Fwhm and return the smallest.
Definition: cr2res_qc.c:710
double cr2res_qc_wave_lamp_effic(const cpl_bivector *spec)
Computes the lamp efficiency.
Definition: cr2res_qc.c:656
double cr2res_qc_wave_central(const cpl_table *tw, int order_idx)
Computes the central WLEN of a given order.
Definition: cr2res_qc.c:452
double cr2res_qc_wave_disp(const cpl_table *tw, int order_idx)
Computes the dispersion of a given order.
Definition: cr2res_qc.c:480
cpl_size cr2res_get_trace_table_index(const cpl_table *trace_wave, int order_idx, int trace_nb)
Get the index in a TRACE_WAVE table.
Definition: cr2res_trace.c:652
cpl_polynomial * cr2res_get_trace_wave_poly(const cpl_table *trace_wave, const char *poly_column, int order_idx, int trace_nb)
Get a polynomial from a TRACE_WAVE table.
Definition: cr2res_trace.c:685
cpl_polynomial * cr2res_convert_array_to_poly(const cpl_array *arr)
Convert an array to polynomial.
cpl_vector * cr2res_threshold_spec(const cpl_vector *invector, int smooth, double thresh)
Find the regions with over-average values in a vector.
cpl_vector * cr2res_polynomial_eval_vector(const cpl_polynomial *poly, const cpl_vector *vec)
Evaluate a polynomial on a vector.
cpl_array * cr2res_convert_poly_to_array(const cpl_polynomial *poly, int size)
Convert a polynomial to array.
int cr2res_plot_wavecal_result(const cpl_bivector *extracted_spec, const cpl_bivector *catalog, const char *title, double wmin, double wmax)
Plot the spectrum with the catalog.
int cr2res_order_real_to_idx(int order_real, int order_zp)
Convert the order_real into order_idx.
Definition: cr2res_utils.c:102
hdrl_image * cr2res_wave_gen_wave_map(const cpl_table *trace_wave)
Compute the wavelength map from the trace_wave table.
Definition: cr2res_wave.c:1608
cpl_polynomial * cr2res_wave_polys_1d_to_2d(cpl_polynomial **poly_1ds, int *orders, int npolys, cpl_size xdegree)
Create a 2D Wavelength Polynomial out of a several 1D ones.
Definition: cr2res_wave.c:1685
cpl_polynomial * cr2res_wave_poly_2d_to_1d(cpl_polynomial *poly_2d, int order)
Create a 1D Wavelength Polynomial out of a 2D one.
Definition: cr2res_wave.c:1776
cpl_polynomial * cr2res_wave_2d(cpl_bivector **spectra, cpl_bivector **spectra_err, cpl_polynomial **wavesol_init, cpl_array **wavesol_init_err, int *orders, int *traces_nb, int ninputs, const char *catalog, cpl_size degree, cpl_size degree_y, double threshold, int n_iterations, int zp_order, int display, cpl_array **wavelength_error, cpl_table **lines_diagnostics)
Compute the 2D wavelength polynomial based on a line spectrum and a reference catalog by finding line...
Definition: cr2res_wave.c:889
int cr2res_wave_apply(cpl_table *tw_in, cpl_table *spectra_tab, const cpl_frame *catalog_frame, int reduce_order, int reduce_trace, cr2res_wavecal_type wavecal_type, int degree, int xdegree, double wl_start, double wl_end, double wl_err, double wl_shift, int log_flag, int fallback_input_wavecal_flag, int keep_higher_degrees_flag, int clean_spectrum, int display, double display_wmin, double display_wmax, int zp_order, int grat1_order, cpl_propertylist **qcs, cpl_table **lines_diagnostics, cpl_table **extracted_out, cpl_table **trace_wave_out)
Apply the Wavelength Calibration.
Definition: cr2res_wave.c:150
cpl_polynomial * cr2res_wave_1d(cpl_bivector *spectrum, cpl_bivector *spectrum_err, cpl_polynomial *wavesol_init, const cpl_array *wave_error_init, int order, int trace_nb, cr2res_wavecal_type wavecal_type, const char *catalog, int degree, int clean_spectrum, int log_flag, int keep_higher_degrees_flag, int display, double display_wmin, double display_wmax, double *best_xcorr, double *lines_resol_fwhm, double *lines_pos, double *lines_resol, double *lines_intens, cpl_array **wavelength_error, cpl_table **lines_diagnostics)
1D Wavelength Calibration
Definition: cr2res_wave.c:689
cpl_array * cr2res_wave_get_estimate(const char *filename, int detector, int order)
Compute the wavelength estimate.
Definition: cr2res_wave.c:1527
cpl_polynomial * cr2res_wave_estimate_compute(double wmin, double wmax)
Compute the polynomial from boundaries.
Definition: cr2res_wave.c:1481
cpl_polynomial * cr2res_wave_etalon(cpl_bivector *spectrum, const cpl_bivector *spectrum_err, cpl_polynomial *wavesol_init, int degree, cpl_array **wavelength_error)
Find solution from etalon.
Definition: cr2res_wave.c:1405
cpl_polynomial * cr2res_wave_xcorr(cpl_bivector *spectrum, cpl_polynomial *wavesol_init, double wl_error, cpl_bivector *lines_list, int degree, int keep_higher_degrees_flag, int cleaning_filter_size, double slit_width, double fwhm, int display, double *best_xcorr, cpl_array **wavelength_error)
Find solution by cross-correlating template spectrum.
Definition: cr2res_wave.c:1221
cpl_image * hdrl_image_get_image(hdrl_image *himg)
get data as cpl image
Definition: hdrl_image.c:105
hdrl_image * hdrl_image_new(cpl_size nx, cpl_size ny)
create new zero filled hdrl image
Definition: hdrl_image.c:311