xref: /petsc/src/mat/impls/dense/seq/dense.c (revision a86d99e19318d085cf9119bf351e5b9153e19789)
1 #ifndef lint
2 static char vcid[] = "$Id: dense.c,v 1.90 1996/02/01 18:52:28 curfman Exp balay $";
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   newi->assembled = PETSC_TRUE;
369   *newmat = newi;
370   return 0;
371 }
372 
373 #include "sysio.h"
374 
375 int MatLoad_SeqDense(Viewer bview,MatType type,Mat *A)
376 {
377   Mat_SeqDense *a;
378   Mat          B;
379   int          *scols, i, j, nz, ierr, fd, header[4], size;
380   int          *rowlengths = 0, M, N, *cols;
381   Scalar       *vals, *svals, *v;
382   PetscObject  vobj = (PetscObject) bview;
383   MPI_Comm     comm = vobj->comm;
384 
385   MPI_Comm_size(comm,&size);
386   if (size > 1) SETERRQ(1,"MatLoad_SeqDense: view must have one processor");
387   ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr);
388   ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr);
389   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqDense:Not matrix object");
390   M = header[1]; N = header[2]; nz = header[3];
391 
392   /* read row lengths */
393   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
394   ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr);
395 
396   /* create our matrix */
397   ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
398   B = *A;
399   a = (Mat_SeqDense *) B->data;
400   v = a->v;
401 
402   /* read column indices and nonzeros */
403   cols = scols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(cols);
404   ierr = SYRead(fd,cols,nz,SYINT); CHKERRQ(ierr);
405   vals = svals = (Scalar *) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals);
406   ierr = SYRead(fd,vals,nz,SYSCALAR); CHKERRQ(ierr);
407 
408   /* insert into matrix */
409   for ( i=0; i<M; i++ ) {
410     for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j];
411     svals += rowlengths[i]; scols += rowlengths[i];
412   }
413   PetscFree(vals); PetscFree(cols); PetscFree(rowlengths);
414 
415   ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
416   ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
417   return 0;
418 }
419 
420 #include "pinclude/pviewer.h"
421 #include "sysio.h"
422 
423 static int MatView_SeqDense_ASCII(Mat A,Viewer viewer)
424 {
425   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
426   int          ierr, i, j, format;
427   FILE         *fd;
428   char         *outputname;
429   Scalar       *v;
430 
431   ierr = ViewerFileGetPointer(viewer,&fd); CHKERRQ(ierr);
432   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
433   ierr = ViewerFileGetFormat_Private(viewer,&format);
434   if (format == FILE_FORMAT_INFO) {
435     ;  /* do nothing for now */
436   }
437   else {
438     for ( i=0; i<a->m; i++ ) {
439       v = a->v + i;
440       for ( j=0; j<a->n; j++ ) {
441 #if defined(PETSC_COMPLEX)
442         fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m;
443 #else
444         fprintf(fd,"%6.4e ",*v); v += a->m;
445 #endif
446       }
447       fprintf(fd,"\n");
448     }
449   }
450   fflush(fd);
451   return 0;
452 }
453 
454 static int MatView_SeqDense_Binary(Mat A,Viewer viewer)
455 {
456   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
457   int          ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n;
458   Scalar       *v, *anonz;
459 
460   ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr);
461   col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens);
462   col_lens[0] = MAT_COOKIE;
463   col_lens[1] = m;
464   col_lens[2] = n;
465   col_lens[3] = nz;
466 
467   /* store lengths of each row and write (including header) to file */
468   for ( i=0; i<m; i++ ) col_lens[4+i] = n;
469   ierr = SYWrite(fd,col_lens,4+m,SYINT,1); CHKERRQ(ierr);
470 
471   /* Possibly should write in smaller increments, not whole matrix at once? */
472  /* store column indices (zero start index) */
473   ict = 0;
474   for ( i=0; i<m; i++ ) {
475     for ( j=0; j<n; j++ ) col_lens[ict++] = j;
476   }
477   ierr = SYWrite(fd,col_lens,nz,SYINT,0); CHKERRQ(ierr);
478   PetscFree(col_lens);
479 
480   /* store nonzero values */
481   anonz = (Scalar *) PetscMalloc((nz)*sizeof(Scalar)); CHKPTRQ(anonz);
482   ict = 0;
483   for ( i=0; i<m; i++ ) {
484     v = a->v + i;
485     for ( j=0; j<n; j++ ) {
486       anonz[ict++] = *v; v += a->m;
487     }
488   }
489   ierr = SYWrite(fd,anonz,nz,SYSCALAR,0); CHKERRQ(ierr);
490   PetscFree(anonz);
491   return 0;
492 }
493 
494 static int MatView_SeqDense(PetscObject obj,Viewer viewer)
495 {
496   Mat          A = (Mat) obj;
497   Mat_SeqDense *a = (Mat_SeqDense*) A->data;
498   PetscObject  vobj = (PetscObject) viewer;
499 
500   if (!viewer) {
501     viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer;
502   }
503   if (vobj->cookie == VIEWER_COOKIE) {
504     if (vobj->type == MATLAB_VIEWER) {
505       return ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v);
506     }
507     else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER) {
508       return MatView_SeqDense_ASCII(A,viewer);
509     }
510     else if (vobj->type == BINARY_FILE_VIEWER) {
511       return MatView_SeqDense_Binary(A,viewer);
512     }
513   }
514   return 0;
515 }
516 
517 static int MatDestroy_SeqDense(PetscObject obj)
518 {
519   Mat          mat = (Mat) obj;
520   Mat_SeqDense *l = (Mat_SeqDense *) mat->data;
521 #if defined(PETSC_LOG)
522   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
523 #endif
524   if (l->pivots) PetscFree(l->pivots);
525   if (!l->user_alloc) PetscFree(l->v);
526   PetscFree(l);
527   PLogObjectDestroy(mat);
528   PetscHeaderDestroy(mat);
529   return 0;
530 }
531 
532 static int MatTranspose_SeqDense(Mat A,Mat *matout)
533 {
534   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
535   int          k, j, m, n;
536   Scalar       *v, tmp;
537 
538   v = mat->v; m = mat->m; n = mat->n;
539   if (matout == PETSC_NULL) { /* in place transpose */
540     if (m != n) SETERRQ(1,"MatTranspose_SeqDense:Supports square matrix only in-place");
541     for ( j=0; j<m; j++ ) {
542       for ( k=0; k<j; k++ ) {
543         tmp = v[j + k*n];
544         v[j + k*n] = v[k + j*n];
545         v[k + j*n] = tmp;
546       }
547     }
548   }
549   else { /* out-of-place transpose */
550     int          ierr;
551     Mat          tmat;
552     Mat_SeqDense *tmatd;
553     Scalar       *v2;
554     ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr);
555     tmatd = (Mat_SeqDense *) tmat->data;
556     v = mat->v; v2 = tmatd->v;
557     for ( j=0; j<n; j++ ) {
558       for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m];
559     }
560     ierr = MatAssemblyBegin(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
561     ierr = MatAssemblyEnd(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
562     *matout = tmat;
563   }
564   return 0;
565 }
566 
567 static int MatEqual_SeqDense(Mat A1,Mat A2, int *flg)
568 {
569   Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data;
570   Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data;
571   int          i;
572   Scalar       *v1 = mat1->v, *v2 = mat2->v;
573 
574   if (A2->type != MATSEQDENSE) SETERRQ(1,"MatEqual_SeqDense:Both matrices should be of type  MATSEQDENSE");
575   if (mat1->m != mat2->m) { *flg = 0; return 0;}
576   if (mat1->n != mat2->n) {*flg =0; return 0;}
577   for ( i=0; i<mat1->m*mat1->n; i++ ) {
578     if (*v1 != *v2) {*flg =0; return 0;}
579     v1++; v2++;
580   }
581   *flg = 1;
582   return 0;
583 }
584 
585 static int MatGetDiagonal_SeqDense(Mat A,Vec v)
586 {
587   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
588   int          i, n;
589   Scalar       *x;
590   VecGetArray(v,&x); VecGetSize(v,&n);
591   if (n != mat->m) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec");
592   for ( i=0; i<mat->m; i++ ) {
593     x[i] = mat->v[i*mat->m + i];
594   }
595   return 0;
596 }
597 
598 static int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
599 {
600   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
601   Scalar       *l,*r,x,*v;
602   int          i,j,m = mat->m, n = mat->n;
603 
604   if (ll) {
605     VecGetArray(ll,&l); VecGetSize(ll,&m);
606     if (m != mat->m) SETERRQ(1,"MatDiagonalScale_SeqDense:Left scaling vec wrong size");
607     PLogFlops(n*m);
608     for ( i=0; i<m; i++ ) {
609       x = l[i];
610       v = mat->v + i;
611       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
612     }
613   }
614   if (rr) {
615     VecGetArray(rr,&r); VecGetSize(rr,&n);
616     if (n != mat->n) SETERRQ(1,"MatDiagonalScale_SeqDense:Right scaling vec wrong size");
617     PLogFlops(n*m);
618     for ( i=0; i<n; i++ ) {
619       x = r[i];
620       v = mat->v + i*m;
621       for ( j=0; j<m; j++ ) { (*v++) *= x;}
622     }
623   }
624   return 0;
625 }
626 
627 static int MatNorm_SeqDense(Mat A,NormType type,double *norm)
628 {
629   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
630   Scalar       *v = mat->v;
631   double       sum = 0.0;
632   int          i, j;
633 
634   if (type == NORM_FROBENIUS) {
635     for (i=0; i<mat->n*mat->m; i++ ) {
636 #if defined(PETSC_COMPLEX)
637       sum += real(conj(*v)*(*v)); v++;
638 #else
639       sum += (*v)*(*v); v++;
640 #endif
641     }
642     *norm = sqrt(sum);
643     PLogFlops(2*mat->n*mat->m);
644   }
645   else if (type == NORM_1) {
646     *norm = 0.0;
647     for ( j=0; j<mat->n; j++ ) {
648       sum = 0.0;
649       for ( i=0; i<mat->m; i++ ) {
650         sum += PetscAbsScalar(*v);  v++;
651       }
652       if (sum > *norm) *norm = sum;
653     }
654     PLogFlops(mat->n*mat->m);
655   }
656   else if (type == NORM_INFINITY) {
657     *norm = 0.0;
658     for ( j=0; j<mat->m; j++ ) {
659       v = mat->v + j;
660       sum = 0.0;
661       for ( i=0; i<mat->n; i++ ) {
662         sum += PetscAbsScalar(*v); v += mat->m;
663       }
664       if (sum > *norm) *norm = sum;
665     }
666     PLogFlops(mat->n*mat->m);
667   }
668   else {
669     SETERRQ(1,"MatNorm_SeqDense:No two norm");
670   }
671   return 0;
672 }
673 
674 static int MatSetOption_SeqDense(Mat A,MatOption op)
675 {
676   Mat_SeqDense *aij = (Mat_SeqDense *) A->data;
677 
678   if (op == ROW_ORIENTED)            aij->roworiented = 1;
679   else if (op == COLUMN_ORIENTED)    aij->roworiented = 0;
680   else if (op == ROWS_SORTED ||
681            op == COLUMNS_SORTED ||
682            op == SYMMETRIC_MATRIX ||
683            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
684            op == NO_NEW_NONZERO_LOCATIONS ||
685            op == YES_NEW_NONZERO_LOCATIONS ||
686            op == NO_NEW_DIAGONALS ||
687            op == YES_NEW_DIAGONALS)
688     PLogInfo((PetscObject)A,"Info:MatSetOption_SeqDense:Option ignored\n");
689   else
690     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqDense:unknown option");}
691   return 0;
692 }
693 
694 static int MatZeroEntries_SeqDense(Mat A)
695 {
696   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
697   PetscMemzero(l->v,l->m*l->n*sizeof(Scalar));
698   return 0;
699 }
700 
701 static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
702 {
703   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
704   int          n = l->n, i, j,ierr,N, *rows;
705   Scalar       *slot;
706 
707   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
708   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
709   for ( i=0; i<N; i++ ) {
710     slot = l->v + rows[i];
711     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
712   }
713   if (diag) {
714     for ( i=0; i<N; i++ ) {
715       slot = l->v + (n+1)*rows[i];
716       *slot = *diag;
717     }
718   }
719   ISRestoreIndices(is,&rows);
720   return 0;
721 }
722 
723 static int MatGetSize_SeqDense(Mat A,int *m,int *n)
724 {
725   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
726   *m = mat->m; *n = mat->n;
727   return 0;
728 }
729 
730 static int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
731 {
732   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
733   *m = 0; *n = mat->m;
734   return 0;
735 }
736 
737 static int MatGetArray_SeqDense(Mat A,Scalar **array)
738 {
739   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
740   *array = mat->v;
741   return 0;
742 }
743 
744 
745 static int MatGetSubMatrixInPlace_SeqDense(Mat A,IS isrow,IS iscol)
746 {
747   SETERRQ(1,"MatGetSubMatrixInPlace_SeqDense:not done");
748 }
749 
750 static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,
751                                     Mat *submat)
752 {
753   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
754   int          nznew, *smap, i, j, ierr, oldcols = mat->n;
755   int          *irow, *icol, nrows, ncols, *cwork;
756   Scalar       *vwork, *val;
757   Mat          newmat;
758 
759   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
760   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
761   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
762   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
763 
764   smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap);
765   cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork);
766   vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
767   PetscMemzero((char*)smap,oldcols*sizeof(int));
768   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
769 
770   /* Create and fill new matrix */
771   ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr);
772   for (i=0; i<nrows; i++) {
773     nznew = 0;
774     val   = mat->v + irow[i];
775     for (j=0; j<oldcols; j++) {
776       if (smap[j]) {
777         cwork[nznew]   = smap[j] - 1;
778         vwork[nznew++] = val[j * mat->m];
779       }
780     }
781     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES);
782            CHKERRQ(ierr);
783   }
784   ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
785   ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
786 
787   /* Free work space */
788   PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
789   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
790   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
791   *submat = newmat;
792   return 0;
793 }
794 
795 static int MatCopy_SeqDense(Mat A, Mat B)
796 {
797   Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data;
798   if (B->type != MATSEQDENSE) return MatCopy_Basic(A,B);
799   if (a->m != b->m || a->n != b->n) SETERRQ(1,"MatCopy_SeqDense:size(B) != size(A)");
800   PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar));
801   return 0;
802 }
803 
804 /* -------------------------------------------------------------------*/
805 static struct _MatOps MatOps = {MatSetValues_SeqDense,
806        MatGetRow_SeqDense, MatRestoreRow_SeqDense,
807        MatMult_SeqDense, MatMultAdd_SeqDense,
808        MatMultTrans_SeqDense, MatMultTransAdd_SeqDense,
809        MatSolve_SeqDense,MatSolveAdd_SeqDense,
810        MatSolveTrans_SeqDense,MatSolveTransAdd_SeqDense,
811        MatLUFactor_SeqDense,MatCholeskyFactor_SeqDense,
812        MatRelax_SeqDense,
813        MatTranspose_SeqDense,
814        MatGetInfo_SeqDense,MatEqual_SeqDense,
815        MatGetDiagonal_SeqDense,MatDiagonalScale_SeqDense,MatNorm_SeqDense,
816        0,0,
817        0, MatSetOption_SeqDense,MatZeroEntries_SeqDense,MatZeroRows_SeqDense,0,
818        MatLUFactorSymbolic_SeqDense,MatLUFactorNumeric_SeqDense,
819        MatCholeskyFactorSymbolic_SeqDense,MatCholeskyFactorNumeric_SeqDense,
820        MatGetSize_SeqDense,MatGetSize_SeqDense,MatGetOwnershipRange_SeqDense,
821        0,0,MatGetArray_SeqDense,0,0,
822        MatGetSubMatrix_SeqDense,MatGetSubMatrixInPlace_SeqDense,
823        MatConvertSameType_SeqDense,0,0,0,0,
824        MatAXPY_SeqDense,0,0,
825        MatGetValues_SeqDense,
826        MatCopy_SeqDense};
827 
828 /*@C
829    MatCreateSeqDense - Creates a sequential dense matrix that
830    is stored in column major order (the usual Fortran 77 manner). Many
831    of the matrix operations use the BLAS and LAPACK routines.
832 
833    Input Parameters:
834 .  comm - MPI communicator, set to MPI_COMM_SELF
835 .  m - number of rows
836 .  n - number of columns
837 .  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
838    to control all matrix memory allocation.
839 
840    Output Parameter:
841 .  newmat - the matrix
842 
843   Notes:
844   The data input variable is intended primarily for Fortran programmers
845   who wish to allocate their own matrix memory space.  Most users should
846   set data=PETSC_NULL.
847 
848 .keywords: dense, matrix, LAPACK, BLAS
849 
850 .seealso: MatCreate(), MatSetValues()
851 @*/
852 int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *newmat)
853 {
854   Mat          mat;
855   Mat_SeqDense *l;
856   int          ierr,flg;
857 
858   *newmat        = 0;
859 
860   PetscHeaderCreate(mat,_Mat,MAT_COOKIE,MATSEQDENSE,comm);
861   PLogObjectCreate(mat);
862   l               = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(l);
863   PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps));
864   mat->destroy    = MatDestroy_SeqDense;
865   mat->view       = MatView_SeqDense;
866   mat->factor     = 0;
867   PLogObjectMemory(mat,sizeof(struct _Mat));
868   mat->data       = (void *) l;
869 
870   l->m            = m;
871   l->n            = n;
872   l->pivots       = 0;
873   l->roworiented  = 1;
874   if (data == PETSC_NULL) {
875     l->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(l->v);
876     PetscMemzero(l->v,m*n*sizeof(Scalar));
877     l->user_alloc = 0;
878     PLogObjectMemory(mat,n*m*sizeof(Scalar));
879   }
880   else { /* user-allocated storage */
881     l->v = data;
882     l->user_alloc = 1;
883   }
884 
885   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
886   if (flg) {
887     ierr = MatPrintHelp(mat); CHKERRQ(ierr);
888   }
889   *newmat = mat;
890   return 0;
891 }
892 
893 int MatCreate_SeqDense(Mat A,Mat *newmat)
894 {
895   Mat_SeqDense *m = (Mat_SeqDense *) A->data;
896   return MatCreateSeqDense(A->comm,m->m,m->n,PETSC_NULL,newmat);
897 }
898