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