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