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