xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 899cda47322a0d0eb8e2428039961ef470104e3e)
1 #define PETSCMAT_DLL
2 
3 /*
4      Defines the basic matrix operations for sequential dense.
5 */
6 
7 #include "src/mat/impls/dense/seq/dense.h"
8 #include "petscblaslapack.h"
9 
10 #undef __FUNCT__
11 #define __FUNCT__ "MatAXPY_SeqDense"
12 PetscErrorCode MatAXPY_SeqDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str)
13 {
14   Mat_SeqDense   *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data;
15   PetscScalar    oalpha = alpha;
16   PetscInt       j;
17   PetscBLASInt   N = (PetscBLASInt)X->rmap.n*X->cmap.n,m=(PetscBLASInt)X->rmap.n,ldax = x->lda,lday=y->lda,one = 1;
18   PetscErrorCode ierr;
19 
20   PetscFunctionBegin;
21   if (ldax>m || lday>m) {
22     for (j=0; j<X->cmap.n; j++) {
23       BLASaxpy_(&m,&oalpha,x->v+j*ldax,&one,y->v+j*lday,&one);
24     }
25   } else {
26     BLASaxpy_(&N,&oalpha,x->v,&one,y->v,&one);
27   }
28   ierr = PetscLogFlops(2*N-1);CHKERRQ(ierr);
29   PetscFunctionReturn(0);
30 }
31 
32 #undef __FUNCT__
33 #define __FUNCT__ "MatGetInfo_SeqDense"
34 PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
35 {
36   PetscInt     N = A->rmap.n*A->cmap.n;
37 
38   PetscFunctionBegin;
39   info->rows_global       = (double)A->rmap.n;
40   info->columns_global    = (double)A->cmap.n;
41   info->rows_local        = (double)A->rmap.n;
42   info->columns_local     = (double)A->cmap.n;
43   info->block_size        = 1.0;
44   info->nz_allocated      = (double)N;
45   info->nz_used           = (double)N;
46   info->nz_unneeded       = (double)0;
47   info->assemblies        = (double)A->num_ass;
48   info->mallocs           = 0;
49   info->memory            = A->mem;
50   info->fill_ratio_given  = 0;
51   info->fill_ratio_needed = 0;
52   info->factor_mallocs    = 0;
53   PetscFunctionReturn(0);
54 }
55 
56 #undef __FUNCT__
57 #define __FUNCT__ "MatScale_SeqDense"
58 PetscErrorCode MatScale_SeqDense(Mat A,PetscScalar alpha)
59 {
60   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
61   PetscScalar    oalpha = alpha;
62   PetscBLASInt   one = 1,lda = a->lda,j,nz;
63   PetscErrorCode ierr;
64 
65   PetscFunctionBegin;
66   if (lda>A->rmap.n) {
67     nz = (PetscBLASInt)A->rmap.n;
68     for (j=0; j<A->cmap.n; j++) {
69       BLASscal_(&nz,&oalpha,a->v+j*lda,&one);
70     }
71   } else {
72     nz = (PetscBLASInt)A->rmap.n*A->cmap.n;
73     BLASscal_(&nz,&oalpha,a->v,&one);
74   }
75   ierr = PetscLogFlops(nz);CHKERRQ(ierr);
76   PetscFunctionReturn(0);
77 }
78 
79 /* ---------------------------------------------------------------*/
80 /* COMMENT: I have chosen to hide row permutation in the pivots,
81    rather than put it in the Mat->row slot.*/
82 #undef __FUNCT__
83 #define __FUNCT__ "MatLUFactor_SeqDense"
84 PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,MatFactorInfo *minfo)
85 {
86 #if defined(PETSC_MISSING_LAPACK_GETRF)
87   PetscFunctionBegin;
88   SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable.");
89 #else
90   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
91   PetscErrorCode ierr;
92   PetscBLASInt   n = (PetscBLASInt)A->cmap.n,m = (PetscBLASInt)A->rmap.n,info;
93 
94   PetscFunctionBegin;
95   if (!mat->pivots) {
96     ierr = PetscMalloc((A->rmap.n+1)*sizeof(PetscBLASInt),&mat->pivots);CHKERRQ(ierr);
97     ierr = PetscLogObjectMemory(A,A->rmap.n*sizeof(PetscBLASInt));CHKERRQ(ierr);
98   }
99   A->factor = FACTOR_LU;
100   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
101   LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info);
102   if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization");
103   if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
104   ierr = PetscLogFlops((2*A->cmap.n*A->cmap.n*A->cmap.n)/3);CHKERRQ(ierr);
105 #endif
106   PetscFunctionReturn(0);
107 }
108 
109 #undef __FUNCT__
110 #define __FUNCT__ "MatDuplicate_SeqDense"
111 PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
112 {
113   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data,*l;
114   PetscErrorCode ierr;
115   PetscInt       lda = (PetscInt)mat->lda,j,m;
116   Mat            newi;
117 
118   PetscFunctionBegin;
119   ierr = MatCreate(A->comm,&newi);CHKERRQ(ierr);
120   ierr = MatSetSizes(newi,A->rmap.n,A->cmap.n,A->rmap.n,A->cmap.n);CHKERRQ(ierr);
121   ierr = MatSetType(newi,A->type_name);CHKERRQ(ierr);
122   ierr = MatSeqDenseSetPreallocation(newi,PETSC_NULL);CHKERRQ(ierr);
123   if (cpvalues == MAT_COPY_VALUES) {
124     l = (Mat_SeqDense*)newi->data;
125     if (lda>A->rmap.n) {
126       m = A->rmap.n;
127       for (j=0; j<A->cmap.n; j++) {
128 	ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
129       }
130     } else {
131       ierr = PetscMemcpy(l->v,mat->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
132     }
133   }
134   newi->assembled = PETSC_TRUE;
135   *newmat = newi;
136   PetscFunctionReturn(0);
137 }
138 
139 #undef __FUNCT__
140 #define __FUNCT__ "MatLUFactorSymbolic_SeqDense"
141 PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,MatFactorInfo *info,Mat *fact)
142 {
143   PetscErrorCode ierr;
144 
145   PetscFunctionBegin;
146   ierr = MatDuplicate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
147   PetscFunctionReturn(0);
148 }
149 
150 #undef __FUNCT__
151 #define __FUNCT__ "MatLUFactorNumeric_SeqDense"
152 PetscErrorCode MatLUFactorNumeric_SeqDense(Mat A,MatFactorInfo *info_dummy,Mat *fact)
153 {
154   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data,*l = (Mat_SeqDense*)(*fact)->data;
155   PetscErrorCode ierr;
156   PetscInt       lda1=mat->lda,lda2=l->lda, m=A->rmap.n,n=A->cmap.n, j;
157   MatFactorInfo  info;
158 
159   PetscFunctionBegin;
160   /* copy the numerical values */
161   if (lda1>m || lda2>m ) {
162     for (j=0; j<n; j++) {
163       ierr = PetscMemcpy(l->v+j*lda2,mat->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
164     }
165   } else {
166     ierr = PetscMemcpy(l->v,mat->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
167   }
168   (*fact)->factor = 0;
169   ierr = MatLUFactor(*fact,0,0,&info);CHKERRQ(ierr);
170   PetscFunctionReturn(0);
171 }
172 
173 #undef __FUNCT__
174 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense"
175 PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,MatFactorInfo *info,Mat *fact)
176 {
177   PetscErrorCode ierr;
178 
179   PetscFunctionBegin;
180   ierr = MatConvert(A,MATSAME,MAT_INITIAL_MATRIX,fact);CHKERRQ(ierr);
181   PetscFunctionReturn(0);
182 }
183 
184 #undef __FUNCT__
185 #define __FUNCT__ "MatCholeskyFactor_SeqDense"
186 PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,MatFactorInfo *factinfo)
187 {
188 #if defined(PETSC_MISSING_LAPACK_POTRF)
189   PetscFunctionBegin;
190   SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable.");
191 #else
192   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
193   PetscErrorCode ierr;
194   PetscBLASInt   n = (PetscBLASInt)A->cmap.n,info;
195 
196   PetscFunctionBegin;
197   if (mat->pivots) {
198     ierr = PetscFree(mat->pivots);CHKERRQ(ierr);
199     ierr = PetscLogObjectMemory(A,-A->rmap.n*sizeof(PetscInt));CHKERRQ(ierr);
200     mat->pivots = 0;
201   }
202   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
203   LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info);
204   if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1);
205   A->factor = FACTOR_CHOLESKY;
206   ierr = PetscLogFlops((A->cmap.n*A->cmap.n*A->cmap.n)/3);CHKERRQ(ierr);
207 #endif
208   PetscFunctionReturn(0);
209 }
210 
211 #undef __FUNCT__
212 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense"
213 PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat A,MatFactorInfo *info_dummy,Mat *fact)
214 {
215   PetscErrorCode ierr;
216   MatFactorInfo  info;
217 
218   PetscFunctionBegin;
219   info.fill = 1.0;
220   ierr = MatCholeskyFactor_SeqDense(*fact,0,&info);CHKERRQ(ierr);
221   PetscFunctionReturn(0);
222 }
223 
224 #undef __FUNCT__
225 #define __FUNCT__ "MatSolve_SeqDense"
226 PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
227 {
228   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
229   PetscErrorCode ierr;
230   PetscBLASInt   m = (PetscBLASInt)A->rmap.n, one = 1,info;
231   PetscScalar    *x,*y;
232 
233   PetscFunctionBegin;
234   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
235   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
236   ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
237   if (A->factor == FACTOR_LU) {
238 #if defined(PETSC_MISSING_LAPACK_GETRS)
239     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
240 #else
241     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
242     if (info) SETERRQ(PETSC_ERR_LIB,"GETRS - Bad solve");
243 #endif
244   } else if (A->factor == FACTOR_CHOLESKY){
245 #if defined(PETSC_MISSING_LAPACK_POTRS)
246     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
247 #else
248     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
249     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS Bad solve");
250 #endif
251   }
252   else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
253   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
254   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
255   ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n - A->cmap.n);CHKERRQ(ierr);
256   PetscFunctionReturn(0);
257 }
258 
259 #undef __FUNCT__
260 #define __FUNCT__ "MatSolveTranspose_SeqDense"
261 PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
262 {
263   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
264   PetscErrorCode ierr;
265   PetscBLASInt   m = (PetscBLASInt) A->rmap.n,one = 1,info;
266   PetscScalar    *x,*y;
267 
268   PetscFunctionBegin;
269   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
270   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
271   ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
272   /* assume if pivots exist then use LU; else Cholesky */
273   if (mat->pivots) {
274 #if defined(PETSC_MISSING_LAPACK_GETRS)
275     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
276 #else
277     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
278     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve");
279 #endif
280   } else {
281 #if defined(PETSC_MISSING_LAPACK_POTRS)
282     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
283 #else
284     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
285     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve");
286 #endif
287   }
288   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
289   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
290   ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n - A->cmap.n);CHKERRQ(ierr);
291   PetscFunctionReturn(0);
292 }
293 
294 #undef __FUNCT__
295 #define __FUNCT__ "MatSolveAdd_SeqDense"
296 PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
297 {
298   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
299   PetscErrorCode ierr;
300   PetscBLASInt   m = (PetscBLASInt)A->rmap.n,one = 1,info;
301   PetscScalar    *x,*y,sone = 1.0;
302   Vec            tmp = 0;
303 
304   PetscFunctionBegin;
305   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
306   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
307   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
308   if (yy == zz) {
309     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
310     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
311     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
312   }
313   ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
314   /* assume if pivots exist then use LU; else Cholesky */
315   if (mat->pivots) {
316 #if defined(PETSC_MISSING_LAPACK_GETRS)
317     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
318 #else
319     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
320     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
321 #endif
322   } else {
323 #if defined(PETSC_MISSING_LAPACK_POTRS)
324     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
325 #else
326     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
327     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
328 #endif
329   }
330   if (tmp) {ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);}
331   else     {ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);}
332   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
333   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
334   ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n);CHKERRQ(ierr);
335   PetscFunctionReturn(0);
336 }
337 
338 #undef __FUNCT__
339 #define __FUNCT__ "MatSolveTransposeAdd_SeqDense"
340 PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
341 {
342   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
343   PetscErrorCode ierr;
344   PetscBLASInt   m = (PetscBLASInt)A->rmap.n,one = 1,info;
345   PetscScalar    *x,*y,sone = 1.0;
346   Vec            tmp;
347 
348   PetscFunctionBegin;
349   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
350   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
351   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
352   if (yy == zz) {
353     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
354     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
355     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
356   }
357   ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
358   /* assume if pivots exist then use LU; else Cholesky */
359   if (mat->pivots) {
360 #if defined(PETSC_MISSING_LAPACK_GETRS)
361     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
362 #else
363     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
364     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
365 #endif
366   } else {
367 #if defined(PETSC_MISSING_LAPACK_POTRS)
368     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
369 #else
370     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
371     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
372 #endif
373   }
374   if (tmp) {
375     ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr);
376     ierr = VecDestroy(tmp);CHKERRQ(ierr);
377   } else {
378     ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);
379   }
380   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
381   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
382   ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n);CHKERRQ(ierr);
383   PetscFunctionReturn(0);
384 }
385 /* ------------------------------------------------------------------*/
386 #undef __FUNCT__
387 #define __FUNCT__ "MatRelax_SeqDense"
388 PetscErrorCode MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx)
389 {
390   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
391   PetscScalar    *x,*b,*v = mat->v,zero = 0.0,xt;
392   PetscErrorCode ierr;
393   PetscInt       m = A->rmap.n,i;
394 #if !defined(PETSC_USE_COMPLEX)
395   PetscBLASInt   bm = (PetscBLASInt)m, o = 1;
396 #endif
397 
398   PetscFunctionBegin;
399   if (flag & SOR_ZERO_INITIAL_GUESS) {
400     /* this is a hack fix, should have another version without the second BLASdot */
401     ierr = VecSet(xx,zero);CHKERRQ(ierr);
402   }
403   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
404   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
405   its  = its*lits;
406   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
407   while (its--) {
408     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
409       for (i=0; i<m; i++) {
410 #if defined(PETSC_USE_COMPLEX)
411         /* cannot use BLAS dot for complex because compiler/linker is
412            not happy about returning a double complex */
413         PetscInt         _i;
414         PetscScalar sum = b[i];
415         for (_i=0; _i<m; _i++) {
416           sum -= PetscConj(v[i+_i*m])*x[_i];
417         }
418         xt = sum;
419 #else
420         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
421 #endif
422         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
423       }
424     }
425     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
426       for (i=m-1; i>=0; i--) {
427 #if defined(PETSC_USE_COMPLEX)
428         /* cannot use BLAS dot for complex because compiler/linker is
429            not happy about returning a double complex */
430         PetscInt         _i;
431         PetscScalar sum = b[i];
432         for (_i=0; _i<m; _i++) {
433           sum -= PetscConj(v[i+_i*m])*x[_i];
434         }
435         xt = sum;
436 #else
437         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
438 #endif
439         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
440       }
441     }
442   }
443   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
444   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
445   PetscFunctionReturn(0);
446 }
447 
448 /* -----------------------------------------------------------------*/
449 #undef __FUNCT__
450 #define __FUNCT__ "MatMultTranspose_SeqDense"
451 PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
452 {
453   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
454   PetscScalar    *v = mat->v,*x,*y;
455   PetscErrorCode ierr;
456   PetscBLASInt   m = (PetscBLASInt)A->rmap.n, n = (PetscBLASInt)A->cmap.n,_One=1;
457   PetscScalar    _DOne=1.0,_DZero=0.0;
458 
459   PetscFunctionBegin;
460   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
461   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
462   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
463   BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One);
464   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
465   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
466   ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n - A->cmap.n);CHKERRQ(ierr);
467   PetscFunctionReturn(0);
468 }
469 
470 #undef __FUNCT__
471 #define __FUNCT__ "MatMult_SeqDense"
472 PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy)
473 {
474   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
475   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0;
476   PetscErrorCode ierr;
477   PetscBLASInt   m = (PetscBLASInt)A->rmap.n, n = (PetscBLASInt)A->cmap.n, _One=1;
478 
479   PetscFunctionBegin;
480   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
481   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
482   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
483   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One);
484   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
485   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
486   ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n - A->rmap.n);CHKERRQ(ierr);
487   PetscFunctionReturn(0);
488 }
489 
490 #undef __FUNCT__
491 #define __FUNCT__ "MatMultAdd_SeqDense"
492 PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
493 {
494   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
495   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0;
496   PetscErrorCode ierr;
497   PetscBLASInt   m = (PetscBLASInt)A->rmap.n, n = (PetscBLASInt)A->cmap.n, _One=1;
498 
499   PetscFunctionBegin;
500   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
501   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
502   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
503   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
504   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
505   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
506   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
507   ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n);CHKERRQ(ierr);
508   PetscFunctionReturn(0);
509 }
510 
511 #undef __FUNCT__
512 #define __FUNCT__ "MatMultTransposeAdd_SeqDense"
513 PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
514 {
515   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
516   PetscScalar    *v = mat->v,*x,*y;
517   PetscErrorCode ierr;
518   PetscBLASInt   m = (PetscBLASInt)A->rmap.n, n = (PetscBLASInt)A->cmap.n, _One=1;
519   PetscScalar    _DOne=1.0;
520 
521   PetscFunctionBegin;
522   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
523   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
524   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
525   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
526   BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
527   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
528   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
529   ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n);CHKERRQ(ierr);
530   PetscFunctionReturn(0);
531 }
532 
533 /* -----------------------------------------------------------------*/
534 #undef __FUNCT__
535 #define __FUNCT__ "MatGetRow_SeqDense"
536 PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
537 {
538   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
539   PetscScalar    *v;
540   PetscErrorCode ierr;
541   PetscInt       i;
542 
543   PetscFunctionBegin;
544   *ncols = A->cmap.n;
545   if (cols) {
546     ierr = PetscMalloc((A->cmap.n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr);
547     for (i=0; i<A->cmap.n; i++) (*cols)[i] = i;
548   }
549   if (vals) {
550     ierr = PetscMalloc((A->cmap.n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr);
551     v    = mat->v + row;
552     for (i=0; i<A->cmap.n; i++) {(*vals)[i] = *v; v += mat->lda;}
553   }
554   PetscFunctionReturn(0);
555 }
556 
557 #undef __FUNCT__
558 #define __FUNCT__ "MatRestoreRow_SeqDense"
559 PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
560 {
561   PetscErrorCode ierr;
562   PetscFunctionBegin;
563   if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);}
564   if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); }
565   PetscFunctionReturn(0);
566 }
567 /* ----------------------------------------------------------------*/
568 #undef __FUNCT__
569 #define __FUNCT__ "MatSetValues_SeqDense"
570 PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv)
571 {
572   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
573   PetscInt     i,j,idx=0;
574 
575   PetscFunctionBegin;
576   if (!mat->roworiented) {
577     if (addv == INSERT_VALUES) {
578       for (j=0; j<n; j++) {
579         if (indexn[j] < 0) {idx += m; continue;}
580 #if defined(PETSC_USE_DEBUG)
581         if (indexn[j] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap.n-1);
582 #endif
583         for (i=0; i<m; i++) {
584           if (indexm[i] < 0) {idx++; continue;}
585 #if defined(PETSC_USE_DEBUG)
586           if (indexm[i] >= A->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap.n-1);
587 #endif
588           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
589         }
590       }
591     } else {
592       for (j=0; j<n; j++) {
593         if (indexn[j] < 0) {idx += m; continue;}
594 #if defined(PETSC_USE_DEBUG)
595         if (indexn[j] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap.n-1);
596 #endif
597         for (i=0; i<m; i++) {
598           if (indexm[i] < 0) {idx++; continue;}
599 #if defined(PETSC_USE_DEBUG)
600           if (indexm[i] >= A->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap.n-1);
601 #endif
602           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
603         }
604       }
605     }
606   } else {
607     if (addv == INSERT_VALUES) {
608       for (i=0; i<m; i++) {
609         if (indexm[i] < 0) { idx += n; continue;}
610 #if defined(PETSC_USE_DEBUG)
611         if (indexm[i] >= A->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap.n-1);
612 #endif
613         for (j=0; j<n; j++) {
614           if (indexn[j] < 0) { idx++; continue;}
615 #if defined(PETSC_USE_DEBUG)
616           if (indexn[j] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap.n-1);
617 #endif
618           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
619         }
620       }
621     } else {
622       for (i=0; i<m; i++) {
623         if (indexm[i] < 0) { idx += n; continue;}
624 #if defined(PETSC_USE_DEBUG)
625         if (indexm[i] >= A->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap.n-1);
626 #endif
627         for (j=0; j<n; j++) {
628           if (indexn[j] < 0) { idx++; continue;}
629 #if defined(PETSC_USE_DEBUG)
630           if (indexn[j] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap.n-1);
631 #endif
632           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
633         }
634       }
635     }
636   }
637   PetscFunctionReturn(0);
638 }
639 
640 #undef __FUNCT__
641 #define __FUNCT__ "MatGetValues_SeqDense"
642 PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[])
643 {
644   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
645   PetscInt     i,j;
646   PetscScalar  *vpt = v;
647 
648   PetscFunctionBegin;
649   /* row-oriented output */
650   for (i=0; i<m; i++) {
651     for (j=0; j<n; j++) {
652       *vpt++ = mat->v[indexn[j]*mat->lda + indexm[i]];
653     }
654   }
655   PetscFunctionReturn(0);
656 }
657 
658 /* -----------------------------------------------------------------*/
659 
660 #include "petscsys.h"
661 
662 #undef __FUNCT__
663 #define __FUNCT__ "MatLoad_SeqDense"
664 PetscErrorCode MatLoad_SeqDense(PetscViewer viewer, MatType type,Mat *A)
665 {
666   Mat_SeqDense   *a;
667   Mat            B;
668   PetscErrorCode ierr;
669   PetscInt       *scols,i,j,nz,header[4];
670   int            fd;
671   PetscMPIInt    size;
672   PetscInt       *rowlengths = 0,M,N,*cols;
673   PetscScalar    *vals,*svals,*v,*w;
674   MPI_Comm       comm = ((PetscObject)viewer)->comm;
675 
676   PetscFunctionBegin;
677   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
678   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
679   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
680   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
681   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
682   M = header[1]; N = header[2]; nz = header[3];
683 
684   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
685     ierr = MatCreate(comm,A);CHKERRQ(ierr);
686     ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr);
687     ierr = MatSetType(*A,type);CHKERRQ(ierr);
688     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
689     B    = *A;
690     a    = (Mat_SeqDense*)B->data;
691     v    = a->v;
692     /* Allocate some temp space to read in the values and then flip them
693        from row major to column major */
694     ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr);
695     /* read in nonzero values */
696     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr);
697     /* now flip the values and store them in the matrix*/
698     for (j=0; j<N; j++) {
699       for (i=0; i<M; i++) {
700         *v++ =w[i*N+j];
701       }
702     }
703     ierr = PetscFree(w);CHKERRQ(ierr);
704     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
705     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
706   } else {
707     /* read row lengths */
708     ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
709     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
710 
711     /* create our matrix */
712     ierr = MatCreate(comm,A);CHKERRQ(ierr);
713     ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr);
714     ierr = MatSetType(*A,type);CHKERRQ(ierr);
715     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
716     B = *A;
717     a = (Mat_SeqDense*)B->data;
718     v = a->v;
719 
720     /* read column indices and nonzeros */
721     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr);
722     cols = scols;
723     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
724     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr);
725     vals = svals;
726     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
727 
728     /* insert into matrix */
729     for (i=0; i<M; i++) {
730       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
731       svals += rowlengths[i]; scols += rowlengths[i];
732     }
733     ierr = PetscFree(vals);CHKERRQ(ierr);
734     ierr = PetscFree(cols);CHKERRQ(ierr);
735     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
736 
737     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
738     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
739   }
740   PetscFunctionReturn(0);
741 }
742 
743 #include "petscsys.h"
744 
745 #undef __FUNCT__
746 #define __FUNCT__ "MatView_SeqDense_ASCII"
747 static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
748 {
749   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
750   PetscErrorCode    ierr;
751   PetscInt          i,j;
752   const char        *name;
753   PetscScalar       *v;
754   PetscViewerFormat format;
755 
756   PetscFunctionBegin;
757   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
758   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
759   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
760     PetscFunctionReturn(0);  /* do nothing for now */
761   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
762     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
763     for (i=0; i<A->rmap.n; i++) {
764       v = a->v + i;
765       ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr);
766       for (j=0; j<A->cmap.n; j++) {
767 #if defined(PETSC_USE_COMPLEX)
768         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
769           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
770         } else if (PetscRealPart(*v)) {
771           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,PetscRealPart(*v));CHKERRQ(ierr);
772         }
773 #else
774         if (*v) {
775           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,*v);CHKERRQ(ierr);
776         }
777 #endif
778         v += a->lda;
779       }
780       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
781     }
782     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
783   } else {
784     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
785 #if defined(PETSC_USE_COMPLEX)
786     PetscTruth allreal = PETSC_TRUE;
787     /* determine if matrix has all real values */
788     v = a->v;
789     for (i=0; i<A->rmap.n*A->cmap.n; i++) {
790 	if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;}
791     }
792 #endif
793     if (format == PETSC_VIEWER_ASCII_MATLAB) {
794       ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
795       ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap.n,A->cmap.n);CHKERRQ(ierr);
796       ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap.n,A->cmap.n);CHKERRQ(ierr);
797       ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr);
798     }
799 
800     for (i=0; i<A->rmap.n; i++) {
801       v = a->v + i;
802       for (j=0; j<A->cmap.n; j++) {
803 #if defined(PETSC_USE_COMPLEX)
804         if (allreal) {
805           ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",PetscRealPart(*v));CHKERRQ(ierr);
806         } else {
807           ierr = PetscViewerASCIIPrintf(viewer,"%6.4e + %6.4e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
808         }
809 #else
810         ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",*v);CHKERRQ(ierr);
811 #endif
812 	v += a->lda;
813       }
814       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
815     }
816     if (format == PETSC_VIEWER_ASCII_MATLAB) {
817       ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr);
818     }
819     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
820   }
821   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
822   PetscFunctionReturn(0);
823 }
824 
825 #undef __FUNCT__
826 #define __FUNCT__ "MatView_SeqDense_Binary"
827 static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
828 {
829   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
830   PetscErrorCode    ierr;
831   int               fd;
832   PetscInt          ict,j,n = A->cmap.n,m = A->rmap.n,i,*col_lens,nz = m*n;
833   PetscScalar       *v,*anonz,*vals;
834   PetscViewerFormat format;
835 
836   PetscFunctionBegin;
837   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
838 
839   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
840   if (format == PETSC_VIEWER_BINARY_NATIVE) {
841     /* store the matrix as a dense matrix */
842     ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
843     col_lens[0] = MAT_FILE_COOKIE;
844     col_lens[1] = m;
845     col_lens[2] = n;
846     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
847     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
848     ierr = PetscFree(col_lens);CHKERRQ(ierr);
849 
850     /* write out matrix, by rows */
851     ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
852     v    = a->v;
853     for (i=0; i<m; i++) {
854       for (j=0; j<n; j++) {
855         vals[i + j*m] = *v++;
856       }
857     }
858     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
859     ierr = PetscFree(vals);CHKERRQ(ierr);
860   } else {
861     ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
862     col_lens[0] = MAT_FILE_COOKIE;
863     col_lens[1] = m;
864     col_lens[2] = n;
865     col_lens[3] = nz;
866 
867     /* store lengths of each row and write (including header) to file */
868     for (i=0; i<m; i++) col_lens[4+i] = n;
869     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
870 
871     /* Possibly should write in smaller increments, not whole matrix at once? */
872     /* store column indices (zero start index) */
873     ict = 0;
874     for (i=0; i<m; i++) {
875       for (j=0; j<n; j++) col_lens[ict++] = j;
876     }
877     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
878     ierr = PetscFree(col_lens);CHKERRQ(ierr);
879 
880     /* store nonzero values */
881     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr);
882     ict  = 0;
883     for (i=0; i<m; i++) {
884       v = a->v + i;
885       for (j=0; j<n; j++) {
886         anonz[ict++] = *v; v += a->lda;
887       }
888     }
889     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
890     ierr = PetscFree(anonz);CHKERRQ(ierr);
891   }
892   PetscFunctionReturn(0);
893 }
894 
895 #undef __FUNCT__
896 #define __FUNCT__ "MatView_SeqDense_Draw_Zoom"
897 PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
898 {
899   Mat               A = (Mat) Aa;
900   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
901   PetscErrorCode    ierr;
902   PetscInt          m = A->rmap.n,n = A->cmap.n,color,i,j;
903   PetscScalar       *v = a->v;
904   PetscViewer       viewer;
905   PetscDraw         popup;
906   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0;
907   PetscViewerFormat format;
908 
909   PetscFunctionBegin;
910 
911   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
912   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
913   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
914 
915   /* Loop over matrix elements drawing boxes */
916   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
917     /* Blue for negative and Red for positive */
918     color = PETSC_DRAW_BLUE;
919     for(j = 0; j < n; j++) {
920       x_l = j;
921       x_r = x_l + 1.0;
922       for(i = 0; i < m; i++) {
923         y_l = m - i - 1.0;
924         y_r = y_l + 1.0;
925 #if defined(PETSC_USE_COMPLEX)
926         if (PetscRealPart(v[j*m+i]) >  0.) {
927           color = PETSC_DRAW_RED;
928         } else if (PetscRealPart(v[j*m+i]) <  0.) {
929           color = PETSC_DRAW_BLUE;
930         } else {
931           continue;
932         }
933 #else
934         if (v[j*m+i] >  0.) {
935           color = PETSC_DRAW_RED;
936         } else if (v[j*m+i] <  0.) {
937           color = PETSC_DRAW_BLUE;
938         } else {
939           continue;
940         }
941 #endif
942         ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
943       }
944     }
945   } else {
946     /* use contour shading to indicate magnitude of values */
947     /* first determine max of all nonzero values */
948     for(i = 0; i < m*n; i++) {
949       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
950     }
951     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
952     ierr  = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
953     if (popup) {ierr  = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);}
954     for(j = 0; j < n; j++) {
955       x_l = j;
956       x_r = x_l + 1.0;
957       for(i = 0; i < m; i++) {
958         y_l   = m - i - 1.0;
959         y_r   = y_l + 1.0;
960         color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i]));
961         ierr  = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
962       }
963     }
964   }
965   PetscFunctionReturn(0);
966 }
967 
968 #undef __FUNCT__
969 #define __FUNCT__ "MatView_SeqDense_Draw"
970 PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
971 {
972   PetscDraw      draw;
973   PetscTruth     isnull;
974   PetscReal      xr,yr,xl,yl,h,w;
975   PetscErrorCode ierr;
976 
977   PetscFunctionBegin;
978   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
979   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
980   if (isnull) PetscFunctionReturn(0);
981 
982   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
983   xr  = A->cmap.n; yr = A->rmap.n; h = yr/10.0; w = xr/10.0;
984   xr += w;    yr += h;  xl = -w;     yl = -h;
985   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
986   ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr);
987   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
988   PetscFunctionReturn(0);
989 }
990 
991 #undef __FUNCT__
992 #define __FUNCT__ "MatView_SeqDense"
993 PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer)
994 {
995   PetscErrorCode ierr;
996   PetscTruth     issocket,iascii,isbinary,isdraw;
997 
998   PetscFunctionBegin;
999   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
1000   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1001   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
1002   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
1003 
1004   if (iascii) {
1005     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
1006 #if defined(PETSC_USE_SOCKET_VIEWER)
1007   } else if (issocket) {
1008     Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1009     if (a->lda>A->rmap.n) SETERRQ(PETSC_ERR_SUP,"Case can not handle LDA");
1010     ierr = PetscViewerSocketPutScalar(viewer,A->rmap.n,A->cmap.n,a->v);CHKERRQ(ierr);
1011 #endif
1012   } else if (isbinary) {
1013     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
1014   } else if (isdraw) {
1015     ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr);
1016   } else {
1017     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name);
1018   }
1019   PetscFunctionReturn(0);
1020 }
1021 
1022 #undef __FUNCT__
1023 #define __FUNCT__ "MatDestroy_SeqDense"
1024 PetscErrorCode MatDestroy_SeqDense(Mat mat)
1025 {
1026   Mat_SeqDense   *l = (Mat_SeqDense*)mat->data;
1027   PetscErrorCode ierr;
1028 
1029   PetscFunctionBegin;
1030 #if defined(PETSC_USE_LOG)
1031   PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap.n,mat->cmap.n);
1032 #endif
1033   if (l->pivots) {ierr = PetscFree(l->pivots);CHKERRQ(ierr);}
1034   if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);}
1035   ierr = PetscFree(l);CHKERRQ(ierr);
1036   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
1037   PetscFunctionReturn(0);
1038 }
1039 
1040 #undef __FUNCT__
1041 #define __FUNCT__ "MatTranspose_SeqDense"
1042 PetscErrorCode MatTranspose_SeqDense(Mat A,Mat *matout)
1043 {
1044   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1045   PetscErrorCode ierr;
1046   PetscInt       k,j,m,n,M;
1047   PetscScalar    *v,tmp;
1048 
1049   PetscFunctionBegin;
1050   v = mat->v; m = A->rmap.n; M = mat->lda; n = A->cmap.n;
1051   if (!matout) { /* in place transpose */
1052     if (m != n) {
1053       SETERRQ(PETSC_ERR_SUP,"Can not transpose non-square matrix in place");
1054     } else {
1055       for (j=0; j<m; j++) {
1056         for (k=0; k<j; k++) {
1057           tmp = v[j + k*M];
1058           v[j + k*M] = v[k + j*M];
1059           v[k + j*M] = tmp;
1060         }
1061       }
1062     }
1063   } else { /* out-of-place transpose */
1064     Mat          tmat;
1065     Mat_SeqDense *tmatd;
1066     PetscScalar  *v2;
1067 
1068     ierr  = MatCreate(A->comm,&tmat);CHKERRQ(ierr);
1069     ierr  = MatSetSizes(tmat,A->cmap.n,A->rmap.n,A->cmap.n,A->rmap.n);CHKERRQ(ierr);
1070     ierr  = MatSetType(tmat,A->type_name);CHKERRQ(ierr);
1071     ierr  = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr);
1072     tmatd = (Mat_SeqDense*)tmat->data;
1073     v = mat->v; v2 = tmatd->v;
1074     for (j=0; j<n; j++) {
1075       for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M];
1076     }
1077     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1078     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1079     *matout = tmat;
1080   }
1081   PetscFunctionReturn(0);
1082 }
1083 
1084 #undef __FUNCT__
1085 #define __FUNCT__ "MatEqual_SeqDense"
1086 PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg)
1087 {
1088   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1089   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
1090   PetscInt     i,j;
1091   PetscScalar  *v1 = mat1->v,*v2 = mat2->v;
1092 
1093   PetscFunctionBegin;
1094   if (A1->rmap.n != A2->rmap.n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1095   if (A1->cmap.n != A2->cmap.n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1096   for (i=0; i<A1->rmap.n; i++) {
1097     v1 = mat1->v+i; v2 = mat2->v+i;
1098     for (j=0; j<A1->cmap.n; j++) {
1099       if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1100       v1 += mat1->lda; v2 += mat2->lda;
1101     }
1102   }
1103   *flg = PETSC_TRUE;
1104   PetscFunctionReturn(0);
1105 }
1106 
1107 #undef __FUNCT__
1108 #define __FUNCT__ "MatGetDiagonal_SeqDense"
1109 PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v)
1110 {
1111   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1112   PetscErrorCode ierr;
1113   PetscInt       i,n,len;
1114   PetscScalar    *x,zero = 0.0;
1115 
1116   PetscFunctionBegin;
1117   ierr = VecSet(v,zero);CHKERRQ(ierr);
1118   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
1119   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1120   len = PetscMin(A->rmap.n,A->cmap.n);
1121   if (n != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
1122   for (i=0; i<len; i++) {
1123     x[i] = mat->v[i*mat->lda + i];
1124   }
1125   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1126   PetscFunctionReturn(0);
1127 }
1128 
1129 #undef __FUNCT__
1130 #define __FUNCT__ "MatDiagonalScale_SeqDense"
1131 PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1132 {
1133   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1134   PetscScalar    *l,*r,x,*v;
1135   PetscErrorCode ierr;
1136   PetscInt       i,j,m = A->rmap.n,n = A->cmap.n;
1137 
1138   PetscFunctionBegin;
1139   if (ll) {
1140     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
1141     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1142     if (m != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1143     for (i=0; i<m; i++) {
1144       x = l[i];
1145       v = mat->v + i;
1146       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1147     }
1148     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1149     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1150   }
1151   if (rr) {
1152     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
1153     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1154     if (n != A->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1155     for (i=0; i<n; i++) {
1156       x = r[i];
1157       v = mat->v + i*m;
1158       for (j=0; j<m; j++) { (*v++) *= x;}
1159     }
1160     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1161     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1162   }
1163   PetscFunctionReturn(0);
1164 }
1165 
1166 #undef __FUNCT__
1167 #define __FUNCT__ "MatNorm_SeqDense"
1168 PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1169 {
1170   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1171   PetscScalar  *v = mat->v;
1172   PetscReal    sum = 0.0;
1173   PetscInt     lda=mat->lda,m=A->rmap.n,i,j;
1174   PetscErrorCode ierr;
1175 
1176   PetscFunctionBegin;
1177   if (type == NORM_FROBENIUS) {
1178     if (lda>m) {
1179       for (j=0; j<A->cmap.n; j++) {
1180 	v = mat->v+j*lda;
1181 	for (i=0; i<m; i++) {
1182 #if defined(PETSC_USE_COMPLEX)
1183 	  sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1184 #else
1185 	  sum += (*v)*(*v); v++;
1186 #endif
1187 	}
1188       }
1189     } else {
1190       for (i=0; i<A->cmap.n*A->rmap.n; i++) {
1191 #if defined(PETSC_USE_COMPLEX)
1192 	sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1193 #else
1194 	sum += (*v)*(*v); v++;
1195 #endif
1196       }
1197     }
1198     *nrm = sqrt(sum);
1199     ierr = PetscLogFlops(2*A->cmap.n*A->rmap.n);CHKERRQ(ierr);
1200   } else if (type == NORM_1) {
1201     *nrm = 0.0;
1202     for (j=0; j<A->cmap.n; j++) {
1203       v = mat->v + j*mat->lda;
1204       sum = 0.0;
1205       for (i=0; i<A->rmap.n; i++) {
1206         sum += PetscAbsScalar(*v);  v++;
1207       }
1208       if (sum > *nrm) *nrm = sum;
1209     }
1210     ierr = PetscLogFlops(A->cmap.n*A->rmap.n);CHKERRQ(ierr);
1211   } else if (type == NORM_INFINITY) {
1212     *nrm = 0.0;
1213     for (j=0; j<A->rmap.n; j++) {
1214       v = mat->v + j;
1215       sum = 0.0;
1216       for (i=0; i<A->cmap.n; i++) {
1217         sum += PetscAbsScalar(*v); v += mat->lda;
1218       }
1219       if (sum > *nrm) *nrm = sum;
1220     }
1221     ierr = PetscLogFlops(A->cmap.n*A->rmap.n);CHKERRQ(ierr);
1222   } else {
1223     SETERRQ(PETSC_ERR_SUP,"No two norm");
1224   }
1225   PetscFunctionReturn(0);
1226 }
1227 
1228 #undef __FUNCT__
1229 #define __FUNCT__ "MatSetOption_SeqDense"
1230 PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op)
1231 {
1232   Mat_SeqDense   *aij = (Mat_SeqDense*)A->data;
1233   PetscErrorCode ierr;
1234 
1235   PetscFunctionBegin;
1236   switch (op) {
1237   case MAT_ROW_ORIENTED:
1238     aij->roworiented = PETSC_TRUE;
1239     break;
1240   case MAT_COLUMN_ORIENTED:
1241     aij->roworiented = PETSC_FALSE;
1242     break;
1243   case MAT_ROWS_SORTED:
1244   case MAT_ROWS_UNSORTED:
1245   case MAT_COLUMNS_SORTED:
1246   case MAT_COLUMNS_UNSORTED:
1247   case MAT_NO_NEW_NONZERO_LOCATIONS:
1248   case MAT_YES_NEW_NONZERO_LOCATIONS:
1249   case MAT_NEW_NONZERO_LOCATION_ERR:
1250   case MAT_NO_NEW_DIAGONALS:
1251   case MAT_YES_NEW_DIAGONALS:
1252   case MAT_IGNORE_OFF_PROC_ENTRIES:
1253   case MAT_USE_HASH_TABLE:
1254     ierr = PetscInfo(A,"Option ignored\n");CHKERRQ(ierr);
1255     break;
1256   case MAT_SYMMETRIC:
1257   case MAT_STRUCTURALLY_SYMMETRIC:
1258   case MAT_NOT_SYMMETRIC:
1259   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
1260   case MAT_HERMITIAN:
1261   case MAT_NOT_HERMITIAN:
1262   case MAT_SYMMETRY_ETERNAL:
1263   case MAT_NOT_SYMMETRY_ETERNAL:
1264     break;
1265   default:
1266     SETERRQ(PETSC_ERR_SUP,"unknown option");
1267   }
1268   PetscFunctionReturn(0);
1269 }
1270 
1271 #undef __FUNCT__
1272 #define __FUNCT__ "MatZeroEntries_SeqDense"
1273 PetscErrorCode MatZeroEntries_SeqDense(Mat A)
1274 {
1275   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
1276   PetscErrorCode ierr;
1277   PetscInt       lda=l->lda,m=A->rmap.n,j;
1278 
1279   PetscFunctionBegin;
1280   if (lda>m) {
1281     for (j=0; j<A->cmap.n; j++) {
1282       ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
1283     }
1284   } else {
1285     ierr = PetscMemzero(l->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
1286   }
1287   PetscFunctionReturn(0);
1288 }
1289 
1290 #undef __FUNCT__
1291 #define __FUNCT__ "MatZeroRows_SeqDense"
1292 PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
1293 {
1294   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
1295   PetscInt       n = A->cmap.n,i,j;
1296   PetscScalar    *slot;
1297 
1298   PetscFunctionBegin;
1299   for (i=0; i<N; i++) {
1300     slot = l->v + rows[i];
1301     for (j=0; j<n; j++) { *slot = 0.0; slot += n;}
1302   }
1303   if (diag != 0.0) {
1304     for (i=0; i<N; i++) {
1305       slot = l->v + (n+1)*rows[i];
1306       *slot = diag;
1307     }
1308   }
1309   PetscFunctionReturn(0);
1310 }
1311 
1312 #undef __FUNCT__
1313 #define __FUNCT__ "MatGetArray_SeqDense"
1314 PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[])
1315 {
1316   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1317 
1318   PetscFunctionBegin;
1319   if (mat->lda != A->rmap.n) SETERRQ(PETSC_ERR_SUP,"Cannot get array for Dense matrices with LDA different from number of rows");
1320   *array = mat->v;
1321   PetscFunctionReturn(0);
1322 }
1323 
1324 #undef __FUNCT__
1325 #define __FUNCT__ "MatRestoreArray_SeqDense"
1326 PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[])
1327 {
1328   PetscFunctionBegin;
1329   *array = 0; /* user cannot accidently use the array later */
1330   PetscFunctionReturn(0);
1331 }
1332 
1333 #undef __FUNCT__
1334 #define __FUNCT__ "MatGetSubMatrix_SeqDense"
1335 static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
1336 {
1337   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1338   PetscErrorCode ierr;
1339   PetscInt       i,j,*irow,*icol,nrows,ncols;
1340   PetscScalar    *av,*bv,*v = mat->v;
1341   Mat            newmat;
1342 
1343   PetscFunctionBegin;
1344   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
1345   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
1346   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
1347   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
1348 
1349   /* Check submatrixcall */
1350   if (scall == MAT_REUSE_MATRIX) {
1351     PetscInt n_cols,n_rows;
1352     ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr);
1353     if (n_rows != nrows || n_cols != ncols) {
1354       /* resize the result result matrix to match number of requested rows/columns */
1355       ierr = MatSetSizes(*B,nrows,nrows,nrows,nrows);CHKERRQ(ierr);
1356     }
1357     newmat = *B;
1358   } else {
1359     /* Create and fill new matrix */
1360     ierr = MatCreate(A->comm,&newmat);CHKERRQ(ierr);
1361     ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
1362     ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr);
1363     ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
1364   }
1365 
1366   /* Now extract the data pointers and do the copy,column at a time */
1367   bv = ((Mat_SeqDense*)newmat->data)->v;
1368 
1369   for (i=0; i<ncols; i++) {
1370     av = v + mat->lda*icol[i];
1371     for (j=0; j<nrows; j++) {
1372       *bv++ = av[irow[j]];
1373     }
1374   }
1375 
1376   /* Assemble the matrices so that the correct flags are set */
1377   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1378   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1379 
1380   /* Free work space */
1381   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
1382   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
1383   *B = newmat;
1384   PetscFunctionReturn(0);
1385 }
1386 
1387 #undef __FUNCT__
1388 #define __FUNCT__ "MatGetSubMatrices_SeqDense"
1389 PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1390 {
1391   PetscErrorCode ierr;
1392   PetscInt       i;
1393 
1394   PetscFunctionBegin;
1395   if (scall == MAT_INITIAL_MATRIX) {
1396     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
1397   }
1398 
1399   for (i=0; i<n; i++) {
1400     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1401   }
1402   PetscFunctionReturn(0);
1403 }
1404 
1405 #undef __FUNCT__
1406 #define __FUNCT__ "MatAssemblyBegin_SeqDense"
1407 PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode)
1408 {
1409   PetscFunctionBegin;
1410   PetscFunctionReturn(0);
1411 }
1412 
1413 #undef __FUNCT__
1414 #define __FUNCT__ "MatAssemblyEnd_SeqDense"
1415 PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode)
1416 {
1417   PetscFunctionBegin;
1418   PetscFunctionReturn(0);
1419 }
1420 
1421 #undef __FUNCT__
1422 #define __FUNCT__ "MatCopy_SeqDense"
1423 PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
1424 {
1425   Mat_SeqDense   *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
1426   PetscErrorCode ierr;
1427   PetscInt       lda1=a->lda,lda2=b->lda, m=A->rmap.n,n=A->cmap.n, j;
1428 
1429   PetscFunctionBegin;
1430   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
1431   if (A->ops->copy != B->ops->copy) {
1432     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
1433     PetscFunctionReturn(0);
1434   }
1435   if (m != B->rmap.n || n != B->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1436   if (lda1>m || lda2>m) {
1437     for (j=0; j<n; j++) {
1438       ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1439     }
1440   } else {
1441     ierr = PetscMemcpy(b->v,a->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
1442   }
1443   PetscFunctionReturn(0);
1444 }
1445 
1446 #undef __FUNCT__
1447 #define __FUNCT__ "MatSetUpPreallocation_SeqDense"
1448 PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A)
1449 {
1450   PetscErrorCode ierr;
1451 
1452   PetscFunctionBegin;
1453   ierr =  MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr);
1454   PetscFunctionReturn(0);
1455 }
1456 
1457 #undef __FUNCT__
1458 #define __FUNCT__ "MatSetSizes_SeqDense"
1459 PetscErrorCode MatSetSizes_SeqDense(Mat A,PetscInt m,PetscInt n,PetscInt M,PetscInt N)
1460 {
1461   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1462   PetscErrorCode ierr;
1463   PetscFunctionBegin;
1464   /* this will not be called before lda, Mmax,  and Nmax have been set */
1465   m = PetscMax(m,M);
1466   n = PetscMax(n,N);
1467   if (m > a->Mmax) SETERRQ2(PETSC_ERR_SUP,"Cannot yet resize number rows of dense matrix larger then its initial size %d, requested %d",a->lda,(int)m);
1468   if (n > a->Nmax) SETERRQ2(PETSC_ERR_SUP,"Cannot yet resize number columns of dense matrix larger then its initial size %d, requested %d",a->Nmax,(int)n);
1469   A->rmap.n = A->rmap.n = m;
1470   A->cmap.n = A->cmap.N = n;
1471   if (a->changelda) a->lda = m;
1472   ierr = PetscMemzero(a->v,m*n*sizeof(PetscScalar));CHKERRQ(ierr);
1473   PetscFunctionReturn(0);
1474 }
1475 
1476 /* ----------------------------------------------------------------*/
1477 #undef __FUNCT__
1478 #define __FUNCT__ "MatMatMult_SeqDense_SeqDense"
1479 PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1480 {
1481   PetscErrorCode ierr;
1482 
1483   PetscFunctionBegin;
1484   if (scall == MAT_INITIAL_MATRIX){
1485     ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1486   }
1487   ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1488   PetscFunctionReturn(0);
1489 }
1490 
1491 #undef __FUNCT__
1492 #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense"
1493 PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1494 {
1495   PetscErrorCode ierr;
1496   PetscInt       m=A->rmap.n,n=B->cmap.n;
1497   Mat            Cmat;
1498 
1499   PetscFunctionBegin;
1500   if (A->cmap.n != B->rmap.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"A->cmap.n %d != B->rmap.n %d\n",A->cmap.n,B->rmap.n);
1501   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1502   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1503   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1504   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1505   Cmat->assembled = PETSC_TRUE;
1506   *C = Cmat;
1507   PetscFunctionReturn(0);
1508 }
1509 
1510 #undef __FUNCT__
1511 #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense"
1512 PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1513 {
1514   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1515   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1516   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
1517   PetscBLASInt   m=(PetscBLASInt)A->rmap.n,n=(PetscBLASInt)B->cmap.n,k=(PetscBLASInt)A->cmap.n;
1518   PetscScalar    _DOne=1.0,_DZero=0.0;
1519 
1520   PetscFunctionBegin;
1521   BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1522   PetscFunctionReturn(0);
1523 }
1524 
1525 #undef __FUNCT__
1526 #define __FUNCT__ "MatMatMultTranspose_SeqDense_SeqDense"
1527 PetscErrorCode MatMatMultTranspose_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1528 {
1529   PetscErrorCode ierr;
1530 
1531   PetscFunctionBegin;
1532   if (scall == MAT_INITIAL_MATRIX){
1533     ierr = MatMatMultTransposeSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1534   }
1535   ierr = MatMatMultTransposeNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1536   PetscFunctionReturn(0);
1537 }
1538 
1539 #undef __FUNCT__
1540 #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqDense_SeqDense"
1541 PetscErrorCode MatMatMultTransposeSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1542 {
1543   PetscErrorCode ierr;
1544   PetscInt       m=A->cmap.n,n=B->cmap.n;
1545   Mat            Cmat;
1546 
1547   PetscFunctionBegin;
1548   if (A->rmap.n != B->rmap.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"A->rmap.n %d != B->rmap.n %d\n",A->rmap.n,B->rmap.n);
1549   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1550   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1551   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1552   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1553   Cmat->assembled = PETSC_TRUE;
1554   *C = Cmat;
1555   PetscFunctionReturn(0);
1556 }
1557 
1558 #undef __FUNCT__
1559 #define __FUNCT__ "MatMatMultTransposeNumeric_SeqDense_SeqDense"
1560 PetscErrorCode MatMatMultTransposeNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1561 {
1562   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1563   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1564   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
1565   PetscBLASInt   m=(PetscBLASInt)A->cmap.n,n=(PetscBLASInt)B->cmap.n,k=(PetscBLASInt)A->rmap.n;
1566   PetscScalar    _DOne=1.0,_DZero=0.0;
1567 
1568   PetscFunctionBegin;
1569   BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1570   PetscFunctionReturn(0);
1571 }
1572 /* -------------------------------------------------------------------*/
1573 static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1574        MatGetRow_SeqDense,
1575        MatRestoreRow_SeqDense,
1576        MatMult_SeqDense,
1577 /* 4*/ MatMultAdd_SeqDense,
1578        MatMultTranspose_SeqDense,
1579        MatMultTransposeAdd_SeqDense,
1580        MatSolve_SeqDense,
1581        MatSolveAdd_SeqDense,
1582        MatSolveTranspose_SeqDense,
1583 /*10*/ MatSolveTransposeAdd_SeqDense,
1584        MatLUFactor_SeqDense,
1585        MatCholeskyFactor_SeqDense,
1586        MatRelax_SeqDense,
1587        MatTranspose_SeqDense,
1588 /*15*/ MatGetInfo_SeqDense,
1589        MatEqual_SeqDense,
1590        MatGetDiagonal_SeqDense,
1591        MatDiagonalScale_SeqDense,
1592        MatNorm_SeqDense,
1593 /*20*/ MatAssemblyBegin_SeqDense,
1594        MatAssemblyEnd_SeqDense,
1595        0,
1596        MatSetOption_SeqDense,
1597        MatZeroEntries_SeqDense,
1598 /*25*/ MatZeroRows_SeqDense,
1599        MatLUFactorSymbolic_SeqDense,
1600        MatLUFactorNumeric_SeqDense,
1601        MatCholeskyFactorSymbolic_SeqDense,
1602        MatCholeskyFactorNumeric_SeqDense,
1603 /*30*/ MatSetUpPreallocation_SeqDense,
1604        0,
1605        0,
1606        MatGetArray_SeqDense,
1607        MatRestoreArray_SeqDense,
1608 /*35*/ MatDuplicate_SeqDense,
1609        0,
1610        0,
1611        0,
1612        0,
1613 /*40*/ MatAXPY_SeqDense,
1614        MatGetSubMatrices_SeqDense,
1615        0,
1616        MatGetValues_SeqDense,
1617        MatCopy_SeqDense,
1618 /*45*/ 0,
1619        MatScale_SeqDense,
1620        0,
1621        0,
1622        0,
1623 /*50*/ 0,
1624        0,
1625        0,
1626        0,
1627        0,
1628 /*55*/ 0,
1629        0,
1630        0,
1631        0,
1632        0,
1633 /*60*/ 0,
1634        MatDestroy_SeqDense,
1635        MatView_SeqDense,
1636        0,
1637        0,
1638 /*65*/ 0,
1639        0,
1640        0,
1641        0,
1642        0,
1643 /*70*/ 0,
1644        0,
1645        0,
1646        0,
1647        0,
1648 /*75*/ 0,
1649        0,
1650        0,
1651        0,
1652        0,
1653 /*80*/ 0,
1654        0,
1655        0,
1656        0,
1657 /*84*/ MatLoad_SeqDense,
1658        0,
1659        0,
1660        0,
1661        0,
1662        0,
1663 /*90*/ MatMatMult_SeqDense_SeqDense,
1664        MatMatMultSymbolic_SeqDense_SeqDense,
1665        MatMatMultNumeric_SeqDense_SeqDense,
1666        0,
1667        0,
1668 /*95*/ 0,
1669        MatMatMultTranspose_SeqDense_SeqDense,
1670        MatMatMultTransposeSymbolic_SeqDense_SeqDense,
1671        MatMatMultTransposeNumeric_SeqDense_SeqDense,
1672        0,
1673 /*100*/0,
1674        0,
1675        0,
1676        0,
1677        MatSetSizes_SeqDense};
1678 
1679 #undef __FUNCT__
1680 #define __FUNCT__ "MatCreateSeqDense"
1681 /*@C
1682    MatCreateSeqDense - Creates a sequential dense matrix that
1683    is stored in column major order (the usual Fortran 77 manner). Many
1684    of the matrix operations use the BLAS and LAPACK routines.
1685 
1686    Collective on MPI_Comm
1687 
1688    Input Parameters:
1689 +  comm - MPI communicator, set to PETSC_COMM_SELF
1690 .  m - number of rows
1691 .  n - number of columns
1692 -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1693    to control all matrix memory allocation.
1694 
1695    Output Parameter:
1696 .  A - the matrix
1697 
1698    Notes:
1699    The data input variable is intended primarily for Fortran programmers
1700    who wish to allocate their own matrix memory space.  Most users should
1701    set data=PETSC_NULL.
1702 
1703    Level: intermediate
1704 
1705 .keywords: dense, matrix, LAPACK, BLAS
1706 
1707 .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
1708 @*/
1709 PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A)
1710 {
1711   PetscErrorCode ierr;
1712 
1713   PetscFunctionBegin;
1714   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1715   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
1716   ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1717   ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1718   PetscFunctionReturn(0);
1719 }
1720 
1721 #undef __FUNCT__
1722 #define __FUNCT__ "MatSeqDensePreallocation"
1723 /*@C
1724    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
1725 
1726    Collective on MPI_Comm
1727 
1728    Input Parameters:
1729 +  A - the matrix
1730 -  data - the array (or PETSC_NULL)
1731 
1732    Notes:
1733    The data input variable is intended primarily for Fortran programmers
1734    who wish to allocate their own matrix memory space.  Most users should
1735    need not call this routine.
1736 
1737    Level: intermediate
1738 
1739 .keywords: dense, matrix, LAPACK, BLAS
1740 
1741 .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
1742 @*/
1743 PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
1744 {
1745   PetscErrorCode ierr,(*f)(Mat,PetscScalar[]);
1746 
1747   PetscFunctionBegin;
1748   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1749   if (f) {
1750     ierr = (*f)(B,data);CHKERRQ(ierr);
1751   }
1752   PetscFunctionReturn(0);
1753 }
1754 
1755 EXTERN_C_BEGIN
1756 #undef __FUNCT__
1757 #define __FUNCT__ "MatSeqDensePreallocation_SeqDense"
1758 PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
1759 {
1760   Mat_SeqDense   *b;
1761   PetscErrorCode ierr;
1762 
1763   PetscFunctionBegin;
1764   B->preallocated = PETSC_TRUE;
1765   b               = (Mat_SeqDense*)B->data;
1766   if (!data) {
1767     ierr          = PetscMalloc((b->lda*b->Nmax+1)*sizeof(PetscScalar),&b->v);CHKERRQ(ierr);
1768     ierr          = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
1769     b->user_alloc = PETSC_FALSE;
1770     ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
1771   } else { /* user-allocated storage */
1772     b->v          = data;
1773     b->user_alloc = PETSC_TRUE;
1774   }
1775   PetscFunctionReturn(0);
1776 }
1777 EXTERN_C_END
1778 
1779 #undef __FUNCT__
1780 #define __FUNCT__ "MatSeqDenseSetLDA"
1781 /*@C
1782   MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array
1783 
1784   Input parameter:
1785 + A - the matrix
1786 - lda - the leading dimension
1787 
1788   Notes:
1789   This routine is to be used in conjunction with MatSeqDenseSetPreallocation;
1790   it asserts that the preallocation has a leading dimension (the LDA parameter
1791   of Blas and Lapack fame) larger than M, the first dimension of the matrix.
1792 
1793   Level: intermediate
1794 
1795 .keywords: dense, matrix, LAPACK, BLAS
1796 
1797 .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize()
1798 @*/
1799 PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat B,PetscInt lda)
1800 {
1801   Mat_SeqDense *b = (Mat_SeqDense*)B->data;
1802 
1803   PetscFunctionBegin;
1804   if (lda < B->rmap.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->rmap.n);
1805   b->lda       = lda;
1806   b->changelda = PETSC_FALSE;
1807   b->Mmax      = PetscMax(b->Mmax,lda);
1808   PetscFunctionReturn(0);
1809 }
1810 
1811 /*MC
1812    MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices.
1813 
1814    Options Database Keys:
1815 . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions()
1816 
1817   Level: beginner
1818 
1819 .seealso: MatCreateSeqDense()
1820 
1821 M*/
1822 
1823 EXTERN_C_BEGIN
1824 #undef __FUNCT__
1825 #define __FUNCT__ "MatCreate_SeqDense"
1826 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqDense(Mat B)
1827 {
1828   Mat_SeqDense   *b;
1829   PetscErrorCode ierr;
1830   PetscMPIInt    size;
1831 
1832   PetscFunctionBegin;
1833   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
1834   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1835 
1836   B->rmap.bs = B->cmap.bs = 1;
1837   ierr = PetscMapInitialize(B->comm,&B->rmap);CHKERRQ(ierr);
1838   ierr = PetscMapInitialize(B->comm,&B->cmap);CHKERRQ(ierr);
1839 
1840   ierr            = PetscNew(Mat_SeqDense,&b);CHKERRQ(ierr);
1841   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1842   B->factor       = 0;
1843   B->mapping      = 0;
1844   ierr = PetscLogObjectMemory(B,sizeof(struct _p_Mat));CHKERRQ(ierr);
1845   B->data         = (void*)b;
1846 
1847 
1848   b->pivots       = 0;
1849   b->roworiented  = PETSC_TRUE;
1850   b->v            = 0;
1851   b->lda          = B->rmap.n;
1852   b->changelda    = PETSC_FALSE;
1853   b->Mmax         = B->rmap.n;
1854   b->Nmax         = B->cmap.n;
1855 
1856   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C",
1857                                     "MatSeqDenseSetPreallocation_SeqDense",
1858                                      MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr);
1859   PetscFunctionReturn(0);
1860 }
1861 
1862 
1863 EXTERN_C_END
1864