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