xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 9d00d63d4b4014b06eeb660baae4e576889960ed)
1 #ifndef lint
2 static char vcid[] = "$Id: dense.c,v 1.33 1995/05/02 23:38:15 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 #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 
23 static int MatGetInfo_Dense(Mat matin,MatInfoType flag,int *nz,
24                             int *nzalloc,int *mem)
25 {
26   Mat_Dense *mat = (Mat_Dense *) matin->data;
27   int    i,N = mat->m*mat->n,count = 0;
28   Scalar *v = mat->v;
29   for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;}
30   *nz = count; *nzalloc = N; *mem = N*sizeof(Scalar);
31   return 0;
32 }
33 
34 /* ---------------------------------------------------------------*/
35 /* COMMENT: I have chosen to hide column permutation in the pivots,
36    rather than put it in the Mat->col slot.*/
37 static int MatLUFactor_Dense(Mat matin,IS row,IS col)
38 {
39   Mat_Dense *mat = (Mat_Dense *) matin->data;
40   int    info;
41   if (!mat->pivots) {
42     mat->pivots = (int *) MALLOC( mat->m*sizeof(int) );
43     CHKPTR(mat->pivots);
44   }
45   LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info);
46   if (info) SETERR(1,"Bad LU factorization");
47   matin->factor = FACTOR_LU;
48   return 0;
49 }
50 static int MatLUFactorSymbolic_Dense(Mat matin,IS row,IS col,Mat *fact)
51 {
52   int ierr;
53   if ((ierr = MatConvert(matin,MATSAME,fact))) SETERR(ierr,0);
54   return 0;
55 }
56 static int MatLUFactorNumeric_Dense(Mat matin,Mat *fact)
57 {
58   return MatLUFactor(*fact,0,0);
59 }
60 static int MatCholeskyFactorSymbolic_Dense(Mat matin,IS row,Mat *fact)
61 {
62   int ierr;
63   if ((ierr = MatConvert(matin,MATSAME,fact))) SETERR(ierr,0);
64   return 0;
65 }
66 static int MatCholeskyFactorNumeric_Dense(Mat matin,Mat *fact)
67 {
68   return MatCholeskyFactor(*fact,0);
69 }
70 static int MatCholeskyFactor_Dense(Mat matin,IS perm)
71 {
72   Mat_Dense    *mat = (Mat_Dense *) matin->data;
73   int       info;
74   if (mat->pivots) {FREE(mat->pivots); mat->pivots = 0;}
75   LApotrf_("L",&mat->n,mat->v,&mat->m,&info);
76   if (info) SETERR(1,"Bad Cholesky factorization");
77   matin->factor = FACTOR_CHOLESKY;
78   return 0;
79 }
80 
81 static int MatSolve_Dense(Mat matin,Vec xx,Vec yy)
82 {
83   Mat_Dense *mat = (Mat_Dense *) matin->data;
84   int    one = 1, info;
85   Scalar *x, *y;
86   VecGetArray(xx,&x); VecGetArray(yy,&y);
87   MEMCPY(y,x,mat->m*sizeof(Scalar));
88   if (matin->factor == FACTOR_LU) {
89     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,
90               y, &mat->m, &info );
91   }
92   else if (matin->factor == FACTOR_CHOLESKY){
93     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
94               y, &mat->m, &info );
95   }
96   else SETERR(1,"Matrix must be factored to solve");
97   if (info) SETERR(1,"Bad solve");
98   return 0;
99 }
100 static int MatSolveTrans_Dense(Mat matin,Vec xx,Vec yy)
101 {
102   Mat_Dense *mat = (Mat_Dense *) matin->data;
103   int    one = 1, info;
104   Scalar *x, *y;
105   VecGetArray(xx,&x); VecGetArray(yy,&y);
106   MEMCPY(y,x,mat->m*sizeof(Scalar));
107   /* assume if pivots exist then LU else Cholesky */
108   if (mat->pivots) {
109     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,
110               y, &mat->m, &info );
111   }
112   else {
113     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
114               y, &mat->m, &info );
115   }
116   if (info) SETERR(1,"Bad solve");
117   return 0;
118 }
119 static int MatSolveAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy)
120 {
121   Mat_Dense *mat = (Mat_Dense *) matin->data;
122   int    one = 1, info,ierr;
123   Scalar *x, *y, sone = 1.0;
124   Vec    tmp = 0;
125   VecGetArray(xx,&x); VecGetArray(yy,&y);
126   if (yy == zz) {
127     ierr = VecCreate(yy,&tmp); CHKERR(ierr);
128     ierr = VecCopy(yy,tmp); CHKERR(ierr);
129   }
130   MEMCPY(y,x,mat->m*sizeof(Scalar));
131   /* assume if pivots exist then LU else Cholesky */
132   if (mat->pivots) {
133     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,
134               y, &mat->m, &info );
135   }
136   else {
137     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
138               y, &mat->m, &info );
139   }
140   if (info) SETERR(1,"Bad solve");
141   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
142   else VecAXPY(&sone,zz,yy);
143   return 0;
144 }
145 static int MatSolveTransAdd_Dense(Mat matin,Vec xx,Vec zz, Vec yy)
146 {
147   Mat_Dense  *mat = (Mat_Dense *) matin->data;
148   int     one = 1, info,ierr;
149   Scalar  *x, *y, sone = 1.0;
150   Vec     tmp;
151   VecGetArray(xx,&x); VecGetArray(yy,&y);
152   if (yy == zz) {
153     ierr = VecCreate(yy,&tmp); CHKERR(ierr);
154     ierr = VecCopy(yy,tmp); CHKERR(ierr);
155   }
156   MEMCPY(y,x,mat->m*sizeof(Scalar));
157   /* assume if pivots exist then LU else Cholesky */
158   if (mat->pivots) {
159     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,
160               y, &mat->m, &info );
161   }
162   else {
163     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
164               y, &mat->m, &info );
165   }
166   if (info) SETERR(1,"Bad solve");
167   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
168   else VecAXPY(&sone,zz,yy);
169   return 0;
170 }
171 /* ------------------------------------------------------------------*/
172 static int MatRelax_Dense(Mat matin,Vec bb,double omega,MatSORType flag,
173                           double shift,int its,Vec xx)
174 {
175   Mat_Dense *mat = (Mat_Dense *) matin->data;
176   Scalar *x, *b, *v = mat->v, zero = 0.0, xt;
177   int    o = 1,ierr, m = mat->m, i;
178 
179   if (flag & SOR_ZERO_INITIAL_GUESS) {
180     /* this is a hack fix, should have another version without
181        the second BLdot */
182     if ((ierr = VecSet(&zero,xx))) SETERR(ierr,0);
183   }
184   VecGetArray(xx,&x); VecGetArray(bb,&b);
185   while (its--) {
186     if (flag & SOR_FORWARD_SWEEP){
187       for ( i=0; i<m; i++ ) {
188         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
189         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
190       }
191     }
192     if (flag & SOR_BACKWARD_SWEEP) {
193       for ( i=m-1; i>=0; i-- ) {
194         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
195         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
196       }
197     }
198   }
199   return 0;
200 }
201 
202 /* -----------------------------------------------------------------*/
203 static int MatMultTrans_Dense(Mat matin,Vec xx,Vec yy)
204 {
205   Mat_Dense *mat = (Mat_Dense *) matin->data;
206   Scalar *v = mat->v, *x, *y;
207   int _One=1;Scalar _DOne=1.0, _DZero=0.0;
208   VecGetArray(xx,&x), VecGetArray(yy,&y);
209   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
210          x, &_One, &_DZero, y, &_One );
211   return 0;
212 }
213 static int MatMult_Dense(Mat matin,Vec xx,Vec yy)
214 {
215   Mat_Dense *mat = (Mat_Dense *) matin->data;
216   Scalar *v = mat->v, *x, *y;
217   int _One=1;Scalar _DOne=1.0, _DZero=0.0;
218   VecGetArray(xx,&x); VecGetArray(yy,&y);
219   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
220          x, &_One, &_DZero, y, &_One );
221   return 0;
222 }
223 static int MatMultAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy)
224 {
225   Mat_Dense *mat = (Mat_Dense *) matin->data;
226   Scalar *v = mat->v, *x, *y, *z;
227   int    _One=1; Scalar _DOne=1.0;
228   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
229   if (zz != yy) MEMCPY(y,z,mat->m*sizeof(Scalar));
230   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
231          x, &_One, &_DOne, y, &_One );
232   return 0;
233 }
234 static int MatMultTransAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy)
235 {
236   Mat_Dense *mat = (Mat_Dense *) matin->data;
237   Scalar *v = mat->v, *x, *y, *z;
238   int    _One=1;
239   Scalar _DOne=1.0;
240   VecGetArray(xx,&x); VecGetArray(yy,&y);
241   VecGetArray(zz,&z);
242   if (zz != yy) MEMCPY(y,z,mat->m*sizeof(Scalar));
243   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
244          x, &_One, &_DOne, y, &_One );
245   return 0;
246 }
247 
248 /* -----------------------------------------------------------------*/
249 static int MatGetRow_Dense(Mat matin,int row,int *ncols,int **cols,
250                         Scalar **vals)
251 {
252   Mat_Dense *mat = (Mat_Dense *) matin->data;
253   Scalar *v;
254   int    i;
255   *ncols = mat->n;
256   if (cols) {
257     *cols = (int *) MALLOC(mat->n*sizeof(int)); CHKPTR(*cols);
258     for ( i=0; i<mat->n; i++ ) *cols[i] = i;
259   }
260   if (vals) {
261     *vals = (Scalar *) MALLOC(mat->n*sizeof(Scalar)); CHKPTR(*vals);
262     v = mat->v + row;
263     for ( i=0; i<mat->n; i++ ) {*vals[i] = *v; v += mat->m;}
264   }
265   return 0;
266 }
267 static int MatRestoreRow_Dense(Mat matin,int row,int *ncols,int **cols,
268                             Scalar **vals)
269 {
270   if (cols) { FREE(*cols); }
271   if (vals) { FREE(*vals); }
272   return 0;
273 }
274 /* ----------------------------------------------------------------*/
275 static int MatInsert_Dense(Mat matin,int m,int *indexm,int n,
276                         int *indexn,Scalar *v,InsertMode addv)
277 {
278   Mat_Dense *mat = (Mat_Dense *) matin->data;
279   int    i,j;
280 
281   if (!mat->roworiented) {
282     if (addv == INSERTVALUES) {
283       for ( j=0; j<n; j++ ) {
284         if (indexn[j] < 0) {*v += m; continue;}
285         for ( i=0; i<m; i++ ) {
286           if (indexm[i] < 0) {*v++; continue;}
287           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
288         }
289       }
290     }
291     else {
292       for ( j=0; j<n; j++ ) {
293         if (indexn[j] < 0) {*v += m; continue;}
294         for ( i=0; i<m; i++ ) {
295           if (indexm[i] < 0) {*v++; continue;}
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         if (indexm[i] < 0) {*v += n; continue;}
305         for ( j=0; j<n; j++ ) {
306           if (indexn[j] < 0) {*v++; continue;}
307           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
308         }
309       }
310     }
311     else {
312       for ( i=0; i<m; i++ ) {
313         if (indexm[i] < 0) {*v += n; continue;}
314         for ( j=0; j<n; j++ ) {
315           if (indexn[j] < 0) {*v++; continue;}
316           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
317         }
318       }
319     }
320   }
321   return 0;
322 }
323 
324 /* -----------------------------------------------------------------*/
325 static int MatCopy_Dense_Private(Mat matin,Mat *newmat)
326 {
327   Mat_Dense *mat = (Mat_Dense *) matin->data;
328   int ierr;
329   Mat newi;
330   Mat_Dense *l;
331   if ((ierr = MatCreateSequentialDense(matin->comm,mat->m,mat->n,&newi)))
332                                                           SETERR(ierr,0);
333   l = (Mat_Dense *) newi->data;
334   MEMCPY(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
335   *newmat = newi;
336   return 0;
337 }
338 #include "viewer.h"
339 
340 int MatView_Dense(PetscObject obj,Viewer ptr)
341 {
342   Mat         matin = (Mat) obj;
343   Mat_Dense   *mat = (Mat_Dense *) matin->data;
344   Scalar      *v;
345   int         i,j;
346   PetscObject vobj = (PetscObject) ptr;
347 
348   if (ptr == 0) {
349     ptr = STDOUT_VIEWER; vobj = (PetscObject) ptr;
350   }
351   if (vobj->cookie == VIEWER_COOKIE && vobj->type == MATLAB_VIEWER) {
352     return ViewerMatlabPutArray_Private(ptr,mat->m,mat->n,mat->v);
353   }
354   else {
355     FILE *fd = ViewerFileGetPointer_Private(ptr);
356     for ( i=0; i<mat->m; i++ ) {
357       v = mat->v + i;
358       for ( j=0; j<mat->n; j++ ) {
359 #if defined(PETSC_COMPLEX)
360         fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m;
361 #else
362         fprintf(fd,"%6.4e ",*v); v += mat->m;
363 #endif
364       }
365       fprintf(fd,"\n");
366     }
367   }
368   return 0;
369 }
370 
371 
372 static int MatDestroy_Dense(PetscObject obj)
373 {
374   Mat    mat = (Mat) obj;
375   Mat_Dense *l = (Mat_Dense *) mat->data;
376 #if defined(PETSC_LOG)
377   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
378 #endif
379   if (l->pivots) FREE(l->pivots);
380   FREE(l);
381   PLogObjectDestroy(mat);
382   PETSCHEADERDESTROY(mat);
383   return 0;
384 }
385 
386 static int MatTrans_Dense(Mat matin)
387 {
388   Mat_Dense *mat = (Mat_Dense *) matin->data;
389   int    k,j;
390   Scalar *v = mat->v, tmp;
391   if (mat->m != mat->n) {
392     SETERR(1,"Cannot transpose rectangular dense matrix");
393   }
394   for ( j=0; j<mat->m; j++ ) {
395     for ( k=0; k<j; k++ ) {
396       tmp = v[j + k*mat->n];
397       v[j + k*mat->n] = v[k + j*mat->n];
398       v[k + j*mat->n] = tmp;
399     }
400   }
401   return 0;
402 }
403 
404 static int MatEqual_Dense(Mat matin1,Mat matin2)
405 {
406   Mat_Dense *mat1 = (Mat_Dense *) matin1->data;
407   Mat_Dense *mat2 = (Mat_Dense *) matin2->data;
408   int    i;
409   Scalar *v1 = mat1->v, *v2 = mat2->v;
410   if (mat1->m != mat2->m) return 0;
411   if (mat1->n != mat2->n) return 0;
412   for ( i=0; i<mat1->m*mat1->n; i++ ) {
413     if (*v1 != *v2) return 0;
414     v1++; v2++;
415   }
416   return 1;
417 }
418 
419 static int MatGetDiagonal_Dense(Mat matin,Vec v)
420 {
421   Mat_Dense *mat = (Mat_Dense *) matin->data;
422   int    i, n;
423   Scalar *x;
424   CHKTYPE(v,SEQVECTOR);
425   VecGetArray(v,&x); VecGetSize(v,&n);
426   if (n != mat->m) SETERR(1,"Nonconforming matrix and vector");
427   for ( i=0; i<mat->m; i++ ) {
428     x[i] = mat->v[i*mat->m + i];
429   }
430   return 0;
431 }
432 
433 static int MatScale_Dense(Mat matin,Vec ll,Vec rr)
434 {
435   Mat_Dense *mat = (Mat_Dense *) matin->data;
436   Scalar *l,*r,x,*v;
437   int    i,j,m = mat->m, n = mat->n;
438   if (ll) {
439     VecGetArray(ll,&l); VecGetSize(ll,&m);
440     if (m != mat->m) SETERR(1,"Left scaling vector wrong length");
441     for ( i=0; i<m; i++ ) {
442       x = l[i];
443       v = mat->v + i;
444       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
445     }
446   }
447   if (rr) {
448     VecGetArray(rr,&r); VecGetSize(rr,&n);
449     if (n != mat->n) SETERR(1,"Right scaling vector wrong length");
450     for ( i=0; i<n; i++ ) {
451       x = r[i];
452       v = mat->v + i*m;
453       for ( j=0; j<m; j++ ) { (*v++) *= x;}
454     }
455   }
456   return 0;
457 }
458 
459 
460 static int MatNorm_Dense(Mat matin,int 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     SETERR(1,"No support for the 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   MEMSET(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); CHKERR(ierr);
533   ierr = ISGetIndices(is,&rows); CHKERR(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   SETERR(1,"MatGetSubMatrixInPlace_Dense not yet 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, jstart;
573   Scalar  *vwork, *val;
574   Mat     newmat;
575 
576   ierr = ISGetIndices(isrow,&irow); CHKERR(ierr);
577   ierr = ISGetIndices(iscol,&icol); CHKERR(ierr);
578   ierr = ISGetSize(isrow,&nrows); CHKERR(ierr);
579   ierr = ISGetSize(iscol,&ncols); CHKERR(ierr);
580 
581   smap = (int *) MALLOC(oldcols*sizeof(int)); CHKPTR(smap);
582   cwork = (int *) MALLOC(ncols*sizeof(int)); CHKPTR(cwork);
583   vwork = (Scalar *) MALLOC(ncols*sizeof(Scalar)); CHKPTR(vwork);
584   memset((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          CHKERR(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            CHKERR(ierr);
601   }
602   ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERR(ierr);
603   ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERR(ierr);
604 
605   /* Free work space */
606   FREE(smap); FREE(cwork); FREE(vwork);
607   ierr = ISRestoreIndices(isrow,&irow); CHKERR(ierr);
608   ierr = ISRestoreIndices(iscol,&icol); CHKERR(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        MatTrans_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        MatCopy_Dense_Private};
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: Mat, dense, matrix, LAPACK, BLAS
648 
649 .seealso: MatCreateSequentialAIJ()
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 *) MALLOC(size); CHKPTR(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   MEMSET(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