xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 3a40ed3dce77c081171d005ae1a6ff4bb9d13b6f) !
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: dense.c,v 1.131 1997/09/26 02:18:55 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)
621           fprintf(fd," %d %g + %g i",j,real(*v),imag(*v));
622         else if (real(*v)) fprintf(fd," %d %g ",j,real(*v));
623         v += a->m;
624 #else
625         if (*v) fprintf(fd," %d %g ",j,*v);
626         v += a->m;
627 #endif
628       }
629       fprintf(fd,"\n");
630     }
631   } else {
632 #if defined(USE_PETSC_COMPLEX)
633     int allreal = 1;
634     /* determine if matrix has all real values */
635     v = a->v;
636     for ( i=0; i<a->m*a->n; i++ ) {
637       if (imag(v[i])) { allreal = 0; break ;}
638     }
639 #endif
640     for ( i=0; i<a->m; i++ ) {
641       v = a->v + i;
642       for ( j=0; j<a->n; j++ ) {
643 #if defined(USE_PETSC_COMPLEX)
644         if (allreal) {
645           fprintf(fd,"%6.4e ",real(*v)); v += a->m;
646         } else {
647           fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m;
648         }
649 #else
650         fprintf(fd,"%6.4e ",*v); v += a->m;
651 #endif
652       }
653       fprintf(fd,"\n");
654     }
655   }
656   fflush(fd);
657   PetscFunctionReturn(0);
658 }
659 
660 #undef __FUNC__
661 #define __FUNC__ "MatView_SeqDense_Binary"
662 static int MatView_SeqDense_Binary(Mat A,Viewer viewer)
663 {
664   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
665   int          ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n;
666   int          format;
667   Scalar       *v, *anonz,*vals;
668 
669   PetscFunctionBegin;
670   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
671 
672   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
673   if (format == VIEWER_FORMAT_BINARY_NATIVE) {
674     /* store the matrix as a dense matrix */
675     col_lens = (int *) PetscMalloc( 4*sizeof(int) ); CHKPTRQ(col_lens);
676     col_lens[0] = MAT_COOKIE;
677     col_lens[1] = m;
678     col_lens[2] = n;
679     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
680     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1); CHKERRQ(ierr);
681     PetscFree(col_lens);
682 
683     /* write out matrix, by rows */
684     vals = (Scalar *) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(vals);
685     v    = a->v;
686     for ( i=0; i<m; i++ ) {
687       for ( j=0; j<n; j++ ) {
688         vals[i + j*m] = *v++;
689       }
690     }
691     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0); CHKERRQ(ierr);
692     PetscFree(vals);
693   } else {
694     col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens);
695     col_lens[0] = MAT_COOKIE;
696     col_lens[1] = m;
697     col_lens[2] = n;
698     col_lens[3] = nz;
699 
700     /* store lengths of each row and write (including header) to file */
701     for ( i=0; i<m; i++ ) col_lens[4+i] = n;
702     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1); CHKERRQ(ierr);
703 
704     /* Possibly should write in smaller increments, not whole matrix at once? */
705     /* store column indices (zero start index) */
706     ict = 0;
707     for ( i=0; i<m; i++ ) {
708       for ( j=0; j<n; j++ ) col_lens[ict++] = j;
709     }
710     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0); CHKERRQ(ierr);
711     PetscFree(col_lens);
712 
713     /* store nonzero values */
714     anonz = (Scalar *) PetscMalloc((nz+1)*sizeof(Scalar)); CHKPTRQ(anonz);
715     ict = 0;
716     for ( i=0; i<m; i++ ) {
717       v = a->v + i;
718       for ( j=0; j<n; j++ ) {
719         anonz[ict++] = *v; v += a->m;
720       }
721     }
722     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0); CHKERRQ(ierr);
723     PetscFree(anonz);
724   }
725   PetscFunctionReturn(0);
726 }
727 
728 #undef __FUNC__
729 #define __FUNC__ "MatView_SeqDense"
730 int MatView_SeqDense(PetscObject obj,Viewer viewer)
731 {
732   Mat          A = (Mat) obj;
733   Mat_SeqDense *a = (Mat_SeqDense*) A->data;
734   ViewerType   vtype;
735   int          ierr;
736 
737   PetscFunctionBegin;
738   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
739 
740   if (vtype == MATLAB_VIEWER) {
741     ierr = ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v); CHKERRQ(ierr);
742   } else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
743     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
744   } else if (vtype == BINARY_FILE_VIEWER) {
745     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
746   }
747   PetscFunctionReturn(0);
748 }
749 
750 #undef __FUNC__
751 #define __FUNC__ "MatDestroy_SeqDense"
752 int MatDestroy_SeqDense(PetscObject obj)
753 {
754   Mat          mat = (Mat) obj;
755   Mat_SeqDense *l = (Mat_SeqDense *) mat->data;
756   int          ierr;
757 
758   PetscFunctionBegin;
759 #if defined(USE_PETSC_LOG)
760   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
761 #endif
762   if (l->pivots) PetscFree(l->pivots);
763   if (!l->user_alloc) PetscFree(l->v);
764   PetscFree(l);
765   if (mat->mapping) {
766     ierr = ISLocalToGlobalMappingDestroy(mat->mapping); CHKERRQ(ierr);
767   }
768   PLogObjectDestroy(mat);
769   PetscHeaderDestroy(mat);
770   PetscFunctionReturn(0);
771 }
772 
773 #undef __FUNC__
774 #define __FUNC__ "MatTranspose_SeqDense"
775 int MatTranspose_SeqDense(Mat A,Mat *matout)
776 {
777   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
778   int          k, j, m, n;
779   Scalar       *v, tmp;
780 
781   PetscFunctionBegin;
782   v = mat->v; m = mat->m; n = mat->n;
783   if (matout == PETSC_NULL) { /* in place transpose */
784     if (m != n) SETERRQ(PETSC_ERR_SUP,0,"Supports square matrix only in-place");
785     for ( j=0; j<m; j++ ) {
786       for ( k=0; k<j; k++ ) {
787         tmp = v[j + k*n];
788         v[j + k*n] = v[k + j*n];
789         v[k + j*n] = tmp;
790       }
791     }
792   } else { /* out-of-place transpose */
793     int          ierr;
794     Mat          tmat;
795     Mat_SeqDense *tmatd;
796     Scalar       *v2;
797     ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr);
798     tmatd = (Mat_SeqDense *) tmat->data;
799     v = mat->v; v2 = tmatd->v;
800     for ( j=0; j<n; j++ ) {
801       for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m];
802     }
803     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
804     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
805     *matout = tmat;
806   }
807   PetscFunctionReturn(0);
808 }
809 
810 #undef __FUNC__
811 #define __FUNC__ "MatEqual_SeqDense"
812 int MatEqual_SeqDense(Mat A1,Mat A2, PetscTruth *flg)
813 {
814   Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data;
815   Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data;
816   int          i;
817   Scalar       *v1 = mat1->v, *v2 = mat2->v;
818 
819   PetscFunctionBegin;
820   if (A2->type != MATSEQDENSE) SETERRQ(PETSC_ERR_SUP,0,"Matrices should be of type  MATSEQDENSE");
821   if (mat1->m != mat2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
822   if (mat1->n != mat2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
823   for ( i=0; i<mat1->m*mat1->n; i++ ) {
824     if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
825     v1++; v2++;
826   }
827   *flg = PETSC_TRUE;
828   PetscFunctionReturn(0);
829 }
830 
831 #undef __FUNC__
832 #define __FUNC__ "MatGetDiagonal_SeqDense"
833 int MatGetDiagonal_SeqDense(Mat A,Vec v)
834 {
835   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
836   int          i, n, len;
837   Scalar       *x, zero = 0.0;
838 
839   PetscFunctionBegin;
840   VecSet(&zero,v);
841   VecGetArray(v,&x); VecGetSize(v,&n);
842   len = PetscMin(mat->m,mat->n);
843   if (n != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming mat and vec");
844   for ( i=0; i<len; i++ ) {
845     x[i] = mat->v[i*mat->m + i];
846   }
847   PetscFunctionReturn(0);
848 }
849 
850 #undef __FUNC__
851 #define __FUNC__ "MatDiagonalScale_SeqDense"
852 int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
853 {
854   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
855   Scalar       *l,*r,x,*v;
856   int          i,j,m = mat->m, n = mat->n;
857 
858   PetscFunctionBegin;
859   if (ll) {
860     VecGetArray(ll,&l); VecGetSize(ll,&m);
861     if (m != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vec wrong size");
862     for ( i=0; i<m; i++ ) {
863       x = l[i];
864       v = mat->v + i;
865       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
866     }
867     PLogFlops(n*m);
868   }
869   if (rr) {
870     VecGetArray(rr,&r); VecGetSize(rr,&n);
871     if (n != mat->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vec wrong size");
872     for ( i=0; i<n; i++ ) {
873       x = r[i];
874       v = mat->v + i*m;
875       for ( j=0; j<m; j++ ) { (*v++) *= x;}
876     }
877     PLogFlops(n*m);
878   }
879   PetscFunctionReturn(0);
880 }
881 
882 #undef __FUNC__
883 #define __FUNC__ "MatNorm_SeqDense"
884 int MatNorm_SeqDense(Mat A,NormType type,double *norm)
885 {
886   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
887   Scalar       *v = mat->v;
888   double       sum = 0.0;
889   int          i, j;
890 
891   PetscFunctionBegin;
892   if (type == NORM_FROBENIUS) {
893     for (i=0; i<mat->n*mat->m; i++ ) {
894 #if defined(USE_PETSC_COMPLEX)
895       sum += real(conj(*v)*(*v)); v++;
896 #else
897       sum += (*v)*(*v); v++;
898 #endif
899     }
900     *norm = sqrt(sum);
901     PLogFlops(2*mat->n*mat->m);
902   } else if (type == NORM_1) {
903     *norm = 0.0;
904     for ( j=0; j<mat->n; j++ ) {
905       sum = 0.0;
906       for ( i=0; i<mat->m; i++ ) {
907         sum += PetscAbsScalar(*v);  v++;
908       }
909       if (sum > *norm) *norm = sum;
910     }
911     PLogFlops(mat->n*mat->m);
912   } else if (type == NORM_INFINITY) {
913     *norm = 0.0;
914     for ( j=0; j<mat->m; j++ ) {
915       v = mat->v + j;
916       sum = 0.0;
917       for ( i=0; i<mat->n; i++ ) {
918         sum += PetscAbsScalar(*v); v += mat->m;
919       }
920       if (sum > *norm) *norm = sum;
921     }
922     PLogFlops(mat->n*mat->m);
923   } else {
924     SETERRQ(PETSC_ERR_SUP,0,"No two norm");
925   }
926   PetscFunctionReturn(0);
927 }
928 
929 #undef __FUNC__
930 #define __FUNC__ "MatSetOption_SeqDense"
931 int MatSetOption_SeqDense(Mat A,MatOption op)
932 {
933   Mat_SeqDense *aij = (Mat_SeqDense *) A->data;
934 
935   PetscFunctionBegin;
936   if (op == MAT_ROW_ORIENTED)            aij->roworiented = 1;
937   else if (op == MAT_COLUMN_ORIENTED)    aij->roworiented = 0;
938   else if (op == MAT_ROWS_SORTED ||
939            op == MAT_ROWS_UNSORTED ||
940            op == MAT_COLUMNS_SORTED ||
941            op == MAT_COLUMNS_UNSORTED ||
942            op == MAT_SYMMETRIC ||
943            op == MAT_STRUCTURALLY_SYMMETRIC ||
944            op == MAT_NO_NEW_NONZERO_LOCATIONS ||
945            op == MAT_YES_NEW_NONZERO_LOCATIONS ||
946            op == MAT_NEW_NONZERO_LOCATION_ERROR ||
947            op == MAT_NO_NEW_DIAGONALS ||
948            op == MAT_YES_NEW_DIAGONALS ||
949            op == MAT_IGNORE_OFF_PROC_ENTRIES)
950     PLogInfo(A,"Info:MatSetOption_SeqDense:Option ignored\n");
951   else {
952     SETERRQ(PETSC_ERR_SUP,0,"unknown option");
953   }
954   PetscFunctionReturn(0);
955 }
956 
957 #undef __FUNC__
958 #define __FUNC__ "MatZeroEntries_SeqDense"
959 int MatZeroEntries_SeqDense(Mat A)
960 {
961   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
962 
963   PetscFunctionBegin;
964   PetscMemzero(l->v,l->m*l->n*sizeof(Scalar));
965   PetscFunctionReturn(0);
966 }
967 
968 #undef __FUNC__
969 #define __FUNC__ "MatGetBlockSize_SeqDense"
970 int MatGetBlockSize_SeqDense(Mat A,int *bs)
971 {
972   PetscFunctionBegin;
973   *bs = 1;
974   PetscFunctionReturn(0);
975 }
976 
977 #undef __FUNC__
978 #define __FUNC__ "MatZeroRows_SeqDense"
979 int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
980 {
981   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
982   int          n = l->n, i, j,ierr,N, *rows;
983   Scalar       *slot;
984 
985   PetscFunctionBegin;
986   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
987   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
988   for ( i=0; i<N; i++ ) {
989     slot = l->v + rows[i];
990     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
991   }
992   if (diag) {
993     for ( i=0; i<N; i++ ) {
994       slot = l->v + (n+1)*rows[i];
995       *slot = *diag;
996     }
997   }
998   ISRestoreIndices(is,&rows);
999   PetscFunctionReturn(0);
1000 }
1001 
1002 #undef __FUNC__
1003 #define __FUNC__ "MatGetSize_SeqDense"
1004 int MatGetSize_SeqDense(Mat A,int *m,int *n)
1005 {
1006   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1007 
1008   PetscFunctionBegin;
1009   *m = mat->m; *n = mat->n;
1010   PetscFunctionReturn(0);
1011 }
1012 
1013 #undef __FUNC__
1014 #define __FUNC__ "MatGetOwnershipRange_SeqDense"
1015 int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
1016 {
1017   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1018 
1019   PetscFunctionBegin;
1020   *m = 0; *n = mat->m;
1021   PetscFunctionReturn(0);
1022 }
1023 
1024 #undef __FUNC__
1025 #define __FUNC__ "MatGetArray_SeqDense"
1026 int MatGetArray_SeqDense(Mat A,Scalar **array)
1027 {
1028   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1029 
1030   PetscFunctionBegin;
1031   *array = mat->v;
1032   PetscFunctionReturn(0);
1033 }
1034 
1035 #undef __FUNC__
1036 #define __FUNC__ "MatRestoreArray_SeqDense"
1037 int MatRestoreArray_SeqDense(Mat A,Scalar **array)
1038 {
1039   PetscFunctionBegin;
1040   PetscFunctionReturn(0);
1041 }
1042 
1043 #undef __FUNC__
1044 #define __FUNC__ "MatGetSubMatrix_SeqDense"
1045 static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,
1046                                     Mat *submat)
1047 {
1048   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1049   int          nznew, *smap, i, j, ierr, oldcols = mat->n;
1050   int          *irow, *icol, nrows, ncols, *cwork;
1051   Scalar       *vwork, *val;
1052   Mat          newmat;
1053 
1054   PetscFunctionBegin;
1055   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
1056   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
1057   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
1058   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
1059 
1060   smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap);
1061   cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork);
1062   vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
1063   PetscMemzero((char*)smap,oldcols*sizeof(int));
1064   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
1065 
1066   /* Create and fill new matrix */
1067   ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr);
1068   for (i=0; i<nrows; i++) {
1069     nznew = 0;
1070     val   = mat->v + irow[i];
1071     for (j=0; j<oldcols; j++) {
1072       if (smap[j]) {
1073         cwork[nznew]   = smap[j] - 1;
1074         vwork[nznew++] = val[j * mat->m];
1075       }
1076     }
1077     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES);CHKERRQ(ierr);
1078   }
1079   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1080   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1081 
1082   /* Free work space */
1083   PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
1084   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1085   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
1086   *submat = newmat;
1087   PetscFunctionReturn(0);
1088 }
1089 
1090 #undef __FUNC__
1091 #define __FUNC__ "MatGetSubMatrices_SeqDense"
1092 int MatGetSubMatrices_SeqDense(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
1093                                     Mat **B)
1094 {
1095   int ierr,i;
1096 
1097   PetscFunctionBegin;
1098   if (scall == MAT_INITIAL_MATRIX) {
1099     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1100   }
1101 
1102   for ( i=0; i<n; i++ ) {
1103     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr);
1104   }
1105   PetscFunctionReturn(0);
1106 }
1107 
1108 #undef __FUNC__
1109 #define __FUNC__ "MatCopy_SeqDense"
1110 int MatCopy_SeqDense(Mat A, Mat B)
1111 {
1112   Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data;
1113   int          ierr;
1114 
1115   PetscFunctionBegin;
1116   if (B->type != MATSEQDENSE) {
1117     ierr = MatCopy_Basic(A,B);CHKERRQ(ierr);
1118     PetscFunctionReturn(0);
1119   }
1120   if (a->m != b->m || a->n != b->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"size(B) != size(A)");
1121   PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar));
1122   PetscFunctionReturn(0);
1123 }
1124 
1125 /* -------------------------------------------------------------------*/
1126 static struct _MatOps MatOps = {MatSetValues_SeqDense,
1127        MatGetRow_SeqDense,
1128        MatRestoreRow_SeqDense,
1129        MatMult_SeqDense,
1130        MatMultAdd_SeqDense,
1131        MatMultTrans_SeqDense,
1132        MatMultTransAdd_SeqDense,
1133        MatSolve_SeqDense,
1134        MatSolveAdd_SeqDense,
1135        MatSolveTrans_SeqDense,
1136        MatSolveTransAdd_SeqDense,
1137        MatLUFactor_SeqDense,
1138        MatCholeskyFactor_SeqDense,
1139        MatRelax_SeqDense,
1140        MatTranspose_SeqDense,
1141        MatGetInfo_SeqDense,
1142        MatEqual_SeqDense,
1143        MatGetDiagonal_SeqDense,
1144        MatDiagonalScale_SeqDense,
1145        MatNorm_SeqDense,
1146        0,
1147        0,
1148        0,
1149        MatSetOption_SeqDense,
1150        MatZeroEntries_SeqDense,
1151        MatZeroRows_SeqDense,
1152        MatLUFactorSymbolic_SeqDense,
1153        MatLUFactorNumeric_SeqDense,
1154        MatCholeskyFactorSymbolic_SeqDense,
1155        MatCholeskyFactorNumeric_SeqDense,
1156        MatGetSize_SeqDense,
1157        MatGetSize_SeqDense,
1158        MatGetOwnershipRange_SeqDense,
1159        0,
1160        0,
1161        MatGetArray_SeqDense,
1162        MatRestoreArray_SeqDense,
1163        MatConvertSameType_SeqDense,0,0,0,0,
1164        MatAXPY_SeqDense,MatGetSubMatrices_SeqDense,0,
1165        MatGetValues_SeqDense,
1166        MatCopy_SeqDense,0,MatScale_SeqDense,
1167        0,0,0,MatGetBlockSize_SeqDense};
1168 
1169 #undef __FUNC__
1170 #define __FUNC__ "MatCreateSeqDense"
1171 /*@C
1172    MatCreateSeqDense - Creates a sequential dense matrix that
1173    is stored in column major order (the usual Fortran 77 manner). Many
1174    of the matrix operations use the BLAS and LAPACK routines.
1175 
1176    Input Parameters:
1177 .  comm - MPI communicator, set to PETSC_COMM_SELF
1178 .  m - number of rows
1179 .  n - number of columns
1180 .  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1181    to control all matrix memory allocation.
1182 
1183    Output Parameter:
1184 .  A - the matrix
1185 
1186   Notes:
1187   The data input variable is intended primarily for Fortran programmers
1188   who wish to allocate their own matrix memory space.  Most users should
1189   set data=PETSC_NULL.
1190 
1191 .keywords: dense, matrix, LAPACK, BLAS
1192 
1193 .seealso: MatCreate(), MatSetValues()
1194 @*/
1195 int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A)
1196 {
1197   Mat          B;
1198   Mat_SeqDense *b;
1199   int          ierr,flg,size;
1200 
1201   PetscFunctionBegin;
1202   MPI_Comm_size(comm,&size);
1203   if (size > 1) SETERRQ(1,0,"Comm must be of size 1");
1204 
1205   *A            = 0;
1206   PetscHeaderCreate(B,_p_Mat,MAT_COOKIE,MATSEQDENSE,comm,MatDestroy,MatView);
1207   PLogObjectCreate(B);
1208   b             = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(b);
1209   PetscMemzero(b,sizeof(Mat_SeqDense));
1210   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1211   B->destroy    = MatDestroy_SeqDense;
1212   B->view       = MatView_SeqDense;
1213   B->factor     = 0;
1214   B->mapping    = 0;
1215   PLogObjectMemory(B,sizeof(struct _p_Mat));
1216   B->data       = (void *) b;
1217 
1218   b->m = m;  B->m = m; B->M = m;
1219   b->n = n;  B->n = n; B->N = n;
1220   b->pivots       = 0;
1221   b->roworiented  = 1;
1222   if (data == PETSC_NULL) {
1223     b->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(b->v);
1224     PetscMemzero(b->v,m*n*sizeof(Scalar));
1225     b->user_alloc = 0;
1226     PLogObjectMemory(B,n*m*sizeof(Scalar));
1227   }
1228   else { /* user-allocated storage */
1229     b->v = data;
1230     b->user_alloc = 1;
1231   }
1232   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
1233   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr);}
1234 
1235   *A = B;
1236   PetscFunctionReturn(0);
1237 }
1238 
1239 
1240 
1241 
1242 
1243 
1244