void PyGSL_multimin_function_fdf(const gsl_vector* x, void* params, double *f, gsl_vector *df) { int flag, i; PyGSL_solver *min_o; FUNC_MESS_BEGIN(); min_o = (PyGSL_solver *) params; assert(PyGSL_solver_check(min_o)); for(i = 0; i<x->size; i++){ DEBUG_MESS(2, "Got a x[%d] of %f", i, gsl_vector_get(x, i)); } assert(min_o->mstatic->n_cbs >= 3); flag = PyGSL_function_wrap_On_O(x, min_o->cbs[2], min_o->args, f, df, x->size, __FUNCTION__); DEBUG_MESS(2, "Got a result of %f", *f); for(i = 0; i<df->size; i++){ DEBUG_MESS(2, "Got df x[%d] of %f", i, gsl_vector_get(df, i)); } if (flag!= GSL_SUCCESS){ *f = gsl_nan(); if(min_o->isset == 1) longjmp(min_o->buffer,flag); } FUNC_MESS_END(); return; }
/* The Callbacks */ double PyGSL_multimin_function_f(const gsl_vector* x, void* params) { double result; int flag; int i; FUNC_MESS_BEGIN(); PyGSL_solver *min_o; min_o = (PyGSL_solver *) params; assert(PyGSL_solver_check(min_o)); for(i = 0; i<x->size; i++){ DEBUG_MESS(2, "Got a x[%d] of %f", i, gsl_vector_get(x, i)); } assert(min_o->mstatic->n_cbs >= 1); flag = PyGSL_function_wrap_On_O(x, min_o->cbs[0], min_o->args, &result, NULL, x->size, __FUNCTION__); if (flag!= GSL_SUCCESS){ result = gsl_nan(); if(min_o->isset == 1) longjmp(min_o->buffer,flag); } DEBUG_MESS(2, "Got a result of %f", result); FUNC_MESS_END(); return result; }
/*Private*/ void sdatio_append_variable(struct sdatio_file * sfile, struct sdatio_variable * svar){ int nvars; struct sdatio_variable ** new_variables; int i; nvars = sfile->n_variables + 1; new_variables = (struct sdatio_variable **) malloc(sizeof(struct sdatio_variable *)*nvars); DEBUG_MESS("Setting variable %s; %d, %d\n", svar->name, nvars, sfile->n_variables); for (i=0; i < nvars-1; i++){ DEBUG_MESS("i %d\n", i); new_variables[i] = sfile->variables[i]; } DEBUG_MESS("Setting new\n"); new_variables[nvars-1] = svar; DEBUG_MESS("About to deallocate old variables\n"); if (sfile->n_variables > 0) free(sfile->variables); sfile->n_variables = nvars; sfile->variables = new_variables; DEBUG_MESS("Deallocated old vars\n"); }
/* Private */ void sdatio_append_dimension(struct sdatio_file * sfile, struct sdatio_dimension * sdim){ int ndims; struct sdatio_dimension ** new_dimensions; int i; ndims = sfile->n_dimensions + 1; new_dimensions = (struct sdatio_dimension **) malloc(sizeof(struct sdatio_dimension *)*ndims); DEBUG_MESS("Setting dimensions; %d, %d\n", ndims, sfile->n_dimensions); for (i=0; i < ndims-1; i++){ DEBUG_MESS("i %d\n", i); new_dimensions[i] = sfile->dimensions[i]; } DEBUG_MESS("Setting new\n"); new_dimensions[ndims-1] = sdim; DEBUG_MESS("About to deallocate old dimensions\n"); if (sfile->n_dimensions > 0) free(sfile->dimensions); sfile->n_dimensions = ndims; sfile->dimensions = new_dimensions; }
static PyArrayObject * PyGSL_copy_gslmatrix_to_pyarray(const gsl_matrix *x) { int i, j; PyGSL_array_index_t dimensions[2]; PyArrayObject *a_array = NULL; double tmp; char *myptr; FUNC_MESS_BEGIN(); dimensions[0] = x->size1; dimensions[1] = x->size2; a_array = (PyArrayObject *) PyGSL_New_Array(2, dimensions, PyArray_DOUBLE); if (a_array == NULL) return NULL; for (i=0;i<dimensions[1];i++){ for (j=0;j<dimensions[0];j++){ myptr = a_array->data + a_array->strides[0] * i + a_array->strides[1] * j; tmp = gsl_matrix_get(x, j, i); *((double *) myptr) = tmp; DEBUG_MESS(3, "\t\ta_array_%d = %f\n", i, tmp); } } FUNC_MESS_END(); return a_array; }
static PyObject * PyGSL_monte_init(PyGSL_solver *self, PyObject *args) { int flag = GSL_EFAILED; monte_csys * csys; FUNC_MESS_BEGIN(); assert(PyGSL_solver_check(self)); csys = (monte_csys *)self->c_sys; assert(csys); switch(csys->type){ case PyGSL_MONTE_plain: flag = gsl_monte_plain_init(self->solver); break; case PyGSL_MONTE_miser: flag = gsl_monte_miser_init(self->solver); break; case PyGSL_MONTE_vegas: flag = gsl_monte_vegas_init(self->solver); break; default: DEBUG_MESS(2, "Monte type %d unknown",flag); PyGSL_ERROR_NULL("Unknown monte type!", GSL_ESANITY); } if(PyGSL_ERROR_FLAG(flag) != GSL_SUCCESS){ PyGSL_add_traceback(module, filename, __FUNCTION__, __LINE__); return NULL; } Py_INCREF(Py_None); FUNC_MESS_END(); return Py_None; }
/* Returns 1 if the given variable exists, 0 otherwise */ int sdatio_variable_exists(struct sdatio_file * sfile, char * variable_name){ int i; DEBUG_MESS("Finding variable in sdatio_variable_exists...\n"); for (i=0;i<sfile->n_variables;i++) if (!strcmp(sfile->variables[i]->name, variable_name)) return 1; return 0; }
static int PyGSL_stride_recalc(PyGSL_array_index_t strides, int basic_type_size, PyGSL_array_index_t * stride_recalc) { int line; FUNC_MESS_BEGIN(); line = __LINE__ + 1; if((strides % basic_type_size) == 0) { *stride_recalc = strides / basic_type_size; DEBUG_MESS(2, "\tRecalculated strides to %ld", (long)*stride_recalc); FUNC_MESS_END(); return GSL_SUCCESS; } DEBUG_MESS(2, "Failed to convert stride. %ld/%d != 0", (long)strides, basic_type_size); pygsl_error("Can not convert the stride to a GSL stride", filename, __LINE__, PyGSL_ESTRIDE); PyGSL_add_traceback(NULL, filename, __FUNCTION__, line); return PyGSL_ESTRIDE; }
void PyGSL_multimin_function_df(const gsl_vector* x, void* params, gsl_vector *df) { int flag, i; PyGSL_solver *min_o; min_o = (PyGSL_solver *) params; FUNC_MESS_BEGIN(); assert(PyGSL_solver_check(min_o)); for(i = 0; i<x->size; i++){ DEBUG_MESS(2, "Got a x[%d] of %f", i, gsl_vector_get(x, i)); } assert(min_o->mstatic->n_cbs >= 2); flag = PyGSL_function_wrap_Op_On(x, df, min_o->cbs[1], min_o->args, x->size, x->size, __FUNCTION__); for(i = 0; i<df->size; i++){ DEBUG_MESS(2, "Got df x[%d] of %f", i, gsl_vector_get(df, i)); } if(flag!=GSL_SUCCESS){ if(min_o->isset == 1) longjmp(min_o->buffer,flag); } FUNC_MESS_END(); }
void sdatio_write_variable_private(struct sdatio_file * sfile, struct sdatio_variable * svar, size_t * counts, size_t * starts, void * address){ int retval; /*if (sfile->is_parallel){}*/ /*else {*/ switch (svar->type){ case (SDATIO_INT): DEBUG_MESS("Writing an integer\n"); if ((retval = nc_put_vara_int(sfile->nc_file_id, svar->nc_id, starts, counts, address))) ERR(retval); break; case (SDATIO_FLOAT): DEBUG_MESS("Writing a float\n"); /*if ((retval = nc_put_var_double(sfile->nc_file_id, svar->nc_id, address))) ERR(retval);*/ if ((retval = nc_put_vara_float(sfile->nc_file_id, svar->nc_id, starts, counts, address))) ERR(retval); break; case (SDATIO_DOUBLE): DEBUG_MESS("Writing a double\n"); /*if ((retval = nc_put_var_double(sfile->nc_file_id, svar->nc_id, address))) ERR(retval);*/ if ((retval = nc_put_vara_double(sfile->nc_file_id, svar->nc_id, starts, counts, address))) ERR(retval); break; } /*}*/ sfile->data_written = 1; }
/* * Get the factors of FFT wavetables.... */ static PyObject * PyGSL_transform_space_get_factors(PyGSL_transform_space *self, PyGSL_transform_space *args) { PyGSL_array_index_t nf, i; long *data=NULL; size_t *cp_data=NULL; PyArrayObject * a_array = NULL; int lineno; FUNC_MESS_BEGIN(); assert(PyGSL_transform_space_check(self)); assert(self->space.v); DEBUG_MESS(2, "Type = %d", self->type); switch(self->type){ case COMPLEX_WAVETABLE: nf = self->space.cwt ->nf; cp_data = self->space.cwt ->factor; break; case REAL_WAVETABLE: nf = self->space.rwt ->nf; cp_data = self->space.rwt ->factor; break; case HALFCOMPLEX_WAVETABLE: nf = self->space.hcwt->nf; cp_data = self->space.hcwt->factor; break; case COMPLEX_WAVETABLE_FLOAT: nf = self->space.cwtf ->nf; cp_data = self->space.cwtf ->factor; break; case REAL_WAVETABLE_FLOAT: nf = self->space.rwtf ->nf; cp_data = self->space.rwtf ->factor; break; case HALFCOMPLEX_WAVETABLE_FLOAT: nf = self->space.hcwtf->nf; cp_data = self->space.hcwtf->factor; break; default: lineno = __LINE__ - 1; pygsl_error("Got unknown switch", filename, lineno, GSL_ESANITY); goto fail; break; } a_array = (PyArrayObject *) PyGSL_New_Array(1, &nf, PyArray_LONG); if(a_array == NULL){ lineno = __LINE__ - 1; goto fail; } data = (long *) a_array->data; for(i=0; i<nf; i++){ data[i] = (long) cp_data[i]; } FUNC_MESS_END(); return (PyObject *) a_array; fail: PyGSL_add_traceback(module, filename, __FUNCTION__, lineno); return NULL; }
struct sdatio_dimension * sdatio_find_dimension(struct sdatio_file * sfile, char * dimension_name){ int i, dimension_number; dimension_number = -1; DEBUG_MESS("Finding dimension...\n"); for (i=0;i<sfile->n_dimensions;i++) if (!strcmp(sfile->dimensions[i]->name, dimension_name)) dimension_number = i; if (dimension_number==-1){ printf("Couldn't find dimension %s\n", dimension_name); abort(); } return sfile->dimensions[dimension_number]; }
struct sdatio_variable * sdatio_find_variable(struct sdatio_file * sfile, char * variable_name){ int i, variable_number; variable_number = -1; DEBUG_MESS("Finding variable...\n"); for (i=0;i<sfile->n_variables;i++) if (!strcmp(sfile->variables[i]->name, variable_name)) variable_number = i; if (variable_number==-1){ printf("Couldn't find variable %s\n", variable_name); abort(); } return sfile->variables[variable_number]; }
static int PyGSL_copy_pyarray_to_gslmatrix(gsl_matrix *f, PyObject *object, PyGSL_array_index_t n, PyGSL_array_index_t p, PyGSL_error_info * info) { PyArrayObject *a_array = NULL; double tmp; char *myptr; int argnum=-1; long i, j; FUNC_MESS_BEGIN(); if (info) argnum = info->argnum; a_array = PyGSL_matrix_check(object, n, p, PyGSL_DARRAY_CINPUT(info->argnum), NULL, NULL, info); if(a_array == NULL){ FUNC_MESS(" PyGSL_PyArray_PREPARE_gsl_matrix_view failed!"); goto fail; } assert(f->size1 == (size_t) n); assert(f->size2 == (size_t) p); for (i=0;i<n;i++){ for (j=0;j<p;j++){ myptr = a_array->data + a_array->strides[0] * i + a_array->strides[1] * j; tmp = *((double *)(myptr)); DEBUG_MESS(3, "\t\ta_array[%ld,%ld] = %f\n", i, j, tmp); gsl_matrix_set(f, i, j, tmp); } } FUNC_MESS_END(); Py_DECREF(a_array); return GSL_SUCCESS; fail: PyGSL_add_traceback(NULL, filename, __FUNCTION__, __LINE__); FUNC_MESS(" Failure"); Py_XDECREF(a_array); return GSL_FAILURE; }
static PyArrayObject * PyGSL_copy_gslvector_to_pyarray(const gsl_vector *x) { PyGSL_array_index_t dimension = -1, i; PyArrayObject *a_array = NULL; double tmp; FUNC_MESS_BEGIN(); dimension = x->size; a_array = (PyArrayObject *) PyGSL_New_Array(1, &dimension, PyArray_DOUBLE); if (a_array == NULL) return NULL; for (i=0;i<dimension;i++){ tmp = gsl_vector_get(x, i); ((double *) a_array->data)[i] = tmp; DEBUG_MESS(3, "\t\ta_array_%ld = %f\n", (long)i, tmp); } FUNC_MESS_END(); return a_array; }
/* * Set a descriptive error. The callback name is listed together with the "GSL Object" * that called it, and a error description. */ static int PyGSL_copy_pyarray_to_gslvector(gsl_vector *f, PyObject *object, PyGSL_array_index_t n, PyGSL_error_info * info) { PyArrayObject *a_array = NULL; double tmp; int i, argnum = -1; FUNC_MESS_BEGIN(); if (info) argnum = info->argnum; a_array = PyGSL_vector_check(object, n, PyGSL_DARRAY_INPUT(argnum), NULL, info); if(a_array == NULL){ FUNC_MESS("PyArray_FromObject failed"); goto fail; } if(DEBUG>2){ fprintf(stderr, "\t\ta_array->dimensions[0] = %d\n", a_array->dimensions[0]); fprintf(stderr, "\t\ta_array->strides[0] = %d\n", a_array->strides[0]); } for (i=0;i<n;i++){ tmp = *((double *) (a_array->data + a_array->strides[0] * i)); gsl_vector_set(f, i, tmp); DEBUG_MESS(3, "\t\ta_array_%d = %f\n", i, tmp); } FUNC_MESS_END(); Py_DECREF(a_array); return GSL_SUCCESS; fail: PyGSL_add_traceback(NULL, filename, __FUNCTION__, __LINE__); FUNC_MESS("Failure"); Py_XDECREF(a_array); return GSL_FAILURE; }
/* make the init look similar to the one */ static void * PyGSL_gsl_monte_alloc(void * type, int n) { enum PyGSL_gsl_monte_type flag = (enum PyGSL_gsl_monte_type) type; void * result = NULL; FUNC_MESS_BEGIN(); switch(flag){ case PyGSL_MONTE_plain: result = gsl_monte_plain_alloc(n); break; case PyGSL_MONTE_miser: result = gsl_monte_miser_alloc(n); break; case PyGSL_MONTE_vegas: result = gsl_monte_vegas_alloc(n); break; default: DEBUG_MESS(2, "Monte type %d unknown",flag); PyGSL_ERROR_NULL("Unknown monte type!", GSL_ESANITY); } FUNC_MESS_END(); return result; }
static PyObject * PyGSL_monte_integrate(PyGSL_solver *self, PyObject *args) { int flag = GSL_EFAILED, line, dim; unsigned int ncalls; monte_csys * csys; gsl_rng * rng; PyObject *func=NULL, *xlo=NULL, *xuo=NULL, *rng_obj=NULL, *ncallso=NULL, *argso=NULL; PyArrayObject *xla=NULL; PyArrayObject *xua=NULL; double result, abserr; gsl_monte_function mfunc; pygsl_monte_fptr_t fptr; FUNC_MESS_BEGIN(); assert(PyGSL_solver_check(self)); if(!PyArg_ParseTuple(args, "OOOOO|O", &func, &xlo, &xuo, &ncallso, &rng_obj, &argso)){ line = __LINE__ - 1; goto fail; } if(!PyCallable_Check(func)){ line = __LINE__ - 1; pygsl_error("The function argument must be callable!", filename, line, GSL_EINVAL); goto fail; } if(PyGSL_PYLONG_TO_UINT(ncallso, &ncalls, NULL) != GSL_SUCCESS){ line = __LINE__ - 1; goto fail; } if(!PyGSL_RNG_Check(rng_obj)){ line = __LINE__ - 1; pygsl_error("The rng object must be a rng!", filename, line, GSL_EINVAL); goto fail; } rng = ((PyGSL_rng *) rng_obj)->rng; xla = PyGSL_vector_check(xlo, -1, PyGSL_DARRAY_CINPUT(2), NULL, NULL); if(xla == NULL){ line = __LINE__ - 2; goto fail; } dim = xla->dimensions[0]; xua = PyGSL_vector_check(xuo, dim, PyGSL_DARRAY_CINPUT(3), NULL, NULL); if(xua == NULL){ line = __LINE__ - 2; goto fail; } if(argso == NULL){ Py_INCREF(Py_None); argso = Py_None; } Py_XDECREF(self->cbs[0]); Py_INCREF(func); self->cbs[0] = func; self->args = argso; csys = (monte_csys *)self->c_sys; switch(csys->type){ case PyGSL_MONTE_plain: fptr = (pygsl_monte_fptr_t)&gsl_monte_plain_integrate; break; case PyGSL_MONTE_miser: fptr = (pygsl_monte_fptr_t)&gsl_monte_miser_integrate; break; case PyGSL_MONTE_vegas: fptr = (pygsl_monte_fptr_t)&gsl_monte_vegas_integrate; break; default: line = __LINE__ - 5; DEBUG_MESS(2, "Monte type %d unknown",flag); pygsl_error("Unknown monte type!", filename, line, GSL_ESANITY); goto fail; } assert(fptr); mfunc.f = &PyGSL_monte_function_wrap; mfunc.f = &PyGSL_g_test; mfunc.dim = dim; mfunc.params = self; if(setjmp(self->buffer) != 0){ self->isset = 0; line = __LINE__ - 2; goto fail; }else{ self->isset = 1; } flag = fptr(&mfunc, (double *) xla->data, (double *)xua->data, dim, ncalls, rng, self->solver, &result, &abserr); self->isset = 0; if(PyGSL_ERROR_FLAG(flag) != GSL_SUCCESS){ line = __LINE__ - 2; return NULL; } Py_DECREF(xla); xla = NULL; Py_DECREF(xua); xua = NULL; Py_DECREF(self->cbs[0]); Py_DECREF(self->args); self->cbs[0] = NULL; self->args = NULL; FUNC_MESS_END(); return Py_BuildValue("dd", result, abserr); fail: FUNC_MESS("FAIL"); Py_XDECREF(xla); Py_XDECREF(xua); Py_XDECREF(self->cbs[0]); Py_XDECREF(self->args); self->cbs[0] = NULL; self->args = NULL; PyGSL_add_traceback(module, filename, __FUNCTION__, line); return NULL; }
void sdatio_get_dimension_ids(struct sdatio_file * sfile, char * dimension_list, struct sdatio_variable * svar){ int ndims; int i,j; /*char dim_name[2];*/ char * dim_name; int * dimension_ids; int dim_name_length; int counter; int sep_size; ndims = svar->ndims; counter = 0; DEBUG_MESS("ndims %d\n", ndims); dimension_ids = (int *) malloc(sizeof(int)*ndims); for (i=0;i<ndims;i++){ /* In the next section we set counter to * the beginning of the next dimension name and * dim_name_length to the size of the name*/ if (sfile->has_long_dim_names){ sep_size = 1; dim_name_length = 0; /*The first condition checks that we haven't reached the end of the string*/ while (dimension_list[counter] && !(dimension_list[counter]==',')){ dim_name_length++; counter++; } counter++; } else { sep_size = 0; dim_name_length = 1; counter++; } if (dim_name_length < 1){ printf("ERROR: zero length dimension name in dimension_list %s\n", dimension_list); abort(); } /*dim_name[0] = dimension_list[i];*/ /*dim_name[1] = dimension_list[ndims];*/ /* Copy the name of the current dimension to the temporary * variable dim_name*/ dim_name = (char *)malloc(sizeof(char)*(dim_name_length+1)); strncpy(dim_name, dimension_list+(counter-dim_name_length-sep_size), dim_name_length); dim_name[dim_name_length] = '\0'; DEBUG_MESS("svar %s: dim_name_length %d, dim_name: %s, counter: %d\n", svar->name, dim_name_length, dim_name, counter); DEBUG_MESS("i %d\n", i); DEBUG_MESS("Getting id for dim %s\n", dim_name); /*Find the id for the dimension named dim_name*/ dimension_ids[i] = -1; for (j=0;j<sfile->n_dimensions;j++){ DEBUG_MESS("j %d\n", j); if (!strcmp(dim_name, sfile->dimensions[j]->name)) dimension_ids[i] = sfile->dimensions[j]->nc_id; } if (dimension_ids[i]==-1){ printf("Dimension %s is undefined!\n", dim_name); abort(); } DEBUG_MESS("Finished loop\n"); DEBUG_MESS("dim %s has id %d \n", dim_name, dimension_ids[i]); free(dim_name); } svar->dimension_ids = dimension_ids; }
static PyObject * PyGSL_diff_generic(PyObject *self, PyObject *args, #ifndef PyGSL_DIFF_MODULE pygsl_deriv_func func #else pygsl_diff_func func #endif ) { PyObject *result=NULL, *myargs=NULL; PyObject *cb=NULL; pygsl_diff_args pargs = {NULL, NULL}; /* Changed to compile using Sun's Compiler */ gsl_function diff_gsl_callback = {NULL, NULL}; double x, value, abserr; int flag; #ifndef PyGSL_DIFF_MODULE double h; if(! PyArg_ParseTuple(args, "Odd|O", &cb, &x, &h, &myargs)){ return NULL; } #else if(! PyArg_ParseTuple(args, "Od|O", &cb, &x, &myargs)){ return NULL; } #endif /* Changed to compile using Sun's Compiler */ diff_gsl_callback.function = diff_callback; diff_gsl_callback.params = (void *) &pargs; if(! PyCallable_Check(cb)) { PyErr_SetString(PyExc_TypeError, "The first parameter must be callable"); return NULL; } Py_INCREF(cb); /* Add a reference to new callback */ pargs.callback = cb; /* Remember new callback */ /* Did I get arguments? If so handle them */ if(NULL == myargs){ Py_INCREF(Py_None); pargs.args= Py_None; }else{ Py_INCREF(myargs); pargs.args= myargs; } if((flag=setjmp(pargs.buffer)) == 0){ /* Jmp buffer set, call the function */ #ifndef PyGSL_DIFF_MODULE flag = func(&diff_gsl_callback, x, h, &value, &abserr); #else flag = func(&diff_gsl_callback, x, &value, &abserr); #endif }else{ DEBUG_MESS(2, "CALLBACK called longjmp! flag =%d", flag); } /* Arguments no longer used */ Py_DECREF(pargs.args); /* Dispose of callback */ Py_DECREF(pargs.callback); if(flag != GSL_SUCCESS){ PyGSL_ERROR_FLAG(flag); return NULL; } result = Py_BuildValue("(dd)", value, abserr); return result; }
static int PyGSL_PyArray_Check(PyArrayObject *a_array, int array_type, int flag, int nd, PyGSL_array_index_t *dimensions, int argnum, PyGSL_error_info * info) { int i; int error_flag = GSL_ESANITY, line = -1; PyGSL_array_index_t dim_temp; FUNC_MESS_BEGIN(); if(!PyArray_Check((PyObject *) a_array)){ pygsl_error("Did not recieve an array!", filename, __LINE__, GSL_ESANITY); line = __LINE__ - 2; error_flag = GSL_ESANITY; goto fail; } if(nd < 1 || nd > 2){ DEBUG_MESS(2, "Got an nd of %d", nd); line = __LINE__ - 2; pygsl_error("nd must either 1 or 2!", filename, __LINE__, GSL_ESANITY); error_flag = GSL_ESANITY; goto fail; } if (a_array->nd != nd){ DEBUG_MESS(3, "array->nd = %d\t nd = %d", a_array->nd, nd); line = __LINE__ - 1; sprintf(pygsl_error_str, "I could not convert argument number % 3d." " I expected a %s, but got an array of % 3d dimensions!\n", argnum, (nd == 1) ? "vector" : "matrix", a_array->nd); if (info){ info->error_description = pygsl_error_str; PyGSL_set_error_string_for_callback(info); } else { pygsl_error(pygsl_error_str, filename, __LINE__, GSL_EBADLEN); } error_flag = GSL_EBADLEN; goto fail; } for(i=0; i<nd; ++i){ if(dimensions[i] == -1 ){ switch(i){ case 0: DEBUG_MESS(2, "\t\t No one cares about its first dimension! %d", 0); break; case 1: DEBUG_MESS(2, "\t\t No one cares about its second dimension! %d", 0); break; default: error_flag = GSL_ESANITY; line = __LINE__ - 3; goto fail; break; } continue; } /* Check to be performed ... */ dim_temp = a_array->dimensions[i]; DEBUG_MESS(9, "Dimension %d has %ld elements", i, (unsigned long) dim_temp); if(dim_temp != (dimensions[i])){ sprintf(pygsl_error_str, "The size of argument % 3d did not match the expected" " size for the %d dimension. I got % 3ld elements but" " expected % 3ld elements!\n", argnum, i, (long)a_array->dimensions[0],(long)dimensions[0]); if (info){ info->error_description = pygsl_error_str; PyGSL_set_error_string_for_callback(info); } else { pygsl_error(pygsl_error_str, filename, __LINE__, GSL_EBADLEN); } error_flag = GSL_EBADLEN; line = __LINE__ - 11; goto fail; } } if(a_array->data == NULL){ pygsl_error("Got an array object were the data was NULL!", filename, __LINE__, GSL_ESANITY); error_flag = GSL_ESANITY; line = __LINE__ - 4; goto fail; } if(a_array->descr->type_num == array_type){ DEBUG_MESS(4, "\t\tArray type matched! %d", 0); }else{ pygsl_error("The array type did not match the spezified one!", filename, __LINE__, GSL_ESANITY); DEBUG_MESS(4, "Found an array type of %d but expected %d", (int) a_array->descr->type_num, array_type); error_flag = GSL_ESANITY; line = __LINE__ - 6; goto fail; } if ((flag & PyGSL_CONTIGUOUS) == 0){ DEBUG_MESS(2, "\t\t Can deal with discontiguous arrays! flag = %d", flag); } else { if(!(a_array->flags & CONTIGUOUS)){ DEBUG_MESS(3, "array->flags %d requested flags %d", a_array->flags, flag); pygsl_error("The array is not contiguous as requested!", filename, __LINE__, GSL_ESANITY); error_flag = GSL_ESANITY; line = __LINE__ - 3; goto fail; } } FUNC_MESS_END(); return GSL_SUCCESS; fail: PyGSL_add_traceback(NULL, filename, __FUNCTION__, line); DEBUG_MESS(4, "common array types: Double %d, CDouble %d", PyArray_DOUBLE, PyArray_CDOUBLE); DEBUG_MESS(4, "integer: Long %d, Int %d, Short %d", PyArray_LONG, PyArray_INT, PyArray_SHORT); /* DEBUG_MESS(8, "Char type %d Byte type %d String type %d", PyArray_CHAR, PyArray_BYTE, PyArray_STRING); */ return error_flag; }
static PyArrayObject * PyGSL_vector_check(PyObject *src, PyGSL_array_index_t size, PyGSL_array_info_t ainfo, PyGSL_array_index_t *stride, PyGSL_error_info * info) { int line=-1, tries; PyArrayObject * a_array = NULL; int array_type, flag, argnum, type_size; FUNC_MESS_BEGIN(); array_type = PyGSL_GET_ARRAYTYPE(ainfo); flag = PyGSL_GET_ARRAYFLAG(ainfo); type_size = PyGSL_GET_TYPESIZE(ainfo); argnum = PyGSL_GET_ARGNUM(ainfo); DEBUG_MESS(2, "Type requests: array_type %d, flag %d, c type_size %d, argnum %d", array_type, flag, type_size, argnum); /* * numpy arrays are which are non contiguous can have a stride which does * not match the basis type. In this case the conversion is repeated and * a contiguous array is demanded. */ for(tries = 0; tries <2; ++tries){ /* try fast conversion first */ #if 0 a_array = PyGSL_VECTOR_CONVERT(src, array_type, flag); if(a_array != NULL && a_array->nd == 1 && (size == -1 || a_array->dimensions[0] == size)) { /* good ... everything fine */ ; #else if(0) { ; #endif } else { /* lets try if that goes okay */ if(a_array != NULL){ /* Array could be allocated previously... */ Py_DECREF(a_array); } a_array = PyGSL_PyArray_prepare_gsl_vector_view(src, array_type, flag, size, argnum, info); if(a_array == NULL){ line = __LINE__ - 2; goto fail; } } if(stride == NULL){ /* no one interested in the stride. so help yourself ... */ break; /* will return array */ } if(PyGSL_STRIDE_RECALC(a_array->strides[0], type_size, stride) == GSL_SUCCESS){ /* everybody happy I hope! */ if((flag & PyGSL_CONTIGUOUS) == 1){ /* * just a check to see ... could be disabled later on when * the code is tested a little */ if(PyGSL_DEBUG_LEVEL() > 0){ if(*stride != 1){ line = __LINE__ - 1; pygsl_error("Stride not one for a contiguous array!", filename, line, GSL_ESANITY); goto fail; } /* stride */ }/*debug level */ }/* contiguous arrays */ break; /* will return array */ } /* Lets try to see if it makes sense to meake a copy */ DEBUG_MESS(2, "Stride recalc failed type size is %ld, array stride[0] is %ld", (long)type_size, (long)a_array->strides[0]); if((flag & PyGSL_CONTIGUOUS) == 1){ line = __LINE__ - 1; pygsl_error("Why does the stride recalc fail for a contigous array?", filename, line, GSL_ESANITY); goto fail; } else { /* keep the flags, but demand contiguous this time */ flag -= (flag & PyGSL_CONTIGUOUS); assert(a_array); Py_DECREF(a_array); a_array = NULL; } }/* number of tries */ DEBUG_MESS(7, "Checking refcount src obj @ %p had %d cts and array @ %p has now %d cts", (void *) src, src->ob_refcnt, (void *)a_array, a_array->ob_refcnt); /* handling failed stride recalc */ FUNC_MESS_END(); return a_array; fail: FUNC_MESS("Fail"); PyGSL_add_traceback(NULL, filename, __FUNCTION__, line); Py_XDECREF(a_array); return NULL; } /* * maximum 0xffffffff * * array_type = flag & 0x000000ff * type_size = flag & 0x0000ff00 * array_flag = flag & 0x00ff0000 * argnum = flag & 0xff000000 */ static PyArrayObject * PyGSL_matrix_check(PyObject *src, PyGSL_array_index_t size1, PyGSL_array_index_t size2, PyGSL_array_info_t ainfo, PyGSL_array_index_t *stride1, PyGSL_array_index_t *stride2, PyGSL_error_info * info) { int line= __LINE__, tries, j; PyArrayObject * a_array = NULL; PyGSL_array_index_t * stride; int array_type, flag, argnum, type_size; FUNC_MESS_BEGIN(); array_type = PyGSL_GET_ARRAYTYPE(ainfo); flag = PyGSL_GET_ARRAYFLAG(ainfo); type_size = PyGSL_GET_TYPESIZE(ainfo); argnum = PyGSL_GET_ARGNUM(ainfo); /* * numpy arrays are which are non contiguous can have a stride which does * not match the basis type. In this case the conversion is repeated and * a contiguous array is demanded. */ for(tries = 0; tries <2; ++tries){ #if 0 a_array = PyGSL_MATRIX_CONVERT(src, array_type, flag); if(a_array != NULL && a_array->nd == 1 && (size1 == -1 || a_array->dimensions[0] == size1) && (size2 == -1 || a_array->dimensions[1] == size2)) #else if(0) #endif { /* good ... everything fine */ ; } else { DEBUG_MESS(4, "PyGSL_MATRIX_CONVERT failed a_array = %p", a_array); /* lets try if that goes okay */ if(a_array != NULL){ /* Array could be allocated previously... */ Py_DECREF(a_array); } a_array = PyGSL_PyArray_prepare_gsl_matrix_view(src, array_type, flag, size1, size2, argnum, info); if(a_array == NULL){ line = __LINE__ - 2; goto fail; } } for(j = 0; j<2; ++j){ switch(j){ case 0: stride = stride1; break; case 1: stride = stride2; break; default: assert(0); } if(stride == NULL){ /* no one interested in the stride. lets check the other one ... */ continue; } if(PyGSL_STRIDE_RECALC(a_array->strides[j], type_size, stride) == GSL_SUCCESS){ /* everybody happy I hope! */ if((flag & PyGSL_CONTIGUOUS) == 1){ /* * just a check to see ... could be disabled later on when * the code is tested a little * * 14. Oktober 2009: but only for the last dimension! */ if(j == 1 && *stride != 1){ line = __LINE__ - 1; DEBUG_MESS(6, "array stride %ld, type size %d, " "found a stride of %ld", (long) a_array->strides[j], type_size, (long) *stride); pygsl_error("Stride not one of a contiguous array!", filename, line, GSL_ESANITY); goto fail; } } } else { /* Lets try to see if it makes sense to meake a copy */ DEBUG_MESS(2, "Stride recalc failed type size is %ld, array stride[0] is %ld", (long)type_size, (long)a_array->strides[j]); if((flag & PyGSL_CONTIGUOUS) == 1){ line = __LINE__ - 1; pygsl_error("Why does the stride recalc fail for a contigous array?", filename, line, GSL_ESANITY); goto fail; } else { /* keep the flags, but demand contiguous this time */ DEBUG_MESS(3, "Matrix %p ot satisfying requests, trying this time contiguous", (void *) a_array); flag -= (flag & PyGSL_CONTIGUOUS); Py_DECREF(a_array); a_array = NULL; break; } } /* recalc stride or try again */ } /* check strides */ }/* number of tries */ /* handling failed stride recalc */ FUNC_MESS_END(); return a_array; fail: PyGSL_add_traceback(NULL, filename, __FUNCTION__, line); Py_XDECREF(a_array); return NULL; }
/* Has to be at the bottom because it uses many * other functions */ void sdatio_open_file(struct sdatio_file * sfile ) { /*printf("called\n");*/ int retval; int ndims, nvars; int i, j; struct sdatio_dimension * sdim; struct sdatio_variable * svar; size_t lengthp; int nunlimdims; int *unlimdims; int is_unlimited; char * name_tmp; int * vardimids; char * dimension_list; nc_type vartype; int vartypeint; int dummy; int *nunlim; DEBUG_MESS("starting sdatio_open_file\n"); retval = 0; if (sfile->is_open){ printf("ERROR: The supplied sdatio_file struct corresponds to an already open file.\n"); abort(); } /* We make the choice to always open the file in readwrite mode*/ sfile->mode = sfile->mode|NC_WRITE; DEBUG_MESS("sdatio_open_file opening file\n"); /* Open the file*/ if (sfile->is_parallel) { #ifdef PARALLEL if ((retval = nc_open_par(sfile->name, sfile->mode, *(sfile->communicator), MPI_INFO_NULL, &(sfile->nc_file_id)))) ERR(retval); #else printf("sdatio was built without --enable-parallel, sdatio_create_file will not work for parallel files\n"); abort(); #endif } else { if ((retval = nc_open(sfile->name, sfile->mode, &(sfile->nc_file_id)))) ERR(retval); } DEBUG_MESS("sdatio_open_file opened file\n"); /*Initialize object file data*/ sfile->is_open = 1; sfile->n_dimensions = 0; sfile->n_variables = 0; sfile->data_written = 0; /* Get number of dimensions in the file*/ if ((retval = nc_inq_ndims(sfile->nc_file_id, &ndims))) ERR(retval); /* Allocate some temp storate*/ name_tmp = (char*)malloc(sizeof(char*)*(NC_MAX_NAME+1)); /* Get a list of unlimited dimensions*/ if ((retval = nc_inq_unlimdims(sfile->nc_file_id, &nunlimdims, NULL))) ERR(retval); unlimdims = (int*)malloc(sizeof(int*)*nunlimdims); if ((retval = nc_inq_unlimdims(sfile->nc_file_id, &nunlimdims, unlimdims))) ERR(retval); /* Add each dimension to the sfile object*/ for (i=0; i<ndims; i++){ if ((retval = nc_inq_dim(sfile->nc_file_id, i, name_tmp, &lengthp))) ERR(retval); sdim = (struct sdatio_dimension *) malloc(sizeof(struct sdatio_dimension)); /*if ((retval = nc_def_dim(sfile->nc_file_id, dimension_name, size, &(sdim->nc_id)))) ERR(retval);*/ /*}*/ /*sdatio_end_definitions(sfile);*/ is_unlimited = 0; for(j=0; j<nunlimdims; j++) if (unlimdims[j] == i) is_unlimited = 1; if (is_unlimited) { sdim->size = SDATIO_UNLIMITED; /* We choose the first write to unlimited variables to be the last * existing record. Users can easily move to the next by * calling sdatio_increment_start before writing anything */ sdim->start = lengthp-1; } else { sdim->size = lengthp; sdim->start = 0; } if (strlen(name_tmp)>1){ sfile->has_long_dim_names = 1; /*printf("Dimension names can only be one character long!\n");*/ /*abort();*/ } sdim->nc_id = i; sdim->name = (char *)malloc(sizeof(char)*(strlen(name_tmp)+1)); strcpy(sdim->name, name_tmp); sdatio_append_dimension(sfile, sdim); } DEBUG_MESS("Finished reading dimensions\n"); /* Get the number of variables in the file*/ if ((retval = nc_inq_nvars(sfile->nc_file_id, &nvars))) ERR(retval); /* Add each variable to the sfile object*/ for (i=0; i<nvars; i++){ if ((retval = nc_inq_varndims(sfile->nc_file_id, i, &ndims))) ERR(retval); vardimids = (int*)malloc(sizeof(int)*ndims); if ((retval = nc_inq_var(sfile->nc_file_id, i, name_tmp, &vartype, &ndims, vardimids, &dummy))) ERR(retval); vartypeint = vartype; vartypeint = sdatio_sdatio_variable_type(vartypeint); svar = (struct sdatio_variable *) malloc(sizeof(struct sdatio_variable)); /*Set variable id*/ svar->nc_id = i; /* Set variable name*/ svar->name = (char *)malloc(sizeof(char)*(strlen(name_tmp)+1)); strcpy(svar->name, name_tmp); /*ndims = strlen(dimension_list);*/ svar->ndims = ndims; DEBUG_MESS("ndims = %d for variable %s\n", ndims, name_tmp); /* Set the dimension_ids*/ svar->dimension_ids = vardimids; /*sdatio_get_dimension_ids(sfile, dimension_list, svar);*/ sdatio_get_dimension_list(sfile, vardimids, svar); /*svar->dimension_ids = dimension_ids;*/ DEBUG_MESS("Setting type for variable %s\n", svar->name); switch (vartypeint){ case SDATIO_INT: svar->type_size = sizeof(int); break; case SDATIO_FLOAT: svar->type_size = sizeof(float); break; case SDATIO_DOUBLE: svar->type_size = sizeof(double); break; case SDATIO_CHAR: svar->type_size = sizeof(char); break; default: printf("Unknown type in sdatio_create_variable\n"); abort(); } svar->type = vartypeint; DEBUG_MESS("Allocating manual starts and counts for variable %s; ndims %d\n", svar->name, ndims); svar->manual_starts=(int*)malloc(sizeof(int)*ndims); svar->manual_counts=(int*)malloc(sizeof(int)*ndims); svar->manual_offsets=(int*)malloc(sizeof(int)*ndims); DEBUG_MESS("Setting manual starts and counts for variable %s\n", svar->name); for (j=0;j<ndims;j++){ svar->manual_starts[j]=-1; svar->manual_counts[j]=-1; svar->manual_offsets[j]=-1; } DEBUG_MESS("Starting sdatio_append_variable\n"); sdatio_append_variable(sfile, svar); DEBUG_MESS("Ending sdatio_append_variable\n"); #ifdef PARALLEL if (sfile->is_parallel){ sdatio_number_of_unlimited_dimensions(sfile, svar->name, &nunlim); if (nunlim > 0) if ((retval = nc_var_par_access(sfile->nc_file_id, svar->nc_id, NC_COLLECTIVE))) ERR(retval); } #endif /*vartypeint = */ } DEBUG_MESS("Finished reading variables\n"); free(unlimdims); free(name_tmp); }
static void /* generic instance destruction */ generic_dealloc (PyObject *self) { DEBUG_MESS(1, " *** generic_dealloc %p\n", (void *) self); PyMem_Free(self); }
void sdatio_get_dimension_list(struct sdatio_file * sfile, int * dimension_ids, struct sdatio_variable * svar){ int ndims; int i,j, loop; /*char dim_name[2];*/ /*char * dim_name;*/ char * dimension_list; struct sdatio_dimension * sdim; int dim_name_length; int counter; int sep_size; ndims = svar->ndims; DEBUG_MESS("ndims %d\n", ndims); /* In the first iteration, we calculate * the length of the dimension list string. * In the second iteration, we copy the names*/ for (loop=0; loop<2; loop++){ if (loop==1) dimension_list = (char *) malloc(sizeof(char)*(counter+1)); counter = 0; for (i=0;i<ndims;i++){ /* First we find the dimension object*/ for (j=0;j<sfile->n_dimensions+1;j++){ if (j==sfile->n_dimensions) { printf("ERROR: dimension with nc_id %d not found for variable %s\n", dimension_ids[i], svar->name); abort(); } if (sfile->dimensions[j]->nc_id==dimension_ids[i]) break; } sdim = sfile->dimensions[j]; /* In the next section we set counter to * the beginning of the dimension name and * dim_name_length to the size of the name*/ if (sfile->has_long_dim_names){ sep_size = 1; dim_name_length = strlen(sdim->name); } else { sep_size = 0; dim_name_length = 1; if (strlen(sdim->name)!=1){ printf("Fatal logic error: dim_name_length!=1 for but not has_long_dim_names\n"); abort(); } } if (dim_name_length < 1){ printf("ERROR: zero length dimension name in sdatio_open_file \n"); abort(); } DEBUG_MESS("In sdatio_get_dimension_list, dim id %d, has name %s, with counter %d\n", i, sdim->name, counter); /*dim_name[0] = dimension_list[i];*/ /*dim_name[1] = dimension_list[ndims];*/ /* Copy the name of the current dimension to the dimension_list*/ if (loop==1) strncpy(dimension_list+counter, sdim->name, dim_name_length); counter = counter + dim_name_length; if (i<(ndims-1)){ if (loop==1 && sep_size==1) strncpy(dimension_list+counter, ",", sep_size); counter = counter + sep_size; } } } dimension_list[counter] = '\0'; /*dim_name[dim_name_length] = '\0';*/ svar->dimension_list = dimension_list; }
void sdatio_create_variable(struct sdatio_file * sfile, int variable_type, char * variable_name, char * dimension_list, char * description, char * units){ int ndims; int nunlim; struct sdatio_variable * svar; int retval; /*int * dimension_ids;*/ /*dimension_ids = (int **)malloc(sizeof(int*));*/ /*printf("dimension_list is %s\n", dimension_list);*/ svar = (struct sdatio_variable *) malloc(sizeof(struct sdatio_variable)); /* Set variable name*/ svar->name = (char *)malloc(sizeof(char)*(strlen(variable_name)+1)); strcpy(svar->name, variable_name); /*ndims = strlen(dimension_list);*/ ndims = sdatio_ndims_from_string(sfile, dimension_list); svar->ndims = ndims; DEBUG_MESS("ndims = %d for variable %s\n", ndims, variable_name); sdatio_get_dimension_ids(sfile, dimension_list, svar); /*svar->dimension_ids = dimension_ids;*/ sdatio_recommence_definitions(sfile); /*if (sfile->is_parallel){}*/ /*else {*/ if ((retval = nc_def_var(sfile->nc_file_id, variable_name, sdatio_netcdf_variable_type(variable_type), ndims, svar->dimension_ids, &(svar->nc_id)))) ERR(retval); if ((retval = nc_put_att_text(sfile->nc_file_id, svar->nc_id, "description", strlen(description), description))) ERR(retval); if ((retval = nc_put_att_text(sfile->nc_file_id, svar->nc_id, "units", strlen(units), units))) ERR(retval); /*}*/ switch (variable_type){ case SDATIO_INT: svar->type_size = sizeof(int); break; case SDATIO_FLOAT: svar->type_size = sizeof(float); break; case SDATIO_DOUBLE: svar->type_size = sizeof(double); break; case SDATIO_CHAR: svar->type_size = sizeof(char); break; default: printf("Unknown type in sdatio_create_variable\n"); abort(); } sdatio_end_definitions(sfile); svar->type = variable_type; svar->dimension_list = (char *)malloc(sizeof(char)*(strlen(dimension_list)+1)); strcpy(svar->dimension_list, dimension_list); svar->manual_starts=(int*)malloc(sizeof(int)*ndims); svar->manual_counts=(int*)malloc(sizeof(int)*ndims); svar->manual_offsets=(int*)malloc(sizeof(int)*ndims); int i; for (i=0;i<ndims;i++){ svar->manual_starts[i]=-1; svar->manual_counts[i]=-1; svar->manual_offsets[i]=-1; } DEBUG_MESS("Starting sdatio_append_variable\n"); sdatio_append_variable(sfile, svar); DEBUG_MESS("Ending sdatio_append_variable\n"); #ifdef PARALLEL if (sfile->is_parallel){ sdatio_number_of_unlimited_dimensions(sfile, variable_name, &nunlim); if (nunlim > 0) if ((retval = nc_var_par_access(sfile->nc_file_id, svar->nc_id, NC_COLLECTIVE))) ERR(retval); } #endif }
static PyObject * PyGSL_odeiv_evolve_apply_array(PyGSL_odeiv_evolve *self, PyObject *args) { PyObject *result = NULL; PyObject *y0_o = NULL; PyObject *t_o = NULL; PyArrayObject *volatile y0 = NULL, *volatile yout = NULL, *volatile ta = NULL; double t=0, h=0, t1 = 0, flag; double *y0_data = NULL, *yout_data = NULL; int dimension, r, dims[2]; int nlen=-1, i, j; assert(PyGSL_ODEIV_EVOLVE_Check(self)); FUNC_MESS_BEGIN(); if(! PyArg_ParseTuple(args, "OdO", &t_o, &h, &y0_o)){ return NULL; } dimension = self->step->system.dimension; ta = PyGSL_PyArray_PREPARE_gsl_vector_view(t_o, PyArray_DOUBLE, 1, -1, 1, NULL); if(ta == NULL) goto fail; nlen = ta->dimensions[0]; y0 = PyGSL_PyArray_PREPARE_gsl_vector_view(y0_o, PyArray_DOUBLE, 1, dimension, 1, NULL); if(y0 == NULL) goto fail; dims[0] = nlen; dims[1] = dimension; yout = (PyArrayObject *) PyArray_FromDims(2, dims, PyArray_DOUBLE); if(yout == NULL) goto fail; if((flag=setjmp(self->step->buffer)) == 0){ FUNC_MESS("\t\t Setting Jmp Buffer"); } else { FUNC_MESS("\t\t Returning from Jmp Buffer"); goto fail; } DEBUG_MESS(5, "\t\t Yout Layout %p, shape=[%d,%d], strides=[%d,%d]", yout->data, yout->dimensions[0], yout->dimensions[1], yout->strides[0], yout->strides[1]); r = GSL_CONTINUE; yout_data = (double *)(yout->data); /* Copy the start vector */ for(j=0; j<dimension; j++){ yout_data[j] = *((double *)(y0->data + y0->strides[0] * j)); } for(i=1; i<nlen; i++){ y0_data = (double *)(yout->data + yout->strides[0]*(i-1)); yout_data = (double *)(yout->data + yout->strides[0]*(i) ); /* Copy the start vector */ for(j=0; j<dimension; j++){ yout_data[j] = y0_data[j]; } t = *((double *) (ta->data + ta->strides[0]*(i-1)) ); t1 = *((double *) (ta->data + ta->strides[0]*(i) ) ); while(t<t1){ r = gsl_odeiv_evolve_apply(self->evolve, self->control->control, self->step->step, &(self->step->system), &t, t1, &h, yout_data); /* All GSL Errors are > 0. */ if (r != GSL_SUCCESS){ goto fail; } } assert(r == GSL_SUCCESS); } result = (PyObject *) yout; /* Deleting the arrays */ Py_DECREF(y0); y0=NULL; Py_DECREF(ta); FUNC_MESS_END(); return result; fail: FUNC_MESS("IN Fail"); PyGSL_add_traceback(module, this_file, odeiv_step_init_err_msg, __LINE__); Py_XDECREF(ta); ta = NULL; Py_XDECREF(y0); y0=NULL; Py_XDECREF(yout); FUNC_MESS("IN Fail End"); return NULL; }
static PyObject* PyGSL_multimin_set_f(PyGSL_solver *self, PyObject *pyargs, PyObject *kw) { PyObject* func = NULL, * args = Py_None, * x = NULL, * steps = NULL; PyArrayObject * xa = NULL, * stepsa = NULL; int n=0, flag=GSL_EFAILED; PyGSL_array_index_t stride_recalc; gsl_vector_view gsl_x; gsl_vector_view gsl_steps; gsl_multimin_function * c_sys = NULL; static const char *kwlist[] = {"f", "x0", "args", "steps", NULL}; FUNC_MESS_BEGIN(); assert(PyGSL_solver_check(self)); if (self->solver == NULL) { pygsl_error("Got a NULL Pointer of min.f", filename, __LINE__ - 3, GSL_EFAULT); return NULL; } assert(pyargs); /* arguments PyFunction, Parameters, start Vector, step Vector */ if (0==PyArg_ParseTupleAndKeywords(pyargs,kw,"OOOO", (char **)kwlist, &func, &x,&args,&steps)) return NULL; if(!PyCallable_Check(func)){ pygsl_error("First argument must be callable", filename, __LINE__ - 3, GSL_EBADFUNC); return NULL; } n=self->problem_dimensions[0]; xa = PyGSL_vector_check(x, n, PyGSL_DARRAY_INPUT(4), &stride_recalc, NULL); if (xa == NULL){ PyGSL_add_traceback(module, filename, __FUNCTION__, __LINE__ - 1); goto fail; } gsl_x = gsl_vector_view_array_with_stride((double *)(xa->data), stride_recalc, xa->dimensions[0]); stepsa = PyGSL_vector_check(steps, n, PyGSL_DARRAY_INPUT(5), &stride_recalc, NULL); if (stepsa == NULL){ PyGSL_add_traceback(module, filename, __FUNCTION__, __LINE__ - 1); goto fail; } gsl_steps = gsl_vector_view_array_with_stride((double *)(stepsa->data), stride_recalc, stepsa->dimensions[0]); if (self->c_sys != NULL) { /* free the previous function and args */ c_sys = self->c_sys; } else { /* allocate function space */ c_sys=calloc(1, sizeof(gsl_multimin_function)); if (c_sys == NULL) { pygsl_error("Could not allocate the object for the minimizer function", filename, __LINE__ - 3, GSL_ENOMEM); goto fail; } } DEBUG_MESS(3, "Everything allocated args = %p", args); /* add new function and parameters */ if(PyGSL_solver_func_set(self, args, func, NULL, NULL) != GSL_SUCCESS) goto fail; /* initialize the function struct */ c_sys->n=n; c_sys->f=PyGSL_multimin_function_f; c_sys->params=(void*)self; DEBUG_MESS(3, "Setting jmp buffer isset = % d", self->isset); if((flag = setjmp(self->buffer)) == 0){ self->isset = 1; flag = gsl_multimin_fminimizer_set(self->solver,c_sys, &gsl_x.vector, &gsl_steps.vector); if(PyGSL_ERROR_FLAG(flag) != GSL_SUCCESS){ goto fail; } } else { goto fail; } DEBUG_MESS(4, "Set evaluated. flag = %d", flag); self->c_sys = c_sys; Py_DECREF(xa); Py_DECREF(stepsa); self->set_called = 1; Py_INCREF(Py_None); self->isset = 0; FUNC_MESS_END(); return Py_None; fail: FUNC_MESS("Fail"); PyGSL_ERROR_FLAG(flag); self->isset = 0; Py_XDECREF(xa); Py_XDECREF(stepsa); return NULL; }