MUSE Pipeline Reference Manual  2.1.1
muse_cplwrappers.h
1 /* -*- Mode: C; tab-width: 2; indent-tabs-mode: nil; c-basic-offset: 2 -*- */
2 /* vim:set sw=2 sts=2 et cin: */
3 /*
4  * This file is part of the MUSE Instrument Pipeline
5  * Copyright (C) 2005-2014 European Southern Observatory
6  *
7  * This program is free software; you can redistribute it and/or modify
8  * it under the terms of the GNU General Public License as published by
9  * the Free Software Foundation; either version 2 of the License, or
10  * (at your option) any later version.
11  *
12  * This program is distributed in the hope that it will be useful,
13  * but WITHOUT ANY WARRANTY; without even the implied warranty of
14  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15  * GNU General Public License for more details.
16  *
17  * You should have received a copy of the GNU General Public License
18  * along with this program; if not, write to the Free Software
19  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
20  */
21 
22 #ifndef MUSE_CPLWRAPPERS_H
23 #define MUSE_CPLWRAPPERS_H
24 
25 /*----------------------------------------------------------------------------*
26  * Includes *
27  *----------------------------------------------------------------------------*/
28 #include <cpl.h>
29 
30 /*----------------------------------------------------------------------------*
31  * Defines *
32  *----------------------------------------------------------------------------*/
33 
34 /*----------------------------------------------------------------------------*
35  * Data structures *
36  *----------------------------------------------------------------------------*/
37 
41 /*----------------------------------------------------------------------------*/
51 /*----------------------------------------------------------------------------*/
52 typedef struct {
54  const char *name;
56  cpl_type type;
58  const char *unit;
60  const char *format;
62  const char *description;
64  const int required;
66 
67 /*----------------------------------------------------------------------------*/
71 /*----------------------------------------------------------------------------*/
72 
73 typedef enum {
83 
84 /*----------------------------------------------------------------------------*/
88 /*----------------------------------------------------------------------------*/
89 typedef int (*muse_cplmatrix_element_compare_func)(double aValue1, double aValue2);
90 
93 /*----------------------------------------------------------------------------*
94  * Function prototypes *
95  *----------------------------------------------------------------------------*/
96 cpl_error_code muse_cplimage_or(cpl_image *, const cpl_image *, unsigned int);
97 cpl_image *muse_cplimage_concat_x(const cpl_image *, const cpl_image *);
98 cpl_image *muse_cplimage_concat_y(const cpl_image *, const cpl_image *);
99 cpl_image *muse_cplimage_filter_median_subtract(cpl_image *, unsigned int, unsigned int);
100 cpl_vector *muse_cplimage_slope_window(const cpl_image *, const cpl_size *);
101 double muse_cplimage_get_percentile(const cpl_image *, double);
102 
103 cpl_image *muse_cplimagelist_collapse_or_create(const cpl_imagelist *);
104 
105 cpl_mask *muse_cplmask_adapt_to_image(const cpl_mask *, const cpl_image *);
106 
107 cpl_matrix *muse_cplmatrix_multiply_create(const cpl_matrix *, const cpl_matrix *);
108 cpl_array *muse_cplmatrix_where(const cpl_matrix *, double, muse_cplmatrix_element_compare_func);
109 cpl_matrix *muse_cplmatrix_extract_selected(const cpl_matrix *, const cpl_array *);
110 
111 double muse_cplvector_get_adev_const(const cpl_vector *, double);
112 double muse_cplvector_get_median_dev(cpl_vector *, double *);
113 double muse_cplvector_get_semiquartile(cpl_vector *);
114 cpl_error_code muse_cplvector_threshold(cpl_vector *, double, double, double, double);
115 cpl_error_code muse_cplvector_erase_element(cpl_vector *, int);
116 cpl_size muse_cplvector_count_unique(const cpl_vector *);
117 cpl_vector *muse_cplvector_get_unique(const cpl_vector *);
118 
119 cpl_table *muse_cpltable_new(const muse_cpltable_def *, cpl_size);
120 cpl_table *muse_cpltable_load(const char *, const char *, const muse_cpltable_def []);
121 cpl_error_code muse_cpltable_check(const cpl_table *, const muse_cpltable_def *);
122 cpl_error_code muse_cpltable_append_file(const cpl_table *, const char *, const char *, const muse_cpltable_def[]);
123 cpl_array *muse_cpltable_extract_column(cpl_table *, const char *);
124 cpl_error_code muse_cpltable_copy_array(cpl_table *, const char *, const cpl_array *);
125 cpl_array *muse_cpltable_get_array_copy(cpl_table *, const char *, cpl_size);
126 cpl_size muse_cpltable_find_sorted(const cpl_table *, const char *, double);
127 
128 cpl_array *muse_cplarray_new_from_image(const cpl_image *);
129 cpl_error_code muse_cplarray_poly1d(cpl_array *, const cpl_array *);
130 double muse_cplarray_poly1d_double(double, const cpl_array *);
131 cpl_array *muse_cplarray_extract(cpl_array *, cpl_size, cpl_size);
132 cpl_error_code muse_cplarray_add_window(cpl_array *a, cpl_size, const cpl_array *);
133 cpl_error_code muse_cplarray_sort(cpl_array *, cpl_boolean order);
134 cpl_bivector *muse_cplarray_histogram(const cpl_array *, double, double, double);
135 cpl_size muse_cplarray_find_sorted(const cpl_array *, double);
136 cpl_boolean muse_cplarray_has_duplicate(const cpl_array *);
137 cpl_error_code muse_cplarray_dump_name(const cpl_array *, const char *);
138 cpl_error_code muse_cplarray_erase_invalid(cpl_array *);
139 cpl_size muse_cplarray_erase_outliers(cpl_array *, const cpl_bivector *, cpl_size, double);
140 cpl_array *muse_cplarray_diff(const cpl_array *, int);
141 cpl_error_code muse_cplarray_erf(cpl_array *);
142 cpl_error_code muse_cplarray_exp(cpl_array *);
143 cpl_array *muse_cplarray_interpolate_linear(const cpl_array *, const cpl_array *, const cpl_array *);
144 cpl_array *muse_cplarray_interpolate_table_linear(const cpl_array *, const cpl_table *, const char *, const char *);
145 cpl_array *muse_cplarray_new_from_delimited_string(const char *, const char *);
146 cpl_array *muse_cplarray_string_to_double(const cpl_array *);
147 
148 cpl_parameter *muse_cplparamerterlist_find_prefix(cpl_parameterlist *, const char *, const char *);
149 cpl_parameterlist *muse_cplparameterlist_from_propertylist(const cpl_propertylist *, int);
150 cpl_parameterlist *muse_cplparameterlist_duplicate(const cpl_parameterlist *);
151 
152 cpl_error_code muse_cplpropertylist_update_long_long(cpl_propertylist *, const char *, cpl_size);
153 
154 cpl_error_code muse_cplframeset_erase_all(cpl_frameset *);
155 cpl_error_code muse_cplframeset_erase_duplicate(cpl_frameset *);
156 
157 void muse_cplerrorstate_dump_some(unsigned, unsigned, unsigned);
159 
160 #endif /* MUSE_CPLWRAPPERS_H */
cpl_error_code muse_cplarray_erase_invalid(cpl_array *)
Erase all invalid values from an array.
cpl_type type
Column type (use CPL_TYPE_POINTER for array columns).
int(* muse_cplmatrix_element_compare_func)(double aValue1, double aValue2)
Matrix element comparison function type definition.
cpl_array * muse_cplarray_interpolate_table_linear(const cpl_array *, const cpl_table *, const char *, const char *)
Linear interpolation of a 1d column.
const int required
Is column required?
cpl_error_code muse_cplarray_dump_name(const cpl_array *, const char *)
Dump a numerical array to stdout with a name prefixed to each line.
cpl_image * muse_cplimage_filter_median_subtract(cpl_image *, unsigned int, unsigned int)
Subtract a median-filtered version of the input image from itself.
cpl_array * muse_cplarray_extract(cpl_array *, cpl_size, cpl_size)
Create an array from a section of another array.
muse_cplframework_type muse_cplframework(void)
Return the CPL framework the recipe is run under.
cpl_array * muse_cplarray_string_to_double(const cpl_array *)
Convert a string array into an array of type double.
cpl_error_code muse_cplarray_sort(cpl_array *, cpl_boolean order)
Sort float, int or double array by quicksort.
const char * unit
Column unit, or NULL.
cpl_error_code muse_cplarray_poly1d(cpl_array *, const cpl_array *)
Apply a polynomial to an array.
double muse_cplvector_get_median_dev(cpl_vector *, double *)
Compute the median and average absolute deviation against the median of a vector. ...
cpl_vector * muse_cplvector_get_unique(const cpl_vector *)
Separate out all unique entries in a given vector into a new one.
double muse_cplvector_get_semiquartile(cpl_vector *)
compute the semi-quartile range of a vector of elements
const char * description
Table column description, or NULL.
cpl_boolean muse_cplarray_has_duplicate(const cpl_array *)
Check, if an array contains duplicate values.
cpl_image * muse_cplimage_concat_y(const cpl_image *, const cpl_image *)
Concatenate two images in y direction.
cpl_error_code muse_cpltable_check(const cpl_table *, const muse_cpltable_def *)
Check whether the table contains the fields of the definition.
cpl_table * muse_cpltable_load(const char *, const char *, const muse_cpltable_def[])
Load a table from disk (and check against definition).
cpl_matrix * muse_cplmatrix_extract_selected(const cpl_matrix *, const cpl_array *)
Extract the elements given by the index array.
cpl_table * muse_cpltable_new(const muse_cpltable_def *, cpl_size)
Create an empty table according to the specified definition.
cpl_error_code muse_cplarray_exp(cpl_array *)
Compute the exponential function of array elements.
cpl_error_code muse_cplimage_or(cpl_image *, const cpl_image *, unsigned int)
Provide an &#39;OR&#39; operation of two integer images.
cpl_array * muse_cplarray_new_from_delimited_string(const char *, const char *)
Convert a delimited string into an array of strings.
cpl_array * muse_cpltable_get_array_copy(cpl_table *, const char *, cpl_size)
Return the copy of an array of a table cell.
cpl_array * muse_cpltable_extract_column(cpl_table *, const char *)
Create an array from a section of a column.
cpl_error_code muse_cplvector_erase_element(cpl_vector *, int)
delete the given element from the input vector
const char * name
Column name.
const char * format
Default print format, or NULL.
cpl_array * muse_cplarray_new_from_image(const cpl_image *)
Copy the image data into an array.
cpl_bivector * muse_cplarray_histogram(const cpl_array *, double, double, double)
Create a histogram for a numerical array.
cpl_array * muse_cplarray_interpolate_linear(const cpl_array *, const cpl_array *, const cpl_array *)
Linear interpolation of a 1d array.
double muse_cplimage_get_percentile(const cpl_image *, double)
Get the percentile of an image.
void muse_cplerrorstate_dump_some(unsigned, unsigned, unsigned)
Dump some CPL errors.
cpl_vector * muse_cplimage_slope_window(const cpl_image *, const cpl_size *)
Compute slopes of an image, both horizontally and vertically.
cpl_error_code muse_cplarray_erf(cpl_array *)
Compute the error function of array elements.
cpl_error_code muse_cplarray_add_window(cpl_array *a, cpl_size, const cpl_array *)
Add the value of an array to a window of a table column.
cpl_error_code muse_cplpropertylist_update_long_long(cpl_propertylist *, const char *, cpl_size)
Update an integer-like property irrespective of the real type.
cpl_size muse_cpltable_find_sorted(const cpl_table *, const char *, double)
Find a row in a table.
cpl_image * muse_cplimage_concat_x(const cpl_image *, const cpl_image *)
Concatenate two images in x direction.
cpl_error_code muse_cpltable_copy_array(cpl_table *, const char *, const cpl_array *)
Copy an array into a table.
cpl_array * muse_cplarray_diff(const cpl_array *, int)
Build the difference of any element and one of the next elements.
cpl_parameterlist * muse_cplparameterlist_duplicate(const cpl_parameterlist *)
Duplicate a CPL parameter list.
cpl_error_code muse_cplframeset_erase_duplicate(cpl_frameset *)
Erase all duplicate frames from a frameset.
cpl_error_code muse_cplframeset_erase_all(cpl_frameset *)
Erase all frames in a frameset.
cpl_size muse_cplarray_erase_outliers(cpl_array *, const cpl_bivector *, cpl_size, double)
Erase outliers from an array using histogram information.
cpl_matrix * muse_cplmatrix_multiply_create(const cpl_matrix *, const cpl_matrix *)
Compute the element-wise product of two matrices.
cpl_array * muse_cplmatrix_where(const cpl_matrix *, double, muse_cplmatrix_element_compare_func)
Select matrix elements according to a condition.
muse_cplframework_type
Type for the framework that called the recipe.
cpl_size muse_cplvector_count_unique(const cpl_vector *)
Count the number of unique entries in a given vector.
Definition of a cpl table structure.
double muse_cplvector_get_adev_const(const cpl_vector *, double)
Compute the average absolute deviation of a (constant) vector.
cpl_parameter * muse_cplparamerterlist_find_prefix(cpl_parameterlist *, const char *, const char *)
Return the full recipe parameter belonging to prefix and shortname.
double muse_cplarray_poly1d_double(double, const cpl_array *)
Apply a polynomial to a double value.
cpl_size muse_cplarray_find_sorted(const cpl_array *, double)
Find a row in an array.
cpl_error_code muse_cplvector_threshold(cpl_vector *, double, double, double, double)
Threshold a vector to a given interval.
cpl_mask * muse_cplmask_adapt_to_image(const cpl_mask *, const cpl_image *)
Adapt mask with masked region in one quadrant to size of an image.
cpl_parameterlist * muse_cplparameterlist_from_propertylist(const cpl_propertylist *, int)
Recreate a cpl_parameterlist from the RECi headers of an output MUSE product.
cpl_image * muse_cplimagelist_collapse_or_create(const cpl_imagelist *)
Compute the OR of an image list to a single image.
cpl_error_code muse_cpltable_append_file(const cpl_table *, const char *, const char *, const muse_cpltable_def[])
Save a table to disk (into a FITS extension)