xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 5a5d4f665b3034485fbcdfbe42345f740a4bd170)
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: dense.c,v 1.132 1997/10/19 03:25:08 bsmith Exp bsmith $";
3 #endif
4 /*
5      Defines the basic matrix operations for sequential dense.
6 */
7 
8 #include "src/mat/impls/dense/seq/dense.h"
9 #include "pinclude/plapack.h"
10 #include "pinclude/pviewer.h"
11 
12 #undef __FUNC__
13 #define __FUNC__ "MatAXPY_SeqDense"
14 int MatAXPY_SeqDense(Scalar *alpha,Mat X,Mat Y)
15 {
16   Mat_SeqDense *x = (Mat_SeqDense*) X->data,*y = (Mat_SeqDense*) Y->data;
17   int          N = x->m*x->n, one = 1;
18 
19   PetscFunctionBegin;
20   BLaxpy_( &N, alpha, x->v, &one, y->v, &one );
21   PLogFlops(2*N-1);
22   PetscFunctionReturn(0);
23 }
24 
25 #undef __FUNC__
26 #define __FUNC__ "MatGetInfo_SeqDense"
27 int MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
28 {
29   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
30   int          i,N = a->m*a->n,count = 0;
31   Scalar       *v = a->v;
32 
33   PetscFunctionBegin;
34   for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;}
35 
36   info->rows_global       = (double)a->m;
37   info->columns_global    = (double)a->n;
38   info->rows_local        = (double)a->m;
39   info->columns_local     = (double)a->n;
40   info->block_size        = 1.0;
41   info->nz_allocated      = (double)N;
42   info->nz_used           = (double)count;
43   info->nz_unneeded       = (double)(N-count);
44   info->assemblies        = (double)A->num_ass;
45   info->mallocs           = 0;
46   info->memory            = A->mem;
47   info->fill_ratio_given  = 0;
48   info->fill_ratio_needed = 0;
49   info->factor_mallocs    = 0;
50 
51   PetscFunctionReturn(0);
52 }
53 
54 #undef __FUNC__
55 #define __FUNC__ "MatScale_SeqDense"
56 int MatScale_SeqDense(Scalar *alpha,Mat inA)
57 {
58   Mat_SeqDense *a = (Mat_SeqDense *) inA->data;
59   int          one = 1, nz;
60 
61   PetscFunctionBegin;
62   nz = a->m*a->n;
63   BLscal_( &nz, alpha, a->v, &one );
64   PLogFlops(nz);
65   PetscFunctionReturn(0);
66 }
67 
68 /* ---------------------------------------------------------------*/
69 /* COMMENT: I have chosen to hide column permutation in the pivots,
70    rather than put it in the Mat->col slot.*/
71 #undef __FUNC__
72 #define __FUNC__ "MatLUFactor_SeqDense"
73 int MatLUFactor_SeqDense(Mat A,IS row,IS col,double f)
74 {
75   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
76   int          info;
77 
78   PetscFunctionBegin;
79   if (!mat->pivots) {
80     mat->pivots = (int *) PetscMalloc((mat->m+1)*sizeof(int));CHKPTRQ(mat->pivots);
81     PLogObjectMemory(A,mat->m*sizeof(int));
82   }
83   LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info);
84   if (info<0) SETERRQ(1,0,"Bad argument to LU factorization");
85   if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,0,"Bad LU factorization");
86   A->factor = FACTOR_LU;
87   PLogFlops((2*mat->n*mat->n*mat->n)/3);
88   PetscFunctionReturn(0);
89 }
90 
91 #undef __FUNC__
92 #define __FUNC__ "MatConvertSameType_SeqDense"
93 int MatConvertSameType_SeqDense(Mat A,Mat *newmat,int cpvalues)
94 {
95   Mat_SeqDense *mat = (Mat_SeqDense *) A->data, *l;
96   int          ierr;
97   Mat          newi;
98 
99   PetscFunctionBegin;
100   ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,PETSC_NULL,&newi); CHKERRQ(ierr);
101   l = (Mat_SeqDense *) newi->data;
102   if (cpvalues == COPY_VALUES) {
103     PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
104   }
105   newi->assembled = PETSC_TRUE;
106   *newmat = newi;
107   PetscFunctionReturn(0);
108 }
109 
110 #undef __FUNC__
111 #define __FUNC__ "MatLUFactorSymbolic_SeqDense"
112 int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,double f,Mat *fact)
113 {
114   int ierr;
115 
116   PetscFunctionBegin;
117   ierr = MatConvertSameType_SeqDense(A,fact,PETSC_FALSE);CHKERRQ(ierr);
118   PetscFunctionReturn(0);
119 }
120 
121 #undef __FUNC__
122 #define __FUNC__ "MatLUFactorNumeric_SeqDense"
123 int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact)
124 {
125   Mat_SeqDense *mat = (Mat_SeqDense*) A->data, *l = (Mat_SeqDense*) (*fact)->data;
126   int          ierr;
127 
128   PetscFunctionBegin;
129   /* copy the numerical values */
130   PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
131   (*fact)->factor = 0;
132   ierr = MatLUFactor(*fact,0,0,1.0);CHKERRQ(ierr);
133   PetscFunctionReturn(0);
134 }
135 
136 #undef __FUNC__
137 #define __FUNC__ "MatCholeskyFactorSymbolic_SeqDense"
138 int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact)
139 {
140   int ierr;
141 
142   PetscFunctionBegin;
143   ierr = MatConvert(A,MATSAME,fact);CHKERRQ(ierr);
144   PetscFunctionReturn(0);
145 }
146 
147 #undef __FUNC__
148 #define __FUNC__ "MatCholeskyFactorNumeric_SeqDense"
149 int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact)
150 {
151   int ierr;
152 
153   PetscFunctionBegin;
154   ierr = MatCholeskyFactor(*fact,0,1.0);CHKERRQ(ierr);
155   PetscFunctionReturn(0);
156 }
157 
158 #undef __FUNC__
159 #define __FUNC__ "MatCholeskyFactor_SeqDense"
160 int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f)
161 {
162   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
163   int           info;
164 
165   PetscFunctionBegin;
166   if (mat->pivots) {
167     PetscFree(mat->pivots);
168     PLogObjectMemory(A,-mat->m*sizeof(int));
169     mat->pivots = 0;
170   }
171   LApotrf_("L",&mat->n,mat->v,&mat->m,&info);
172   if (info) SETERRQ(PETSC_ERR_MAT_CH_ZRPVT,0,"Bad factorization");
173   A->factor = FACTOR_CHOLESKY;
174   PLogFlops((mat->n*mat->n*mat->n)/3);
175   PetscFunctionReturn(0);
176 }
177 
178 #undef __FUNC__
179 #define __FUNC__ "MatSolve_SeqDense"
180 int MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
181 {
182   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
183   int          one = 1, info, ierr;
184   Scalar       *x, *y;
185 
186   PetscFunctionBegin;
187   ierr = VecGetArray(xx,&x); CHKERRQ(ierr); ierr = VecGetArray(yy,&y); CHKERRQ(ierr);
188   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
189   if (A->factor == FACTOR_LU) {
190     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
191   }
192   else if (A->factor == FACTOR_CHOLESKY){
193     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
194   }
195   else SETERRQ(1,0,"Matrix must be factored to solve");
196   if (info) SETERRQ(1,0,"MBad solve");
197   PLogFlops(mat->n*mat->n - mat->n);
198   PetscFunctionReturn(0);
199 }
200 
201 #undef __FUNC__
202 #define __FUNC__ "MatSolveTrans_SeqDense"
203 int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy)
204 {
205   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
206   int          one = 1, info;
207   Scalar       *x, *y;
208 
209   PetscFunctionBegin;
210   VecGetArray(xx,&x); VecGetArray(yy,&y);
211   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
212   /* assume if pivots exist then use LU; else Cholesky */
213   if (mat->pivots) {
214     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
215   }
216   else {
217     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
218   }
219   if (info) SETERRQ(1,0,"Bad solve");
220   PLogFlops(mat->n*mat->n - mat->n);
221   PetscFunctionReturn(0);
222 }
223 
224 #undef __FUNC__
225 #define __FUNC__ "MatSolveAdd_SeqDense"
226 int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
227 {
228   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
229   int          one = 1, info,ierr;
230   Scalar       *x, *y, sone = 1.0;
231   Vec          tmp = 0;
232 
233   PetscFunctionBegin;
234   VecGetArray(xx,&x); VecGetArray(yy,&y);
235   if (yy == zz) {
236     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
237     PLogObjectParent(A,tmp);
238     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
239   }
240   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
241   /* assume if pivots exist then use LU; else Cholesky */
242   if (mat->pivots) {
243     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
244   }
245   else {
246     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
247   }
248   if (info) SETERRQ(1,0,"Bad solve");
249   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
250   else VecAXPY(&sone,zz,yy);
251   PLogFlops(mat->n*mat->n - mat->n);
252   PetscFunctionReturn(0);
253 }
254 
255 #undef __FUNC__
256 #define __FUNC__ "MatSolveTransAdd_SeqDense"
257 int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy)
258 {
259   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
260   int           one = 1, info,ierr;
261   Scalar        *x, *y, sone = 1.0;
262   Vec           tmp;
263 
264   PetscFunctionBegin;
265   VecGetArray(xx,&x); VecGetArray(yy,&y);
266   if (yy == zz) {
267     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
268     PLogObjectParent(A,tmp);
269     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
270   }
271   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
272   /* assume if pivots exist then use LU; else Cholesky */
273   if (mat->pivots) {
274     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
275   } else {
276     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
277   }
278   if (info) SETERRQ(1,0,"Bad solve");
279   if (tmp) {
280     ierr = VecAXPY(&sone,tmp,yy);  CHKERRQ(ierr);
281     ierr = VecDestroy(tmp); CHKERRQ(ierr);
282   } else {
283     ierr = VecAXPY(&sone,zz,yy); CHKERRQ(ierr);
284   }
285   PLogFlops(mat->n*mat->n - mat->n);
286   PetscFunctionReturn(0);
287 }
288 /* ------------------------------------------------------------------*/
289 #undef __FUNC__
290 #define __FUNC__ "MatRelax_SeqDense"
291 int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag,
292                           double shift,int its,Vec xx)
293 {
294   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
295   Scalar       *x, *b, *v = mat->v, zero = 0.0, xt;
296   int          o = 1,ierr, m = mat->m, i;
297 
298   PetscFunctionBegin;
299   if (flag & SOR_ZERO_INITIAL_GUESS) {
300     /* this is a hack fix, should have another version without the second BLdot */
301     ierr = VecSet(&zero,xx); CHKERRQ(ierr);
302   }
303   VecGetArray(xx,&x); VecGetArray(bb,&b);
304   while (its--) {
305     if (flag & SOR_FORWARD_SWEEP){
306       for ( i=0; i<m; i++ ) {
307 #if defined(USE_PETSC_COMPLEX)
308         /* cannot use BLAS dot for complex because compiler/linker is
309            not happy about returning a double complex */
310         int    _i;
311         Scalar sum = b[i];
312         for ( _i=0; _i<m; _i++ ) {
313           sum -= conj(v[i+_i*m])*x[_i];
314         }
315         xt = sum;
316 #else
317         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
318 #endif
319         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
320       }
321     }
322     if (flag & SOR_BACKWARD_SWEEP) {
323       for ( i=m-1; i>=0; i-- ) {
324 #if defined(USE_PETSC_COMPLEX)
325         /* cannot use BLAS dot for complex because compiler/linker is
326            not happy about returning a double complex */
327         int    _i;
328         Scalar sum = b[i];
329         for ( _i=0; _i<m; _i++ ) {
330           sum -= conj(v[i+_i*m])*x[_i];
331         }
332         xt = sum;
333 #else
334         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
335 #endif
336         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
337       }
338     }
339   }
340   PetscFunctionReturn(0);
341 }
342 
343 /* -----------------------------------------------------------------*/
344 #undef __FUNC__
345 #define __FUNC__ "MatMultTrans_SeqDense"
346 int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy)
347 {
348   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
349   Scalar       *v = mat->v, *x, *y;
350   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
351 
352   PetscFunctionBegin;
353   VecGetArray(xx,&x), VecGetArray(yy,&y);
354   LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One);
355   PLogFlops(2*mat->m*mat->n - mat->n);
356   PetscFunctionReturn(0);
357 }
358 
359 #undef __FUNC__
360 #define __FUNC__ "MatMult_SeqDense"
361 int MatMult_SeqDense(Mat A,Vec xx,Vec yy)
362 {
363   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
364   Scalar       *v = mat->v, *x, *y;
365   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
366 
367   PetscFunctionBegin;
368   VecGetArray(xx,&x); VecGetArray(yy,&y);
369   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One);
370   PLogFlops(2*mat->m*mat->n - mat->m);
371   PetscFunctionReturn(0);
372 }
373 
374 #undef __FUNC__
375 #define __FUNC__ "MatMultAdd_SeqDense"
376 int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
377 {
378   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
379   Scalar       *v = mat->v, *x, *y, *z;
380   int          _One=1; Scalar _DOne=1.0;
381 
382   PetscFunctionBegin;
383   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
384   if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar));
385   LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One);
386   PLogFlops(2*mat->m*mat->n);
387   PetscFunctionReturn(0);
388 }
389 
390 #undef __FUNC__
391 #define __FUNC__ "MatMultTransAdd_SeqDense"
392 int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
393 {
394   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
395   Scalar       *v = mat->v, *x, *y, *z;
396   int          _One=1;
397   Scalar       _DOne=1.0;
398 
399   PetscFunctionBegin;
400   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
401   if (zz != yy) PetscMemcpy(y,z,mat->n*sizeof(Scalar));
402   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One);
403   PLogFlops(2*mat->m*mat->n);
404   PetscFunctionReturn(0);
405 }
406 
407 /* -----------------------------------------------------------------*/
408 #undef __FUNC__
409 #define __FUNC__ "MatGetRow_SeqDense"
410 int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
411 {
412   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
413   Scalar       *v;
414   int          i;
415 
416   PetscFunctionBegin;
417   *ncols = mat->n;
418   if (cols) {
419     *cols = (int *) PetscMalloc((mat->n+1)*sizeof(int)); CHKPTRQ(*cols);
420     for ( i=0; i<mat->n; i++ ) (*cols)[i] = i;
421   }
422   if (vals) {
423     *vals = (Scalar *) PetscMalloc((mat->n+1)*sizeof(Scalar)); CHKPTRQ(*vals);
424     v = mat->v + row;
425     for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;}
426   }
427   PetscFunctionReturn(0);
428 }
429 
430 #undef __FUNC__
431 #define __FUNC__ "MatRestoreRow_SeqDense"
432 int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
433 {
434   if (cols) { PetscFree(*cols); }
435   if (vals) { PetscFree(*vals); }
436   PetscFunctionReturn(0);
437 }
438 /* ----------------------------------------------------------------*/
439 #undef __FUNC__
440 #define __FUNC__ "MatSetValues_SeqDense"
441 int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n,
442                                     int *indexn,Scalar *v,InsertMode addv)
443 {
444   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
445   int          i,j;
446 
447   PetscFunctionBegin;
448   if (!mat->roworiented) {
449     if (addv == INSERT_VALUES) {
450       for ( j=0; j<n; j++ ) {
451 #if defined(USE_PETSC_BOPT_g)
452         if (indexn[j] < 0) SETERRQ(1,0,"Negative column");
453         if (indexn[j] >= A->n) SETERRQ(1,0,"Column too large");
454 #endif
455         for ( i=0; i<m; i++ ) {
456 #if defined(USE_PETSC_BOPT_g)
457           if (indexm[i] < 0) SETERRQ(1,0,"Negative row");
458           if (indexm[i] >= A->m) SETERRQ(1,0,"Row too large");
459 #endif
460           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
461         }
462       }
463     } else {
464       for ( j=0; j<n; j++ ) {
465 #if defined(USE_PETSC_BOPT_g)
466         if (indexn[j] < 0) SETERRQ(1,0,"Negative column");
467         if (indexn[j] >= A->n) SETERRQ(1,0,"Column too large");
468 #endif
469         for ( i=0; i<m; i++ ) {
470 #if defined(USE_PETSC_BOPT_g)
471           if (indexm[i] < 0) SETERRQ(1,0,"Negative row");
472           if (indexm[i] >= A->m) SETERRQ(1,0,"Row too large");
473 #endif
474           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
475         }
476       }
477     }
478   } else {
479     if (addv == INSERT_VALUES) {
480       for ( i=0; i<m; i++ ) {
481 #if defined(USE_PETSC_BOPT_g)
482         if (indexm[i] < 0) SETERRQ(1,0,"Negative row");
483         if (indexm[i] >= A->m) SETERRQ(1,0,"Row too large");
484 #endif
485         for ( j=0; j<n; j++ ) {
486 #if defined(USE_PETSC_BOPT_g)
487           if (indexn[j] < 0) SETERRQ(1,0,"Negative column");
488           if (indexn[j] >= A->n) SETERRQ(1,0,"Column too large");
489 #endif
490           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
491         }
492       }
493     } else {
494       for ( i=0; i<m; i++ ) {
495 #if defined(USE_PETSC_BOPT_g)
496         if (indexm[i] < 0) SETERRQ(1,0,"Negative row");
497         if (indexm[i] >= A->m) SETERRQ(1,0,"Row too large");
498 #endif
499         for ( j=0; j<n; j++ ) {
500 #if defined(USE_PETSC_BOPT_g)
501           if (indexn[j] < 0) SETERRQ(1,0,"Negative column");
502           if (indexn[j] >= A->n) SETERRQ(1,0,"Column too large");
503 #endif
504           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
505         }
506       }
507     }
508   }
509   PetscFunctionReturn(0);
510 }
511 
512 #undef __FUNC__
513 #define __FUNC__ "MatGetValues_SeqDense"
514 int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v)
515 {
516   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
517   int          i, j;
518   Scalar       *vpt = v;
519 
520   PetscFunctionBegin;
521   /* row-oriented output */
522   for ( i=0; i<m; i++ ) {
523     for ( j=0; j<n; j++ ) {
524       *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]];
525     }
526   }
527   PetscFunctionReturn(0);
528 }
529 
530 /* -----------------------------------------------------------------*/
531 
532 #include "sys.h"
533 
534 #undef __FUNC__
535 #define __FUNC__ "MatLoad_SeqDense"
536 int MatLoad_SeqDense(Viewer viewer,MatType type,Mat *A)
537 {
538   Mat_SeqDense *a;
539   Mat          B;
540   int          *scols, i, j, nz, ierr, fd, header[4], size;
541   int          *rowlengths = 0, M, N, *cols;
542   Scalar       *vals, *svals, *v;
543   MPI_Comm     comm = ((PetscObject)viewer)->comm;
544 
545   PetscFunctionBegin;
546   MPI_Comm_size(comm,&size);
547   if (size > 1) SETERRQ(1,0,"view must have one processor");
548   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
549   ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr);
550   if (header[0] != MAT_COOKIE) SETERRQ(1,0,"Not matrix object");
551   M = header[1]; N = header[2]; nz = header[3];
552 
553   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
554     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
555     B    = *A;
556     a    = (Mat_SeqDense *) B->data;
557 
558     /* read in nonzero values */
559     ierr = PetscBinaryRead(fd,a->v,M*N,PETSC_SCALAR); CHKERRQ(ierr);
560 
561     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
562     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
563     /* ierr = MatTranspose(B,PETSC_NULL); CHKERRQ(ierr); */
564   } else {
565     /* read row lengths */
566     rowlengths = (int*) PetscMalloc( (M+1)*sizeof(int) ); CHKPTRQ(rowlengths);
567     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr);
568 
569     /* create our matrix */
570     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
571     B = *A;
572     a = (Mat_SeqDense *) B->data;
573     v = a->v;
574 
575     /* read column indices and nonzeros */
576     cols = scols = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(cols);
577     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT); CHKERRQ(ierr);
578     vals = svals = (Scalar *) PetscMalloc( (nz+1)*sizeof(Scalar) ); CHKPTRQ(vals);
579     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr);
580 
581     /* insert into matrix */
582     for ( i=0; i<M; i++ ) {
583       for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j];
584       svals += rowlengths[i]; scols += rowlengths[i];
585     }
586     PetscFree(vals); PetscFree(cols); PetscFree(rowlengths);
587 
588     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
589     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
590   }
591   PetscFunctionReturn(0);
592 }
593 
594 #include "pinclude/pviewer.h"
595 #include "sys.h"
596 
597 #undef __FUNC__
598 #define __FUNC__ "MatView_SeqDense_ASCII"
599 static int MatView_SeqDense_ASCII(Mat A,Viewer viewer)
600 {
601   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
602   int          ierr, i, j, format;
603   FILE         *fd;
604   char         *outputname;
605   Scalar       *v;
606 
607   PetscFunctionBegin;
608   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
609   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
610   ierr = ViewerGetFormat(viewer,&format);
611   if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) {
612     PetscFunctionReturn(0);  /* do nothing for now */
613   }
614   else if (format == VIEWER_FORMAT_ASCII_COMMON) {
615     for ( i=0; i<a->m; i++ ) {
616       v = a->v + i;
617       fprintf(fd,"row %d:",i);
618       for ( j=0; j<a->n; j++ ) {
619 #if defined(USE_PETSC_COMPLEX)
620         if (real(*v) != 0.0 && imag(*v) != 0.0) fprintf(fd," %d %g + %g i",j,real(*v),imag(*v));
621         else if (real(*v)) fprintf(fd," %d %g ",j,real(*v));
622 #else
623         if (*v) fprintf(fd," %d %g ",j,*v);
624 #endif
625         v += a->m;
626       }
627       fprintf(fd,"\n");
628     }
629   } else {
630 #if defined(USE_PETSC_COMPLEX)
631     int allreal = 1;
632     /* determine if matrix has all real values */
633     v = a->v;
634     for ( i=0; i<a->m*a->n; i++ ) {
635       if (imag(v[i])) { allreal = 0; break ;}
636     }
637 #endif
638     for ( i=0; i<a->m; i++ ) {
639       v = a->v + i;
640       for ( j=0; j<a->n; j++ ) {
641 #if defined(USE_PETSC_COMPLEX)
642         if (allreal) {
643           fprintf(fd,"%6.4e ",real(*v)); v += a->m;
644         } else {
645           fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m;
646         }
647 #else
648         fprintf(fd,"%6.4e ",*v); v += a->m;
649 #endif
650       }
651       fprintf(fd,"\n");
652     }
653   }
654   fflush(fd);
655   PetscFunctionReturn(0);
656 }
657 
658 #undef __FUNC__
659 #define __FUNC__ "MatView_SeqDense_Binary"
660 static int MatView_SeqDense_Binary(Mat A,Viewer viewer)
661 {
662   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
663   int          ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n;
664   int          format;
665   Scalar       *v, *anonz,*vals;
666 
667   PetscFunctionBegin;
668   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
669 
670   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
671   if (format == VIEWER_FORMAT_BINARY_NATIVE) {
672     /* store the matrix as a dense matrix */
673     col_lens = (int *) PetscMalloc( 4*sizeof(int) ); CHKPTRQ(col_lens);
674     col_lens[0] = MAT_COOKIE;
675     col_lens[1] = m;
676     col_lens[2] = n;
677     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
678     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1); CHKERRQ(ierr);
679     PetscFree(col_lens);
680 
681     /* write out matrix, by rows */
682     vals = (Scalar *) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(vals);
683     v    = a->v;
684     for ( i=0; i<m; i++ ) {
685       for ( j=0; j<n; j++ ) {
686         vals[i + j*m] = *v++;
687       }
688     }
689     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0); CHKERRQ(ierr);
690     PetscFree(vals);
691   } else {
692     col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens);
693     col_lens[0] = MAT_COOKIE;
694     col_lens[1] = m;
695     col_lens[2] = n;
696     col_lens[3] = nz;
697 
698     /* store lengths of each row and write (including header) to file */
699     for ( i=0; i<m; i++ ) col_lens[4+i] = n;
700     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1); CHKERRQ(ierr);
701 
702     /* Possibly should write in smaller increments, not whole matrix at once? */
703     /* store column indices (zero start index) */
704     ict = 0;
705     for ( i=0; i<m; i++ ) {
706       for ( j=0; j<n; j++ ) col_lens[ict++] = j;
707     }
708     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0); CHKERRQ(ierr);
709     PetscFree(col_lens);
710 
711     /* store nonzero values */
712     anonz = (Scalar *) PetscMalloc((nz+1)*sizeof(Scalar)); CHKPTRQ(anonz);
713     ict = 0;
714     for ( i=0; i<m; i++ ) {
715       v = a->v + i;
716       for ( j=0; j<n; j++ ) {
717         anonz[ict++] = *v; v += a->m;
718       }
719     }
720     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0); CHKERRQ(ierr);
721     PetscFree(anonz);
722   }
723   PetscFunctionReturn(0);
724 }
725 
726 #undef __FUNC__
727 #define __FUNC__ "MatView_SeqDense"
728 int MatView_SeqDense(PetscObject obj,Viewer viewer)
729 {
730   Mat          A = (Mat) obj;
731   Mat_SeqDense *a = (Mat_SeqDense*) A->data;
732   ViewerType   vtype;
733   int          ierr;
734 
735   PetscFunctionBegin;
736   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
737 
738   if (vtype == MATLAB_VIEWER) {
739     ierr = ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v); CHKERRQ(ierr);
740   } else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
741     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
742   } else if (vtype == BINARY_FILE_VIEWER) {
743     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
744   }
745   PetscFunctionReturn(0);
746 }
747 
748 #undef __FUNC__
749 #define __FUNC__ "MatDestroy_SeqDense"
750 int MatDestroy_SeqDense(PetscObject obj)
751 {
752   Mat          mat = (Mat) obj;
753   Mat_SeqDense *l = (Mat_SeqDense *) mat->data;
754   int          ierr;
755 
756   PetscFunctionBegin;
757 #if defined(USE_PETSC_LOG)
758   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
759 #endif
760   if (l->pivots) PetscFree(l->pivots);
761   if (!l->user_alloc) PetscFree(l->v);
762   PetscFree(l);
763   if (mat->mapping) {
764     ierr = ISLocalToGlobalMappingDestroy(mat->mapping); CHKERRQ(ierr);
765   }
766   PLogObjectDestroy(mat);
767   PetscHeaderDestroy(mat);
768   PetscFunctionReturn(0);
769 }
770 
771 #undef __FUNC__
772 #define __FUNC__ "MatTranspose_SeqDense"
773 int MatTranspose_SeqDense(Mat A,Mat *matout)
774 {
775   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
776   int          k, j, m, n;
777   Scalar       *v, tmp;
778 
779   PetscFunctionBegin;
780   v = mat->v; m = mat->m; n = mat->n;
781   if (matout == PETSC_NULL) { /* in place transpose */
782     if (m != n) { /* malloc temp to hold transpose */
783       Scalar *w = (Scalar *) PetscMalloc((m+1)*(n+1)*sizeof(Scalar));CHKPTRQ(w);
784       for ( j=0; j<m; j++ ) {
785         for ( k=0; k<n; k++ ) {
786           w[k + j*n] = v[j + k*m];
787         }
788       }
789       PetscMemcpy(v,w,m*n*sizeof(Scalar));
790       PetscFree(w);
791     } else {
792       for ( j=0; j<m; j++ ) {
793         for ( k=0; k<j; k++ ) {
794           tmp = v[j + k*n];
795           v[j + k*n] = v[k + j*n];
796           v[k + j*n] = tmp;
797         }
798       }
799     }
800   } else { /* out-of-place transpose */
801     int          ierr;
802     Mat          tmat;
803     Mat_SeqDense *tmatd;
804     Scalar       *v2;
805     ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr);
806     tmatd = (Mat_SeqDense *) tmat->data;
807     v = mat->v; v2 = tmatd->v;
808     for ( j=0; j<n; j++ ) {
809       for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m];
810     }
811     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
812     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
813     *matout = tmat;
814   }
815   PetscFunctionReturn(0);
816 }
817 
818 #undef __FUNC__
819 #define __FUNC__ "MatEqual_SeqDense"
820 int MatEqual_SeqDense(Mat A1,Mat A2, PetscTruth *flg)
821 {
822   Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data;
823   Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data;
824   int          i;
825   Scalar       *v1 = mat1->v, *v2 = mat2->v;
826 
827   PetscFunctionBegin;
828   if (A2->type != MATSEQDENSE) SETERRQ(PETSC_ERR_SUP,0,"Matrices should be of type  MATSEQDENSE");
829   if (mat1->m != mat2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
830   if (mat1->n != mat2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
831   for ( i=0; i<mat1->m*mat1->n; i++ ) {
832     if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
833     v1++; v2++;
834   }
835   *flg = PETSC_TRUE;
836   PetscFunctionReturn(0);
837 }
838 
839 #undef __FUNC__
840 #define __FUNC__ "MatGetDiagonal_SeqDense"
841 int MatGetDiagonal_SeqDense(Mat A,Vec v)
842 {
843   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
844   int          i, n, len;
845   Scalar       *x, zero = 0.0;
846 
847   PetscFunctionBegin;
848   VecSet(&zero,v);
849   VecGetArray(v,&x); VecGetSize(v,&n);
850   len = PetscMin(mat->m,mat->n);
851   if (n != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming mat and vec");
852   for ( i=0; i<len; i++ ) {
853     x[i] = mat->v[i*mat->m + i];
854   }
855   PetscFunctionReturn(0);
856 }
857 
858 #undef __FUNC__
859 #define __FUNC__ "MatDiagonalScale_SeqDense"
860 int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
861 {
862   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
863   Scalar       *l,*r,x,*v;
864   int          i,j,m = mat->m, n = mat->n;
865 
866   PetscFunctionBegin;
867   if (ll) {
868     VecGetArray(ll,&l); VecGetSize(ll,&m);
869     if (m != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vec wrong size");
870     for ( i=0; i<m; i++ ) {
871       x = l[i];
872       v = mat->v + i;
873       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
874     }
875     PLogFlops(n*m);
876   }
877   if (rr) {
878     VecGetArray(rr,&r); VecGetSize(rr,&n);
879     if (n != mat->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vec wrong size");
880     for ( i=0; i<n; i++ ) {
881       x = r[i];
882       v = mat->v + i*m;
883       for ( j=0; j<m; j++ ) { (*v++) *= x;}
884     }
885     PLogFlops(n*m);
886   }
887   PetscFunctionReturn(0);
888 }
889 
890 #undef __FUNC__
891 #define __FUNC__ "MatNorm_SeqDense"
892 int MatNorm_SeqDense(Mat A,NormType type,double *norm)
893 {
894   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
895   Scalar       *v = mat->v;
896   double       sum = 0.0;
897   int          i, j;
898 
899   PetscFunctionBegin;
900   if (type == NORM_FROBENIUS) {
901     for (i=0; i<mat->n*mat->m; i++ ) {
902 #if defined(USE_PETSC_COMPLEX)
903       sum += real(conj(*v)*(*v)); v++;
904 #else
905       sum += (*v)*(*v); v++;
906 #endif
907     }
908     *norm = sqrt(sum);
909     PLogFlops(2*mat->n*mat->m);
910   } else if (type == NORM_1) {
911     *norm = 0.0;
912     for ( j=0; j<mat->n; j++ ) {
913       sum = 0.0;
914       for ( i=0; i<mat->m; i++ ) {
915         sum += PetscAbsScalar(*v);  v++;
916       }
917       if (sum > *norm) *norm = sum;
918     }
919     PLogFlops(mat->n*mat->m);
920   } else if (type == NORM_INFINITY) {
921     *norm = 0.0;
922     for ( j=0; j<mat->m; j++ ) {
923       v = mat->v + j;
924       sum = 0.0;
925       for ( i=0; i<mat->n; i++ ) {
926         sum += PetscAbsScalar(*v); v += mat->m;
927       }
928       if (sum > *norm) *norm = sum;
929     }
930     PLogFlops(mat->n*mat->m);
931   } else {
932     SETERRQ(PETSC_ERR_SUP,0,"No two norm");
933   }
934   PetscFunctionReturn(0);
935 }
936 
937 #undef __FUNC__
938 #define __FUNC__ "MatSetOption_SeqDense"
939 int MatSetOption_SeqDense(Mat A,MatOption op)
940 {
941   Mat_SeqDense *aij = (Mat_SeqDense *) A->data;
942 
943   PetscFunctionBegin;
944   if (op == MAT_ROW_ORIENTED)            aij->roworiented = 1;
945   else if (op == MAT_COLUMN_ORIENTED)    aij->roworiented = 0;
946   else if (op == MAT_ROWS_SORTED ||
947            op == MAT_ROWS_UNSORTED ||
948            op == MAT_COLUMNS_SORTED ||
949            op == MAT_COLUMNS_UNSORTED ||
950            op == MAT_SYMMETRIC ||
951            op == MAT_STRUCTURALLY_SYMMETRIC ||
952            op == MAT_NO_NEW_NONZERO_LOCATIONS ||
953            op == MAT_YES_NEW_NONZERO_LOCATIONS ||
954            op == MAT_NEW_NONZERO_LOCATION_ERROR ||
955            op == MAT_NO_NEW_DIAGONALS ||
956            op == MAT_YES_NEW_DIAGONALS ||
957            op == MAT_IGNORE_OFF_PROC_ENTRIES)
958     PLogInfo(A,"Info:MatSetOption_SeqDense:Option ignored\n");
959   else {
960     SETERRQ(PETSC_ERR_SUP,0,"unknown option");
961   }
962   PetscFunctionReturn(0);
963 }
964 
965 #undef __FUNC__
966 #define __FUNC__ "MatZeroEntries_SeqDense"
967 int MatZeroEntries_SeqDense(Mat A)
968 {
969   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
970 
971   PetscFunctionBegin;
972   PetscMemzero(l->v,l->m*l->n*sizeof(Scalar));
973   PetscFunctionReturn(0);
974 }
975 
976 #undef __FUNC__
977 #define __FUNC__ "MatGetBlockSize_SeqDense"
978 int MatGetBlockSize_SeqDense(Mat A,int *bs)
979 {
980   PetscFunctionBegin;
981   *bs = 1;
982   PetscFunctionReturn(0);
983 }
984 
985 #undef __FUNC__
986 #define __FUNC__ "MatZeroRows_SeqDense"
987 int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
988 {
989   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
990   int          n = l->n, i, j,ierr,N, *rows;
991   Scalar       *slot;
992 
993   PetscFunctionBegin;
994   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
995   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
996   for ( i=0; i<N; i++ ) {
997     slot = l->v + rows[i];
998     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
999   }
1000   if (diag) {
1001     for ( i=0; i<N; i++ ) {
1002       slot = l->v + (n+1)*rows[i];
1003       *slot = *diag;
1004     }
1005   }
1006   ISRestoreIndices(is,&rows);
1007   PetscFunctionReturn(0);
1008 }
1009 
1010 #undef __FUNC__
1011 #define __FUNC__ "MatGetSize_SeqDense"
1012 int MatGetSize_SeqDense(Mat A,int *m,int *n)
1013 {
1014   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1015 
1016   PetscFunctionBegin;
1017   *m = mat->m; *n = mat->n;
1018   PetscFunctionReturn(0);
1019 }
1020 
1021 #undef __FUNC__
1022 #define __FUNC__ "MatGetOwnershipRange_SeqDense"
1023 int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
1024 {
1025   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1026 
1027   PetscFunctionBegin;
1028   *m = 0; *n = mat->m;
1029   PetscFunctionReturn(0);
1030 }
1031 
1032 #undef __FUNC__
1033 #define __FUNC__ "MatGetArray_SeqDense"
1034 int MatGetArray_SeqDense(Mat A,Scalar **array)
1035 {
1036   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1037 
1038   PetscFunctionBegin;
1039   *array = mat->v;
1040   PetscFunctionReturn(0);
1041 }
1042 
1043 #undef __FUNC__
1044 #define __FUNC__ "MatRestoreArray_SeqDense"
1045 int MatRestoreArray_SeqDense(Mat A,Scalar **array)
1046 {
1047   PetscFunctionBegin;
1048   PetscFunctionReturn(0);
1049 }
1050 
1051 #undef __FUNC__
1052 #define __FUNC__ "MatGetSubMatrix_SeqDense"
1053 static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,
1054                                     Mat *submat)
1055 {
1056   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1057   int          nznew, *smap, i, j, ierr, oldcols = mat->n;
1058   int          *irow, *icol, nrows, ncols, *cwork;
1059   Scalar       *vwork, *val;
1060   Mat          newmat;
1061 
1062   PetscFunctionBegin;
1063   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
1064   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
1065   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
1066   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
1067 
1068   smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap);
1069   cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork);
1070   vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
1071   PetscMemzero((char*)smap,oldcols*sizeof(int));
1072   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
1073 
1074   /* Create and fill new matrix */
1075   ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr);
1076   for (i=0; i<nrows; i++) {
1077     nznew = 0;
1078     val   = mat->v + irow[i];
1079     for (j=0; j<oldcols; j++) {
1080       if (smap[j]) {
1081         cwork[nznew]   = smap[j] - 1;
1082         vwork[nznew++] = val[j * mat->m];
1083       }
1084     }
1085     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES);CHKERRQ(ierr);
1086   }
1087   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1088   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1089 
1090   /* Free work space */
1091   PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
1092   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1093   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
1094   *submat = newmat;
1095   PetscFunctionReturn(0);
1096 }
1097 
1098 #undef __FUNC__
1099 #define __FUNC__ "MatGetSubMatrices_SeqDense"
1100 int MatGetSubMatrices_SeqDense(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
1101                                     Mat **B)
1102 {
1103   int ierr,i;
1104 
1105   PetscFunctionBegin;
1106   if (scall == MAT_INITIAL_MATRIX) {
1107     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1108   }
1109 
1110   for ( i=0; i<n; i++ ) {
1111     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr);
1112   }
1113   PetscFunctionReturn(0);
1114 }
1115 
1116 #undef __FUNC__
1117 #define __FUNC__ "MatCopy_SeqDense"
1118 int MatCopy_SeqDense(Mat A, Mat B)
1119 {
1120   Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data;
1121   int          ierr;
1122 
1123   PetscFunctionBegin;
1124   if (B->type != MATSEQDENSE) {
1125     ierr = MatCopy_Basic(A,B);CHKERRQ(ierr);
1126     PetscFunctionReturn(0);
1127   }
1128   if (a->m != b->m || a->n != b->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"size(B) != size(A)");
1129   PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar));
1130   PetscFunctionReturn(0);
1131 }
1132 
1133 /* -------------------------------------------------------------------*/
1134 static struct _MatOps MatOps = {MatSetValues_SeqDense,
1135        MatGetRow_SeqDense,
1136        MatRestoreRow_SeqDense,
1137        MatMult_SeqDense,
1138        MatMultAdd_SeqDense,
1139        MatMultTrans_SeqDense,
1140        MatMultTransAdd_SeqDense,
1141        MatSolve_SeqDense,
1142        MatSolveAdd_SeqDense,
1143        MatSolveTrans_SeqDense,
1144        MatSolveTransAdd_SeqDense,
1145        MatLUFactor_SeqDense,
1146        MatCholeskyFactor_SeqDense,
1147        MatRelax_SeqDense,
1148        MatTranspose_SeqDense,
1149        MatGetInfo_SeqDense,
1150        MatEqual_SeqDense,
1151        MatGetDiagonal_SeqDense,
1152        MatDiagonalScale_SeqDense,
1153        MatNorm_SeqDense,
1154        0,
1155        0,
1156        0,
1157        MatSetOption_SeqDense,
1158        MatZeroEntries_SeqDense,
1159        MatZeroRows_SeqDense,
1160        MatLUFactorSymbolic_SeqDense,
1161        MatLUFactorNumeric_SeqDense,
1162        MatCholeskyFactorSymbolic_SeqDense,
1163        MatCholeskyFactorNumeric_SeqDense,
1164        MatGetSize_SeqDense,
1165        MatGetSize_SeqDense,
1166        MatGetOwnershipRange_SeqDense,
1167        0,
1168        0,
1169        MatGetArray_SeqDense,
1170        MatRestoreArray_SeqDense,
1171        MatConvertSameType_SeqDense,0,0,0,0,
1172        MatAXPY_SeqDense,MatGetSubMatrices_SeqDense,0,
1173        MatGetValues_SeqDense,
1174        MatCopy_SeqDense,0,MatScale_SeqDense,
1175        0,0,0,MatGetBlockSize_SeqDense};
1176 
1177 #undef __FUNC__
1178 #define __FUNC__ "MatCreateSeqDense"
1179 /*@C
1180    MatCreateSeqDense - Creates a sequential dense matrix that
1181    is stored in column major order (the usual Fortran 77 manner). Many
1182    of the matrix operations use the BLAS and LAPACK routines.
1183 
1184    Input Parameters:
1185 .  comm - MPI communicator, set to PETSC_COMM_SELF
1186 .  m - number of rows
1187 .  n - number of columns
1188 .  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1189    to control all matrix memory allocation.
1190 
1191    Output Parameter:
1192 .  A - the matrix
1193 
1194   Notes:
1195   The data input variable is intended primarily for Fortran programmers
1196   who wish to allocate their own matrix memory space.  Most users should
1197   set data=PETSC_NULL.
1198 
1199 .keywords: dense, matrix, LAPACK, BLAS
1200 
1201 .seealso: MatCreate(), MatSetValues()
1202 @*/
1203 int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A)
1204 {
1205   Mat          B;
1206   Mat_SeqDense *b;
1207   int          ierr,flg,size;
1208 
1209   PetscFunctionBegin;
1210   MPI_Comm_size(comm,&size);
1211   if (size > 1) SETERRQ(1,0,"Comm must be of size 1");
1212 
1213   *A            = 0;
1214   PetscHeaderCreate(B,_p_Mat,MAT_COOKIE,MATSEQDENSE,comm,MatDestroy,MatView);
1215   PLogObjectCreate(B);
1216   b             = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(b);
1217   PetscMemzero(b,sizeof(Mat_SeqDense));
1218   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1219   B->destroy    = MatDestroy_SeqDense;
1220   B->view       = MatView_SeqDense;
1221   B->factor     = 0;
1222   B->mapping    = 0;
1223   PLogObjectMemory(B,sizeof(struct _p_Mat));
1224   B->data       = (void *) b;
1225 
1226   b->m = m;  B->m = m; B->M = m;
1227   b->n = n;  B->n = n; B->N = n;
1228   b->pivots       = 0;
1229   b->roworiented  = 1;
1230   if (data == PETSC_NULL) {
1231     b->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(b->v);
1232     PetscMemzero(b->v,m*n*sizeof(Scalar));
1233     b->user_alloc = 0;
1234     PLogObjectMemory(B,n*m*sizeof(Scalar));
1235   }
1236   else { /* user-allocated storage */
1237     b->v = data;
1238     b->user_alloc = 1;
1239   }
1240   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
1241   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr);}
1242 
1243   *A = B;
1244   PetscFunctionReturn(0);
1245 }
1246 
1247 
1248 
1249 
1250 
1251 
1252