X-shooter Pipeline Reference Manual 3.8.15
xsh_data_wavesol.h
Go to the documentation of this file.
1/* *
2 * This file is part of the ESO X-shooter Pipeline *
3 * Copyright (C) 2006 European Southern Observatory *
4 * *
5 * This library is/ free software; you can redistribute it and/or modify *
6 * it under the terms of the GNU General Public License as published by *
7 * the Free Software Foundation; either version 2 of the License, or *
8 * (at your option) any later version. *
9 * *
10 * This program is distributed in the hope that it will be useful, *
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
13 * GNU General Public License for more details. *
14 * *
15 * You should have received a copy of the GNU General Public License *
16 * along with this program; if not, write to the Free Software *
17 * Foundation, 51 Franklin St, Fifth Floor, Boston, MA 02111-1307 USA *
18 * */
19
20/*
21 * $Author: amodigli $
22 * $Date: 2011-12-09 10:05:11 $
23 * $Revision: 1.24 $
24 * $Name: not supported by cvs2svn $
25 */
26#ifndef XSH_DATA_WAVESOL_H
27#define XSH_DATA_WAVESOL_H
28
29#include <cpl.h>
30#include <xsh_cpl_size.h>
31#include <xsh_parameters.h>
32
33#define XSH_SLIT_RANGE 1.2
34
35#define XSH_WAVESOL_TABLE_NB_COL 4
36#define XSH_WAVESOL_TABLE_NB_ROWS 2
37#define XSH_WAVESOL_TABLE_COLNAME_AXIS "AXIS"
38#define XSH_WAVESOL_TABLE_COLNAME_DEGLAMBDA "DEGLAMBDA"
39#define XSH_WAVESOL_TABLE_COLNAME_DEGORDER "DEGORDER"
40#define XSH_WAVESOL_TABLE_COLNAME_DEGSLIT "DEGSLIT"
41
46};
47
48typedef struct{
50 int bin_x;
51 int bin_y;
52 cpl_polynomial* polx;
53 cpl_polynomial* poly;
54 cpl_propertylist* header;
55 cpl_vector* dim;
56 int * coefs;
61 double min_lambda;
62 double max_lambda;
63 double min_order;
64 double max_order;
65 double min_slit;
66 double max_slit;
67 double min_x;
68 double max_x;
69 double min_y;
70 double max_y;
71
73
74//xsh_wavesol* xsh_wavesol_create(cpl_frame* frame);
76 cpl_frame* spectral_format_frame, xsh_detect_arclines_param* p,
78
80void xsh_wavesol_add_poly( xsh_wavesol * to, xsh_wavesol * from ) ;
81void xsh_wavesol_dump( xsh_wavesol * wsol, const char * fname, int nb ) ;
82
83void xsh_wavesol_set_type(xsh_wavesol * wsol, enum wavesol_type type);
85
86cpl_polynomial* xsh_wavesol_get_poly(xsh_wavesol* sol);
87cpl_polynomial* xsh_wavesol_get_polx(xsh_wavesol* sol);
88cpl_propertylist* xsh_wavesol_get_header(xsh_wavesol* sol);
89double xsh_wavesol_eval_polx(xsh_wavesol* sol, double lambda, double order,
90 double slit);
91double xsh_wavesol_eval_poly(xsh_wavesol* sol, double lambda, double order,
92 double slit);
94 double* posdata, double *minpos, double *maxpos, double* lambda, double* order,
95 double* slit, cpl_polynomial* res);
97 double* new_pos, double* lambda,
98 double* order, double* slit,
99 cpl_polynomial* result, char axis) ;
100
101cpl_frame*
102xsh_wavesol_save(xsh_wavesol *w, cpl_table* trace, const char* filename,const char* tag);
104xsh_wavesol * xsh_wavesol_load( cpl_frame * frame,
106cpl_table*
107xsh_wavesol_trace( xsh_wavesol * wsol, double* lambda,
108 double* order, double* slit,int size);
109void xsh_wavesol_set_bin_x( xsh_wavesol * wsol, int bin ) ;
110void xsh_wavesol_set_bin_y( xsh_wavesol * wsol, int bin ) ;
111void xsh_wavesol_apply_shift( xsh_wavesol *wsol, float shift_x, float shift_y);
112
113#endif /* XSH_WAVESOL_H */
static xsh_instrument * instrument
void xsh_wavesol_apply_shift(xsh_wavesol *wsol, float shift_x, float shift_y)
Apply a shift on X and Y to wave solution.
void xsh_wavesol_set_bin_y(xsh_wavesol *wsol, int bin)
Set the bin of wave table in y.
void xsh_wavesol_set_bin_x(xsh_wavesol *wsol, int bin)
Set the bin of wave table in x.
cpl_table * xsh_wavesol_trace(xsh_wavesol *wsol, double *lambda, double *order, double *slit, int size)
void xsh_wavesol_compute(xsh_wavesol *sol, int size, double *posdata, double *minpos, double *maxpos, double *lambda, double *order, double *slit, cpl_polynomial *res)
compute a wavelength solution
cpl_polynomial * xsh_wavesol_get_polx(xsh_wavesol *sol)
get the solution in X
void xsh_wavesol_dump(xsh_wavesol *wsol, const char *fname, int nb)
xsh_wavesol * xsh_wavesol_load(cpl_frame *frame, xsh_instrument *instrument)
load a wavelength solution
double xsh_wavesol_eval_poly(xsh_wavesol *sol, double lambda, double order, double slit)
eval the polynomial solution in Y
void xsh_wavesol_residual(xsh_wavesol *sol, xsh_wavesol *adj, int size, double *new_pos, double *lambda, double *order, double *slit, cpl_polynomial *result, char axis)
xsh_wavesol * xsh_wavesol_create(cpl_frame *spectral_format_frame, xsh_detect_arclines_param *p, xsh_instrument *instrument)
Create a new wavelength solution structure.
cpl_polynomial * xsh_wavesol_get_poly(xsh_wavesol *sol)
get the solution in Y
void xsh_wavesol_free(xsh_wavesol **w)
free wavelength solution structure
double xsh_wavesol_eval_polx(xsh_wavesol *sol, double lambda, double order, double slit)
eval the polynomial solution in X
enum wavesol_type xsh_wavesol_get_type(xsh_wavesol *wsol)
get the type of the wave table
void xsh_wavesol_set_type(xsh_wavesol *wsol, enum wavesol_type type)
set the type of the wave table
cpl_propertylist * xsh_wavesol_get_header(xsh_wavesol *sol)
get header of the table
xsh_wavesol * xsh_wavesol_duplicate(xsh_wavesol *org)
duplicate a wavelength solution structure
void xsh_wavesol_add_poly(xsh_wavesol *to, xsh_wavesol *from)
cpl_frame * xsh_wavesol_save(xsh_wavesol *w, cpl_table *trace, const char *filename, const char *tag)
save a wavelength solution
int size
cpl_polynomial * polx
cpl_propertylist * header
cpl_vector * dim
enum wavesol_type type
cpl_polynomial * poly
wavesol_type
@ XSH_WAVESOL_GUESS
@ XSH_WAVESOL_2D
@ XSH_WAVESOL_UNDEFINED
int order
Definition: xsh_detmon_lg.c:80