![]() |
![]() |
![]() |
PHOEBE Reference Manual | |
---|---|---|---|---|
enum bool; enum PHOEBE_type; char* phoebe_type_get_name (PHOEBE_type type); PHOEBE_vector; PHOEBE_vector* phoebe_vector_new (); PHOEBE_vector* phoebe_vector_new_from_qualifier (char *qualifier); PHOEBE_vector* phoebe_vector_new_from_column (char *filename, int col); PHOEBE_vector* phoebe_vector_new_from_range (int dim, double start, double end); PHOEBE_vector* phoebe_vector_duplicate (PHOEBE_vector *vec); int phoebe_vector_alloc (PHOEBE_vector *vec, int dimension); int phoebe_vector_realloc (PHOEBE_vector *vec, int dimension); int phoebe_vector_pad (PHOEBE_vector *vec, double value); int phoebe_vector_free (PHOEBE_vector *vec); int phoebe_vector_add (PHOEBE_vector *result, PHOEBE_vector *fac1, PHOEBE_vector *fac2); int phoebe_vector_subtract (PHOEBE_vector *result, PHOEBE_vector *fac1, PHOEBE_vector *fac2); int phoebe_vector_multiply (PHOEBE_vector *result, PHOEBE_vector *fac1, PHOEBE_vector *fac2); int phoebe_vector_divide (PHOEBE_vector *result, PHOEBE_vector *fac1, PHOEBE_vector *fac2); int phoebe_vector_raise (PHOEBE_vector *result, PHOEBE_vector *fac1, PHOEBE_vector *fac2); int phoebe_vector_offset (PHOEBE_vector *vec, double offset); int phoebe_vector_sum (PHOEBE_vector *vec, double *sum); int phoebe_vector_mean (PHOEBE_vector *vec, double *mean); int phoebe_vector_median (PHOEBE_vector *vec, double *median); int phoebe_vector_standard_deviation (PHOEBE_vector *vec, double *sigma); int phoebe_vector_multiply_by (PHOEBE_vector *vec, double factor); int phoebe_vector_dot_product (double *result, PHOEBE_vector *fac1, PHOEBE_vector *fac2); int phoebe_vector_vec_product (PHOEBE_vector *result, PHOEBE_vector *fac1, PHOEBE_vector *fac2); int phoebe_vector_submit (PHOEBE_vector *result, PHOEBE_vector *vec); int phoebe_vector_norm (double *result, PHOEBE_vector *vec); int phoebe_vector_dim (int *result, PHOEBE_vector *vec); int phoebe_vector_randomize (PHOEBE_vector *result, double limit); int phoebe_vector_min_max (PHOEBE_vector *vec, double *min, double *max); int phoebe_vector_min_index (PHOEBE_vector *vec, int *index); int phoebe_vector_max_index (PHOEBE_vector *vec, int *index); int phoebe_vector_rescale (PHOEBE_vector *vec, double ll, double ul); bool phoebe_vector_compare (PHOEBE_vector *vec1, PHOEBE_vector *vec2); int phoebe_vector_less_than (bool *result, PHOEBE_vector *vec1, PHOEBE_vector *vec2); int phoebe_vector_leq_than (bool *result, PHOEBE_vector *vec1, PHOEBE_vector *vec2); int phoebe_vector_greater_than (bool *result, PHOEBE_vector *vec1, PHOEBE_vector *vec2); int phoebe_vector_geq_than (bool *result, PHOEBE_vector *vec1, PHOEBE_vector *vec2); int phoebe_vector_append_element (PHOEBE_vector *vec, double val); int phoebe_vector_remove_element (PHOEBE_vector *vec, int index); PHOEBE_matrix; PHOEBE_matrix* phoebe_matrix_new (); int phoebe_matrix_alloc (PHOEBE_matrix *matrix, int cols, int rows); int phoebe_matrix_free (PHOEBE_matrix *matrix); int phoebe_matrix_get_row (PHOEBE_vector *vec, PHOEBE_matrix *matrix, int row); int phoebe_matrix_set_row (PHOEBE_matrix *matrix, PHOEBE_vector *vec, int row); PHOEBE_array; PHOEBE_array* phoebe_array_new (PHOEBE_type type); PHOEBE_array* phoebe_array_new_from_qualifier (char *qualifier); PHOEBE_array* phoebe_array_new_from_column (char *filename, int col); int phoebe_array_alloc (PHOEBE_array *array, int dimension); int phoebe_array_realloc (PHOEBE_array *array, int dimension); PHOEBE_array* phoebe_array_duplicate (PHOEBE_array *array); bool phoebe_array_compare (PHOEBE_array *array1, PHOEBE_array *array2); int phoebe_array_free (PHOEBE_array *array); PHOEBE_vector* phoebe_vector_new_from_array (PHOEBE_array *array); PHOEBE_hist; enum PHOEBE_hist_rebin_type; PHOEBE_hist* phoebe_hist_new (); PHOEBE_hist* phoebe_hist_new_from_arrays (int bins, double *binarray, double *valarray); PHOEBE_hist* phoebe_hist_new_from_file (char *filename); PHOEBE_hist* phoebe_hist_duplicate (PHOEBE_hist *hist); int phoebe_hist_alloc (PHOEBE_hist *hist, int bins); int phoebe_hist_realloc (PHOEBE_hist *hist, int bins); int phoebe_hist_free (PHOEBE_hist *hist); int phoebe_hist_set_ranges (PHOEBE_hist *hist, PHOEBE_vector *bin_centers); int phoebe_hist_set_values (PHOEBE_hist *hist, PHOEBE_vector *values); int phoebe_hist_get_bin_centers (PHOEBE_hist *hist, PHOEBE_vector *bin_centers); int phoebe_hist_get_bin (int *bin, PHOEBE_hist *hist, double r); int phoebe_hist_evaluate (double *y, PHOEBE_hist *hist, double x); int phoebe_hist_integrate (double *integral, PHOEBE_hist *hist, double ll, double ul); int phoebe_hist_shift (PHOEBE_hist *hist, double shift); int phoebe_hist_correlate (double *cfval, PHOEBE_hist *h1, PHOEBE_hist *h2, double sigma1, double sigma2, double ll, double ul, double xi); int phoebe_hist_pad (PHOEBE_hist *hist, double val); int phoebe_hist_crop (PHOEBE_hist *hist, double ll, double ul); bool phoebe_hist_compare (PHOEBE_hist *hist1, PHOEBE_hist *hist2); int phoebe_hist_resample (PHOEBE_hist *out, PHOEBE_hist *in, PHOEBE_hist_rebin_type type); int phoebe_hist_rebin (PHOEBE_hist *out, PHOEBE_hist *in, PHOEBE_hist_rebin_type type); PHOEBE_ld; PHOEBE_passband; enum PHOEBE_curve_type; int phoebe_curve_type_get_name (PHOEBE_curve_type ctype, char **name); enum PHOEBE_column_type; int phoebe_column_type_get_name (PHOEBE_column_type ctype, char **name); int phoebe_column_get_type (PHOEBE_column_type *type, const char *string); enum PHOEBE_data_flag; PHOEBE_curve; PHOEBE_curve* phoebe_curve_new (); PHOEBE_curve* phoebe_curve_new_from_file (char *filename); PHOEBE_curve* phoebe_curve_new_from_pars (PHOEBE_curve_type ctype, int index); PHOEBE_curve* phoebe_curve_duplicate (PHOEBE_curve *curve); int phoebe_curve_alloc (PHOEBE_curve *curve, int dim); int phoebe_curve_realloc (PHOEBE_curve *curve, int dim); int phoebe_curve_compute (PHOEBE_curve *curve, PHOEBE_vector *nodes, int index, PHOEBE_column_type itype, PHOEBE_column_type dtype); int phoebe_curve_transform (PHOEBE_curve *curve, PHOEBE_column_type itype, PHOEBE_column_type dtype, PHOEBE_column_type wtype); int phoebe_curve_alias (PHOEBE_curve *curve, double phmin, double phmax); int phoebe_curve_set_properties (PHOEBE_curve *curve, PHOEBE_curve_type type, char *filename, PHOEBE_passband *passband, PHOEBE_column_type itype, PHOEBE_column_type dtype, PHOEBE_column_type wtype, double sigma); int phoebe_curve_free (PHOEBE_curve *curve); enum PHOEBE_spectrum_dispersion; PHOEBE_spectrum; enum PHOEBE_minimizer_type; int phoebe_minimizer_type_get_name (PHOEBE_minimizer_type minimizer, char **name); PHOEBE_minimizer_feedback; PHOEBE_minimizer_feedback* phoebe_minimizer_feedback_new (); PHOEBE_minimizer_feedback* phoebe_minimizer_feedback_duplicate (PHOEBE_minimizer_feedback *feedback); int phoebe_minimizer_feedback_alloc (PHOEBE_minimizer_feedback *feedback, int tba, int cno); int phoebe_minimizer_feedback_accept (PHOEBE_minimizer_feedback *feedback); int phoebe_minimizer_feedback_free (PHOEBE_minimizer_feedback *feedback); PHOEBE_spectrum* phoebe_spectrum_duplicate (PHOEBE_spectrum *spectrum); PHOEBE_value phoebe_value_duplicate (PHOEBE_type type, PHOEBE_value val);
typedef enum PHOEBE_type { TYPE_INT, TYPE_BOOL, TYPE_DOUBLE, TYPE_STRING, TYPE_INT_ARRAY, TYPE_BOOL_ARRAY, TYPE_DOUBLE_ARRAY, TYPE_STRING_ARRAY, TYPE_CURVE, TYPE_SPECTRUM, TYPE_MINIMIZER_FEEDBACK, TYPE_ANY } PHOEBE_type;
Various data types used in Phoebe.
The integer type. | |
The boolean type. | |
The double precision floating-point type. | |
The string type. | |
Array of integers. | |
Array of booleans. | |
Array of doubles. | |
Array of strings. | |
The PHOEBE_curve type. | |
The PHOEBE_spectrum type. | |
The PHOEBE_minimizer_feedback type. | |
Can be any PHOEBE_type. |
typedef struct { int dim; double *val; } PHOEBE_vector;
int |
Vector dimension. |
double * |
An array of vector values. |
PHOEBE_vector* phoebe_vector_new ();
Initializes a vector for allocation.
Returns : |
A PHOEBE_vector. |
PHOEBE_vector* phoebe_vector_new_from_qualifier (char *qualifier);
Allocates a new PHOEBE_vector and assigns it elements
taken from the array parameter represented by qualifier
.
|
The array parameter. |
Returns : |
A PHOEBE_vector, or NULL if an error occured. |
PHOEBE_vector* phoebe_vector_new_from_column (char *filename, int col);
Reads in the col
-th column from file filename
,
parses it and stores it into the returned vector.
|
The full path to the file to be read. |
|
The column to be read. |
Returns : |
A PHOEBE_vector, or NULL if an error occured. |
PHOEBE_vector* phoebe_vector_new_from_range (int dim, double start, double end);
Creates a vector that starts at start
, ends at end
, and has dim
equidistant steps.
|
vector dimension |
|
first vector element |
|
last vector element |
Returns : |
PHOEBE_vector* phoebe_vector_duplicate (PHOEBE_vector *vec);
Makes a duplicate copy of PHOEBE_vector vec
.
|
The PHOEBE_vector to copy. |
Returns : |
A PHOEBE_vector. |
int phoebe_vector_alloc (PHOEBE_vector *vec, int dimension);
Allocates storage memory for a PHOEBE_vector of dimension
.
|
The PHOEBE_vector to store. |
|
The size of the new PHOEBE_vector. |
Returns : |
A PHOEBE_error_code indicating the success of the operation. Possible errors: ERROR_VECTOR_ALREADY_ALLOCATED and ERROR_VECTOR_INVALID_DIMENSION. |
int phoebe_vector_realloc (PHOEBE_vector *vec, int dimension);
Reallocates storage memory for a PHOEBE_vector of dimension
.
|
PHOEBE_vector to reallocate. |
|
new size for vec .
|
Returns : |
A PHOEBE_error_code indicating the success of the operation. Possible errors: ERROR_VECTOR_INVALID_DIMENSION. |
int phoebe_vector_pad (PHOEBE_vector *vec, double value);
Pads all components of vec
with value
.
|
The PHOEBE_vector to pad with value .
|
|
The new value for all elements of vec .
|
Returns : |
A PHOEBE_error_code indicating the success of the operation. |
int phoebe_vector_free (PHOEBE_vector *vec);
Frees the storage memory allocated for PHOEBE_vector vec
.
|
The PHOEBE_vector to free. |
Returns : |
A PHOEBE_error_code indicating the success of the operation. |
int phoebe_vector_add (PHOEBE_vector *result, PHOEBE_vector *fac1, PHOEBE_vector *fac2);
Adds PHOEBE_vector fac1
to fac2
and returns the sum PHOEBE_vector result
.
|
The placeholder for a return value. |
|
Vector 1. |
|
Vector 2. |
Returns : |
PHOEBE_error_code |
int phoebe_vector_subtract (PHOEBE_vector *result, PHOEBE_vector *fac1, PHOEBE_vector *fac2);
Subtracts PHOEBE_vector fac2
from fac1
and returns the difference
PHOEBE_vector result
.
|
The placeholder for a return value. |
|
Vector 1. |
|
Vector 2. |
Returns : |
A PHOEBE_error_code indicating the success of the operation. Possible errors: ERROR_VECTOR_DIMENSIONS_MISMATCH. |
int phoebe_vector_multiply (PHOEBE_vector *result, PHOEBE_vector *fac1, PHOEBE_vector *fac2);
Multiplies PHOEBE_vector fac1
with fac2
and returns the product
PHOEBE_vector result
.
|
The placeholder for a return value. |
|
Vector 1. |
|
Vector 2. |
Returns : |
A PHOEBE_error_code indicating the success of the operation. Possible errors: ERROR_VECTOR_DIMENSIONS_MISMATCH. |
int phoebe_vector_divide (PHOEBE_vector *result, PHOEBE_vector *fac1, PHOEBE_vector *fac2);
Divides PHOEBE_vector fac1
with fac2
and returns the quotient
PHOEBE_vector result
.
|
The placeholder for a return value. |
|
Vector 1. |
|
Vector 2. |
Returns : |
A PHOEBE_error_code indicating the success of the operation. Possible errors: ERROR_VECTOR_DIMENSIONS_MISMATCH. |
int phoebe_vector_raise (PHOEBE_vector *result, PHOEBE_vector *fac1, PHOEBE_vector *fac2);
Raises all elements of PHOEBE_vector fac1
to the corresponding element
of PHOEBE_vector fac2
, basically result[i] = fac1[i]^fac2[i].
|
The placeholder for a return value. |
|
Vector 1. |
|
Vector 2. |
Returns : |
A PHOEBE_error_code indicating the success of the operation. Possible errors: ERROR_VECTOR_DIMENSIONS_MISMATCH. |
int phoebe_vector_offset (PHOEBE_vector *vec, double offset);
Adds offset
to all elements of PHOEBE_vector vec
.
|
PHOEBE_vector to be offset |
|
offset value. |
Returns : |
PHOEBE_error_code. |
int phoebe_vector_sum (PHOEBE_vector *vec, double *sum);
Computes a sum of all vector elements.
|
PHOEBE_vector for which a sum is computed |
|
|
Returns : |
PHOEBE_error_code. |
int phoebe_vector_mean (PHOEBE_vector *vec, double *mean);
Computes the mean (arithmetic average) of all vector elements.
|
PHOEBE_vector for which a mean is computed |
|
|
Returns : |
PHOEBE_error_code. |
int phoebe_vector_median (PHOEBE_vector *vec, double *median);
Sorts the vector elements and returns a median value. The passed
vector vec
is not changed.
|
PHOEBE_vector for which a median is computed |
|
placeholder for the median value |
Returns : |
PHOEBE_error_code. |
int phoebe_vector_standard_deviation (PHOEBE_vector *vec, double *sigma);
Computes standard deviation (sigma
) of vector elements.
|
PHOEBE_vector for which standard deviation is computed |
|
placeholder for standard deviation value |
Returns : |
PHOEBE_error_code. |
int phoebe_vector_multiply_by (PHOEBE_vector *vec, double factor);
Multiplies all elements of the PHOEBE_vector vec
with the scalar
value factor
.
|
PHOEBE_vector to be modified |
|
scaling factor |
Returns : |
PHOEBE_error_code. |
int phoebe_vector_dot_product (double *result, PHOEBE_vector *fac1, PHOEBE_vector *fac2);
Returns the scalar (dot) product of two PHOEBE_vectors.
|
The placeholder for a return value. |
|
Vector 1. |
|
Vector 2. |
Returns : |
A PHOEBE_error_code indicating the success of the operation. Possible errors: ERROR_VECTOR_DIMENSIONS_MISMATCH. |
int phoebe_vector_vec_product (PHOEBE_vector *result, PHOEBE_vector *fac1, PHOEBE_vector *fac2);
Returns the vector product of the two PHOEBE_vectors. Both fac1
and fac2
need to have exactly 3 elements.
|
The placeholder for a return value. |
|
Vector 1. |
|
Vector 2. |
Returns : |
A PHOEBE_error_code indicating the success of the operation. Possible errors: ERROR_VECTOR_DIMENSION_NOT_THREE, and ERROR_VECTOR_DIMENSIONS_MISMATCH. |
int phoebe_vector_submit (PHOEBE_vector *result, PHOEBE_vector *vec);
Calculates the functional value of func
for each element of vec
.
|
The placeholder for a return value. |
|
The PHOEBE_vector to submit to func .
|
Returns : |
A PHOEBE_error_code indicating the success of the operation. |
int phoebe_vector_norm (double *result, PHOEBE_vector *vec);
Returns the norm of PHOEBE_vector vec
.
|
The placeholder for a return value. |
|
The PHOEBE_vector to norm. |
Returns : |
A PHOEBE_error_code indicating the success of the operation. |
int phoebe_vector_dim (int *result, PHOEBE_vector *vec);
Returns the dimension of vec
.
|
The placeholder for a return value. |
|
The PHOEBE_vector to examine. |
Returns : |
A PHOEBE_error_code indicating the success of the operation. |
int phoebe_vector_randomize (PHOEBE_vector *result, double limit);
Fills all vector elements with random numbers from the
interval [0, limit
]. limit
may also be negative, then it's [limit
, 0].
The vector must be allocated prior to calling this function.
|
The placeholder for a return value. |
|
The upper limit for generated numbers. |
Returns : |
A PHOEBE_error_code indicating the success of the operation. Possible errors: ERROR_VECTOR_IS_EMPTY. |
int phoebe_vector_min_max (PHOEBE_vector *vec, double *min, double *max);
Scans through the dataset of and assigns the minimal and the
maximal value encountered to min
and max
variables.
|
The PHOEBE_vector to scan for minimum and maximum. |
|
The placeholder for the minimal value in vec .
|
|
The placeholder for the maximal value in vec .
|
Returns : |
A PHOEBE_error_code indicating the success of the operation. |
int phoebe_vector_min_index (PHOEBE_vector *vec, int *index);
Scans through the PHOEBE_vector vec
and assigns the index of
the lowest value in vec
to index
.
|
The PHOEBE_vector to scan for minimum. |
|
The placeholder for the index of the minimal value in vec .
|
Returns : |
A PHOEBE_error_code indicating the success of the operation. |
int phoebe_vector_max_index (PHOEBE_vector *vec, int *index);
Scans through the PHOEBE_vector vec
and assigns the index of
the highest value in vec
to index
.
|
The PHOEBE_vector to scan for maximum. |
|
The placeholder for the index of the maximal value in vec .
|
Returns : |
A PHOEBE_error_code indicating the success of the operation. |
int phoebe_vector_rescale (PHOEBE_vector *vec, double ll, double ul);
Linearly rescales the values of elements in the vector vec
to the
[ll
, ul
] interval.
|
The PHOEBE_vector to rescale. |
|
The lower limit for rescaling. |
|
The upper limit for rescaling. |
Returns : |
PHOEBE_error_code. |
bool phoebe_vector_compare (PHOEBE_vector *vec1, PHOEBE_vector *vec2);
Compares two passed PHOEBE_vectors. It returns TRUE if all vector elements are the same; it returns FALSE otherwise. The comparison is done by comparing the difference of both elements to PHOEBE_NUMERICAL_ACCURACY.
|
Vector 1. |
|
Vector 2. |
Returns : |
A boolean indicating whether the two vectors have identical elements. |
int phoebe_vector_less_than (bool *result, PHOEBE_vector *vec1, PHOEBE_vector *vec2);
Tests whether *all* vector elements of vec1
are less
than their respective counterparts of vec2
. If so, TRUE is returned,
otherwise FALSE is returned.
|
The placeholder for the result. |
|
Vector 1. |
|
Vector 2. |
Returns : |
A PHOEBE_error_code indicating the success of the operation. Possible errors: ERROR_VECTOR_NOT_INITIALIZED and ERROR_VECTOR_DIMENSIONS_MISMATCH. |
int phoebe_vector_leq_than (bool *result, PHOEBE_vector *vec1, PHOEBE_vector *vec2);
Tests whether *all* vector elements of vec1
are less or
equal to their respective counterparts of vec2
. If so, TRUE is returned,
otherwise FALSE is returned.
|
The placeholder for the result. |
|
Vector 1. |
|
Vector 2. |
Returns : |
A PHOEBE_error_code indicating the success of the operation. Possible errors: ERROR_VECTOR_NOT_INITIALIZED and ERROR_VECTOR_DIMENSIONS_MISMATCH. |
int phoebe_vector_greater_than (bool *result, PHOEBE_vector *vec1, PHOEBE_vector *vec2);
Tests whether *all* vector elements of vec1
are greater
than their respective counterparts of vec2
. If so, TRUE is returned,
otherwise FALSE is returned.
|
The placeholder for the result. |
|
Vector 1. |
|
Vector 2. |
Returns : |
A PHOEBE_error_code indicating the success of the operation. Possible errors: ERROR_VECTOR_NOT_INITIALIZED and ERROR_VECTOR_DIMENSIONS_MISMATCH. |
int phoebe_vector_geq_than (bool *result, PHOEBE_vector *vec1, PHOEBE_vector *vec2);
Tests whether *all* vector elements of vec1
are greater
or equal to their respective counterparts of vec2
. If so, TRUE is
returned, otherwise FALSE is returned.
|
The placeholder for the result. |
|
Vector 1. |
|
Vector 2. |
Returns : |
A PHOEBE_error_code indicating the success of the operation. Possible errors: ERROR_VECTOR_NOT_INITIALIZED and ERROR_VECTOR_DIMENSIONS_MISMATCH. |
int phoebe_vector_append_element (PHOEBE_vector *vec, double val);
Appends an element with value val
to the PHOEBE_vector vec
.
|
The vector to extend. |
|
The value to append. |
Returns : |
A PHOEBE_error_code indicating the success of the operation. |
int phoebe_vector_remove_element (PHOEBE_vector *vec, int index);
Removes the index
-th element from PHOEBE_vector vec
.
|
The vector to shorten. |
|
The index of the element to be removed. |
Returns : |
A PHOEBE_error_code indicating the success of the operation. Possible errors: ERROR_INDEX_OUT_OF_RANGE. |
typedef struct { int rows; int cols; double **val; } PHOEBE_matrix;
int |
The horizontal dimension of the matrix. |
int |
The vertical dimension of the matrix. |
double ** |
The elements of the matrix. |
PHOEBE_matrix* phoebe_matrix_new ();
Initializes an unallocated PHOEBE_matrix.
Returns : |
An empty PHOEBE_matrix. |
int phoebe_matrix_alloc (PHOEBE_matrix *matrix, int cols, int rows);
Allocates storage memory for matrix
. The
subscripts are such that the elements of matrix
should be accessed
by matrix
[row
][col
].
|
A PHOEBE_matrix to allocate. |
|
The number of columns for matrix .
|
|
The number of rows for matrix .
|
Returns : |
A PHOEBE_error_code indicating the success of the operation. Possible errors: ERROR_MATRIX_ALREADY_ALLOCATED and ERROR_MATRIX_INVALID_DIMENSION. |
int phoebe_matrix_free (PHOEBE_matrix *matrix);
Frees the storage memory allocated for matrix
.
|
The PHOEBE_matrix to be freed. |
Returns : |
A PHOEBE_error_code indicating the success of the operation. |
int phoebe_matrix_get_row (PHOEBE_vector *vec, PHOEBE_matrix *matrix, int row);
Copies the contents of row
of matrix
to vec
, assuming that vec
is already
allocated.
|
The PHOEBE_vector to store the values from row .
|
|
The PHOEBE_matrix to read row from.
|
|
The row of matrix to be copied to vec .
|
Returns : |
A PHOEBE_error_code indicating the success of the operation. |
int phoebe_matrix_set_row (PHOEBE_matrix *matrix, PHOEBE_vector *vec, int row);
Sets the elements of row
of matrix
to the values of elements in vec
.
|
A PHOEBE_matrix to modify. |
|
The PHOEBE_vector to be placed in matrix .
|
|
The row of matrix to be replaced with vec .
|
Returns : |
A PHOEBE_error_code indicating the success of the operation. |
typedef struct { int dim; PHOEBE_type type; union { int *iarray; double *darray; bool *barray; char **strarray; } val; } PHOEBE_array;
int |
The size of the array. |
PHOEBE_type |
The type of the array. |
PHOEBE_array* phoebe_array_new_from_qualifier (char *qualifier);
|
|
Returns : |
PHOEBE_array* phoebe_array_new_from_column (char *filename, int col);
Reads in the col
-th column from file filename
, parses it and stores
it into the returned array. The first element determines the type of
the array.
|
absolute path to the file to be read. |
|
column to be read. |
Returns : |
PHOEBE_array, or NULL if an error occured. |
int phoebe_array_alloc (PHOEBE_array *array, int dimension);
|
|
|
|
Returns : |
int phoebe_array_realloc (PHOEBE_array *array, int dimension);
|
|
|
|
Returns : |
PHOEBE_array* phoebe_array_duplicate (PHOEBE_array *array);
|
|
Returns : |
bool phoebe_array_compare (PHOEBE_array *array1, PHOEBE_array *array2);
Compares two passed PHOEBE_arrays. It returns TRUE if all array elements are the same; it returns FALSE otherwise. In cases of double precision numbers the comparison is done by comparing the difference of both elements to PHOEBE_NUMERICAL_ACCURACY.
|
Array 1. |
|
Array 2. |
Returns : |
Boolean indicating whether the two arrays have identical elements. |
PHOEBE_vector* phoebe_vector_new_from_array (PHOEBE_array *array);
Converts the array of doubles into PHOEBE_vector.
|
PHOEBE_array of type TYPE_INT or TYPE_DOUBLE |
Returns : |
PHOEBE_vector, or NULL in case of failure. |
typedef struct { int bins; double *range; double *val; } PHOEBE_hist;
int |
Number of histogram bins. |
double * |
A vector of histogram ranges. |
double * |
A vector of histogram values. |
typedef enum PHOEBE_hist_rebin_type { PHOEBE_HIST_CONSERVE_VALUES = 20, PHOEBE_HIST_CONSERVE_DENSITY } PHOEBE_hist_rebin_type;
There are two ways to rebin a histogram:
1) conserve the values and 2) conserve the value densities.
The first option is better if we are degrading the histogram, and the second option is better if we are oversampling the histogram.
PHOEBE_hist* phoebe_hist_new_from_arrays (int bins, double *binarray, double *valarray);
|
|
|
|
|
|
Returns : |
PHOEBE_hist* phoebe_hist_new_from_file (char *filename);
|
|
Returns : |
int phoebe_hist_realloc (PHOEBE_hist *hist, int bins);
|
|
|
|
Returns : |
int phoebe_hist_set_ranges (PHOEBE_hist *hist, PHOEBE_vector *bin_centers);
|
|
|
|
Returns : |
int phoebe_hist_set_values (PHOEBE_hist *hist, PHOEBE_vector *values);
Sets bin values to the passed vector elements. The histogram hist
must be initialized and allocated, and the dimension of values
must
match the number of bins.
|
histogram to be modified |
|
a vector of values to be copied to the histogram |
Returns : |
PHOEBE_error_code. |
int phoebe_hist_get_bin_centers (PHOEBE_hist *hist, PHOEBE_vector *bin_centers);
|
|
|
|
Returns : |
int phoebe_hist_get_bin (int *bin, PHOEBE_hist *hist, double r);
|
|
|
|
|
|
Returns : |
int phoebe_hist_evaluate (double *y, PHOEBE_hist *hist, double x);
|
|
|
|
|
|
Returns : |
int phoebe_hist_integrate (double *integral, PHOEBE_hist *hist, double ll, double ul);
This function integrates the passed histogram hist
from ll
to ul
by a simple rectangular scheme.
|
placeholder for the integral value |
|
histogram to be integrated |
|
lower integral limit |
|
upper integral limit |
Returns : |
PHOEBE_error_code. |
int phoebe_hist_shift (PHOEBE_hist *hist, double shift);
Shifts the contents of the histogram hist
in pixel-space by the passed
shift
. If shift
is positive, the contents are shifted to the right;
if it is negative, they are shifted to the left. The bin structure is
retained, the bins outside the shifted range are padded with 0.
|
histogram to be shifted |
|
pixel-space shift |
Returns : |
PHOEBE_error_code. |
int phoebe_hist_correlate (double *cfval, PHOEBE_hist *h1, PHOEBE_hist *h2, double sigma1, double sigma2, double ll, double ul, double xi);
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
Returns : |
int phoebe_hist_crop (PHOEBE_hist *hist, double ll, double ul);
|
|
|
|
|
|
Returns : |
bool phoebe_hist_compare (PHOEBE_hist *hist1, PHOEBE_hist *hist2);
Compares the contents of histograms hist1
and hist2
. All ranges and
values are compared by evaluating the absolute difference between the
elements and comparing that to PHOEBE_NUMERICAL_ACCURACY.
int phoebe_hist_resample (PHOEBE_hist *out, PHOEBE_hist *in, PHOEBE_hist_rebin_type type);
|
|
|
|
|
|
Returns : |
int phoebe_hist_rebin (PHOEBE_hist *out, PHOEBE_hist *in, PHOEBE_hist_rebin_type type);
|
|
|
|
|
|
Returns : |
typedef struct { char *set; char *name; char *reftable; PHOEBE_vector *lin_x; PHOEBE_vector *log_x; PHOEBE_vector *log_y; PHOEBE_vector *sqrt_x; PHOEBE_vector *sqrt_y; } PHOEBE_ld;
typedef struct { int id; char *set; char *name; double effwl; PHOEBE_hist *tf; PHOEBE_ld *ld; } PHOEBE_passband;
int |
ID number of the passband |
char * |
Filter-set of the passband, i.e. "Cousins" |
char * |
Passband identifier, i.e. "Rc" |
double |
Effective wavelength of the passband. |
PHOEBE_hist * |
Passband transmission function |
PHOEBE_ld * |
Limb darkening table (attached optionally) |
typedef enum PHOEBE_curve_type { PHOEBE_CURVE_UNDEFINED, PHOEBE_CURVE_LC, PHOEBE_CURVE_RV } PHOEBE_curve_type;
int phoebe_curve_type_get_name (PHOEBE_curve_type ctype, char **name);
|
|
|
|
Returns : |
typedef enum PHOEBE_column_type { PHOEBE_COLUMN_UNDEFINED, PHOEBE_COLUMN_HJD, PHOEBE_COLUMN_PHASE, PHOEBE_COLUMN_MAGNITUDE, PHOEBE_COLUMN_FLUX, PHOEBE_COLUMN_PRIMARY_RV, PHOEBE_COLUMN_SECONDARY_RV, PHOEBE_COLUMN_SIGMA, PHOEBE_COLUMN_WEIGHT, PHOEBE_COLUMN_INVALID } PHOEBE_column_type;
Various sorts of data that can be found in files usable by Phoebe.
int phoebe_column_type_get_name (PHOEBE_column_type ctype, char **name);
|
|
|
|
Returns : |
int phoebe_column_get_type (PHOEBE_column_type *type, const char *string);
Parses the passed string
and converts it to the enumerated type of
the column (PHOEBE_column_type).
|
placeholder for the enumerated column type |
|
string representation of the column type |
Returns : |
PHOEBE_error_code. |
typedef enum PHOEBE_data_flag { PHOEBE_DATA_REGULAR, PHOEBE_DATA_ALIASED, PHOEBE_DATA_DELETED, PHOEBE_DATA_OMITTED } PHOEBE_data_flag;
typedef struct { PHOEBE_curve_type type; PHOEBE_passband *passband; PHOEBE_vector *indep; PHOEBE_vector *dep; PHOEBE_vector *weight; PHOEBE_array *flag; PHOEBE_column_type itype; PHOEBE_column_type dtype; PHOEBE_column_type wtype; char *filename; double sigma; } PHOEBE_curve;
PHOEBE_curve_type |
Type of the curve. |
PHOEBE_passband * |
Passband of the curve. |
PHOEBE_vector * |
Elements of the independant variable vector. |
PHOEBE_vector * |
Elements of the dependant variable vector. |
PHOEBE_vector * |
Elements of the weight vector. |
PHOEBE_array * |
data flag of the enumerated PHOEBE_data_flag type |
PHOEBE_column_type |
Column type of the independant variable. |
PHOEBE_column_type |
Column type of the dependant variable. |
PHOEBE_column_type |
Column type of the weights. |
char * |
Absolute path to the file containing the curve. |
double |
Sigma value of the curve. |
PHOEBE_curve* phoebe_curve_new_from_file (char *filename);
|
|
Returns : |
PHOEBE_curve* phoebe_curve_new_from_pars (PHOEBE_curve_type ctype, int index);
|
|
|
|
Returns : |
PHOEBE_curve* phoebe_curve_duplicate (PHOEBE_curve *curve);
|
|
Returns : |
int phoebe_curve_alloc (PHOEBE_curve *curve, int dim);
|
|
|
|
Returns : |
int phoebe_curve_realloc (PHOEBE_curve *curve, int dim);
Reallocates storage memory for PHOEBE_curve curve
.
|
PHOEBE_curve to reallocate. |
|
the new size of curve .
|
Returns : |
PHOEBE_error_code. Error codes: ERROR_CURVE_INVALID_DIMENSION SUCCESS |
int phoebe_curve_compute (PHOEBE_curve *curve, PHOEBE_vector *nodes, int index, PHOEBE_column_type itype, PHOEBE_column_type dtype);
Computes the index
-th model light curve or RV curve in nodes
.
The computation is governed by the enumerated choices of itype
and dtype
.
|
a pointer to the initialized PHOEBE_curve |
|
a vector of nodes in which the curve should be computed |
|
curve index |
|
requested independent data type (see PHOEBE_column_type) |
|
requested dependent data type (see PHOEBE_column_type) |
Returns : |
PHOEBE_error_code. |
int phoebe_curve_transform (PHOEBE_curve *curve, PHOEBE_column_type itype, PHOEBE_column_type dtype, PHOEBE_column_type wtype);
Transforms curve columns to requested types.
|
PHOEBE_curve to be transformed |
|
requested independent variable (time or phase) |
|
requested dependent variable (flux, magnitude, RV) |
|
requested weight variable (standard weight, standard deviation, none) |
Returns : |
PHOEBE_error_code. |
int phoebe_curve_alias (PHOEBE_curve *curve, double phmin, double phmax);
This function redimensiones the array of data phases by aliasing points to outside the [-0.5, 0.5] range. If the new interval is narrower, the points are omitted, otherwise they are aliased.
|
curve to be aliased |
|
start phase |
|
end phase |
Returns : |
PHOEBE_error_code. |
int phoebe_curve_set_properties (PHOEBE_curve *curve, PHOEBE_curve_type type, char *filename, PHOEBE_passband *passband, PHOEBE_column_type itype, PHOEBE_column_type dtype, PHOEBE_column_type wtype, double sigma);
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
Returns : |
typedef enum PHOEBE_spectrum_dispersion { PHOEBE_SPECTRUM_DISPERSION_LINEAR, PHOEBE_SPECTRUM_DISPERSION_LOG, PHOEBE_SPECTRUM_DISPERSION_NONE } PHOEBE_spectrum_dispersion;
Spectrum dispersion tells us how dx is connected to dlambda. If the dispersion is linear, then dx is proportional to dlambda. If it is log, dx is proportional to dlambda/lambda. If there is no dispersion function, that means that there is no simple transformation from wavelength space to pixel space. In that case everything must be done according to histogram ranges.
typedef struct { double R; double Rs; PHOEBE_spectrum_dispersion disp; PHOEBE_hist *data; } PHOEBE_spectrum;
The spectrum structure is defined here, but the manipulation functions reside in phoebe_spectra files.
typedef enum PHOEBE_minimizer_type { PHOEBE_MINIMIZER_NMS, PHOEBE_MINIMIZER_DC } PHOEBE_minimizer_type;
The available minimizers.
int phoebe_minimizer_type_get_name (PHOEBE_minimizer_type minimizer, char **name);
|
|
|
|
Returns : |
typedef struct { PHOEBE_minimizer_type algorithm; bool converged; double cputime; int iters; double cfval; PHOEBE_array *qualifiers; PHOEBE_vector *initvals; PHOEBE_vector *newvals; PHOEBE_vector *ferrors; PHOEBE_vector *chi2s; PHOEBE_vector *wchi2s; PHOEBE_matrix *cormat; } PHOEBE_minimizer_feedback;
PHOEBE_minimizer_type |
Minimizer (algorithm) type |
bool |
|
double |
CPU time required for algorithm execution |
int |
Number of performed iterations |
double |
Cost function value (combined chi2) |
PHOEBE_array * |
A list of TBA qualifiers |
PHOEBE_vector * |
A list of initial parameter values |
PHOEBE_vector * |
A list of new parameter values |
PHOEBE_vector * |
A list of formal error estimates |
PHOEBE_vector * |
A list of passband chi2 values |
PHOEBE_vector * |
A list of weighted passband chi2 values |
PHOEBE_matrix * |
Correlation matrix |
PHOEBE_minimizer_feedback* phoebe_minimizer_feedback_new ();
Returns : |
PHOEBE_minimizer_feedback* phoebe_minimizer_feedback_duplicate (PHOEBE_minimizer_feedback *feedback);
|
|
Returns : |
int phoebe_minimizer_feedback_alloc (PHOEBE_minimizer_feedback *feedback, int tba, int cno);
|
|
|
|
|
|
Returns : |
int phoebe_minimizer_feedback_accept (PHOEBE_minimizer_feedback *feedback);
Traverses through all the parameters stored in the feedback structure and copies the values to the currently active parameter table. After all the values have been updated, the function satisfies all constraints as well.
|
minimizer feedback with new values of parameters. |
Returns : |
PHOEBE_error_code. |
int phoebe_minimizer_feedback_free (PHOEBE_minimizer_feedback *feedback);
|
|
Returns : |
PHOEBE_spectrum* phoebe_spectrum_duplicate (PHOEBE_spectrum *spectrum);
Makes a duplicate copy of spectrum
.
|
|
Returns : |
PHOEBE_value phoebe_value_duplicate (PHOEBE_type type, PHOEBE_value val);
|
|
|
|
Returns : |