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