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