62#define COLUMN_ORDER1 "Order1"
63#define COLUMN_ORDER2 "Order2"
64#define COLUMN_COEFF "Coeff"
110 cpl_matrix * product;
111 const double * ai = cpl_matrix_get_data_const(self);
114 const size_t m = cpl_matrix_get_nrow(self);
115 const size_t n = cpl_matrix_get_ncol(self);
119 cpl_ensure(self != NULL, CPL_ERROR_NULL_INPUT, NULL);
125 product = cpl_matrix_new(
m,
m);
126 bwrite = cpl_matrix_get_data(product);
128 bwrite = (
double *) cpl_malloc((
size_t)
m * (size_t)
m *
sizeof(
double));
129 product = cpl_matrix_wrap(
m,
m, bwrite);
133 for (i = 0; i <
m; i++, bwrite +=
m, ai +=
n) {
135 for (j = i; j <
m; j++, aj +=
n) {
137 for (k = 0; k <
n; k++) {
138 sum += (
long double)ai[k] * (
long double)aj[k];
153 const cpl_matrix *rhs)
156 cpl_matrix * solution;
160 cpl_ensure(coeff != NULL, CPL_ERROR_NULL_INPUT, NULL);
161 cpl_ensure(rhs != NULL, CPL_ERROR_NULL_INPUT, NULL);
164 At = cpl_matrix_transpose_create(coeff);
165 solution = cpl_matrix_product_create(At, rhs);
169 cpl_matrix_delete(At);
175 if (cpl_matrix_decomp_chol(AtA) != CPL_ERROR_NONE ||
176 cpl_matrix_solve_chol(AtA, solution) != CPL_ERROR_NONE) {
177 cpl_matrix_delete(solution);
180 (void)cpl_error_set_where(cpl_func);
183 cpl_matrix_delete(AtA);
213 assure(pol != NULL, CPL_ERROR_ILLEGAL_INPUT,
"Null polynomial");
219 check_msg( p->
dimension = cpl_polynomial_get_dimension(pol),
"Error reading dimension");
235 check_msg( p->
pol = cpl_polynomial_duplicate(pol),
"Error copying polynomial");
238 if (cpl_error_get_code() != CPL_ERROR_NONE)
257 cpl_polynomial *p = NULL;
259 assure( dim >= 1, CPL_ERROR_ILLEGAL_INPUT,
"Illegal dimension: %d", dim);
261 p = cpl_polynomial_new(dim);
298 if (*p == NULL)
return;
299 cpl_polynomial_delete((*p)->pol);
300 cpl_vector_delete((*p)->vec);
301 cpl_free((*p)->shift);
302 cpl_free((*p)->scale);
318 assure( p != NULL, CPL_ERROR_NULL_INPUT,
"Null polynomial");
320 result = cpl_polynomial_get_degree(p->
pol);
340 assure( p != NULL, CPL_ERROR_NULL_INPUT,
"Null polynomial");
344 "Error allocating polynomial");
353 if (cpl_error_get_code() != CPL_ERROR_NONE)
383 assure( p != NULL, CPL_ERROR_NULL_INPUT,
"Null polynomial");
385 CPL_ERROR_ILLEGAL_INPUT,
"Polynomial must be 2D");
387 degree = cpl_polynomial_get_degree(p->
pol);
394 cpl_table_new_column(t,
COLUMN_COEFF , CPL_TYPE_DOUBLE);
425 for (i = 0; i <=
degree; i++){
426 for (j = 0; j+i <=
degree; j++){
432 coeff = cpl_polynomial_get_coeff(p->
pol, power);
459 cpl_polynomial *pol = NULL;
464 check_msg( pol = cpl_polynomial_new(2),
"Error initializing polynomial");
467 assure(t != NULL, CPL_ERROR_NULL_INPUT,
"Null table");
476 assure(type == CPL_TYPE_INT , CPL_ERROR_INVALID_TYPE,
481 assure(type == CPL_TYPE_INT , CPL_ERROR_INVALID_TYPE,
486 assure(type == CPL_TYPE_DOUBLE, CPL_ERROR_INVALID_TYPE,
487 "Column '%s' has type %s. Double expected",
COLUMN_COEFF ,
490 assure(cpl_table_get_nrow(t) > 1 + 2 + 1 + 2, CPL_ERROR_ILLEGAL_INPUT,
491 "Table must contain at least one coefficient");
494 for(i = 3 + 3; i < cpl_table_get_nrow(t); i++) {
500 coeff = cpl_table_get_double(t,
COLUMN_COEFF , i, NULL)),
501 "Error reading table row %d", i);
503 xsh_msg_debug(
"Pol.coeff.(%" CPL_SIZE_FORMAT
", %" CPL_SIZE_FORMAT
") = %e", power[0], power[1], coeff);
505 check_msg( cpl_polynomial_set_coeff(pol, power, coeff),
"Error creating polynomial");
519 if (cpl_error_get_code() != CPL_ERROR_NONE)
537 assure(p != NULL, CPL_ERROR_ILLEGAL_INPUT,
"Null polynomial");
558 fprintf(stream,
"Null polynomial\n");
561 cpl_polynomial_dump(p->
pol, stream);
562 fprintf(stream,
"shift_y \t= %f \tscale_y \t= %f\n", p->
shift[0], p->
scale[0]);
565 fprintf(stream,
"shift_x%d \t= %f \tscale_x%d \t= %f\n",
590 assure(p != NULL, CPL_ERROR_NULL_INPUT,
"Null polynomial");
592 CPL_ERROR_ILLEGAL_INPUT,
"Illegal variable number: %d", varno);
606 p->
shift[varno] *= scale;
607 p->
scale[varno] *= scale;
610 return cpl_error_get_code();
631 assure(p != NULL, CPL_ERROR_NULL_INPUT,
"Null polynomial");
633 CPL_ERROR_ILLEGAL_INPUT,
"Illegal variable number: %d", varno);
644 p->
shift[varno] += shift;
647 return cpl_error_get_code();
665 assure(p != NULL, CPL_ERROR_NULL_INPUT,
"Null polynomial");
667 CPL_ERROR_ILLEGAL_INPUT,
"Polynomial must be 1d");
670 cpl_polynomial_eval_1d(p->
pol, (
x - p->
shift[1])/p->
scale[1], NULL)
672 "Could not evaluate polynomial");
696 assure(p != NULL, CPL_ERROR_NULL_INPUT,
"Null polynomial");
698 "Polynomial must be 2d. It's %dd", p->
dimension);
700 double scale = p->
scale[0];
701 double shift = p->
shift[0];
708 result = cpl_polynomial_eval(p->
pol, p->
vec) * scale + shift;
736 assure(p != NULL, CPL_ERROR_NULL_INPUT,
"Null polynomial");
738 "Polynomial must be 1d");
749 check_msg(( coeff0 = cpl_polynomial_get_coeff(p->
pol, power),
750 cpl_polynomial_set_coeff(p->
pol, power, coeff0 + (p->
shift[0] - value)/p->
scale[0])),
751 "Error setting coefficient");
754 &result, multiplicity),
"Could not find root");
756 cpl_polynomial_set_coeff(p->
pol, power, coeff0);
785 int multiplicity,
int varno,
double x_value)
790 assure( 1 <= varno && varno <= 2, CPL_ERROR_ILLEGAL_INPUT,
791 "Illegal variable number: %d", varno);
794 "Could not collapse polynomial");
797 "Could not find root");
820 assure (1 <= varno && varno <= 2, CPL_ERROR_ILLEGAL_INPUT,
821 "Illegal variable number (%d)", varno);
823 assure(p != NULL, CPL_ERROR_NULL_INPUT,
"Null polynomial");
840 int degree = cpl_polynomial_get_degree(p->
pol);
847 yj *= (varno == 1) ? x2 : x1)
860 for (i =
degree; i >= 1; i--)
864 power[0] = (varno == 1) ? i : j;
865 power[1] = (varno == 1) ? j : i;
867 c_ij = cpl_polynomial_get_coeff(p->
pol, power);
870 if (i >= 2) sum *= (varno == 1) ? x1 : x2;
878 result *= p->
scale[0];
908 assure(p != NULL, CPL_ERROR_NULL_INPUT,
"Null polynomial");
910 CPL_ERROR_ILLEGAL_INPUT,
"Polynomial must be 1d");
913 "Error evaluating derivative");
931 cpl_polynomial *pol = NULL;
933 assure(p1 != NULL, CPL_ERROR_NULL_INPUT,
"Null polynomial");
934 assure(p2 != NULL, CPL_ERROR_NULL_INPUT,
"Null polynomial");
936 CPL_ERROR_ILLEGAL_INPUT,
"Polynomial must be 2d");
938 CPL_ERROR_ILLEGAL_INPUT,
"Polynomial must be 2d");
955 pol = cpl_polynomial_new(2);
956 for (i = 0; i <=
degree; i++)
957 for (j = 0; j <=
degree; j++) {
958 double coeff1, coeff2;
967 cpl_polynomial_set_coeff(pol, power, coeff1 + coeff2);
999 assure(p != NULL, CPL_ERROR_NULL_INPUT,
"Null polynomial");
1000 dimension = cpl_polynomial_get_dimension(p);
1001 degree = cpl_polynomial_get_degree(p);
1005 "Illegal variable number: %d", varno);
1010 for(i = 0; i <=
degree; i++)
1016 coeff = cpl_polynomial_get_coeff(p, power);
1019 cpl_polynomial_set_coeff(p, power, (i+1) * coeff);
1026 for(i = 0; i <=
degree; i++)
1028 for(j = 0; i + j <=
degree; j++)
1031 power[varno - 1] = i+1;
1032 power[2 - varno] = j;
1034 coeff = cpl_polynomial_get_coeff(p, power);
1036 power[varno - 1] = i;
1038 cpl_polynomial_set_coeff(p, power, (i+1) * coeff);
1044 return cpl_error_get_code();
1063 assure( p != NULL, CPL_ERROR_NULL_INPUT,
"Null polynomial");
1066 "Illegal variable number: %d", varno);
1084 "Error calculating derivative of CPL-polynomial");
1087 return cpl_error_get_code();
1109 assure( p != NULL, CPL_ERROR_NULL_INPUT,
"Null polynomial");
1112 assure( 0 <= degree1, CPL_ERROR_ILLEGAL_INPUT,
"Illegal degree: %d", degree1);
1113 assure( 0 <= degree2, CPL_ERROR_ILLEGAL_INPUT,
"Illegal degree: %d", degree2);
1127 factorial *= degree1;
1135 factorial *= degree2;
1140 "Error evaluating polynomial");
1165 assure( p != NULL, CPL_ERROR_NULL_INPUT,
"Null polynomial");
1187 "Error evaluating polynomial");
1216 cpl_polynomial *pol = NULL;
1217 cpl_size *power = NULL;
1222 assure(p != NULL, CPL_ERROR_NULL_INPUT,
"Null polynomial");
1225 "Polynomial has non-positive dimension: %d",
dimension);
1227 "Don't collapse a 1d polynomial. Evaluate it!");
1232 assure(
dimension == 2, CPL_ERROR_ILLEGAL_INPUT,
"Polynomial must be 2d");
1235 "Wrong variable number");
1236 value = (value - p->
shift[varno]) / p->
scale[varno];
1239 degree = cpl_polynomial_get_degree(p->
pol);
1240 pol = cpl_polynomial_new(
dimension - 1);
1241 power = cpl_malloc(
sizeof(cpl_size) *
dimension);
1243 for (i = 0; i <=
degree; i++)
1251 for (j =
degree - i; j >= 0; j--)
1254 coeff += cpl_polynomial_get_coeff(p->
pol, power);
1255 if (j > 0) coeff *= value;
1259 cpl_polynomial_set_coeff(pol, power, coeff);
1283 assure(cpl_error_get_code() == CPL_ERROR_NONE, cpl_error_get_code(),
1284 "Error collapsing polynomial");
1287 cpl_free(power); power = NULL;
1289 if (cpl_error_get_code() != CPL_ERROR_NONE)
1320 const cpl_vector * x_pos,
1321 const cpl_vector * values,
1322 const cpl_vector * sigmas,
1328 cpl_matrix * ma = NULL;
1329 cpl_matrix * mb = NULL;
1330 cpl_matrix * mx = NULL;
1331 const double * x_pos_data ;
1332 const double * values_data ;
1333 const double * sigmas_data = NULL;
1334 double mean_x, mean_z;
1336 cpl_polynomial * out ;
1337 cpl_vector * x_val = NULL;
1341 assure_nomsg( x_pos != NULL && values != NULL, CPL_ERROR_NULL_INPUT);
1342 assure( poly_deg >= 0, CPL_ERROR_ILLEGAL_INPUT,
1343 "Polynomial degree is %d. Must be non-negative", poly_deg);
1344 np = cpl_vector_get_size(x_pos) ;
1347 assure( np >= nc, CPL_ERROR_ILLEGAL_INPUT,
1348 "Not enough points (%d) to fit %d-order polynomial. %d point(s) needed",
1355 ma = cpl_matrix_new(np, nc) ;
1356 mb = cpl_matrix_new(np, 1) ;
1359 mean_x = cpl_vector_get_mean(x_pos);
1360 mean_z = cpl_vector_get_mean(values);
1363 x_pos_data = cpl_vector_get_data_const(x_pos) ;
1364 values_data = cpl_vector_get_data_const(values) ;
1367 sigmas_data = cpl_vector_get_data_const(sigmas) ;
1372 for (i=0 ; i<np ; i++)
1375 if (sigmas_data[i] == 0)
1379 assure(
false, CPL_ERROR_DIVISION_BY_ZERO,
1380 "Sigmas must be non-zero");
1382 for (j=0 ; j<nc ; j++)
1384 cpl_matrix_set(ma, i, j,
1389 cpl_matrix_set(mb, i, 0, (values_data[i] - mean_z) / sigmas_data[i]);
1394 for (i=0 ; i<np ; i++)
1396 for (j=0 ; j<nc ; j++)
1398 cpl_matrix_set(ma, i, j,
1402 cpl_matrix_set(mb, i, 0, (values_data[i] - mean_z) / 1) ;
1408 "Could not invert matrix");
1413 out = cpl_polynomial_new(1) ;
1415 for (deg=0 ; deg<nc ; deg++) {
1416 cpl_polynomial_set_coeff(out, °, cpl_matrix_get(mx, deg, 0)) ;
1423 x_val = cpl_vector_new(1) ;
1424 for (i=0 ; i<np ; i++)
1427 cpl_vector_set(x_val, 0, x_pos_data[i] - mean_x) ;
1429 residual = (values_data[i] - mean_z) - cpl_polynomial_eval(out, x_val);
1430 *mse += residual*residual;
1434 *mse /= (double)np ;
1500 const cpl_vector * sigmas,
int poly_deg1,
int poly_deg2,
double * mse,
1501 double * red_chisq,
polynomial ** variance) {
1511 cpl_matrix * mat_ma;
1512 cpl_matrix * cov = NULL;
1513 const double * xy_pos_data_x;
1514 const double * xy_pos_data_y;
1515 const double * values_data;
1516 const double * sigmas_data = NULL;
1517 const cpl_vector* xy_pos_x;
1518 const cpl_vector* xy_pos_y;
1519 double mean_x, mean_y, mean_z;
1520 cpl_polynomial * out;
1521 cpl_polynomial * variance_cpl;
1527 assure(xy_pos && values, CPL_ERROR_NULL_INPUT,
"Null input");
1528 assure(poly_deg1 >= 0, CPL_ERROR_ILLEGAL_INPUT,
1529 "Polynomial degree1 is %d", poly_deg1);
1530 assure(poly_deg2 >= 0, CPL_ERROR_ILLEGAL_INPUT,
1531 "Polynomial degree2 is %d", poly_deg2);
1532 np = cpl_bivector_get_size(xy_pos);
1535 assure( (variance == NULL && red_chisq == NULL) || sigmas != NULL,
1536 CPL_ERROR_ILLEGAL_INPUT,
1537 "Cannot calculate variance or chi_sq without knowing");
1540 nc = (1 + poly_deg1) * (1 + poly_deg2);
1542 assure(np >= nc, CPL_ERROR_SINGULAR_MATRIX,
1543 "%d coefficients. Only %d points", nc, np);
1548 assure(red_chisq == NULL || np > nc, CPL_ERROR_ILLEGAL_INPUT,
1549 "%d coefficients. %d points. Cannot calculate chi square", nc, np);
1551 degx_tab = cpl_malloc(nc *
sizeof(
int));
1554 degy_tab = cpl_malloc(nc *
sizeof(
int));
1555 if (degy_tab == NULL) {
1562 for (degy = 0; degy <= poly_deg2; degy++) {
1563 for (degx = 0; degx <= poly_deg1; degx++) {
1575 ma = cpl_matrix_new(np, nc);
1576 mb = cpl_matrix_new(np, 1);
1579 xy_pos_x = cpl_bivector_get_x_const(xy_pos);
1580 xy_pos_y = cpl_bivector_get_y_const(xy_pos);
1582 mean_x = cpl_vector_get_mean(xy_pos_x);
1583 mean_y = cpl_vector_get_mean(xy_pos_y);
1584 mean_z = cpl_vector_get_mean(values);
1588 xy_pos_data_x = cpl_vector_get_data_const(xy_pos_x);
1589 xy_pos_data_y = cpl_vector_get_data_const(xy_pos_y);
1590 values_data = cpl_vector_get_data_const(values);
1591 if (sigmas != NULL) {
1592 sigmas_data = cpl_vector_get_data_const(sigmas);
1595 if (sigmas != NULL) {
1597 for (i = 0; i < np; i++) {
1598 double *ma_data = cpl_matrix_get_data(ma);
1599 double *mb_data = cpl_matrix_get_data(mb);
1605 if (sigmas_data[i] == 0) {
1610 assure(
false, CPL_ERROR_DIVISION_BY_ZERO,
1611 "Sigmas must be non-zero. sigma[%d] is %f", i, sigmas_data[i]);
1614 for (degy = 0; degy <= poly_deg2; degy++) {
1616 for (degx = 0; degx <= poly_deg1; degx++) {
1617 ma_data[j + i_nc] = valx * valy / sigmas_data[i];
1618 valx *= (xy_pos_data_x[i] - mean_x);
1621 valy *= (xy_pos_data_y[i] - mean_y);
1626 mb_data[i] = (values_data[i] - mean_z) / sigmas_data[i];
1631 for (i = 0; i < np; i++) {
1632 double *ma_data = cpl_matrix_get_data(ma);
1633 double *mb_data = cpl_matrix_get_data(mb);
1638 for (degy = 0; degy <= poly_deg2; degy++) {
1640 for (degx = 0; degx <= poly_deg1; degx++) {
1641 ma_data[j + i_nc] = valx * valy / 1;
1642 valx *= (xy_pos_data_x[i] - mean_x);
1645 valy *= (xy_pos_data_y[i] - mean_y);
1650 mb_data[i] = values_data[i] - mean_z;
1656 if (variance != NULL) {
1657 mat = cpl_matrix_transpose_create(ma);
1659 mat_ma = cpl_matrix_product_create(
mat, ma);
1660 if (mat_ma != NULL) {
1661 cov = cpl_matrix_invert_create(mat_ma);
1666 variance_cpl = cpl_polynomial_new(2);
1682 assure(
false, CPL_ERROR_ILLEGAL_OUTPUT,
"Matrix inversion failed");
1686 out = cpl_polynomial_new(2);
1687 powers = cpl_malloc(2 *
sizeof(cpl_size));
1688 if (powers == NULL) {
1699 for (i = 0; i < nc; i++) {
1700 powers[0] = degx_tab[i];
1701 powers[1] = degy_tab[i];
1702 cpl_polynomial_set_coeff(out, powers, cpl_matrix_get(mx, i, 0));
1705 if (variance != NULL &&
1706 cov != NULL && variance_cpl != NULL
1709 for (j = 0; j < nc; j++) {
1717 powers[0] = degx_tab[i] + degx_tab[j];
1718 powers[1] = degy_tab[i] + degy_tab[j];
1720 coeff = cpl_polynomial_get_coeff(variance_cpl, powers);
1721 cpl_polynomial_set_coeff(variance_cpl, powers,
1722 coeff + cpl_matrix_get(cov, i, j));
1742 if (variance != NULL) {
1757 if (mse != NULL || red_chisq != NULL) {
1762 if (red_chisq != NULL)
1764 for (i = 0; i < np; i++) {
1767 double residual = values_data[i] - regress;
1770 *mse += residual * residual;
1772 if (red_chisq != NULL) {
1773 *red_chisq += (residual / sigmas_data[i]) * (residual / sigmas_data[i]) ;
1778 *mse /= (double) np;
1780 if (red_chisq != NULL) {
1781 passure( np > nc,
"%d %d", np, nc);
1783 *red_chisq /= (double) (np - nc);
1787 cleanup:
return result;
const char * xsh_tostring_cpl_type(cpl_type t)
Convert a CPL type to a string.
#define assure(CONDITION, ERROR_CODE,...)
#define assure_nomsg(BOOL, CODE)
#define check_msg(COMMAND,...)
#define passure(CONDITION,...)
#define xsh_msg_debug(...)
Print a debug message.
polynomial * xsh_polynomial_new_zero(int dim)
Create a zero polynomial.
polynomial * xsh_polynomial_fit_2d(const cpl_bivector *xy_pos, const cpl_vector *values, const cpl_vector *sigmas, int poly_deg1, int poly_deg2, double *mse, double *red_chisq, polynomial **variance)
Fit a 2d surface with a polynomial in x and y.
cpl_matrix * xsh_matrix_solve_normal(const cpl_matrix *coeff, const cpl_matrix *rhs)
void xsh_polynomial_delete(polynomial **p)
Delete a polynomial.
double xsh_polynomial_evaluate_1d(const polynomial *p, double x)
Evaluate a 1d polynomial.
double xsh_polynomial_solve_1d(const polynomial *p, double value, double guess, int multiplicity)
Solve p(x) = value.
static cpl_error_code derivative_cpl_polynomial(cpl_polynomial *p, int varno)
Calculate the partial derivative of a CPL-polynomial.
polynomial * xsh_polynomial_collapse(const polynomial *p, int varno, double value)
Collapse a polynomial by fixing one variable to a constant.
void xsh_polynomial_delete_const(const polynomial **p)
Delete a const polynomial.
double xsh_polynomial_get_coeff_1d(const polynomial *p, int degree)
Get a coefficient of a 1D polynomial.
double xsh_polynomial_solve_2d(const polynomial *p, double value, double guess, int multiplicity, int varno, double x_value)
Solve p(x1, x2) = value.
polynomial * xsh_polynomial_fit_1d(const cpl_vector *x_pos, const cpl_vector *values, const cpl_vector *sigmas, int poly_deg, double *mse)
Fit a 1d function with a polynomial.
cpl_error_code xsh_polynomial_rescale(polynomial *p, int varno, double scale)
Rescale a polynomial.
double xsh_polynomial_derivative_1d(const polynomial *p, double x)
Evaluate the derivative of a 1d polynomial.
polynomial * xsh_polynomial_new(const cpl_polynomial *pol)
Create a polynomial.
polynomial * xsh_polynomial_add_2d(const polynomial *p1, const polynomial *p2)
Add two polynomials.
int xsh_polynomial_get_dimension(const polynomial *p)
Get the dimension of a polynomial.
int xsh_polynomial_get_degree(const polynomial *p)
Get degree.
polynomial * xsh_polynomial_convert_from_table(cpl_table *t)
Convert a table to a polynomial.
double xsh_polynomial_get_coeff_2d(const polynomial *p, int degree1, int degree2)
Get a coefficient of a 2D polynomial.
double xsh_polynomial_derivative_2d(const polynomial *p, double x1, double x2, int varno)
Evaluate the partial derivative of a 2d polynomial.
void xsh_polynomial_dump(const polynomial *p, FILE *stream)
Print a polynomial.
cpl_error_code xsh_polynomial_derivative(polynomial *p, int varno)
Calculate the partial derivative of a polynomial.
cpl_error_code xsh_polynomial_shift(polynomial *p, int varno, double shift)
Shift a polynomial.
cpl_matrix * xsh_matrix_product_normal_create(const cpl_matrix *self)
cpl_table * xsh_polynomial_convert_to_table(const polynomial *p)
Convert a polynomial to a table.
double xsh_polynomial_evaluate_2d(const polynomial *p, double x1, double x2)
Evaluate a 2d polynomial.
polynomial * xsh_polynomial_duplicate(const polynomial *p)
Copy a polynomial.
int xsh_max_int(int x, int y)
Maximum of two numbers.
void xsh_free_polynomial(cpl_polynomial **p)
Deallocate a polynomial and set the pointer to NULL.
void xsh_free_vector(cpl_vector **v)
Deallocate a vector and set the pointer to NULL.
void xsh_free_matrix(cpl_matrix **m)
Deallocate a matrix and set the pointer to NULL.
double xsh_pow_int(double x, int y)
Computes x^y.
void xsh_free(const void *mem)
Deallocate memory.