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