xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 5e42470ac62cb63ed26e96b335b06ee82e25fc39)
1 #ifndef lint
2 static char vcid[] = "$Id: dense.c,v 1.19 1995/03/25 01:11:56 bsmith Exp curfman $";
3 #endif
4 
5 /*
6     Standard Fortran style matrices
7 */
8 #include "ptscimpl.h"
9 #include "plapack.h"
10 #include "matimpl.h"
11 #include "math.h"
12 #include "vec/vecimpl.h"
13 
14 typedef struct {
15   Scalar *v;
16   int    roworiented;
17   int    m,n,pad;
18   int    *pivots;   /* pivots in LU factorization */
19 } Mat_Dense;
20 
21 
22 static int MatGetInfo_Dense(Mat matin,int flag,int *nz,int *nzalloc,int *mem)
23 {
24   Mat_Dense *mat = (Mat_Dense *) matin->data;
25   int    i,N = mat->m*mat->n,count = 0;
26   Scalar *v = mat->v;
27   for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;}
28   *nz = count; *nzalloc = N; *mem = N*sizeof(Scalar);
29   return 0;
30 }
31 
32 /* ---------------------------------------------------------------*/
33 /* COMMENT: I have chosen to hide column permutation in the pivots,
34    rather than put it in the Mat->col slot.*/
35 static int MatLUFactor_Dense(Mat matin,IS row,IS col)
36 {
37   Mat_Dense *mat = (Mat_Dense *) matin->data;
38   int    info;
39   if (!mat->pivots) {
40     mat->pivots = (int *) MALLOC( mat->m*sizeof(int) );
41     CHKPTR(mat->pivots);
42   }
43   LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info);
44   if (info) SETERR(1,"Bad LU factorization");
45   matin->factor = FACTOR_LU;
46   return 0;
47 }
48 static int MatLUFactorSymbolic_Dense(Mat matin,IS row,IS col,Mat *fact)
49 {
50   int ierr;
51   if ((ierr = MatCopy(matin,fact))) SETERR(ierr,0);
52   return 0;
53 }
54 static int MatLUFactorNumeric_Dense(Mat matin,Mat *fact)
55 {
56   return MatLUFactor(*fact,0,0);
57 }
58 static int MatChFactorSymbolic_Dense(Mat matin,IS row,Mat *fact)
59 {
60   int ierr;
61   if ((ierr = MatCopy(matin,fact))) SETERR(ierr,0);
62   return 0;
63 }
64 static int MatChFactorNumeric_Dense(Mat matin,Mat *fact)
65 {
66   return MatCholeskyFactor(*fact,0);
67 }
68 static int MatChFactor_Dense(Mat matin,IS perm)
69 {
70   Mat_Dense    *mat = (Mat_Dense *) matin->data;
71   int       info;
72   if (mat->pivots) {FREE(mat->pivots); mat->pivots = 0;}
73   LApotrf_("L",&mat->n,mat->v,&mat->m,&info);
74   if (info) SETERR(1,"Bad Cholesky factorization");
75   matin->factor = FACTOR_CHOLESKY;
76   return 0;
77 }
78 
79 static int MatSolve_Dense(Mat matin,Vec xx,Vec yy)
80 {
81   Mat_Dense *mat = (Mat_Dense *) matin->data;
82   int    one = 1, info;
83   Scalar *x, *y;
84   VecGetArray(xx,&x); VecGetArray(yy,&y);
85   MEMCPY(y,x,mat->m*sizeof(Scalar));
86   if (matin->factor == FACTOR_LU) {
87     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,
88               y, &mat->m, &info );
89   }
90   else if (matin->factor == FACTOR_CHOLESKY){
91     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
92               y, &mat->m, &info );
93   }
94   else SETERR(1,"Matrix must be factored to solve");
95   if (info) SETERR(1,"Bad solve");
96   return 0;
97 }
98 static int MatSolveTrans_Dense(Mat matin,Vec xx,Vec yy)
99 {
100   Mat_Dense *mat = (Mat_Dense *) matin->data;
101   int    one = 1, info;
102   Scalar *x, *y;
103   VecGetArray(xx,&x); VecGetArray(yy,&y);
104   MEMCPY(y,x,mat->m*sizeof(Scalar));
105   /* assume if pivots exist then LU else Cholesky */
106   if (mat->pivots) {
107     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,
108               y, &mat->m, &info );
109   }
110   else {
111     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
112               y, &mat->m, &info );
113   }
114   if (info) SETERR(1,"Bad solve");
115   return 0;
116 }
117 static int MatSolveAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy)
118 {
119   Mat_Dense *mat = (Mat_Dense *) matin->data;
120   int    one = 1, info,ierr;
121   Scalar *x, *y, sone = 1.0;
122   Vec    tmp = 0;
123   VecGetArray(xx,&x); VecGetArray(yy,&y);
124   if (yy == zz) {
125     ierr = VecCreate(yy,&tmp); CHKERR(ierr);
126     ierr = VecCopy(yy,tmp); CHKERR(ierr);
127   }
128   MEMCPY(y,x,mat->m*sizeof(Scalar));
129   /* assume if pivots exist then LU else Cholesky */
130   if (mat->pivots) {
131     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,
132               y, &mat->m, &info );
133   }
134   else {
135     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
136               y, &mat->m, &info );
137   }
138   if (info) SETERR(1,"Bad solve");
139   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
140   else VecAXPY(&sone,zz,yy);
141   return 0;
142 }
143 static int MatSolveTransAdd_Dense(Mat matin,Vec xx,Vec zz, Vec yy)
144 {
145   Mat_Dense  *mat = (Mat_Dense *) matin->data;
146   int     one = 1, info,ierr;
147   Scalar  *x, *y, sone = 1.0;
148   Vec     tmp;
149   VecGetArray(xx,&x); VecGetArray(yy,&y);
150   if (yy == zz) {
151     ierr = VecCreate(yy,&tmp); CHKERR(ierr);
152     ierr = VecCopy(yy,tmp); CHKERR(ierr);
153   }
154   MEMCPY(y,x,mat->m*sizeof(Scalar));
155   /* assume if pivots exist then LU else Cholesky */
156   if (mat->pivots) {
157     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,
158               y, &mat->m, &info );
159   }
160   else {
161     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
162               y, &mat->m, &info );
163   }
164   if (info) SETERR(1,"Bad solve");
165   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
166   else VecAXPY(&sone,zz,yy);
167   return 0;
168 }
169 /* ------------------------------------------------------------------*/
170 static int MatRelax_Dense(Mat matin,Vec bb,double omega,int flag,double shift,
171                        int its,Vec xx)
172 {
173   Mat_Dense *mat = (Mat_Dense *) matin->data;
174   Scalar *x, *b, *v = mat->v, zero = 0.0, xt;
175   int    o = 1,ierr, m = mat->m, i;
176 
177   if (flag & SOR_ZERO_INITIAL_GUESS) {
178     /* this is a hack fix, should have another version without
179        the second BLdot */
180     if ((ierr = VecSet(&zero,xx))) SETERR(ierr,0);
181   }
182   VecGetArray(xx,&x); VecGetArray(bb,&b);
183   while (its--) {
184     if (flag & SOR_FORWARD_SWEEP){
185       for ( i=0; i<m; i++ ) {
186         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
187         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
188       }
189     }
190     if (flag & SOR_BACKWARD_SWEEP) {
191       for ( i=m-1; i>=0; i-- ) {
192         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
193         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
194       }
195     }
196   }
197   return 0;
198 }
199 
200 /* -----------------------------------------------------------------*/
201 static int MatMultTrans_Dense(Mat matin,Vec xx,Vec yy)
202 {
203   Mat_Dense *mat = (Mat_Dense *) matin->data;
204   Scalar *v = mat->v, *x, *y;
205   int _One=1;Scalar _DOne=1.0, _DZero=0.0;
206   VecGetArray(xx,&x), VecGetArray(yy,&y);
207   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
208          x, &_One, &_DZero, y, &_One );
209   return 0;
210 }
211 static int MatMult_Dense(Mat matin,Vec xx,Vec yy)
212 {
213   Mat_Dense *mat = (Mat_Dense *) matin->data;
214   Scalar *v = mat->v, *x, *y;
215   int _One=1;Scalar _DOne=1.0, _DZero=0.0;
216   VecGetArray(xx,&x); VecGetArray(yy,&y);
217   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
218          x, &_One, &_DZero, y, &_One );
219   return 0;
220 }
221 static int MatMultAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy)
222 {
223   Mat_Dense *mat = (Mat_Dense *) matin->data;
224   Scalar *v = mat->v, *x, *y, *z;
225   int    _One=1; Scalar _DOne=1.0;
226   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
227   if (zz != yy) MEMCPY(y,z,mat->m*sizeof(Scalar));
228   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
229          x, &_One, &_DOne, y, &_One );
230   return 0;
231 }
232 static int MatMultTransAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy)
233 {
234   Mat_Dense *mat = (Mat_Dense *) matin->data;
235   Scalar *v = mat->v, *x, *y, *z;
236   int    _One=1;
237   Scalar _DOne=1.0;
238   VecGetArray(xx,&x); VecGetArray(yy,&y);
239   VecGetArray(zz,&z);
240   if (zz != yy) MEMCPY(y,z,mat->m*sizeof(Scalar));
241   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
242          x, &_One, &_DOne, y, &_One );
243   return 0;
244 }
245 
246 /* -----------------------------------------------------------------*/
247 static int MatGetRow_Dense(Mat matin,int row,int *ncols,int **cols,
248                         Scalar **vals)
249 {
250   Mat_Dense *mat = (Mat_Dense *) matin->data;
251   Scalar *v;
252   int    i;
253   *ncols = mat->n;
254   if (cols) {
255     *cols = (int *) MALLOC(mat->n*sizeof(int)); CHKPTR(*cols);
256     for ( i=0; i<mat->n; i++ ) *cols[i] = i;
257   }
258   if (vals) {
259     *vals = (Scalar *) MALLOC(mat->n*sizeof(Scalar)); CHKPTR(*vals);
260     v = mat->v + row;
261     for ( i=0; i<mat->n; i++ ) {*vals[i] = *v; v += mat->m;}
262   }
263   return 0;
264 }
265 static int MatRestoreRow_Dense(Mat matin,int row,int *ncols,int **cols,
266                             Scalar **vals)
267 {
268   if (cols) { FREE(*cols); }
269   if (vals) { FREE(*vals); }
270   return 0;
271 }
272 /* ----------------------------------------------------------------*/
273 static int MatInsert_Dense(Mat matin,int m,int *indexm,int n,
274                         int *indexn,Scalar *v,InsertMode addv)
275 {
276   Mat_Dense *mat = (Mat_Dense *) matin->data;
277   int    i,j;
278 
279   if (!mat->roworiented) {
280     if (addv == InsertValues) {
281       for ( j=0; j<n; j++ ) {
282         if (indexn[j] < 0) {*v += m; continue;}
283         for ( i=0; i<m; i++ ) {
284           if (indexm[i] < 0) {*v++; continue;}
285           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
286         }
287       }
288     }
289     else {
290       for ( j=0; j<n; j++ ) {
291         if (indexn[j] < 0) {*v += m; continue;}
292         for ( i=0; i<m; i++ ) {
293           if (indexm[i] < 0) {*v++; continue;}
294           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
295         }
296       }
297     }
298   }
299   else {
300     if (addv == InsertValues) {
301       for ( i=0; i<m; i++ ) {
302         if (indexm[i] < 0) {*v += n; continue;}
303         for ( j=0; j<n; j++ ) {
304           if (indexn[j] < 0) {*v++; continue;}
305           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
306         }
307       }
308     }
309     else {
310       for ( i=0; i<m; i++ ) {
311         if (indexm[i] < 0) {*v += n; continue;}
312         for ( j=0; j<n; j++ ) {
313           if (indexn[j] < 0) {*v++; continue;}
314           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
315         }
316       }
317     }
318   }
319   return 0;
320 }
321 
322 /* -----------------------------------------------------------------*/
323 static int MatCopy_Dense(Mat matin,Mat *newmat)
324 {
325   Mat_Dense *mat = (Mat_Dense *) matin->data;
326   int ierr;
327   Mat newi;
328   Mat_Dense *l;
329   if ((ierr = MatCreateSequentialDense(mat->m,mat->n,&newi))) SETERR(ierr,0);
330   l = (Mat_Dense *) newi->data;
331   MEMCPY(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
332   *newmat = newi;
333   return 0;
334 }
335 #include "viewer.h"
336 
337 int MatView_Dense(PetscObject obj,Viewer ptr)
338 {
339   Mat         matin = (Mat) obj;
340   Mat_Dense      *mat = (Mat_Dense *) matin->data;
341   Scalar      *v;
342   int         i,j;
343   PetscObject ojb = (PetscObject) ptr;
344 
345   if (ojb && ojb->cookie == VIEWER_COOKIE && ojb->type == MATLAB_VIEWER) {
346     return ViewerMatlabPutArray(ptr,mat->m,mat->n,mat->v);
347   }
348   else {
349     FILE *fd = ViewerFileGetPointer(ptr);
350     for ( i=0; i<mat->m; i++ ) {
351       v = mat->v + i;
352       for ( j=0; j<mat->n; j++ ) {
353 #if defined(PETSC_COMPLEX)
354         fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m;
355 #else
356         fprintf(fd,"%6.4e ",*v); v += mat->m;
357 #endif
358       }
359       fprintf(fd,"\n");
360     }
361   }
362   return 0;
363 }
364 
365 
366 static int MatDestroy_Dense(PetscObject obj)
367 {
368   Mat    mat = (Mat) obj;
369   Mat_Dense *l = (Mat_Dense *) mat->data;
370 #if defined(PETSC_LOG)
371   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
372 #endif
373   if (l->pivots) FREE(l->pivots);
374   FREE(l);
375   PLogObjectDestroy(mat);
376   PETSCHEADERDESTROY(mat);
377   return 0;
378 }
379 
380 static int MatTrans_Dense(Mat matin)
381 {
382   Mat_Dense *mat = (Mat_Dense *) matin->data;
383   int    k,j;
384   Scalar *v = mat->v, tmp;
385   if (mat->m != mat->n) {
386     SETERR(1,"Cannot transpose rectangular dense matrix");
387   }
388   for ( j=0; j<mat->m; j++ ) {
389     for ( k=0; k<j; k++ ) {
390       tmp = v[j + k*mat->n];
391       v[j + k*mat->n] = v[k + j*mat->n];
392       v[k + j*mat->n] = tmp;
393     }
394   }
395   return 0;
396 }
397 
398 static int MatEqual_Dense(Mat matin1,Mat matin2)
399 {
400   Mat_Dense *mat1 = (Mat_Dense *) matin1->data;
401   Mat_Dense *mat2 = (Mat_Dense *) matin2->data;
402   int    i;
403   Scalar *v1 = mat1->v, *v2 = mat2->v;
404   if (mat1->m != mat2->m) return 0;
405   if (mat1->n != mat2->n) return 0;
406   for ( i=0; i<mat1->m*mat1->n; i++ ) {
407     if (*v1 != *v2) return 0;
408     v1++; v2++;
409   }
410   return 1;
411 }
412 
413 static int MatGetDiag_Dense(Mat matin,Vec v)
414 {
415   Mat_Dense *mat = (Mat_Dense *) matin->data;
416   int    i, n;
417   Scalar *x;
418   CHKTYPE(v,SEQVECTOR);
419   VecGetArray(v,&x); VecGetSize(v,&n);
420   if (n != mat->m) SETERR(1,"Nonconforming matrix and vector");
421   for ( i=0; i<mat->m; i++ ) {
422     x[i] = mat->v[i*mat->m + i];
423   }
424   return 0;
425 }
426 
427 static int MatScale_Dense(Mat matin,Vec ll,Vec rr)
428 {
429   Mat_Dense *mat = (Mat_Dense *) matin->data;
430   Scalar *l,*r,x,*v;
431   int    i,j,m = mat->m, n = mat->n;
432   if (ll) {
433     VecGetArray(ll,&l); VecGetSize(ll,&m);
434     if (m != mat->m) SETERR(1,"Left scaling vector wrong length");
435     for ( i=0; i<m; i++ ) {
436       x = l[i];
437       v = mat->v + i;
438       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
439     }
440   }
441   if (rr) {
442     VecGetArray(rr,&r); VecGetSize(rr,&n);
443     if (n != mat->n) SETERR(1,"Right scaling vector wrong length");
444     for ( i=0; i<n; i++ ) {
445       x = r[i];
446       v = mat->v + i*m;
447       for ( j=0; j<m; j++ ) { (*v++) *= x;}
448     }
449   }
450   return 0;
451 }
452 
453 
454 static int MatNorm_Dense(Mat matin,int type,double *norm)
455 {
456   Mat_Dense *mat = (Mat_Dense *) matin->data;
457   Scalar *v = mat->v;
458   double sum = 0.0;
459   int    i, j;
460   if (type == NORM_FROBENIUS) {
461     for (i=0; i<mat->n*mat->m; i++ ) {
462 #if defined(PETSC_COMPLEX)
463       sum += real(conj(*v)*(*v)); v++;
464 #else
465       sum += (*v)*(*v); v++;
466 #endif
467     }
468     *norm = sqrt(sum);
469   }
470   else if (type == NORM_1) {
471     *norm = 0.0;
472     for ( j=0; j<mat->n; j++ ) {
473       sum = 0.0;
474       for ( i=0; i<mat->m; i++ ) {
475 #if defined(PETSC_COMPLEX)
476         sum += abs(*v++);
477 #else
478         sum += fabs(*v++);
479 #endif
480       }
481       if (sum > *norm) *norm = sum;
482     }
483   }
484   else if (type == NORM_INFINITY) {
485     *norm = 0.0;
486     for ( j=0; j<mat->m; j++ ) {
487       v = mat->v + j;
488       sum = 0.0;
489       for ( i=0; i<mat->n; i++ ) {
490 #if defined(PETSC_COMPLEX)
491         sum += abs(*v); v += mat->m;
492 #else
493         sum += fabs(*v); v += mat->m;
494 #endif
495       }
496       if (sum > *norm) *norm = sum;
497     }
498   }
499   else {
500     SETERR(1,"No support for the two norm yet");
501   }
502   return 0;
503 }
504 
505 static int MatInsOpt_Dense(Mat aijin,int op)
506 {
507   Mat_Dense *aij = (Mat_Dense *) aijin->data;
508   if (op == ROW_ORIENTED)            aij->roworiented = 1;
509   else if (op == COLUMN_ORIENTED)    aij->roworiented = 0;
510   /* doesn't care about sorted rows or columns */
511   return 0;
512 }
513 
514 static int MatZero_Dense(Mat A)
515 {
516   Mat_Dense *l = (Mat_Dense *) A->data;
517   MEMSET(l->v,0,l->m*l->n*sizeof(Scalar));
518   return 0;
519 }
520 
521 static int MatZeroRows_Dense(Mat A,IS is,Scalar *diag)
522 {
523   Mat_Dense *l = (Mat_Dense *) A->data;
524   int    n = l->n, i, j,ierr,N, *rows;
525   Scalar *slot;
526   ierr = ISGetLocalSize(is,&N); CHKERR(ierr);
527   ierr = ISGetIndices(is,&rows); CHKERR(ierr);
528   for ( i=0; i<N; i++ ) {
529     slot = l->v + rows[i];
530     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
531   }
532   if (diag) {
533     for ( i=0; i<N; i++ ) {
534       slot = l->v + (n+1)*rows[i];
535       *slot = *diag;
536     }
537   }
538   ISRestoreIndices(is,&rows);
539   return 0;
540 }
541 
542 static int MatGetSize_Dense(Mat matin,int *m,int *n)
543 {
544   Mat_Dense *mat = (Mat_Dense *) matin->data;
545   *m = mat->m; *n = mat->n;
546   return 0;
547 }
548 
549 static int MatGetArray_Dense(Mat matin,Scalar **array)
550 {
551   Mat_Dense *mat = (Mat_Dense *) matin->data;
552   *array = mat->v;
553   return 0;
554 }
555 /* -------------------------------------------------------------------*/
556 static struct _MatOps MatOps = {MatInsert_Dense,
557        MatGetRow_Dense, MatRestoreRow_Dense,
558        MatMult_Dense, MatMultAdd_Dense,
559        MatMultTrans_Dense, MatMultTransAdd_Dense,
560        MatSolve_Dense,MatSolveAdd_Dense,
561        MatSolveTrans_Dense,MatSolveTransAdd_Dense,
562        MatLUFactor_Dense,MatChFactor_Dense,
563        MatRelax_Dense,
564        MatTrans_Dense,
565        MatGetInfo_Dense,MatEqual_Dense,
566        MatCopy_Dense,
567        MatGetDiag_Dense,MatScale_Dense,MatNorm_Dense,
568        0,0,
569        0, MatInsOpt_Dense,MatZero_Dense,MatZeroRows_Dense,0,
570        MatLUFactorSymbolic_Dense,MatLUFactorNumeric_Dense,
571        MatChFactorSymbolic_Dense,MatChFactorNumeric_Dense,
572        MatGetSize_Dense,MatGetSize_Dense,0,
573        0,0,MatGetArray_Dense
574 };
575 /*@
576     MatCreateSequentialDense - Creates a sequential dense matrix that
577         is stored in the usual Fortran 77 manner. Many of the matrix
578         operations use the BLAS and LAPACK routines.
579 
580   Input Parameters:
581 .   m, n - the number of rows and columns in the matrix.
582 
583   Output Parameter:
584 .  newmat - the matrix created.
585 
586   Keywords: dense matrix, lapack, blas
587 @*/
588 int MatCreateSequentialDense(int m,int n,Mat *newmat)
589 {
590   int       size = sizeof(Mat_Dense) + m*n*sizeof(Scalar);
591   Mat mat;
592   Mat_Dense    *l;
593   *newmat        = 0;
594   PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATDENSE,MPI_COMM_SELF);
595   PLogObjectCreate(mat);
596   l              = (Mat_Dense *) MALLOC(size); CHKPTR(l);
597   mat->ops       = &MatOps;
598   mat->destroy   = MatDestroy_Dense;
599   mat->view      = MatView_Dense;
600   mat->data      = (void *) l;
601   mat->factor    = 0;
602   mat->col       = 0;
603   mat->row       = 0;
604 
605   l->m           = m;
606   l->n           = n;
607   l->v           = (Scalar *) (l + 1);
608   l->pivots      = 0;
609   l->roworiented = 1;
610 
611   MEMSET(l->v,0,m*n*sizeof(Scalar));
612   *newmat = mat;
613   return 0;
614 }
615 
616 int MatCreate_Dense(Mat matin,Mat *newmat)
617 {
618   Mat_Dense *m = (Mat_Dense *) matin->data;
619   return MatCreateSequentialDense(m->m,m->n,newmat);
620 }
621