1 | |
---|
2 | #include "Python.h" |
---|
3 | |
---|
4 | #include <stdio.h> /* gets */ |
---|
5 | #include <stdlib.h> /* atoi, malloc */ |
---|
6 | #include <string.h> /* strcpy */ |
---|
7 | #include "numpy/arrayobject.h" |
---|
8 | #include "math.h" |
---|
9 | |
---|
10 | #include "util_ext.h" /* in utilities */ |
---|
11 | #include "sparse_dok.h" /*in utilities */ |
---|
12 | #include "quad_tree.h" |
---|
13 | #include "omp.h" |
---|
14 | |
---|
15 | // Errors defined for netcdf reading |
---|
16 | #define ERRCODE 2 |
---|
17 | #define ERR(e) {printf("Error: %s\n", nc_strerror(e)); exit(ERRCODE);} |
---|
18 | |
---|
19 | #include "patchlevel.h" |
---|
20 | |
---|
21 | // PYVERSION273 used to check python version for use of PyCapsule |
---|
22 | #if PY_MAJOR_VERSION>=2 && PY_MINOR_VERSION>=7 && PY_MICRO_VERSION>=3 |
---|
23 | #define PYVERSION273 |
---|
24 | #endif |
---|
25 | |
---|
26 | //-------------------------- QUANTITY FITTING ------------------------------ |
---|
27 | |
---|
28 | |
---|
29 | // Builds the matrix D used to smooth the interpolation |
---|
30 | // of a variables from scattered data points to a mesh. See fit.py for more details.s |
---|
31 | int _build_smoothing_matrix(int n, |
---|
32 | long* triangles, |
---|
33 | double* areas, |
---|
34 | double* vertex_coordinates, |
---|
35 | int* strides, |
---|
36 | sparse_dok * smoothing_mat) |
---|
37 | { |
---|
38 | |
---|
39 | |
---|
40 | int k; |
---|
41 | int k3,k6; |
---|
42 | int err = 0; |
---|
43 | edge_key_t key; |
---|
44 | |
---|
45 | double det,area,x0,x1,x2,y0,y1,y2; |
---|
46 | double a0,b0,a1,b1,a2,b2,e01,e12,e20; |
---|
47 | int v0,v1,v2; |
---|
48 | double smoothing_val; |
---|
49 | |
---|
50 | |
---|
51 | for(k=0; k<n; k++) { |
---|
52 | // multiple k by 3 to give index for triangles which store in 3-blocks |
---|
53 | k3=k*3; |
---|
54 | k6=k*6; |
---|
55 | // store the area for the current triangle |
---|
56 | area = areas[k]; |
---|
57 | // store current triangles global vertex indicies |
---|
58 | v0 = triangles[k3]; |
---|
59 | v1 = triangles[k3+1]; |
---|
60 | v2 = triangles[k3+2]; |
---|
61 | // store the locations of the three verticies |
---|
62 | x0 = vertex_coordinates[k6]; |
---|
63 | y0 = vertex_coordinates[k6+1]; |
---|
64 | x1 = vertex_coordinates[k6+2]; |
---|
65 | y1 = vertex_coordinates[k6+3]; |
---|
66 | x2 = vertex_coordinates[k6+4]; |
---|
67 | y2 = vertex_coordinates[k6+5]; |
---|
68 | |
---|
69 | // calculate gradients (move to external function?) |
---|
70 | |
---|
71 | det = (y2-y0)*(x1-x0) - (y1-y0)*(x2-x0); |
---|
72 | a0 = (y2-y0)*(0-1) - (y1-y0)*(0-1); |
---|
73 | a0 /= det; |
---|
74 | |
---|
75 | b0 = (x1-x0)*(0-1) - (x2-x0)*(0-1); |
---|
76 | b0 /= det; |
---|
77 | |
---|
78 | a1 = (y2-y0)*(1-0) - (y1-y0)*(0-0); |
---|
79 | a1 /= det; |
---|
80 | |
---|
81 | b1 = (x1-x0)*(0-0) - (x2-x0)*(1-0); |
---|
82 | b1 /= det; |
---|
83 | |
---|
84 | a2 = (y2-y0)*(0-0) - (y1-y0)*(1-0); |
---|
85 | a2 /= det; |
---|
86 | |
---|
87 | b2 = (x1-x0)*(1-0) - (x2-x0)*(0-0); |
---|
88 | b2 /= det; |
---|
89 | |
---|
90 | // insert diagonal contributions |
---|
91 | |
---|
92 | // v0,v0 |
---|
93 | key.i = v0; |
---|
94 | key.j = v0; |
---|
95 | smoothing_val = (a0*a0 + b0*b0)*area; |
---|
96 | add_dok_entry(smoothing_mat,key,smoothing_val); |
---|
97 | |
---|
98 | // v1,v1 |
---|
99 | key.i = v1; |
---|
100 | key.j = v1; |
---|
101 | smoothing_val = (a1*a1 + b1*b1)*area; |
---|
102 | add_dok_entry(smoothing_mat,key,smoothing_val); |
---|
103 | |
---|
104 | // v2,v2 |
---|
105 | key.i = v2; |
---|
106 | key.j = v2; |
---|
107 | smoothing_val = (a2*a2 + b2*b2)*area; |
---|
108 | add_dok_entry(smoothing_mat,key,smoothing_val); |
---|
109 | |
---|
110 | |
---|
111 | // insert off diagonal contributions |
---|
112 | e01 = (a0*a1 + b0*b1)*area; |
---|
113 | // v0,v1 (v1,v0) |
---|
114 | key.i = v0; |
---|
115 | key.j = v1; |
---|
116 | add_dok_entry(smoothing_mat,key,e01); |
---|
117 | key.i = v1; |
---|
118 | key.j = v0; |
---|
119 | add_dok_entry(smoothing_mat,key,e01); |
---|
120 | |
---|
121 | e12 = (a1*a2 + b1*b2)*area; |
---|
122 | // v1,v2 (v2,v1) |
---|
123 | key.i = v1; |
---|
124 | key.j = v2; |
---|
125 | add_dok_entry(smoothing_mat,key,e12); |
---|
126 | key.i = v2; |
---|
127 | key.j = v1; |
---|
128 | add_dok_entry(smoothing_mat,key,e12); |
---|
129 | |
---|
130 | e20 = (a2*a0 + b2*b0)*area; |
---|
131 | // v2,v0 (v0,v2) |
---|
132 | key.i = v2; |
---|
133 | key.j = v0; |
---|
134 | add_dok_entry(smoothing_mat,key,e20); |
---|
135 | key.i = v0; |
---|
136 | key.j = v2; |
---|
137 | add_dok_entry(smoothing_mat,key,e20); |
---|
138 | } |
---|
139 | |
---|
140 | return err; |
---|
141 | |
---|
142 | } |
---|
143 | |
---|
144 | // Builds a quad tree out of a list of triangles for quick |
---|
145 | // searching. |
---|
146 | quad_tree * _build_quad_tree(int n, |
---|
147 | long* triangles, |
---|
148 | double* vertex_coordinates, |
---|
149 | double* extents) |
---|
150 | { |
---|
151 | |
---|
152 | int k,k6; |
---|
153 | double x0,y0,x1,y1,x2,y2; |
---|
154 | |
---|
155 | // set up quad tree and allocate memory |
---|
156 | quad_tree * tree = new_quad_tree(extents[0],extents[1],extents[2],extents[3]); |
---|
157 | |
---|
158 | // iterate through triangles |
---|
159 | for(k=0; k<n; k++) { |
---|
160 | // multiple k by 3 to give index for triangles which store in 3-blocks |
---|
161 | k6=k*6; |
---|
162 | // store the locations of the three verticies |
---|
163 | x0 = vertex_coordinates[k6]; |
---|
164 | y0 = vertex_coordinates[k6 + 1]; |
---|
165 | x1 = vertex_coordinates[k6 + 2]; |
---|
166 | y1 = vertex_coordinates[k6 + 3]; |
---|
167 | x2 = vertex_coordinates[k6 + 4]; |
---|
168 | y2 = vertex_coordinates[k6 + 5]; |
---|
169 | triangle * T = new_triangle(k,x0,y0,x1,y1,x2,y2); |
---|
170 | quad_tree_insert_triangle(tree,T); |
---|
171 | } |
---|
172 | |
---|
173 | // return pointer to new tree struct |
---|
174 | return tree; |
---|
175 | |
---|
176 | } |
---|
177 | |
---|
178 | // Builds the AtA and Atz interpolation matrix |
---|
179 | // and residual. Uses a quad_tree for fast access to the triangles of the mesh. |
---|
180 | // This function takes a list of point coordinates, and associated point values |
---|
181 | // (for any number of attributes). |
---|
182 | int _build_matrix_AtA_Atz_points(int N, long * triangles, |
---|
183 | double * point_coordinates, double * point_values, |
---|
184 | int zdims, int npts, |
---|
185 | sparse_dok * AtA, |
---|
186 | double ** Atz,quad_tree * quadtree) |
---|
187 | { |
---|
188 | |
---|
189 | |
---|
190 | |
---|
191 | int k; |
---|
192 | |
---|
193 | |
---|
194 | |
---|
195 | int i,w; |
---|
196 | for(w=0;w<zdims;w++){ |
---|
197 | for(i=0;i<N;i++){ |
---|
198 | Atz[w][i]=0; |
---|
199 | } |
---|
200 | } |
---|
201 | |
---|
202 | edge_key_t key; |
---|
203 | |
---|
204 | |
---|
205 | |
---|
206 | |
---|
207 | #pragma omp parallel for private(k,i,key,w) |
---|
208 | for(k=0;k<npts;k++){ |
---|
209 | |
---|
210 | |
---|
211 | double x = point_coordinates[2*k]; |
---|
212 | double y = point_coordinates[2*k+1]; |
---|
213 | triangle * T = search(quadtree,x,y); |
---|
214 | |
---|
215 | if(T!=NULL){ |
---|
216 | double * sigma = calculate_sigma(T,x,y); |
---|
217 | int js[3]; |
---|
218 | for(i=0;i<3;i++){ |
---|
219 | js[i]=triangles[3*(T->index)+i]; |
---|
220 | } |
---|
221 | |
---|
222 | #pragma omp critical |
---|
223 | { |
---|
224 | for(i=0;i<3;i++){ |
---|
225 | |
---|
226 | for(w=0;w<zdims;w++){ |
---|
227 | Atz[w][js[i]] += sigma[i]*point_values[zdims*k+w]; |
---|
228 | } |
---|
229 | |
---|
230 | for(w=0;w<3;w++){ |
---|
231 | |
---|
232 | key.i=js[i]; |
---|
233 | key.j=js[w]; |
---|
234 | |
---|
235 | add_dok_entry(AtA,key,sigma[i]*sigma[w]); |
---|
236 | } |
---|
237 | } |
---|
238 | } |
---|
239 | free(sigma); |
---|
240 | sigma=NULL; |
---|
241 | |
---|
242 | } |
---|
243 | } |
---|
244 | |
---|
245 | return 0; |
---|
246 | } |
---|
247 | |
---|
248 | // Combines two sparse_dok matricies and two vectors of doubles. |
---|
249 | void _combine_partial_AtA_Atz(sparse_dok * dok_AtA1,sparse_dok * dok_AtA2, |
---|
250 | double* Atz1, |
---|
251 | double* Atz2, |
---|
252 | int n, int zdim){ |
---|
253 | |
---|
254 | add_sparse_dok(dok_AtA1,1,dok_AtA2,1); |
---|
255 | |
---|
256 | int i; |
---|
257 | for(i=0;i<n*zdim;i++){ |
---|
258 | Atz1[i]+=Atz2[i]; |
---|
259 | } |
---|
260 | } |
---|
261 | |
---|
262 | // --------------------------- Utilities ----------------------------------- |
---|
263 | |
---|
264 | // Converts a double array into a PyList object for return to python. |
---|
265 | static PyObject *c_double_array_to_list(double * mat,int cols){ |
---|
266 | |
---|
267 | int j; |
---|
268 | |
---|
269 | PyObject *lst; |
---|
270 | |
---|
271 | lst = PyList_New(cols); |
---|
272 | if (!lst) return NULL; |
---|
273 | for (j = 0; j < cols; j++) { |
---|
274 | PyObject *num = PyFloat_FromDouble(mat[j]); |
---|
275 | if (!num) { |
---|
276 | Py_DECREF(lst); |
---|
277 | return NULL; |
---|
278 | } |
---|
279 | PyList_SET_ITEM(lst, j, num); // reference to num stolen |
---|
280 | } |
---|
281 | |
---|
282 | return lst; |
---|
283 | } |
---|
284 | |
---|
285 | // Converts a int array into a PyList object for return to python. |
---|
286 | static PyObject *c_int_array_to_list(int * mat,int cols){ |
---|
287 | |
---|
288 | int j; |
---|
289 | |
---|
290 | PyObject *lst; |
---|
291 | |
---|
292 | lst = PyList_New(cols); |
---|
293 | if (!lst) return NULL; |
---|
294 | for (j = 0; j < cols; j++) { |
---|
295 | PyObject *num = PyInt_FromLong(mat[j]); |
---|
296 | if (!num) { |
---|
297 | Py_DECREF(lst); |
---|
298 | return NULL; |
---|
299 | } |
---|
300 | PyList_SET_ITEM(lst, j, num); // reference to num stolen |
---|
301 | } |
---|
302 | |
---|
303 | return lst; |
---|
304 | } |
---|
305 | |
---|
306 | // ---------------------------------------------------------------------------- |
---|
307 | |
---|
308 | // If using python 2.7.3 or later, build with PyCapsules |
---|
309 | |
---|
310 | // Delete capsule containing a quad tree - name of capsule must be exactly |
---|
311 | // "quad tree". |
---|
312 | |
---|
313 | #ifdef PYVERSION273 |
---|
314 | |
---|
315 | |
---|
316 | void delete_quad_tree_cap(PyObject * cap){ |
---|
317 | quad_tree * kill = (quad_tree*) PyCapsule_GetPointer(cap,"quad tree"); |
---|
318 | if(kill!=NULL){ |
---|
319 | delete_quad_tree(kill); |
---|
320 | } else{ |
---|
321 | } |
---|
322 | } |
---|
323 | |
---|
324 | // Delete capsule containing a sparse_dok - name of capsule must be exactly |
---|
325 | // "sparse dok". |
---|
326 | void delete_dok_cap(PyObject *cap){ |
---|
327 | |
---|
328 | sparse_dok * kill = (sparse_dok*) PyCapsule_GetPointer(cap,"sparse dok"); |
---|
329 | if(kill!=NULL){ |
---|
330 | delete_dok_matrix(kill); |
---|
331 | } |
---|
332 | |
---|
333 | } |
---|
334 | |
---|
335 | #else |
---|
336 | |
---|
337 | // If using python earlier version, build with PyCObject |
---|
338 | |
---|
339 | // Delete cobj containing a quad tree |
---|
340 | void delete_quad_tree_cobj(void * cobj){ |
---|
341 | quad_tree * kill = (quad_tree*) cobj; |
---|
342 | if(kill!=NULL){ |
---|
343 | delete_quad_tree(kill); |
---|
344 | } |
---|
345 | } |
---|
346 | |
---|
347 | // Delete cobj containing a sparse_dok |
---|
348 | void delete_dok_cobj(void *cobj){ |
---|
349 | |
---|
350 | sparse_dok * kill = (sparse_dok*) cobj; |
---|
351 | if(kill!=NULL){ |
---|
352 | delete_dok_matrix(kill); |
---|
353 | } |
---|
354 | |
---|
355 | } |
---|
356 | |
---|
357 | #endif |
---|
358 | |
---|
359 | //----------------------- PYTHON WRAPPER FUNCTION ----------------------------- |
---|
360 | |
---|
361 | // Parses python information about triangles |
---|
362 | // and vertex coordiantes to build a quad tree for efficient location of points in |
---|
363 | // the triangle mesh. |
---|
364 | PyObject *build_quad_tree(PyObject *self, PyObject *args) { |
---|
365 | |
---|
366 | |
---|
367 | |
---|
368 | PyArrayObject *triangles; |
---|
369 | PyArrayObject *vertex_coordinates; |
---|
370 | PyArrayObject *extents; |
---|
371 | |
---|
372 | |
---|
373 | // Convert Python arguments to C |
---|
374 | if (!PyArg_ParseTuple(args, "OOO", |
---|
375 | &triangles, |
---|
376 | &vertex_coordinates, |
---|
377 | &extents |
---|
378 | )) { |
---|
379 | PyErr_SetString(PyExc_RuntimeError, |
---|
380 | "fitsmooth.c: could not parse input"); |
---|
381 | return NULL; |
---|
382 | } |
---|
383 | |
---|
384 | CHECK_C_CONTIG(triangles); |
---|
385 | CHECK_C_CONTIG(vertex_coordinates); |
---|
386 | CHECK_C_CONTIG(extents); |
---|
387 | |
---|
388 | int n = triangles->dimensions[0]; |
---|
389 | |
---|
390 | #ifdef PYVERSION273 |
---|
391 | |
---|
392 | return PyCapsule_New((void*) _build_quad_tree(n, |
---|
393 | (long*) triangles -> data, |
---|
394 | (double*) vertex_coordinates -> data, |
---|
395 | (double*) extents -> data), |
---|
396 | "quad tree", |
---|
397 | &delete_quad_tree_cap); |
---|
398 | |
---|
399 | #else |
---|
400 | |
---|
401 | return PyCObject_FromVoidPtr((void*) _build_quad_tree(n, |
---|
402 | (long*) triangles -> data, |
---|
403 | (double*) vertex_coordinates -> data, |
---|
404 | (double*) extents -> data), |
---|
405 | &delete_quad_tree_cobj); |
---|
406 | #endif |
---|
407 | } |
---|
408 | |
---|
409 | // Builds the smoothing matrix D. Parses |
---|
410 | // input mesh information from python and then builds a sparse_dok, returning |
---|
411 | // a capsule object wrapping a pointer to this struct. |
---|
412 | PyObject *build_smoothing_matrix(PyObject *self, PyObject *args) { |
---|
413 | |
---|
414 | |
---|
415 | PyArrayObject *triangles; |
---|
416 | PyArrayObject *areas; |
---|
417 | PyArrayObject *vertex_coordinates; |
---|
418 | int err; |
---|
419 | |
---|
420 | // Convert Python arguments to C |
---|
421 | if (!PyArg_ParseTuple(args, "OOO" , |
---|
422 | &triangles, |
---|
423 | &areas, |
---|
424 | &vertex_coordinates |
---|
425 | )) { |
---|
426 | PyErr_SetString(PyExc_RuntimeError, |
---|
427 | "fitsmooth.c: could not parse input"); |
---|
428 | return NULL; |
---|
429 | } |
---|
430 | |
---|
431 | CHECK_C_CONTIG(triangles); |
---|
432 | CHECK_C_CONTIG(areas); |
---|
433 | CHECK_C_CONTIG(vertex_coordinates); |
---|
434 | |
---|
435 | int n = triangles->dimensions[0]; |
---|
436 | |
---|
437 | |
---|
438 | sparse_dok * smoothing_mat; // Should be an input argument? |
---|
439 | smoothing_mat = make_dok(); |
---|
440 | |
---|
441 | err = _build_smoothing_matrix(n, |
---|
442 | (long*) triangles -> data, |
---|
443 | (double*) areas -> data, |
---|
444 | (double*) vertex_coordinates -> data, |
---|
445 | (int *) vertex_coordinates -> strides, |
---|
446 | smoothing_mat); |
---|
447 | |
---|
448 | |
---|
449 | if (err != 0) { |
---|
450 | PyErr_SetString(PyExc_RuntimeError, |
---|
451 | "Unknown Error"); |
---|
452 | return NULL; |
---|
453 | } |
---|
454 | |
---|
455 | #ifdef PYVERSION273 |
---|
456 | |
---|
457 | return PyCapsule_New((void*) smoothing_mat, |
---|
458 | "sparse dok", |
---|
459 | &delete_dok_cap); |
---|
460 | |
---|
461 | #else |
---|
462 | |
---|
463 | return PyCObject_FromVoidPtr((void*) smoothing_mat, |
---|
464 | &delete_dok_cobj); |
---|
465 | |
---|
466 | #endif |
---|
467 | |
---|
468 | } |
---|
469 | |
---|
470 | // Builds AtA and Atz directly from a list of |
---|
471 | // points and values Returns a pointer to the sparse dok matrix AtA wrapped in |
---|
472 | // a capsule object and a python list for the array Atz. |
---|
473 | PyObject *build_matrix_AtA_Atz_points(PyObject *self, PyObject *args) { |
---|
474 | |
---|
475 | |
---|
476 | PyArrayObject *triangles; |
---|
477 | PyArrayObject *point_coordinates; |
---|
478 | PyArrayObject *z; |
---|
479 | int N; // Number of triangles |
---|
480 | int err; |
---|
481 | int npts; |
---|
482 | int zdims; |
---|
483 | PyObject *tree; |
---|
484 | |
---|
485 | // Convert Python arguments to C |
---|
486 | if (!PyArg_ParseTuple(args, "OiOOOii",&tree, &N, |
---|
487 | &triangles, |
---|
488 | &point_coordinates, |
---|
489 | &z, |
---|
490 | &zdims, |
---|
491 | &npts |
---|
492 | )) { |
---|
493 | PyErr_SetString(PyExc_RuntimeError, |
---|
494 | "fitsmooth.c: could not parse input"); |
---|
495 | return NULL; |
---|
496 | } |
---|
497 | |
---|
498 | |
---|
499 | CHECK_C_CONTIG(triangles); |
---|
500 | CHECK_C_CONTIG(point_coordinates); |
---|
501 | CHECK_C_CONTIG(z); |
---|
502 | |
---|
503 | #ifdef PYVERSION273 |
---|
504 | quad_tree * quadtree = (quad_tree*) PyCapsule_GetPointer(tree,"quad tree"); |
---|
505 | #else |
---|
506 | quad_tree * quadtree = (quad_tree*) PyCObject_AsVoidPtr(tree); |
---|
507 | #endif |
---|
508 | |
---|
509 | sparse_dok * dok_AtA; // Should be an input argument? |
---|
510 | dok_AtA = make_dok(); |
---|
511 | double ** Atz = malloc(sizeof(double*)*zdims); |
---|
512 | int i; |
---|
513 | for(i=0;i<zdims;i++){ |
---|
514 | Atz[i] = malloc(sizeof(double)*N); |
---|
515 | } |
---|
516 | |
---|
517 | err = _build_matrix_AtA_Atz_points(N,(long*) triangles->data, |
---|
518 | (double*) point_coordinates->data, |
---|
519 | (double*) z->data, |
---|
520 | zdims, |
---|
521 | npts, |
---|
522 | dok_AtA, |
---|
523 | Atz, |
---|
524 | quadtree); |
---|
525 | |
---|
526 | |
---|
527 | if (err != 0) { |
---|
528 | PyErr_SetString(PyExc_RuntimeError, |
---|
529 | "Unknown Error"); |
---|
530 | return NULL; |
---|
531 | } |
---|
532 | |
---|
533 | #ifdef PYVERSION273 |
---|
534 | PyObject * AtA_cap = PyCapsule_New((void*) dok_AtA, |
---|
535 | "sparse dok", |
---|
536 | &delete_dok_cap); |
---|
537 | #else |
---|
538 | PyObject * AtA_cap = PyCObject_FromVoidPtr((void*) dok_AtA, |
---|
539 | &delete_dok_cobj); |
---|
540 | #endif |
---|
541 | PyObject * Atz_ret = PyList_New(zdims); |
---|
542 | for(i=0;i<zdims;i++){ |
---|
543 | PyList_SET_ITEM(Atz_ret,i,c_double_array_to_list(Atz[i],N)); |
---|
544 | free(Atz[i]); |
---|
545 | } |
---|
546 | free(Atz); |
---|
547 | |
---|
548 | PyObject *lst = PyList_New(2); |
---|
549 | PyList_SET_ITEM(lst, 0, AtA_cap); |
---|
550 | PyList_SET_ITEM(lst, 1, Atz_ret); |
---|
551 | return lst; |
---|
552 | |
---|
553 | } |
---|
554 | |
---|
555 | // Searches a quad tree struct for the triangle containing a given point, |
---|
556 | // returns the sigma values produced by this point and the triangle found, |
---|
557 | // and the triangle index. Found is returned as 0 if no triangle is found |
---|
558 | // and 1 otherwise. |
---|
559 | PyObject *individual_tree_search(PyObject *self, PyObject *args) { |
---|
560 | |
---|
561 | // Setting up variables to parse input |
---|
562 | PyObject *tree; |
---|
563 | PyArrayObject *point; |
---|
564 | |
---|
565 | // Convert Python arguments to C |
---|
566 | if (!PyArg_ParseTuple(args, "OO",&tree, &point |
---|
567 | )) { |
---|
568 | PyErr_SetString(PyExc_RuntimeError, |
---|
569 | "fitsmooth.c: could not parse input"); |
---|
570 | return NULL; |
---|
571 | } |
---|
572 | |
---|
573 | CHECK_C_CONTIG(point); |
---|
574 | |
---|
575 | #ifdef PYVERSION273 |
---|
576 | quad_tree * quadtree = (quad_tree*) PyCapsule_GetPointer(tree,"quad tree"); |
---|
577 | #else |
---|
578 | quad_tree * quadtree = (quad_tree*) PyCObject_AsVoidPtr(tree); |
---|
579 | #endif |
---|
580 | |
---|
581 | double xp,yp; |
---|
582 | if(PyArray_TYPE(point)==7){ |
---|
583 | int *pointd = (int*)point->data; |
---|
584 | xp = (double)(pointd[0]); |
---|
585 | yp = (double)(pointd[1]); |
---|
586 | } |
---|
587 | else{ |
---|
588 | double *pointd = (double*)point->data; |
---|
589 | xp = pointd[0]; |
---|
590 | yp = pointd[1]; |
---|
591 | } |
---|
592 | |
---|
593 | triangle * T = search(quadtree,xp,yp); |
---|
594 | PyObject * sigmalist = PyList_New(3); |
---|
595 | long found; |
---|
596 | long index; |
---|
597 | |
---|
598 | if(T!=NULL){ |
---|
599 | double * sigma = calculate_sigma(T,xp,yp); |
---|
600 | sigmalist = c_double_array_to_list(sigma,3); |
---|
601 | free(sigma); |
---|
602 | found = 1; |
---|
603 | index = (long)T->index; |
---|
604 | }else{ |
---|
605 | double sigma[3]; |
---|
606 | sigma[0]=sigma[1]=sigma[2]=-1; |
---|
607 | sigmalist = c_double_array_to_list(sigma,3); |
---|
608 | index = -10; |
---|
609 | found = 0; |
---|
610 | } |
---|
611 | PyObject * retlist = PyList_New(3); |
---|
612 | PyList_SET_ITEM(retlist,0,PyInt_FromLong(found)); |
---|
613 | PyList_SET_ITEM(retlist,1,sigmalist); |
---|
614 | PyList_SET_ITEM(retlist,2,PyInt_FromLong(index)); |
---|
615 | return retlist; |
---|
616 | } |
---|
617 | |
---|
618 | // Returns the total number of triangles stored in quad_tree. |
---|
619 | // Takes a capsule object holding a pointer to the quad_tree as input. |
---|
620 | // |
---|
621 | // PADARN NOTE: This function was only added for the purpose of passing unit |
---|
622 | // tests. Perhaps this, and other tree functions, should be moved to another |
---|
623 | // file to provide a more comprehensive python interface to the tree structure. |
---|
624 | PyObject *items_in_tree(PyObject *self, PyObject *args) { |
---|
625 | |
---|
626 | // Setting up variables to parse input |
---|
627 | PyObject *tree; // capsule to hold quad_tree pointer |
---|
628 | |
---|
629 | // Convert Python arguments to C |
---|
630 | if (!PyArg_ParseTuple(args, "O",&tree |
---|
631 | )) { |
---|
632 | PyErr_SetString(PyExc_RuntimeError, |
---|
633 | "fitsmooth.c: could not parse input"); |
---|
634 | return NULL; |
---|
635 | } |
---|
636 | |
---|
637 | // Extract quad_tree pointer from capsule |
---|
638 | #ifdef PYVERSION273 |
---|
639 | quad_tree * quadtree = (quad_tree*) PyCapsule_GetPointer(tree,"quad tree"); |
---|
640 | #else |
---|
641 | quad_tree * quadtree = (quad_tree*) PyCObject_AsVoidPtr(tree); |
---|
642 | #endif |
---|
643 | |
---|
644 | // Return the number of elements in the tree (stored in struct) |
---|
645 | return PyInt_FromLong((long)quadtree->count); |
---|
646 | |
---|
647 | } |
---|
648 | |
---|
649 | // Returns the total number of nodes stored in quad_tree. |
---|
650 | // Takes a capsule object holding a pointer to the quad_tree as input. |
---|
651 | PyObject *nodes_in_tree(PyObject *self, PyObject *args) { |
---|
652 | |
---|
653 | // Setting up variables to parse input |
---|
654 | PyObject *tree; // capsule to hold quad_tree pointer |
---|
655 | |
---|
656 | // Convert Python arguments to C |
---|
657 | if (!PyArg_ParseTuple(args, "O",&tree |
---|
658 | )) { |
---|
659 | PyErr_SetString(PyExc_RuntimeError, |
---|
660 | "fitsmooth.c: could not parse input"); |
---|
661 | return NULL; |
---|
662 | } |
---|
663 | |
---|
664 | // Extract quad_tree pointer from capsule |
---|
665 | #ifdef PYVERSION273 |
---|
666 | quad_tree * quadtree = (quad_tree*) PyCapsule_GetPointer(tree,"quad tree"); |
---|
667 | #else |
---|
668 | quad_tree * quadtree = (quad_tree*) PyCObject_AsVoidPtr(tree); |
---|
669 | #endif |
---|
670 | |
---|
671 | // Return the number of elements in the tree (stored in struct) |
---|
672 | return PyInt_FromLong((long)quad_tree_node_count(quadtree)); |
---|
673 | |
---|
674 | } |
---|
675 | |
---|
676 | // Combines two sparse_dok and two double arrays together |
---|
677 | // which represent partial completion of the AtA and Atz matrices, when they are build |
---|
678 | // in parts due to points being read in blocks in python. Result is stored in the first |
---|
679 | // sparse_dok and double array. |
---|
680 | // Function takes as arguments two capsule objects holding pointers to the sparse_dok |
---|
681 | // structs, and two numpy double arrays. Aslo the size of the double arrays, and the number |
---|
682 | // of columns (variables) for the array Atz are passed as arguments. |
---|
683 | // |
---|
684 | // PADARN NOTE: Blocking the points in python is far slower than reading them directly |
---|
685 | // in c and bypassing the need for this. |
---|
686 | PyObject *combine_partial_AtA_Atz(PyObject *self, PyObject *args) { |
---|
687 | |
---|
688 | // Setting up variables to parse input |
---|
689 | PyObject *AtA_cap1, *AtA_cap2; |
---|
690 | PyArrayObject *Atz1, *Atz2; |
---|
691 | int n,zdim; // number of nodes |
---|
692 | |
---|
693 | // Convert Python arguments to C |
---|
694 | if (!PyArg_ParseTuple(args, "OOOOii",&AtA_cap1, &AtA_cap2, |
---|
695 | &Atz1, &Atz2, &zdim, &n |
---|
696 | )) { |
---|
697 | PyErr_SetString(PyExc_RuntimeError, |
---|
698 | "fitsmooth.c: could not parse input"); |
---|
699 | return NULL; |
---|
700 | } |
---|
701 | |
---|
702 | // Get pointers to sparse_dok objects from capsules |
---|
703 | #ifdef PYVERSION273 |
---|
704 | sparse_dok * dok_AtA1 = (sparse_dok*) PyCapsule_GetPointer(AtA_cap1,"sparse dok"); |
---|
705 | sparse_dok * dok_AtA2 = (sparse_dok*) PyCapsule_GetPointer(AtA_cap2,"sparse dok"); |
---|
706 | #else |
---|
707 | sparse_dok * dok_AtA1 = (sparse_dok*) PyCObject_AsVoidPtr(AtA_cap1); |
---|
708 | sparse_dok * dok_AtA2 = (sparse_dok*) PyCObject_AsVoidPtr(AtA_cap2); |
---|
709 | #endif |
---|
710 | |
---|
711 | // Combine the partial AtA and Atz |
---|
712 | _combine_partial_AtA_Atz(dok_AtA1,dok_AtA2, |
---|
713 | (double*) Atz1->data, |
---|
714 | (double*) Atz2->data, |
---|
715 | n, zdim); |
---|
716 | |
---|
717 | // Return nothing interesting |
---|
718 | return Py_BuildValue(""); |
---|
719 | } |
---|
720 | |
---|
721 | // Converts a sparse_dok matrix to a full non-compressed matrix expressed |
---|
722 | // as a list of lists (python). Takes as input a capsule object containing a pointer to the |
---|
723 | // sparse_dok object. Also takes an integer n as input, specifying the (n x n) size of the |
---|
724 | // matrix. Matrix is assumed square. |
---|
725 | // |
---|
726 | // PADARN NOTE: This function does not seem particularly sensible, but was required for the |
---|
727 | // purposes of passing unit tests. |
---|
728 | PyObject *return_full_D(PyObject *self, PyObject *args) { |
---|
729 | |
---|
730 | // Setting up variables to parse input |
---|
731 | PyObject *D_cap; // capsule object holding sparse_dok pointer |
---|
732 | int n; // number of matrix columns/rows |
---|
733 | |
---|
734 | // Convert Python arguments to C |
---|
735 | if (!PyArg_ParseTuple(args, "Oi",&D_cap, &n |
---|
736 | )) { |
---|
737 | PyErr_SetString(PyExc_RuntimeError, |
---|
738 | "fitsmooth.return_full_D: could not parse input"); |
---|
739 | return NULL; |
---|
740 | } |
---|
741 | |
---|
742 | // Get pointer to spars_dok struct |
---|
743 | #ifdef PYVERSION273 |
---|
744 | sparse_dok * D_mat = (sparse_dok*) PyCapsule_GetPointer(D_cap,"sparse dok"); |
---|
745 | #else |
---|
746 | sparse_dok * D_mat = (sparse_dok*) PyCObject_AsVoidPtr(D_cap); |
---|
747 | #endif |
---|
748 | |
---|
749 | // Check to make sure that the specified size of the matrix is at least as big |
---|
750 | // as the entries stored in the sparse_dok. |
---|
751 | if (D_mat->num_rows>n || get_dok_rows(D_mat)>n){ |
---|
752 | PyErr_SetString(PyExc_RuntimeError, |
---|
753 | "fitsmooth.return_full_D: sparse_dok is bigger than size specified for return."); |
---|
754 | return NULL; |
---|
755 | } |
---|
756 | |
---|
757 | // Build new python list containing the full matrix to return |
---|
758 | PyObject *ret_D = PyList_New(n); |
---|
759 | int i,j; |
---|
760 | edge_key_t key; |
---|
761 | edge_t *s; |
---|
762 | for(i=0;i<n;i++) |
---|
763 | { |
---|
764 | PyObject *temp = PyList_New(n); |
---|
765 | for(j=0;j<n;j++) |
---|
766 | { |
---|
767 | key.i=i; |
---|
768 | key.j=j; |
---|
769 | s = find_dok_entry(D_mat,key); |
---|
770 | if (s){ |
---|
771 | PyList_SET_ITEM(temp,j,PyFloat_FromDouble(s->entry)); |
---|
772 | }else{ |
---|
773 | PyList_SET_ITEM(temp,j,PyFloat_FromDouble(0)); |
---|
774 | } |
---|
775 | PyList_SET_ITEM(ret_D,i,temp); |
---|
776 | } |
---|
777 | } |
---|
778 | |
---|
779 | return ret_D; |
---|
780 | } |
---|
781 | |
---|
782 | // Takes as input two capsule objects corresponding to the two matricies |
---|
783 | // D and AtA, along with a regularization coefficient alpha. |
---|
784 | // Returns three matricies corresponding to the CSR format storage of the matrix resulting |
---|
785 | // from the computation B = AtA + alpha * D. |
---|
786 | // |
---|
787 | // Capsule objects are not freed and are left alone to be freed as they naturally go out of |
---|
788 | // scope in Python. The temporary CSR matrix created is cleaned up. |
---|
789 | PyObject *build_matrix_B(PyObject *self, PyObject *args) { |
---|
790 | |
---|
791 | |
---|
792 | // Setting up variables to parse input |
---|
793 | double alpha; // Regularization coefficient |
---|
794 | PyObject *smoothing_mat_cap; // Capsule storing D pointer (sparse_dok struct) |
---|
795 | PyObject *AtA_cap; // Capsule storing AtA pointer (sparse_dok struct) |
---|
796 | |
---|
797 | // Convert Python arguments to C |
---|
798 | if (!PyArg_ParseTuple(args, "OOd",&smoothing_mat_cap, &AtA_cap, &alpha |
---|
799 | )) { |
---|
800 | PyErr_SetString(PyExc_RuntimeError, |
---|
801 | "fitsmooth.build_matrix_B: could not parse input"); |
---|
802 | return NULL; |
---|
803 | } |
---|
804 | |
---|
805 | // Extract pointers to c structs from capsules |
---|
806 | #ifdef PYVERSION273 |
---|
807 | sparse_dok * smoothing_mat = (sparse_dok*) PyCapsule_GetPointer(smoothing_mat_cap,"sparse dok"); |
---|
808 | sparse_dok * dok_AtA = (sparse_dok*) PyCapsule_GetPointer(AtA_cap,"sparse dok"); |
---|
809 | #else |
---|
810 | sparse_dok * smoothing_mat = (sparse_dok*) PyCObject_AsVoidPtr(smoothing_mat_cap); |
---|
811 | sparse_dok * dok_AtA = (sparse_dok*) PyCObject_AsVoidPtr(AtA_cap); |
---|
812 | #endif |
---|
813 | |
---|
814 | // Add two sparse_dok matrices |
---|
815 | add_sparse_dok(smoothing_mat,alpha,dok_AtA,1); |
---|
816 | |
---|
817 | // Create sparse_csr matrix and convert result to this format |
---|
818 | sparse_csr * B; |
---|
819 | B = make_csr(); |
---|
820 | convert_to_csr_ptr(B,smoothing_mat); |
---|
821 | |
---|
822 | // Extract the sparse_csr data to be returned as python lists |
---|
823 | PyObject *data = c_double_array_to_list(B->data, |
---|
824 | B->num_entries); |
---|
825 | PyObject *colind = c_int_array_to_list(B->colind, |
---|
826 | B->num_entries); |
---|
827 | PyObject *row_ptr = c_int_array_to_list(B->row_ptr, |
---|
828 | B->num_rows); |
---|
829 | |
---|
830 | // Build python list for return |
---|
831 | PyObject *lst = PyList_New(3); |
---|
832 | PyList_SET_ITEM(lst, 0, data); |
---|
833 | PyList_SET_ITEM(lst, 1, colind); |
---|
834 | PyList_SET_ITEM(lst, 2, row_ptr); |
---|
835 | |
---|
836 | // Clean up |
---|
837 | delete_csr_matrix(B); |
---|
838 | |
---|
839 | return lst; |
---|
840 | |
---|
841 | |
---|
842 | |
---|
843 | } |
---|
844 | //------------------------------------------------------------------------------ |
---|
845 | |
---|
846 | |
---|
847 | // ------------------------------ PYTHON GLUE ---------------------------------- |
---|
848 | |
---|
849 | //============================================================================== |
---|
850 | // Structures to allow calling from python |
---|
851 | //============================================================================== |
---|
852 | |
---|
853 | // Method table for python module |
---|
854 | static struct PyMethodDef MethodTable[] = { |
---|
855 | {"items_in_tree",items_in_tree, METH_VARARGS, "Print out"}, |
---|
856 | {"nodes_in_tree",nodes_in_tree, METH_VARARGS, "Print out"}, |
---|
857 | {"return_full_D",return_full_D, METH_VARARGS, "Print out"}, |
---|
858 | {"build_matrix_B",build_matrix_B, METH_VARARGS, "Print out"}, |
---|
859 | {"build_quad_tree",build_quad_tree, METH_VARARGS, "Print out"}, |
---|
860 | {"build_smoothing_matrix",build_smoothing_matrix, METH_VARARGS, "Print out"}, |
---|
861 | {"build_matrix_AtA_Atz_points",build_matrix_AtA_Atz_points, METH_VARARGS, "Print out"}, |
---|
862 | {"combine_partial_AtA_Atz",combine_partial_AtA_Atz, METH_VARARGS, "Print out"}, |
---|
863 | {"individual_tree_search",individual_tree_search, METH_VARARGS, "Print out"}, |
---|
864 | {NULL, NULL, 0, NULL} // sentinel |
---|
865 | }; |
---|
866 | |
---|
867 | |
---|
868 | // Module initialisation |
---|
869 | void initfitsmooth(void){ |
---|
870 | Py_InitModule("fitsmooth", MethodTable); |
---|
871 | import_array(); // Necessary for handling of NumPY structures |
---|
872 | |
---|
873 | } |
---|
874 | |
---|
875 | // -------------------------------------------------------------------------------- |
---|