1 #include <pygsl/solver.h>
2 #include <gsl/gsl_multifit_nlin.h>
3 
4 const char  * filename = __FILE__;
5 PyObject *module = NULL;
6 static const char multifit_f_type_name[] = "F-MultiFitSolver";
7 static const char multifit_fdf_type_name[] = "FdF-MultiFitSolver";
8 
9 
10 static int
PyGSL_multifit_function_wrap(const gsl_vector * x,void * params,gsl_vector * f)11 PyGSL_multifit_function_wrap(const gsl_vector *x, void *params, gsl_vector *f)
12 {
13      PyGSL_solver *self = (PyGSL_solver *) params;
14      return PyGSL_function_wrap_Op_On(x, f, self->cbs[0], self->args, x->size, f->size,
15 				      __FUNCTION__);
16 }
17 
18 
19 static int
PyGSL_multifit_function_wrap_df(const gsl_vector * x,void * params,gsl_matrix * df)20 PyGSL_multifit_function_wrap_df(const gsl_vector *x, void *params, gsl_matrix *df)
21 {
22      PyGSL_solver *self = (PyGSL_solver *) params;
23      /* size 1 or size 2 from matrix ? */
24      return PyGSL_function_wrap_Op_Opn(x, df, self->cbs[1], self->args, df->size1, x->size,
25 				       __FUNCTION__);
26 }
27 
28 static int
PyGSL_multifit_function_wrap_fdf(const gsl_vector * x,void * params,gsl_vector * f,gsl_matrix * df)29 PyGSL_multifit_function_wrap_fdf(const gsl_vector *x, void *params, gsl_vector *f, gsl_matrix *df)
30 {
31      PyGSL_solver *self = (PyGSL_solver *) params;
32      /* size 1 or size 2 from matrix ? */
33      return PyGSL_function_wrap_Op_On_Opn(x, f, df, self->cbs[2], self->args, f->size, x->size,
34 					  __FUNCTION__);
35 }
36 
37 static PyObject *
PyGSL_multifit_fdfsolver_set(PyGSL_solver * self,PyObject * pyargs,PyObject * kw)38 PyGSL_multifit_fdfsolver_set(PyGSL_solver *self, PyObject *pyargs, PyObject *kw)
39 {
40 
41      gsl_multifit_function_fdf * c_sys;
42      struct pygsl_solver_n_set info = {1, NULL, (set_m_t)gsl_multifit_fdfsolver_set};
43      PyObject * tmp;
44 
45      FUNC_MESS_BEGIN();
46      if(self->c_sys == NULL){
47 	  if((c_sys = calloc(1, sizeof(gsl_multifit_function_fdf))) == NULL){
48 	       PyGSL_ERROR_NULL("Could not allocate the memory for the c_sys", GSL_ENOMEM);
49 	  }
50 	  c_sys->n=self->problem_dimensions[1];
51 	  c_sys->p=self->problem_dimensions[0];
52 	  c_sys->f  = PyGSL_multifit_function_wrap;
53 	  c_sys->df = PyGSL_multifit_function_wrap_df;
54 	  c_sys->fdf = PyGSL_multifit_function_wrap_fdf;
55 	  c_sys->params=(void*)self;
56 	  info.c_sys = c_sys;
57      }else{
58 	  info.c_sys = self->c_sys;
59      }
60      tmp =  PyGSL_solver_n_set(self, pyargs, kw, &info);
61      if(tmp == NULL){
62 	  PyGSL_add_traceback(module, __FILE__, __FUNCTION__, __LINE__ - 2);
63      }
64      FUNC_MESS_END();
65      return tmp;
66 }
67 
68 static PyObject*
PyGSL_multifit_fdfsolver_position(PyGSL_solver * self,PyObject * args)69 PyGSL_multifit_fdfsolver_position(PyGSL_solver *self, PyObject *args)
70 {
71      return PyGSL_solver_ret_vec(self, args, (ret_vec)gsl_multifit_fdfsolver_position);
72 }
73 static gsl_multifit_fdfsolver *
PyGSL_get_multifit_solver(PyGSL_solver * self)74 PyGSL_get_multifit_solver(PyGSL_solver *self)
75 {
76      FUNC_MESS_BEGIN();
77      assert(PyGSL_solver_check(self));
78      FUNC_MESS_END();
79      return (gsl_multifit_fdfsolver *) (self->solver);
80 }
81 
82 static PyObject*
PyGSL_multifit_fdfsolver_x(PyGSL_solver * self,PyObject * args)83 PyGSL_multifit_fdfsolver_x(PyGSL_solver *self, PyObject *args)
84 {
85      return (PyObject *) PyGSL_copy_gslvector_to_pyarray(PyGSL_get_multifit_solver(self)->x);
86 }
87 
88 static PyObject*
PyGSL_multifit_fdfsolver_dx(PyGSL_solver * self,PyObject * args)89 PyGSL_multifit_fdfsolver_dx(PyGSL_solver *self, PyObject *args)
90 {
91      return (PyObject *) PyGSL_copy_gslvector_to_pyarray(PyGSL_get_multifit_solver(self)->dx);
92 }
93 
94 static PyObject*
PyGSL_multifit_fdfsolver_f(PyGSL_solver * self,PyObject * args)95 PyGSL_multifit_fdfsolver_f(PyGSL_solver *self, PyObject *args)
96 {
97      return (PyObject *) PyGSL_copy_gslvector_to_pyarray(PyGSL_get_multifit_solver(self)->f);
98 }
99 
100 static PyObject*
PyGSL_multifit_fdfsolver_J(PyGSL_solver * self,PyObject * args)101 PyGSL_multifit_fdfsolver_J(PyGSL_solver *self, PyObject *args)
102 {
103      return (PyObject *) PyGSL_copy_gslmatrix_to_pyarray(PyGSL_get_multifit_solver(self)->J);
104 }
105 
106 static PyObject*
PyGSL_multifit_fdfsolver_test_delta(PyGSL_solver * self,PyObject * args)107 PyGSL_multifit_fdfsolver_test_delta(PyGSL_solver *self, PyObject *args)
108 {
109      int flag;
110      double epsabs, epsrel;
111      gsl_multifit_fdfsolver *s = self->solver;
112      if(!PyArg_ParseTuple(args, "dd", &epsabs, &epsrel))
113 	  return NULL;
114      flag = gsl_multifit_test_delta(s->dx, s->x, epsabs, epsrel);
115      return PyGSL_ERROR_FLAG_TO_PYINT(flag);
116 }
117 
118 static PyObject*
PyGSL_multifit_fdfsolver_test_gradient(PyGSL_solver * self,PyObject * args)119 PyGSL_multifit_fdfsolver_test_gradient(PyGSL_solver *self, PyObject *args)
120 {
121      int flag;
122      double epsabs;
123      gsl_vector *g = NULL;
124      gsl_multifit_fdfsolver *s = self->solver;
125 
126      if(!PyArg_ParseTuple(args, "d", &epsabs))
127 	  return NULL;
128      flag = gsl_multifit_gradient(s->J, s->f, g);
129      if(PyGSL_ERROR_FLAG(flag) != GSL_SUCCESS)
130 	  return NULL;
131      flag = gsl_multifit_test_gradient(g, epsabs);
132      return PyGSL_ERROR_FLAG_TO_PYINT(flag);
133 }
134 
135 static PyMethodDef PyGSL_multifit_fmethods[] = {
136      {NULL, NULL, 0, NULL}
137 };
138 
139 static PyMethodDef PyGSL_multifit_fdfmethods[] = {
140      {"J",    (PyCFunction)PyGSL_multifit_fdfsolver_J,   METH_NOARGS, NULL},
141      {"dx",   (PyCFunction)PyGSL_multifit_fdfsolver_dx,   METH_NOARGS, NULL},
142      {"x",    (PyCFunction)PyGSL_multifit_fdfsolver_x,   METH_NOARGS, NULL},
143      {"position",    (PyCFunction)PyGSL_multifit_fdfsolver_position,   METH_NOARGS, NULL},
144      {"f",    (PyCFunction)PyGSL_multifit_fdfsolver_f,    METH_NOARGS, NULL},
145      {"set",  (PyCFunction)PyGSL_multifit_fdfsolver_set,  METH_VARARGS|METH_KEYWORDS, NULL},
146      {"test_delta",  (PyCFunction)PyGSL_multifit_fdfsolver_test_delta,  METH_VARARGS, NULL},
147      {"test_gradient",  (PyCFunction)PyGSL_multifit_fdfsolver_test_gradient,  METH_VARARGS, NULL},
148      {NULL, NULL, 0, NULL}
149 };
150 
151 
152 
153 const struct _SolverStatic
154 multifit_solver_f   = {{ (void_m_t) gsl_multifit_fsolver_free,
155 			 /* gsl_multifit_fsolver_restart */  (void_m_t) NULL,
156 			 (name_m_t) gsl_multifit_fsolver_name,
157 			 (int_m_t) gsl_multifit_fsolver_iterate},
158 		       1, PyGSL_multifit_fmethods, multifit_f_type_name},
159 multifit_solver_fdf = {{(void_m_t) gsl_multifit_fdfsolver_free,
160 		     /* gsl_multifit_fdfsolver_restart (void_m_t) */ NULL,
161 		     (name_m_t) gsl_multifit_fdfsolver_name,
162 		     (int_m_t)  gsl_multifit_fdfsolver_iterate},
163 		       3, PyGSL_multifit_fdfmethods, multifit_fdf_type_name};
164 
165 static PyObject*
PyGSL_multifit_f_init(PyObject * self,PyObject * args,const gsl_multifit_fsolver_type * type)166 PyGSL_multifit_f_init(PyObject *self, PyObject *args,
167 		      const gsl_multifit_fsolver_type * type)
168 {
169 
170      PyObject *tmp=NULL;
171      solver_alloc_struct s = {type, (void_an_t) gsl_multifit_fsolver_alloc,
172 			      &multifit_solver_f};
173      FUNC_MESS_BEGIN();
174      tmp = PyGSL_solver_dn_init(self, args, &s, 2);
175      FUNC_MESS_END();
176      return tmp;
177 }
178 
179 static PyObject*
PyGSL_multifit_fdf_init(PyObject * self,PyObject * args,const gsl_multifit_fdfsolver_type * type)180 PyGSL_multifit_fdf_init(PyObject *self, PyObject *args,
181 		      const gsl_multifit_fdfsolver_type * type)
182 {
183 
184      PyObject *tmp=NULL;
185      solver_alloc_struct s = {type, (void_an_t) gsl_multifit_fdfsolver_alloc,
186 			      &multifit_solver_fdf};
187      FUNC_MESS_BEGIN();
188      tmp = PyGSL_solver_dn_init(self, args, &s, 2);
189      FUNC_MESS_END();
190      return tmp;
191 }
192 
193 #define AMFIT_FDF(name)                                                  \
194 static PyObject* PyGSL_multifit_init_ ## name (PyObject *self, PyObject *args)\
195 {                                                                             \
196      PyObject *tmp = NULL;                                                    \
197      FUNC_MESS_BEGIN();                                                       \
198      tmp = PyGSL_multifit_fdf_init(self, args,  gsl_multifit_fdfsolver_ ## name); \
199      if (tmp == NULL){                                                        \
200 	  PyGSL_add_traceback(module, __FILE__, __FUNCTION__, __LINE__); \
201      }                                                                        \
202      FUNC_MESS_END();                                                         \
203      return tmp;                                                              \
204 }
205 
206 AMFIT_FDF(lmsder)
AMFIT_FDF(lmder)207 AMFIT_FDF(lmder)
208 
209 PyObject *
210 PyGSL_multifit_gradient(PyObject *self, PyObject *args)
211 {
212   PyArrayObject *J_a = NULL, *f_a = NULL, *g_a = NULL;
213   PyObject *J_o = NULL, *f_o = NULL;
214   gsl_vector_view f;
215   gsl_vector_view g;
216   gsl_matrix_view J;
217 
218   PyGSL_array_index_t stride_recalc, dimension;
219   int flag;
220 
221   if(!PyArg_ParseTuple(args, "OO:gsl_multifit_gradient", &J_o, &f_o)){
222        return NULL;
223   }
224 
225   J_a = PyGSL_matrix_check(J_o, -1, -1, PyGSL_DARRAY_CINPUT(1), NULL, NULL, NULL);
226   if(J_a == NULL) goto fail;
227 
228   dimension = PyArray_DIM(J_a, 0);
229   /* Numpy calculates strides in bytes, gsl in basis type */
230   f_a = PyGSL_vector_check(f_o, dimension, PyGSL_DARRAY_INPUT(2), &stride_recalc, NULL);
231   if(f_a == NULL) goto fail;
232 
233   dimension = PyArray_DIM(J_a, 1);
234   g_a = (PyArrayObject *) PyGSL_New_Array(1, &dimension, NPY_DOUBLE);
235   if(g_a == NULL) goto fail;
236 
237   J = gsl_matrix_view_array((double *) PyArray_DATA(J_a), PyArray_DIM(J_a, 0), PyArray_DIM(J_a, 1));
238   f = gsl_vector_view_array_with_stride((double *) PyArray_DATA(f_a), stride_recalc, PyArray_DIM(f_a, 0));
239   g = gsl_vector_view_array((double *) PyArray_DATA(g_a), dimension);
240   flag = gsl_multifit_gradient(&J.matrix, &f.vector, &g.vector);
241 
242   Py_DECREF(J_a);
243   Py_DECREF(f_a);
244 
245   if((PyGSL_ERROR_FLAG(flag)) != GSL_SUCCESS)
246        goto fail;
247 
248   return (PyObject * )g_a;
249 
250   fail :
251     Py_XDECREF(J_a);
252     Py_XDECREF(f_a);
253     Py_XDECREF(g_a);
254     return NULL;
255 }
256 
257 PyObject *
PyGSL_multifit_covar(PyObject * self,PyObject * args)258 PyGSL_multifit_covar(PyObject *self, PyObject *args)
259 {
260   PyArrayObject *J_a = NULL, *C_a = NULL;
261   PyObject *J_o = NULL;
262   gsl_matrix_view J, C;
263   PyGSL_array_index_t dimensions[2];
264   int flag;
265   double epsrel;
266 
267 
268   if(!PyArg_ParseTuple(args, "Od:gsl_multifit_covar", &J_o, &epsrel)){
269     return NULL;
270   }
271 
272   J_a = PyGSL_matrix_check(J_o, -1, -1, PyGSL_DARRAY_CINPUT(1), NULL, NULL, NULL);
273   if(J_a == NULL) goto fail;
274 
275   dimensions[0] = dimensions[1] = PyArray_DIM(J_a, 1);
276   C_a = (PyArrayObject *) PyGSL_New_Array(2, dimensions, NPY_DOUBLE);
277   if(C_a == NULL) goto fail;
278 
279   J = gsl_matrix_view_array((double *) PyArray_DATA(J_a), PyArray_DIM(J_a, 0), PyArray_DIM(J_a, 1));
280   C = gsl_matrix_view_array((double *) PyArray_DATA(C_a), PyArray_DIM(C_a, 0), PyArray_DIM(C_a, 1));
281 
282   flag = gsl_multifit_covar(&J.matrix, epsrel, &C.matrix);
283 
284   Py_DECREF(J_a);
285   if((PyGSL_ERROR_FLAG(flag)) != GSL_SUCCESS)
286        goto fail;
287   return (PyObject * )C_a;
288 
289   fail :
290     Py_XDECREF(J_a);
291     Py_XDECREF(C_a);
292     return NULL;
293 }
294 
295 static PyObject *
PyGSL_multifit_test_delta(PyObject * self,PyObject * args)296 PyGSL_multifit_test_delta(PyObject * self, PyObject * args)
297 {
298      return PyGSL_solver_vvdd_i(self, args, gsl_multifit_test_delta);
299 }
300 
301 static PyObject *
PyGSL_multifit_test_gradient(PyObject * self,PyObject * args)302 PyGSL_multifit_test_gradient(PyObject * self, PyObject * args)
303 {
304      return PyGSL_solver_vd_i(self, args, gsl_multifit_test_gradient);
305 }
306 
307 static PyMethodDef mMethods[] = {
308      /* multifit solvers */
309      {"lmder",          PyGSL_multifit_init_lmder,  METH_VARARGS, NULL},
310      {"lmsder",          PyGSL_multifit_init_lmsder,  METH_VARARGS, NULL},
311      /* multifit funcs */
312      {"fit_test_delta",    PyGSL_multifit_test_delta,     METH_VARARGS, NULL},
313      {"fit_test_gradient", PyGSL_multifit_test_gradient,  METH_VARARGS, NULL},
314      {"gradient",          PyGSL_multifit_gradient,  METH_VARARGS, NULL},
315      {"covar",             PyGSL_multifit_covar,  METH_VARARGS, NULL},
316      {NULL, NULL, 0, NULL}
317 };
318 
319 static const char PyGSL_multifit_nlin_module_doc[] = "XXX Missing \n";
320 void
initmultifit_nlin(void)321 initmultifit_nlin(void)
322 {
323      PyObject* m, *dict, *item;
324      FUNC_MESS_BEGIN();
325 
326      m=Py_InitModule("multifit_nlin", mMethods);
327      module = m;
328      assert(m);
329      dict = PyModule_GetDict(m);
330      if(!dict)
331 	  goto fail;
332 
333      init_pygsl()
334      import_pygsl_solver();
335      assert(PyGSL_API);
336 
337 
338      if (!(item = PyString_FromString((char*)PyGSL_multifit_nlin_module_doc))){
339 	  PyErr_SetString(PyExc_ImportError,
340 			  "I could not generate module doc string!");
341 	  goto fail;
342      }
343      if (PyDict_SetItemString(dict, "__doc__", item) != 0){
344 	  PyErr_SetString(PyExc_ImportError,
345 			  "I could not init doc string!");
346 	  goto fail;
347      }
348 
349      FUNC_MESS_END();
350      return;
351 
352  fail:
353      FUNC_MESS("FAIL");
354      return;
355 }
356