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