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