xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 4e627d3a034d8dbe6c0b444d105a274c665112f3)
1 #ifndef lint
2 static char vcid[] = "$Id: dense.c,v 1.102 1996/04/07 22:44:50 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 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 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 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 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); VecGetArray(zz,&z);
280   if (zz != yy) PetscMemcpy(y,z,mat->n*sizeof(Scalar));
281   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One);
282   PLogFlops(2*mat->m*mat->n);
283   return 0;
284 }
285 
286 /* -----------------------------------------------------------------*/
287 static int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
288 {
289   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
290   Scalar       *v;
291   int          i;
292 
293   *ncols = mat->n;
294   if (cols) {
295     *cols = (int *) PetscMalloc(mat->n*sizeof(int)); CHKPTRQ(*cols);
296     for ( i=0; i<mat->n; i++ ) (*cols)[i] = i;
297   }
298   if (vals) {
299     *vals = (Scalar *) PetscMalloc(mat->n*sizeof(Scalar)); CHKPTRQ(*vals);
300     v = mat->v + row;
301     for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;}
302   }
303   return 0;
304 }
305 static int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
306 {
307   if (cols) { PetscFree(*cols); }
308   if (vals) { PetscFree(*vals); }
309   return 0;
310 }
311 /* ----------------------------------------------------------------*/
312 static int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n,
313                                     int *indexn,Scalar *v,InsertMode addv)
314 {
315   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
316   int          i,j;
317 
318   if (!mat->roworiented) {
319     if (addv == INSERT_VALUES) {
320       for ( j=0; j<n; j++ ) {
321         for ( i=0; i<m; i++ ) {
322           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
323         }
324       }
325     }
326     else {
327       for ( j=0; j<n; j++ ) {
328         for ( i=0; i<m; i++ ) {
329           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
330         }
331       }
332     }
333   }
334   else {
335     if (addv == INSERT_VALUES) {
336       for ( i=0; i<m; i++ ) {
337         for ( j=0; j<n; j++ ) {
338           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
339         }
340       }
341     }
342     else {
343       for ( i=0; i<m; i++ ) {
344         for ( j=0; j<n; j++ ) {
345           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
346         }
347       }
348     }
349   }
350   return 0;
351 }
352 
353 static int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v)
354 {
355   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
356   int          i, j;
357   Scalar       *vpt = v;
358 
359   /* row-oriented output */
360   for ( i=0; i<m; i++ ) {
361     for ( j=0; j<n; j++ ) {
362       *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]];
363     }
364   }
365   return 0;
366 }
367 
368 /* -----------------------------------------------------------------*/
369 static int MatConvertSameType_SeqDense(Mat A,Mat *newmat,int cpvalues)
370 {
371   Mat_SeqDense *mat = (Mat_SeqDense *) A->data, *l;
372   int          ierr;
373   Mat          newi;
374 
375   ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,PETSC_NULL,&newi); CHKERRQ(ierr);
376   l = (Mat_SeqDense *) newi->data;
377   if (cpvalues == COPY_VALUES) {
378     PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
379   }
380   newi->assembled = PETSC_TRUE;
381   *newmat = newi;
382   return 0;
383 }
384 
385 #include "sys.h"
386 
387 int MatLoad_SeqDense(Viewer viewer,MatType type,Mat *A)
388 {
389   Mat_SeqDense *a;
390   Mat          B;
391   int          *scols, i, j, nz, ierr, fd, header[4], size;
392   int          *rowlengths = 0, M, N, *cols;
393   Scalar       *vals, *svals, *v;
394   MPI_Comm     comm = ((PetscObject)viewer)->comm;
395 
396   MPI_Comm_size(comm,&size);
397   if (size > 1) SETERRQ(1,"MatLoad_SeqDense: view must have one processor");
398   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
399   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
400   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqDense:Not matrix object");
401   M = header[1]; N = header[2]; nz = header[3];
402 
403   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
404     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
405     B = *A;
406     a = (Mat_SeqDense *) B->data;
407 
408     /* read in nonzero values */
409     ierr = PetscBinaryRead(fd,a->v,M*N,BINARY_SCALAR); CHKERRQ(ierr);
410 
411     ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
412     ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
413     ierr = MatTranspose(B,PETSC_NULL); CHKERRQ(ierr);
414   } else {
415     /* read row lengths */
416     rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
417     ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
418 
419     /* create our matrix */
420     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
421     B = *A;
422     a = (Mat_SeqDense *) B->data;
423     v = a->v;
424 
425     /* read column indices and nonzeros */
426     cols = scols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(cols);
427     ierr = PetscBinaryRead(fd,cols,nz,BINARY_INT); CHKERRQ(ierr);
428     vals = svals = (Scalar *) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals);
429     ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr);
430 
431     /* insert into matrix */
432     for ( i=0; i<M; i++ ) {
433       for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j];
434       svals += rowlengths[i]; scols += rowlengths[i];
435     }
436     PetscFree(vals); PetscFree(cols); PetscFree(rowlengths);
437 
438     ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
439     ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
440   }
441   return 0;
442 }
443 
444 #include "pinclude/pviewer.h"
445 #include "sys.h"
446 
447 static int MatView_SeqDense_ASCII(Mat A,Viewer viewer)
448 {
449   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
450   int          ierr, i, j, format;
451   FILE         *fd;
452   char         *outputname;
453   Scalar       *v;
454 
455   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
456   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
457   ierr = ViewerGetFormat(viewer,&format);
458   if (format == ASCII_FORMAT_INFO || format == ASCII_FORMAT_INFO_DETAILED) {
459     return 0;  /* do nothing for now */
460   }
461   else if (format == ASCII_FORMAT_COMMON) {
462     for ( i=0; i<a->m; i++ ) {
463       v = a->v + i;
464       fprintf(fd,"row %d:",i);
465       for ( j=0; j<a->n; j++ ) {
466 #if defined(PETSC_COMPLEX)
467         if (real(*v) != 0.0 && imag(*v) != 0.0)
468           fprintf(fd," %d %g + %g i",j,real(*v),imag(*v));
469         else if (real(*v)) fprintf(fd," %d %g ",j,real(*v));
470         v += a->m;
471 #else
472         if (*v) fprintf(fd," %d %g ",j,*v);
473         v += a->m;
474 #endif
475       }
476       fprintf(fd,"\n");
477     }
478   }
479   else {
480 #if defined(PETSC_COMPLEX)
481     int allreal = 1;
482     /* determine if matrix has all real values */
483     v = a->v;
484     for ( i=0; i<a->m*a->n; i++ ) {
485       if (imag(v[i])) { allreal = 0; break ;}
486     }
487 #endif
488     for ( i=0; i<a->m; i++ ) {
489       v = a->v + i;
490       for ( j=0; j<a->n; j++ ) {
491 #if defined(PETSC_COMPLEX)
492         if (allreal) {
493           fprintf(fd,"%6.4e ",real(*v)); v += a->m;
494         } else {
495           fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m;
496         }
497 #else
498         fprintf(fd,"%6.4e ",*v); v += a->m;
499 #endif
500       }
501       fprintf(fd,"\n");
502     }
503   }
504   fflush(fd);
505   return 0;
506 }
507 
508 static int MatView_SeqDense_Binary(Mat A,Viewer viewer)
509 {
510   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
511   int          ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n;
512   int          format;
513   Scalar       *v, *anonz,*vals;
514 
515   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
516 
517   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
518   if (format == BINARY_FORMAT_NATIVE) {
519     /* store the matrix as a dense matrix */
520     col_lens = (int *) PetscMalloc( 4*sizeof(int) ); CHKPTRQ(col_lens);
521     col_lens[0] = MAT_COOKIE;
522     col_lens[1] = m;
523     col_lens[2] = n;
524     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
525     ierr = PetscBinaryWrite(fd,col_lens,4,BINARY_INT,1); CHKERRQ(ierr);
526     PetscFree(col_lens);
527 
528     /* write out matrix, by rows */
529     vals = (Scalar *) PetscMalloc(m*n*sizeof(Scalar)); CHKPTRQ(vals);
530     v    = a->v;
531     for ( i=0; i<m; i++ ) {
532       for ( j=0; j<n; j++ ) {
533         vals[i + j*m] = *v++;
534       }
535     }
536     ierr = PetscBinaryWrite(fd,vals,n*m,BINARY_SCALAR,0); CHKERRQ(ierr);
537     PetscFree(vals);
538   } else {
539     col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens);
540     col_lens[0] = MAT_COOKIE;
541     col_lens[1] = m;
542     col_lens[2] = n;
543     col_lens[3] = nz;
544 
545     /* store lengths of each row and write (including header) to file */
546     for ( i=0; i<m; i++ ) col_lens[4+i] = n;
547     ierr = PetscBinaryWrite(fd,col_lens,4+m,BINARY_INT,1); CHKERRQ(ierr);
548 
549     /* Possibly should write in smaller increments, not whole matrix at once? */
550     /* store column indices (zero start index) */
551     ict = 0;
552     for ( i=0; i<m; i++ ) {
553       for ( j=0; j<n; j++ ) col_lens[ict++] = j;
554     }
555     ierr = PetscBinaryWrite(fd,col_lens,nz,BINARY_INT,0); CHKERRQ(ierr);
556     PetscFree(col_lens);
557 
558     /* store nonzero values */
559     anonz = (Scalar *) PetscMalloc((nz)*sizeof(Scalar)); CHKPTRQ(anonz);
560     ict = 0;
561     for ( i=0; i<m; i++ ) {
562       v = a->v + i;
563       for ( j=0; j<n; j++ ) {
564         anonz[ict++] = *v; v += a->m;
565       }
566     }
567     ierr = PetscBinaryWrite(fd,anonz,nz,BINARY_SCALAR,0); CHKERRQ(ierr);
568     PetscFree(anonz);
569   }
570   return 0;
571 }
572 
573 static int MatView_SeqDense(PetscObject obj,Viewer viewer)
574 {
575   Mat          A = (Mat) obj;
576   Mat_SeqDense *a = (Mat_SeqDense*) A->data;
577   ViewerType   vtype;
578   int          ierr;
579 
580   if (!viewer) {
581     viewer = STDOUT_VIEWER_SELF;
582   }
583   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
584 
585   if (vtype == MATLAB_VIEWER) {
586     return ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v);
587   }
588   else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
589     return MatView_SeqDense_ASCII(A,viewer);
590   }
591   else if (vtype == BINARY_FILE_VIEWER) {
592     return MatView_SeqDense_Binary(A,viewer);
593   }
594   return 0;
595 }
596 
597 static int MatDestroy_SeqDense(PetscObject obj)
598 {
599   Mat          mat = (Mat) obj;
600   Mat_SeqDense *l = (Mat_SeqDense *) mat->data;
601 #if defined(PETSC_LOG)
602   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
603 #endif
604   if (l->pivots) PetscFree(l->pivots);
605   if (!l->user_alloc) PetscFree(l->v);
606   PetscFree(l);
607   PLogObjectDestroy(mat);
608   PetscHeaderDestroy(mat);
609   return 0;
610 }
611 
612 static int MatTranspose_SeqDense(Mat A,Mat *matout)
613 {
614   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
615   int          k, j, m, n;
616   Scalar       *v, tmp;
617 
618   v = mat->v; m = mat->m; n = mat->n;
619   if (matout == PETSC_NULL) { /* in place transpose */
620     if (m != n) SETERRQ(1,"MatTranspose_SeqDense:Supports square matrix only in-place");
621     for ( j=0; j<m; j++ ) {
622       for ( k=0; k<j; k++ ) {
623         tmp = v[j + k*n];
624         v[j + k*n] = v[k + j*n];
625         v[k + j*n] = tmp;
626       }
627     }
628   }
629   else { /* out-of-place transpose */
630     int          ierr;
631     Mat          tmat;
632     Mat_SeqDense *tmatd;
633     Scalar       *v2;
634     ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr);
635     tmatd = (Mat_SeqDense *) tmat->data;
636     v = mat->v; v2 = tmatd->v;
637     for ( j=0; j<n; j++ ) {
638       for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m];
639     }
640     ierr = MatAssemblyBegin(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
641     ierr = MatAssemblyEnd(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
642     *matout = tmat;
643   }
644   return 0;
645 }
646 
647 static int MatEqual_SeqDense(Mat A1,Mat A2, PetscTruth *flg)
648 {
649   Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data;
650   Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data;
651   int          i;
652   Scalar       *v1 = mat1->v, *v2 = mat2->v;
653 
654   if (A2->type != MATSEQDENSE) SETERRQ(1,"MatEqual_SeqDense:Both matrices should be of type  MATSEQDENSE");
655   if (mat1->m != mat2->m) {*flg = PETSC_FALSE; return 0;}
656   if (mat1->n != mat2->n) {*flg = PETSC_FALSE; return 0;}
657   for ( i=0; i<mat1->m*mat1->n; i++ ) {
658     if (*v1 != *v2) {*flg = PETSC_FALSE; return 0;}
659     v1++; v2++;
660   }
661   *flg = PETSC_TRUE;
662   return 0;
663 }
664 
665 static int MatGetDiagonal_SeqDense(Mat A,Vec v)
666 {
667   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
668   int          i, n, len;
669   Scalar       *x, zero = 0.0;
670 
671   VecSet(&zero,v);
672   VecGetArray(v,&x); VecGetSize(v,&n);
673   len = PetscMin(mat->m,mat->n);
674   if (n != mat->m) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec");
675   for ( i=0; i<len; i++ ) {
676     x[i] = mat->v[i*mat->m + i];
677   }
678   return 0;
679 }
680 
681 static int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
682 {
683   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
684   Scalar       *l,*r,x,*v;
685   int          i,j,m = mat->m, n = mat->n;
686 
687   if (ll) {
688     VecGetArray(ll,&l); VecGetSize(ll,&m);
689     if (m != mat->m) SETERRQ(1,"MatDiagonalScale_SeqDense:Left scaling vec wrong size");
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     PLogFlops(n*m);
696   }
697   if (rr) {
698     VecGetArray(rr,&r); VecGetSize(rr,&n);
699     if (n != mat->n) SETERRQ(1,"MatDiagonalScale_SeqDense:Right scaling vec wrong size");
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     PLogFlops(n*m);
706   }
707   return 0;
708 }
709 
710 static int MatNorm_SeqDense(Mat A,NormType type,double *norm)
711 {
712   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
713   Scalar       *v = mat->v;
714   double       sum = 0.0;
715   int          i, j;
716 
717   if (type == NORM_FROBENIUS) {
718     for (i=0; i<mat->n*mat->m; i++ ) {
719 #if defined(PETSC_COMPLEX)
720       sum += real(conj(*v)*(*v)); v++;
721 #else
722       sum += (*v)*(*v); v++;
723 #endif
724     }
725     *norm = sqrt(sum);
726     PLogFlops(2*mat->n*mat->m);
727   }
728   else if (type == NORM_1) {
729     *norm = 0.0;
730     for ( j=0; j<mat->n; j++ ) {
731       sum = 0.0;
732       for ( i=0; i<mat->m; i++ ) {
733         sum += PetscAbsScalar(*v);  v++;
734       }
735       if (sum > *norm) *norm = sum;
736     }
737     PLogFlops(mat->n*mat->m);
738   }
739   else if (type == NORM_INFINITY) {
740     *norm = 0.0;
741     for ( j=0; j<mat->m; j++ ) {
742       v = mat->v + j;
743       sum = 0.0;
744       for ( i=0; i<mat->n; i++ ) {
745         sum += PetscAbsScalar(*v); v += mat->m;
746       }
747       if (sum > *norm) *norm = sum;
748     }
749     PLogFlops(mat->n*mat->m);
750   }
751   else {
752     SETERRQ(1,"MatNorm_SeqDense:No two norm");
753   }
754   return 0;
755 }
756 
757 static int MatSetOption_SeqDense(Mat A,MatOption op)
758 {
759   Mat_SeqDense *aij = (Mat_SeqDense *) A->data;
760 
761   if (op == ROW_ORIENTED)            aij->roworiented = 1;
762   else if (op == COLUMN_ORIENTED)    aij->roworiented = 0;
763   else if (op == ROWS_SORTED ||
764            op == COLUMNS_SORTED ||
765            op == SYMMETRIC_MATRIX ||
766            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
767            op == NO_NEW_NONZERO_LOCATIONS ||
768            op == YES_NEW_NONZERO_LOCATIONS ||
769            op == NO_NEW_DIAGONALS ||
770            op == YES_NEW_DIAGONALS)
771     PLogInfo(A,"Info:MatSetOption_SeqDense:Option ignored\n");
772   else
773     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqDense:unknown option");}
774   return 0;
775 }
776 
777 static int MatZeroEntries_SeqDense(Mat A)
778 {
779   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
780   PetscMemzero(l->v,l->m*l->n*sizeof(Scalar));
781   return 0;
782 }
783 
784 static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
785 {
786   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
787   int          n = l->n, i, j,ierr,N, *rows;
788   Scalar       *slot;
789 
790   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
791   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
792   for ( i=0; i<N; i++ ) {
793     slot = l->v + rows[i];
794     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
795   }
796   if (diag) {
797     for ( i=0; i<N; i++ ) {
798       slot = l->v + (n+1)*rows[i];
799       *slot = *diag;
800     }
801   }
802   ISRestoreIndices(is,&rows);
803   return 0;
804 }
805 
806 static int MatGetSize_SeqDense(Mat A,int *m,int *n)
807 {
808   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
809   *m = mat->m; *n = mat->n;
810   return 0;
811 }
812 
813 static int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
814 {
815   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
816   *m = 0; *n = mat->m;
817   return 0;
818 }
819 
820 static int MatGetArray_SeqDense(Mat A,Scalar **array)
821 {
822   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
823   *array = mat->v;
824   return 0;
825 }
826 
827 static int MatRestoreArray_SeqDense(Mat A,Scalar **array)
828 {
829   return 0;
830 }
831 
832 static int MatGetSubMatrixInPlace_SeqDense(Mat A,IS isrow,IS iscol)
833 {
834   SETERRQ(1,"MatGetSubMatrixInPlace_SeqDense:not done");
835 }
836 
837 static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,
838                                     Mat *submat)
839 {
840   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
841   int          nznew, *smap, i, j, ierr, oldcols = mat->n;
842   int          *irow, *icol, nrows, ncols, *cwork;
843   Scalar       *vwork, *val;
844   Mat          newmat;
845 
846   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
847   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
848   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
849   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
850 
851   smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap);
852   cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork);
853   vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
854   PetscMemzero((char*)smap,oldcols*sizeof(int));
855   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
856 
857   /* Create and fill new matrix */
858   ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr);
859   for (i=0; i<nrows; i++) {
860     nznew = 0;
861     val   = mat->v + irow[i];
862     for (j=0; j<oldcols; j++) {
863       if (smap[j]) {
864         cwork[nznew]   = smap[j] - 1;
865         vwork[nznew++] = val[j * mat->m];
866       }
867     }
868     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES);
869            CHKERRQ(ierr);
870   }
871   ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
872   ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
873 
874   /* Free work space */
875   PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
876   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
877   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
878   *submat = newmat;
879   return 0;
880 }
881 
882 static int MatCopy_SeqDense(Mat A, Mat B)
883 {
884   Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data;
885   if (B->type != MATSEQDENSE) return MatCopy_Basic(A,B);
886   if (a->m != b->m || a->n != b->n) SETERRQ(1,"MatCopy_SeqDense:size(B) != size(A)");
887   PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar));
888   return 0;
889 }
890 
891 /* -------------------------------------------------------------------*/
892 static struct _MatOps MatOps = {MatSetValues_SeqDense,
893        MatGetRow_SeqDense, MatRestoreRow_SeqDense,
894        MatMult_SeqDense, MatMultAdd_SeqDense,
895        MatMultTrans_SeqDense, MatMultTransAdd_SeqDense,
896        MatSolve_SeqDense,MatSolveAdd_SeqDense,
897        MatSolveTrans_SeqDense,MatSolveTransAdd_SeqDense,
898        MatLUFactor_SeqDense,MatCholeskyFactor_SeqDense,
899        MatRelax_SeqDense,
900        MatTranspose_SeqDense,
901        MatGetInfo_SeqDense,MatEqual_SeqDense,
902        MatGetDiagonal_SeqDense,MatDiagonalScale_SeqDense,MatNorm_SeqDense,
903        0,0,
904        0, MatSetOption_SeqDense,MatZeroEntries_SeqDense,MatZeroRows_SeqDense,0,
905        MatLUFactorSymbolic_SeqDense,MatLUFactorNumeric_SeqDense,
906        MatCholeskyFactorSymbolic_SeqDense,MatCholeskyFactorNumeric_SeqDense,
907        MatGetSize_SeqDense,MatGetSize_SeqDense,MatGetOwnershipRange_SeqDense,
908        0,0, MatGetArray_SeqDense, MatRestoreArray_SeqDense,0,
909        MatGetSubMatrix_SeqDense,MatGetSubMatrixInPlace_SeqDense,
910        MatConvertSameType_SeqDense,0,0,0,0,
911        MatAXPY_SeqDense,0,0,
912        MatGetValues_SeqDense,
913        MatCopy_SeqDense,0,MatScale_SeqDense};
914 
915 
916 /*@C
917    MatCreateSeqDense - Creates a sequential dense matrix that
918    is stored in column major order (the usual Fortran 77 manner). Many
919    of the matrix operations use the BLAS and LAPACK routines.
920 
921    Input Parameters:
922 .  comm - MPI communicator, set to MPI_COMM_SELF
923 .  m - number of rows
924 .  n - number of columns
925 .  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
926    to control all matrix memory allocation.
927 
928    Output Parameter:
929 .  A - the matrix
930 
931   Notes:
932   The data input variable is intended primarily for Fortran programmers
933   who wish to allocate their own matrix memory space.  Most users should
934   set data=PETSC_NULL.
935 
936 .keywords: dense, matrix, LAPACK, BLAS
937 
938 .seealso: MatCreate(), MatSetValues()
939 @*/
940 int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A)
941 {
942   Mat          B;
943   Mat_SeqDense *b;
944   int          ierr,flg;
945 
946   *A            = 0;
947   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQDENSE,comm);
948   PLogObjectCreate(B);
949   b             = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(b);
950   PetscMemzero(b,sizeof(Mat_SeqDense));
951   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
952   B->destroy    = MatDestroy_SeqDense;
953   B->view       = MatView_SeqDense;
954   B->factor     = 0;
955   PLogObjectMemory(B,sizeof(struct _Mat));
956   B->data       = (void *) b;
957 
958   b->m = m;  B->m = m; B->M = m;
959   b->n = n;  B->n = n; B->N = n;
960   b->pivots       = 0;
961   b->roworiented  = 1;
962   if (data == PETSC_NULL) {
963     b->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(b->v);
964     PetscMemzero(b->v,m*n*sizeof(Scalar));
965     b->user_alloc = 0;
966     PLogObjectMemory(B,n*m*sizeof(Scalar));
967   }
968   else { /* user-allocated storage */
969     b->v = data;
970     b->user_alloc = 1;
971   }
972   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
973   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr);}
974   *A = B;
975   return 0;
976 }
977 
978 int MatCreate_SeqDense(Mat A,Mat *newmat)
979 {
980   Mat_SeqDense *m = (Mat_SeqDense *) A->data;
981   return MatCreateSeqDense(A->comm,m->m,m->n,PETSC_NULL,newmat);
982 }
983