xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 70f55243aafb320636e2a54ff30cab5d1e8d3d7b)
1 #ifndef lint
2 static char vcid[] = "$Id: dense.c,v 1.107 1996/08/06 04:02:12 bsmith Exp bsmith $";
3 #endif
4 /*
5      Defines the basic matrix operations for sequential dense.
6 */
7 
8 #include "src/mat/impls/dense/seq/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<0) SETERRQ(1,"MatLUFactor_SeqDense:Bad argument to LU factorization");
58   if (info>0) SETERRQ(1,"MatLUFactor_SeqDense:Bad LU factorization");
59   A->factor = FACTOR_LU;
60   PLogFlops((2*mat->n*mat->n*mat->n)/3);
61   return 0;
62 }
63 static int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,double f,Mat *fact)
64 {
65   return MatConvert(A,MATSAME,fact);
66 }
67 static int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact)
68 {
69   return MatLUFactor(*fact,0,0,1.0);
70 }
71 static int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact)
72 {
73   return MatConvert(A,MATSAME,fact);
74 }
75 static int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact)
76 {
77   return MatCholeskyFactor(*fact,0,1.0);
78 }
79 static int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f)
80 {
81   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
82   int           info;
83 
84   if (mat->pivots) {
85     PetscFree(mat->pivots);
86     PLogObjectMemory(A,-mat->m*sizeof(int));
87     mat->pivots = 0;
88   }
89   LApotrf_("L",&mat->n,mat->v,&mat->m,&info);
90   if (info) SETERRQ(1,"MatCholeskyFactor_SeqDense:Bad factorization");
91   A->factor = FACTOR_CHOLESKY;
92   PLogFlops((mat->n*mat->n*mat->n)/3);
93   return 0;
94 }
95 
96 static int MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
97 {
98   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
99   int          one = 1, info, ierr;
100   Scalar       *x, *y;
101 
102   ierr = VecGetArray(xx,&x); CHKERRQ(ierr); ierr = VecGetArray(yy,&y); CHKERRQ(ierr);
103   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
104   if (A->factor == FACTOR_LU) {
105     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
106   }
107   else if (A->factor == FACTOR_CHOLESKY){
108     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
109   }
110   else SETERRQ(1,"MatSolve_SeqDense:Matrix must be factored to solve");
111   if (info) SETERRQ(1,"MatSolve_SeqDense:Bad solve");
112   PLogFlops(mat->n*mat->n - mat->n);
113   return 0;
114 }
115 static int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy)
116 {
117   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
118   int          one = 1, info;
119   Scalar       *x, *y;
120 
121   VecGetArray(xx,&x); VecGetArray(yy,&y);
122   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
123   /* assume if pivots exist then use LU; else Cholesky */
124   if (mat->pivots) {
125     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
126   }
127   else {
128     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
129   }
130   if (info) SETERRQ(1,"MatSolveTrans_SeqDense:Bad solve");
131   PLogFlops(mat->n*mat->n - mat->n);
132   return 0;
133 }
134 static int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
135 {
136   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
137   int          one = 1, info,ierr;
138   Scalar       *x, *y, sone = 1.0;
139   Vec          tmp = 0;
140 
141   VecGetArray(xx,&x); VecGetArray(yy,&y);
142   if (yy == zz) {
143     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
144     PLogObjectParent(A,tmp);
145     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
146   }
147   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
148   /* assume if pivots exist then use LU; else Cholesky */
149   if (mat->pivots) {
150     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
151   }
152   else {
153     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
154   }
155   if (info) SETERRQ(1,"MatSolveAdd_SeqDense:Bad solve");
156   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
157   else VecAXPY(&sone,zz,yy);
158   PLogFlops(mat->n*mat->n - mat->n);
159   return 0;
160 }
161 
162 static int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy)
163 {
164   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
165   int           one = 1, info,ierr;
166   Scalar        *x, *y, sone = 1.0;
167   Vec           tmp;
168 
169   VecGetArray(xx,&x); VecGetArray(yy,&y);
170   if (yy == zz) {
171     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
172     PLogObjectParent(A,tmp);
173     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
174   }
175   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
176   /* assume if pivots exist then use LU; else Cholesky */
177   if (mat->pivots) {
178     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
179   }
180   else {
181     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
182   }
183   if (info) SETERRQ(1,"MatSolveTransAdd_SeqDense:Bad solve");
184   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
185   else VecAXPY(&sone,zz,yy);
186   PLogFlops(mat->n*mat->n - mat->n);
187   return 0;
188 }
189 /* ------------------------------------------------------------------*/
190 static int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag,
191                           double shift,int its,Vec xx)
192 {
193   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
194   Scalar       *x, *b, *v = mat->v, zero = 0.0, xt;
195   int          o = 1,ierr, m = mat->m, i;
196 
197   if (flag & SOR_ZERO_INITIAL_GUESS) {
198     /* this is a hack fix, should have another version without
199        the second BLdot */
200     ierr = VecSet(&zero,xx); CHKERRQ(ierr);
201   }
202   VecGetArray(xx,&x); VecGetArray(bb,&b);
203   while (its--) {
204     if (flag & SOR_FORWARD_SWEEP){
205       for ( i=0; i<m; i++ ) {
206 #if defined(PETSC_COMPLEX)
207         /* cannot use BLAS dot for complex because compiler/linker is
208            not happy about returning a double complex */
209         int    _i;
210         Scalar sum = b[i];
211         for ( _i=0; _i<m; _i++ ) {
212           sum -= conj(v[i+_i*m])*x[_i];
213         }
214         xt = sum;
215 #else
216         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
217 #endif
218         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
219       }
220     }
221     if (flag & SOR_BACKWARD_SWEEP) {
222       for ( i=m-1; i>=0; i-- ) {
223 #if defined(PETSC_COMPLEX)
224         /* cannot use BLAS dot for complex because compiler/linker is
225            not happy about returning a double complex */
226         int    _i;
227         Scalar sum = b[i];
228         for ( _i=0; _i<m; _i++ ) {
229           sum -= conj(v[i+_i*m])*x[_i];
230         }
231         xt = sum;
232 #else
233         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
234 #endif
235         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
236       }
237     }
238   }
239   return 0;
240 }
241 
242 /* -----------------------------------------------------------------*/
243 int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy)
244 {
245   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
246   Scalar       *v = mat->v, *x, *y;
247   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
248   VecGetArray(xx,&x), VecGetArray(yy,&y);
249   LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One);
250   PLogFlops(2*mat->m*mat->n - mat->n);
251   return 0;
252 }
253 int MatMult_SeqDense(Mat A,Vec xx,Vec yy)
254 {
255   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
256   Scalar       *v = mat->v, *x, *y;
257   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
258   VecGetArray(xx,&x); VecGetArray(yy,&y);
259   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One);
260   PLogFlops(2*mat->m*mat->n - mat->m);
261   return 0;
262 }
263 int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
264 {
265   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
266   Scalar       *v = mat->v, *x, *y, *z;
267   int          _One=1; Scalar _DOne=1.0;
268   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
269   if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar));
270   LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One);
271   PLogFlops(2*mat->m*mat->n);
272   return 0;
273 }
274 int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
275 {
276   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
277   Scalar       *v = mat->v, *x, *y, *z;
278   int          _One=1;
279   Scalar       _DOne=1.0;
280   VecGetArray(xx,&x); VecGetArray(yy,&y); 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,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
413     ierr = MatAssemblyEnd(B,MAT_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,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
440     ierr = MatAssemblyEnd(B,MAT_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 || format == ASCII_FORMAT_INFO_DETAILED) {
460     return 0;  /* do nothing for now */
461   }
462   else if (format == ASCII_FORMAT_COMMON) {
463     for ( i=0; i<a->m; i++ ) {
464       v = a->v + i;
465       fprintf(fd,"row %d:",i);
466       for ( j=0; j<a->n; j++ ) {
467 #if defined(PETSC_COMPLEX)
468         if (real(*v) != 0.0 && imag(*v) != 0.0)
469           fprintf(fd," %d %g + %g i",j,real(*v),imag(*v));
470         else if (real(*v)) fprintf(fd," %d %g ",j,real(*v));
471         v += a->m;
472 #else
473         if (*v) fprintf(fd," %d %g ",j,*v);
474         v += a->m;
475 #endif
476       }
477       fprintf(fd,"\n");
478     }
479   }
480   else {
481 #if defined(PETSC_COMPLEX)
482     int allreal = 1;
483     /* determine if matrix has all real values */
484     v = a->v;
485     for ( i=0; i<a->m*a->n; i++ ) {
486       if (imag(v[i])) { allreal = 0; break ;}
487     }
488 #endif
489     for ( i=0; i<a->m; i++ ) {
490       v = a->v + i;
491       for ( j=0; j<a->n; j++ ) {
492 #if defined(PETSC_COMPLEX)
493         if (allreal) {
494           fprintf(fd,"%6.4e ",real(*v)); v += a->m;
495         } else {
496           fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m;
497         }
498 #else
499         fprintf(fd,"%6.4e ",*v); v += a->m;
500 #endif
501       }
502       fprintf(fd,"\n");
503     }
504   }
505   fflush(fd);
506   return 0;
507 }
508 
509 static int MatView_SeqDense_Binary(Mat A,Viewer viewer)
510 {
511   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
512   int          ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n;
513   int          format;
514   Scalar       *v, *anonz,*vals;
515 
516   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
517 
518   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
519   if (format == BINARY_FORMAT_NATIVE) {
520     /* store the matrix as a dense matrix */
521     col_lens = (int *) PetscMalloc( 4*sizeof(int) ); CHKPTRQ(col_lens);
522     col_lens[0] = MAT_COOKIE;
523     col_lens[1] = m;
524     col_lens[2] = n;
525     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
526     ierr = PetscBinaryWrite(fd,col_lens,4,BINARY_INT,1); CHKERRQ(ierr);
527     PetscFree(col_lens);
528 
529     /* write out matrix, by rows */
530     vals = (Scalar *) PetscMalloc(m*n*sizeof(Scalar)); CHKPTRQ(vals);
531     v    = a->v;
532     for ( i=0; i<m; i++ ) {
533       for ( j=0; j<n; j++ ) {
534         vals[i + j*m] = *v++;
535       }
536     }
537     ierr = PetscBinaryWrite(fd,vals,n*m,BINARY_SCALAR,0); CHKERRQ(ierr);
538     PetscFree(vals);
539   } else {
540     col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens);
541     col_lens[0] = MAT_COOKIE;
542     col_lens[1] = m;
543     col_lens[2] = n;
544     col_lens[3] = nz;
545 
546     /* store lengths of each row and write (including header) to file */
547     for ( i=0; i<m; i++ ) col_lens[4+i] = n;
548     ierr = PetscBinaryWrite(fd,col_lens,4+m,BINARY_INT,1); CHKERRQ(ierr);
549 
550     /* Possibly should write in smaller increments, not whole matrix at once? */
551     /* store column indices (zero start index) */
552     ict = 0;
553     for ( i=0; i<m; i++ ) {
554       for ( j=0; j<n; j++ ) col_lens[ict++] = j;
555     }
556     ierr = PetscBinaryWrite(fd,col_lens,nz,BINARY_INT,0); CHKERRQ(ierr);
557     PetscFree(col_lens);
558 
559     /* store nonzero values */
560     anonz = (Scalar *) PetscMalloc((nz)*sizeof(Scalar)); CHKPTRQ(anonz);
561     ict = 0;
562     for ( i=0; i<m; i++ ) {
563       v = a->v + i;
564       for ( j=0; j<n; j++ ) {
565         anonz[ict++] = *v; v += a->m;
566       }
567     }
568     ierr = PetscBinaryWrite(fd,anonz,nz,BINARY_SCALAR,0); CHKERRQ(ierr);
569     PetscFree(anonz);
570   }
571   return 0;
572 }
573 
574 static int MatView_SeqDense(PetscObject obj,Viewer viewer)
575 {
576   Mat          A = (Mat) obj;
577   Mat_SeqDense *a = (Mat_SeqDense*) A->data;
578   ViewerType   vtype;
579   int          ierr;
580 
581   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
582 
583   if (vtype == MATLAB_VIEWER) {
584     return ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v);
585   }
586   else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
587     return MatView_SeqDense_ASCII(A,viewer);
588   }
589   else if (vtype == BINARY_FILE_VIEWER) {
590     return MatView_SeqDense_Binary(A,viewer);
591   }
592   return 0;
593 }
594 
595 static int MatDestroy_SeqDense(PetscObject obj)
596 {
597   Mat          mat = (Mat) obj;
598   Mat_SeqDense *l = (Mat_SeqDense *) mat->data;
599 #if defined(PETSC_LOG)
600   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
601 #endif
602   if (l->pivots) PetscFree(l->pivots);
603   if (!l->user_alloc) PetscFree(l->v);
604   PetscFree(l);
605   PLogObjectDestroy(mat);
606   PetscHeaderDestroy(mat);
607   return 0;
608 }
609 
610 static int MatTranspose_SeqDense(Mat A,Mat *matout)
611 {
612   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
613   int          k, j, m, n;
614   Scalar       *v, tmp;
615 
616   v = mat->v; m = mat->m; n = mat->n;
617   if (matout == PETSC_NULL) { /* in place transpose */
618     if (m != n) SETERRQ(1,"MatTranspose_SeqDense:Supports square matrix only in-place");
619     for ( j=0; j<m; j++ ) {
620       for ( k=0; k<j; k++ ) {
621         tmp = v[j + k*n];
622         v[j + k*n] = v[k + j*n];
623         v[k + j*n] = tmp;
624       }
625     }
626   }
627   else { /* out-of-place transpose */
628     int          ierr;
629     Mat          tmat;
630     Mat_SeqDense *tmatd;
631     Scalar       *v2;
632     ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr);
633     tmatd = (Mat_SeqDense *) tmat->data;
634     v = mat->v; v2 = tmatd->v;
635     for ( j=0; j<n; j++ ) {
636       for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m];
637     }
638     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
639     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
640     *matout = tmat;
641   }
642   return 0;
643 }
644 
645 static int MatEqual_SeqDense(Mat A1,Mat A2, PetscTruth *flg)
646 {
647   Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data;
648   Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data;
649   int          i;
650   Scalar       *v1 = mat1->v, *v2 = mat2->v;
651 
652   if (A2->type != MATSEQDENSE) SETERRQ(1,"MatEqual_SeqDense:Both matrices should be of type  MATSEQDENSE");
653   if (mat1->m != mat2->m) {*flg = PETSC_FALSE; return 0;}
654   if (mat1->n != mat2->n) {*flg = PETSC_FALSE; return 0;}
655   for ( i=0; i<mat1->m*mat1->n; i++ ) {
656     if (*v1 != *v2) {*flg = PETSC_FALSE; return 0;}
657     v1++; v2++;
658   }
659   *flg = PETSC_TRUE;
660   return 0;
661 }
662 
663 static int MatGetDiagonal_SeqDense(Mat A,Vec v)
664 {
665   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
666   int          i, n, len;
667   Scalar       *x, zero = 0.0;
668 
669   VecSet(&zero,v);
670   VecGetArray(v,&x); VecGetSize(v,&n);
671   len = PetscMin(mat->m,mat->n);
672   if (n != mat->m) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec");
673   for ( i=0; i<len; i++ ) {
674     x[i] = mat->v[i*mat->m + i];
675   }
676   return 0;
677 }
678 
679 static int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
680 {
681   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
682   Scalar       *l,*r,x,*v;
683   int          i,j,m = mat->m, n = mat->n;
684 
685   if (ll) {
686     VecGetArray(ll,&l); VecGetSize(ll,&m);
687     if (m != mat->m) SETERRQ(1,"MatDiagonalScale_SeqDense:Left scaling vec wrong size");
688     for ( i=0; i<m; i++ ) {
689       x = l[i];
690       v = mat->v + i;
691       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
692     }
693     PLogFlops(n*m);
694   }
695   if (rr) {
696     VecGetArray(rr,&r); VecGetSize(rr,&n);
697     if (n != mat->n) SETERRQ(1,"MatDiagonalScale_SeqDense:Right scaling vec wrong size");
698     for ( i=0; i<n; i++ ) {
699       x = r[i];
700       v = mat->v + i*m;
701       for ( j=0; j<m; j++ ) { (*v++) *= x;}
702     }
703     PLogFlops(n*m);
704   }
705   return 0;
706 }
707 
708 static int MatNorm_SeqDense(Mat A,NormType type,double *norm)
709 {
710   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
711   Scalar       *v = mat->v;
712   double       sum = 0.0;
713   int          i, j;
714 
715   if (type == NORM_FROBENIUS) {
716     for (i=0; i<mat->n*mat->m; i++ ) {
717 #if defined(PETSC_COMPLEX)
718       sum += real(conj(*v)*(*v)); v++;
719 #else
720       sum += (*v)*(*v); v++;
721 #endif
722     }
723     *norm = sqrt(sum);
724     PLogFlops(2*mat->n*mat->m);
725   }
726   else if (type == NORM_1) {
727     *norm = 0.0;
728     for ( j=0; j<mat->n; j++ ) {
729       sum = 0.0;
730       for ( i=0; i<mat->m; i++ ) {
731         sum += PetscAbsScalar(*v);  v++;
732       }
733       if (sum > *norm) *norm = sum;
734     }
735     PLogFlops(mat->n*mat->m);
736   }
737   else if (type == NORM_INFINITY) {
738     *norm = 0.0;
739     for ( j=0; j<mat->m; j++ ) {
740       v = mat->v + j;
741       sum = 0.0;
742       for ( i=0; i<mat->n; i++ ) {
743         sum += PetscAbsScalar(*v); v += mat->m;
744       }
745       if (sum > *norm) *norm = sum;
746     }
747     PLogFlops(mat->n*mat->m);
748   }
749   else {
750     SETERRQ(1,"MatNorm_SeqDense:No two norm");
751   }
752   return 0;
753 }
754 
755 static int MatSetOption_SeqDense(Mat A,MatOption op)
756 {
757   Mat_SeqDense *aij = (Mat_SeqDense *) A->data;
758 
759   if (op == MAT_ROW_ORIENTED)            aij->roworiented = 1;
760   else if (op == MAT_COLUMN_ORIENTED)    aij->roworiented = 0;
761   else if (op == MAT_ROWS_SORTED ||
762            op == MAT_COLUMNS_SORTED ||
763            op == MAT_SYMMETRIC ||
764            op == MAT_STRUCTURALLY_SYMMETRIC ||
765            op == MAT_NO_NEW_NONZERO_LOCATIONS ||
766            op == MAT_YES_NEW_NONZERO_LOCATIONS ||
767            op == MAT_NO_NEW_DIAGONALS ||
768            op == MAT_YES_NEW_DIAGONALS)
769     PLogInfo(A,"Info:MatSetOption_SeqDense:Option ignored\n");
770   else
771     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqDense:unknown option");}
772   return 0;
773 }
774 
775 static int MatZeroEntries_SeqDense(Mat A)
776 {
777   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
778   PetscMemzero(l->v,l->m*l->n*sizeof(Scalar));
779   return 0;
780 }
781 
782 static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
783 {
784   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
785   int          n = l->n, i, j,ierr,N, *rows;
786   Scalar       *slot;
787 
788   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
789   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
790   for ( i=0; i<N; i++ ) {
791     slot = l->v + rows[i];
792     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
793   }
794   if (diag) {
795     for ( i=0; i<N; i++ ) {
796       slot = l->v + (n+1)*rows[i];
797       *slot = *diag;
798     }
799   }
800   ISRestoreIndices(is,&rows);
801   return 0;
802 }
803 
804 static int MatGetSize_SeqDense(Mat A,int *m,int *n)
805 {
806   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
807   *m = mat->m; *n = mat->n;
808   return 0;
809 }
810 
811 static int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
812 {
813   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
814   *m = 0; *n = mat->m;
815   return 0;
816 }
817 
818 static int MatGetArray_SeqDense(Mat A,Scalar **array)
819 {
820   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
821   *array = mat->v;
822   return 0;
823 }
824 
825 static int MatRestoreArray_SeqDense(Mat A,Scalar **array)
826 {
827   return 0;
828 }
829 
830 static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,
831                                     Mat *submat)
832 {
833   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
834   int          nznew, *smap, i, j, ierr, oldcols = mat->n;
835   int          *irow, *icol, nrows, ncols, *cwork;
836   Scalar       *vwork, *val;
837   Mat          newmat;
838 
839   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
840   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
841   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
842   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
843 
844   smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap);
845   cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork);
846   vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
847   PetscMemzero((char*)smap,oldcols*sizeof(int));
848   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
849 
850   /* Create and fill new matrix */
851   ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr);
852   for (i=0; i<nrows; i++) {
853     nznew = 0;
854     val   = mat->v + irow[i];
855     for (j=0; j<oldcols; j++) {
856       if (smap[j]) {
857         cwork[nznew]   = smap[j] - 1;
858         vwork[nznew++] = val[j * mat->m];
859       }
860     }
861     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES);
862            CHKERRQ(ierr);
863   }
864   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
865   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
866 
867   /* Free work space */
868   PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
869   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
870   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
871   *submat = newmat;
872   return 0;
873 }
874 
875 static int MatGetSubMatrices_SeqDense(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
876                                     Mat **B)
877 {
878   int ierr,i;
879 
880   if (scall == MAT_INITIAL_MATRIX) {
881     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
882   }
883 
884   for ( i=0; i<n; i++ ) {
885     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr);
886   }
887   return 0;
888 }
889 
890 static int MatCopy_SeqDense(Mat A, Mat B)
891 {
892   Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data;
893   if (B->type != MATSEQDENSE) return MatCopy_Basic(A,B);
894   if (a->m != b->m || a->n != b->n) SETERRQ(1,"MatCopy_SeqDense:size(B) != size(A)");
895   PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar));
896   return 0;
897 }
898 
899 /* -------------------------------------------------------------------*/
900 static struct _MatOps MatOps = {MatSetValues_SeqDense,
901        MatGetRow_SeqDense,
902        MatRestoreRow_SeqDense,
903        MatMult_SeqDense,
904        MatMultAdd_SeqDense,
905        MatMultTrans_SeqDense,
906        MatMultTransAdd_SeqDense,
907        MatSolve_SeqDense,
908        MatSolveAdd_SeqDense,
909        MatSolveTrans_SeqDense,
910        MatSolveTransAdd_SeqDense,
911        MatLUFactor_SeqDense,
912        MatCholeskyFactor_SeqDense,
913        MatRelax_SeqDense,
914        MatTranspose_SeqDense,
915        MatGetInfo_SeqDense,
916        MatEqual_SeqDense,
917        MatGetDiagonal_SeqDense,
918        MatDiagonalScale_SeqDense,
919        MatNorm_SeqDense,
920        0,
921        0,
922        0,
923        MatSetOption_SeqDense,
924        MatZeroEntries_SeqDense,
925        MatZeroRows_SeqDense,
926        0,
927        MatLUFactorSymbolic_SeqDense,
928        MatLUFactorNumeric_SeqDense,
929        MatCholeskyFactorSymbolic_SeqDense,
930        MatCholeskyFactorNumeric_SeqDense,
931        MatGetSize_SeqDense,
932        MatGetSize_SeqDense,
933        MatGetOwnershipRange_SeqDense,
934        0,
935        0,
936        MatGetArray_SeqDense,
937        MatRestoreArray_SeqDense,
938        0,
939        MatConvertSameType_SeqDense,0,0,0,0,
940        MatAXPY_SeqDense,MatGetSubMatrices_SeqDense,0,
941        MatGetValues_SeqDense,
942        MatCopy_SeqDense,0,MatScale_SeqDense};
943 
944 
945 /*@C
946    MatCreateSeqDense - Creates a sequential dense matrix that
947    is stored in column major order (the usual Fortran 77 manner). Many
948    of the matrix operations use the BLAS and LAPACK routines.
949 
950    Input Parameters:
951 .  comm - MPI communicator, set to MPI_COMM_SELF
952 .  m - number of rows
953 .  n - number of columns
954 .  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
955    to control all matrix memory allocation.
956 
957    Output Parameter:
958 .  A - the matrix
959 
960   Notes:
961   The data input variable is intended primarily for Fortran programmers
962   who wish to allocate their own matrix memory space.  Most users should
963   set data=PETSC_NULL.
964 
965 .keywords: dense, matrix, LAPACK, BLAS
966 
967 .seealso: MatCreate(), MatSetValues()
968 @*/
969 int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A)
970 {
971   Mat          B;
972   Mat_SeqDense *b;
973   int          ierr,flg;
974 
975   *A            = 0;
976   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQDENSE,comm);
977   PLogObjectCreate(B);
978   b             = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(b);
979   PetscMemzero(b,sizeof(Mat_SeqDense));
980   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
981   B->destroy    = MatDestroy_SeqDense;
982   B->view       = MatView_SeqDense;
983   B->factor     = 0;
984   PLogObjectMemory(B,sizeof(struct _Mat));
985   B->data       = (void *) b;
986 
987   b->m = m;  B->m = m; B->M = m;
988   b->n = n;  B->n = n; B->N = n;
989   b->pivots       = 0;
990   b->roworiented  = 1;
991   if (data == PETSC_NULL) {
992     b->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(b->v);
993     PetscMemzero(b->v,m*n*sizeof(Scalar));
994     b->user_alloc = 0;
995     PLogObjectMemory(B,n*m*sizeof(Scalar));
996   }
997   else { /* user-allocated storage */
998     b->v = data;
999     b->user_alloc = 1;
1000   }
1001   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
1002   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr);}
1003   *A = B;
1004   return 0;
1005 }
1006 
1007 int MatCreate_SeqDense(Mat A,Mat *newmat)
1008 {
1009   Mat_SeqDense *m = (Mat_SeqDense *) A->data;
1010   return MatCreateSeqDense(A->comm,m->m,m->n,PETSC_NULL,newmat);
1011 }
1012