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