xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 28adb3f739167aae2a3fdf9bf501dbd0150de1de)
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->factor == 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->factor == 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->factor = 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->factor = 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)->factor = ftype;
436   PetscFunctionReturn(0);
437 }
438 EXTERN_C_END
439 
440 /* ------------------------------------------------------------------*/
441 #undef __FUNCT__
442 #define __FUNCT__ "MatRelax_SeqDense"
443 PetscErrorCode MatRelax_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 (!mat->roworiented) {
640     if (addv == INSERT_VALUES) {
641       for (j=0; j<n; j++) {
642         if (indexn[j] < 0) {idx += m; continue;}
643 #if defined(PETSC_USE_DEBUG)
644         if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
645 #endif
646         for (i=0; i<m; i++) {
647           if (indexm[i] < 0) {idx++; continue;}
648 #if defined(PETSC_USE_DEBUG)
649           if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
650 #endif
651           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
652         }
653       }
654     } else {
655       for (j=0; j<n; j++) {
656         if (indexn[j] < 0) {idx += m; continue;}
657 #if defined(PETSC_USE_DEBUG)
658         if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
659 #endif
660         for (i=0; i<m; i++) {
661           if (indexm[i] < 0) {idx++; continue;}
662 #if defined(PETSC_USE_DEBUG)
663           if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
664 #endif
665           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
666         }
667       }
668     }
669   } else {
670     if (addv == INSERT_VALUES) {
671       for (i=0; i<m; i++) {
672         if (indexm[i] < 0) { idx += n; continue;}
673 #if defined(PETSC_USE_DEBUG)
674         if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
675 #endif
676         for (j=0; j<n; j++) {
677           if (indexn[j] < 0) { idx++; continue;}
678 #if defined(PETSC_USE_DEBUG)
679           if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
680 #endif
681           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
682         }
683       }
684     } else {
685       for (i=0; i<m; i++) {
686         if (indexm[i] < 0) { idx += n; continue;}
687 #if defined(PETSC_USE_DEBUG)
688         if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
689 #endif
690         for (j=0; j<n; j++) {
691           if (indexn[j] < 0) { idx++; continue;}
692 #if defined(PETSC_USE_DEBUG)
693           if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
694 #endif
695           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
696         }
697       }
698     }
699   }
700   PetscFunctionReturn(0);
701 }
702 
703 #undef __FUNCT__
704 #define __FUNCT__ "MatGetValues_SeqDense"
705 PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[])
706 {
707   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
708   PetscInt     i,j;
709 
710   PetscFunctionBegin;
711   /* row-oriented output */
712   for (i=0; i<m; i++) {
713     if (indexm[i] < 0) {v += n;continue;}
714     if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row %D requested larger than number rows %D",indexm[i],A->rmap->n);
715     for (j=0; j<n; j++) {
716       if (indexn[j] < 0) {v++; continue;}
717       if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column %D requested larger than number columns %D",indexn[j],A->cmap->n);
718       *v++ = mat->v[indexn[j]*mat->lda + indexm[i]];
719     }
720   }
721   PetscFunctionReturn(0);
722 }
723 
724 /* -----------------------------------------------------------------*/
725 
726 #include "petscsys.h"
727 
728 #undef __FUNCT__
729 #define __FUNCT__ "MatLoad_SeqDense"
730 PetscErrorCode MatLoad_SeqDense(PetscViewer viewer, const MatType type,Mat *A)
731 {
732   Mat_SeqDense   *a;
733   Mat            B;
734   PetscErrorCode ierr;
735   PetscInt       *scols,i,j,nz,header[4];
736   int            fd;
737   PetscMPIInt    size;
738   PetscInt       *rowlengths = 0,M,N,*cols;
739   PetscScalar    *vals,*svals,*v,*w;
740   MPI_Comm       comm = ((PetscObject)viewer)->comm;
741 
742   PetscFunctionBegin;
743   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
744   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
745   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
746   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
747   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
748   M = header[1]; N = header[2]; nz = header[3];
749 
750   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
751     ierr = MatCreate(comm,A);CHKERRQ(ierr);
752     ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr);
753     ierr = MatSetType(*A,type);CHKERRQ(ierr);
754     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
755     B    = *A;
756     a    = (Mat_SeqDense*)B->data;
757     v    = a->v;
758     /* Allocate some temp space to read in the values and then flip them
759        from row major to column major */
760     ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr);
761     /* read in nonzero values */
762     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr);
763     /* now flip the values and store them in the matrix*/
764     for (j=0; j<N; j++) {
765       for (i=0; i<M; i++) {
766         *v++ =w[i*N+j];
767       }
768     }
769     ierr = PetscFree(w);CHKERRQ(ierr);
770     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
771     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
772   } else {
773     /* read row lengths */
774     ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
775     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
776 
777     /* create our matrix */
778     ierr = MatCreate(comm,A);CHKERRQ(ierr);
779     ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr);
780     ierr = MatSetType(*A,type);CHKERRQ(ierr);
781     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
782     B = *A;
783     a = (Mat_SeqDense*)B->data;
784     v = a->v;
785 
786     /* read column indices and nonzeros */
787     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr);
788     cols = scols;
789     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
790     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr);
791     vals = svals;
792     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
793 
794     /* insert into matrix */
795     for (i=0; i<M; i++) {
796       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
797       svals += rowlengths[i]; scols += rowlengths[i];
798     }
799     ierr = PetscFree(vals);CHKERRQ(ierr);
800     ierr = PetscFree(cols);CHKERRQ(ierr);
801     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
802 
803     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
804     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
805   }
806   PetscFunctionReturn(0);
807 }
808 
809 #include "petscsys.h"
810 
811 #undef __FUNCT__
812 #define __FUNCT__ "MatView_SeqDense_ASCII"
813 static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
814 {
815   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
816   PetscErrorCode    ierr;
817   PetscInt          i,j;
818   const char        *name;
819   PetscScalar       *v;
820   PetscViewerFormat format;
821 #if defined(PETSC_USE_COMPLEX)
822   PetscTruth        allreal = PETSC_TRUE;
823 #endif
824 
825   PetscFunctionBegin;
826   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
827   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
828   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
829     PetscFunctionReturn(0);  /* do nothing for now */
830   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
831     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
832     for (i=0; i<A->rmap->n; i++) {
833       v = a->v + i;
834       ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr);
835       for (j=0; j<A->cmap->n; j++) {
836 #if defined(PETSC_USE_COMPLEX)
837         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
838           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
839         } else if (PetscRealPart(*v)) {
840           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,PetscRealPart(*v));CHKERRQ(ierr);
841         }
842 #else
843         if (*v) {
844           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,*v);CHKERRQ(ierr);
845         }
846 #endif
847         v += a->lda;
848       }
849       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
850     }
851     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
852   } else {
853     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
854 #if defined(PETSC_USE_COMPLEX)
855     /* determine if matrix has all real values */
856     v = a->v;
857     for (i=0; i<A->rmap->n*A->cmap->n; i++) {
858 	if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;}
859     }
860 #endif
861     if (format == PETSC_VIEWER_ASCII_MATLAB) {
862       ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
863       ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap->n,A->cmap->n);CHKERRQ(ierr);
864       ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap->n,A->cmap->n);CHKERRQ(ierr);
865       ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr);
866     }
867 
868     for (i=0; i<A->rmap->n; i++) {
869       v = a->v + i;
870       for (j=0; j<A->cmap->n; j++) {
871 #if defined(PETSC_USE_COMPLEX)
872         if (allreal) {
873           ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",PetscRealPart(*v));CHKERRQ(ierr);
874         } else {
875           ierr = PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
876         }
877 #else
878         ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",*v);CHKERRQ(ierr);
879 #endif
880 	v += a->lda;
881       }
882       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
883     }
884     if (format == PETSC_VIEWER_ASCII_MATLAB) {
885       ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr);
886     }
887     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
888   }
889   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
890   PetscFunctionReturn(0);
891 }
892 
893 #undef __FUNCT__
894 #define __FUNCT__ "MatView_SeqDense_Binary"
895 static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
896 {
897   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
898   PetscErrorCode    ierr;
899   int               fd;
900   PetscInt          ict,j,n = A->cmap->n,m = A->rmap->n,i,*col_lens,nz = m*n;
901   PetscScalar       *v,*anonz,*vals;
902   PetscViewerFormat format;
903 
904   PetscFunctionBegin;
905   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
906 
907   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
908   if (format == PETSC_VIEWER_NATIVE) {
909     /* store the matrix as a dense matrix */
910     ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
911     col_lens[0] = MAT_FILE_COOKIE;
912     col_lens[1] = m;
913     col_lens[2] = n;
914     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
915     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
916     ierr = PetscFree(col_lens);CHKERRQ(ierr);
917 
918     /* write out matrix, by rows */
919     ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
920     v    = a->v;
921     for (j=0; j<n; j++) {
922       for (i=0; i<m; i++) {
923         vals[j + i*n] = *v++;
924       }
925     }
926     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
927     ierr = PetscFree(vals);CHKERRQ(ierr);
928   } else {
929     ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
930     col_lens[0] = MAT_FILE_COOKIE;
931     col_lens[1] = m;
932     col_lens[2] = n;
933     col_lens[3] = nz;
934 
935     /* store lengths of each row and write (including header) to file */
936     for (i=0; i<m; i++) col_lens[4+i] = n;
937     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
938 
939     /* Possibly should write in smaller increments, not whole matrix at once? */
940     /* store column indices (zero start index) */
941     ict = 0;
942     for (i=0; i<m; i++) {
943       for (j=0; j<n; j++) col_lens[ict++] = j;
944     }
945     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
946     ierr = PetscFree(col_lens);CHKERRQ(ierr);
947 
948     /* store nonzero values */
949     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr);
950     ict  = 0;
951     for (i=0; i<m; i++) {
952       v = a->v + i;
953       for (j=0; j<n; j++) {
954         anonz[ict++] = *v; v += a->lda;
955       }
956     }
957     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
958     ierr = PetscFree(anonz);CHKERRQ(ierr);
959   }
960   PetscFunctionReturn(0);
961 }
962 
963 #undef __FUNCT__
964 #define __FUNCT__ "MatView_SeqDense_Draw_Zoom"
965 PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
966 {
967   Mat               A = (Mat) Aa;
968   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
969   PetscErrorCode    ierr;
970   PetscInt          m = A->rmap->n,n = A->cmap->n,color,i,j;
971   PetscScalar       *v = a->v;
972   PetscViewer       viewer;
973   PetscDraw         popup;
974   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0;
975   PetscViewerFormat format;
976 
977   PetscFunctionBegin;
978 
979   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
980   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
981   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
982 
983   /* Loop over matrix elements drawing boxes */
984   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
985     /* Blue for negative and Red for positive */
986     color = PETSC_DRAW_BLUE;
987     for(j = 0; j < n; j++) {
988       x_l = j;
989       x_r = x_l + 1.0;
990       for(i = 0; i < m; i++) {
991         y_l = m - i - 1.0;
992         y_r = y_l + 1.0;
993 #if defined(PETSC_USE_COMPLEX)
994         if (PetscRealPart(v[j*m+i]) >  0.) {
995           color = PETSC_DRAW_RED;
996         } else if (PetscRealPart(v[j*m+i]) <  0.) {
997           color = PETSC_DRAW_BLUE;
998         } else {
999           continue;
1000         }
1001 #else
1002         if (v[j*m+i] >  0.) {
1003           color = PETSC_DRAW_RED;
1004         } else if (v[j*m+i] <  0.) {
1005           color = PETSC_DRAW_BLUE;
1006         } else {
1007           continue;
1008         }
1009 #endif
1010         ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
1011       }
1012     }
1013   } else {
1014     /* use contour shading to indicate magnitude of values */
1015     /* first determine max of all nonzero values */
1016     for(i = 0; i < m*n; i++) {
1017       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
1018     }
1019     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
1020     ierr  = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
1021     if (popup) {ierr  = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);}
1022     for(j = 0; j < n; j++) {
1023       x_l = j;
1024       x_r = x_l + 1.0;
1025       for(i = 0; i < m; i++) {
1026         y_l   = m - i - 1.0;
1027         y_r   = y_l + 1.0;
1028         color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i]));
1029         ierr  = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
1030       }
1031     }
1032   }
1033   PetscFunctionReturn(0);
1034 }
1035 
1036 #undef __FUNCT__
1037 #define __FUNCT__ "MatView_SeqDense_Draw"
1038 PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
1039 {
1040   PetscDraw      draw;
1041   PetscTruth     isnull;
1042   PetscReal      xr,yr,xl,yl,h,w;
1043   PetscErrorCode ierr;
1044 
1045   PetscFunctionBegin;
1046   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
1047   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
1048   if (isnull) PetscFunctionReturn(0);
1049 
1050   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
1051   xr  = A->cmap->n; yr = A->rmap->n; h = yr/10.0; w = xr/10.0;
1052   xr += w;    yr += h;  xl = -w;     yl = -h;
1053   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
1054   ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr);
1055   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
1056   PetscFunctionReturn(0);
1057 }
1058 
1059 #undef __FUNCT__
1060 #define __FUNCT__ "MatView_SeqDense"
1061 PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer)
1062 {
1063   PetscErrorCode ierr;
1064   PetscTruth     iascii,isbinary,isdraw;
1065 
1066   PetscFunctionBegin;
1067   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1068   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
1069   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
1070 
1071   if (iascii) {
1072     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
1073   } else if (isbinary) {
1074     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
1075   } else if (isdraw) {
1076     ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr);
1077   } else {
1078     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name);
1079   }
1080   PetscFunctionReturn(0);
1081 }
1082 
1083 #undef __FUNCT__
1084 #define __FUNCT__ "MatDestroy_SeqDense"
1085 PetscErrorCode MatDestroy_SeqDense(Mat mat)
1086 {
1087   Mat_SeqDense   *l = (Mat_SeqDense*)mat->data;
1088   PetscErrorCode ierr;
1089 
1090   PetscFunctionBegin;
1091 #if defined(PETSC_USE_LOG)
1092   PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap->n,mat->cmap->n);
1093 #endif
1094   ierr = PetscFree(l->pivots);CHKERRQ(ierr);
1095   if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);}
1096   ierr = PetscFree(l);CHKERRQ(ierr);
1097 
1098   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
1099   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
1100   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
1101   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
1102   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
1103   PetscFunctionReturn(0);
1104 }
1105 
1106 #undef __FUNCT__
1107 #define __FUNCT__ "MatTranspose_SeqDense"
1108 PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout)
1109 {
1110   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1111   PetscErrorCode ierr;
1112   PetscInt       k,j,m,n,M;
1113   PetscScalar    *v,tmp;
1114 
1115   PetscFunctionBegin;
1116   v = mat->v; m = A->rmap->n; M = mat->lda; n = A->cmap->n;
1117   if (reuse == MAT_REUSE_MATRIX && *matout == A) { /* in place transpose */
1118     if (m != n) {
1119       SETERRQ(PETSC_ERR_SUP,"Can not transpose non-square matrix in place");
1120     } else {
1121       for (j=0; j<m; j++) {
1122         for (k=0; k<j; k++) {
1123           tmp = v[j + k*M];
1124           v[j + k*M] = v[k + j*M];
1125           v[k + j*M] = tmp;
1126         }
1127       }
1128     }
1129   } else { /* out-of-place transpose */
1130     Mat          tmat;
1131     Mat_SeqDense *tmatd;
1132     PetscScalar  *v2;
1133 
1134     if (reuse == MAT_INITIAL_MATRIX) {
1135       ierr  = MatCreate(((PetscObject)A)->comm,&tmat);CHKERRQ(ierr);
1136       ierr  = MatSetSizes(tmat,A->cmap->n,A->rmap->n,A->cmap->n,A->rmap->n);CHKERRQ(ierr);
1137       ierr  = MatSetType(tmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
1138       ierr  = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr);
1139     } else {
1140       tmat = *matout;
1141     }
1142     tmatd = (Mat_SeqDense*)tmat->data;
1143     v = mat->v; v2 = tmatd->v;
1144     for (j=0; j<n; j++) {
1145       for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M];
1146     }
1147     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1148     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1149     *matout = tmat;
1150   }
1151   PetscFunctionReturn(0);
1152 }
1153 
1154 #undef __FUNCT__
1155 #define __FUNCT__ "MatEqual_SeqDense"
1156 PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg)
1157 {
1158   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1159   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
1160   PetscInt     i,j;
1161   PetscScalar  *v1 = mat1->v,*v2 = mat2->v;
1162 
1163   PetscFunctionBegin;
1164   if (A1->rmap->n != A2->rmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1165   if (A1->cmap->n != A2->cmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1166   for (i=0; i<A1->rmap->n; i++) {
1167     v1 = mat1->v+i; v2 = mat2->v+i;
1168     for (j=0; j<A1->cmap->n; j++) {
1169       if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1170       v1 += mat1->lda; v2 += mat2->lda;
1171     }
1172   }
1173   *flg = PETSC_TRUE;
1174   PetscFunctionReturn(0);
1175 }
1176 
1177 #undef __FUNCT__
1178 #define __FUNCT__ "MatGetDiagonal_SeqDense"
1179 PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v)
1180 {
1181   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1182   PetscErrorCode ierr;
1183   PetscInt       i,n,len;
1184   PetscScalar    *x,zero = 0.0;
1185 
1186   PetscFunctionBegin;
1187   ierr = VecSet(v,zero);CHKERRQ(ierr);
1188   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
1189   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1190   len = PetscMin(A->rmap->n,A->cmap->n);
1191   if (n != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
1192   for (i=0; i<len; i++) {
1193     x[i] = mat->v[i*mat->lda + i];
1194   }
1195   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1196   PetscFunctionReturn(0);
1197 }
1198 
1199 #undef __FUNCT__
1200 #define __FUNCT__ "MatDiagonalScale_SeqDense"
1201 PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1202 {
1203   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1204   PetscScalar    *l,*r,x,*v;
1205   PetscErrorCode ierr;
1206   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n;
1207 
1208   PetscFunctionBegin;
1209   if (ll) {
1210     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
1211     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1212     if (m != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1213     for (i=0; i<m; i++) {
1214       x = l[i];
1215       v = mat->v + i;
1216       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1217     }
1218     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1219     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1220   }
1221   if (rr) {
1222     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
1223     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1224     if (n != A->cmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1225     for (i=0; i<n; i++) {
1226       x = r[i];
1227       v = mat->v + i*m;
1228       for (j=0; j<m; j++) { (*v++) *= x;}
1229     }
1230     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1231     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1232   }
1233   PetscFunctionReturn(0);
1234 }
1235 
1236 #undef __FUNCT__
1237 #define __FUNCT__ "MatNorm_SeqDense"
1238 PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1239 {
1240   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1241   PetscScalar  *v = mat->v;
1242   PetscReal    sum = 0.0;
1243   PetscInt     lda=mat->lda,m=A->rmap->n,i,j;
1244   PetscErrorCode ierr;
1245 
1246   PetscFunctionBegin;
1247   if (type == NORM_FROBENIUS) {
1248     if (lda>m) {
1249       for (j=0; j<A->cmap->n; j++) {
1250 	v = mat->v+j*lda;
1251 	for (i=0; i<m; i++) {
1252 #if defined(PETSC_USE_COMPLEX)
1253 	  sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1254 #else
1255 	  sum += (*v)*(*v); v++;
1256 #endif
1257 	}
1258       }
1259     } else {
1260       for (i=0; i<A->cmap->n*A->rmap->n; i++) {
1261 #if defined(PETSC_USE_COMPLEX)
1262 	sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1263 #else
1264 	sum += (*v)*(*v); v++;
1265 #endif
1266       }
1267     }
1268     *nrm = sqrt(sum);
1269     ierr = PetscLogFlops(2.0*A->cmap->n*A->rmap->n);CHKERRQ(ierr);
1270   } else if (type == NORM_1) {
1271     *nrm = 0.0;
1272     for (j=0; j<A->cmap->n; j++) {
1273       v = mat->v + j*mat->lda;
1274       sum = 0.0;
1275       for (i=0; i<A->rmap->n; i++) {
1276         sum += PetscAbsScalar(*v);  v++;
1277       }
1278       if (sum > *nrm) *nrm = sum;
1279     }
1280     ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr);
1281   } else if (type == NORM_INFINITY) {
1282     *nrm = 0.0;
1283     for (j=0; j<A->rmap->n; j++) {
1284       v = mat->v + j;
1285       sum = 0.0;
1286       for (i=0; i<A->cmap->n; i++) {
1287         sum += PetscAbsScalar(*v); v += mat->lda;
1288       }
1289       if (sum > *nrm) *nrm = sum;
1290     }
1291     ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr);
1292   } else {
1293     SETERRQ(PETSC_ERR_SUP,"No two norm");
1294   }
1295   PetscFunctionReturn(0);
1296 }
1297 
1298 #undef __FUNCT__
1299 #define __FUNCT__ "MatSetOption_SeqDense"
1300 PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscTruth flg)
1301 {
1302   Mat_SeqDense   *aij = (Mat_SeqDense*)A->data;
1303   PetscErrorCode ierr;
1304 
1305   PetscFunctionBegin;
1306   switch (op) {
1307   case MAT_ROW_ORIENTED:
1308     aij->roworiented = flg;
1309     break;
1310   case MAT_NEW_NONZERO_LOCATIONS:
1311   case MAT_NEW_NONZERO_LOCATION_ERR:
1312   case MAT_NEW_NONZERO_ALLOCATION_ERR:
1313   case MAT_NEW_DIAGONALS:
1314   case MAT_IGNORE_OFF_PROC_ENTRIES:
1315   case MAT_USE_HASH_TABLE:
1316   case MAT_SYMMETRIC:
1317   case MAT_STRUCTURALLY_SYMMETRIC:
1318   case MAT_HERMITIAN:
1319   case MAT_SYMMETRY_ETERNAL:
1320   case MAT_IGNORE_LOWER_TRIANGULAR:
1321     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
1322     break;
1323   default:
1324     SETERRQ1(PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
1325   }
1326   PetscFunctionReturn(0);
1327 }
1328 
1329 #undef __FUNCT__
1330 #define __FUNCT__ "MatZeroEntries_SeqDense"
1331 PetscErrorCode MatZeroEntries_SeqDense(Mat A)
1332 {
1333   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
1334   PetscErrorCode ierr;
1335   PetscInt       lda=l->lda,m=A->rmap->n,j;
1336 
1337   PetscFunctionBegin;
1338   if (lda>m) {
1339     for (j=0; j<A->cmap->n; j++) {
1340       ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
1341     }
1342   } else {
1343     ierr = PetscMemzero(l->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
1344   }
1345   PetscFunctionReturn(0);
1346 }
1347 
1348 #undef __FUNCT__
1349 #define __FUNCT__ "MatZeroRows_SeqDense"
1350 PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
1351 {
1352   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
1353   PetscInt       n = A->cmap->n,i,j;
1354   PetscScalar    *slot;
1355 
1356   PetscFunctionBegin;
1357   for (i=0; i<N; i++) {
1358     slot = l->v + rows[i];
1359     for (j=0; j<n; j++) { *slot = 0.0; slot += n;}
1360   }
1361   if (diag != 0.0) {
1362     for (i=0; i<N; i++) {
1363       slot = l->v + (n+1)*rows[i];
1364       *slot = diag;
1365     }
1366   }
1367   PetscFunctionReturn(0);
1368 }
1369 
1370 #undef __FUNCT__
1371 #define __FUNCT__ "MatGetArray_SeqDense"
1372 PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[])
1373 {
1374   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1375 
1376   PetscFunctionBegin;
1377   if (mat->lda != A->rmap->n) SETERRQ(PETSC_ERR_SUP,"Cannot get array for Dense matrices with LDA different from number of rows");
1378   *array = mat->v;
1379   PetscFunctionReturn(0);
1380 }
1381 
1382 #undef __FUNCT__
1383 #define __FUNCT__ "MatRestoreArray_SeqDense"
1384 PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[])
1385 {
1386   PetscFunctionBegin;
1387   *array = 0; /* user cannot accidently use the array later */
1388   PetscFunctionReturn(0);
1389 }
1390 
1391 #undef __FUNCT__
1392 #define __FUNCT__ "MatGetSubMatrix_SeqDense"
1393 static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
1394 {
1395   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1396   PetscErrorCode ierr;
1397   PetscInt       i,j,nrows,ncols;
1398   const PetscInt *irow,*icol;
1399   PetscScalar    *av,*bv,*v = mat->v;
1400   Mat            newmat;
1401 
1402   PetscFunctionBegin;
1403   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
1404   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
1405   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
1406   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
1407 
1408   /* Check submatrixcall */
1409   if (scall == MAT_REUSE_MATRIX) {
1410     PetscInt n_cols,n_rows;
1411     ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr);
1412     if (n_rows != nrows || n_cols != ncols) {
1413       /* resize the result result matrix to match number of requested rows/columns */
1414       ierr = MatSetSizes(*B,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
1415     }
1416     newmat = *B;
1417   } else {
1418     /* Create and fill new matrix */
1419     ierr = MatCreate(((PetscObject)A)->comm,&newmat);CHKERRQ(ierr);
1420     ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
1421     ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
1422     ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
1423   }
1424 
1425   /* Now extract the data pointers and do the copy,column at a time */
1426   bv = ((Mat_SeqDense*)newmat->data)->v;
1427 
1428   for (i=0; i<ncols; i++) {
1429     av = v + mat->lda*icol[i];
1430     for (j=0; j<nrows; j++) {
1431       *bv++ = av[irow[j]];
1432     }
1433   }
1434 
1435   /* Assemble the matrices so that the correct flags are set */
1436   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1437   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1438 
1439   /* Free work space */
1440   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
1441   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
1442   *B = newmat;
1443   PetscFunctionReturn(0);
1444 }
1445 
1446 #undef __FUNCT__
1447 #define __FUNCT__ "MatGetSubMatrices_SeqDense"
1448 PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1449 {
1450   PetscErrorCode ierr;
1451   PetscInt       i;
1452 
1453   PetscFunctionBegin;
1454   if (scall == MAT_INITIAL_MATRIX) {
1455     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
1456   }
1457 
1458   for (i=0; i<n; i++) {
1459     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1460   }
1461   PetscFunctionReturn(0);
1462 }
1463 
1464 #undef __FUNCT__
1465 #define __FUNCT__ "MatAssemblyBegin_SeqDense"
1466 PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode)
1467 {
1468   PetscFunctionBegin;
1469   PetscFunctionReturn(0);
1470 }
1471 
1472 #undef __FUNCT__
1473 #define __FUNCT__ "MatAssemblyEnd_SeqDense"
1474 PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode)
1475 {
1476   PetscFunctionBegin;
1477   PetscFunctionReturn(0);
1478 }
1479 
1480 #undef __FUNCT__
1481 #define __FUNCT__ "MatCopy_SeqDense"
1482 PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
1483 {
1484   Mat_SeqDense   *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
1485   PetscErrorCode ierr;
1486   PetscInt       lda1=a->lda,lda2=b->lda, m=A->rmap->n,n=A->cmap->n, j;
1487 
1488   PetscFunctionBegin;
1489   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
1490   if (A->ops->copy != B->ops->copy) {
1491     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
1492     PetscFunctionReturn(0);
1493   }
1494   if (m != B->rmap->n || n != B->cmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1495   if (lda1>m || lda2>m) {
1496     for (j=0; j<n; j++) {
1497       ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1498     }
1499   } else {
1500     ierr = PetscMemcpy(b->v,a->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
1501   }
1502   PetscFunctionReturn(0);
1503 }
1504 
1505 #undef __FUNCT__
1506 #define __FUNCT__ "MatSetUpPreallocation_SeqDense"
1507 PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A)
1508 {
1509   PetscErrorCode ierr;
1510 
1511   PetscFunctionBegin;
1512   ierr =  MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr);
1513   PetscFunctionReturn(0);
1514 }
1515 
1516 #undef __FUNCT__
1517 #define __FUNCT__ "MatSetSizes_SeqDense"
1518 PetscErrorCode MatSetSizes_SeqDense(Mat A,PetscInt m,PetscInt n,PetscInt M,PetscInt N)
1519 {
1520   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1521   PetscFunctionBegin;
1522   /* this will not be called before lda, Mmax,  and Nmax have been set */
1523   m = PetscMax(m,M);
1524   n = PetscMax(n,N);
1525   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);
1526   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);
1527   A->rmap->n = A->rmap->n = m;
1528   A->cmap->n = A->cmap->N = n;
1529   if (a->changelda) a->lda = m;
1530   PetscFunctionReturn(0);
1531 }
1532 
1533 
1534 /* ----------------------------------------------------------------*/
1535 #undef __FUNCT__
1536 #define __FUNCT__ "MatMatMult_SeqDense_SeqDense"
1537 PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1538 {
1539   PetscErrorCode ierr;
1540 
1541   PetscFunctionBegin;
1542   if (scall == MAT_INITIAL_MATRIX){
1543     ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1544   }
1545   ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1546   PetscFunctionReturn(0);
1547 }
1548 
1549 #undef __FUNCT__
1550 #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense"
1551 PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1552 {
1553   PetscErrorCode ierr;
1554   PetscInt       m=A->rmap->n,n=B->cmap->n;
1555   Mat            Cmat;
1556 
1557   PetscFunctionBegin;
1558   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);
1559   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1560   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1561   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1562   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1563   Cmat->assembled = PETSC_TRUE;
1564   *C = Cmat;
1565   PetscFunctionReturn(0);
1566 }
1567 
1568 #undef __FUNCT__
1569 #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense"
1570 PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1571 {
1572   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1573   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1574   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
1575   PetscBLASInt   m,n,k;
1576   PetscScalar    _DOne=1.0,_DZero=0.0;
1577 
1578   PetscFunctionBegin;
1579   m = PetscBLASIntCast(A->rmap->n);
1580   n = PetscBLASIntCast(B->cmap->n);
1581   k = PetscBLASIntCast(A->cmap->n);
1582   BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1583   PetscFunctionReturn(0);
1584 }
1585 
1586 #undef __FUNCT__
1587 #define __FUNCT__ "MatMatMultTranspose_SeqDense_SeqDense"
1588 PetscErrorCode MatMatMultTranspose_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1589 {
1590   PetscErrorCode ierr;
1591 
1592   PetscFunctionBegin;
1593   if (scall == MAT_INITIAL_MATRIX){
1594     ierr = MatMatMultTransposeSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1595   }
1596   ierr = MatMatMultTransposeNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1597   PetscFunctionReturn(0);
1598 }
1599 
1600 #undef __FUNCT__
1601 #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqDense_SeqDense"
1602 PetscErrorCode MatMatMultTransposeSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1603 {
1604   PetscErrorCode ierr;
1605   PetscInt       m=A->cmap->n,n=B->cmap->n;
1606   Mat            Cmat;
1607 
1608   PetscFunctionBegin;
1609   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);
1610   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1611   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1612   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1613   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1614   Cmat->assembled = PETSC_TRUE;
1615   *C = Cmat;
1616   PetscFunctionReturn(0);
1617 }
1618 
1619 #undef __FUNCT__
1620 #define __FUNCT__ "MatMatMultTransposeNumeric_SeqDense_SeqDense"
1621 PetscErrorCode MatMatMultTransposeNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1622 {
1623   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1624   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1625   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
1626   PetscBLASInt   m,n,k;
1627   PetscScalar    _DOne=1.0,_DZero=0.0;
1628 
1629   PetscFunctionBegin;
1630   m = PetscBLASIntCast(A->cmap->n);
1631   n = PetscBLASIntCast(B->cmap->n);
1632   k = PetscBLASIntCast(A->rmap->n);
1633   /*
1634      Note the m and n arguments below are the number rows and columns of A', not A!
1635   */
1636   BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1637   PetscFunctionReturn(0);
1638 }
1639 
1640 #undef __FUNCT__
1641 #define __FUNCT__ "MatGetRowMax_SeqDense"
1642 PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[])
1643 {
1644   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1645   PetscErrorCode ierr;
1646   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1647   PetscScalar    *x;
1648   MatScalar      *aa = a->v;
1649 
1650   PetscFunctionBegin;
1651   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1652 
1653   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1654   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1655   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1656   if (p != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1657   for (i=0; i<m; i++) {
1658     x[i] = aa[i]; if (idx) idx[i] = 0;
1659     for (j=1; j<n; j++){
1660       if (PetscRealPart(x[i]) < PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
1661     }
1662   }
1663   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1664   PetscFunctionReturn(0);
1665 }
1666 
1667 #undef __FUNCT__
1668 #define __FUNCT__ "MatGetRowMaxAbs_SeqDense"
1669 PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[])
1670 {
1671   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1672   PetscErrorCode ierr;
1673   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1674   PetscScalar    *x;
1675   PetscReal      atmp;
1676   MatScalar      *aa = a->v;
1677 
1678   PetscFunctionBegin;
1679   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1680 
1681   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1682   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1683   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1684   if (p != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1685   for (i=0; i<m; i++) {
1686     x[i] = PetscAbsScalar(aa[i]);
1687     for (j=1; j<n; j++){
1688       atmp = PetscAbsScalar(aa[i+m*j]);
1689       if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;}
1690     }
1691   }
1692   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1693   PetscFunctionReturn(0);
1694 }
1695 
1696 #undef __FUNCT__
1697 #define __FUNCT__ "MatGetRowMin_SeqDense"
1698 PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[])
1699 {
1700   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1701   PetscErrorCode ierr;
1702   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1703   PetscScalar    *x;
1704   MatScalar      *aa = a->v;
1705 
1706   PetscFunctionBegin;
1707   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1708 
1709   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1710   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1711   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1712   if (p != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1713   for (i=0; i<m; i++) {
1714     x[i] = aa[i]; if (idx) idx[i] = 0;
1715     for (j=1; j<n; j++){
1716       if (PetscRealPart(x[i]) > PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
1717     }
1718   }
1719   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1720   PetscFunctionReturn(0);
1721 }
1722 
1723 #undef __FUNCT__
1724 #define __FUNCT__ "MatGetColumnVector_SeqDense"
1725 PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col)
1726 {
1727   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1728   PetscErrorCode ierr;
1729   PetscScalar    *x;
1730 
1731   PetscFunctionBegin;
1732   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1733 
1734   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1735   ierr = PetscMemcpy(x,a->v+col*a->lda,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
1736   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1737   PetscFunctionReturn(0);
1738 }
1739 
1740 /* -------------------------------------------------------------------*/
1741 static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1742        MatGetRow_SeqDense,
1743        MatRestoreRow_SeqDense,
1744        MatMult_SeqDense,
1745 /* 4*/ MatMultAdd_SeqDense,
1746        MatMultTranspose_SeqDense,
1747        MatMultTransposeAdd_SeqDense,
1748        0,
1749        0,
1750        0,
1751 /*10*/ 0,
1752        MatLUFactor_SeqDense,
1753        MatCholeskyFactor_SeqDense,
1754        MatRelax_SeqDense,
1755        MatTranspose_SeqDense,
1756 /*15*/ MatGetInfo_SeqDense,
1757        MatEqual_SeqDense,
1758        MatGetDiagonal_SeqDense,
1759        MatDiagonalScale_SeqDense,
1760        MatNorm_SeqDense,
1761 /*20*/ MatAssemblyBegin_SeqDense,
1762        MatAssemblyEnd_SeqDense,
1763        MatSetOption_SeqDense,
1764        MatZeroEntries_SeqDense,
1765 /*24*/ MatZeroRows_SeqDense,
1766        0,
1767        0,
1768        0,
1769        0,
1770 /*29*/ MatSetUpPreallocation_SeqDense,
1771        0,
1772        0,
1773        MatGetArray_SeqDense,
1774        MatRestoreArray_SeqDense,
1775 /*34*/ MatDuplicate_SeqDense,
1776        0,
1777        0,
1778        0,
1779        0,
1780 /*39*/ MatAXPY_SeqDense,
1781        MatGetSubMatrices_SeqDense,
1782        0,
1783        MatGetValues_SeqDense,
1784        MatCopy_SeqDense,
1785 /*44*/ MatGetRowMax_SeqDense,
1786        MatScale_SeqDense,
1787        0,
1788        0,
1789        0,
1790 /*49*/ 0,
1791        0,
1792        0,
1793        0,
1794        0,
1795 /*54*/ 0,
1796        0,
1797        0,
1798        0,
1799        0,
1800 /*59*/ 0,
1801        MatDestroy_SeqDense,
1802        MatView_SeqDense,
1803        0,
1804        0,
1805 /*64*/ 0,
1806        0,
1807        0,
1808        0,
1809        0,
1810 /*69*/ MatGetRowMaxAbs_SeqDense,
1811        0,
1812        0,
1813        0,
1814        0,
1815 /*74*/ 0,
1816        0,
1817        0,
1818        0,
1819        0,
1820 /*79*/ 0,
1821        0,
1822        0,
1823        0,
1824 /*83*/ MatLoad_SeqDense,
1825        0,
1826        MatIsHermitian_SeqDense,
1827        0,
1828        0,
1829        0,
1830 /*89*/ MatMatMult_SeqDense_SeqDense,
1831        MatMatMultSymbolic_SeqDense_SeqDense,
1832        MatMatMultNumeric_SeqDense_SeqDense,
1833        0,
1834        0,
1835 /*94*/ 0,
1836        MatMatMultTranspose_SeqDense_SeqDense,
1837        MatMatMultTransposeSymbolic_SeqDense_SeqDense,
1838        MatMatMultTransposeNumeric_SeqDense_SeqDense,
1839        0,
1840 /*99*/ 0,
1841        0,
1842        0,
1843        0,
1844        MatSetSizes_SeqDense,
1845        0,
1846        0,
1847        0,
1848        0,
1849        0,
1850 /*109*/0,
1851        0,
1852        MatGetRowMin_SeqDense,
1853        MatGetColumnVector_SeqDense
1854 };
1855 
1856 #undef __FUNCT__
1857 #define __FUNCT__ "MatCreateSeqDense"
1858 /*@C
1859    MatCreateSeqDense - Creates a sequential dense matrix that
1860    is stored in column major order (the usual Fortran 77 manner). Many
1861    of the matrix operations use the BLAS and LAPACK routines.
1862 
1863    Collective on MPI_Comm
1864 
1865    Input Parameters:
1866 +  comm - MPI communicator, set to PETSC_COMM_SELF
1867 .  m - number of rows
1868 .  n - number of columns
1869 -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1870    to control all matrix memory allocation.
1871 
1872    Output Parameter:
1873 .  A - the matrix
1874 
1875    Notes:
1876    The data input variable is intended primarily for Fortran programmers
1877    who wish to allocate their own matrix memory space.  Most users should
1878    set data=PETSC_NULL.
1879 
1880    Level: intermediate
1881 
1882 .keywords: dense, matrix, LAPACK, BLAS
1883 
1884 .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
1885 @*/
1886 PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A)
1887 {
1888   PetscErrorCode ierr;
1889 
1890   PetscFunctionBegin;
1891   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1892   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
1893   ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1894   ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1895   PetscFunctionReturn(0);
1896 }
1897 
1898 #undef __FUNCT__
1899 #define __FUNCT__ "MatSeqDenseSetPreallocation"
1900 /*@C
1901    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
1902 
1903    Collective on MPI_Comm
1904 
1905    Input Parameters:
1906 +  A - the matrix
1907 -  data - the array (or PETSC_NULL)
1908 
1909    Notes:
1910    The data input variable is intended primarily for Fortran programmers
1911    who wish to allocate their own matrix memory space.  Most users should
1912    need not call this routine.
1913 
1914    Level: intermediate
1915 
1916 .keywords: dense, matrix, LAPACK, BLAS
1917 
1918 .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues(), MatSeqDenseSetLDA()
1919 
1920 @*/
1921 PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
1922 {
1923   PetscErrorCode ierr,(*f)(Mat,PetscScalar[]);
1924 
1925   PetscFunctionBegin;
1926   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1927   if (f) {
1928     ierr = (*f)(B,data);CHKERRQ(ierr);
1929   }
1930   PetscFunctionReturn(0);
1931 }
1932 
1933 EXTERN_C_BEGIN
1934 #undef __FUNCT__
1935 #define __FUNCT__ "MatSeqDenseSetPreallocation_SeqDense"
1936 PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
1937 {
1938   Mat_SeqDense   *b;
1939   PetscErrorCode ierr;
1940 
1941   PetscFunctionBegin;
1942   B->preallocated = PETSC_TRUE;
1943   b               = (Mat_SeqDense*)B->data;
1944   if (b->lda <= 0) b->lda = B->rmap->n;
1945   if (!data) { /* petsc-allocated storage */
1946     if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); }
1947     ierr = PetscMalloc(b->lda*b->Nmax*sizeof(PetscScalar),&b->v);CHKERRQ(ierr);
1948     ierr = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
1949     ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
1950     b->user_alloc = PETSC_FALSE;
1951   } else { /* user-allocated storage */
1952     if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); }
1953     b->v          = data;
1954     b->user_alloc = PETSC_TRUE;
1955   }
1956   B->assembled = PETSC_TRUE;
1957   PetscFunctionReturn(0);
1958 }
1959 EXTERN_C_END
1960 
1961 #undef __FUNCT__
1962 #define __FUNCT__ "MatSeqDenseSetLDA"
1963 /*@C
1964   MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array
1965 
1966   Input parameter:
1967 + A - the matrix
1968 - lda - the leading dimension
1969 
1970   Notes:
1971   This routine is to be used in conjunction with MatSeqDenseSetPreallocation();
1972   it asserts that the preallocation has a leading dimension (the LDA parameter
1973   of Blas and Lapack fame) larger than M, the first dimension of the matrix.
1974 
1975   Level: intermediate
1976 
1977 .keywords: dense, matrix, LAPACK, BLAS
1978 
1979 .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize()
1980 
1981 @*/
1982 PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat B,PetscInt lda)
1983 {
1984   Mat_SeqDense *b = (Mat_SeqDense*)B->data;
1985 
1986   PetscFunctionBegin;
1987   if (lda < B->rmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->rmap->n);
1988   b->lda       = lda;
1989   b->changelda = PETSC_FALSE;
1990   b->Mmax      = PetscMax(b->Mmax,lda);
1991   PetscFunctionReturn(0);
1992 }
1993 
1994 /*MC
1995    MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices.
1996 
1997    Options Database Keys:
1998 . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions()
1999 
2000   Level: beginner
2001 
2002 .seealso: MatCreateSeqDense()
2003 
2004 M*/
2005 
2006 EXTERN_C_BEGIN
2007 #undef __FUNCT__
2008 #define __FUNCT__ "MatCreate_SeqDense"
2009 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqDense(Mat B)
2010 {
2011   Mat_SeqDense   *b;
2012   PetscErrorCode ierr;
2013   PetscMPIInt    size;
2014 
2015   PetscFunctionBegin;
2016   ierr = MPI_Comm_size(((PetscObject)B)->comm,&size);CHKERRQ(ierr);
2017   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
2018 
2019   ierr = PetscMapSetBlockSize(B->rmap,1);CHKERRQ(ierr);
2020   ierr = PetscMapSetBlockSize(B->cmap,1);CHKERRQ(ierr);
2021   ierr = PetscMapSetUp(B->rmap);CHKERRQ(ierr);
2022   ierr = PetscMapSetUp(B->cmap);CHKERRQ(ierr);
2023 
2024   ierr            = PetscNewLog(B,Mat_SeqDense,&b);CHKERRQ(ierr);
2025   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
2026   B->mapping      = 0;
2027   B->data         = (void*)b;
2028 
2029 
2030   b->pivots       = 0;
2031   b->roworiented  = PETSC_TRUE;
2032   b->v            = 0;
2033   b->lda          = B->rmap->n;
2034   b->changelda    = PETSC_FALSE;
2035   b->Mmax         = B->rmap->n;
2036   b->Nmax         = B->cmap->n;
2037 
2038 
2039   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_seqdense_petsc_C",
2040                                      "MatGetFactor_seqdense_petsc",
2041                                       MatGetFactor_seqdense_petsc);CHKERRQ(ierr);
2042   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C",
2043                                     "MatSeqDenseSetPreallocation_SeqDense",
2044                                      MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr);
2045   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMult_seqaij_seqdense_C",
2046                                      "MatMatMult_SeqAIJ_SeqDense",
2047                                       MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr);
2048   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C",
2049                                      "MatMatMultSymbolic_SeqAIJ_SeqDense",
2050                                       MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr);
2051   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C",
2052                                      "MatMatMultNumeric_SeqAIJ_SeqDense",
2053                                       MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr);
2054   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);CHKERRQ(ierr);
2055   PetscFunctionReturn(0);
2056 }
2057 
2058 
2059 EXTERN_C_END
2060