X-shooter Pipeline Reference Manual 3.8.15
xsh_eqwidth_lib.c
Go to the documentation of this file.
1
2#ifdef HAVE_CONFIG_H
3#include <config.h>
4#endif
5
6/*-----------------------------------------------------------------------------
7 Includes
8 -----------------------------------------------------------------------------*/
9
10#include "xsh_eqwidth_lib.h"
11#include <xsh_msg.h>
12#include <gsl/gsl_version.h>
13
16/*----------------------------------------------------------------------------*/
22/*----------------------------------------------------------------------------*/
23
24
25/*----------------------------------------------------------------------------*/
26/*
27 @brief get index close to a given wavelenght in a table spectrum
28 @param spec_1d spectrum table
29 @param lambda wavelegth
30 @return return index (-1) if not valid????
31
32
33 This function returns a the index of the closest wavelenght greatest than
34 lambda
35-------------------------------------------------------------------------------*/
36static cpl_size get_index_from_spec(cpl_table* spec_total, double lambda){
37 //int* error=NULL;
38 cpl_size index;
39
40 /*
41 double cdelta1=cpl_table_get(spec_total, "WAVEL",1, error) -
42 cpl_table_get(spec_total, "WAVEL",0, error);
43 double crval1=cpl_table_get(spec_total, "WAVEL",0, error);
44 */
45 // Obsolete
46 //index = (long) (1./cdelta1*lambda-crval1/cdelta1) + 1;
47 // Correct formulae for velocity-binned spectra (TO BE USED FOR EFFICIENCY!?)
48 //double dv = cdelta1 / crval1;
49 //index = (long)ceil(log(lambda / crval1) / dv);
50
51 // Alternative instructions, using CPL functions
52 cpl_table_unselect_all(spec_total);
53 cpl_table_or_selected_double(spec_total, "WAVEL", CPL_NOT_GREATER_THAN,
54 lambda);
55 index = (int)cpl_table_count_selected(spec_total);
56
57 if (index <= cpl_table_get_nrow(spec_total))
58 { return index; }
59 else { return -1;}
60}
61
62
63#if 0
64/* NOT USED ! */
65static cpl_size get_index_from_spec_alt(cpl_table* spec_total, double lambda){
66 int* error=NULL;
67 cpl_size index;
68 double cdelta1, crval1;
69
70 cdelta1=cpl_table_get(spec_total, "WAVEL",1, error) -
71 cpl_table_get(spec_total, "WAVEL",0, error);
72
73 crval1=cpl_table_get(spec_total, "WAVEL",0, error);
74 index = (long) (1./cdelta1*lambda-crval1/cdelta1) + 1;
75
76// if spectra is binned in velocity
77 double dv = (cdelta1)/crval1;
78 double r = 1.+dv;
79 index = (long) log(lambda*r/cdelta1)/log(r);
80
81// Python: math.log(ls*(1+dv/c)/ll[0])/math.log(r)
82
83 if (index <= cpl_table_get_nrow(spec_total))
84 { return index; }
85 else { return -1;}
86}
87#endif
88
89
90
91/*----------------------------------------------------------------------------*/
102/*----------------------------------------------------------------------------*/
103
104static cpl_size get_pixel_to_nm_scale(cpl_table* spec_total, double nm){
105 int* error=NULL;
106 double cdelta1= cpl_table_get(spec_total, "WAVEL",1, error) -
107 cpl_table_get(spec_total, "WAVEL",0, error);
108 return nm/cdelta1;
109}
110
111/*----------------------------------------------------------------------------*/
123/*----------------------------------------------------------------------------*/
124
125
126cpl_error_code select_local_spec(cpl_table* spec_total,
127 double ew_space,
128 double lambda,
129 cpl_table** spec_region)
130{
131 cpl_errorstate prestate = cpl_errorstate_get();
132 //int* error = NULL;
133
134
135
136 cpl_size index_line = get_index_from_spec(spec_total, lambda);
137
138// 2* space is just a definition...
139 cpl_size count = get_pixel_to_nm_scale(spec_total, 2 * ew_space);
140 cpl_size begin = index_line - count/2;
141
142
143 *spec_region = cpl_table_extract(spec_total,
144 begin,
145 count);
146
147
148 if (!cpl_errorstate_is_equal(prestate)) {
149 return (int)cpl_error_set_message(cpl_func, cpl_error_get_code(),
150 "Unable to Get region of the spectrum");
151 }
152
153
154 return CPL_ERROR_NONE;
155}
156
157
158
159
160
161
162void find_left_right_continuum_pos(int* x1, int* x2, cpl_table* spec_region,double cont_rejt, double line){
163
164
165
166// Interface from CPL tables to C arrays
167 int i,n = cpl_table_get_nrow(spec_region);
168 double x[n], y[n];
169 for(i=0; i < n; i++){
170 x[i] = cpl_table_get_double(spec_region,"WAVEL",i,NULL);
171 y[i] = cpl_table_get_double(spec_region,"FLUX",i,NULL);
172 }
173
174
175// Code from ARES (could be improved
176// It starts from the borders to the inside and check the continuum closest
177// location to the line (Confuse but it works)
178 int xind1=0,xind2=n-1,hjk;
179 double klo=0.01; //Nm instead of Angstroms
180 for (hjk=0; hjk < n; hjk++){
181 if ((y[hjk] > cont_rejt) &&
182 (x[hjk] - (line-klo) > x[xind1] - (line-klo)) &&
183 (x[hjk] - (line-klo) < 0))
184 xind1=hjk;
185 if ((y[hjk] > cont_rejt) &&
186 (x[hjk] - (line+klo) < x[xind2] - (line+klo)) &&
187 (x[hjk] - (line+klo) > 0) )
188 xind2=hjk;
189 }
190 *x1=xind1;
191 *x2=xind2;
192
193}
194
195static cpl_table* esp_spec_deriv(cpl_table* spec_region){
196
197// Interface from CPL tables to C arrays
198 int i,n = cpl_table_get_nrow(spec_region);
199 double x[n], y[n], dy[n];
200
201 for(i=0; i < n; i++){
202 x[i] = cpl_table_get_double(spec_region,"WAVEL",i,NULL);
203 y[i] = cpl_table_get_double(spec_region,"FLUX",i,NULL);
204 }
205
206 cpl_table* table_deriv = NULL;
207 table_deriv = cpl_table_duplicate(spec_region);
208
209 deriv(x,y,dy,n);
210
211// Interface from c arrays to CPL Tables:
212// Updating the spectral region
213 for(i=0; i < n; i++){
214 cpl_table_set_double(table_deriv,"FLUX",i,dy[i]);
215 }
216
217 return table_deriv;
218}
219
220
221static void esp_spec_smooth(cpl_table* spec_region, int smwidth){
222
223// Interface from CPL tables to C arrays
224 int i,n = cpl_table_get_nrow(spec_region);
225 double y[n], sy[n];
226 for(i=0; i < n; i++){
227 y[i] = cpl_table_get_double(spec_region,"FLUX",i,NULL);
228 }
229
230 smooth(y, n, smwidth, sy);
231
232// Interface from c arrays to CPL Tables:
233// Updating the spectral region
234 for(i=0; i < n; i++){
235 cpl_table_set_double(spec_region,"FLUX",i,sy[i]);
236 }
237
238}
239
240
241
242
243void smooth(double vec[], long n, int w, double svec[]) {
244 int i,j;
245 double soma;
246
247 if (w%2 != 1)
248 w++;
249 for (i=0; (i < (w-1)/2);i++)
250 svec[i]=vec[i];
251 for (i=(w-1)/2;i<n-((w-1)/2);i++) {
252 soma=0.;
253 for (j=i-((w-1)/2); j<=i+((w-1)/2);j++)
254 soma+=vec[j];
255 svec[i]=soma/w;
256 }
257 for (i=n-((w-1)/2); i<n;i++)
258 svec[i]=vec[i];
259}
260
261
262
263static void zeroscenterfind(double iy[], double y[], double dy[], double ddy[], long n, long center[], long *ncenter,double det_line_thres) {
264 double maxdy;
265 long ctot=0, i, centertot[n];
266 int signalc=0, signalc_ant=0;
267 if (ddy[0] == abs(ddy[0]))
268 signalc=1;
269 signalc_ant=signalc;
270 maxdy=maxele_vec(dy,n);
271
272 for (i=0; i<n; i++) {
273 signalc=0;
274 if ( (float) ddy[i] == fabs( (float) ddy[i]) )
275 signalc=1;
276 // EN: when the signal changes, the local maximum of the 2nd derivative
277 // is below the the noise, and the 3rd derivative is negative enough (due
278 // to noise)
279 // EN: for the 0.98 the idea was to have the tree/rejt value, but it does not
280 // work so well. It identifies to many lines in cases of good S/N. This way
281 // we only accept identified lines that have a depth of at least 0.98.
282 // EN : The 3rd condition was hard codded as iy[i] < 0.98
283 if ((signalc != signalc_ant) &&
284 (dy[i] > 0.01*maxdy) &&
285 (iy[i] < 1. - det_line_thres) &&
286 (ddy[i] < -0.1)) {
287 centertot[ctot]=i;
288 ctot++;
289 }
290 signalc_ant=signalc;
291 }
292
293
294 if (ctot != 0) {
295 *ncenter=ctot;
296 for (i=0;i<ctot;i++) center[i]=centertot[i];
297 } else {
298 center[0]=-1;
299 *ncenter=0;
300 }
301}
302
303
304double maxele_vec(double vec[], long nvec) {
305 long i;
306 double maxi=vec[1+nvec/20];
307 for (i=1+nvec/20; i<nvec-nvec/20; i++)
308 maxi = max(maxi,vec[i]);
309 return maxi;
310}
311
312
313// NEEDS GSL:::
314
315
316/*----------------------------------------------------------------------------*/
330/*----------------------------------------------------------------------------*/
331
332cpl_error_code esp_fit_ngauss(cpl_table* spec_cont_region,
333 cpl_table* line_table,
334 double fit_ngauss_width)
335{
336
337 cpl_errorstate prestate = cpl_errorstate_get();
338
339
340// Interface from CPL tables to C arrays
341 int i,ns = cpl_table_get_nrow(spec_cont_region);
342 int nl = cpl_table_get_nrow(line_table);
343
344 double xfit[ns], yfit[ns], sigma[ns];
345
346 for(i=0; i < ns; i++){
347 xfit[i] = cpl_table_get_double(spec_cont_region,"WAVEL",i,NULL);
348 yfit[i] = cpl_table_get_double(spec_cont_region,"FLUX",i,NULL) - 1.0;
349// sigma[i] = cpl_table_get_double(spec_cont_region,"FLUXERR",i,NULL);
350// One good approximation for testing purposes is 1.-cont_rejt
351 sigma[i] = 1.- 0.997573;
352 }
353
354 int npara=3*nl;
355 int igauss=0;
356 double acoef[npara], acoef_er[npara];
357 for (i=0;i<nl;i++) {
358 acoef[3*igauss]=cpl_table_get_double(line_table,"PEAK",i,NULL) - 1.0;
359// CAREFUL here with this guess value
360 acoef[3*igauss+1]=fit_ngauss_width;
361 acoef[3*igauss+2]=cpl_table_get_double(line_table,"WAVEL",i,NULL);
362 igauss++;
363 }
364
365 int status2;
366
367 fitngauss(xfit,yfit,sigma,ns,acoef,acoef_er,npara,&status2);
368
369// We probably need to control here the sucess of the fitting...
370
371// Interface from c arrays to CPL Tables:
372// Updating the spectral region
373 for(i=0; i < nl; i++){
374 cpl_table_set_double(line_table,"PEAK",i,-acoef[3*i]);
375 cpl_table_set_double(line_table,"PEAK_ERR",i,acoef_er[3*i]);
376 cpl_table_set_double(line_table,"WAVEL",i,acoef[3*i+2]);
377 cpl_table_set_double(line_table,"WAVEL_ERR",i,acoef_er[3*i+2]);
378// FWHM for the defined Gaussian (ARES): F(X)=Aexp(-Lambda(x-c)^2) => FWHM=2*sqrt(ln(2)/lambda)
379 cpl_table_set_double(line_table,"FWHM",i,2.*sqrt(log(2)/acoef[3*i+1]));
380 double ind_ew = acoef[3*i]*sqrt(CPL_MATH_PI/acoef[3*i+1])*-1.;
381 cpl_table_set_double(line_table,"EW",i,ind_ew*10000.);
382 cpl_table_set_double(line_table,"SIGMA",i,sqrt(0.5/(acoef[3*i+1])));
383 cpl_table_set_double(line_table,"SIGMA_ERR",i,0.5*CPL_MATH_SQRT1_2*(acoef_er[3*i+1])*pow(acoef[3*i+1],-3.*0.5));
384// using error equation propagation
385// medida_er_square+=medida*medida * ( acoef_er[3*hj]*acoef_er[3*hj]/acoef[3*hj]/acoef[3*hj] + (0.5*0.5*acoef_er[3*hj+1]*acoef_er[3*hj+1]/acoef[3*hj+1]/acoef[3*hj+1]));
386// printf("medida: %15.10lf\n", ind_ew);
387 double ind_ew_err_square = ind_ew*ind_ew * ( acoef_er[3*i]*acoef_er[3*i]/acoef[3*i]/acoef[3*i] + (0.5*0.5*acoef_er[3*i+1]*acoef_er[3*i+1]/acoef[3*i+1]/acoef[3*i+1]));
388 cpl_table_set_double(line_table,"EW_ERR",i,sqrt(ind_ew_err_square)*10000.);
389 }
390
391
392 if (!cpl_errorstate_is_equal(prestate)) {
393 return (int)cpl_error_set_message(cpl_func, cpl_error_get_code(),
394 "Unable to Get region of the spectrum");
395 }
396
397
398 return CPL_ERROR_NONE;
399
400}
401
402
403double check_ew(cpl_table* line_table,
404 double line,
405 double det_line_resol,
406 int* index_line,
407 int* n_lines,
408 double* ew_error){
409
410 int nl = cpl_table_get_nrow(line_table);
411
412 double ew=0.;
413 double ew_er=0.;
414 int i;
415
416 *index_line = 0;
417 *n_lines=0;
418
419 for (i=0; i < nl; i++) {
420 if ( fabs( line - cpl_table_get_double(line_table,"WAVEL",i,NULL) ) < det_line_resol ) {
421 ew+=cpl_table_get_double(line_table,"EW",i,NULL);
422 ew_er+=cpl_table_get_double(line_table,"EW_ERR",i,NULL);
423 (*n_lines)++;
424 *index_line=i;
425 }
426 }
427 *ew_error=ew_er;
428
429 return ew;
430}
431
432
433
434static void
435poly_fitn(double xvec[], double yvec[], double err[], long n, long ord, double coefs[]) {
436 int i, j, k;
437 double xi, yi, ei, chisq,xi2;
438 gsl_matrix *X, *cov;
439 gsl_vector *y, *w, *c;
440 ord++;
441 X = gsl_matrix_alloc (n, ord);
442 y = gsl_vector_alloc (n);
443 w = gsl_vector_alloc (n);
444 c = gsl_vector_alloc (ord);
445 cov = gsl_matrix_alloc (ord, ord);
446 for (i = 0; i < n; i++) {
447 xi=xvec[i];
448 yi=yvec[i];
449 ei=err[i];
450 for (j = 0; j < ord; j++) {
451 xi2=1.0;
452 for (k=0; k<j; k++) xi2*=xi;
453 gsl_matrix_set (X, i, j, xi2);
454 }
455 gsl_vector_set (y, i, yi);
456 gsl_vector_set (w, i, 1.0/(ei*ei));
457 }
458
459 gsl_multifit_linear_workspace * work = gsl_multifit_linear_alloc (n, ord);
460 gsl_multifit_wlinear (X, w, y, c, cov, &chisq, work);
461 gsl_multifit_linear_free (work);
462
463 #define C(i) (gsl_vector_get(c,(i)))
464 #define COV(i,j) (gsl_matrix_get(cov,(i),(j)))
465
466 for (j = 0; j < ord; j++)
467 coefs[j]=C(j);
468 gsl_vector_free (y);
469 gsl_vector_free (w);
470 gsl_vector_free (c);
471 gsl_matrix_free (X);
472 gsl_matrix_free (cov);
473
474}
475
476void deriv(double x[], double y[], double dy[], long n) {
477 int i;
478 gsl_interp_accel *acc = gsl_interp_accel_alloc ();
479 gsl_interp *interp = gsl_interp_alloc (gsl_interp_cspline, n);
480 //gsl_interp *interp = gsl_interp_alloc (gsl_interp_akima, n);
481
482 gsl_interp_init (interp, x, y, n);
483 for (i=0; i<n; i++)
484 dy[i]=gsl_interp_eval_deriv (interp, x, y,x[i],acc);
485 gsl_interp_free (interp);
486 gsl_interp_accel_free (acc);
487}
488
489
490void fitngauss(double t[], double y[], double sigma[], long nvec,
491 double acoef[], double acoef_er[], int para, int *status2)
492{
493 const gsl_multifit_fdfsolver_type *T;
494 gsl_multifit_fdfsolver *s;
495
496 int status;
497 size_t i, iter = 0;
498 long N=nvec;
499 const size_t n = N;
500 const size_t p = para;
501
502 gsl_matrix *covar = gsl_matrix_alloc (p, p);
503#if defined GSL_MAJOR_VERSION && GSL_MAJOR_VERSION >= 2
504 gsl_matrix *J;
505#endif
506// double dy[N];
507 struct data d = { n, para, t, y, sigma};
508 gsl_multifit_function_fdf f;
509
510 double x_init[para];
511 for (i=0; i< (size_t) para; i++)
512 x_init[i]=acoef[i];
513
514 f.f = &expb_f;
515 f.df = &expb_df;
516 f.fdf = &expb_fdf;
517 f.n = n;
518 f.p = p;
519 f.params = &d;
520
521 gsl_vector_view x = gsl_vector_view_array (x_init, p);
522
523 T = gsl_multifit_fdfsolver_lmder;
524
525 s = gsl_multifit_fdfsolver_alloc (T, n, p);
526
527 gsl_multifit_fdfsolver_set (s, &f, &x.vector);
528 do
529 {
530 iter++;
531 status = gsl_multifit_fdfsolver_iterate (s);
532
533// int i=0;
534 if (status)
535 break;
536
537 status = gsl_multifit_test_delta (s->dx, s->x,
538 1e-6, 1e-6);
539 }
540 while (status == GSL_CONTINUE && iter < 5000);
541#if defined GSL_MAJOR_VERSION && GSL_MAJOR_VERSION >= 2
542 J = gsl_matrix_alloc(n, p);
543 gsl_multifit_fdfsolver_jac(s, J);
544 gsl_multifit_covar (J, 0.0, covar);
545 gsl_matrix_free (J);
546#else
547 gsl_multifit_covar (s->J, 0.0, covar);
548#endif
549
550#define FIT(i) gsl_vector_get(s->x, i)
551#define ERR(i) sqrt(gsl_matrix_get(covar,i,i))
552
553 {
554
555 double chi = gsl_blas_dnrm2(s->f);
556 double dof = n - p;
557 double c = GSL_MAX_DBL(1, chi / sqrt(dof));
558
559 for (i=0; i < (size_t) para; i++)
560 {
561 acoef[i]=FIT(i);
562 acoef_er[i]=c*ERR(i);
563 }
564
565 *status2=status;
566 *status2=0; //sometimes we have bad fits but the result is perfectably acceptable
567 }
568 gsl_multifit_fdfsolver_free (s);
569 gsl_matrix_free (covar);
570}
571
572
573// TO BE Located in another file:
574
575
576/*----------------------------------------------------------------------------*/
586/*----------------------------------------------------------------------------*/
587
588cpl_error_code espda_create_line_table(cpl_table **tab,
589 cpl_size size)
590
591{
592 cpl_errorstate prestate = cpl_errorstate_get();
593 cpl_size ceil;
594
595/*
596 * Table of given row size is created. Columns are those specified for the
597 * structure SPEC.
598 */
599
600 *tab = cpl_table_new(size);
601 cpl_table_new_column(*tab, "WAVEL", CPL_TYPE_DOUBLE);
602 cpl_table_new_column(*tab, "WAVEL_ERR", CPL_TYPE_DOUBLE);
603 cpl_table_new_column(*tab, "PEAK", CPL_TYPE_DOUBLE);
604 cpl_table_new_column(*tab, "PEAK_ERR", CPL_TYPE_DOUBLE);
605 cpl_table_new_column(*tab, "MU", CPL_TYPE_DOUBLE);
606 cpl_table_new_column(*tab, "MU_ERR", CPL_TYPE_DOUBLE);
607 cpl_table_new_column(*tab, "SIGMA", CPL_TYPE_DOUBLE);
608 cpl_table_new_column(*tab, "SIGMA_ERR", CPL_TYPE_DOUBLE);
609 cpl_table_new_column(*tab, "EW", CPL_TYPE_DOUBLE);
610 cpl_table_new_column(*tab, "EW_ERR", CPL_TYPE_DOUBLE);
611 cpl_table_new_column(*tab, "FWHM", CPL_TYPE_DOUBLE);
612 cpl_table_new_column(*tab, "FWHM_ERR", CPL_TYPE_DOUBLE);
613
614 if (!cpl_errorstate_is_equal(prestate)) {
615 return (int)cpl_error_set_message(cpl_func, cpl_error_get_code(),
616 "Unable to create table.");
617 }
618
619/*
620 * Columns are initialized with default values.
621 */
622
623 if (size > 0) {
624 ceil = size;
625 }
626 else {
627 ceil = 0;
628 }
629 cpl_table_fill_column_window_double(*tab, "WAVEL", 0, ceil, -9999.0);
630 cpl_table_fill_column_window_double(*tab, "WAVEL_ERR", 0, ceil, -9999.0);
631 cpl_table_fill_column_window_double(*tab, "PEAK", 0, ceil, -9999.0);
632 cpl_table_fill_column_window_double(*tab, "PEAK_ERR", 0, ceil, -9999.0);
633 cpl_table_fill_column_window_double(*tab, "MU", 0, ceil, -9999.0);
634 cpl_table_fill_column_window_double(*tab, "MU_ERR", 0, ceil, -9999.0);
635 cpl_table_fill_column_window_double(*tab, "SIGMA", 0, ceil, -9999.0);
636 cpl_table_fill_column_window_double(*tab, "SIGMA_ERR", 0, ceil, -9999.0);
637 cpl_table_fill_column_window_double(*tab, "EW", 0, ceil, -9999.0);
638 cpl_table_fill_column_window_double(*tab, "EW_ERR", 0, ceil, -9999.0);
639 cpl_table_fill_column_window_double(*tab, "FWHM", 0, ceil, -9999.0);
640 cpl_table_fill_column_window_double(*tab, "FWHM_ERR", 0, ceil, -9999.0);
641
642 if (!cpl_errorstate_is_equal(prestate)) {
643 return (int)cpl_error_set_message(cpl_func, cpl_error_get_code(),
644 "Unable to initialize table.");
645 }
646
647/*
648 * The new table is passed along.
649 */
650
651 return CPL_ERROR_NONE;
652}
653
654/*----------------------------------------------------------------------------*/
668/*----------------------------------------------------------------------------*/
669
670
671cpl_error_code esp_fit_lcont(cpl_table* spec_region,
672 double cont_rejt,
673 int cont_iter){
674
675 cpl_errorstate prestate = cpl_errorstate_get();
676
677 int n = cpl_table_get_nrow(spec_region);
678
679
680 int order,i,j;
681 order=2;
682 double coefs[order+1];
683 double x[n], y[n], err[n], ynorm[n];
684 long nvec;
685
686
687// Interface from CPL tables to C arrays
688 for(i=0; i < n; i++){
689 x[i] = cpl_table_get_double(spec_region,"WAVEL",i,NULL);
690 y[i] = cpl_table_get_double(spec_region,"FLUX",i,NULL);
691// err[i] = cpl_table_get_double(spec_region,"FLUXERR",i,NULL);
692 err[i] = 1.0;
693 }
694
695
696// Calling ARES routines and code (continuum_det5)
697
698 poly_fitn(x,y,err,n,order,coefs);
699
700// Creation of the array with the values of the polynomial
701 double xi=1.;
702 for(i=0; i < n; i++) {
703 ynorm[i]=0.;
704 xi=1.;
705 for (j=0; j < order + 1; j++) {
706 ynorm[i]+=coefs[j]*xi;
707 xi*=x[i];
708 }
709 }
710
711 double vecx[n],vecy[n];
712 int jk;
713 for (jk=0; jk < cont_iter; jk++) {
714 nvec=0;
715// Selection of the points for the next fit
716 for (i=0; i<n-1; i++) {
717//This was to try to avoid cosmic rays, although does not works well always.
718// original comment: test were made with 0.01, there should not be any problem to put it larger. I choosed 0.1 in ARES
719
720 if (y[i] > ynorm[i]*cont_rejt && fabs(y[i]-y[i+1]) < 0.1*y[i]){
721 vecx[nvec]=x[i];
722 vecy[nvec]=y[i];
723 nvec++;
724 }
725 }
726 poly_fitn(vecx,vecy,err,nvec,order,coefs);
727
728 for(i=0; i < n; i++) {
729 ynorm[i]=0.;
730 xi=1.;
731 for (j=0;j<order+1;j++) {
732 ynorm[i]+=coefs[j]*xi;
733 xi*=x[i];
734 }
735 }
736 }
737
738// normalization
739 for (i=0; i < n; i++)
740 ynorm[i]=y[i]/(coefs[0]+coefs[1]*x[i]+coefs[2]*x[i]*x[i]);
741
742// Interface from c arrays to CPL Tables:
743// Updating the spectral region
744 for(i=0; i < n; i++){
745 cpl_table_set_double(spec_region,"FLUX",i,ynorm[i]);
746 }
747
748 if (!cpl_errorstate_is_equal(prestate)) {
749 return (int)cpl_error_set_message(cpl_func, cpl_error_get_code(),
750 "Unable to Get region of the spectrum");
751 }
752
753
754 return CPL_ERROR_NONE;
755}
756
757
758
759/*----------------------------------------------------------------------------*/
773/*----------------------------------------------------------------------------*/
774
775cpl_error_code esp_det_line(cpl_table* spec_cont_region,
776 double det_line_thres,
777 double det_line_resol,
778 int det_line_smwidth,
779 cpl_table** line_table
780 ){
781
782 cpl_errorstate prestate = cpl_errorstate_get();
783
784
785 cpl_table* spec_cont_region_der1 = NULL;
786 cpl_table* spec_cont_region_der2 = NULL;
787 cpl_table* spec_cont_region_der3 = NULL;
788
789 spec_cont_region_der1 = esp_spec_deriv(spec_cont_region);
790
791 esp_spec_smooth(spec_cont_region_der1, det_line_smwidth);
792
793 spec_cont_region_der2 = esp_spec_deriv(spec_cont_region_der1);
794 esp_spec_smooth(spec_cont_region_der2, det_line_smwidth);
795
796 spec_cont_region_der3 = esp_spec_deriv(spec_cont_region_der2);
797 esp_spec_smooth(spec_cont_region_der3, det_line_smwidth);
798
799
800// Interface from CPL tables to C arrays
801 int i,n = cpl_table_get_nrow(spec_cont_region);
802 double x[n], iy[n], y[n], dy[n], ddy[n];
803 for(i=0; i < n; i++){
804 x[i] = cpl_table_get_double(spec_cont_region,"WAVEL",i,NULL);
805 iy[i] = cpl_table_get_double(spec_cont_region,"FLUX",i,NULL);
806 y[i] = cpl_table_get_double(spec_cont_region_der1,"FLUX",i,NULL);
807 dy[i] = cpl_table_get_double(spec_cont_region_der2,"FLUX",i,NULL);
808 ddy[i] = cpl_table_get_double(spec_cont_region_der3,"FLUX",i,NULL);
809 }
810
811 cpl_table_delete(spec_cont_region_der1);
812 cpl_table_delete(spec_cont_region_der2);
813 cpl_table_delete(spec_cont_region_der3);
814
815 long ncenter=n, center[n];
816
817// this obtains the indexes close to the lines center
818 zeroscenterfind(iy, y, dy, ddy, n, center, &ncenter,det_line_thres);
819
820// calculating/interpolating the position of the center of the lines in the spectra:
821// Maybe this is not necessary...
822
823// Only if lines were detected
824 if (center[0] != -1 && ncenter != 0) {
825 double xlinhas[ncenter], ylinhas[ncenter];
826 int i1,i1m;
827 double den, diff_y;
828
829 for (i=0; i<ncenter; i++) {
830 i1=center[i];
831 i1m=i1-1;
832 den=1./(x[i1]-x[i1-1]);
833 diff_y=(ddy[i1] - ddy[i1m]);
834 xlinhas[i]= ( -ddy[i1m] + diff_y*den * x[i1] ) / ( diff_y*den );
835 ylinhas[i]= diff_y*den * xlinhas[i] + iy[i1m] - ( iy[i1]- iy[i1m] )*den * x[i1] ;
836 }
837
838
839//RESAMPLING, Negleting "lines" that are to close to each other (noisy lines)
840 double xvec2[ncenter], yvec2[ncenter];
841 int nvec2,j;
842 xvec2[0]=xlinhas[0];
843 yvec2[0]=ylinhas[0];
844 j=0;
845 for(i=1;i<ncenter;i++) {
846 if (fabs(xvec2[j]-xlinhas[i]) < det_line_resol ) {
847 xvec2[j]=(xvec2[j]+xlinhas[i])*0.5;
848 yvec2[j]=(yvec2[j]+ylinhas[i])*0.5;
849 } else {
850 j++;
851 xvec2[j]=xlinhas[i];
852 yvec2[j]=ylinhas[i];
853 }
854 }
855 nvec2=j+1;
856
857
858// Interface C arrays to CPL Table
859// Creating and filling the line table
860
861 cpl_ensure_code(espda_create_line_table(line_table,
862 nvec2)
863 == CPL_ERROR_NONE, cpl_error_get_code());
864
865 for (i=0; i<nvec2; i++) {
866 cpl_table_set_double(*line_table,"WAVEL",i,xvec2[i]);
867 cpl_table_set_double(*line_table,"PEAK",i,yvec2[i]);
868 }
869// If no lines are detected then we create a line table without rows
870 } else {
871 int nvec2=0;
872 cpl_ensure_code(espda_create_line_table(line_table,
873 nvec2)
874 == CPL_ERROR_NONE, cpl_error_get_code());
875
876 }
877
878
879 if (!cpl_errorstate_is_equal(prestate)) {
880 return (int)cpl_error_set_message(cpl_func, cpl_error_get_code(),
881 "Unable to Get region of the spectrum");
882 }
883
884
885 return CPL_ERROR_NONE;
886}
887
888
889
890
891
892
893
894
#define N
static double sigma
int size
int * y
int * x
static SimAnneal s
Definition: xsh_model_sa.c:99
double * x
Definition: xsh_model_sa.c:94
double * t
Definition: xsh_ngaussfdf.h:13
int para
Definition: xsh_ngaussfdf.h:12
int n
Definition: xsh_detmon_lg.c:92
int order
Definition: xsh_detmon_lg.c:80
cpl_error_code esp_det_line(cpl_table *spec_cont_region, double det_line_thres, double det_line_resol, int det_line_smwidth, cpl_table **line_table)
Performs a local normalization of the spectrum region.
double maxele_vec(double vec[], long nvec)
cpl_error_code esp_fit_ngauss(cpl_table *spec_cont_region, cpl_table *line_table, double fit_ngauss_width)
Fit a list of absorption lines with n Gaussian profiles and compute the equivalent width of each line...
void smooth(double vec[], long n, int w, double svec[])
void find_left_right_continuum_pos(int *x1, int *x2, cpl_table *spec_region, double cont_rejt, double line)
static void esp_spec_smooth(cpl_table *spec_region, int smwidth)
#define C(i)
#define FIT(i)
void fitngauss(double t[], double y[], double sigma[], long nvec, double acoef[], double acoef_er[], int para, int *status2)
double check_ew(cpl_table *line_table, double line, double det_line_resol, int *index_line, int *n_lines, double *ew_error)
static void zeroscenterfind(double iy[], double y[], double dy[], double ddy[], long n, long center[], long *ncenter, double det_line_thres)
#define ERR(i)
cpl_error_code espda_create_line_table(cpl_table **tab, cpl_size size)
Create a LINE table with a given row size.
static cpl_size get_pixel_to_nm_scale(cpl_table *spec_total, double nm)
get scale from nm to pixels
static void poly_fitn(double xvec[], double yvec[], double err[], long n, long ord, double coefs[])
static cpl_size get_index_from_spec(cpl_table *spec_total, double lambda)
cpl_error_code esp_fit_lcont(cpl_table *spec_region, double cont_rejt, int cont_iter)
Performs a local normalization of the spectrum region.
void deriv(double x[], double y[], double dy[], long n)
static cpl_table * esp_spec_deriv(cpl_table *spec_region)
cpl_error_code select_local_spec(cpl_table *spec_total, double ew_space, double lambda, cpl_table **spec_region)
Select_local_spec.
#define max(a, b)
DOUBLE vec[4]
int expb_fdf(const gsl_vector *x, void *data, gsl_vector *f, gsl_matrix *J)
Definition: xsh_ngaussfdf.h:92
int expb_f(const gsl_vector *x, void *data, gsl_vector *f)
Definition: xsh_ngaussfdf.h:19
int expb_df(const gsl_vector *x, void *data, gsl_matrix *J)
Definition: xsh_ngaussfdf.h:53