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