File: _superluobject.c

package info (click to toggle)
python-scipy 0.7.2%2Bdfsg1-1
  • links: PTS, VCS
  • area: main
  • in suites: squeeze
  • size: 28,500 kB
  • ctags: 36,081
  • sloc: cpp: 216,880; fortran: 76,016; python: 71,576; ansic: 62,118; makefile: 243; sh: 17
file content (416 lines) | stat: -rw-r--r-- 11,910 bytes parent folder | download | duplicates (4)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416

#include "Python.h"
#include "SuperLU/SRC/zsp_defs.h"
#define NO_IMPORT_ARRAY
#include "_superluobject.h"
#include <setjmp.h>

extern jmp_buf _superlu_py_jmpbuf;

/*********************************************************************** 
 * SciPyLUObject methods
 */

static char solve_doc[] = "x = self.solve(b, trans)\n\
\n\
solves linear system of equations with one or sereral right hand sides.\n\
\n\
parameters\n\
----------\n\
\n\
b        array, right hand side(s) of equation\n\
x        array, solution vector(s)\n\
trans    'N': solve A   * x == b\n\
         'T': solve A^T * x == b\n\
         'H': solve A^H * x == b (not yet implemented)\n\
         (optional, default value 'N')\n\
";

static PyObject *
SciPyLU_solve(SciPyLUObject *self, PyObject *args, PyObject *kwds) {
  PyArrayObject *b, *x=NULL;
  SuperMatrix B;
  char itrans = 'N';
  int info;
  trans_t trans;
  SuperLUStat_t stat;

  static char *kwlist[] = {"rhs","trans",NULL};

  if (!PyArg_ParseTupleAndKeywords(args, kwds, "O!|c", kwlist,
                                   &PyArray_Type, &b, 
                                   &itrans))
    return NULL;

  /* solve transposed system: matrix was passed row-wise instead of column-wise */
  if (itrans == 'n' || itrans == 'N')
      trans = NOTRANS;
  else if (itrans == 't' || itrans == 'T')
      trans = TRANS;
  else if (itrans == 'h' || itrans == 'H')
      trans = CONJ;
  else {
    PyErr_SetString(PyExc_ValueError, "trans must be N, T, or H");
    return NULL;
  }

  if ((x = (PyArrayObject *) \
       PyArray_CopyFromObject((PyObject *)b,self->type,1,2))==NULL) return NULL;

  if (b->dimensions[0] != self->n) goto fail;


  if (setjmp(_superlu_py_jmpbuf)) goto fail; 

  if (DenseSuper_from_Numeric(&B, (PyObject *)x)) goto fail;

  StatInit(&stat);

  /* Solve the system, overwriting vector x. */
  switch(self->type) {
  case PyArray_FLOAT:
      sgstrs(trans, &self->L, &self->U, self->perm_c, self->perm_r, &B, &stat, &info);      
      break;
  case PyArray_DOUBLE:
      dgstrs(trans, &self->L, &self->U, self->perm_c, self->perm_r, &B, &stat, &info);      
      break;
  case PyArray_CFLOAT:
      cgstrs(trans, &self->L, &self->U, self->perm_c, self->perm_r, &B, &stat, &info);      
      break;
  case PyArray_CDOUBLE:
      zgstrs(trans, &self->L, &self->U, self->perm_c, self->perm_r, &B, &stat, &info); 
      break;
  default:
      PyErr_SetString(PyExc_TypeError, "Invalid type for array.");
      goto fail;
  }

  if (info) { 
      PyErr_SetString(PyExc_SystemError, "gstrs was called with invalid arguments");
      goto fail;
  }
  
  /* free memory */
  Destroy_SuperMatrix_Store(&B);
  StatFree(&stat);
  return (PyObject *)x;

 fail:
  Destroy_SuperMatrix_Store(&B);  
  StatFree(&stat);
  Py_XDECREF(x);
  return NULL;
}

/** table of object methods
 */
PyMethodDef SciPyLU_methods[] = {
  {"solve", (PyCFunction)SciPyLU_solve, METH_VARARGS|METH_KEYWORDS, solve_doc},
  {NULL, NULL}			/* sentinel */
};


/*********************************************************************** 
 * SciPySuperLUType methods
 */

static void
SciPyLU_dealloc(SciPyLUObject *self)
{
  SUPERLU_FREE(self->perm_r);
  SUPERLU_FREE(self->perm_c);
  Destroy_SuperNode_Matrix(&self->L);
  Destroy_CompCol_Matrix(&self->U);
  PyObject_Del(self);
}

static PyObject *
SciPyLU_getattr(SciPyLUObject *self, char *name)
{
  if (strcmp(name, "shape") == 0)
    return Py_BuildValue("(i,i)", self->m, self->n);
  if (strcmp(name, "nnz") == 0)
    return Py_BuildValue("i", ((SCformat *)self->L.Store)->nnz + ((SCformat *)self->U.Store)->nnz);
  if (strcmp(name, "__members__") == 0) {
    char *members[] = {"shape", "nnz"};
    int i;

    PyObject *list = PyList_New(sizeof(members)/sizeof(char *));
    if (list != NULL) {
      for (i = 0; i < sizeof(members)/sizeof(char *); i ++)
	PyList_SetItem(list, i, PyString_FromString(members[i]));
      if (PyErr_Occurred()) {
	Py_DECREF(list);
	list = NULL;
      }
    }
    return list;
  }
  return Py_FindMethod(SciPyLU_methods, (PyObject *)self, name);
}


