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