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