xref: /petsc/src/mat/impls/dense/seq/dense.c (revision c5df96a59d3f7e354f9770fb35f00818af2dc9fc)
1 
2 /*
3      Defines the basic matrix operations for sequential dense.
4 */
5 
6 #include <../src/mat/impls/dense/seq/dense.h> /*I "petscmat.h" I*/
7 #include <petscblaslapack.h>
8 
9 #include <../src/mat/impls/aij/seq/aij.h>
10 EXTERN_C_BEGIN
11 #undef __FUNCT__
12 #define __FUNCT__ "MatConvert_SeqDense_SeqAIJ"
13 PetscErrorCode  MatConvert_SeqDense_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
14 {
15   Mat            B;
16   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
17   PetscErrorCode ierr;
18   PetscInt       i;
19   PetscInt       *rows;
20   MatScalar      *aa = a->v;
21 
22   PetscFunctionBegin;
23   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
24   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
25   ierr = MatSetType(B,MATSEQAIJ);CHKERRQ(ierr);
26   ierr = MatSeqAIJSetPreallocation(B,A->cmap->n,PETSC_NULL);CHKERRQ(ierr);
27 
28   ierr = PetscMalloc(A->rmap->n*sizeof(PetscInt),&rows);CHKERRQ(ierr);
29   for (i=0; i<A->rmap->n; i++) rows[i] = i;
30 
31   for (i=0; i<A->cmap->n; i++) {
32     ierr  = MatSetValues(B,A->rmap->n,rows,1,&i,aa,INSERT_VALUES);CHKERRQ(ierr);
33     aa   += a->lda;
34   }
35   ierr = PetscFree(rows);CHKERRQ(ierr);
36   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
37   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
38 
39   if (reuse == MAT_REUSE_MATRIX) {
40     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
41   } else {
42     *newmat = B;
43   }
44   PetscFunctionReturn(0);
45 }
46 EXTERN_C_END
47 
48 #undef __FUNCT__
49 #define __FUNCT__ "MatAXPY_SeqDense"
50 PetscErrorCode MatAXPY_SeqDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str)
51 {
52   Mat_SeqDense   *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data;
53   PetscScalar    oalpha = alpha;
54   PetscInt       j;
55   PetscBLASInt   N,m,ldax,lday,one = 1;
56   PetscErrorCode ierr;
57 
58   PetscFunctionBegin;
59   ierr = PetscBLASIntCast(X->rmap->n*X->cmap->n,&N);CHKERRQ(ierr);
60   ierr = PetscBLASIntCast(X->rmap->n,&m);CHKERRQ(ierr);
61   ierr = PetscBLASIntCast(x->lda,&ldax);CHKERRQ(ierr);
62   ierr = PetscBLASIntCast(y->lda,&lday);CHKERRQ(ierr);
63   if (ldax>m || lday>m) {
64     for (j=0; j<X->cmap->n; j++) {
65       BLASaxpy_(&m,&oalpha,x->v+j*ldax,&one,y->v+j*lday,&one);
66     }
67   } else {
68     BLASaxpy_(&N,&oalpha,x->v,&one,y->v,&one);
69   }
70   ierr = PetscLogFlops(PetscMax(2*N-1,0));CHKERRQ(ierr);
71   PetscFunctionReturn(0);
72 }
73 
74 #undef __FUNCT__
75 #define __FUNCT__ "MatGetInfo_SeqDense"
76 PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
77 {
78   PetscInt     N = A->rmap->n*A->cmap->n;
79 
80   PetscFunctionBegin;
81   info->block_size        = 1.0;
82   info->nz_allocated      = (double)N;
83   info->nz_used           = (double)N;
84   info->nz_unneeded       = (double)0;
85   info->assemblies        = (double)A->num_ass;
86   info->mallocs           = 0;
87   info->memory            = ((PetscObject)A)->mem;
88   info->fill_ratio_given  = 0;
89   info->fill_ratio_needed = 0;
90   info->factor_mallocs    = 0;
91   PetscFunctionReturn(0);
92 }
93 
94 #undef __FUNCT__
95 #define __FUNCT__ "MatScale_SeqDense"
96 PetscErrorCode MatScale_SeqDense(Mat A,PetscScalar alpha)
97 {
98   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
99   PetscScalar    oalpha = alpha;
100   PetscErrorCode ierr;
101   PetscBLASInt   one = 1,j,nz,lda;
102 
103   PetscFunctionBegin;
104   ierr = PetscBLASIntCast(a->lda,&lda);CHKERRQ(ierr);
105   if (lda>A->rmap->n) {
106     ierr = PetscBLASIntCast(A->rmap->n,&nz);CHKERRQ(ierr);
107     for (j=0; j<A->cmap->n; j++) {
108       BLASscal_(&nz,&oalpha,a->v+j*lda,&one);
109     }
110   } else {
111     ierr = PetscBLASIntCast(A->rmap->n*A->cmap->n,&nz);CHKERRQ(ierr);
112     BLASscal_(&nz,&oalpha,a->v,&one);
113   }
114   ierr = PetscLogFlops(nz);CHKERRQ(ierr);
115   PetscFunctionReturn(0);
116 }
117 
118 #undef __FUNCT__
119 #define __FUNCT__ "MatIsHermitian_SeqDense"
120 PetscErrorCode MatIsHermitian_SeqDense(Mat A,PetscReal rtol,PetscBool  *fl)
121 {
122   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
123   PetscInt       i,j,m = A->rmap->n,N;
124   PetscScalar    *v = a->v;
125 
126   PetscFunctionBegin;
127   *fl = PETSC_FALSE;
128   if (A->rmap->n != A->cmap->n) PetscFunctionReturn(0);
129   N = a->lda;
130 
131   for (i=0; i<m; i++) {
132     for (j=i+1; j<m; j++) {
133       if (PetscAbsScalar(v[i+j*N] - PetscConj(v[j+i*N])) > rtol) PetscFunctionReturn(0);
134     }
135   }
136   *fl = PETSC_TRUE;
137   PetscFunctionReturn(0);
138 }
139 
140 #undef __FUNCT__
141 #define __FUNCT__ "MatDuplicateNoCreate_SeqDense"
142 PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat newi,Mat A,MatDuplicateOption cpvalues)
143 {
144   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data,*l;
145   PetscErrorCode ierr;
146   PetscInt       lda = (PetscInt)mat->lda,j,m;
147 
148   PetscFunctionBegin;
149   ierr = PetscLayoutReference(A->rmap,&newi->rmap);CHKERRQ(ierr);
150   ierr = PetscLayoutReference(A->cmap,&newi->cmap);CHKERRQ(ierr);
151   ierr = MatSeqDenseSetPreallocation(newi,PETSC_NULL);CHKERRQ(ierr);
152   if (cpvalues == MAT_COPY_VALUES) {
153     l = (Mat_SeqDense*)newi->data;
154     if (lda>A->rmap->n) {
155       m = A->rmap->n;
156       for (j=0; j<A->cmap->n; j++) {
157         ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
158       }
159     } else {
160       ierr = PetscMemcpy(l->v,mat->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
161     }
162   }
163   newi->assembled = PETSC_TRUE;
164   PetscFunctionReturn(0);
165 }
166 
167 #undef __FUNCT__
168 #define __FUNCT__ "MatDuplicate_SeqDense"
169 PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
170 {
171   PetscErrorCode ierr;
172 
173   PetscFunctionBegin;
174   ierr = MatCreate(((PetscObject)A)->comm,newmat);CHKERRQ(ierr);
175   ierr = MatSetSizes(*newmat,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr);
176   ierr = MatSetType(*newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
177   ierr = MatDuplicateNoCreate_SeqDense(*newmat,A,cpvalues);CHKERRQ(ierr);
178   PetscFunctionReturn(0);
179 }
180 
181 
182 extern PetscErrorCode MatLUFactor_SeqDense(Mat,IS,IS,const MatFactorInfo*);
183 
184 #undef __FUNCT__
185 #define __FUNCT__ "MatLUFactorNumeric_SeqDense"
186 PetscErrorCode MatLUFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy)
187 {
188   MatFactorInfo  info;
189   PetscErrorCode ierr;
190 
191   PetscFunctionBegin;
192   ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr);
193   ierr = MatLUFactor_SeqDense(fact,0,0,&info);CHKERRQ(ierr);
194   PetscFunctionReturn(0);
195 }
196 
197 #undef __FUNCT__
198 #define __FUNCT__ "MatSolve_SeqDense"
199 PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
200 {
201   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
202   PetscErrorCode ierr;
203   PetscScalar    *x,*y;
204   PetscBLASInt   one = 1,info,m;
205 
206   PetscFunctionBegin;
207   ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr);
208   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
209   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
210   ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
211   if (A->factortype == MAT_FACTOR_LU) {
212 #if defined(PETSC_MISSING_LAPACK_GETRS)
213     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
214 #else
215     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
216     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve");
217 #endif
218   } else if (A->factortype == MAT_FACTOR_CHOLESKY) {
219 #if defined(PETSC_MISSING_LAPACK_POTRS)
220     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
221 #else
222     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
223     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve");
224 #endif
225   }
226   else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
227   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
228   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
229   ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr);
230   PetscFunctionReturn(0);
231 }
232 
233 #undef __FUNCT__
234 #define __FUNCT__ "MatMatSolve_SeqDense"
235 PetscErrorCode MatMatSolve_SeqDense(Mat A,Mat B,Mat X)
236 {
237   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
238   PetscErrorCode ierr;
239   PetscScalar    *b,*x;
240   PetscInt       n;
241   PetscBLASInt   nrhs,info,m;
242   PetscBool      flg;
243 
244   PetscFunctionBegin;
245   ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr);
246   ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,PETSC_NULL);CHKERRQ(ierr);
247   if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
248   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,PETSC_NULL);CHKERRQ(ierr);
249   if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
250 
251   ierr = MatGetSize(B,PETSC_NULL,&n);CHKERRQ(ierr);
252   ierr = PetscBLASIntCast(n,&nrhs);CHKERRQ(ierr);
253   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
254   ierr = MatDenseGetArray(X,&x);CHKERRQ(ierr);
255 
256   ierr = PetscMemcpy(x,b,m*nrhs*sizeof(PetscScalar));CHKERRQ(ierr);
257 
258   if (A->factortype == MAT_FACTOR_LU) {
259 #if defined(PETSC_MISSING_LAPACK_GETRS)
260     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
261 #else
262     LAPACKgetrs_("N",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info);
263     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve");
264 #endif
265   } else if (A->factortype == MAT_FACTOR_CHOLESKY) {
266 #if defined(PETSC_MISSING_LAPACK_POTRS)
267     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
268 #else
269     LAPACKpotrs_("L",&m,&nrhs,mat->v,&mat->lda,x,&m,&info);
270     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve");
271 #endif
272   }
273   else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
274 
275   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
276   ierr = MatDenseRestoreArray(X,&x);CHKERRQ(ierr);
277   ierr = PetscLogFlops(nrhs*(2.0*m*m - m));CHKERRQ(ierr);
278   PetscFunctionReturn(0);
279 }
280 
281 #undef __FUNCT__
282 #define __FUNCT__ "MatSolveTranspose_SeqDense"
283 PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
284 {
285   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
286   PetscErrorCode ierr;
287   PetscScalar    *x,*y;
288   PetscBLASInt   one = 1,info,m;
289 
290   PetscFunctionBegin;
291   ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr);
292   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
293   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
294   ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
295   /* assume if pivots exist then use LU; else Cholesky */
296   if (mat->pivots) {
297 #if defined(PETSC_MISSING_LAPACK_GETRS)
298     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
299 #else
300     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
301     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve");
302 #endif
303   } else {
304 #if defined(PETSC_MISSING_LAPACK_POTRS)
305     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
306 #else
307     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
308     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve");
309 #endif
310   }
311   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
312   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
313   ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr);
314   PetscFunctionReturn(0);
315 }
316 
317 #undef __FUNCT__
318 #define __FUNCT__ "MatSolveAdd_SeqDense"
319 PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
320 {
321   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
322   PetscErrorCode ierr;
323   PetscScalar    *x,*y,sone = 1.0;
324   Vec            tmp = 0;
325   PetscBLASInt   one = 1,info,m;
326 
327   PetscFunctionBegin;
328   ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr);
329   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
330   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
331   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
332   if (yy == zz) {
333     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
334     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
335     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
336   }
337   ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
338   /* assume if pivots exist then use LU; else Cholesky */
339   if (mat->pivots) {
340 #if defined(PETSC_MISSING_LAPACK_GETRS)
341     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
342 #else
343     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
344     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve");
345 #endif
346   } else {
347 #if defined(PETSC_MISSING_LAPACK_POTRS)
348     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
349 #else
350     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
351     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve");
352 #endif
353   }
354   if (tmp) {
355     ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr);
356     ierr = VecDestroy(&tmp);CHKERRQ(ierr);
357   } else {
358     ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);
359   }
360   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
361   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
362   ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n);CHKERRQ(ierr);
363   PetscFunctionReturn(0);
364 }
365 
366 #undef __FUNCT__
367 #define __FUNCT__ "MatSolveTransposeAdd_SeqDense"
368 PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
369 {
370   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
371   PetscErrorCode ierr;
372   PetscScalar    *x,*y,sone = 1.0;
373   Vec            tmp;
374   PetscBLASInt   one = 1,info,m;
375 
376   PetscFunctionBegin;
377   ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr);
378   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
379   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
380   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
381   if (yy == zz) {
382     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
383     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
384     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
385   }
386   ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
387   /* assume if pivots exist then use LU; else Cholesky */
388   if (mat->pivots) {
389 #if defined(PETSC_MISSING_LAPACK_GETRS)
390     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
391 #else
392     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
393     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve");
394 #endif
395   } else {
396 #if defined(PETSC_MISSING_LAPACK_POTRS)
397     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
398 #else
399     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
400     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve");
401 #endif
402   }
403   if (tmp) {
404     ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr);
405     ierr = VecDestroy(&tmp);CHKERRQ(ierr);
406   } else {
407     ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);
408   }
409   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
410   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
411   ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n);CHKERRQ(ierr);
412   PetscFunctionReturn(0);
413 }
414 
415 /* ---------------------------------------------------------------*/
416 /* COMMENT: I have chosen to hide row permutation in the pivots,
417    rather than put it in the Mat->row slot.*/
418 #undef __FUNCT__
419 #define __FUNCT__ "MatLUFactor_SeqDense"
420 PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,const MatFactorInfo *minfo)
421 {
422 #if defined(PETSC_MISSING_LAPACK_GETRF)
423   PetscFunctionBegin;
424   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable.");
425 #else
426   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
427   PetscErrorCode ierr;
428   PetscBLASInt   n,m,info;
429 
430   PetscFunctionBegin;
431   ierr = PetscBLASIntCast(A->cmap->n,&n);CHKERRQ(ierr);
432   ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr);
433   if (!mat->pivots) {
434     ierr = PetscMalloc((A->rmap->n+1)*sizeof(PetscBLASInt),&mat->pivots);CHKERRQ(ierr);
435     ierr = PetscLogObjectMemory(A,A->rmap->n*sizeof(PetscBLASInt));CHKERRQ(ierr);
436   }
437   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
438   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
439   LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info);
440   ierr = PetscFPTrapPop();CHKERRQ(ierr);
441 
442   if (info<0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to LU factorization");
443   if (info>0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
444   A->ops->solve             = MatSolve_SeqDense;
445   A->ops->solvetranspose    = MatSolveTranspose_SeqDense;
446   A->ops->solveadd          = MatSolveAdd_SeqDense;
447   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense;
448   A->factortype             = MAT_FACTOR_LU;
449 
450   ierr = PetscLogFlops((2.0*A->cmap->n*A->cmap->n*A->cmap->n)/3);CHKERRQ(ierr);
451 #endif
452   PetscFunctionReturn(0);
453 }
454 
455 #undef __FUNCT__
456 #define __FUNCT__ "MatCholeskyFactor_SeqDense"
457 PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,const MatFactorInfo *factinfo)
458 {
459 #if defined(PETSC_MISSING_LAPACK_POTRF)
460   PetscFunctionBegin;
461   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable.");
462 #else
463   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
464   PetscErrorCode ierr;
465   PetscBLASInt   info,n;
466 
467   PetscFunctionBegin;
468   ierr = PetscBLASIntCast(A->cmap->n,&n);CHKERRQ(ierr);
469   ierr = PetscFree(mat->pivots);CHKERRQ(ierr);
470 
471   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
472   LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info);
473   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1);
474   A->ops->solve             = MatSolve_SeqDense;
475   A->ops->solvetranspose    = MatSolveTranspose_SeqDense;
476   A->ops->solveadd          = MatSolveAdd_SeqDense;
477   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense;
478   A->factortype             = MAT_FACTOR_CHOLESKY;
479   ierr = PetscLogFlops((A->cmap->n*A->cmap->n*A->cmap->n)/3.0);CHKERRQ(ierr);
480 #endif
481   PetscFunctionReturn(0);
482 }
483 
484 
485 #undef __FUNCT__
486 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense"
487 PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy)
488 {
489   PetscErrorCode ierr;
490   MatFactorInfo  info;
491 
492   PetscFunctionBegin;
493   info.fill = 1.0;
494   ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr);
495   ierr = MatCholeskyFactor_SeqDense(fact,0,&info);CHKERRQ(ierr);
496   PetscFunctionReturn(0);
497 }
498 
499 #undef __FUNCT__
500 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense"
501 PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,const MatFactorInfo *info)
502 {
503   PetscFunctionBegin;
504   fact->assembled                  = PETSC_TRUE;
505   fact->preallocated               = PETSC_TRUE;
506   fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqDense;
507   PetscFunctionReturn(0);
508 }
509 
510 #undef __FUNCT__
511 #define __FUNCT__ "MatLUFactorSymbolic_SeqDense"
512 PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info)
513 {
514   PetscFunctionBegin;
515   fact->preallocated         = PETSC_TRUE;
516   fact->assembled            = PETSC_TRUE;
517   fact->ops->lufactornumeric = MatLUFactorNumeric_SeqDense;
518   PetscFunctionReturn(0);
519 }
520 
521 EXTERN_C_BEGIN
522 #undef __FUNCT__
523 #define __FUNCT__ "MatGetFactor_seqdense_petsc"
524 PetscErrorCode MatGetFactor_seqdense_petsc(Mat A,MatFactorType ftype,Mat *fact)
525 {
526   PetscErrorCode ierr;
527 
528   PetscFunctionBegin;
529   ierr = MatCreate(((PetscObject)A)->comm,fact);CHKERRQ(ierr);
530   ierr = MatSetSizes(*fact,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr);
531   ierr = MatSetType(*fact,((PetscObject)A)->type_name);CHKERRQ(ierr);
532   if (ftype == MAT_FACTOR_LU) {
533     (*fact)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqDense;
534   } else {
535     (*fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqDense;
536   }
537   (*fact)->factortype = ftype;
538   PetscFunctionReturn(0);
539 }
540 EXTERN_C_END
541 
542 /* ------------------------------------------------------------------*/
543 #undef __FUNCT__
544 #define __FUNCT__ "MatSOR_SeqDense"
545 PetscErrorCode MatSOR_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx)
546 {
547   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
548   PetscScalar    *x,*b,*v = mat->v,zero = 0.0,xt;
549   PetscErrorCode ierr;
550   PetscInt       m = A->rmap->n,i;
551   PetscBLASInt   o = 1,bm;
552 
553   PetscFunctionBegin;
554   ierr = PetscBLASIntCast(m,&bm);CHKERRQ(ierr);
555   if (flag & SOR_ZERO_INITIAL_GUESS) {
556     /* this is a hack fix, should have another version without the second BLASdot */
557     ierr = VecSet(xx,zero);CHKERRQ(ierr);
558   }
559   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
560   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
561   its  = its*lits;
562   if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
563   while (its--) {
564     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
565       for (i=0; i<m; i++) {
566         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
567         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
568       }
569     }
570     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
571       for (i=m-1; i>=0; i--) {
572         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
573         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
574       }
575     }
576   }
577   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
578   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
579   PetscFunctionReturn(0);
580 }
581 
582 /* -----------------------------------------------------------------*/
583 #undef __FUNCT__
584 #define __FUNCT__ "MatMultTranspose_SeqDense"
585 PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
586 {
587   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
588   PetscScalar    *v = mat->v,*x,*y;
589   PetscErrorCode ierr;
590   PetscBLASInt   m, n,_One=1;
591   PetscScalar    _DOne=1.0,_DZero=0.0;
592 
593   PetscFunctionBegin;
594   ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr);
595   ierr = PetscBLASIntCast(A->cmap->n,&n);CHKERRQ(ierr);
596   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
597   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
598   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
599   BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One);
600   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
601   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
602   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr);
603   PetscFunctionReturn(0);
604 }
605 
606 #undef __FUNCT__
607 #define __FUNCT__ "MatMult_SeqDense"
608 PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy)
609 {
610   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
611   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0;
612   PetscErrorCode ierr;
613   PetscBLASInt   m, n, _One=1;
614 
615   PetscFunctionBegin;
616   ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr);
617   ierr = PetscBLASIntCast(A->cmap->n,&n);CHKERRQ(ierr);
618   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
619   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
620   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
621   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One);
622   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
623   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
624   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->rmap->n);CHKERRQ(ierr);
625   PetscFunctionReturn(0);
626 }
627 
628 #undef __FUNCT__
629 #define __FUNCT__ "MatMultAdd_SeqDense"
630 PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
631 {
632   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
633   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0;
634   PetscErrorCode ierr;
635   PetscBLASInt   m, n, _One=1;
636 
637   PetscFunctionBegin;
638   ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr);
639   ierr = PetscBLASIntCast(A->cmap->n,&n);CHKERRQ(ierr);
640   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
641   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
642   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
643   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
644   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
645   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
646   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
647   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr);
648   PetscFunctionReturn(0);
649 }
650 
651 #undef __FUNCT__
652 #define __FUNCT__ "MatMultTransposeAdd_SeqDense"
653 PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
654 {
655   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
656   PetscScalar    *v = mat->v,*x,*y;
657   PetscErrorCode ierr;
658   PetscBLASInt   m, n, _One=1;
659   PetscScalar    _DOne=1.0;
660 
661   PetscFunctionBegin;
662   ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr);
663   ierr = PetscBLASIntCast(A->cmap->n,&n);CHKERRQ(ierr);
664   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
665   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
666   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
667   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
668   BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
669   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
670   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
671   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr);
672   PetscFunctionReturn(0);
673 }
674 
675 /* -----------------------------------------------------------------*/
676 #undef __FUNCT__
677 #define __FUNCT__ "MatGetRow_SeqDense"
678 PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
679 {
680   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
681   PetscScalar    *v;
682   PetscErrorCode ierr;
683   PetscInt       i;
684 
685   PetscFunctionBegin;
686   *ncols = A->cmap->n;
687   if (cols) {
688     ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr);
689     for (i=0; i<A->cmap->n; i++) (*cols)[i] = i;
690   }
691   if (vals) {
692     ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr);
693     v    = mat->v + row;
694     for (i=0; i<A->cmap->n; i++) {(*vals)[i] = *v; v += mat->lda;}
695   }
696   PetscFunctionReturn(0);
697 }
698 
699 #undef __FUNCT__
700 #define __FUNCT__ "MatRestoreRow_SeqDense"
701 PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
702 {
703   PetscErrorCode ierr;
704 
705   PetscFunctionBegin;
706   if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);}
707   if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); }
708   PetscFunctionReturn(0);
709 }
710 /* ----------------------------------------------------------------*/
711 #undef __FUNCT__
712 #define __FUNCT__ "MatSetValues_SeqDense"
713 PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv)
714 {
715   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
716   PetscInt     i,j,idx=0;
717 
718   PetscFunctionBegin;
719   if (v) PetscValidScalarPointer(v,6);
720   if (!mat->roworiented) {
721     if (addv == INSERT_VALUES) {
722       for (j=0; j<n; j++) {
723         if (indexn[j] < 0) {idx += m; continue;}
724 #if defined(PETSC_USE_DEBUG)
725         if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
726 #endif
727         for (i=0; i<m; i++) {
728           if (indexm[i] < 0) {idx++; continue;}
729 #if defined(PETSC_USE_DEBUG)
730           if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
731 #endif
732           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
733         }
734       }
735     } else {
736       for (j=0; j<n; j++) {
737         if (indexn[j] < 0) {idx += m; continue;}
738 #if defined(PETSC_USE_DEBUG)
739         if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
740 #endif
741         for (i=0; i<m; i++) {
742           if (indexm[i] < 0) {idx++; continue;}
743 #if defined(PETSC_USE_DEBUG)
744           if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
745 #endif
746           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
747         }
748       }
749     }
750   } else {
751     if (addv == INSERT_VALUES) {
752       for (i=0; i<m; i++) {
753         if (indexm[i] < 0) { idx += n; continue;}
754 #if defined(PETSC_USE_DEBUG)
755         if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
756 #endif
757         for (j=0; j<n; j++) {
758           if (indexn[j] < 0) { idx++; continue;}
759 #if defined(PETSC_USE_DEBUG)
760           if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
761 #endif
762           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
763         }
764       }
765     } else {
766       for (i=0; i<m; i++) {
767         if (indexm[i] < 0) { idx += n; continue;}
768 #if defined(PETSC_USE_DEBUG)
769         if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
770 #endif
771         for (j=0; j<n; j++) {
772           if (indexn[j] < 0) { idx++; continue;}
773 #if defined(PETSC_USE_DEBUG)
774           if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
775 #endif
776           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
777         }
778       }
779     }
780   }
781   PetscFunctionReturn(0);
782 }
783 
784 #undef __FUNCT__
785 #define __FUNCT__ "MatGetValues_SeqDense"
786 PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[])
787 {
788   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
789   PetscInt     i,j;
790 
791   PetscFunctionBegin;
792   /* row-oriented output */
793   for (i=0; i<m; i++) {
794     if (indexm[i] < 0) {v += n;continue;}
795     if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %D requested larger than number rows %D",indexm[i],A->rmap->n);
796     for (j=0; j<n; j++) {
797       if (indexn[j] < 0) {v++; continue;}
798       if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column %D requested larger than number columns %D",indexn[j],A->cmap->n);
799       *v++ = mat->v[indexn[j]*mat->lda + indexm[i]];
800     }
801   }
802   PetscFunctionReturn(0);
803 }
804 
805 /* -----------------------------------------------------------------*/
806 
807 #undef __FUNCT__
808 #define __FUNCT__ "MatLoad_SeqDense"
809 PetscErrorCode MatLoad_SeqDense(Mat newmat,PetscViewer viewer)
810 {
811   Mat_SeqDense   *a;
812   PetscErrorCode ierr;
813   PetscInt       *scols,i,j,nz,header[4];
814   int            fd;
815   PetscMPIInt    size;
816   PetscInt       *rowlengths = 0,M,N,*cols,grows,gcols;
817   PetscScalar    *vals,*svals,*v,*w;
818   MPI_Comm       comm = ((PetscObject)viewer)->comm;
819 
820   PetscFunctionBegin;
821   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
822   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"view must have one processor");
823   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
824   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
825   if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
826   M = header[1]; N = header[2]; nz = header[3];
827 
828   /* set global size if not set already*/
829   if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) {
830     ierr = MatSetSizes(newmat,M,N,M,N);CHKERRQ(ierr);
831   } else {
832     /* if sizes and type are already set, check if the vector global sizes are correct */
833     ierr = MatGetSize(newmat,&grows,&gcols);CHKERRQ(ierr);
834     if (M != grows ||  N != gcols) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Matrix in file of different length (%d, %d) than the input matrix (%d, %d)",M,N,grows,gcols);
835   }
836   ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
837 
838   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
839     a    = (Mat_SeqDense*)newmat->data;
840     v    = a->v;
841     /* Allocate some temp space to read in the values and then flip them
842        from row major to column major */
843     ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr);
844     /* read in nonzero values */
845     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr);
846     /* now flip the values and store them in the matrix*/
847     for (j=0; j<N; j++) {
848       for (i=0; i<M; i++) {
849         *v++ =w[i*N+j];
850       }
851     }
852     ierr = PetscFree(w);CHKERRQ(ierr);
853   } else {
854     /* read row lengths */
855     ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
856     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
857 
858     a = (Mat_SeqDense*)newmat->data;
859     v = a->v;
860 
861     /* read column indices and nonzeros */
862     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr);
863     cols = scols;
864     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
865     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr);
866     vals = svals;
867     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
868 
869     /* insert into matrix */
870     for (i=0; i<M; i++) {
871       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
872       svals += rowlengths[i]; scols += rowlengths[i];
873     }
874     ierr = PetscFree(vals);CHKERRQ(ierr);
875     ierr = PetscFree(cols);CHKERRQ(ierr);
876     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
877   }
878   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
879   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
880 
881   PetscFunctionReturn(0);
882 }
883 
884 #undef __FUNCT__
885 #define __FUNCT__ "MatView_SeqDense_ASCII"
886 static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
887 {
888   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
889   PetscErrorCode    ierr;
890   PetscInt          i,j;
891   const char        *name;
892   PetscScalar       *v;
893   PetscViewerFormat format;
894 #if defined(PETSC_USE_COMPLEX)
895   PetscBool         allreal = PETSC_TRUE;
896 #endif
897 
898   PetscFunctionBegin;
899   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
900   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
901     PetscFunctionReturn(0);  /* do nothing for now */
902   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
903     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
904     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)A,viewer,"Matrix Object");CHKERRQ(ierr);
905     for (i=0; i<A->rmap->n; i++) {
906       v = a->v + i;
907       ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr);
908       for (j=0; j<A->cmap->n; j++) {
909 #if defined(PETSC_USE_COMPLEX)
910         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
911           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
912         } else if (PetscRealPart(*v)) {
913           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,PetscRealPart(*v));CHKERRQ(ierr);
914         }
915 #else
916         if (*v) {
917           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,*v);CHKERRQ(ierr);
918         }
919 #endif
920         v += a->lda;
921       }
922       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
923     }
924     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
925   } else {
926     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
927 #if defined(PETSC_USE_COMPLEX)
928     /* determine if matrix has all real values */
929     v = a->v;
930     for (i=0; i<A->rmap->n*A->cmap->n; i++) {
931         if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;}
932     }
933 #endif
934     if (format == PETSC_VIEWER_ASCII_MATLAB) {
935       ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
936       ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap->n,A->cmap->n);CHKERRQ(ierr);
937       ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap->n,A->cmap->n);CHKERRQ(ierr);
938       ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr);
939     } else {
940       ierr = PetscObjectPrintClassNamePrefixType((PetscObject)A,viewer,"Matrix Object");CHKERRQ(ierr);
941     }
942 
943     for (i=0; i<A->rmap->n; i++) {
944       v = a->v + i;
945       for (j=0; j<A->cmap->n; j++) {
946 #if defined(PETSC_USE_COMPLEX)
947         if (allreal) {
948           ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)PetscRealPart(*v));CHKERRQ(ierr);
949         } else {
950           ierr = PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16ei ",(double)PetscRealPart(*v),(double)PetscImaginaryPart(*v));CHKERRQ(ierr);
951         }
952 #else
953         ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)*v);CHKERRQ(ierr);
954 #endif
955         v += a->lda;
956       }
957       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
958     }
959     if (format == PETSC_VIEWER_ASCII_MATLAB) {
960       ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr);
961     }
962     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
963   }
964   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
965   PetscFunctionReturn(0);
966 }
967 
968 #undef __FUNCT__
969 #define __FUNCT__ "MatView_SeqDense_Binary"
970 static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
971 {
972   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
973   PetscErrorCode    ierr;
974   int               fd;
975   PetscInt          ict,j,n = A->cmap->n,m = A->rmap->n,i,*col_lens,nz = m*n;
976   PetscScalar       *v,*anonz,*vals;
977   PetscViewerFormat format;
978 
979   PetscFunctionBegin;
980   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
981 
982   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
983   if (format == PETSC_VIEWER_NATIVE) {
984     /* store the matrix as a dense matrix */
985     ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
986     col_lens[0] = MAT_FILE_CLASSID;
987     col_lens[1] = m;
988     col_lens[2] = n;
989     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
990     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
991     ierr = PetscFree(col_lens);CHKERRQ(ierr);
992 
993     /* write out matrix, by rows */
994     ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
995     v    = a->v;
996     for (j=0; j<n; j++) {
997       for (i=0; i<m; i++) {
998         vals[j + i*n] = *v++;
999       }
1000     }
1001     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
1002     ierr = PetscFree(vals);CHKERRQ(ierr);
1003   } else {
1004     ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
1005     col_lens[0] = MAT_FILE_CLASSID;
1006     col_lens[1] = m;
1007     col_lens[2] = n;
1008     col_lens[3] = nz;
1009 
1010     /* store lengths of each row and write (including header) to file */
1011     for (i=0; i<m; i++) col_lens[4+i] = n;
1012     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
1013 
1014     /* Possibly should write in smaller increments, not whole matrix at once? */
1015     /* store column indices (zero start index) */
1016     ict = 0;
1017     for (i=0; i<m; i++) {
1018       for (j=0; j<n; j++) col_lens[ict++] = j;
1019     }
1020     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
1021     ierr = PetscFree(col_lens);CHKERRQ(ierr);
1022 
1023     /* store nonzero values */
1024     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr);
1025     ict  = 0;
1026     for (i=0; i<m; i++) {
1027       v = a->v + i;
1028       for (j=0; j<n; j++) {
1029         anonz[ict++] = *v; v += a->lda;
1030       }
1031     }
1032     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
1033     ierr = PetscFree(anonz);CHKERRQ(ierr);
1034   }
1035   PetscFunctionReturn(0);
1036 }
1037 
1038 #undef __FUNCT__
1039 #define __FUNCT__ "MatView_SeqDense_Draw_Zoom"
1040 PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
1041 {
1042   Mat               A = (Mat) Aa;
1043   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
1044   PetscErrorCode    ierr;
1045   PetscInt          m = A->rmap->n,n = A->cmap->n,color,i,j;
1046   PetscScalar       *v = a->v;
1047   PetscViewer       viewer;
1048   PetscDraw         popup;
1049   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0;
1050   PetscViewerFormat format;
1051 
1052   PetscFunctionBegin;
1053   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
1054   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1055   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
1056 
1057   /* Loop over matrix elements drawing boxes */
1058   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
1059     /* Blue for negative and Red for positive */
1060     color = PETSC_DRAW_BLUE;
1061     for (j = 0; j < n; j++) {
1062       x_l = j;
1063       x_r = x_l + 1.0;
1064       for (i = 0; i < m; i++) {
1065         y_l = m - i - 1.0;
1066         y_r = y_l + 1.0;
1067         if (PetscRealPart(v[j*m+i]) >  0.) {
1068           color = PETSC_DRAW_RED;
1069         } else if (PetscRealPart(v[j*m+i]) <  0.) {
1070           color = PETSC_DRAW_BLUE;
1071         } else {
1072           continue;
1073         }
1074         ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
1075       }
1076     }
1077   } else {
1078     /* use contour shading to indicate magnitude of values */
1079     /* first determine max of all nonzero values */
1080     for (i = 0; i < m*n; i++) {
1081       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
1082     }
1083     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
1084     ierr  = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
1085     if (popup) {ierr  = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);}
1086     for (j = 0; j < n; j++) {
1087       x_l = j;
1088       x_r = x_l + 1.0;
1089       for (i = 0; i < m; i++) {
1090         y_l   = m - i - 1.0;
1091         y_r   = y_l + 1.0;
1092         color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i]));
1093         ierr  = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
1094       }
1095     }
1096   }
1097   PetscFunctionReturn(0);
1098 }
1099 
1100 #undef __FUNCT__
1101 #define __FUNCT__ "MatView_SeqDense_Draw"
1102 PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
1103 {
1104   PetscDraw      draw;
1105   PetscBool      isnull;
1106   PetscReal      xr,yr,xl,yl,h,w;
1107   PetscErrorCode ierr;
1108 
1109   PetscFunctionBegin;
1110   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
1111   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
1112   if (isnull) PetscFunctionReturn(0);
1113 
1114   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
1115   xr  = A->cmap->n; yr = A->rmap->n; h = yr/10.0; w = xr/10.0;
1116   xr += w;    yr += h;  xl = -w;     yl = -h;
1117   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
1118   ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr);
1119   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
1120   PetscFunctionReturn(0);
1121 }
1122 
1123 #undef __FUNCT__
1124 #define __FUNCT__ "MatView_SeqDense"
1125 PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer)
1126 {
1127   PetscErrorCode ierr;
1128   PetscBool      iascii,isbinary,isdraw;
1129 
1130   PetscFunctionBegin;
1131   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1132   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
1133   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
1134 
1135   if (iascii) {
1136     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
1137   } else if (isbinary) {
1138     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
1139   } else if (isdraw) {
1140     ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr);
1141   }
1142   PetscFunctionReturn(0);
1143 }
1144 
1145 #undef __FUNCT__
1146 #define __FUNCT__ "MatDestroy_SeqDense"
1147 PetscErrorCode MatDestroy_SeqDense(Mat mat)
1148 {
1149   Mat_SeqDense   *l = (Mat_SeqDense*)mat->data;
1150   PetscErrorCode ierr;
1151 
1152   PetscFunctionBegin;
1153 #if defined(PETSC_USE_LOG)
1154   PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap->n,mat->cmap->n);
1155 #endif
1156   ierr = PetscFree(l->pivots);CHKERRQ(ierr);
1157   if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);}
1158   ierr = PetscFree(mat->data);CHKERRQ(ierr);
1159 
1160   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
1161   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatDenseGetArray_C","",PETSC_NULL);CHKERRQ(ierr);
1162   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatDenseRestoreArray_C","",PETSC_NULL);CHKERRQ(ierr);
1163   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
1164   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
1165   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
1166   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
1167   PetscFunctionReturn(0);
1168 }
1169 
1170 #undef __FUNCT__
1171 #define __FUNCT__ "MatTranspose_SeqDense"
1172 PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout)
1173 {
1174   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1175   PetscErrorCode ierr;
1176   PetscInt       k,j,m,n,M;
1177   PetscScalar    *v,tmp;
1178 
1179   PetscFunctionBegin;
1180   v = mat->v; m = A->rmap->n; M = mat->lda; n = A->cmap->n;
1181   if (reuse == MAT_REUSE_MATRIX && *matout == A) { /* in place transpose */
1182     if (m != n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can not transpose non-square matrix in place");
1183     else {
1184       for (j=0; j<m; j++) {
1185         for (k=0; k<j; k++) {
1186           tmp = v[j + k*M];
1187           v[j + k*M] = v[k + j*M];
1188           v[k + j*M] = tmp;
1189         }
1190       }
1191     }
1192   } else { /* out-of-place transpose */
1193     Mat          tmat;
1194     Mat_SeqDense *tmatd;
1195     PetscScalar  *v2;
1196 
1197     if (reuse == MAT_INITIAL_MATRIX) {
1198       ierr  = MatCreate(((PetscObject)A)->comm,&tmat);CHKERRQ(ierr);
1199       ierr  = MatSetSizes(tmat,A->cmap->n,A->rmap->n,A->cmap->n,A->rmap->n);CHKERRQ(ierr);
1200       ierr  = MatSetType(tmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
1201       ierr  = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr);
1202     } else {
1203       tmat = *matout;
1204     }
1205     tmatd = (Mat_SeqDense*)tmat->data;
1206     v = mat->v; v2 = tmatd->v;
1207     for (j=0; j<n; j++) {
1208       for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M];
1209     }
1210     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1211     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1212     *matout = tmat;
1213   }
1214   PetscFunctionReturn(0);
1215 }
1216 
1217 #undef __FUNCT__
1218 #define __FUNCT__ "MatEqual_SeqDense"
1219 PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscBool  *flg)
1220 {
1221   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1222   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
1223   PetscInt     i,j;
1224   PetscScalar  *v1,*v2;
1225 
1226   PetscFunctionBegin;
1227   if (A1->rmap->n != A2->rmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1228   if (A1->cmap->n != A2->cmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1229   for (i=0; i<A1->rmap->n; i++) {
1230     v1 = mat1->v+i; v2 = mat2->v+i;
1231     for (j=0; j<A1->cmap->n; j++) {
1232       if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1233       v1 += mat1->lda; v2 += mat2->lda;
1234     }
1235   }
1236   *flg = PETSC_TRUE;
1237   PetscFunctionReturn(0);
1238 }
1239 
1240 #undef __FUNCT__
1241 #define __FUNCT__ "MatGetDiagonal_SeqDense"
1242 PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v)
1243 {
1244   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1245   PetscErrorCode ierr;
1246   PetscInt       i,n,len;
1247   PetscScalar    *x,zero = 0.0;
1248 
1249   PetscFunctionBegin;
1250   ierr = VecSet(v,zero);CHKERRQ(ierr);
1251   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
1252   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1253   len = PetscMin(A->rmap->n,A->cmap->n);
1254   if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
1255   for (i=0; i<len; i++) {
1256     x[i] = mat->v[i*mat->lda + i];
1257   }
1258   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1259   PetscFunctionReturn(0);
1260 }
1261 
1262 #undef __FUNCT__
1263 #define __FUNCT__ "MatDiagonalScale_SeqDense"
1264 PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1265 {
1266   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1267   PetscScalar    *l,*r,x,*v;
1268   PetscErrorCode ierr;
1269   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n;
1270 
1271   PetscFunctionBegin;
1272   if (ll) {
1273     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
1274     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1275     if (m != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1276     for (i=0; i<m; i++) {
1277       x = l[i];
1278       v = mat->v + i;
1279       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1280     }
1281     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1282     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1283   }
1284   if (rr) {
1285     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
1286     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1287     if (n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1288     for (i=0; i<n; i++) {
1289       x = r[i];
1290       v = mat->v + i*m;
1291       for (j=0; j<m; j++) { (*v++) *= x;}
1292     }
1293     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1294     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1295   }
1296   PetscFunctionReturn(0);
1297 }
1298 
1299 #undef __FUNCT__
1300 #define __FUNCT__ "MatNorm_SeqDense"
1301 PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1302 {
1303   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1304   PetscScalar  *v = mat->v;
1305   PetscReal    sum = 0.0;
1306   PetscInt     lda=mat->lda,m=A->rmap->n,i,j;
1307   PetscErrorCode ierr;
1308 
1309   PetscFunctionBegin;
1310   if (type == NORM_FROBENIUS) {
1311     if (lda>m) {
1312       for (j=0; j<A->cmap->n; j++) {
1313         v = mat->v+j*lda;
1314         for (i=0; i<m; i++) {
1315           sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1316         }
1317       }
1318     } else {
1319       for (i=0; i<A->cmap->n*A->rmap->n; i++) {
1320         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1321       }
1322     }
1323     *nrm = PetscSqrtReal(sum);
1324     ierr = PetscLogFlops(2.0*A->cmap->n*A->rmap->n);CHKERRQ(ierr);
1325   } else if (type == NORM_1) {
1326     *nrm = 0.0;
1327     for (j=0; j<A->cmap->n; j++) {
1328       v = mat->v + j*mat->lda;
1329       sum = 0.0;
1330       for (i=0; i<A->rmap->n; i++) {
1331         sum += PetscAbsScalar(*v);  v++;
1332       }
1333       if (sum > *nrm) *nrm = sum;
1334     }
1335     ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr);
1336   } else if (type == NORM_INFINITY) {
1337     *nrm = 0.0;
1338     for (j=0; j<A->rmap->n; j++) {
1339       v = mat->v + j;
1340       sum = 0.0;
1341       for (i=0; i<A->cmap->n; i++) {
1342         sum += PetscAbsScalar(*v); v += mat->lda;
1343       }
1344       if (sum > *nrm) *nrm = sum;
1345     }
1346     ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr);
1347   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No two norm");
1348   PetscFunctionReturn(0);
1349 }
1350 
1351 #undef __FUNCT__
1352 #define __FUNCT__ "MatSetOption_SeqDense"
1353 PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscBool  flg)
1354 {
1355   Mat_SeqDense   *aij = (Mat_SeqDense*)A->data;
1356   PetscErrorCode ierr;
1357 
1358   PetscFunctionBegin;
1359   switch (op) {
1360   case MAT_ROW_ORIENTED:
1361     aij->roworiented = flg;
1362     break;
1363   case MAT_NEW_NONZERO_LOCATIONS:
1364   case MAT_NEW_NONZERO_LOCATION_ERR:
1365   case MAT_NEW_NONZERO_ALLOCATION_ERR:
1366   case MAT_NEW_DIAGONALS:
1367   case MAT_KEEP_NONZERO_PATTERN:
1368   case MAT_IGNORE_OFF_PROC_ENTRIES:
1369   case MAT_USE_HASH_TABLE:
1370   case MAT_IGNORE_LOWER_TRIANGULAR:
1371     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
1372     break;
1373   case MAT_SPD:
1374   case MAT_SYMMETRIC:
1375   case MAT_STRUCTURALLY_SYMMETRIC:
1376   case MAT_HERMITIAN:
1377   case MAT_SYMMETRY_ETERNAL:
1378     /* These options are handled directly by MatSetOption() */
1379     break;
1380   default:
1381     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
1382   }
1383   PetscFunctionReturn(0);
1384 }
1385 
1386 #undef __FUNCT__
1387 #define __FUNCT__ "MatZeroEntries_SeqDense"
1388 PetscErrorCode MatZeroEntries_SeqDense(Mat A)
1389 {
1390   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
1391   PetscErrorCode ierr;
1392   PetscInt       lda=l->lda,m=A->rmap->n,j;
1393 
1394   PetscFunctionBegin;
1395   if (lda>m) {
1396     for (j=0; j<A->cmap->n; j++) {
1397       ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
1398     }
1399   } else {
1400     ierr = PetscMemzero(l->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
1401   }
1402   PetscFunctionReturn(0);
1403 }
1404 
1405 #undef __FUNCT__
1406 #define __FUNCT__ "MatZeroRows_SeqDense"
1407 PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
1408 {
1409   PetscErrorCode    ierr;
1410   Mat_SeqDense      *l = (Mat_SeqDense*)A->data;
1411   PetscInt          m = l->lda, n = A->cmap->n, i,j;
1412   PetscScalar       *slot,*bb;
1413   const PetscScalar *xx;
1414 
1415   PetscFunctionBegin;
1416 #if defined(PETSC_USE_DEBUG)
1417   for (i=0; i<N; i++) {
1418     if (rows[i] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row requested to be zeroed");
1419     if (rows[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %D requested to be zeroed greater than or equal number of rows %D",rows[i],A->rmap->n);
1420   }
1421 #endif
1422 
1423   /* fix right hand side if needed */
1424   if (x && b) {
1425     ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
1426     ierr = VecGetArray(b,&bb);CHKERRQ(ierr);
1427     for (i=0; i<N; i++) {
1428       bb[rows[i]] = diag*xx[rows[i]];
1429     }
1430     ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
1431     ierr = VecRestoreArray(b,&bb);CHKERRQ(ierr);
1432   }
1433 
1434   for (i=0; i<N; i++) {
1435     slot = l->v + rows[i];
1436     for (j=0; j<n; j++) { *slot = 0.0; slot += m;}
1437   }
1438   if (diag != 0.0) {
1439     if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only coded for square matrices");
1440     for (i=0; i<N; i++) {
1441       slot = l->v + (m+1)*rows[i];
1442       *slot = diag;
1443     }
1444   }
1445   PetscFunctionReturn(0);
1446 }
1447 
1448 #undef __FUNCT__
1449 #define __FUNCT__ "MatDenseGetArray_SeqDense"
1450 PetscErrorCode MatDenseGetArray_SeqDense(Mat A,PetscScalar *array[])
1451 {
1452   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1453 
1454   PetscFunctionBegin;
1455   if (mat->lda != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot get array for Dense matrices with LDA different from number of rows");
1456   *array = mat->v;
1457   PetscFunctionReturn(0);
1458 }
1459 
1460 #undef __FUNCT__
1461 #define __FUNCT__ "MatDenseRestoreArray_SeqDense"
1462 PetscErrorCode MatDenseRestoreArray_SeqDense(Mat A,PetscScalar *array[])
1463 {
1464   PetscFunctionBegin;
1465   *array = 0; /* user cannot accidently use the array later */
1466   PetscFunctionReturn(0);
1467 }
1468 
1469 #undef __FUNCT__
1470 #define __FUNCT__ "MatDenseGetArray"
1471 /*@C
1472    MatDenseGetArray - gives access to the array where the data for a SeqDense matrix is stored
1473 
1474    Not Collective
1475 
1476    Input Parameter:
1477 .  mat - a MATSEQDENSE matrix
1478 
1479    Output Parameter:
1480 .   array - pointer to the data
1481 
1482    Level: intermediate
1483 
1484 .seealso: MatDenseRestoreArray()
1485 @*/
1486 PetscErrorCode  MatDenseGetArray(Mat A,PetscScalar **array)
1487 {
1488   PetscErrorCode ierr;
1489 
1490   PetscFunctionBegin;
1491   ierr = PetscUseMethod(A,"MatDenseGetArray_C",(Mat,PetscScalar**),(A,array));CHKERRQ(ierr);
1492   PetscFunctionReturn(0);
1493 }
1494 
1495 #undef __FUNCT__
1496 #define __FUNCT__ "MatDenseRestoreArray"
1497 /*@C
1498    MatDenseRestoreArray - returns access to the array where the data for a SeqDense matrix is stored obtained by MatDenseGetArray()
1499 
1500    Not Collective
1501 
1502    Input Parameters:
1503 .  mat - a MATSEQDENSE matrix
1504 .  array - pointer to the data
1505 
1506    Level: intermediate
1507 
1508 .seealso: MatDenseGetArray()
1509 @*/
1510 PetscErrorCode  MatDenseRestoreArray(Mat A,PetscScalar **array)
1511 {
1512   PetscErrorCode ierr;
1513 
1514   PetscFunctionBegin;
1515   ierr = PetscUseMethod(A,"MatDenseRestoreArray_C",(Mat,PetscScalar**),(A,array));CHKERRQ(ierr);
1516   PetscFunctionReturn(0);
1517 }
1518 
1519 #undef __FUNCT__
1520 #define __FUNCT__ "MatGetSubMatrix_SeqDense"
1521 static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
1522 {
1523   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1524   PetscErrorCode ierr;
1525   PetscInt       i,j,nrows,ncols;
1526   const PetscInt *irow,*icol;
1527   PetscScalar    *av,*bv,*v = mat->v;
1528   Mat            newmat;
1529 
1530   PetscFunctionBegin;
1531   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
1532   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
1533   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
1534   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
1535 
1536   /* Check submatrixcall */
1537   if (scall == MAT_REUSE_MATRIX) {
1538     PetscInt n_cols,n_rows;
1539     ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr);
1540     if (n_rows != nrows || n_cols != ncols) {
1541       /* resize the result matrix to match number of requested rows/columns */
1542       ierr = MatSetSizes(*B,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
1543     }
1544     newmat = *B;
1545   } else {
1546     /* Create and fill new matrix */
1547     ierr = MatCreate(((PetscObject)A)->comm,&newmat);CHKERRQ(ierr);
1548     ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
1549     ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
1550     ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
1551   }
1552 
1553   /* Now extract the data pointers and do the copy,column at a time */
1554   bv = ((Mat_SeqDense*)newmat->data)->v;
1555 
1556   for (i=0; i<ncols; i++) {
1557     av = v + mat->lda*icol[i];
1558     for (j=0; j<nrows; j++) {
1559       *bv++ = av[irow[j]];
1560     }
1561   }
1562 
1563   /* Assemble the matrices so that the correct flags are set */
1564   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1565   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1566 
1567   /* Free work space */
1568   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
1569   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
1570   *B = newmat;
1571   PetscFunctionReturn(0);
1572 }
1573 
1574 #undef __FUNCT__
1575 #define __FUNCT__ "MatGetSubMatrices_SeqDense"
1576 PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1577 {
1578   PetscErrorCode ierr;
1579   PetscInt       i;
1580 
1581   PetscFunctionBegin;
1582   if (scall == MAT_INITIAL_MATRIX) {
1583     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
1584   }
1585 
1586   for (i=0; i<n; i++) {
1587     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1588   }
1589   PetscFunctionReturn(0);
1590 }
1591 
1592 #undef __FUNCT__
1593 #define __FUNCT__ "MatAssemblyBegin_SeqDense"
1594 PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode)
1595 {
1596   PetscFunctionBegin;
1597   PetscFunctionReturn(0);
1598 }
1599 
1600 #undef __FUNCT__
1601 #define __FUNCT__ "MatAssemblyEnd_SeqDense"
1602 PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode)
1603 {
1604   PetscFunctionBegin;
1605   PetscFunctionReturn(0);
1606 }
1607 
1608 #undef __FUNCT__
1609 #define __FUNCT__ "MatCopy_SeqDense"
1610 PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
1611 {
1612   Mat_SeqDense   *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
1613   PetscErrorCode ierr;
1614   PetscInt       lda1=a->lda,lda2=b->lda, m=A->rmap->n,n=A->cmap->n, j;
1615 
1616   PetscFunctionBegin;
1617   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
1618   if (A->ops->copy != B->ops->copy) {
1619     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
1620     PetscFunctionReturn(0);
1621   }
1622   if (m != B->rmap->n || n != B->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1623   if (lda1>m || lda2>m) {
1624     for (j=0; j<n; j++) {
1625       ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1626     }
1627   } else {
1628     ierr = PetscMemcpy(b->v,a->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
1629   }
1630   PetscFunctionReturn(0);
1631 }
1632 
1633 #undef __FUNCT__
1634 #define __FUNCT__ "MatSetUp_SeqDense"
1635 PetscErrorCode MatSetUp_SeqDense(Mat A)
1636 {
1637   PetscErrorCode ierr;
1638 
1639   PetscFunctionBegin;
1640   ierr =  MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr);
1641   PetscFunctionReturn(0);
1642 }
1643 
1644 #undef __FUNCT__
1645 #define __FUNCT__ "MatConjugate_SeqDense"
1646 static PetscErrorCode MatConjugate_SeqDense(Mat A)
1647 {
1648   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1649   PetscInt       i,nz = A->rmap->n*A->cmap->n;
1650   PetscScalar    *aa = a->v;
1651 
1652   PetscFunctionBegin;
1653   for (i=0; i<nz; i++) aa[i] = PetscConj(aa[i]);
1654   PetscFunctionReturn(0);
1655 }
1656 
1657 #undef __FUNCT__
1658 #define __FUNCT__ "MatRealPart_SeqDense"
1659 static PetscErrorCode MatRealPart_SeqDense(Mat A)
1660 {
1661   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1662   PetscInt       i,nz = A->rmap->n*A->cmap->n;
1663   PetscScalar    *aa = a->v;
1664 
1665   PetscFunctionBegin;
1666   for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
1667   PetscFunctionReturn(0);
1668 }
1669 
1670 #undef __FUNCT__
1671 #define __FUNCT__ "MatImaginaryPart_SeqDense"
1672 static PetscErrorCode MatImaginaryPart_SeqDense(Mat A)
1673 {
1674   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1675   PetscInt       i,nz = A->rmap->n*A->cmap->n;
1676   PetscScalar    *aa = a->v;
1677 
1678   PetscFunctionBegin;
1679   for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
1680   PetscFunctionReturn(0);
1681 }
1682 
1683 /* ----------------------------------------------------------------*/
1684 #undef __FUNCT__
1685 #define __FUNCT__ "MatMatMult_SeqDense_SeqDense"
1686 PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1687 {
1688   PetscErrorCode ierr;
1689 
1690   PetscFunctionBegin;
1691   if (scall == MAT_INITIAL_MATRIX) {
1692     ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1693   }
1694   ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1695   PetscFunctionReturn(0);
1696 }
1697 
1698 #undef __FUNCT__
1699 #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense"
1700 PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1701 {
1702   PetscErrorCode ierr;
1703   PetscInt       m=A->rmap->n,n=B->cmap->n;
1704   Mat            Cmat;
1705 
1706   PetscFunctionBegin;
1707   if (A->cmap->n != B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"A->cmap->n %d != B->rmap->n %d\n",A->cmap->n,B->rmap->n);
1708   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1709   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1710   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1711   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1712 
1713   *C = Cmat;
1714   PetscFunctionReturn(0);
1715 }
1716 
1717 #undef __FUNCT__
1718 #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense"
1719 PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1720 {
1721   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1722   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1723   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
1724   PetscBLASInt   m,n,k;
1725   PetscScalar    _DOne=1.0,_DZero=0.0;
1726   PetscErrorCode ierr;
1727 
1728   PetscFunctionBegin;
1729   ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr);
1730   ierr = PetscBLASIntCast(B->cmap->n,&n);CHKERRQ(ierr);
1731   ierr = PetscBLASIntCast(A->cmap->n,&k);CHKERRQ(ierr);
1732   BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1733   PetscFunctionReturn(0);
1734 }
1735 
1736 #undef __FUNCT__
1737 #define __FUNCT__ "MatTransposeMatMult_SeqDense_SeqDense"
1738 PetscErrorCode MatTransposeMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1739 {
1740   PetscErrorCode ierr;
1741 
1742   PetscFunctionBegin;
1743   if (scall == MAT_INITIAL_MATRIX) {
1744     ierr = MatTransposeMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1745   }
1746   ierr = MatTransposeMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1747   PetscFunctionReturn(0);
1748 }
1749 
1750 #undef __FUNCT__
1751 #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqDense_SeqDense"
1752 PetscErrorCode MatTransposeMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1753 {
1754   PetscErrorCode ierr;
1755   PetscInt       m=A->cmap->n,n=B->cmap->n;
1756   Mat            Cmat;
1757 
1758   PetscFunctionBegin;
1759   if (A->rmap->n != B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"A->rmap->n %d != B->rmap->n %d\n",A->rmap->n,B->rmap->n);
1760   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1761   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1762   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1763   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1764   Cmat->assembled = PETSC_TRUE;
1765   *C = Cmat;
1766   PetscFunctionReturn(0);
1767 }
1768 
1769 #undef __FUNCT__
1770 #define __FUNCT__ "MatTransposeMatMultNumeric_SeqDense_SeqDense"
1771 PetscErrorCode MatTransposeMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1772 {
1773   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1774   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1775   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
1776   PetscBLASInt   m,n,k;
1777   PetscScalar    _DOne=1.0,_DZero=0.0;
1778   PetscErrorCode ierr;
1779 
1780   PetscFunctionBegin;
1781   ierr = PetscBLASIntCast(A->cmap->n,&m);CHKERRQ(ierr);
1782   ierr = PetscBLASIntCast(B->cmap->n,&n);CHKERRQ(ierr);
1783   ierr = PetscBLASIntCast(A->rmap->n,&k);CHKERRQ(ierr);
1784   /*
1785      Note the m and n arguments below are the number rows and columns of A', not A!
1786   */
1787   BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1788   PetscFunctionReturn(0);
1789 }
1790 
1791 #undef __FUNCT__
1792 #define __FUNCT__ "MatGetRowMax_SeqDense"
1793 PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[])
1794 {
1795   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1796   PetscErrorCode ierr;
1797   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1798   PetscScalar    *x;
1799   MatScalar      *aa = a->v;
1800 
1801   PetscFunctionBegin;
1802   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1803 
1804   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1805   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1806   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1807   if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1808   for (i=0; i<m; i++) {
1809     x[i] = aa[i]; if (idx) idx[i] = 0;
1810     for (j=1; j<n; j++) {
1811       if (PetscRealPart(x[i]) < PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
1812     }
1813   }
1814   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1815   PetscFunctionReturn(0);
1816 }
1817 
1818 #undef __FUNCT__
1819 #define __FUNCT__ "MatGetRowMaxAbs_SeqDense"
1820 PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[])
1821 {
1822   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1823   PetscErrorCode ierr;
1824   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1825   PetscScalar    *x;
1826   PetscReal      atmp;
1827   MatScalar      *aa = a->v;
1828 
1829   PetscFunctionBegin;
1830   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1831 
1832   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1833   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1834   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1835   if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1836   for (i=0; i<m; i++) {
1837     x[i] = PetscAbsScalar(aa[i]);
1838     for (j=1; j<n; j++) {
1839       atmp = PetscAbsScalar(aa[i+m*j]);
1840       if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;}
1841     }
1842   }
1843   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1844   PetscFunctionReturn(0);
1845 }
1846 
1847 #undef __FUNCT__
1848 #define __FUNCT__ "MatGetRowMin_SeqDense"
1849 PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[])
1850 {
1851   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1852   PetscErrorCode ierr;
1853   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1854   PetscScalar    *x;
1855   MatScalar      *aa = a->v;
1856 
1857   PetscFunctionBegin;
1858   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1859 
1860   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1861   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1862   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1863   if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1864   for (i=0; i<m; i++) {
1865     x[i] = aa[i]; if (idx) idx[i] = 0;
1866     for (j=1; j<n; j++) {
1867       if (PetscRealPart(x[i]) > PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
1868     }
1869   }
1870   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1871   PetscFunctionReturn(0);
1872 }
1873 
1874 #undef __FUNCT__
1875 #define __FUNCT__ "MatGetColumnVector_SeqDense"
1876 PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col)
1877 {
1878   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1879   PetscErrorCode ierr;
1880   PetscScalar    *x;
1881 
1882   PetscFunctionBegin;
1883   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1884 
1885   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1886   ierr = PetscMemcpy(x,a->v+col*a->lda,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
1887   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1888   PetscFunctionReturn(0);
1889 }
1890 
1891 
1892 #undef __FUNCT__
1893 #define __FUNCT__ "MatGetColumnNorms_SeqDense"
1894 PetscErrorCode MatGetColumnNorms_SeqDense(Mat A,NormType type,PetscReal *norms)
1895 {
1896   PetscErrorCode ierr;
1897   PetscInt       i,j,m,n;
1898   PetscScalar    *a;
1899 
1900   PetscFunctionBegin;
1901   ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr);
1902   ierr = PetscMemzero(norms,n*sizeof(PetscReal));CHKERRQ(ierr);
1903   ierr = MatDenseGetArray(A,&a);CHKERRQ(ierr);
1904   if (type == NORM_2) {
1905     for (i=0; i<n; i++) {
1906       for (j=0; j<m; j++) {
1907         norms[i] += PetscAbsScalar(a[j]*a[j]);
1908       }
1909       a += m;
1910     }
1911   } else if (type == NORM_1) {
1912     for (i=0; i<n; i++) {
1913       for (j=0; j<m; j++) {
1914         norms[i] += PetscAbsScalar(a[j]);
1915       }
1916       a += m;
1917     }
1918   } else if (type == NORM_INFINITY) {
1919     for (i=0; i<n; i++) {
1920       for (j=0; j<m; j++) {
1921         norms[i] = PetscMax(PetscAbsScalar(a[j]),norms[i]);
1922       }
1923       a += m;
1924     }
1925   } else SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Unknown NormType");
1926   ierr = MatDenseRestoreArray(A,&a);CHKERRQ(ierr);
1927   if (type == NORM_2) {
1928     for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]);
1929   }
1930   PetscFunctionReturn(0);
1931 }
1932 
1933 #undef __FUNCT__
1934 #define __FUNCT__ "MatSetRandom_SeqDense"
1935 static PetscErrorCode  MatSetRandom_SeqDense(Mat x,PetscRandom rctx)
1936 {
1937   PetscErrorCode ierr;
1938   PetscScalar    *a;
1939   PetscInt       m,n,i;
1940 
1941   PetscFunctionBegin;
1942   ierr = MatGetSize(x,&m,&n);CHKERRQ(ierr);
1943   ierr = MatDenseGetArray(x,&a);CHKERRQ(ierr);
1944   for (i=0; i<m*n; i++) {
1945     ierr = PetscRandomGetValue(rctx,a+i);CHKERRQ(ierr);
1946   }
1947   ierr = MatDenseRestoreArray(x,&a);CHKERRQ(ierr);
1948   PetscFunctionReturn(0);
1949 }
1950 
1951 
1952 /* -------------------------------------------------------------------*/
1953 static struct _MatOps MatOps_Values = { MatSetValues_SeqDense,
1954                                         MatGetRow_SeqDense,
1955                                         MatRestoreRow_SeqDense,
1956                                         MatMult_SeqDense,
1957                                 /*  4*/ MatMultAdd_SeqDense,
1958                                         MatMultTranspose_SeqDense,
1959                                         MatMultTransposeAdd_SeqDense,
1960                                         0,
1961                                         0,
1962                                         0,
1963                                 /* 10*/ 0,
1964                                         MatLUFactor_SeqDense,
1965                                         MatCholeskyFactor_SeqDense,
1966                                         MatSOR_SeqDense,
1967                                         MatTranspose_SeqDense,
1968                                 /* 15*/ MatGetInfo_SeqDense,
1969                                         MatEqual_SeqDense,
1970                                         MatGetDiagonal_SeqDense,
1971                                         MatDiagonalScale_SeqDense,
1972                                         MatNorm_SeqDense,
1973                                 /* 20*/ MatAssemblyBegin_SeqDense,
1974                                         MatAssemblyEnd_SeqDense,
1975                                         MatSetOption_SeqDense,
1976                                         MatZeroEntries_SeqDense,
1977                                 /* 24*/ MatZeroRows_SeqDense,
1978                                         0,
1979                                         0,
1980                                         0,
1981                                         0,
1982                                 /* 29*/ MatSetUp_SeqDense,
1983                                         0,
1984                                         0,
1985                                         0,
1986                                         0,
1987                                 /* 34*/ MatDuplicate_SeqDense,
1988                                         0,
1989                                         0,
1990                                         0,
1991                                         0,
1992                                 /* 39*/ MatAXPY_SeqDense,
1993                                         MatGetSubMatrices_SeqDense,
1994                                         0,
1995                                         MatGetValues_SeqDense,
1996                                         MatCopy_SeqDense,
1997                                 /* 44*/ MatGetRowMax_SeqDense,
1998                                         MatScale_SeqDense,
1999                                         0,
2000                                         0,
2001                                         0,
2002                                 /* 49*/ MatSetRandom_SeqDense,
2003                                         0,
2004                                         0,
2005                                         0,
2006                                         0,
2007                                 /* 54*/ 0,
2008                                         0,
2009                                         0,
2010                                         0,
2011                                         0,
2012                                 /* 59*/ 0,
2013                                         MatDestroy_SeqDense,
2014                                         MatView_SeqDense,
2015                                         0,
2016                                         0,
2017                                 /* 64*/ 0,
2018                                         0,
2019                                         0,
2020                                         0,
2021                                         0,
2022                                 /* 69*/ MatGetRowMaxAbs_SeqDense,
2023                                         0,
2024                                         0,
2025                                         0,
2026                                         0,
2027                                 /* 74*/ 0,
2028                                         0,
2029                                         0,
2030                                         0,
2031                                         0,
2032                                 /* 79*/ 0,
2033                                         0,
2034                                         0,
2035                                         0,
2036                                 /* 83*/ MatLoad_SeqDense,
2037                                         0,
2038                                         MatIsHermitian_SeqDense,
2039                                         0,
2040                                         0,
2041                                         0,
2042                                 /* 89*/ MatMatMult_SeqDense_SeqDense,
2043                                         MatMatMultSymbolic_SeqDense_SeqDense,
2044                                         MatMatMultNumeric_SeqDense_SeqDense,
2045                                         0,
2046                                         0,
2047                                 /* 94*/ 0,
2048                                         0,
2049                                         0,
2050                                         0,
2051                                         0,
2052                                 /* 99*/ 0,
2053                                         0,
2054                                         0,
2055                                         MatConjugate_SeqDense,
2056                                         0,
2057                                 /*104*/ 0,
2058                                         MatRealPart_SeqDense,
2059                                         MatImaginaryPart_SeqDense,
2060                                         0,
2061                                         0,
2062                                 /*109*/ MatMatSolve_SeqDense,
2063                                         0,
2064                                         MatGetRowMin_SeqDense,
2065                                         MatGetColumnVector_SeqDense,
2066                                         0,
2067                                 /*114*/ 0,
2068                                         0,
2069                                         0,
2070                                         0,
2071                                         0,
2072                                 /*119*/ 0,
2073                                         0,
2074                                         0,
2075                                         0,
2076                                         0,
2077                                 /*124*/ 0,
2078                                         MatGetColumnNorms_SeqDense,
2079                                         0,
2080                                         0,
2081                                         0,
2082                                 /*129*/ 0,
2083                                         MatTransposeMatMult_SeqDense_SeqDense,
2084                                         MatTransposeMatMultSymbolic_SeqDense_SeqDense,
2085                                         MatTransposeMatMultNumeric_SeqDense_SeqDense,
2086 };
2087 
2088 #undef __FUNCT__
2089 #define __FUNCT__ "MatCreateSeqDense"
2090 /*@C
2091    MatCreateSeqDense - Creates a sequential dense matrix that
2092    is stored in column major order (the usual Fortran 77 manner). Many
2093    of the matrix operations use the BLAS and LAPACK routines.
2094 
2095    Collective on MPI_Comm
2096 
2097    Input Parameters:
2098 +  comm - MPI communicator, set to PETSC_COMM_SELF
2099 .  m - number of rows
2100 .  n - number of columns
2101 -  data - optional location of matrix data in column major order.  Set data=PETSC_NULL for PETSc
2102    to control all matrix memory allocation.
2103 
2104    Output Parameter:
2105 .  A - the matrix
2106 
2107    Notes:
2108    The data input variable is intended primarily for Fortran programmers
2109    who wish to allocate their own matrix memory space.  Most users should
2110    set data=PETSC_NULL.
2111 
2112    Level: intermediate
2113 
2114 .keywords: dense, matrix, LAPACK, BLAS
2115 
2116 .seealso: MatCreate(), MatCreateDense(), MatSetValues()
2117 @*/
2118 PetscErrorCode  MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A)
2119 {
2120   PetscErrorCode ierr;
2121 
2122   PetscFunctionBegin;
2123   ierr = MatCreate(comm,A);CHKERRQ(ierr);
2124   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
2125   ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
2126   ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
2127   PetscFunctionReturn(0);
2128 }
2129 
2130 #undef __FUNCT__
2131 #define __FUNCT__ "MatSeqDenseSetPreallocation"
2132 /*@C
2133    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
2134 
2135    Collective on MPI_Comm
2136 
2137    Input Parameters:
2138 +  A - the matrix
2139 -  data - the array (or PETSC_NULL)
2140 
2141    Notes:
2142    The data input variable is intended primarily for Fortran programmers
2143    who wish to allocate their own matrix memory space.  Most users should
2144    need not call this routine.
2145 
2146    Level: intermediate
2147 
2148 .keywords: dense, matrix, LAPACK, BLAS
2149 
2150 .seealso: MatCreate(), MatCreateDense(), MatSetValues(), MatSeqDenseSetLDA()
2151 
2152 @*/
2153 PetscErrorCode  MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
2154 {
2155   PetscErrorCode ierr;
2156 
2157   PetscFunctionBegin;
2158   ierr = PetscTryMethod(B,"MatSeqDenseSetPreallocation_C",(Mat,PetscScalar[]),(B,data));CHKERRQ(ierr);
2159   PetscFunctionReturn(0);
2160 }
2161 
2162 EXTERN_C_BEGIN
2163 #undef __FUNCT__
2164 #define __FUNCT__ "MatSeqDenseSetPreallocation_SeqDense"
2165 PetscErrorCode  MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
2166 {
2167   Mat_SeqDense   *b;
2168   PetscErrorCode ierr;
2169 
2170   PetscFunctionBegin;
2171   B->preallocated = PETSC_TRUE;
2172 
2173   ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr);
2174   ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr);
2175 
2176   b       = (Mat_SeqDense*)B->data;
2177   b->Mmax = B->rmap->n;
2178   b->Nmax = B->cmap->n;
2179   if (b->lda <= 0 || b->changelda) b->lda = B->rmap->n;
2180 
2181   if (!data) { /* petsc-allocated storage */
2182     if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); }
2183     ierr = PetscMalloc(b->lda*b->Nmax*sizeof(PetscScalar),&b->v);CHKERRQ(ierr);
2184     ierr = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
2185     ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
2186     b->user_alloc = PETSC_FALSE;
2187   } else { /* user-allocated storage */
2188     if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); }
2189     b->v          = data;
2190     b->user_alloc = PETSC_TRUE;
2191   }
2192   B->assembled = PETSC_TRUE;
2193   PetscFunctionReturn(0);
2194 }
2195 EXTERN_C_END
2196 
2197 #undef __FUNCT__
2198 #define __FUNCT__ "MatSeqDenseSetLDA"
2199 /*@C
2200   MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array
2201 
2202   Input parameter:
2203 + A - the matrix
2204 - lda - the leading dimension
2205 
2206   Notes:
2207   This routine is to be used in conjunction with MatSeqDenseSetPreallocation();
2208   it asserts that the preallocation has a leading dimension (the LDA parameter
2209   of Blas and Lapack fame) larger than M, the first dimension of the matrix.
2210 
2211   Level: intermediate
2212 
2213 .keywords: dense, matrix, LAPACK, BLAS
2214 
2215 .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize()
2216 
2217 @*/
2218 PetscErrorCode  MatSeqDenseSetLDA(Mat B,PetscInt lda)
2219 {
2220   Mat_SeqDense *b = (Mat_SeqDense*)B->data;
2221 
2222   PetscFunctionBegin;
2223   if (lda < B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->rmap->n);
2224   b->lda       = lda;
2225   b->changelda = PETSC_FALSE;
2226   b->Mmax      = PetscMax(b->Mmax,lda);
2227   PetscFunctionReturn(0);
2228 }
2229 
2230 /*MC
2231    MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices.
2232 
2233    Options Database Keys:
2234 . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions()
2235 
2236   Level: beginner
2237 
2238 .seealso: MatCreateSeqDense()
2239 
2240 M*/
2241 
2242 EXTERN_C_BEGIN
2243 #undef __FUNCT__
2244 #define __FUNCT__ "MatCreate_SeqDense"
2245 PetscErrorCode  MatCreate_SeqDense(Mat B)
2246 {
2247   Mat_SeqDense   *b;
2248   PetscErrorCode ierr;
2249   PetscMPIInt    size;
2250 
2251   PetscFunctionBegin;
2252   ierr = MPI_Comm_size(((PetscObject)B)->comm,&size);CHKERRQ(ierr);
2253   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
2254 
2255   ierr            = PetscNewLog(B,Mat_SeqDense,&b);CHKERRQ(ierr);
2256   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
2257   B->data         = (void*)b;
2258 
2259   b->pivots       = 0;
2260   b->roworiented  = PETSC_TRUE;
2261   b->v            = 0;
2262   b->changelda    = PETSC_FALSE;
2263 
2264   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatDenseGetArray_C","MatDenseGetArray_SeqDense",MatDenseGetArray_SeqDense);CHKERRQ(ierr);
2265   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatDenseRestoreArray_C","MatDenseRestoreArray_SeqDense",MatDenseRestoreArray_SeqDense);CHKERRQ(ierr);
2266 
2267   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqdense_seqaij_C","MatConvert_SeqDense_SeqAIJ",MatConvert_SeqDense_SeqAIJ);CHKERRQ(ierr);
2268   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_petsc_C",
2269                                            "MatGetFactor_seqdense_petsc",
2270                                             MatGetFactor_seqdense_petsc);CHKERRQ(ierr);
2271   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C",
2272                                            "MatSeqDenseSetPreallocation_SeqDense",
2273                                             MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr);
2274 
2275   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMult_seqaij_seqdense_C",
2276                                            "MatMatMult_SeqAIJ_SeqDense",
2277                                             MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr);
2278 
2279 
2280   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C",
2281                                            "MatMatMultSymbolic_SeqAIJ_SeqDense",
2282                                             MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr);
2283   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C",
2284                                            "MatMatMultNumeric_SeqAIJ_SeqDense",
2285                                             MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr);
2286   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);CHKERRQ(ierr);
2287   PetscFunctionReturn(0);
2288 }
2289 EXTERN_C_END
2290