|
X-shooter Pipeline Reference Manual 3.8.15
|
Go to the source code of this file.
Macros | |
| #define | C(i) (gsl_vector_get(c,(i))) |
| #define | COV(i, j) (gsl_matrix_get(cov,(i),(j))) |
| #define | FIT(i) gsl_vector_get(s->x, i) |
| #define | ERR(i) sqrt(gsl_matrix_get(covar,i,i)) |
Functions | |
| static cpl_size | get_index_from_spec (cpl_table *spec_total, double lambda) |
| static cpl_size | get_pixel_to_nm_scale (cpl_table *spec_total, double nm) |
| get scale from nm to pixels | |
| cpl_error_code | select_local_spec (cpl_table *spec_total, double ew_space, double lambda, cpl_table **spec_region) |
| Select_local_spec. | |
| void | find_left_right_continuum_pos (int *x1, int *x2, cpl_table *spec_region, double cont_rejt, double line) |
| static cpl_table * | esp_spec_deriv (cpl_table *spec_region) |
| static void | esp_spec_smooth (cpl_table *spec_region, int smwidth) |
| void | smooth (double vec[], long n, int w, double svec[]) |
| static void | zeroscenterfind (double iy[], double y[], double dy[], double ddy[], long n, long center[], long *ncenter, double det_line_thres) |
| 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 from the fit. | |
| double | check_ew (cpl_table *line_table, double line, double det_line_resol, int *index_line, int *n_lines, double *ew_error) |
| static void | poly_fitn (double xvec[], double yvec[], double err[], long n, long ord, double coefs[]) |
| void | deriv (double x[], double y[], double dy[], long n) |
| void | fitngauss (double t[], double y[], double sigma[], long nvec, double acoef[], double acoef_er[], int para, int *status2) |
| cpl_error_code | espda_create_line_table (cpl_table **tab, cpl_size size) |
| Create a LINE table with a given row size. | |
| cpl_error_code | esp_fit_lcont (cpl_table *spec_region, double cont_rejt, int cont_iter) |
| Performs a local normalization of the spectrum region. | |
| cpl_error_code | 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. | |
| #define C | ( | i | ) | (gsl_vector_get(c,(i))) |
| #define COV | ( | i, | |
| j | |||
| ) | (gsl_matrix_get(cov,(i),(j))) |
| #define ERR | ( | i | ) | sqrt(gsl_matrix_get(covar,i,i)) |
| double check_ew | ( | cpl_table * | line_table, |
| double | line, | ||
| double | det_line_resol, | ||
| int * | index_line, | ||
| int * | n_lines, | ||
| double * | ew_error | ||
| ) |
Definition at line 403 of file xsh_eqwidth_lib.c.
| void deriv | ( | double | x[], |
| double | y[], | ||
| double | dy[], | ||
| long | n | ||
| ) |
Definition at line 476 of file xsh_eqwidth_lib.c.
Referenced by esp_spec_deriv().
| 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.
| spec_cont_region | spectrum region table (normalized) |
| det_line_thres | Normalized threshold for line acceptance The default value is 0.02 |
| det_line_resol | Minimum accepted separation between two lines (in Angstrom) |
| det_line_smwidth | Number of pixels in the boxcar to be averaged |
This function detects absorption lines in a spectrum after normalizing it.
Definition at line 775 of file xsh_eqwidth_lib.c.
References esp_spec_deriv(), esp_spec_smooth(), espda_create_line_table(), n, x, y, and zeroscenterfind().
Referenced by test_spectrum_detect_lines().
| cpl_error_code esp_fit_lcont | ( | cpl_table * | spec_region, |
| double | cont_rejt, | ||
| int | cont_iter | ||
| ) |
Performs a local normalization of the spectrum region.
| spec_region | spectrum region table |
| cont-rejt | Rejection threshold for the continuum determination. The default value depends on the signal-to-noise parameter for the interval around the line |
| cont-iter | Number of iterations to fit the continuum |
This function updates the region of the spectra for the calculations with the normalization
Definition at line 671 of file xsh_eqwidth_lib.c.
References n, order, poly_fitn(), x, and y.
Referenced by test_spectrum_detect_lines().
| 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 from the fit.
| spec_cont_region | spectrum region table (normalized) |
| line_table | table with the lines to be fitted |
| fit_ngauss_width | Width of the interval used for fitting (nm)??? Currently being used as the guess sigma for the gaussian fitting |
This function detects absorption lines in a spectrum after normalizing it.
Definition at line 332 of file xsh_eqwidth_lib.c.
References fitngauss(), and sigma.
|
static |
Definition at line 195 of file xsh_eqwidth_lib.c.
References deriv(), n, x, and y.
Referenced by esp_det_line().
|
static |
Definition at line 221 of file xsh_eqwidth_lib.c.
References n, smooth(), and y.
Referenced by esp_det_line().
| cpl_error_code espda_create_line_table | ( | cpl_table ** | tab, |
| cpl_size | size | ||
| ) |
Create a LINE table with a given row size.
| tab | New table |
| size | Row size |
This functions creates a CPL table with structure LINE and a given row size.
Definition at line 588 of file xsh_eqwidth_lib.c.
References size.
Referenced by esp_det_line().
| void find_left_right_continuum_pos | ( | int * | x1, |
| int * | x2, | ||
| cpl_table * | spec_region, | ||
| double | cont_rejt, | ||
| double | line | ||
| ) |
Definition at line 162 of file xsh_eqwidth_lib.c.
| void fitngauss | ( | double | t[], |
| double | y[], | ||
| double | sigma[], | ||
| long | nvec, | ||
| double | acoef[], | ||
| double | acoef_er[], | ||
| int | para, | ||
| int * | status2 | ||
| ) |
Definition at line 490 of file xsh_eqwidth_lib.c.
References ERR, expb_df(), expb_f(), expb_fdf(), FIT, N, n, data::para, s, sigma, data::t, x, SimAnneal::x, and y.
Referenced by esp_fit_ngauss().
|
static |
Definition at line 36 of file xsh_eqwidth_lib.c.
Referenced by select_local_spec().
|
static |
get scale from nm to pixels
| spec_1d | spectrum table |
| nm | space scale |
This function returns a the index of the closest wavelenght greatest than lambda
Definition at line 104 of file xsh_eqwidth_lib.c.
Referenced by select_local_spec().
| double maxele_vec | ( | double | vec[], |
| long | nvec | ||
| ) |
|
static |
Definition at line 435 of file xsh_eqwidth_lib.c.
Referenced by esp_fit_lcont().
| cpl_error_code select_local_spec | ( | cpl_table * | spec_total, |
| double | ew_space, | ||
| double | lambda, | ||
| cpl_table ** | spec_region | ||
| ) |
Select_local_spec.
| spec_1d | spectrum table |
| ew_space | parameter for the interval around the line |
| line | line center for the spectral region |
| spec_1d_out | Number of frames in the new frameset |
This function returns a region of the spectra for the calculations
Definition at line 126 of file xsh_eqwidth_lib.c.
References get_index_from_spec(), and get_pixel_to_nm_scale().
Referenced by test_spectrum_detect_lines().
| void smooth | ( | double | vec[], |
| long | n, | ||
| int | w, | ||
| double | svec[] | ||
| ) |
Definition at line 243 of file xsh_eqwidth_lib.c.
References n.
Referenced by esp_spec_smooth(), and xsh_atrous().
|
static |
Definition at line 263 of file xsh_eqwidth_lib.c.
References maxele_vec(), and n.
Referenced by esp_det_line().