/***********************************************************************
 * SciPySuperLUType structure
 */

PyTypeObject SciPySuperLUType = {
  PyObject_HEAD_INIT(NULL)
  0,
  "factored_lu",
  sizeof(SciPyLUObject),
  0,
  (destructor)SciPyLU_dealloc,   /* tp_dealloc */
  0,				/* tp_print */
  (getattrfunc)SciPyLU_getattr,  /* tp_getattr */
  0,				/* tp_setattr */
  0,				/* tp_compare */
  0,				/* tp_repr */
  0,				/* tp_as_number*/
  0,				/* tp_as_sequence*/
  0,				/* tp_as_mapping*/
  0,				/* tp_hash */
};


int DenseSuper_from_Numeric(SuperMatrix *X, PyObject *PyX)
{
  int m, n, ldx, nd;
  PyArrayObject *aX;
  
  if (!PyArray_Check(PyX)) {
    PyErr_SetString(PyExc_TypeError, "dgssv: Second argument is not an array.");
    return -1;
  }

  aX = (PyArrayObject *)PyX;
  nd = aX->nd;

  if (nd == 1) {
    m = aX->dimensions[0];
    n = 1;
    ldx = m;
  }
  else {  /* nd == 2 */
    m = aX->dimensions[1];
    n = aX->dimensions[0];
    ldx = m;
  }
  
  if (setjmp(_superlu_py_jmpbuf)) return -1;
  else 
      switch (aX->descr->type_num) {
      case PyArray_FLOAT:
          sCreate_Dense_Matrix(X, m, n, (float *)aX->data, ldx, SLU_DN, SLU_S, SLU_GE);
          break;
      case PyArray_DOUBLE:
          dCreate_Dense_Matrix(X, m, n, (double *)aX->data, ldx, SLU_DN, SLU_D, SLU_GE);
          break;
      case PyArray_CFLOAT:
          cCreate_Dense_Matrix(X, m, n, (complex *)aX->data, ldx, SLU_DN, SLU_C, SLU_GE);
          break;
      case PyArray_CDOUBLE:
          zCreate_Dense_Matrix(X, m, n, (doublecomplex *)aX->data, ldx, SLU_DN, SLU_Z, SLU_GE);
          break;
      default:
          PyErr_SetString(PyExc_TypeError, "Invalid type for Numeric array.");
          return -1;  
      }
  
  return 0;
}

/* Natively handles Compressed Sparse Row and CSC */

int NRFormat_from_spMatrix(SuperMatrix *A, int m, int n, int nnz, PyArrayObject *nzvals, PyArrayObject *colind, PyArrayObject *rowptr, int typenum)
{
  int err = 0;
    
  err = (nzvals->descr->type_num != typenum);
  err += (nzvals->nd != 1);
  err += (nnz > nzvals->dimensions[0]);
  if (err) {
    PyErr_SetString(PyExc_TypeError, "Fourth argument must be a 1-D array at least as big as third argument.");
    return -1;
  }

  if (setjmp(_superlu_py_jmpbuf)) return -1;
  else 
      switch (nzvals->descr->type_num) {
      case PyArray_FLOAT:
          sCreate_CompRow_Matrix(A, m, n, nnz, (float *)nzvals->data, (int *)colind->data, \
                                 (int *)rowptr->data, SLU_NR, SLU_S, SLU_GE);
          break;
      case PyArray_DOUBLE:
          dCreate_CompRow_Matrix(A, m, n, nnz, (double *)nzvals->data, (int *)colind->data, \
                                 (int *)rowptr->data, SLU_NR, SLU_D, SLU_GE);
          break;
      case PyArray_CFLOAT:
          cCreate_CompRow_Matrix(A, m, n, nnz, (complex *)nzvals->data, (int *)colind->data, \
                                 (int *)rowptr->data, SLU_NR, SLU_C, SLU_GE);
          break;
      case PyArray_CDOUBLE:
          zCreate_CompRow_Matrix(A, m, n, nnz, (doublecomplex *)nzvals->data, (int *)colind->data, \
                                 (int *)rowptr->data, SLU_NR, SLU_Z, SLU_GE);
          break;
      default:
          PyErr_SetString(PyExc_TypeError, "Invalid type for array.");
          return -1;  
      }

  return 0;
}

