xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 80cd9d9306eb90d00e6fd99dae46147660087c21) !
1 #ifndef lint
2 static char vcid[] = "$Id: dense.c,v 1.100 1996/04/05 19:41:25 curfman Exp curfman $";
3 #endif
4 /*
5      Defines the basic matrix operations for sequential dense.
6 */
7 
8 #include "dense.h"
9 #include "pinclude/plapack.h"
10 #include "pinclude/pviewer.h"
11 
12 int MatAXPY_SeqDense(Scalar *alpha,Mat X,Mat Y)
13 {
14   Mat_SeqDense *x = (Mat_SeqDense*) X->data,*y = (Mat_SeqDense*) Y->data;
15   int          N = x->m*x->n, one = 1;
16   BLaxpy_( &N, alpha, x->v, &one, y->v, &one );
17   PLogFlops(2*N-1);
18   return 0;
19 }
20 
21 static int MatGetInfo_SeqDense(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem)
22 {
23   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
24   int          i,N = mat->m*mat->n,count = 0;
25   Scalar       *v = mat->v;
26   for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;}
27   if (nz)      *nz      = count;
28   if (nzalloc) *nzalloc = N;
29   if (mem)     *mem     = (int)A->mem;
30   return 0;
31 }
32 
33 static int MatScale_SeqDense(Scalar *alpha,Mat inA)
34 {
35   Mat_SeqDense *a = (Mat_SeqDense *) inA->data;
36   int          one = 1, nz;
37 
38   nz = a->m*a->n;
39   BLscal_( &nz, alpha, a->v, &one );
40   PLogFlops(nz);
41   return 0;
42 }
43 
44 /* ---------------------------------------------------------------*/
45 /* COMMENT: I have chosen to hide column permutation in the pivots,
46    rather than put it in the Mat->col slot.*/
47 static int MatLUFactor_SeqDense(Mat A,IS row,IS col,double f)
48 {
49   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
50   int          info;
51 
52   if (!mat->pivots) {
53     mat->pivots = (int *) PetscMalloc(mat->m*sizeof(int));CHKPTRQ(mat->pivots);
54     PLogObjectMemory(A,mat->m*sizeof(int));
55   }
56   LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info);
57   if (info) SETERRQ(1,"MatLUFactor_SeqDense:Bad LU factorization");
58   A->factor = FACTOR_LU;
59   PLogFlops((2*mat->n*mat->n*mat->n)/3);
60   return 0;
61 }
62 static int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,double f,Mat *fact)
63 {
64   return MatConvert(A,MATSAME,fact);
65 }
66 static int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact)
67 {
68   return MatLUFactor(*fact,0,0,1.0);
69 }
70 static int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact)
71 {
72   return MatConvert(A,MATSAME,fact);
73 }
74 static int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact)
75 {
76   return MatCholeskyFactor(*fact,0,1.0);
77 }
78 static int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f)
79 {
80   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
81   int           info;
82 
83   if (mat->pivots) {
84     PetscFree(mat->pivots);
85     PLogObjectMemory(A,-mat->m*sizeof(int));
86     mat->pivots = 0;
87   }
88   LApotrf_("L",&mat->n,mat->v,&mat->m,&info);
89   if (info) SETERRQ(1,"MatCholeskyFactor_SeqDense:Bad factorization");
90   A->factor = FACTOR_CHOLESKY;
91   PLogFlops((mat->n*mat->n*mat->n)/3);
92   return 0;
93 }
94 
95 static int MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
96 {
97   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
98   int          one = 1, info, ierr;
99   Scalar       *x, *y;
100 
101   ierr = VecGetArray(xx,&x); CHKERRQ(ierr); ierr = VecGetArray(yy,&y); CHKERRQ(ierr);
102   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
103   if (A->factor == FACTOR_LU) {
104     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
105   }
106   else if (A->factor == FACTOR_CHOLESKY){
107     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
108   }
109   else SETERRQ(1,"MatSolve_SeqDense:Matrix must be factored to solve");
110   if (info) SETERRQ(1,"MatSolve_SeqDense:Bad solve");
111   PLogFlops(mat->n*mat->n - mat->n);
112   return 0;
113 }
114 static int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy)
115 {
116   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
117   int          one = 1, info;
118   Scalar       *x, *y;
119 
120   VecGetArray(xx,&x); VecGetArray(yy,&y);
121   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
122   /* assume if pivots exist then use LU; else Cholesky */
123   if (mat->pivots) {
124     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
125   }
126   else {
127     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
128   }
129   if (info) SETERRQ(1,"MatSolveTrans_SeqDense:Bad solve");
130   PLogFlops(mat->n*mat->n - mat->n);
131   return 0;
132 }
133 static int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
134 {
135   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
136   int          one = 1, info,ierr;
137   Scalar       *x, *y, sone = 1.0;
138   Vec          tmp = 0;
139 
140   VecGetArray(xx,&x); VecGetArray(yy,&y);
141   if (yy == zz) {
142     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
143     PLogObjectParent(A,tmp);
144     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
145   }
146   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
147   /* assume if pivots exist then use LU; else Cholesky */
148   if (mat->pivots) {
149     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
150   }
151   else {
152     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
153   }
154   if (info) SETERRQ(1,"MatSolveAdd_SeqDense:Bad solve");
155   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
156   else VecAXPY(&sone,zz,yy);
157   PLogFlops(mat->n*mat->n - mat->n);
158   return 0;
159 }
160 
161 static int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy)
162 {
163   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
164   int           one = 1, info,ierr;
165   Scalar        *x, *y, sone = 1.0;
166   Vec           tmp;
167 
168   VecGetArray(xx,&x); VecGetArray(yy,&y);
169   if (yy == zz) {
170     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
171     PLogObjectParent(A,tmp);
172     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
173   }
174   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
175   /* assume if pivots exist then use LU; else Cholesky */
176   if (mat->pivots) {
177     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
178   }
179   else {
180     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
181   }
182   if (info) SETERRQ(1,"MatSolveTransAdd_SeqDense:Bad solve");
183   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
184   else VecAXPY(&sone,zz,yy);
185   PLogFlops(mat->n*mat->n - mat->n);
186   return 0;
187 }
188 /* ------------------------------------------------------------------*/
189 static int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag,
190                           double shift,int its,Vec xx)
191 {
192   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
193   Scalar       *x, *b, *v = mat->v, zero = 0.0, xt;
194   int          o = 1,ierr, m = mat->m, i;
195 
196   if (flag & SOR_ZERO_INITIAL_GUESS) {
197     /* this is a hack fix, should have another version without
198        the second BLdot */
199     ierr = VecSet(&zero,xx); CHKERRQ(ierr);
200   }
201   VecGetArray(xx,&x); VecGetArray(bb,&b);
202   while (its--) {
203     if (flag & SOR_FORWARD_SWEEP){
204       for ( i=0; i<m; i++ ) {
205 #if defined(PETSC_COMPLEX)
206         /* cannot use BLAS dot for complex because compiler/linker is
207            not happy about returning a double complex */
208         int    _i;
209         Scalar sum = b[i];
210         for ( _i=0; _i<m; _i++ ) {
211           sum -= conj(v[i+_i*m])*x[_i];
212         }
213         xt = sum;
214 #else
215         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
216 #endif
217         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
218       }
219     }
220     if (flag & SOR_BACKWARD_SWEEP) {
221       for ( i=m-1; i>=0; i-- ) {
222 #if defined(PETSC_COMPLEX)
223         /* cannot use BLAS dot for complex because compiler/linker is
224            not happy about returning a double complex */
225         int    _i;
226         Scalar sum = b[i];
227         for ( _i=0; _i<m; _i++ ) {
228           sum -= conj(v[i+_i*m])*x[_i];
229         }
230         xt = sum;
231 #else
232         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
233 #endif
234         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
235       }
236     }
237   }
238   return 0;
239 }
240 
241 /* -----------------------------------------------------------------*/
242 static int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy)
243 {
244   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
245   Scalar       *v = mat->v, *x, *y;
246   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
247   VecGetArray(xx,&x), VecGetArray(yy,&y);
248   LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One);
249   PLogFlops(2*mat->m*mat->n - mat->n);
250   return 0;
251 }
252 static int MatMult_SeqDense(Mat A,Vec xx,Vec yy)
253 {
254   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
255   Scalar       *v = mat->v, *x, *y;
256   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
257   VecGetArray(xx,&x); VecGetArray(yy,&y);
258   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One);
259   PLogFlops(2*mat->m*mat->n - mat->m);
260   return 0;
261 }
262 static int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
263 {
264   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
265   Scalar       *v = mat->v, *x, *y, *z;
266   int          _One=1; Scalar _DOne=1.0;
267   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
268   if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar));
269   LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One);
270   PLogFlops(2*mat->m*mat->n);
271   return 0;
272 }
273 static int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
274 {
275   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
276   Scalar       *v = mat->v, *x, *y, *z;
277   int          _One=1;
278   Scalar       _DOne=1.0;
279   VecGetArray(xx,&x); VecGetArray(yy,&y);
280   VecGetArray(zz,&z);
281   if (zz != yy) PetscMemcpy(y,z,mat->n*sizeof(Scalar));
282   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One);
283   PLogFlops(2*mat->m*mat->n);
284   return 0;
285 }
286 
287 /* -----------------------------------------------------------------*/
288 static int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
289 {
290   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
291   Scalar       *v;
292   int          i;
293 
294   *ncols = mat->n;
295   if (cols) {
296     *cols = (int *) PetscMalloc(mat->n*sizeof(int)); CHKPTRQ(*cols);
297     for ( i=0; i<mat->n; i++ ) (*cols)[i] = i;
298   }
299   if (vals) {
300     *vals = (Scalar *) PetscMalloc(mat->n*sizeof(Scalar)); CHKPTRQ(*vals);
301     v = mat->v + row;
302     for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;}
303   }
304   return 0;
305 }
306 static int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
307 {
308   if (cols) { PetscFree(*cols); }
309   if (vals) { PetscFree(*vals); }
310   return 0;
311 }
312 /* ----------------------------------------------------------------*/
313 static int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n,
314                                     int *indexn,Scalar *v,InsertMode addv)
315 {
316   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
317   int          i,j;
318 
319   if (!mat->roworiented) {
320     if (addv == INSERT_VALUES) {
321       for ( j=0; j<n; j++ ) {
322         for ( i=0; i<m; i++ ) {
323           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
324         }
325       }
326     }
327     else {
328       for ( j=0; j<n; j++ ) {
329         for ( i=0; i<m; i++ ) {
330           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
331         }
332       }
333     }
334   }
335   else {
336     if (addv == INSERT_VALUES) {
337       for ( i=0; i<m; i++ ) {
338         for ( j=0; j<n; j++ ) {
339           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
340         }
341       }
342     }
343     else {
344       for ( i=0; i<m; i++ ) {
345         for ( j=0; j<n; j++ ) {
346           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
347         }
348       }
349     }
350   }
351   return 0;
352 }
353 
354 static int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v)
355 {
356   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
357   int          i, j;
358   Scalar       *vpt = v;
359 
360   /* row-oriented output */
361   for ( i=0; i<m; i++ ) {
362     for ( j=0; j<n; j++ ) {
363       *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]];
364     }
365   }
366   return 0;
367 }
368 
369 /* -----------------------------------------------------------------*/
370 static int MatConvertSameType_SeqDense(Mat A,Mat *newmat,int cpvalues)
371 {
372   Mat_SeqDense *mat = (Mat_SeqDense *) A->data, *l;
373   int          ierr;
374   Mat          newi;
375 
376   ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,PETSC_NULL,&newi); CHKERRQ(ierr);
377   l = (Mat_SeqDense *) newi->data;
378   if (cpvalues == COPY_VALUES) {
379     PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
380   }
381   newi->assembled = PETSC_TRUE;
382   *newmat = newi;
383   return 0;
384 }
385 
386 #include "sys.h"
387 
388 int MatLoad_SeqDense(Viewer viewer,MatType type,Mat *A)
389 {
390   Mat_SeqDense *a;
391   Mat          B;
392   int          *scols, i, j, nz, ierr, fd, header[4], size;
393   int          *rowlengths = 0, M, N, *cols;
394   Scalar       *vals, *svals, *v;
395   MPI_Comm     comm = ((PetscObject)viewer)->comm;
396 
397   MPI_Comm_size(comm,&size);
398   if (size > 1) SETERRQ(1,"MatLoad_SeqDense: view must have one processor");
399   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
400   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
401   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqDense:Not matrix object");
402   M = header[1]; N = header[2]; nz = header[3];
403 
404   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
405     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
406     B = *A;
407     a = (Mat_SeqDense *) B->data;
408 
409     /* read in nonzero values */
410     ierr = PetscBinaryRead(fd,a->v,M*N,BINARY_SCALAR); CHKERRQ(ierr);
411 
412     ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
413     ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
414     ierr = MatTranspose(B,PETSC_NULL); CHKERRQ(ierr);
415   } else {
416     /* read row lengths */
417     rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
418     ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
419 
420     /* create our matrix */
421     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
422     B = *A;
423     a = (Mat_SeqDense *) B->data;
424     v = a->v;
425 
426     /* read column indices and nonzeros */
427     cols = scols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(cols);
428     ierr = PetscBinaryRead(fd,cols,nz,BINARY_INT); CHKERRQ(ierr);
429     vals = svals = (Scalar *) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals);
430     ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr);
431 
432     /* insert into matrix */
433     for ( i=0; i<M; i++ ) {
434       for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j];
435       svals += rowlengths[i]; scols += rowlengths[i];
436     }
437     PetscFree(vals); PetscFree(cols); PetscFree(rowlengths);
438 
439     ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
440     ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
441   }
442   return 0;
443 }
444 
445 #include "pinclude/pviewer.h"
446 #include "sys.h"
447 
448 static int MatView_SeqDense_ASCII(Mat A,Viewer viewer)
449 {
450   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
451   int          ierr, i, j, format;
452   FILE         *fd;
453   char         *outputname;
454   Scalar       *v;
455 
456   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
457   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
458   ierr = ViewerGetFormat(viewer,&format);
459   if (format == ASCII_FORMAT_INFO) {
460     ;  /* do nothing for now */
461   }
462   /* We retain dense format as default; use ASCII_FORMAT_IMPL to get
463      sparse format output ... maybe should switch this? */
464   else if (format == ASCII_FORMAT_IMPL) {
465     for ( i=0; i<a->m; i++ ) {
466       fprintf(fd,"row %d:",i);
467       for ( j=0; j<a->n; j++ ) {
468 #if defined(PETSC_COMPLEX)
469         if (real(*v) != 0.0 && imag(*v) != 0.0)
470           fprintf(fd," %d %g + %g i",j,real(*v),imag(*v));
471         else if (real(*v)) fprintf(fd," %d %g ",j,real(*v));
472         v += a->m;
473 #else
474         if (*v) fprintf(fd," %d %g ",j,*v);
475         v += a->m;
476 #endif
477       }
478       fprintf(fd,"\n");
479     }
480   }
481   else {
482 #if defined(PETSC_COMPLEX)
483     int allreal = 1;
484     /* determine if matrix has all real values */
485     v = a->v;
486     for ( i=0; i<a->m*a->n; i++ ) {
487       if (imag(v[i])) { allreal = 0; break ;}
488     }
489 #endif
490     for ( i=0; i<a->m; i++ ) {
491       v = a->v + i;
492       for ( j=0; j<a->n; j++ ) {
493 #if defined(PETSC_COMPLEX)
494         if (allreal) {
495           fprintf(fd,"%6.4e ",real(*v)); v += a->m;
496         } else {
497           fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m;
498         }
499 #else
500         fprintf(fd,"%6.4e ",*v); v += a->m;
501 #endif
502       }
503       fprintf(fd,"\n");
504     }
505   }
506   fflush(fd);
507   return 0;
508 }
509 
510 static int MatView_SeqDense_Binary(Mat A,Viewer viewer)
511 {
512   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
513   int          ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n;
514   int          format;
515   Scalar       *v, *anonz,*vals;
516 
517   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
518 
519   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
520   if (format == BINARY_FORMAT_NATIVE) {
521     /* store the matrix as a dense matrix */
522     col_lens = (int *) PetscMalloc( 4*sizeof(int) ); CHKPTRQ(col_lens);
523     col_lens[0] = MAT_COOKIE;
524     col_lens[1] = m;
525     col_lens[2] = n;
526     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
527     ierr = PetscBinaryWrite(fd,col_lens,4,BINARY_INT,1); CHKERRQ(ierr);
528     PetscFree(col_lens);
529 
530     /* write out matrix, by rows */
531     vals = (Scalar *) PetscMalloc(m*n*sizeof(Scalar)); CHKPTRQ(vals);
532     v    = a->v;
533     for ( i=0; i<m; i++ ) {
534       for ( j=0; j<n; j++ ) {
535         vals[i + j*m] = *v++;
536       }
537     }
538     ierr = PetscBinaryWrite(fd,vals,n*m,BINARY_SCALAR,0); CHKERRQ(ierr);
539     PetscFree(vals);
540   } else {
541     col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens);
542     col_lens[0] = MAT_COOKIE;
543     col_lens[1] = m;
544     col_lens[2] = n;
545     col_lens[3] = nz;
546 
547     /* store lengths of each row and write (including header) to file */
548     for ( i=0; i<m; i++ ) col_lens[4+i] = n;
549     ierr = PetscBinaryWrite(fd,col_lens,4+m,BINARY_INT,1); CHKERRQ(ierr);
550 
551     /* Possibly should write in smaller increments, not whole matrix at once? */
552     /* store column indices (zero start index) */
553     ict = 0;
554     for ( i=0; i<m; i++ ) {
555       for ( j=0; j<n; j++ ) col_lens[ict++] = j;
556     }
557     ierr = PetscBinaryWrite(fd,col_lens,nz,BINARY_INT,0); CHKERRQ(ierr);
558     PetscFree(col_lens);
559 
560     /* store nonzero values */
561     anonz = (Scalar *) PetscMalloc((nz)*sizeof(Scalar)); CHKPTRQ(anonz);
562     ict = 0;
563     for ( i=0; i<m; i++ ) {
564       v = a->v + i;
565       for ( j=0; j<n; j++ ) {
566         anonz[ict++] = *v; v += a->m;
567       }
568     }
569     ierr = PetscBinaryWrite(fd,anonz,nz,BINARY_SCALAR,0); CHKERRQ(ierr);
570     PetscFree(anonz);
571   }
572   return 0;
573 }
574 
575 static int MatView_SeqDense(PetscObject obj,Viewer viewer)
576 {
577   Mat          A = (Mat) obj;
578   Mat_SeqDense *a = (Mat_SeqDense*) A->data;
579   ViewerType   vtype;
580   int          ierr;
581 
582   if (!viewer) {
583     viewer = STDOUT_VIEWER_SELF;
584   }
585   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
586 
587   if (vtype == MATLAB_VIEWER) {
588     return ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v);
589   }
590   else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
591     return MatView_SeqDense_ASCII(A,viewer);
592   }
593   else if (vtype == BINARY_FILE_VIEWER) {
594     return MatView_SeqDense_Binary(A,viewer);
595   }
596   return 0;
597 }
598 
599 static int MatDestroy_SeqDense(PetscObject obj)
600 {
601   Mat          mat = (Mat) obj;
602   Mat_SeqDense *l = (Mat_SeqDense *) mat->data;
603 #if defined(PETSC_LOG)
604   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
605 #endif
606   if (l->pivots) PetscFree(l->pivots);
607   if (!l->user_alloc) PetscFree(l->v);
608   PetscFree(l);
609   PLogObjectDestroy(mat);
610   PetscHeaderDestroy(mat);
611   return 0;
612 }
613 
614 static int MatTranspose_SeqDense(Mat A,Mat *matout)
615 {
616   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
617   int          k, j, m, n;
618   Scalar       *v, tmp;
619 
620   v = mat->v; m = mat->m; n = mat->n;
621   if (matout == PETSC_NULL) { /* in place transpose */
622     if (m != n) SETERRQ(1,"MatTranspose_SeqDense:Supports square matrix only in-place");
623     for ( j=0; j<m; j++ ) {
624       for ( k=0; k<j; k++ ) {
625         tmp = v[j + k*n];
626         v[j + k*n] = v[k + j*n];
627         v[k + j*n] = tmp;
628       }
629     }
630   }
631   else { /* out-of-place transpose */
632     int          ierr;
633     Mat          tmat;
634     Mat_SeqDense *tmatd;
635     Scalar       *v2;
636     ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr);
637     tmatd = (Mat_SeqDense *) tmat->data;
638     v = mat->v; v2 = tmatd->v;
639     for ( j=0; j<n; j++ ) {
640       for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m];
641     }
642     ierr = MatAssemblyBegin(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
643     ierr = MatAssemblyEnd(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
644     *matout = tmat;
645   }
646   return 0;
647 }
648 
649 static int MatEqual_SeqDense(Mat A1,Mat A2, PetscTruth *flg)
650 {
651   Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data;
652   Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data;
653   int          i;
654   Scalar       *v1 = mat1->v, *v2 = mat2->v;
655 
656   if (A2->type != MATSEQDENSE) SETERRQ(1,"MatEqual_SeqDense:Both matrices should be of type  MATSEQDENSE");
657   if (mat1->m != mat2->m) {*flg = PETSC_FALSE; return 0;}
658   if (mat1->n != mat2->n) {*flg = PETSC_FALSE; return 0;}
659   for ( i=0; i<mat1->m*mat1->n; i++ ) {
660     if (*v1 != *v2) {*flg = PETSC_FALSE; return 0;}
661     v1++; v2++;
662   }
663   *flg = PETSC_TRUE;
664   return 0;
665 }
666 
667 static int MatGetDiagonal_SeqDense(Mat A,Vec v)
668 {
669   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
670   int          i, n;
671   Scalar       *x;
672   VecGetArray(v,&x); VecGetSize(v,&n);
673   if (n != mat->m) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec");
674   for ( i=0; i<mat->m; i++ ) {
675     x[i] = mat->v[i*mat->m + i];
676   }
677   return 0;
678 }
679 
680 static int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
681 {
682   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
683   Scalar       *l,*r,x,*v;
684   int          i,j,m = mat->m, n = mat->n;
685 
686   if (ll) {
687     VecGetArray(ll,&l); VecGetSize(ll,&m);
688     if (m != mat->m) SETERRQ(1,"MatDiagonalScale_SeqDense:Left scaling vec wrong size");
689     PLogFlops(n*m);
690     for ( i=0; i<m; i++ ) {
691       x = l[i];
692       v = mat->v + i;
693       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
694     }
695   }
696   if (rr) {
697     VecGetArray(rr,&r); VecGetSize(rr,&n);
698     if (n != mat->n) SETERRQ(1,"MatDiagonalScale_SeqDense:Right scaling vec wrong size");
699     PLogFlops(n*m);
700     for ( i=0; i<n; i++ ) {
701       x = r[i];
702       v = mat->v + i*m;
703       for ( j=0; j<m; j++ ) { (*v++) *= x;}
704     }
705   }
706   return 0;
707 }
708 
709 static int MatNorm_SeqDense(Mat A,NormType type,double *norm)
710 {
711   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
712   Scalar       *v = mat->v;
713   double       sum = 0.0;
714   int          i, j;
715 
716   if (type == NORM_FROBENIUS) {
717     for (i=0; i<mat->n*mat->m; i++ ) {
718 #if defined(PETSC_COMPLEX)
719       sum += real(conj(*v)*(*v)); v++;
720 #else
721       sum += (*v)*(*v); v++;
722 #endif
723     }
724     *norm = sqrt(sum);
725     PLogFlops(2*mat->n*mat->m);
726   }
727   else if (type == NORM_1) {
728     *norm = 0.0;
729     for ( j=0; j<mat->n; j++ ) {
730       sum = 0.0;
731       for ( i=0; i<mat->m; i++ ) {
732         sum += PetscAbsScalar(*v);  v++;
733       }
734       if (sum > *norm) *norm = sum;
735     }
736     PLogFlops(mat->n*mat->m);
737   }
738   else if (type == NORM_INFINITY) {
739     *norm = 0.0;
740     for ( j=0; j<mat->m; j++ ) {
741       v = mat->v + j;
742       sum = 0.0;
743       for ( i=0; i<mat->n; i++ ) {
744         sum += PetscAbsScalar(*v); v += mat->m;
745       }
746       if (sum > *norm) *norm = sum;
747     }
748     PLogFlops(mat->n*mat->m);
749   }
750   else {
751     SETERRQ(1,"MatNorm_SeqDense:No two norm");
752   }
753   return 0;
754 }
755 
756 static int MatSetOption_SeqDense(Mat A,MatOption op)
757 {
758   Mat_SeqDense *aij = (Mat_SeqDense *) A->data;
759 
760   if (op == ROW_ORIENTED)            aij->roworiented = 1;
761   else if (op == COLUMN_ORIENTED)    aij->roworiented = 0;
762   else if (op == ROWS_SORTED ||
763            op == COLUMNS_SORTED ||
764            op == SYMMETRIC_MATRIX ||
765            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
766            op == NO_NEW_NONZERO_LOCATIONS ||
767            op == YES_NEW_NONZERO_LOCATIONS ||
768            op == NO_NEW_DIAGONALS ||
769            op == YES_NEW_DIAGONALS)
770     PLogInfo(A,"Info:MatSetOption_SeqDense:Option ignored\n");
771   else
772     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqDense:unknown option");}
773   return 0;
774 }
775 
776 static int MatZeroEntries_SeqDense(Mat A)
777 {
778   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
779   PetscMemzero(l->v,l->m*l->n*sizeof(Scalar));
780   return 0;
781 }
782 
783 static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
784 {
785   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
786   int          n = l->n, i, j,ierr,N, *rows;
787   Scalar       *slot;
788 
789   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
790   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
791   for ( i=0; i<N; i++ ) {
792     slot = l->v + rows[i];
793     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
794   }
795   if (diag) {
796     for ( i=0; i<N; i++ ) {
797       slot = l->v + (n+1)*rows[i];
798       *slot = *diag;
799     }
800   }
801   ISRestoreIndices(is,&rows);
802   return 0;
803 }
804 
805 static int MatGetSize_SeqDense(Mat A,int *m,int *n)
806 {
807   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
808   *m = mat->m; *n = mat->n;
809   return 0;
810 }
811 
812 static int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
813 {
814   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
815   *m = 0; *n = mat->m;
816   return 0;
817 }
818 
819 static int MatGetArray_SeqDense(Mat A,Scalar **array)
820 {
821   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
822   *array = mat->v;
823   return 0;
824 }
825 
826 static int MatRestoreArray_SeqDense(Mat A,Scalar **array)
827 {
828   return 0;
829 }
830 
831 static int MatGetSubMatrixInPlace_SeqDense(Mat A,IS isrow,IS iscol)
832 {
833   SETERRQ(1,"MatGetSubMatrixInPlace_SeqDense:not done");
834 }
835 
836 static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,
837                                     Mat *submat)
838 {
839   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
840   int          nznew, *smap, i, j, ierr, oldcols = mat->n;
841   int          *irow, *icol, nrows, ncols, *cwork;
842   Scalar       *vwork, *val;
843   Mat          newmat;
844 
845   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
846   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
847   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
848   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
849 
850   smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap);
851   cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork);
852   vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
853   PetscMemzero((char*)smap,oldcols*sizeof(int));
854   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
855 
856   /* Create and fill new matrix */
857   ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr);
858   for (i=0; i<nrows; i++) {
859     nznew = 0;
860     val   = mat->v + irow[i];
861     for (j=0; j<oldcols; j++) {
862       if (smap[j]) {
863         cwork[nznew]   = smap[j] - 1;
864         vwork[nznew++] = val[j * mat->m];
865       }
866     }
867     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES);
868            CHKERRQ(ierr);
869   }
870   ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
871   ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
872 
873   /* Free work space */
874   PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
875   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
876   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
877   *submat = newmat;
878   return 0;
879 }
880 
881 static int MatCopy_SeqDense(Mat A, Mat B)
882 {
883   Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data;
884   if (B->type != MATSEQDENSE) return MatCopy_Basic(A,B);
885   if (a->m != b->m || a->n != b->n) SETERRQ(1,"MatCopy_SeqDense:size(B) != size(A)");
886   PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar));
887   return 0;
888 }
889 
890 /* -------------------------------------------------------------------*/
891 static struct _MatOps MatOps = {MatSetValues_SeqDense,
892        MatGetRow_SeqDense, MatRestoreRow_SeqDense,
893        MatMult_SeqDense, MatMultAdd_SeqDense,
894        MatMultTrans_SeqDense, MatMultTransAdd_SeqDense,
895        MatSolve_SeqDense,MatSolveAdd_SeqDense,
896        MatSolveTrans_SeqDense,MatSolveTransAdd_SeqDense,
897        MatLUFactor_SeqDense,MatCholeskyFactor_SeqDense,
898        MatRelax_SeqDense,
899        MatTranspose_SeqDense,
900        MatGetInfo_SeqDense,MatEqual_SeqDense,
901        MatGetDiagonal_SeqDense,MatDiagonalScale_SeqDense,MatNorm_SeqDense,
902        0,0,
903        0, MatSetOption_SeqDense,MatZeroEntries_SeqDense,MatZeroRows_SeqDense,0,
904        MatLUFactorSymbolic_SeqDense,MatLUFactorNumeric_SeqDense,
905        MatCholeskyFactorSymbolic_SeqDense,MatCholeskyFactorNumeric_SeqDense,
906        MatGetSize_SeqDense,MatGetSize_SeqDense,MatGetOwnershipRange_SeqDense,
907        0,0, MatGetArray_SeqDense, MatRestoreArray_SeqDense,0,
908        MatGetSubMatrix_SeqDense,MatGetSubMatrixInPlace_SeqDense,
909        MatConvertSameType_SeqDense,0,0,0,0,
910        MatAXPY_SeqDense,0,0,
911        MatGetValues_SeqDense,
912        MatCopy_SeqDense,0,MatScale_SeqDense};
913 
914 
915 /*@C
916    MatCreateSeqDense - Creates a sequential dense matrix that
917    is stored in column major order (the usual Fortran 77 manner). Many
918    of the matrix operations use the BLAS and LAPACK routines.
919 
920    Input Parameters:
921 .  comm - MPI communicator, set to MPI_COMM_SELF
922 .  m - number of rows
923 .  n - number of columns
924 .  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
925    to control all matrix memory allocation.
926 
927    Output Parameter:
928 .  newmat - the matrix
929 
930   Notes:
931   The data input variable is intended primarily for Fortran programmers
932   who wish to allocate their own matrix memory space.  Most users should
933   set data=PETSC_NULL.
934 
935 .keywords: dense, matrix, LAPACK, BLAS
936 
937 .seealso: MatCreate(), MatSetValues()
938 @*/
939 int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *newmat)
940 {
941   Mat          mat;
942   Mat_SeqDense *l;
943   int          ierr,flg;
944 
945   *newmat        = 0;
946 
947   PetscHeaderCreate(mat,_Mat,MAT_COOKIE,MATSEQDENSE,comm);
948   PLogObjectCreate(mat);
949   l               = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(l);
950   PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps));
951   mat->destroy    = MatDestroy_SeqDense;
952   mat->view       = MatView_SeqDense;
953   mat->factor     = 0;
954   PLogObjectMemory(mat,sizeof(struct _Mat));
955   mat->data       = (void *) l;
956 
957   l->m            = m;
958   l->n            = n;
959   l->pivots       = 0;
960   l->roworiented  = 1;
961   if (data == PETSC_NULL) {
962     l->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(l->v);
963     PetscMemzero(l->v,m*n*sizeof(Scalar));
964     l->user_alloc = 0;
965     PLogObjectMemory(mat,n*m*sizeof(Scalar));
966   }
967   else { /* user-allocated storage */
968     l->v = data;
969     l->user_alloc = 1;
970   }
971 
972   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
973   if (flg) {
974     ierr = MatPrintHelp(mat); CHKERRQ(ierr);
975   }
976   *newmat = mat;
977   return 0;
978 }
979 
980 int MatCreate_SeqDense(Mat A,Mat *newmat)
981 {
982   Mat_SeqDense *m = (Mat_SeqDense *) A->data;
983   return MatCreateSeqDense(A->comm,m->m,m->n,PETSC_NULL,newmat);
984 }
985