X-shooter Pipeline Reference Manual 3.8.15
xsh_eqwidth_lib.h
Go to the documentation of this file.
1
2#ifndef XSH_EQWIDTH_LIB_H
3#define XSH_EQWIDTH_LIB_H
4
5/*-----------------------------------------------------------------------------
6 Includes
7 -----------------------------------------------------------------------------*/
8
9#include <cpl.h>
10#include <gsl/gsl_multifit.h>
11#include <gsl/gsl_interp.h>
12#include <gsl/gsl_rng.h>
13#include <gsl/gsl_randist.h>
14#include <gsl/gsl_vector.h>
15#include <gsl/gsl_blas.h>
16#include <gsl/gsl_multifit_nlin.h>
17#include <gsl/gsl_interp.h>
18#include <gsl/gsl_sf_bessel.h>
19#include "xsh_ngaussfdf.h"
20//#include "espda_util_das.h"
21
22/*-----------------------------------------------------------------------------
23 Macros -----------------------------------------------------------------------------*/
24#define max(a,b) (((a)>(b))?(a):(b))
25
26/*-----------------------------------------------------------------------------
27 Functions prototypes
28 -----------------------------------------------------------------------------*/
29cpl_error_code select_local_spec(cpl_table* spec_total,
30 double ew_space,
31 double lambda,
32 cpl_table** spec_region);
33
34//cpl_size get_index_from_spec(cpl_table* spec_total, double lambda);
35//cpl_size get_index_from_spec_alt(cpl_table* spec_total, double lambda);
36
37cpl_size get_nm_to_index_scale(cpl_table* spec_total, double nm);
38
39cpl_error_code esp_fit_lcont(cpl_table* spec_region,
40 double cont_rejt,
41 int cont_iter);
42
43cpl_error_code esp_det_line(cpl_table* spec_region,
44 double det_line_thres,
45 double det_line_resol,
46 int det_line_smwidth,
47 cpl_table** line_table);
48
49
50//cpl_table* esp_spec_deriv(cpl_table* spec_region);
51
52//void esp_spec_smooth(cpl_table* spec_region, int smwidth);
53
54
55void find_left_right_continuum_pos(int* xind1, int* xind2,
56 cpl_table* spec_region,
57 double cont_rejt, double line);
58
59
60void smooth(double vec[], long n, int w, double svec[]);
61
62/*
63void zeroscenterfind(double iy[], double y[],
64 double dy[], double ddy[],
65 long n, long center[], long *ncenter,double det_line_thres);
66*/
67
68double maxele_vec(double vec[], long nvec);
69
70cpl_error_code esp_fit_ngauss(cpl_table* spec_cont_region,
71 cpl_table* line_table,
72 double fit_ngauss_width);
73
74double check_ew(cpl_table* line_table,
75 double line,
76 double det_line_resol,
77 int* index_line,
78 int* n_lines,
79 double* ew_error);
80/*
81void poly_fitn(double xvec[], double yvec[], double err[],
82 long n, long ord, double coefs[]);
83 */
84
85void deriv(double x[], double y[], double dy[], long n);
86
87void fitngauss(double t[], double y[], double sigma[], long nvec,
88 double acoef[], double acoef_er[], int para, int *status2);
89
90cpl_error_code espda_create_line_table(cpl_table **tab,
91 cpl_size size);
92
93
94#endif
static double sigma
int size
int * y
int * x
int n
Definition: xsh_detmon_lg.c:92
cpl_error_code esp_det_line(cpl_table *spec_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[])
cpl_size get_nm_to_index_scale(cpl_table *spec_total, double nm)
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)
cpl_error_code espda_create_line_table(cpl_table **tab, cpl_size size)
Create a LINE table with a given row size.
void find_left_right_continuum_pos(int *xind1, int *xind2, cpl_table *spec_region, double cont_rejt, double line)
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)
cpl_error_code select_local_spec(cpl_table *spec_total, double ew_space, double lambda, cpl_table **spec_region)
Select_local_spec.
DOUBLE vec[4]