int NCFormat_from_spMatrix(SuperMatrix *A, int m, int n, int nnz, PyArrayObject *nzvals, PyArrayObject *rowind, PyArrayObject *colptr, int typenum)
{
  int err=0;

  err = (nzvals->descr->type_num != typenum);
  err += (nzvals->nd != 1);
  err += (nnz > nzvals->dimensions[0]);
  if (err) {
    PyErr_SetString(PyExc_TypeError, "Fifth argument must be a 1-D array at least as big as fourth argument.");
    return -1;
  }


  if (setjmp(_superlu_py_jmpbuf)) return -1;
  else 
      switch (nzvals->descr->type_num) {
      case PyArray_FLOAT:
          sCreate_CompCol_Matrix(A, m, n, nnz, (float *)nzvals->data, (int *)rowind->data, \
                                 (int *)colptr->data, SLU_NC, SLU_S, SLU_GE);
          break;
      case PyArray_DOUBLE:
          dCreate_CompCol_Matrix(A, m, n, nnz, (double *)nzvals->data, (int *)rowind->data, \
                                 (int *)colptr->data, SLU_NC, SLU_D, SLU_GE);
          break;
      case PyArray_CFLOAT:
          cCreate_CompCol_Matrix(A, m, n, nnz, (complex *)nzvals->data, (int *)rowind->data, \
                                 (int *)colptr->data, SLU_NC, SLU_C, SLU_GE);
          break;
      case PyArray_CDOUBLE:
          zCreate_CompCol_Matrix(A, m, n, nnz, (doublecomplex *)nzvals->data, (int *)rowind->data, \
                                 (int *)colptr->data, SLU_NC, SLU_Z, SLU_GE);
          break;
      default:
          PyErr_SetString(PyExc_TypeError, "Invalid type for array.");
          return -1;  
      }

  return 0;
}

colperm_t superlu_module_getpermc(int permc_spec)
{
  switch(permc_spec) {
  case 0:
    return NATURAL;
  case 1:
    return MMD_ATA;
  case 2:
    return MMD_AT_PLUS_A;
  case 3:
    return COLAMD;
  }
  ABORT("Invalid input for permc_spec.");
  return NATURAL; /* compiler complains... */
}

PyObject *
newSciPyLUObject(SuperMatrix *A, double diag_pivot_thresh,
		 double drop_tol, int relax, int panel_size, int permc_spec,
                 int intype)
{

   /* A must be in SLU_NC format used by the factorization routine. */
  SciPyLUObject *self;
  SuperMatrix AC;     /* Matrix postmultiplied by Pc */
  int lwork = 0;
  int *etree=NULL;
  int info;
  int n;
  superlu_options_t options;
  SuperLUStat_t stat;
  
  n = A->ncol;

  /* Create SciPyLUObject */
  self = PyObject_New(SciPyLUObject, &SciPySuperLUType);
  if (self == NULL)
    return PyErr_NoMemory();
  self->m = A->nrow;
  self->n = n;
  self->perm_r = NULL;
  self->perm_c = NULL;
  self->type = intype;

  if (setjmp(_superlu_py_jmpbuf)) goto fail;
  
  /* Calculate and apply minimum degree ordering*/
  etree = intMalloc(n);
  self->perm_r = intMalloc(n);
  self->perm_c = intMalloc(n);
  
  set_default_options(&options);
  options.ColPerm=superlu_module_getpermc(permc_spec);
  options.DiagPivotThresh = diag_pivot_thresh;
  StatInit(&stat);
  
  get_perm_c(permc_spec, A, self->perm_c); /* calc column permutation */
  sp_preorder(&options, A, self->perm_c, etree, &AC); /* apply column permutation */
  

  /* Perform factorization */
  switch (A->Dtype) {
  case SLU_S:
      sgstrf(&options, &AC, (float) drop_tol, relax, panel_size,
             etree, NULL, lwork, self->perm_c, self->perm_r,
             &self->L, &self->U, &stat, &info);
      break;
  case SLU_D:
      dgstrf(&options, &AC, drop_tol, relax, panel_size,
             etree, NULL, lwork, self->perm_c, self->perm_r,
             &self->L, &self->U, &stat, &info);
      break;
  case SLU_C:          
      cgstrf(&options, &AC, (float) drop_tol, relax, panel_size,
             etree, NULL, lwork, self->perm_c, self->perm_r,
             &self->L, &self->U, &stat, &info);
      break;
  case SLU_Z:          
      zgstrf(&options, &AC, drop_tol, relax, panel_size,
             etree, NULL, lwork, self->perm_c, self->perm_r,
             &self->L, &self->U, &stat, &info);
      break;
  default:
      PyErr_SetString(PyExc_ValueError, "Invalid type in SuperMatrix.");
      goto fail;
  }
  
  if (info) {
    if (info < 0)
        PyErr_SetString(PyExc_SystemError, "dgstrf was called with invalid arguments");
    else {
        if (info <= n) 
            PyErr_SetString(PyExc_RuntimeError, "Factor is exactly singular");
        else
            PyErr_NoMemory();
    }
    goto fail;
  }
  
  /* free memory */
  SUPERLU_FREE(etree);
  Destroy_CompCol_Permuted(&AC);
  StatFree(&stat);
  
  return (PyObject *)self;

 fail:
  SUPERLU_FREE(etree);
  Destroy_CompCol_Permuted(&AC);
  StatFree(&stat);
  SciPyLU_dealloc(self);
  return NULL;
}