NEURON
ocmatrix.cpp
Go to the documentation of this file.
1 #include <../../nrnconf.h>
2 #include <vector>
3 #include <InterViews/resource.h>
4 
5 #define v_elem(v, i) (*(vector_vec(v) + i))
6 
7 #include "ivocvect.h"
8 #include "oc2iv.h"
9 #include <OS/math.h>
10 
11 #undef error
12 
13 extern "C" {
14 #undef OUT /* /usr/x86_64-w64-mingw32/sys-root/mingw/include/windef.h */
15 #include "matrix.h" //meschach
16 #include "matrix2.h"
17 #include "sparse.h"
18 #include "sparse2.h"
19 extern MAT* m_get(int, int);
20 } // extern "C"
21 
22 int nrn_matrix_dim(void*, int);
23 
24 #include "ocmatrix.h"
25 using std::vector;
26 
27 int nrn_matrix_dim(void* vm, int d) {
28  OcMatrix* m = (OcMatrix*) vm;
29  return d ? m->ncol() : m->nrow();
30 }
31 
32 static void Vect2VEC(Vect* v1, VEC& v2) {
33 #ifdef WIN32
34  v2.ve = vector_vec(v1);
35  v2.dim = vector_capacity(v1);
36  v2.max_dim = vector_buffer_size(v1);
37 #else
38  v2.ve = v1->data();
39  v2.dim = v1->size();
40  v2.max_dim = v1->buffer_size();
41 #endif
42 }
43 
45  obj_ = NULL;
46  type_ = type;
47 }
49 
50 OcMatrix* OcMatrix::instance(int nrow, int ncol, int type) {
51  switch (type) {
52  default:
53  case MFULL:
54  return new OcFullMatrix(nrow, ncol);
55  case MSPARSE:
56  return new OcSparseMatrix(nrow, ncol);
57  }
58 }
59 
61  hoc_execerror("Matrix method not implemented for this type matrix", 0);
62 }
63 
64 void OcMatrix::nonzeros(vector<int>& m, vector<int>& n) {
65  m.clear();
66  n.clear();
67  for (int i = 0; i < nrow(); i++) {
68  for (int j = 0; j < ncol(); j++) {
69  if (getval(i, j)) {
70  m.push_back(i);
71  n.push_back(j);
72  }
73  }
74  }
75 }
76 
77 void OcSparseMatrix::nonzeros(vector<int>& m, vector<int>& n) {
78  m.clear();
79  n.clear();
80  for (int i = 0; i < m_->m; i++) {
81  SPROW* const r = m_->row + i;
82  row_elt* r_elt = r->elt;
83  for (int k = 0; k < r->len; k++) {
84  int j = r_elt[k].col;
85  m.push_back(i);
86  n.push_back(j);
87  }
88  }
89 }
90 
91 
93  if (type_ != MFULL) { // could clone one maybe
94  hoc_execerror("Matrix is not a FULL matrix (type 1)", 0);
95  }
96  return (OcFullMatrix*) this;
97 }
98 
99 OcFullMatrix::OcFullMatrix(int nrow, int ncol)
100  : OcMatrix(MFULL) {
101  lu_factor_ = NULL;
102  lu_pivot_ = NULL;
103  m_ = m_get(nrow, ncol);
104 }
106  if (lu_factor_) {
109  }
110  M_FREE(m_);
111 }
112 double* OcFullMatrix::mep(int i, int j) {
113  return &m_->me[i][j];
114 }
115 double OcFullMatrix::getval(int i, int j) {
116  return m_->me[i][j];
117 }
119  return m_->m;
120 }
122  return m_->n;
123 }
124 
125 void OcFullMatrix::resize(int i, int j) {
126  m_resize(m_, i, j);
127 }
128 
129 void OcFullMatrix::mulv(Vect* vin, Vect* vout) {
130  VEC v1, v2;
131  Vect2VEC(vin, v1);
132  Vect2VEC(vout, v2);
133  mv_mlt(m_, &v1, &v2);
134 }
135 
137  m_mlt(m_, in->full()->m_, out->full()->m_);
138 }
139 
140 void OcFullMatrix::muls(double s, Matrix* out) {
141  sm_mlt(s, m_, out->full()->m_);
142 }
143 
144 void OcFullMatrix::add(Matrix* in, Matrix* out) {
145  m_add(m_, in->full()->m_, out->full()->m_);
146 }
147 
149  m_copy(m_, out->full()->m_);
150 }
151 
152 void OcFullMatrix::bcopy(Matrix* out, int i0, int j0, int n0, int m0, int i1, int j1) {
153  m_move(m_, i0, j0, n0, m0, out->full()->m_, i1, j1);
154 }
155 
157  m_transp(m_, out->full()->m_);
158 }
159 
160 void OcFullMatrix::symmeigen(Matrix* mout, Vect* vout) {
161  VEC v1;
162  Vect2VEC(vout, v1);
163  symmeig(m_, mout->full()->m_, &v1);
164 }
165 
167  VEC v1;
168  Vect2VEC(d, v1);
169  svd(m_, u ? u->full()->m_ : NULL, v ? v->full()->m_ : NULL, &v1);
170 }
171 
172 void OcFullMatrix::getrow(int k, Vect* out) {
173  VEC v1;
174  Vect2VEC(out, v1);
175  get_row(m_, k, &v1);
176 }
177 
178 void OcFullMatrix::getcol(int k, Vect* out) {
179  VEC v1;
180  Vect2VEC(out, v1);
181  get_col(m_, k, &v1);
182 }
183 
184 void OcFullMatrix::getdiag(int k, Vect* out) {
185  int i, j, row, col;
186  row = nrow();
187  col = ncol();
188  if (k >= 0) {
189  for (i = 0, j = k; i < row && j < col; ++i, ++j) {
190 #ifdef WIN32
191  v_elem(out, i) = m_entry(m_, i, j);
192 #else
193  out->elem(i) = m_entry(m_, i, j);
194 #endif
195  }
196  } else {
197  for (i = -k, j = 0; i < row && j < col; ++i, ++j) {
198 #ifdef WIN32
199  v_elem(out, i) = m_entry(m_, i, j);
200 #else
201  out->elem(i) = m_entry(m_, i, j);
202 #endif
203  }
204  }
205 }
206 
207 void OcFullMatrix::setrow(int k, Vect* in) {
208  VEC v1;
209  Vect2VEC(in, v1);
210  set_row(m_, k, &v1);
211 }
212 
213 void OcFullMatrix::setcol(int k, Vect* in) {
214  VEC v1;
215  Vect2VEC(in, v1);
216  set_col(m_, k, &v1);
217 }
218 
219 void OcFullMatrix::setdiag(int k, Vect* in) {
220  int i, j, row, col;
221  row = nrow();
222  col = ncol();
223  if (k >= 0) {
224  for (i = 0, j = k; i < row && j < col; ++i, ++j) {
225 #ifdef WIN32
226  m_set_val(m_, i, j, v_elem(in, i));
227 #else
228  m_set_val(m_, i, j, in->elem(i));
229 #endif
230  }
231  } else {
232  for (i = -k, j = 0; i < row && j < col; ++i, ++j) {
233 #ifdef WIN32
234  m_set_val(m_, i, j, v_elem(in, i));
235 #else
236  m_set_val(m_, i, j, in->elem(i));
237 #endif
238  }
239  }
240 }
241 
242 void OcFullMatrix::setrow(int k, double in) {
243  int i, col = ncol();
244  for (i = 0; i < col; ++i) {
245  m_set_val(m_, k, i, in);
246  }
247 }
248 
249 void OcFullMatrix::setcol(int k, double in) {
250  int i, row = nrow();
251  for (i = 0; i < row; ++i) {
252  m_set_val(m_, i, k, in);
253  }
254 }
255 
256 void OcFullMatrix::setdiag(int k, double in) {
257  int i, j, row, col;
258  row = nrow();
259  col = ncol();
260  if (k >= 0) {
261  for (i = 0, j = k; i < row && j < col; ++i, ++j) {
262  m_set_val(m_, i, j, in);
263  }
264  } else {
265  for (i = -k, j = 0; i < row && j < col; ++i, ++j) {
266  m_set_val(m_, i, j, in);
267  }
268  }
269 }
270 
272  m_zero(m_);
273 }
274 
276  m_ident(m_);
277 }
278 
280  m_exp(m_, 0., out->full()->m_);
281 }
282 
283 void OcFullMatrix::pow(int i, Matrix* out) {
284  m_pow(m_, i, out->full()->m_);
285 }
286 
288  m_inverse(m_, out->full()->m_);
289 }
290 
291 void OcFullMatrix::solv(Vect* in, Vect* out, bool use_lu) {
292  bool call_lufac = true;
293  if (!lu_factor_) {
294  lu_factor_ = m_get(nrow(), nrow());
295  lu_pivot_ = px_get(nrow());
296  } else if (use_lu && lu_factor_->m == nrow()) {
297  call_lufac = false;
298  }
299  VEC v1, v2;
300  Vect2VEC(in, v1);
301  Vect2VEC(out, v2);
302  if (call_lufac) {
303  m_resize(lu_factor_, nrow(), nrow());
304  m_copy(m_, lu_factor_);
307  }
308  LUsolve(lu_factor_, lu_pivot_, &v1, &v2);
309 }
310 
311 double OcFullMatrix::det(int* e) {
312  int n = nrow();
313  MAT* lu = m_get(n, n);
314  PERM* piv = px_get(n);
315  m_copy(m_, lu);
316  LUfactor(lu, piv);
317 #if 0
318 printf("LU\n");
319 for (int i = 0; i < n; ++i) {
320  for (int j = 0; j < n; ++j) {
321  printf(" %g", lu->me[i][j]);
322  }
323  printf("\t%d\n", piv->pe[i]);
324 }
325 #endif
326  double m = 1.0;
327  *e = 0;
328  for (int i = 0; i < n; ++i) {
329  m *= lu->me[i][i];
330  if (m == 0.0) {
331  break;
332  }
333  while (Math::abs(m) >= 1e12) {
334  m *= 1e-12;
335  *e += 12;
336  }
337  while (Math::abs(m) < 1e-12) {
338  m *= 1e12;
339  *e -= 12;
340  }
341  }
342  if (m) {
343  while (Math::abs(m) >= 10.0) {
344  m *= 0.1;
345  *e += 1;
346  }
347  while (Math::abs(m) < 1.0) {
348  m *= 10.0;
349  *e -= 1;
350  }
351  }
352  m *= double(px_sign(piv));
353  M_FREE(lu);
354  PX_FREE(piv);
355  return m;
356 }
357 
358 //--------------------------
359 
361  : OcMatrix(MSPARSE) {
362  /* sp_get -- get sparse matrix
363  -- len is number of elements available for each row without
364  allocating further memory */
365 
366  int len = 4;
367  m_ = sp_get(nrow, ncol, len);
368  lu_factor_ = NULL;
369  lu_pivot_ = NULL;
370 }
372  if (lu_factor_) {
375  }
376  SP_FREE(m_);
377 }
378 
379 // returns pointer to sparse element. NULL if it does not exist.
380 double* OcSparseMatrix::pelm(int i, int j) {
381  SPROW* r = m_->row + i;
382  int idx = sprow_idx(r, j);
383  if (idx >= 0) {
384  return &r->elt[idx].val;
385  } else {
386  return NULL;
387  }
388 }
389 
390 double* OcSparseMatrix::mep(int i, int j) {
391  SPROW* r = m_->row + i;
392  int idx = sprow_idx(r, j);
393  if (idx >= 0) {
394  return &r->elt[idx].val;
395  }
396  // does not exist so create it with a value of 0
397  sp_set_val(m_, i, j, 0.);
398  // and try again
399  idx = sprow_idx(r, j);
400  return &r->elt[idx].val;
401 }
402 
404  sp_zero(m_);
405 }
406 
407 double OcSparseMatrix::getval(int i, int j) {
408  return sp_get_val(m_, i, j);
409 }
411  return m_->m;
412 }
414  return m_->n;
415 }
416 void OcSparseMatrix::mulv(Vect* vin, Vect* vout) {
417  VEC v1, v2;
418  Vect2VEC(vin, v1);
419  Vect2VEC(vout, v2);
420  sp_mv_mlt(m_, &v1, &v2);
421 }
422 
423 void OcSparseMatrix::solv(Vect* in, Vect* out, bool use_lu) {
424  bool call_lufac = true;
425  if (!lu_factor_) {
426  lu_factor_ = sp_get(nrow(), nrow(), 4);
427  lu_pivot_ = px_get(nrow());
428  } else if (use_lu && lu_factor_->m == nrow()) {
429  call_lufac = false;
430  }
431  VEC v1, v2;
432  Vect2VEC(in, v1);
433  Vect2VEC(out, v2);
434  if (call_lufac) {
435  sp_resize(lu_factor_, nrow(), nrow());
439  }
440  spLUsolve(lu_factor_, lu_pivot_, &v1, &v2);
441 }
442 
443 void OcSparseMatrix::setrow(int k, Vect* in) {
444  VEC v1;
445  Vect2VEC(in, v1);
446  int i, n = ncol();
447  double* p;
448  for (i = 0; i < n; ++i) {
449  if ((p = pelm(k, i)) != NULL) {
450 #ifdef WIN32
451  *p = v_elem(in, i);
452  } else if (v_elem(in, i)) {
453  sp_set_val(m_, k, i, v_elem(in, i));
454 #else
455  *p = in->elem(i);
456  } else if (in->elem(i)) {
457  sp_set_val(m_, k, i, in->elem(i));
458 #endif
459  }
460  }
461 }
462 
463 void OcSparseMatrix::setcol(int k, Vect* in) {
464  VEC v1;
465  Vect2VEC(in, v1);
466  int i, n = nrow();
467  double* p;
468  for (i = 0; i < n; ++i) {
469  if ((p = pelm(i, k)) != NULL) {
470 #ifdef WIN32
471  *p = v_elem(in, i);
472  } else if (v_elem(in, i)) {
473  sp_set_val(m_, i, k, v_elem(in, i));
474 #else
475  *p = in->elem(i);
476  } else if (in->elem(i)) {
477  sp_set_val(m_, i, k, in->elem(i));
478 #endif
479  }
480  }
481 }
482 
484  int i, j, row, col;
485  row = nrow();
486  col = ncol();
487  double* p;
488  if (k >= 0) {
489  for (i = 0, j = k; i < row && j < col; ++i, ++j) {
490  if ((p = pelm(i, j)) != NULL) {
491 #ifdef WIN32
492  *p = v_elem(in, i);
493  } else if (v_elem(in, i)) {
494  sp_set_val(m_, i, j, v_elem(in, i));
495 #else
496  *p = in->elem(i);
497  } else if (in->elem(i)) {
498  sp_set_val(m_, i, j, in->elem(i));
499 #endif
500  }
501  }
502  } else {
503  for (i = -k, j = 0; i < row && j < col; ++i, ++j) {
504  if ((p = pelm(i, j)) != NULL) {
505 #ifdef WIN32
506  *p = v_elem(in, i);
507  } else if (v_elem(in, i)) {
508  sp_set_val(m_, i, j, v_elem(in, i));
509 #else
510  *p = in->elem(i);
511  } else if (in->elem(i)) {
512  sp_set_val(m_, i, j, in->elem(i));
513 #endif
514  }
515  }
516  }
517 }
518 
519 void OcSparseMatrix::setrow(int k, double in) {
520  int i, col = ncol();
521  for (i = 0; i < col; ++i) {
522  sp_set_val(m_, k, i, in);
523  }
524 }
525 
526 void OcSparseMatrix::setcol(int k, double in) {
527  int i, row = nrow();
528  for (i = 0; i < row; ++i) {
529  sp_set_val(m_, i, k, in);
530  }
531 }
532 
534  setdiag(0, 1);
535 }
536 
537 void OcSparseMatrix::setdiag(int k, double in) {
538  int i, j, row, col;
539  row = nrow();
540  col = ncol();
541  if (k >= 0) {
542  for (i = 0, j = k; i < row && j < col; ++i, ++j) {
543  sp_set_val(m_, i, j, in);
544  }
545  } else {
546  for (i = -k, j = 0; i < row && j < col; ++i, ++j) {
547  sp_set_val(m_, i, j, in);
548  }
549  }
550 }
551 
553  return m_->row[i].len;
554 }
555 
556 double OcSparseMatrix::spgetrowval(int i, int jindx, int* j) {
557  *j = m_->row[i].elt[jindx].col;
558  return m_->row[i].elt[jindx].val;
559 }
short type
Definition: cabvars.h:9
static int abs(int)
Definition: math.cpp:43
virtual void ident()
Definition: ocmatrix.cpp:275
PERM * lu_pivot_
Definition: ocmatrix.h:196
virtual int ncol()
Definition: ocmatrix.cpp:121
virtual void getdiag(int, Vect *out)
Definition: ocmatrix.cpp:184
virtual void exp(Matrix *out)
Definition: ocmatrix.cpp:279
virtual void pow(int, Matrix *out)
Definition: ocmatrix.cpp:283
virtual void setdiag(int, Vect *in)
Definition: ocmatrix.cpp:219
virtual void getcol(int, Vect *out)
Definition: ocmatrix.cpp:178
MAT * lu_factor_
Definition: ocmatrix.h:195
virtual void setrow(int, Vect *in)
Definition: ocmatrix.cpp:207
virtual void bcopy(Matrix *mout, int i0, int j0, int n0, int m0, int i1, int j1)
Definition: ocmatrix.cpp:152
virtual void setcol(int, Vect *in)
Definition: ocmatrix.cpp:213
MAT * m_
Definition: ocmatrix.h:194
virtual double det(int *exponent)
Definition: ocmatrix.cpp:311
virtual double * mep(int, int)
Definition: ocmatrix.cpp:112
virtual int nrow()
Definition: ocmatrix.cpp:118
virtual void transpose(Matrix *out)
Definition: ocmatrix.cpp:156
virtual double getval(int i, int j)
Definition: ocmatrix.cpp:115
virtual void solv(Vect *vin, Vect *vout, bool use_lu)
Definition: ocmatrix.cpp:291
virtual void getrow(int, Vect *out)
Definition: ocmatrix.cpp:172
virtual void copy(Matrix *out)
Definition: ocmatrix.cpp:148
virtual void zero()
Definition: ocmatrix.cpp:271
virtual void inverse(Matrix *out)
Definition: ocmatrix.cpp:287
virtual void svd1(Matrix *u, Matrix *v, Vect *d)
Definition: ocmatrix.cpp:166
virtual void add(Matrix *, Matrix *out)
Definition: ocmatrix.cpp:144
virtual void symmeigen(Matrix *mout, Vect *vout)
Definition: ocmatrix.cpp:160
virtual void muls(double, Matrix *out)
Definition: ocmatrix.cpp:140
virtual void resize(int, int)
Definition: ocmatrix.cpp:125
virtual void mulv(Vect *in, Vect *out)
Definition: ocmatrix.cpp:129
virtual void mulm(Matrix *in, Matrix *out)
Definition: ocmatrix.cpp:136
virtual ~OcFullMatrix()
Definition: ocmatrix.cpp:105
OcFullMatrix(int, int)
Definition: ocmatrix.cpp:99
Object * obj_
Definition: ocmatrix.h:148
virtual double getval(int i, int j)
Definition: ocmatrix.h:33
virtual void nonzeros(vector< int > &m, vector< int > &n)
Definition: ocmatrix.cpp:64
virtual int ncol()
Definition: ocmatrix.h:41
@ MFULL
Definition: ocmatrix.h:21
@ MSPARSE
Definition: ocmatrix.h:21
int type_
Definition: ocmatrix.h:151
static OcMatrix * instance(int nrow, int ncol, int type=MFULL)
Definition: ocmatrix.cpp:50
OcFullMatrix * full()
Definition: ocmatrix.cpp:92
OcMatrix(int type)
Definition: ocmatrix.cpp:44
void unimp()
Definition: ocmatrix.cpp:60
virtual ~OcMatrix()
Definition: ocmatrix.cpp:48
virtual int nrow()
Definition: ocmatrix.h:37
virtual void setrow(int, Vect *in)
Definition: ocmatrix.cpp:443
virtual double * mep(int, int)
Definition: ocmatrix.cpp:390
virtual double getval(int, int)
Definition: ocmatrix.cpp:407
virtual void zero()
Definition: ocmatrix.cpp:403
virtual void nonzeros(vector< int > &m, vector< int > &n)
Definition: ocmatrix.cpp:77
virtual double spgetrowval(int i, int jindx, int *j)
Definition: ocmatrix.cpp:556
virtual void mulv(Vect *in, Vect *out)
Definition: ocmatrix.cpp:416
virtual int sprowlen(int)
Definition: ocmatrix.cpp:552
SPMAT * lu_factor_
Definition: ocmatrix.h:229
virtual void solv(Vect *vin, Vect *vout, bool use_lu)
Definition: ocmatrix.cpp:423
virtual int ncol()
Definition: ocmatrix.cpp:413
OcSparseMatrix(int, int)
Definition: ocmatrix.cpp:360
SPMAT * m_
Definition: ocmatrix.h:228
virtual int nrow()
Definition: ocmatrix.cpp:410
virtual void setcol(int, Vect *in)
Definition: ocmatrix.cpp:463
virtual void setdiag(int, Vect *in)
Definition: ocmatrix.cpp:483
virtual double * pelm(int, int)
Definition: ocmatrix.cpp:380
virtual void ident(void)
Definition: ocmatrix.cpp:533
virtual ~OcSparseMatrix()
Definition: ocmatrix.cpp:371
PERM * lu_pivot_
Definition: ocmatrix.h:230
void hoc_execerror(const char *, const char *)
Definition: hoc.cpp:754
int vector_capacity(Vect *v)
Definition: ivocvect.cpp:319
int vector_buffer_size(Vect *v)
Definition: ivocvect.cpp:316
double * vector_vec(Vect *v)
Definition: ivocvect.cpp:328
#define Vect
Definition: ivocvect.h:14
MAT * LUfactor(MAT *A, PERM *pivot)
Definition: lufactor.c:46
VEC * LUsolve(MAT *A, PERM *pivot, VEC *b, VEC *x)
Definition: lufactor.c:140
MAT * sm_mlt(double scalar, MAT *matrix, MAT *out)
Definition: matop.c:241
VEC * mv_mlt(MAT *A, VEC *b, VEC *out)
Definition: matop.c:205
MAT * m_transp(MAT *in, MAT *out)
Definition: matop.c:300
MAT * m_mlt(MAT *A, MAT *B, MAT *OUT)
Definition: matop.c:88
VEC * symmeig(MAT *A, MAT *Q, VEC *out)
Definition: symmeig.c:174
VEC * svd(MAT *A, MAT *U, MAT *V, VEC *out)
Definition: svd.c:350
static Object ** m_pow(void *v)
Definition: matrix.cpp:531
static Object ** m_zero(void *v)
Definition: matrix.cpp:512
static Object ** m_exp(void *v)
Definition: matrix.cpp:524
static Object ** m_ident(void *v)
Definition: matrix.cpp:518
static Object ** m_resize(void *v)
Definition: matrix.cpp:191
static Object ** m_add(void *v)
Definition: matrix.cpp:262
static Object ** m_inverse(void *v)
Definition: matrix.cpp:539
#define m_set_val(A, i, j, val)
Definition: matrix.h:277
#define M_FREE(mat)
Definition: matrix.h:213
#define m_copy(in, out)
Definition: matrix.h:403
VEC * get_row(MAT *, u_int, VEC *)
PERM * px_get(int)
#define set_col(mat, col, vec)
Definition: matrix.h:564
PERM * px_resize(PERM *, int)
Definition: memory.c:399
#define m_entry(A, i, j)
Definition: matrix.h:274
#define PX_FREE(px)
Definition: matrix.h:215
#define m_move
Definition: matrix.h:50
#define set_row(mat, row, vec)
Definition: matrix.h:562
int px_sign(PERM *)
Definition: pxop.c:269
VEC * get_col(MAT *, u_int, VEC *)
#define v
Definition: md1redef.h:4
#define i
Definition: md1redef.h:12
#define printf
Definition: mwprefix.h:26
static unsigned row
Definition: nonlin.cpp:85
int const size_t const size_t n
Definition: nrngsl.h:11
size_t p
size_t j
static philox4x32_key_t k
Definition: nrnran123.cpp:11
static void Vect2VEC(Vect *v1, VEC &v2)
Definition: ocmatrix.cpp:32
#define v_elem(v, i)
Definition: ocmatrix.cpp:5
int nrn_matrix_dim(void *, int)
Definition: ocmatrix.cpp:27
MAT * m_get(int, int)
Definition: memory.c:36
#define e
Definition: passive0.cpp:22
VEC * spLUsolve(SPMAT *, PERM *, VEC *, VEC *)
SPMAT * spLUfactor(SPMAT *, PERM *, double)
Definition: splufctr.c:49
SPMAT * sp_zero(SPMAT *A)
Definition: sparse.c:633
SPMAT * sp_resize(SPMAT *A, int m, int n)
Definition: sparse.c:727
double sp_set_val(SPMAT *A, int i, int j, double val)
Definition: sparse.c:66
double sp_get_val(SPMAT *A, int i, int j)
Definition: sparse.c:45
SPMAT * sp_get(int m, int n, int maxlen)
Definition: sparse.c:198
VEC * sp_mv_mlt(SPMAT *A, VEC *x, VEC *out)
Definition: sparse.c:123
SPMAT * sp_copy2(SPMAT *A, SPMAT *OUT)
Definition: sparse.c:655
#define SP_FREE(A)
Definition: sparse.h:183
int sprow_idx(SPROW *, int)
Definition: sprow.c:75
#define NULL
Definition: sptree.h:16
MatrixPtr Matrix
Definition: sputils.c:601
Definition: matrix.h:73
u_int n
Definition: matrix.h:74
u_int m
Definition: matrix.h:74
Real ** me
Definition: matrix.h:76
Definition: matrix.h:87
u_int * pe
Definition: matrix.h:88
SPROW * row
Definition: sparse.h:57
int n
Definition: sparse.h:55
int m
Definition: sparse.h:55
Definition: sparse.h:49
int len
Definition: sparse.h:50
row_elt * elt
Definition: sparse.h:51
Definition: matrix.h:67
u_int dim
Definition: matrix.h:68
u_int max_dim
Definition: matrix.h:68
Real * ve
Definition: matrix.h:69
Definition: sparse.h:44
int col
Definition: sparse.h:45
Real val
Definition: sparse.h:46