xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 0b6271093bf45363169cd314acc5fc2d6b1e3e10)
1 #ifndef lint
2 static char vcid[] = "$Id: dense.c,v 1.36 1995/05/06 20:11:13 curfman 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         if (indexn[j] < 0) {*v += m; continue;}
284         for ( i=0; i<m; i++ ) {
285           if (indexm[i] < 0) {*v++; continue;}
286           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
287         }
288       }
289     }
290     else {
291       for ( j=0; j<n; j++ ) {
292         if (indexn[j] < 0) {*v += m; continue;}
293         for ( i=0; i<m; i++ ) {
294           if (indexm[i] < 0) {*v++; continue;}
295           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
296         }
297       }
298     }
299   }
300   else {
301     if (addv == INSERTVALUES) {
302       for ( i=0; i<m; i++ ) {
303         if (indexm[i] < 0) {*v += n; continue;}
304         for ( j=0; j<n; j++ ) {
305           if (indexn[j] < 0) {*v++; continue;}
306           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
307         }
308       }
309     }
310     else {
311       for ( i=0; i<m; i++ ) {
312         if (indexm[i] < 0) {*v += n; continue;}
313         for ( j=0; j<n; j++ ) {
314           if (indexn[j] < 0) {*v++; continue;}
315           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
316         }
317       }
318     }
319   }
320   return 0;
321 }
322 
323 /* -----------------------------------------------------------------*/
324 static int MatCopy_Dense_Private(Mat matin,Mat *newmat)
325 {
326   Mat_Dense *mat = (Mat_Dense *) matin->data;
327   int ierr;
328   Mat newi;
329   Mat_Dense *l;
330   if ((ierr = MatCreateSequentialDense(matin->comm,mat->m,mat->n,&newi)))
331                                                           SETERR(ierr,0);
332   l = (Mat_Dense *) newi->data;
333   MEMCPY(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
334   *newmat = newi;
335   return 0;
336 }
337 #include "viewer.h"
338 
339 int MatView_Dense(PetscObject obj,Viewer ptr)
340 {
341   Mat         matin = (Mat) obj;
342   Mat_Dense   *mat = (Mat_Dense *) matin->data;
343   Scalar      *v;
344   int         i,j;
345   PetscObject vobj = (PetscObject) ptr;
346 
347   if (ptr == 0) {
348     ptr = STDOUT_VIEWER; vobj = (PetscObject) ptr;
349   }
350   if (vobj->cookie == VIEWER_COOKIE && vobj->type == MATLAB_VIEWER) {
351     return ViewerMatlabPutArray_Private(ptr,mat->m,mat->n,mat->v);
352   }
353   else {
354     FILE *fd = ViewerFileGetPointer_Private(ptr);
355     for ( i=0; i<mat->m; i++ ) {
356       v = mat->v + i;
357       for ( j=0; j<mat->n; j++ ) {
358 #if defined(PETSC_COMPLEX)
359         fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m;
360 #else
361         fprintf(fd,"%6.4e ",*v); v += mat->m;
362 #endif
363       }
364       fprintf(fd,"\n");
365     }
366   }
367   return 0;
368 }
369 
370 
371 static int MatDestroy_Dense(PetscObject obj)
372 {
373   Mat    mat = (Mat) obj;
374   Mat_Dense *l = (Mat_Dense *) mat->data;
375 #if defined(PETSC_LOG)
376   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
377 #endif
378   if (l->pivots) FREE(l->pivots);
379   FREE(l);
380   PLogObjectDestroy(mat);
381   PETSCHEADERDESTROY(mat);
382   return 0;
383 }
384 
385 static int MatTrans_Dense(Mat matin)
386 {
387   Mat_Dense *mat = (Mat_Dense *) matin->data;
388   int    k,j;
389   Scalar *v = mat->v, tmp;
390   if (mat->m != mat->n) {
391     SETERR(1,"Cannot transpose rectangular dense matrix");
392   }
393   for ( j=0; j<mat->m; j++ ) {
394     for ( k=0; k<j; k++ ) {
395       tmp = v[j + k*mat->n];
396       v[j + k*mat->n] = v[k + j*mat->n];
397       v[k + j*mat->n] = tmp;
398     }
399   }
400   return 0;
401 }
402 
403 static int MatEqual_Dense(Mat matin1,Mat matin2)
404 {
405   Mat_Dense *mat1 = (Mat_Dense *) matin1->data;
406   Mat_Dense *mat2 = (Mat_Dense *) matin2->data;
407   int    i;
408   Scalar *v1 = mat1->v, *v2 = mat2->v;
409   if (mat1->m != mat2->m) return 0;
410   if (mat1->n != mat2->n) return 0;
411   for ( i=0; i<mat1->m*mat1->n; i++ ) {
412     if (*v1 != *v2) return 0;
413     v1++; v2++;
414   }
415   return 1;
416 }
417 
418 static int MatGetDiagonal_Dense(Mat matin,Vec v)
419 {
420   Mat_Dense *mat = (Mat_Dense *) matin->data;
421   int    i, n;
422   Scalar *x;
423   CHKTYPE(v,SEQVECTOR);
424   VecGetArray(v,&x); VecGetSize(v,&n);
425   if (n != mat->m) SETERR(1,"Nonconforming matrix and vector");
426   for ( i=0; i<mat->m; i++ ) {
427     x[i] = mat->v[i*mat->m + i];
428   }
429   return 0;
430 }
431 
432 static int MatScale_Dense(Mat matin,Vec ll,Vec rr)
433 {
434   Mat_Dense *mat = (Mat_Dense *) matin->data;
435   Scalar *l,*r,x,*v;
436   int    i,j,m = mat->m, n = mat->n;
437   if (ll) {
438     VecGetArray(ll,&l); VecGetSize(ll,&m);
439     if (m != mat->m) SETERR(1,"Left scaling vector wrong length");
440     for ( i=0; i<m; i++ ) {
441       x = l[i];
442       v = mat->v + i;
443       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
444     }
445   }
446   if (rr) {
447     VecGetArray(rr,&r); VecGetSize(rr,&n);
448     if (n != mat->n) SETERR(1,"Right scaling vector wrong length");
449     for ( i=0; i<n; i++ ) {
450       x = r[i];
451       v = mat->v + i*m;
452       for ( j=0; j<m; j++ ) { (*v++) *= x;}
453     }
454   }
455   return 0;
456 }
457 
458 
459 static int MatNorm_Dense(Mat matin,int type,double *norm)
460 {
461   Mat_Dense *mat = (Mat_Dense *) matin->data;
462   Scalar *v = mat->v;
463   double sum = 0.0;
464   int    i, j;
465   if (type == NORM_FROBENIUS) {
466     for (i=0; i<mat->n*mat->m; i++ ) {
467 #if defined(PETSC_COMPLEX)
468       sum += real(conj(*v)*(*v)); v++;
469 #else
470       sum += (*v)*(*v); v++;
471 #endif
472     }
473     *norm = sqrt(sum);
474   }
475   else if (type == NORM_1) {
476     *norm = 0.0;
477     for ( j=0; j<mat->n; j++ ) {
478       sum = 0.0;
479       for ( i=0; i<mat->m; i++ ) {
480 #if defined(PETSC_COMPLEX)
481         sum += abs(*v++);
482 #else
483         sum += fabs(*v++);
484 #endif
485       }
486       if (sum > *norm) *norm = sum;
487     }
488   }
489   else if (type == NORM_INFINITY) {
490     *norm = 0.0;
491     for ( j=0; j<mat->m; j++ ) {
492       v = mat->v + j;
493       sum = 0.0;
494       for ( i=0; i<mat->n; i++ ) {
495 #if defined(PETSC_COMPLEX)
496         sum += abs(*v); v += mat->m;
497 #else
498         sum += fabs(*v); v += mat->m;
499 #endif
500       }
501       if (sum > *norm) *norm = sum;
502     }
503   }
504   else {
505     SETERR(1,"No support for the two norm yet");
506   }
507   return 0;
508 }
509 
510 static int MatSetOption_Dense(Mat aijin,MatOption op)
511 {
512   Mat_Dense *aij = (Mat_Dense *) aijin->data;
513   if (op == ROW_ORIENTED)            aij->roworiented = 1;
514   else if (op == COLUMN_ORIENTED)    aij->roworiented = 0;
515   /* doesn't care about sorted rows or columns */
516   return 0;
517 }
518 
519 static int MatZero_Dense(Mat A)
520 {
521   Mat_Dense *l = (Mat_Dense *) A->data;
522   MEMSET(l->v,0,l->m*l->n*sizeof(Scalar));
523   return 0;
524 }
525 
526 static int MatZeroRows_Dense(Mat A,IS is,Scalar *diag)
527 {
528   Mat_Dense *l = (Mat_Dense *) A->data;
529   int    n = l->n, i, j,ierr,N, *rows;
530   Scalar *slot;
531   ierr = ISGetLocalSize(is,&N); CHKERR(ierr);
532   ierr = ISGetIndices(is,&rows); CHKERR(ierr);
533   for ( i=0; i<N; i++ ) {
534     slot = l->v + rows[i];
535     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
536   }
537   if (diag) {
538     for ( i=0; i<N; i++ ) {
539       slot = l->v + (n+1)*rows[i];
540       *slot = *diag;
541     }
542   }
543   ISRestoreIndices(is,&rows);
544   return 0;
545 }
546 
547 static int MatGetSize_Dense(Mat matin,int *m,int *n)
548 {
549   Mat_Dense *mat = (Mat_Dense *) matin->data;
550   *m = mat->m; *n = mat->n;
551   return 0;
552 }
553 
554 static int MatGetArray_Dense(Mat matin,Scalar **array)
555 {
556   Mat_Dense *mat = (Mat_Dense *) matin->data;
557   *array = mat->v;
558   return 0;
559 }
560 
561 
562 static int MatGetSubMatrixInPlace_Dense(Mat matin,IS isrow,IS iscol)
563 {
564   SETERR(1,"MatGetSubMatrixInPlace_Dense not yet done.");
565 }
566 
567 static int MatGetSubMatrix_Dense(Mat matin,IS isrow,IS iscol,Mat *submat)
568 {
569   Mat_Dense *mat = (Mat_Dense *) matin->data;
570   int     nznew, *smap, i, j, ierr, oldcols = mat->n;
571   int     *irow, *icol, nrows, ncols, *cwork;
572   Scalar  *vwork, *val;
573   Mat     newmat;
574 
575   ierr = ISGetIndices(isrow,&irow); CHKERR(ierr);
576   ierr = ISGetIndices(iscol,&icol); CHKERR(ierr);
577   ierr = ISGetSize(isrow,&nrows); CHKERR(ierr);
578   ierr = ISGetSize(iscol,&ncols); CHKERR(ierr);
579 
580   smap = (int *) MALLOC(oldcols*sizeof(int)); CHKPTR(smap);
581   cwork = (int *) MALLOC(ncols*sizeof(int)); CHKPTR(cwork);
582   vwork = (Scalar *) MALLOC(ncols*sizeof(Scalar)); CHKPTR(vwork);
583   memset((char*)smap,0,oldcols*sizeof(int));
584   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
585 
586   /* Create and fill new matrix */
587   ierr = MatCreateSequentialDense(matin->comm,nrows,ncols,&newmat);
588          CHKERR(ierr);
589   for (i=0; i<nrows; i++) {
590     nznew = 0;
591     val   = mat->v + irow[i];
592     for (j=0; j<oldcols; j++) {
593       if (smap[j]) {
594         cwork[nznew]   = smap[j] - 1;
595         vwork[nznew++] = val[j * mat->m];
596       }
597     }
598     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERTVALUES);
599            CHKERR(ierr);
600   }
601   ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERR(ierr);
602   ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERR(ierr);
603 
604   /* Free work space */
605   FREE(smap); FREE(cwork); FREE(vwork);
606   ierr = ISRestoreIndices(isrow,&irow); CHKERR(ierr);
607   ierr = ISRestoreIndices(iscol,&icol); CHKERR(ierr);
608   *submat = newmat;
609   return 0;
610 }
611 
612 /* -------------------------------------------------------------------*/
613 static struct _MatOps MatOps = {MatInsert_Dense,
614        MatGetRow_Dense, MatRestoreRow_Dense,
615        MatMult_Dense, MatMultAdd_Dense,
616        MatMultTrans_Dense, MatMultTransAdd_Dense,
617        MatSolve_Dense,MatSolveAdd_Dense,
618        MatSolveTrans_Dense,MatSolveTransAdd_Dense,
619        MatLUFactor_Dense,MatCholeskyFactor_Dense,
620        MatRelax_Dense,
621        MatTrans_Dense,
622        MatGetInfo_Dense,MatEqual_Dense,
623        MatGetDiagonal_Dense,MatScale_Dense,MatNorm_Dense,
624        0,0,
625        0, MatSetOption_Dense,MatZero_Dense,MatZeroRows_Dense,0,
626        MatLUFactorSymbolic_Dense,MatLUFactorNumeric_Dense,
627        MatCholeskyFactorSymbolic_Dense,MatCholeskyFactorNumeric_Dense,
628        MatGetSize_Dense,MatGetSize_Dense,0,
629        0,0,MatGetArray_Dense,0,0,
630        MatGetSubMatrix_Dense,MatGetSubMatrixInPlace_Dense,
631        MatCopy_Dense_Private};
632 
633 /*@
634    MatCreateSequentialDense - Creates a sequential dense matrix that
635    is stored in column major order (the usual Fortran 77 manner). Many
636    of the matrix operations use the BLAS and LAPACK routines.
637 
638    Input Parameters:
639 .  comm - MPI communicator, set to MPI_COMM_SELF
640 .  m - number of rows
641 .  n - number of column
642 
643    Output Parameter:
644 .  newmat - the matrix
645 
646 .keywords: dense, matrix, LAPACK, BLAS
647 
648 .seealso: MatCreate(), MatSetValues()
649 @*/
650 int MatCreateSequentialDense(MPI_Comm comm,int m,int n,Mat *newmat)
651 {
652   int       size = sizeof(Mat_Dense) + m*n*sizeof(Scalar);
653   Mat mat;
654   Mat_Dense    *l;
655   *newmat        = 0;
656   PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATDENSE,comm);
657   PLogObjectCreate(mat);
658   l              = (Mat_Dense *) MALLOC(size); CHKPTR(l);
659   mat->ops       = &MatOps;
660   mat->destroy   = MatDestroy_Dense;
661   mat->view      = MatView_Dense;
662   mat->data      = (void *) l;
663   mat->factor    = 0;
664 
665   l->m           = m;
666   l->n           = n;
667   l->v           = (Scalar *) (l + 1);
668   l->pivots      = 0;
669   l->roworiented = 1;
670 
671   MEMSET(l->v,0,m*n*sizeof(Scalar));
672   *newmat = mat;
673   return 0;
674 }
675 
676 int MatCreate_Dense(Mat matin,Mat *newmat)
677 {
678   Mat_Dense *m = (Mat_Dense *) matin->data;
679   return MatCreateSequentialDense(matin->comm,m->m,m->n,newmat);
680 }
681