X-shooter Pipeline Reference Manual 3.8.15
xsh_utils.h
Go to the documentation of this file.
1/* a
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: 2013-07-15 08:10:37 $
23 * $Revision: 1.114 $
24 * $Name: not supported by cvs2svn $
25 */
26#ifndef XSH_UTILS_H
27#define XSH_UTILS_H
28
29/*----------------------------------------------------------------------------
30 Includes
31 ----------------------------------------------------------------------------*/
32
33#include <cpl.h>
34#include <sys/time.h>
35#include <sys/resource.h>
36#include <xsh_globals.h>
37#include <xsh_data_grid.h>
39#include <xsh_data_instrument.h>
40#include <xsh_parameters.h>
41
42#ifndef M_PI
43#define M_PI 3.1415926535897932384626433832795
44#endif
45
46#define XSH_MAX(A,B)\
47 A > B ? A : B
48
49#define XSH_MALLOC( POINTER, TYPE, SIZE) \
50 assure(POINTER == NULL, CPL_ERROR_ILLEGAL_OUTPUT,\
51 "Try to allocate non NULL pointer");\
52 POINTER = (TYPE*)(cpl_malloc(SIZE*sizeof(TYPE)));\
53 assure (POINTER != NULL, CPL_ERROR_ILLEGAL_OUTPUT,\
54 "Memory allocation failed!")
55
56#define XSH_CALLOC( POINTER, TYPE, SIZE) \
57 assure(POINTER == NULL, CPL_ERROR_ILLEGAL_OUTPUT,\
58 "Try to allocate non NULL pointer");\
59 POINTER = (TYPE*)(cpl_calloc(SIZE,sizeof(TYPE)));\
60 assure (POINTER != NULL, CPL_ERROR_ILLEGAL_OUTPUT,\
61 "Memory allocation failed!")
62
63#define XSH_REALLOC( POINTER, TYPE, SIZE ) \
64 assure(POINTER != NULL, CPL_ERROR_ILLEGAL_INPUT,\
65 "Try to re-allocate NULL pointer") ;\
66 POINTER = (TYPE *)cpl_realloc(POINTER,SIZE*sizeof(TYPE)));\
67 assure( POINTER != NULL, CPL_ERROR_ILLEGAL_OUTPUT,\
68 "Memory re-allocation failed!")
69
70#define XSH_NEW_PROPERTYLIST( POINTER) \
71 assure(POINTER == NULL, CPL_ERROR_ILLEGAL_OUTPUT,\
72 "Try to allocate non NULL pointer");\
73 POINTER = cpl_propertylist_new();\
74 assure (POINTER != NULL, CPL_ERROR_ILLEGAL_OUTPUT,\
75 "Memory allocation for propertylist failed!")
76
77#define XSH_NEW_FRAME( POINTER) \
78 assure(POINTER == NULL, CPL_ERROR_ILLEGAL_OUTPUT,\
79 "Try to allocate non NULL pointer");\
80 POINTER = cpl_frame_new();\
81 assure (POINTER != NULL, CPL_ERROR_ILLEGAL_OUTPUT,\
82 "Memory allocation for frame failed!")
83
84#define XSH_NEW_FRAMESET( POINTER) \
85 assure(POINTER == NULL, CPL_ERROR_ILLEGAL_OUTPUT,\
86 "Try to allocate non NULL pointer");\
87 POINTER = cpl_frameset_new();\
88 assure (POINTER != NULL, CPL_ERROR_ILLEGAL_OUTPUT,\
89 "Memory allocation for frameset failed!")
90
91
92#define XSH_FREE(POINTER)\
93 if(POINTER!=NULL) cpl_free(POINTER);\
94 POINTER = NULL
95
96#define XSH_PREFIX(prefix,name,instr) \
97 XSH_FREE(prefix);\
98 prefix = xsh_stringcat_any(name,"_",\
99 xsh_instrument_arm_tostring(instr ),\
100 "" ) ;\
101 XSH_ASSURE_NOT_NULL(prefix)
102
103#define XSH_MODE_PREFIX(prefix,name,instr) \
104 XSH_FREE(prefix);\
105 prefix = xsh_stringcat_any(name,"_",\
106 xsh_instrument_mode_tostring(instr ),\
107 "_",\
108 xsh_instrument_arm_tostring(instr ),\
109 "" ) ;\
110 XSH_ASSURE_NOT_NULL(prefix)
111
112
113#define XSH_TABLE_NEW_COL(TABLE, NAME, UNIT, TYPE) \
114 check( cpl_table_new_column(TABLE, NAME, TYPE));\
115 check( cpl_table_set_column_unit( TABLE, NAME, UNIT))
116
117#define BOOLEAN_TO_STRING( boolean) \
118 boolean == 0 ? "false" : "true"
119typedef struct{
120 void* data;
121 int idx;
123
124
125/*
126 Structure to return gaussian fit results
127*/
128typedef struct {
129 double peakpos,
135
136enum {
139} ;
140
141
142/*----------------------------------------------------------------------------
143 INLINE FONCTIONS (could not be used outside the library)
144----------------------------------------------------------------------------*/
145long xsh_round_double(double x);
146double xsh_max_double(double x, double y);
147double xsh_pow_int(double x, int y);
148
149
150
151/*----------------------------------------------------------------------------
152 Prototypes
153----------------------------------------------------------------------------*/
154void xsh_random_init(void);
155int
156xsh_get_random_int_window(const int v1, const int v2);
157double
158xsh_get_random_double_window(const double v1, const double v2);
159
160cpl_frame* xsh_frame_inv( cpl_frame* in, const char *filename,
161 xsh_instrument* instr);
162
163cpl_frame* xsh_frame_abs( cpl_frame* in, xsh_instrument* instr,
164 cpl_frame** sign);
165
166cpl_frame* xsh_frame_mult( cpl_frame* in, xsh_instrument* instr,
167 cpl_frame* sign);
168
169cpl_parameterlist*
170xsh_parameterlist_duplicate(const cpl_parameterlist* pin);
171
172void
173xsh_plist_dump(cpl_propertylist *plist);
174
175cpl_error_code
176xsh_frameset_dump(cpl_frameset* set);
177
178char * xsh_get_basename(const char *filename);
179
180const char * xsh_get_license(void) ;
181void xsh_init(void);
182
183int xsh_min_int(int x, int y);
184int xsh_max_int(int x, int y);
185
186void xsh_free(const void *mem);
187
188void xsh_free_temporary_files(void);
189cpl_error_code xsh_end(const char *recipe_id, cpl_frameset *frames,
190 cpl_parameterlist * list );
191cpl_error_code xsh_begin(cpl_frameset *frames,
192 const cpl_parameterlist *parameters,
193 xsh_instrument ** instr,
194 cpl_frameset ** raws, cpl_frameset ** calib,
195 const char * tag_list[],
196 int tag_list_size,
197 const char *recipe_id,
198 unsigned int binary_version, const char *short_descr);
199void xsh_add_temporary_file( const char *name ) ;
200
201cpl_error_code xsh_get_property_value(const cpl_propertylist *plist,
202 const char *keyword, cpl_type keywordtype, void *result);
203
204char *xsh_sdate_utc( time_t * t ) ;
205
206char *xsh_stringdup (const char *s1);
207char *xsh_stringcat (const char *s1, const char *s2);
208char *xsh_stringcat_3(const char *s1, const char *s2, const char *s3);
209char *xsh_stringcat_4(const char *s1, const char *s2, const char *s3,
210 const char *s4);
211char *xsh_stringcat_5(const char *s1, const char *s2, const char *s3,
212 const char *s4, const char *s5);
213char *xsh_stringcat_6(const char *s1, const char *s2, const char *s3,
214 const char *s4, const char *s5, const char *s6);
215char *xsh_stringcat_any( const char *s, ...) ;
216void xsh_reindex(double* data, int* idx, int size);
217void xsh_reindex_float( float * data, int* idx, int size);
218void xsh_reindex_int( int * data, int* idx, int size);
219int* xsh_sort(void* base, size_t nmemb, size_t size,
220 int (*compar)(const void *, const void *));
221void xsh_tools_min_max(int size, double *tab, double* min, double* max);
222void xsh_tools_get_statistics(double* tab, int size, double* median,
223 double* mean, double* stdev);
224
225void xsh_free_table(cpl_table **t);
226void xsh_free_image(cpl_image **i);
227void xsh_free_mask(cpl_mask **m);
228void xsh_free_imagelist(cpl_imagelist **i);
229void xsh_free_propertylist(cpl_propertylist **p);
230void xsh_free_polynomial(cpl_polynomial **p);
231void xsh_free_matrix(cpl_matrix **m);
232void xsh_free_array(cpl_array **v);
233void xsh_free_vector(cpl_vector **v);
234void xsh_free_stats(cpl_stats **s);
235void xsh_unwrap_image( cpl_image **i);
236void xsh_unwrap_vector(cpl_vector **v);
237void xsh_unwrap_array(cpl_array **a);
238void xsh_unwrap_bivector_vectors(cpl_bivector **b);
239void xsh_free_parameterlist(cpl_parameterlist **p);
240void xsh_free_parameter(cpl_parameter **p);
241void xsh_free_frameset(cpl_frameset **f);
242void xsh_free_frame(cpl_frame **f);
243
244void xsh_show_time( const char * comment ) ;
245
246cpl_error_code xsh_tools_sort_double( double * pix_arr, int size ) ;
247cpl_error_code xsh_tools_sort_float( float * pix_arr, int size ) ;
248cpl_error_code xsh_tools_sort_int( int * pix_arr, int size ) ;
249void xsh_tools_tchebitchev_transform_tab(int size, double* pos, double min,
250 double max, double* tcheb_pos);
251double xsh_tools_tchebitchev_transform(double pos, double min,
252 double max);
253double xsh_tools_tchebitchev_reverse_transform(double pos, double min,
254 double max);
255
256cpl_vector* xsh_tools_tchebitchev_poly_eval( int n, double X);
257double xsh_tools_get_median_double( double *array, int size ) ;
258int xsh_tools_running_median_1d_get_max( double * tab, int size, int wsize ) ;
259
260void xsh_image_fit_spline(cpl_image* img, xsh_grid* grid);
261
262void xsh_vector_fit_gaussian( cpl_vector * x, cpl_vector * y,
263 XSH_GAUSSIAN_FIT * result ) ;
264double xsh_vector_get_err_median( cpl_vector *vect);
265double xsh_vector_get_err_mean( cpl_vector *vect);
266
267int xsh_debug_level_set( int level ) ;
268int xsh_debug_level_get( void ) ;
269const char * xsh_debug_level_tostring( void ) ;
270
271int xsh_time_stamp_set( int ts ) ;
272int xsh_time_stamp_get( void ) ;
273
274void xsh_mem_dump( const char * prompt) ;
276 const cpl_imagelist * imlist,
277 double sigma_low,
278 double sigma_upp,
279 const int niter);
280
281double convert_bin_to_data( double bin_data, int binning);
282double convert_data_to_bin( double data, int binning);
283
284cpl_frameset * xsh_order_frameset_by_date( cpl_frameset * frameset ) ;
285
286cpl_error_code xsh_set_cd_matrix(cpl_propertylist* plist);
287cpl_error_code xsh_set_cd_matrix1d(cpl_propertylist* plist);
288cpl_error_code xsh_set_cd_matrix2d(cpl_propertylist* plist);
289cpl_error_code xsh_set_cd_matrix3d(cpl_propertylist* plist);
290int xsh_erase_table_rows(cpl_table *t, const char *column,
291 cpl_table_select_operator operator,
292 double value);
293int xsh_select_table_rows(cpl_table *t, const char *column,
294 cpl_table_select_operator operator,
295 double value);
296
297
300 const char *X1, const char *X2, const char *Y,
301 const char *sigmaY,
302 int degree1, int degree2,
303 const char *polynomial_fit, const char *residual_square,
304 const char *variance_fit,
305 double *mse, double *red_chisq,
306 polynomial **variance, double kappa,
307 double min_reject);
308cpl_error_code
309xsh_check_input_is_unbinned(cpl_frame* in);
310
311cpl_error_code
312xsh_update_pheader_in_image_multi(cpl_frame *frame,
313 const cpl_propertylist* pheader);
314
315
316cpl_error_code
317xsh_monitor_flux(cpl_frame* frm_ima,const cpl_frame* frm_tab,
318 xsh_instrument* instrument, const char* qc_key_prefix);
319cpl_error_code
320xsh_frameset_dump_nod_info(cpl_frameset* set);
321void
322xsh_frame_spectrum_save(cpl_frame* frm,const char* name_o);
323void
324xsh_frame_image_save(cpl_frame* frm,const char* name_o);
325
326void
327xsh_frame_table_save(cpl_frame* frm,const char* name_o);
328
329char*
330xsh_set_recipe_file_prefix(cpl_frameset* raw,const char* recipe);
331const char*
332xsh_set_recipe_sky_file_prefix(char* rec_prefix);
333cpl_frame*
334xsh_frameset_average(cpl_frameset *set, const char* tag);
335cpl_frame* xsh_frameset_add( cpl_frameset *set, xsh_instrument *instr,const int decode_bp);
336
337int xsh_fileutils_move (const char *srcpath, const char *dstpath);
338int xsh_fileutils_copy (const char * srcpath, const char * dstpath);
339void xsh_add_product_file( const char *name);
340void xsh_free_product_files( void);
341const char* xsh_string_tolower(char* s);
342const char* xsh_string_toupper(char* s);
343double
344xsh_spline_hermite_table( double xp, const cpl_table *t, const char *column_x,
345 const char *column_y, int *istart );
346double
347xsh_spline_hermite( double xp, const double *x, const double *y, int n, int *istart );
348cpl_frame*
349xsh_util_multiply_by_response(cpl_frame* merged_sci, cpl_frame* response,
350 const char* tag);
351
352cpl_frame*
353xsh_util_multiply_by_response_ord(cpl_frame* merged_sci, cpl_frame* response,
354 const char* tag);
355
356cpl_frame*
357xsh_util_frameset_collapse_mean(cpl_frameset* set,
359
360cpl_frame*
361xsh_spectrum_resample(cpl_frame* frame_inp,
362 const double wstep,
363 const double wmin,
364 const double wmax,
365 xsh_instrument* instr);
366
367cpl_frame*
368xsh_spectrum_resample2(cpl_frame* frame_inp,
369 const double wstep,
370 const double wmin,
371 const double wmax,
372 xsh_instrument* instr);
373
374
375cpl_frame*
376xsh_spectrum_interpolate(cpl_frame* table_frame,
377 const double wstep,
378 const double wmin,
379 const double wmax);
380
381cpl_frame*
382xsh_spectrum_interpolate_linear(cpl_frame* table_frame,
383 const double wstep,
384 const double wmin,
385 const double wmax);
386
387cpl_image*
388xsh_vector_to_image(const cpl_vector* vector,cpl_type type);
389cpl_vector *
390xsh_image_to_vector( cpl_image * spectrum );
391
392cpl_image *
393xsh_normalize_spectrum_image(const cpl_image *spectrum,
394 const cpl_image *spectrum_error,
395 const cpl_propertylist *spectrum_header,
396 const int binx,
397 const double gain,
398 const double exptime,
399 const double airmass,
400 const int n_traces,
401 const cpl_table *atm_extinction,
402 cpl_image **scaled_error);
403
404cpl_frame *
405xsh_normalize_spectrum(const cpl_frame *obj_frame,
406 const cpl_frame *atm_ext_frame,
407 cpl_boolean correct_binning,
409 const char* tag);
410
411cpl_frame *
412xsh_normalize_spectrum_ord(const cpl_frame *obj_frame,
413 const cpl_frame *atm_ext_frame,
414 cpl_boolean correct_binning,
416 const char* tag);
417
418void xsh_array_clip_mean( cpl_array *array, double kappa, int niter,
419 double frac_min, double *mean, double *stdev);
420
421void xsh_array_clip_median( cpl_array *array, double kappa, int niter,
422 double frac_min, double *median, double *stdev);
423
424void xsh_array_clip_poly1d( cpl_vector *pos_array, cpl_vector *val_array,
425 double kappa, int niter, double frac_min, int deg, cpl_polynomial **poly,
426 double *chisq, int **flags);
427
428cpl_error_code
429xsh_rectify_params_set_defaults(cpl_parameterlist* pars,
430 const char* rec_id,
431 xsh_instrument* inst,
432 xsh_rectify_param * rectify_par);
433
434cpl_error_code
435xsh_remove_crh_single_params_set_defaults(cpl_parameterlist* pars,
436 const char* rec_id,
437 xsh_instrument* inst,
438 xsh_remove_crh_single_param * rectify_par);
439
440void xsh_gsl_init_gaussian_fit( cpl_vector *xpos_vect, cpl_vector *ypos_vect,
441 double *init_par);
442void xsh_gsl_fit_gaussian( cpl_vector *xpos_vect, cpl_vector *ypos_vect, int deg,
443 double *params, double *errs, int *status);
444
445double xsh_hms2deg(const double hms);
446double xsh_sess2deg(const double sess);
447
448double*
450 double * line_i,
451 int width_i,
452 double * line_t,
453 int width_t,
454 int half_search,
455 int normalise,
456 double * xcorr_max,
457 double * delta
458 );
459
460int raw_mjd_frame_compare(const cpl_frame* f1,const cpl_frame* f2);
461cpl_error_code get_average_qc_from_raws(cpl_frameset* fset,xsh_instrument* instrument,cpl_propertylist* qclist);
462cpl_error_code get_chromatic_eff(cpl_frame* eff,cpl_propertylist* qclist);
463cpl_error_code calc_curve_qc(cpl_image* im,xsh_instrument* instrument,cpl_propertylist* qclist);
464
465cpl_error_code calc_flat_slit_qc(cpl_image* im,int xa1, int xa2, cpl_propertylist* qclist);
466cpl_error_code calc_resp_qc(cpl_table* mtab, cpl_table* rtab,xsh_instrument* instrument, int lower,int upper,const char* label, cpl_propertylist* qclist);
467
468#endif
static double exptime
static xsh_instrument * instrument
int binx
int size
int * y
int * x
static SimAnneal s
Definition: xsh_model_sa.c:99
polynomial * xsh_polynomial_regression_2d(cpl_table *t, const char *X1, const char *X2, const char *Y, const char *sigmaY, int degree1, int degree2, const char *polynomial_fit, const char *residual_square, const char *variance_fit, double *mse, double *red_chisq, polynomial **variance, double kappa, double min_reject)
Fit a 2d polynomial to three table columns.
Definition: xsh_utils.c:3423
cpl_frame * xsh_normalize_spectrum(const cpl_frame *obj_frame, const cpl_frame *atm_ext_frame, cpl_boolean correct_binning, xsh_instrument *instrument, const char *tag)
Normalize a spectrum.
Definition: xsh_utils.c:5807
int xsh_max_int(int x, int y)
Maximum of two numbers.
Definition: xsh_utils.c:4415
void xsh_unwrap_vector(cpl_vector **v)
Unwrap a vector and set the pointer to NULL.
Definition: xsh_utils.c:2345
cpl_frame * xsh_frameset_average(cpl_frameset *set, const char *tag)
Dump propertylist.
Definition: xsh_utils.c:875
void xsh_vector_fit_gaussian(cpl_vector *x, cpl_vector *y, XSH_GAUSSIAN_FIT *result)
set debug level
Definition: xsh_utils.c:3101
void xsh_image_fit_spline(cpl_image *img, xsh_grid *grid)
perform spline fit
Definition: xsh_utils.c:2980
void xsh_free_polynomial(cpl_polynomial **p)
Deallocate a polynomial and set the pointer to NULL.
Definition: xsh_utils.c:2194
void xsh_free_vector(cpl_vector **v)
Deallocate a vector and set the pointer to NULL.
Definition: xsh_utils.c:2284
void xsh_free_parameterlist(cpl_parameterlist **p)
Deallocate a parameter list and set the pointer to NULL.
Definition: xsh_utils.c:2224
const char * xsh_set_recipe_sky_file_prefix(char *rec_prefix)
Set recipe sky frames prefix.
Definition: xsh_utils.c:576
cpl_error_code xsh_check_input_is_unbinned(cpl_frame *in)
Check if an input frame is not binned.
Definition: xsh_utils.c:375
cpl_image * xsh_normalize_spectrum_image(const cpl_image *spectrum, const cpl_image *spectrum_error, const cpl_propertylist *spectrum_header, const int binx, const double gain, const double exptime, const double airmass, const int n_traces, const cpl_table *atm_extinction, cpl_image **scaled_error)
Normalize a spectrum.
Definition: xsh_utils.c:5996
int xsh_time_stamp_set(int ts)
set timestamp
Definition: xsh_utils.c:3179
cpl_error_code get_average_qc_from_raws(cpl_frameset *fset, xsh_instrument *instrument, cpl_propertylist *qclist)
Definition: xsh_utils.c:7216
int xsh_fileutils_copy(const char *srcpath, const char *dstpath)
Definition: xsh_utils.c:429
char * xsh_stringcat_4(const char *s1, const char *s2, const char *s3, const char *s4)
String concatenation.
Definition: xsh_utils.c:1791
const char * xsh_string_tolower(char *s)
Convert all uppercase characters in a string into lowercase characters.
Definition: xsh_utils.c:4508
int xsh_select_table_rows(cpl_table *t, const char *column, cpl_table_select_operator operator, double value)
Select table rows.
Definition: xsh_utils.c:3784
char * xsh_stringcat_5(const char *s1, const char *s2, const char *s3, const char *s4, const char *s5)
String concatenation.
Definition: xsh_utils.c:1831
void xsh_tools_min_max(int size, double *tab, double *min, double *max)
computes min & max in ab array
Definition: xsh_utils.c:2799
cpl_error_code xsh_frameset_dump_nod_info(cpl_frameset *set)
Dump frameset nod info.
Definition: xsh_utils.c:1088
const char * xsh_string_toupper(char *s)
Convert all lowercase characters in a string into uppercase characters.
Definition: xsh_utils.c:4542
void xsh_free_image(cpl_image **i)
Deallocate an image and set the pointer to NULL.
Definition: xsh_utils.c:2116
int xsh_debug_level_get(void)
get debug level
Definition: xsh_utils.c:3142
int xsh_fileutils_move(const char *srcpath, const char *dstpath)
Definition: xsh_utils.c:534
cpl_error_code xsh_update_pheader_in_image_multi(cpl_frame *frame, const cpl_propertylist *pheader)
Update FITS header.
Definition: xsh_utils.c:4258
double xsh_spline_hermite(double xp, const double *x, const double *y, int n, int *istart)
Spline interpolation based on Hermite polynomials.
Definition: xsh_utils.c:4577
void xsh_unwrap_bivector_vectors(cpl_bivector **b)
Unwrap a bi-vector and set the pointer to NULL.
Definition: xsh_utils.c:2376
double convert_data_to_bin(double data, int binning)
Definition: xsh_utils.c:3246
cpl_error_code xsh_get_property_value(const cpl_propertylist *plist, const char *keyword, cpl_type keywordtype, void *result)
Read a property value from a property list.
Definition: xsh_utils.c:1600
int xsh_min_int(int x, int y)
Minimum of two numbers.
Definition: xsh_utils.c:4399
cpl_frame * xsh_spectrum_interpolate_linear(cpl_frame *table_frame, const double wstep, const double wmin, const double wmax)
spectra interpolation
Definition: xsh_utils.c:6213
void xsh_free_frame(cpl_frame **f)
Deallocate a frame and set the pointer to NULL.
Definition: xsh_utils.c:2269
double xsh_vector_get_err_median(cpl_vector *vect)
Computes median error on a data set.
Definition: xsh_utils.c:4319
void xsh_free_frameset(cpl_frameset **f)
Deallocate a frame set and set the pointer to NULL.
Definition: xsh_utils.c:2254
int xsh_debug_level_set(int level)
set debug level
Definition: xsh_utils.c:3125
char * xsh_stringcat_6(const char *s1, const char *s2, const char *s3, const char *s4, const char *s5, const char *s6)
String concatenation.
Definition: xsh_utils.c:1875
char * xsh_set_recipe_file_prefix(cpl_frameset *raw, const char *recipe)
Set recipe frames prefix.
Definition: xsh_utils.c:601
char * xsh_stringcat_3(const char *s1, const char *s2, const char *s3)
String concatenation.
Definition: xsh_utils.c:1753
void xsh_free_temporary_files(void)
Free temprary files list.
Definition: xsh_utils.c:1451
cpl_frame * xsh_normalize_spectrum_ord(const cpl_frame *obj_frame, const cpl_frame *atm_ext_frame, cpl_boolean correct_binning, xsh_instrument *instrument, const char *tag)
Normalize a spectrum.
Definition: xsh_utils.c:5893
void xsh_free_array(cpl_array **v)
Deallocate an array and set the pointer to NULL.
Definition: xsh_utils.c:2299
double * xsh_function1d_xcorrelate(double *line_i, int width_i, double *line_t, int width_t, int half_search, int normalise, double *xcorr_max, double *delta)
Definition: xsh_utils.c:7101
double xsh_tools_tchebitchev_transform(double pos, double min, double max)
computes Tchebitchev transformation
Definition: xsh_utils.c:2917
cpl_frame * xsh_util_multiply_by_response_ord(cpl_frame *merged_sci, cpl_frame *response, const char *tag)
Multiply input frame by response frame.
Definition: xsh_utils.c:5016
double xsh_vector_get_err_mean(cpl_vector *vect)
Computes mean error on a data set.
Definition: xsh_utils.c:4354
const char * xsh_debug_level_tostring(void)
set debug level
Definition: xsh_utils.c:3155
void xsh_tools_get_statistics(double *tab, int size, double *median, double *mean, double *stdev)
Compute median, stdev and mean for the tab.
Definition: xsh_utils.c:2421
void xsh_mem_dump(const char *prompt)
Definition: xsh_utils.c:3208
void xsh_free_stats(cpl_stats **s)
Deallocate a stats object and set the pointer to NULL.
Definition: xsh_utils.c:2314
void xsh_reindex_float(float *data, int *idx, int size)
TO BE DESCRIBED.
Definition: xsh_utils.c:2038
void xsh_frame_table_save(cpl_frame *frm, const char *name_o)
Save a table frame.
Definition: xsh_utils.c:4017
cpl_error_code xsh_monitor_flux(cpl_frame *frm_ima, const cpl_frame *frm_tab, xsh_instrument *instrument, const char *qc_key_prefix)
Monitor Flux level along the orders traces given by an input table
Definition: xsh_utils.c:4107
cpl_frame * xsh_frame_mult(cpl_frame *in, xsh_instrument *instr, cpl_frame *sign)
Computes product of two input frames.
Definition: xsh_utils.c:4061
void xsh_random_init(void)
Definition: xsh_utils.c:109
cpl_error_code calc_flat_slit_qc(cpl_image *im, int xa1, int xa2, cpl_propertylist *qclist)
Definition: xsh_utils.c:7429
cpl_vector * xsh_tools_tchebitchev_poly_eval(int n, double X)
Compute tchebitchev Tn(X) first coefficient for tchebitchev polynomial.
Definition: xsh_utils.c:2836
void xsh_plist_dump(cpl_propertylist *plist)
Dump propertylist.
Definition: xsh_utils.c:1026
int xsh_get_random_int_window(const int v1, const int v2)
generates random integer values in range [v1,v2]
Definition: xsh_utils.c:128
double xsh_tools_tchebitchev_reverse_transform(double pos, double min, double max)
computes reverse Tchebitchev transformation
Definition: xsh_utils.c:2954
void xsh_gsl_init_gaussian_fit(cpl_vector *xpos_vect, cpl_vector *ypos_vect, double *init_par)
Definition: xsh_utils.c:6862
cpl_error_code xsh_set_cd_matrix3d(cpl_propertylist *plist)
Set CD matrix.
Definition: xsh_utils.c:744
cpl_frame * xsh_frame_abs(cpl_frame *in, xsh_instrument *instr, cpl_frame **sign)
Computes absolute value of a frame.
Definition: xsh_utils.c:3911
void xsh_free_mask(cpl_mask **m)
Deallocate an image mask and set the pointer to NULL.
Definition: xsh_utils.c:2149
char * xsh_get_basename(const char *filename)
Return base filename.
Definition: xsh_utils.c:1175
int xsh_tools_running_median_1d_get_max(double *tab, int size, int wsize)
get max of a list of doubles after running median
Definition: xsh_utils.c:2748
cpl_error_code calc_curve_qc(cpl_image *im, xsh_instrument *instrument, cpl_propertylist *qclist)
Definition: xsh_utils.c:7322
void xsh_free_product_files(void)
Free temprary files list.
Definition: xsh_utils.c:1491
int raw_mjd_frame_compare(const cpl_frame *f1, const cpl_frame *f2)
Definition: xsh_utils.c:7196
const char * xsh_get_license(void)
Get the pipeline copyright and license.
Definition: xsh_utils.c:1193
char * xsh_stringcat_any(const char *s,...)
Concatenate an arbitrary number of strings.
Definition: xsh_utils.c:1925
void xsh_reindex_int(int *data, int *idx, int size)
TO BE DESCRIBED.
Definition: xsh_utils.c:2070
cpl_error_code xsh_tools_sort_double(double *pix_arr, int size)
Sort a double array.
Definition: xsh_utils.c:2461
cpl_error_code xsh_frameset_dump(cpl_frameset *set)
Dump frameset.
Definition: xsh_utils.c:1054
cpl_error_code xsh_set_cd_matrix2d(cpl_propertylist *plist)
Set CD matrix.
Definition: xsh_utils.c:718
double xsh_hms2deg(const double hms)
Convert a double from hours minute seconds to deg:
Definition: xsh_utils.c:312
int xsh_time_stamp_get(void)
get timestamp
Definition: xsh_utils.c:3196
double xsh_get_random_double_window(const double v1, const double v2)
generates random integer values in range [v1,v2]
Definition: xsh_utils.c:150
cpl_frame * xsh_spectrum_resample(cpl_frame *frame_inp, const double wstep, const double wmin, const double wmax, xsh_instrument *instr)
resample a spectrum
Definition: xsh_utils.c:5336
void xsh_array_clip_poly1d(cpl_vector *pos_array, cpl_vector *val_array, double kappa, int niter, double frac_min, int deg, cpl_polynomial **poly, double *chisq, int **flags)
clip outliers from a 1D poly fit
Definition: xsh_utils.c:6568
char * xsh_stringdup(const char *s1)
String duplication.
Definition: xsh_utils.c:1658
int xsh_erase_table_rows(cpl_table *t, const char *column, cpl_table_select_operator operator, double value)
Erase table rows.
Definition: xsh_utils.c:3849
cpl_error_code xsh_tools_sort_int(int *pix_arr, int size)
Sort an integer array.
Definition: xsh_utils.c:2633
cpl_error_code get_chromatic_eff(cpl_frame *eff, cpl_propertylist *qclist)
Definition: xsh_utils.c:7306
void xsh_frame_image_save(cpl_frame *frm, const char *name_o)
save an image frame
Definition: xsh_utils.c:3983
cpl_frame * xsh_util_frameset_collapse_mean(cpl_frameset *set, xsh_instrument *instrument)
Compute mean frame from a list of (IMAGE) framesets.
Definition: xsh_utils.c:5538
cpl_vector * xsh_image_to_vector(cpl_image *spectrum)
Definition: xsh_utils.c:4680
cpl_frame * xsh_frame_inv(cpl_frame *in, const char *filename, xsh_instrument *instr)
Inverse the flux of a PRE frame.
Definition: xsh_utils.c:3878
void xsh_free_parameter(cpl_parameter **p)
Deallocate a parameter and set the pointer to NULL.
Definition: xsh_utils.c:2239
long xsh_round_double(double x)
Computes round(x)
Definition: xsh_utils.c:4383
cpl_error_code xsh_set_cd_matrix(cpl_propertylist *plist)
Set CD matrix.
Definition: xsh_utils.c:676
cpl_frame * xsh_util_multiply_by_response(cpl_frame *merged_sci, cpl_frame *response, const char *tag)
Multiply input frame by response frame.
Definition: xsh_utils.c:4798
int * xsh_sort(void *base, size_t nmemb, size_t size, int(*compar)(const void *, const void *))
Sort an array and give is index table.
Definition: xsh_utils.c:1966
void xsh_reindex(double *data, int *idx, int size)
TO BE DESCRIBED.
Definition: xsh_utils.c:2008
void xsh_array_clip_median(cpl_array *array, double kappa, int niter, double frac_min, double *median, double *stdev)
median clip of an array
Definition: xsh_utils.c:6493
void xsh_init(void)
Reset library state.
Definition: xsh_utils.c:1160
double xsh_tools_get_median_double(double *array, int size)
Calculates the median value of an array of double.
Definition: xsh_utils.c:2721
cpl_frameset * xsh_order_frameset_by_date(cpl_frameset *frameset)
Order frameset by date.
Definition: xsh_utils.c:3306
char * xsh_sdate_utc(time_t *t)
Definition: xsh_utils.c:1691
void xsh_unwrap_image(cpl_image **i)
Unwrap an image and set the pointer to NULL.
Definition: xsh_utils.c:2330
cpl_image * xsh_imagelist_collapse_sigclip_iter_create(const cpl_imagelist *imlist, double sigma_low, double sigma_upp, const int niter)
Average with sigma-clipping rejection an imagelist to a single image.
Definition: xsh_utils.c:193
void xsh_array_clip_mean(cpl_array *array, double kappa, int niter, double frac_min, double *mean, double *stdev)
mean clip of an array
Definition: xsh_utils.c:6415
cpl_error_code xsh_tools_sort_float(float *pix_arr, int size)
Sort a float array.
Definition: xsh_utils.c:2548
double xsh_sess2deg(const double sess)
Convert a double from ssessagesimal to deg: 203049.197= 20:30:49.197 = 20.5136658333.
Definition: xsh_utils.c:345
void xsh_free_matrix(cpl_matrix **m)
Deallocate a matrix and set the pointer to NULL.
Definition: xsh_utils.c:2209
char * xsh_stringcat(const char *s1, const char *s2)
String concatenation.
Definition: xsh_utils.c:1719
void xsh_gsl_fit_gaussian(cpl_vector *xpos_vect, cpl_vector *ypos_vect, int deg, double *params, double *errs, int *status)
Definition: xsh_utils.c:6936
cpl_error_code xsh_remove_crh_single_params_set_defaults(cpl_parameterlist *pars, const char *rec_id, xsh_instrument *inst, xsh_remove_crh_single_param *rectify_par)
Definition: xsh_utils.c:6723
cpl_error_code xsh_begin(cpl_frameset *frames, const cpl_parameterlist *parameters, xsh_instrument **instr, cpl_frameset **raws, cpl_frameset **calib, const char *tag_list[], int tag_list_size, const char *recipe_id, unsigned int binary_version, const char *short_descr)
Recipe initialization.
Definition: xsh_utils.c:1244
cpl_error_code xsh_rectify_params_set_defaults(cpl_parameterlist *pars, const char *rec_id, xsh_instrument *inst, xsh_rectify_param *rectify_par)
Definition: xsh_utils.c:6681
double convert_bin_to_data(double bin_data, int binning)
Definition: xsh_utils.c:3228
void xsh_free_table(cpl_table **t)
Deallocate a table and set the pointer to NULL.
Definition: xsh_utils.c:2133
cpl_error_code calc_resp_qc(cpl_table *mtab, cpl_table *rtab, xsh_instrument *instrument, int lower, int upper, const char *label, cpl_propertylist *qclist)
Definition: xsh_utils.c:7504
void xsh_free_propertylist(cpl_propertylist **p)
Deallocate a property list and set the pointer to NULL.
Definition: xsh_utils.c:2179
cpl_error_code xsh_end(const char *recipe_id, cpl_frameset *frames, cpl_parameterlist *list)
Recipe termination.
Definition: xsh_utils.c:1519
cpl_image * xsh_vector_to_image(const cpl_vector *vector, cpl_type type)
Convert a vector to a 1d image.
Definition: xsh_utils.c:4745
void xsh_add_product_file(const char *name)
Add temporary file to temprary files list.
Definition: xsh_utils.c:1472
cpl_frame * xsh_spectrum_interpolate(cpl_frame *table_frame, const double wstep, const double wmin, const double wmax)
spectra interpolation
Definition: xsh_utils.c:6306
void xsh_frame_spectrum_save(cpl_frame *frm, const char *name_o)
save an spectrum frame
Definition: xsh_utils.c:3959
cpl_frame * xsh_frameset_add(cpl_frameset *set, xsh_instrument *instr, const int decode_bp)
coadd frames in a frameset
Definition: xsh_utils.c:921
double xsh_spline_hermite_table(double xp, const cpl_table *t, const char *column_x, const char *column_y, int *istart)
Spline interpolation based on Hermite polynomials.
Definition: xsh_utils.c:4648
void xsh_tools_tchebitchev_transform_tab(int size, double *pos, double min, double max, double *tcheb_pos)
computes Tchebitchev transformation
Definition: xsh_utils.c:2881
void xsh_free_imagelist(cpl_imagelist **i)
Deallocate an image list and set the pointer to NULL.
Definition: xsh_utils.c:2164
double xsh_pow_int(double x, int y)
Computes x^y.
Definition: xsh_utils.c:4463
double xsh_max_double(double x, double y)
Maximum of two numbers.
Definition: xsh_utils.c:4449
void xsh_free(const void *mem)
Deallocate memory.
Definition: xsh_utils.c:2102
void xsh_unwrap_array(cpl_array **a)
Unwrap an array and set the pointer to NULL.
Definition: xsh_utils.c:2361
cpl_parameterlist * xsh_parameterlist_duplicate(const cpl_parameterlist *pin)
Extract frames with given tag from frameset.
Definition: xsh_utils.c:771
void xsh_add_temporary_file(const char *name)
Add temporary file to temprary files list.
Definition: xsh_utils.c:1432
cpl_error_code xsh_set_cd_matrix1d(cpl_propertylist *plist)
Set CD matrix.
Definition: xsh_utils.c:701
void xsh_show_time(const char *comment)
show time
Definition: xsh_utils.c:2389
void * data
Definition: xsh_utils.h:120
double kappa
Definition: xsh_detmon_lg.c:81
int niter
Definition: xsh_detmon_lg.c:82
int m
Definition: xsh_detmon_lg.c:91
int n
Definition: xsh_detmon_lg.c:92
#define max(a, b)
@ XSH_DEBUG_LEVEL_HIGH
Definition: xsh_utils.h:138
@ XSH_DEBUG_LEVEL_NONE
Definition: xsh_utils.h:137
@ XSH_DEBUG_LEVEL_LOW
Definition: xsh_utils.h:137
@ XSH_DEBUG_LEVEL_MEDIUM
Definition: xsh_utils.h:138
cpl_frame * xsh_spectrum_resample2(cpl_frame *frame_inp, const double wstep, const double wmin, const double wmax, xsh_instrument *instr)