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