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