Actual source code: dense.c

petsc-master 2015-03-28
Report Typos and Errors
  2: /*
  3:      Defines the basic matrix operations for sequential dense.
  4: */

  6: #include <../src/mat/impls/dense/seq/dense.h> /*I "petscmat.h" I*/
  7: #include <petscblaslapack.h>

  9: #include <../src/mat/impls/aij/seq/aij.h>

 13: PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqDense(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
 14: {
 15:   Mat            B;
 16:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
 17:   Mat_SeqDense   *b;
 19:   PetscInt       *ai=a->i,*aj=a->j,m=A->rmap->N,n=A->cmap->N,i;
 20:   MatScalar      *av=a->a;

 23:   MatCreate(PetscObjectComm((PetscObject)A),&B);
 24:   MatSetSizes(B,m,n,m,n);
 25:   MatSetType(B,MATSEQDENSE);
 26:   MatSeqDenseSetPreallocation(B,NULL);

 28:   b  = (Mat_SeqDense*)(B->data);
 29:   for (i=0; i<m; i++) {
 30:     PetscInt j;
 31:     for (j=0;j<ai[1]-ai[0];j++) {
 32:       b->v[*aj*m+i] = *av;
 33:       aj++;
 34:       av++;
 35:     }
 36:     ai++;
 37:   }
 38:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 39:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

 41:   if (reuse == MAT_REUSE_MATRIX) {
 42:     MatHeaderReplace(A,B);
 43:   } else {
 44:     *newmat = B;
 45:   }
 46:   return(0);
 47: }

 51: PETSC_EXTERN PetscErrorCode MatConvert_SeqDense_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
 52: {
 53:   Mat            B;
 54:   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
 56:   PetscInt       i, j;
 57:   PetscInt       *rows, *nnz;
 58:   MatScalar      *aa = a->v, *vals;

 61:   MatCreate(PetscObjectComm((PetscObject)A),&B);
 62:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
 63:   MatSetType(B,MATSEQAIJ);
 64:   PetscCalloc3(A->rmap->n,&rows,A->rmap->n,&nnz,A->rmap->n,&vals);
 65:   for (j=0; j<A->cmap->n; j++) {
 66:     for (i=0; i<A->rmap->n; i++) if (aa[i] != 0.0 || i == j) ++nnz[i];
 67:     aa += a->lda;
 68:   }
 69:   MatSeqAIJSetPreallocation(B,PETSC_DETERMINE,nnz);
 70:   aa = a->v;
 71:   for (j=0; j<A->cmap->n; j++) {
 72:     PetscInt numRows = 0;
 73:     for (i=0; i<A->rmap->n; i++) if (aa[i] != 0.0 || i == j) {rows[numRows] = i; vals[numRows++] = aa[i];}
 74:     MatSetValues(B,numRows,rows,1,&j,vals,INSERT_VALUES);
 75:     aa  += a->lda;
 76:   }
 77:   PetscFree3(rows,nnz,vals);
 78:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 79:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

 81:   if (reuse == MAT_REUSE_MATRIX) {
 82:     MatHeaderReplace(A,B);
 83:   } else {
 84:     *newmat = B;
 85:   }
 86:   return(0);
 87: }

 91: PetscErrorCode MatAXPY_SeqDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str)
 92: {
 93:   Mat_SeqDense   *x     = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data;
 94:   PetscScalar    oalpha = alpha;
 95:   PetscInt       j;
 96:   PetscBLASInt   N,m,ldax,lday,one = 1;

100:   PetscBLASIntCast(X->rmap->n*X->cmap->n,&N);
101:   PetscBLASIntCast(X->rmap->n,&m);
102:   PetscBLASIntCast(x->lda,&ldax);
103:   PetscBLASIntCast(y->lda,&lday);
104:   if (ldax>m || lday>m) {
105:     for (j=0; j<X->cmap->n; j++) {
106:       PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&m,&oalpha,x->v+j*ldax,&one,y->v+j*lday,&one));
107:     }
108:   } else {
109:     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&N,&oalpha,x->v,&one,y->v,&one));
110:   }
111:   PetscObjectStateIncrease((PetscObject)Y);
112:   PetscLogFlops(PetscMax(2*N-1,0));
113:   return(0);
114: }

118: PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
119: {
120:   PetscInt N = A->rmap->n*A->cmap->n;

123:   info->block_size        = 1.0;
124:   info->nz_allocated      = (double)N;
125:   info->nz_used           = (double)N;
126:   info->nz_unneeded       = (double)0;
127:   info->assemblies        = (double)A->num_ass;
128:   info->mallocs           = 0;
129:   info->memory            = ((PetscObject)A)->mem;
130:   info->fill_ratio_given  = 0;
131:   info->fill_ratio_needed = 0;
132:   info->factor_mallocs    = 0;
133:   return(0);
134: }

138: PetscErrorCode MatScale_SeqDense(Mat A,PetscScalar alpha)
139: {
140:   Mat_SeqDense   *a     = (Mat_SeqDense*)A->data;
141:   PetscScalar    oalpha = alpha;
143:   PetscBLASInt   one = 1,j,nz,lda;

146:   PetscBLASIntCast(a->lda,&lda);
147:   if (lda>A->rmap->n) {
148:     PetscBLASIntCast(A->rmap->n,&nz);
149:     for (j=0; j<A->cmap->n; j++) {
150:       PetscStackCallBLAS("BLASscal",BLASscal_(&nz,&oalpha,a->v+j*lda,&one));
151:     }
152:   } else {
153:     PetscBLASIntCast(A->rmap->n*A->cmap->n,&nz);
154:     PetscStackCallBLAS("BLASscal",BLASscal_(&nz,&oalpha,a->v,&one));
155:   }
156:   PetscLogFlops(nz);
157:   return(0);
158: }

162: PetscErrorCode MatIsHermitian_SeqDense(Mat A,PetscReal rtol,PetscBool  *fl)
163: {
164:   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
165:   PetscInt     i,j,m = A->rmap->n,N;
166:   PetscScalar  *v = a->v;

169:   *fl = PETSC_FALSE;
170:   if (A->rmap->n != A->cmap->n) return(0);
171:   N = a->lda;

173:   for (i=0; i<m; i++) {
174:     for (j=i+1; j<m; j++) {
175:       if (PetscAbsScalar(v[i+j*N] - PetscConj(v[j+i*N])) > rtol) return(0);
176:     }
177:   }
178:   *fl = PETSC_TRUE;
179:   return(0);
180: }

184: PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat newi,Mat A,MatDuplicateOption cpvalues)
185: {
186:   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data,*l;
188:   PetscInt       lda = (PetscInt)mat->lda,j,m;

191:   PetscLayoutReference(A->rmap,&newi->rmap);
192:   PetscLayoutReference(A->cmap,&newi->cmap);
193:   MatSeqDenseSetPreallocation(newi,NULL);
194:   if (cpvalues == MAT_COPY_VALUES) {
195:     l = (Mat_SeqDense*)newi->data;
196:     if (lda>A->rmap->n) {
197:       m = A->rmap->n;
198:       for (j=0; j<A->cmap->n; j++) {
199:         PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));
200:       }
201:     } else {
202:       PetscMemcpy(l->v,mat->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));
203:     }
204:   }
205:   newi->assembled = PETSC_TRUE;
206:   return(0);
207: }

211: PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
212: {

216:   MatCreate(PetscObjectComm((PetscObject)A),newmat);
217:   MatSetSizes(*newmat,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);
218:   MatSetType(*newmat,((PetscObject)A)->type_name);
219:   MatDuplicateNoCreate_SeqDense(*newmat,A,cpvalues);
220:   return(0);
221: }


224: extern PetscErrorCode MatLUFactor_SeqDense(Mat,IS,IS,const MatFactorInfo*);

228: PetscErrorCode MatLUFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy)
229: {
230:   MatFactorInfo  info;

234:   MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);
235:   MatLUFactor_SeqDense(fact,0,0,&info);
236:   return(0);
237: }

241: PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
242: {
243:   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
244:   PetscErrorCode    ierr;
245:   const PetscScalar *x;
246:   PetscScalar       *y;
247:   PetscBLASInt      one = 1,info,m;

250:   PetscBLASIntCast(A->rmap->n,&m);
251:   VecGetArrayRead(xx,&x);
252:   VecGetArray(yy,&y);
253:   PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));
254:   if (A->factortype == MAT_FACTOR_LU) {
255: #if defined(PETSC_MISSING_LAPACK_GETRS)
256:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
257: #else
258:     PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
259:     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve");
260: #endif
261:   } else if (A->factortype == MAT_FACTOR_CHOLESKY) {
262: #if defined(PETSC_MISSING_LAPACK_POTRS)
263:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
264: #else
265:     PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info));
266:     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve");
267: #endif
268:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
269:   VecRestoreArrayRead(xx,&x);
270:   VecRestoreArray(yy,&y);
271:   PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);
272:   return(0);
273: }

277: PetscErrorCode MatMatSolve_SeqDense(Mat A,Mat B,Mat X)
278: {
279:   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
281:   PetscScalar    *b,*x;
282:   PetscInt       n;
283:   PetscBLASInt   nrhs,info,m;
284:   PetscBool      flg;

287:   PetscBLASIntCast(A->rmap->n,&m);
288:   PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
289:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
290:   PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
291:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");

293:   MatGetSize(B,NULL,&n);
294:   PetscBLASIntCast(n,&nrhs);
295:   MatDenseGetArray(B,&b);
296:   MatDenseGetArray(X,&x);

298:   PetscMemcpy(x,b,m*nrhs*sizeof(PetscScalar));

300:   if (A->factortype == MAT_FACTOR_LU) {
301: #if defined(PETSC_MISSING_LAPACK_GETRS)
302:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
303: #else
304:     PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info));
305:     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve");
306: #endif
307:   } else if (A->factortype == MAT_FACTOR_CHOLESKY) {
308: #if defined(PETSC_MISSING_LAPACK_POTRS)
309:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
310: #else
311:     PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&m,&nrhs,mat->v,&mat->lda,x,&m,&info));
312:     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve");
313: #endif
314:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");

316:   MatDenseRestoreArray(B,&b);
317:   MatDenseRestoreArray(X,&x);
318:   PetscLogFlops(nrhs*(2.0*m*m - m));
319:   return(0);
320: }

324: PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
325: {
326:   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
327:   PetscErrorCode    ierr;
328:   const PetscScalar *x;
329:   PetscScalar       *y;
330:   PetscBLASInt      one = 1,info,m;

333:   PetscBLASIntCast(A->rmap->n,&m);
334:   VecGetArrayRead(xx,&x);
335:   VecGetArray(yy,&y);
336:   PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));
337:   /* assume if pivots exist then use LU; else Cholesky */
338:   if (mat->pivots) {
339: #if defined(PETSC_MISSING_LAPACK_GETRS)
340:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
341: #else
342:     PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
343:     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve");
344: #endif
345:   } else {
346: #if defined(PETSC_MISSING_LAPACK_POTRS)
347:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
348: #else
349:     PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info));
350:     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve");
351: #endif
352:   }
353:   VecRestoreArrayRead(xx,&x);
354:   VecRestoreArray(yy,&y);
355:   PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);
356:   return(0);
357: }

361: PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
362: {
363:   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
364:   PetscErrorCode    ierr;
365:   const PetscScalar *x;
366:   PetscScalar       *y,sone = 1.0;
367:   Vec               tmp = 0;
368:   PetscBLASInt      one = 1,info,m;

371:   PetscBLASIntCast(A->rmap->n,&m);
372:   VecGetArrayRead(xx,&x);
373:   VecGetArray(yy,&y);
374:   if (!A->rmap->n || !A->cmap->n) return(0);
375:   if (yy == zz) {
376:     VecDuplicate(yy,&tmp);
377:     PetscLogObjectParent((PetscObject)A,(PetscObject)tmp);
378:     VecCopy(yy,tmp);
379:   }
380:   PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));
381:   /* assume if pivots exist then use LU; else Cholesky */
382:   if (mat->pivots) {
383: #if defined(PETSC_MISSING_LAPACK_GETRS)
384:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
385: #else
386:     PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
387:     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve");
388: #endif
389:   } else {
390: #if defined(PETSC_MISSING_LAPACK_POTRS)
391:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
392: #else
393:     PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info));
394:     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve");
395: #endif
396:   }
397:   if (tmp) {
398:     VecAXPY(yy,sone,tmp);
399:     VecDestroy(&tmp);
400:   } else {
401:     VecAXPY(yy,sone,zz);
402:   }
403:   VecRestoreArrayRead(xx,&x);
404:   VecRestoreArray(yy,&y);
405:   PetscLogFlops(2.0*A->cmap->n*A->cmap->n);
406:   return(0);
407: }

411: PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
412: {
413:   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
414:   PetscErrorCode    ierr;
415:   const PetscScalar *x;
416:   PetscScalar       *y,sone = 1.0;
417:   Vec               tmp;
418:   PetscBLASInt      one = 1,info,m;

421:   PetscBLASIntCast(A->rmap->n,&m);
422:   if (!A->rmap->n || !A->cmap->n) return(0);
423:   VecGetArrayRead(xx,&x);
424:   VecGetArray(yy,&y);
425:   if (yy == zz) {
426:     VecDuplicate(yy,&tmp);
427:     PetscLogObjectParent((PetscObject)A,(PetscObject)tmp);
428:     VecCopy(yy,tmp);
429:   }
430:   PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));
431:   /* assume if pivots exist then use LU; else Cholesky */
432:   if (mat->pivots) {
433: #if defined(PETSC_MISSING_LAPACK_GETRS)
434:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
435: #else
436:     PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
437:     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve");
438: #endif
439:   } else {
440: #if defined(PETSC_MISSING_LAPACK_POTRS)
441:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
442: #else
443:     PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info));
444:     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve");
445: #endif
446:   }
447:   if (tmp) {
448:     VecAXPY(yy,sone,tmp);
449:     VecDestroy(&tmp);
450:   } else {
451:     VecAXPY(yy,sone,zz);
452:   }
453:   VecRestoreArrayRead(xx,&x);
454:   VecRestoreArray(yy,&y);
455:   PetscLogFlops(2.0*A->cmap->n*A->cmap->n);
456:   return(0);
457: }

459: /* ---------------------------------------------------------------*/
460: /* COMMENT: I have chosen to hide row permutation in the pivots,
461:    rather than put it in the Mat->row slot.*/
464: PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,const MatFactorInfo *minfo)
465: {
466: #if defined(PETSC_MISSING_LAPACK_GETRF)
468:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable.");
469: #else
470:   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
472:   PetscBLASInt   n,m,info;

475:   PetscBLASIntCast(A->cmap->n,&n);
476:   PetscBLASIntCast(A->rmap->n,&m);
477:   if (!mat->pivots) {
478:     PetscMalloc1(A->rmap->n+1,&mat->pivots);
479:     PetscLogObjectMemory((PetscObject)A,A->rmap->n*sizeof(PetscBLASInt));
480:   }
481:   if (!A->rmap->n || !A->cmap->n) return(0);
482:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
483:   PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info));
484:   PetscFPTrapPop();

486:   if (info<0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to LU factorization");
487:   if (info>0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
488:   A->ops->solve             = MatSolve_SeqDense;
489:   A->ops->solvetranspose    = MatSolveTranspose_SeqDense;
490:   A->ops->solveadd          = MatSolveAdd_SeqDense;
491:   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense;
492:   A->factortype             = MAT_FACTOR_LU;

494:   PetscLogFlops((2.0*A->cmap->n*A->cmap->n*A->cmap->n)/3);
495: #endif
496:   return(0);
497: }

501: PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,const MatFactorInfo *factinfo)
502: {
503: #if defined(PETSC_MISSING_LAPACK_POTRF)
505:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable.");
506: #else
507:   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
509:   PetscBLASInt   info,n;

512:   PetscBLASIntCast(A->cmap->n,&n);
513:   PetscFree(mat->pivots);

515:   if (!A->rmap->n || !A->cmap->n) return(0);
516:   PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info));
517:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1);
518:   A->ops->solve             = MatSolve_SeqDense;
519:   A->ops->solvetranspose    = MatSolveTranspose_SeqDense;
520:   A->ops->solveadd          = MatSolveAdd_SeqDense;
521:   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense;
522:   A->factortype             = MAT_FACTOR_CHOLESKY;

524:   PetscLogFlops((1.0*A->cmap->n*A->cmap->n*A->cmap->n)/3.0);
525: #endif
526:   return(0);
527: }


532: PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy)
533: {
535:   MatFactorInfo  info;

538:   info.fill = 1.0;

540:   MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);
541:   MatCholeskyFactor_SeqDense(fact,0,&info);
542:   return(0);
543: }

547: PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,const MatFactorInfo *info)
548: {
550:   fact->assembled                  = PETSC_TRUE;
551:   fact->preallocated               = PETSC_TRUE;
552:   fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqDense;
553:   return(0);
554: }

558: PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info)
559: {
561:   fact->preallocated         = PETSC_TRUE;
562:   fact->assembled            = PETSC_TRUE;
563:   fact->ops->lufactornumeric = MatLUFactorNumeric_SeqDense;
564:   return(0);
565: }

569: PETSC_EXTERN PetscErrorCode MatGetFactor_seqdense_petsc(Mat A,MatFactorType ftype,Mat *fact)
570: {

574:   MatCreate(PetscObjectComm((PetscObject)A),fact);
575:   MatSetSizes(*fact,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);
576:   MatSetType(*fact,((PetscObject)A)->type_name);
577:   if (ftype == MAT_FACTOR_LU) {
578:     (*fact)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqDense;
579:   } else {
580:     (*fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqDense;
581:   }
582:   (*fact)->factortype = ftype;
583:   return(0);
584: }

586: /* ------------------------------------------------------------------*/
589: PetscErrorCode MatSOR_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx)
590: {
591:   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
592:   PetscScalar       *x,*v = mat->v,zero = 0.0,xt;
593:   const PetscScalar *b;
594:   PetscErrorCode    ierr;
595:   PetscInt          m = A->rmap->n,i;
596:   PetscBLASInt      o = 1,bm;

599:   PetscBLASIntCast(m,&bm);
600:   if (flag & SOR_ZERO_INITIAL_GUESS) {
601:     /* this is a hack fix, should have another version without the second BLASdot */
602:     VecSet(xx,zero);
603:   }
604:   VecGetArray(xx,&x);
605:   VecGetArrayRead(bb,&b);
606:   its  = its*lits;
607:   if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
608:   while (its--) {
609:     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
610:       for (i=0; i<m; i++) {
611:         PetscStackCallBLAS("BLASdot",xt   = b[i] - BLASdot_(&bm,v+i,&bm,x,&o));
612:         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
613:       }
614:     }
615:     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
616:       for (i=m-1; i>=0; i--) {
617:         PetscStackCallBLAS("BLASdot",xt   = b[i] - BLASdot_(&bm,v+i,&bm,x,&o));
618:         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
619:       }
620:     }
621:   }
622:   VecRestoreArrayRead(bb,&b);
623:   VecRestoreArray(xx,&x);
624:   return(0);
625: }

627: /* -----------------------------------------------------------------*/
630: PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
631: {
632:   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
633:   const PetscScalar *v   = mat->v,*x;
634:   PetscScalar       *y;
635:   PetscErrorCode    ierr;
636:   PetscBLASInt      m, n,_One=1;
637:   PetscScalar       _DOne=1.0,_DZero=0.0;

640:   PetscBLASIntCast(A->rmap->n,&m);
641:   PetscBLASIntCast(A->cmap->n,&n);
642:   if (!A->rmap->n || !A->cmap->n) return(0);
643:   VecGetArrayRead(xx,&x);
644:   VecGetArray(yy,&y);
645:   PetscStackCallBLAS("BLASgemv",BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One));
646:   VecRestoreArrayRead(xx,&x);
647:   VecRestoreArray(yy,&y);
648:   PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->cmap->n);
649:   return(0);
650: }

654: PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy)
655: {
656:   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
657:   PetscScalar       *y,_DOne=1.0,_DZero=0.0;
658:   PetscErrorCode    ierr;
659:   PetscBLASInt      m, n, _One=1;
660:   const PetscScalar *v = mat->v,*x;

663:   PetscBLASIntCast(A->rmap->n,&m);
664:   PetscBLASIntCast(A->cmap->n,&n);
665:   if (!A->rmap->n || !A->cmap->n) return(0);
666:   VecGetArrayRead(xx,&x);
667:   VecGetArray(yy,&y);
668:   PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One));
669:   VecRestoreArrayRead(xx,&x);
670:   VecRestoreArray(yy,&y);
671:   PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->rmap->n);
672:   return(0);
673: }

677: PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
678: {
679:   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
680:   const PetscScalar *v = mat->v,*x;
681:   PetscScalar       *y,_DOne=1.0;
682:   PetscErrorCode    ierr;
683:   PetscBLASInt      m, n, _One=1;

686:   PetscBLASIntCast(A->rmap->n,&m);
687:   PetscBLASIntCast(A->cmap->n,&n);
688:   if (!A->rmap->n || !A->cmap->n) return(0);
689:   if (zz != yy) {VecCopy(zz,yy);}
690:   VecGetArrayRead(xx,&x);
691:   VecGetArray(yy,&y);
692:   PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One));
693:   VecRestoreArrayRead(xx,&x);
694:   VecRestoreArray(yy,&y);
695:   PetscLogFlops(2.0*A->rmap->n*A->cmap->n);
696:   return(0);
697: }

701: PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
702: {
703:   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
704:   const PetscScalar *v = mat->v,*x;
705:   PetscScalar       *y;
706:   PetscErrorCode    ierr;
707:   PetscBLASInt      m, n, _One=1;
708:   PetscScalar       _DOne=1.0;

711:   PetscBLASIntCast(A->rmap->n,&m);
712:   PetscBLASIntCast(A->cmap->n,&n);
713:   if (!A->rmap->n || !A->cmap->n) return(0);
714:   if (zz != yy) {VecCopy(zz,yy);}
715:   VecGetArrayRead(xx,&x);
716:   VecGetArray(yy,&y);
717:   PetscStackCallBLAS("BLASgemv",BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One));
718:   VecRestoreArrayRead(xx,&x);
719:   VecRestoreArray(yy,&y);
720:   PetscLogFlops(2.0*A->rmap->n*A->cmap->n);
721:   return(0);
722: }

724: /* -----------------------------------------------------------------*/
727: PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
728: {
729:   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
730:   PetscScalar    *v;
732:   PetscInt       i;

735:   *ncols = A->cmap->n;
736:   if (cols) {
737:     PetscMalloc1(A->cmap->n+1,cols);
738:     for (i=0; i<A->cmap->n; i++) (*cols)[i] = i;
739:   }
740:   if (vals) {
741:     PetscMalloc1(A->cmap->n+1,vals);
742:     v    = mat->v + row;
743:     for (i=0; i<A->cmap->n; i++) {(*vals)[i] = *v; v += mat->lda;}
744:   }
745:   return(0);
746: }

750: PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
751: {

755:   if (cols) {PetscFree(*cols);}
756:   if (vals) {PetscFree(*vals); }
757:   return(0);
758: }
759: /* ----------------------------------------------------------------*/
762: PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv)
763: {
764:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
765:   PetscInt     i,j,idx=0;

768:   if (!mat->roworiented) {
769:     if (addv == INSERT_VALUES) {
770:       for (j=0; j<n; j++) {
771:         if (indexn[j] < 0) {idx += m; continue;}
772: #if defined(PETSC_USE_DEBUG)
773:         if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
774: #endif
775:         for (i=0; i<m; i++) {
776:           if (indexm[i] < 0) {idx++; continue;}
777: #if defined(PETSC_USE_DEBUG)
778:           if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
779: #endif
780:           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
781:         }
782:       }
783:     } else {
784:       for (j=0; j<n; j++) {
785:         if (indexn[j] < 0) {idx += m; continue;}
786: #if defined(PETSC_USE_DEBUG)
787:         if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
788: #endif
789:         for (i=0; i<m; i++) {
790:           if (indexm[i] < 0) {idx++; continue;}
791: #if defined(PETSC_USE_DEBUG)
792:           if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
793: #endif
794:           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
795:         }
796:       }
797:     }
798:   } else {
799:     if (addv == INSERT_VALUES) {
800:       for (i=0; i<m; i++) {
801:         if (indexm[i] < 0) { idx += n; continue;}
802: #if defined(PETSC_USE_DEBUG)
803:         if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
804: #endif
805:         for (j=0; j<n; j++) {
806:           if (indexn[j] < 0) { idx++; continue;}
807: #if defined(PETSC_USE_DEBUG)
808:           if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
809: #endif
810:           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
811:         }
812:       }
813:     } else {
814:       for (i=0; i<m; i++) {
815:         if (indexm[i] < 0) { idx += n; continue;}
816: #if defined(PETSC_USE_DEBUG)
817:         if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
818: #endif
819:         for (j=0; j<n; j++) {
820:           if (indexn[j] < 0) { idx++; continue;}
821: #if defined(PETSC_USE_DEBUG)
822:           if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
823: #endif
824:           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
825:         }
826:       }
827:     }
828:   }
829:   return(0);
830: }

834: PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[])
835: {
836:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
837:   PetscInt     i,j;

840:   /* row-oriented output */
841:   for (i=0; i<m; i++) {
842:     if (indexm[i] < 0) {v += n;continue;}
843:     if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %D requested larger than number rows %D",indexm[i],A->rmap->n);
844:     for (j=0; j<n; j++) {
845:       if (indexn[j] < 0) {v++; continue;}
846:       if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column %D requested larger than number columns %D",indexn[j],A->cmap->n);
847:       *v++ = mat->v[indexn[j]*mat->lda + indexm[i]];
848:     }
849:   }
850:   return(0);
851: }

853: /* -----------------------------------------------------------------*/

857: PetscErrorCode MatLoad_SeqDense(Mat newmat,PetscViewer viewer)
858: {
859:   Mat_SeqDense   *a;
861:   PetscInt       *scols,i,j,nz,header[4];
862:   int            fd;
863:   PetscMPIInt    size;
864:   PetscInt       *rowlengths = 0,M,N,*cols,grows,gcols;
865:   PetscScalar    *vals,*svals,*v,*w;
866:   MPI_Comm       comm;

869:   /* force binary viewer to load .info file if it has not yet done so */
870:   PetscViewerSetUp(viewer);
871:   PetscObjectGetComm((PetscObject)viewer,&comm);
872:   MPI_Comm_size(comm,&size);
873:   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"view must have one processor");
874:   PetscViewerBinaryGetDescriptor(viewer,&fd);
875:   PetscBinaryRead(fd,header,4,PETSC_INT);
876:   if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
877:   M = header[1]; N = header[2]; nz = header[3];

879:   /* set global size if not set already*/
880:   if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) {
881:     MatSetSizes(newmat,M,N,M,N);
882:   } else {
883:     /* if sizes and type are already set, check if the vector global sizes are correct */
884:     MatGetSize(newmat,&grows,&gcols);
885:     if (M != grows ||  N != gcols) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Matrix in file of different length (%d, %d) than the input matrix (%d, %d)",M,N,grows,gcols);
886:   }
887:   MatSeqDenseSetPreallocation(newmat,NULL);

889:   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
890:     a = (Mat_SeqDense*)newmat->data;
891:     v = a->v;
892:     /* Allocate some temp space to read in the values and then flip them
893:        from row major to column major */
894:     PetscMalloc1(M*N > 0 ? M*N : 1,&w);
895:     /* read in nonzero values */
896:     PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);
897:     /* now flip the values and store them in the matrix*/
898:     for (j=0; j<N; j++) {
899:       for (i=0; i<M; i++) {
900:         *v++ =w[i*N+j];
901:       }
902:     }
903:     PetscFree(w);
904:   } else {
905:     /* read row lengths */
906:     PetscMalloc1(M+1,&rowlengths);
907:     PetscBinaryRead(fd,rowlengths,M,PETSC_INT);

909:     a = (Mat_SeqDense*)newmat->data;
910:     v = a->v;

912:     /* read column indices and nonzeros */
913:     PetscMalloc1(nz+1,&scols);
914:     cols = scols;
915:     PetscBinaryRead(fd,cols,nz,PETSC_INT);
916:     PetscMalloc1(nz+1,&svals);
917:     vals = svals;
918:     PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);

920:     /* insert into matrix */
921:     for (i=0; i<M; i++) {
922:       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
923:       svals += rowlengths[i]; scols += rowlengths[i];
924:     }
925:     PetscFree(vals);
926:     PetscFree(cols);
927:     PetscFree(rowlengths);
928:   }
929:   MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
930:   MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);
931:   return(0);
932: }

936: static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
937: {
938:   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
939:   PetscErrorCode    ierr;
940:   PetscInt          i,j;
941:   const char        *name;
942:   PetscScalar       *v;
943:   PetscViewerFormat format;
944: #if defined(PETSC_USE_COMPLEX)
945:   PetscBool         allreal = PETSC_TRUE;
946: #endif

949:   PetscViewerGetFormat(viewer,&format);
950:   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
951:     return(0);  /* do nothing for now */
952:   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
953:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
954:     for (i=0; i<A->rmap->n; i++) {
955:       v    = a->v + i;
956:       PetscViewerASCIIPrintf(viewer,"row %D:",i);
957:       for (j=0; j<A->cmap->n; j++) {
958: #if defined(PETSC_USE_COMPLEX)
959:         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
960:           PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",j,(double)PetscRealPart(*v),(double)PetscImaginaryPart(*v));
961:         } else if (PetscRealPart(*v)) {
962:           PetscViewerASCIIPrintf(viewer," (%D, %g) ",j,(double)PetscRealPart(*v));
963:         }
964: #else
965:         if (*v) {
966:           PetscViewerASCIIPrintf(viewer," (%D, %g) ",j,(double)*v);
967:         }
968: #endif
969:         v += a->lda;
970:       }
971:       PetscViewerASCIIPrintf(viewer,"\n");
972:     }
973:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
974:   } else {
975:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
976: #if defined(PETSC_USE_COMPLEX)
977:     /* determine if matrix has all real values */
978:     v = a->v;
979:     for (i=0; i<A->rmap->n*A->cmap->n; i++) {
980:       if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break;}
981:     }
982: #endif
983:     if (format == PETSC_VIEWER_ASCII_MATLAB) {
984:       PetscObjectGetName((PetscObject)A,&name);
985:       PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap->n,A->cmap->n);
986:       PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap->n,A->cmap->n);
987:       PetscViewerASCIIPrintf(viewer,"%s = [\n",name);
988:     }

990:     for (i=0; i<A->rmap->n; i++) {
991:       v = a->v + i;
992:       for (j=0; j<A->cmap->n; j++) {
993: #if defined(PETSC_USE_COMPLEX)
994:         if (allreal) {
995:           PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)PetscRealPart(*v));
996:         } else {
997:           PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16ei ",(double)PetscRealPart(*v),(double)PetscImaginaryPart(*v));
998:         }
999: #else
1000:         PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)*v);
1001: #endif
1002:         v += a->lda;
1003:       }
1004:       PetscViewerASCIIPrintf(viewer,"\n");
1005:     }
1006:     if (format == PETSC_VIEWER_ASCII_MATLAB) {
1007:       PetscViewerASCIIPrintf(viewer,"];\n");
1008:     }
1009:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
1010:   }
1011:   PetscViewerFlush(viewer);
1012:   return(0);
1013: }

1017: static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
1018: {
1019:   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
1020:   PetscErrorCode    ierr;
1021:   int               fd;
1022:   PetscInt          ict,j,n = A->cmap->n,m = A->rmap->n,i,*col_lens,nz = m*n;
1023:   PetscScalar       *v,*anonz,*vals;
1024:   PetscViewerFormat format;

1027:   PetscViewerBinaryGetDescriptor(viewer,&fd);

1029:   PetscViewerGetFormat(viewer,&format);
1030:   if (format == PETSC_VIEWER_NATIVE) {
1031:     /* store the matrix as a dense matrix */
1032:     PetscMalloc1(4,&col_lens);

1034:     col_lens[0] = MAT_FILE_CLASSID;
1035:     col_lens[1] = m;
1036:     col_lens[2] = n;
1037:     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;

1039:     PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);
1040:     PetscFree(col_lens);

1042:     /* write out matrix, by rows */
1043:     PetscMalloc1(m*n+1,&vals);
1044:     v    = a->v;
1045:     for (j=0; j<n; j++) {
1046:       for (i=0; i<m; i++) {
1047:         vals[j + i*n] = *v++;
1048:       }
1049:     }
1050:     PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);
1051:     PetscFree(vals);
1052:   } else {
1053:     PetscMalloc1(4+nz,&col_lens);

1055:     col_lens[0] = MAT_FILE_CLASSID;
1056:     col_lens[1] = m;
1057:     col_lens[2] = n;
1058:     col_lens[3] = nz;

1060:     /* store lengths of each row and write (including header) to file */
1061:     for (i=0; i<m; i++) col_lens[4+i] = n;
1062:     PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);

1064:     /* Possibly should write in smaller increments, not whole matrix at once? */
1065:     /* store column indices (zero start index) */
1066:     ict = 0;
1067:     for (i=0; i<m; i++) {
1068:       for (j=0; j<n; j++) col_lens[ict++] = j;
1069:     }
1070:     PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);
1071:     PetscFree(col_lens);

1073:     /* store nonzero values */
1074:     PetscMalloc1(nz+1,&anonz);
1075:     ict  = 0;
1076:     for (i=0; i<m; i++) {
1077:       v = a->v + i;
1078:       for (j=0; j<n; j++) {
1079:         anonz[ict++] = *v; v += a->lda;
1080:       }
1081:     }
1082:     PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);
1083:     PetscFree(anonz);
1084:   }
1085:   return(0);
1086: }

1088: #include <petscdraw.h>
1091: PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
1092: {
1093:   Mat               A  = (Mat) Aa;
1094:   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
1095:   PetscErrorCode    ierr;
1096:   PetscInt          m  = A->rmap->n,n = A->cmap->n,color,i,j;
1097:   PetscScalar       *v = a->v;
1098:   PetscViewer       viewer;
1099:   PetscDraw         popup;
1100:   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0;
1101:   PetscViewerFormat format;

1104:   PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);
1105:   PetscViewerGetFormat(viewer,&format);
1106:   PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);

1108:   /* Loop over matrix elements drawing boxes */
1109:   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
1110:     /* Blue for negative and Red for positive */
1111:     color = PETSC_DRAW_BLUE;
1112:     for (j = 0; j < n; j++) {
1113:       x_l = j;
1114:       x_r = x_l + 1.0;
1115:       for (i = 0; i < m; i++) {
1116:         y_l = m - i - 1.0;
1117:         y_r = y_l + 1.0;
1118:         if (PetscRealPart(v[j*m+i]) >  0.) {
1119:           color = PETSC_DRAW_RED;
1120:         } else if (PetscRealPart(v[j*m+i]) <  0.) {
1121:           color = PETSC_DRAW_BLUE;
1122:         } else {
1123:           continue;
1124:         }
1125:         PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
1126:       }
1127:     }
1128:   } else {
1129:     /* use contour shading to indicate magnitude of values */
1130:     /* first determine max of all nonzero values */
1131:     for (i = 0; i < m*n; i++) {
1132:       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
1133:     }
1134:     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
1135:     PetscDrawGetPopup(draw,&popup);
1136:     if (popup) {PetscDrawScalePopup(popup,0.0,maxv);}
1137:     for (j = 0; j < n; j++) {
1138:       x_l = j;
1139:       x_r = x_l + 1.0;
1140:       for (i = 0; i < m; i++) {
1141:         y_l   = m - i - 1.0;
1142:         y_r   = y_l + 1.0;
1143:         color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i]));
1144:         PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
1145:       }
1146:     }
1147:   }
1148:   return(0);
1149: }

1153: PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
1154: {
1155:   PetscDraw      draw;
1156:   PetscBool      isnull;
1157:   PetscReal      xr,yr,xl,yl,h,w;

1161:   PetscViewerDrawGetDraw(viewer,0,&draw);
1162:   PetscDrawIsNull(draw,&isnull);
1163:   if (isnull) return(0);

1165:   PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
1166:   xr   = A->cmap->n; yr = A->rmap->n; h = yr/10.0; w = xr/10.0;
1167:   xr  += w;    yr += h;  xl = -w;     yl = -h;
1168:   PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
1169:   PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);
1170:   PetscObjectCompose((PetscObject)A,"Zoomviewer",NULL);
1171:   return(0);
1172: }

1176: PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer)
1177: {
1179:   PetscBool      iascii,isbinary,isdraw;

1182:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1183:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
1184:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);

1186:   if (iascii) {
1187:     MatView_SeqDense_ASCII(A,viewer);
1188:   } else if (isbinary) {
1189:     MatView_SeqDense_Binary(A,viewer);
1190:   } else if (isdraw) {
1191:     MatView_SeqDense_Draw(A,viewer);
1192:   }
1193:   return(0);
1194: }

1198: PetscErrorCode MatDestroy_SeqDense(Mat mat)
1199: {
1200:   Mat_SeqDense   *l = (Mat_SeqDense*)mat->data;

1204: #if defined(PETSC_USE_LOG)
1205:   PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap->n,mat->cmap->n);
1206: #endif
1207:   PetscFree(l->pivots);
1208:   if (!l->user_alloc) {PetscFree(l->v);}
1209:   PetscFree(mat->data);

1211:   PetscObjectChangeTypeName((PetscObject)mat,0);
1212:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArray_C",NULL);
1213:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArray_C",NULL);
1214:   PetscObjectComposeFunction((PetscObject)mat,"MatConvert_seqdense_seqaij_C",NULL);
1215: #if defined(PETSC_HAVE_ELEMENTAL)
1216:   PetscObjectComposeFunction((PetscObject)mat,"MatConvert_seqdense_elemental_C",NULL);
1217: #endif
1218:   PetscObjectComposeFunction((PetscObject)mat,"MatSeqDenseSetPreallocation_C",NULL);
1219:   PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_seqaij_seqdense_C",NULL);
1220:   PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C",NULL);
1221:   PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C",NULL);
1222:   PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMult_seqaij_seqdense_C",NULL);
1223:   PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultSymbolic_seqaij_seqdense_C",NULL);
1224:   PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultNumeric_seqaij_seqdense_C",NULL);
1225:   return(0);
1226: }

1230: PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout)
1231: {
1232:   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1234:   PetscInt       k,j,m,n,M;
1235:   PetscScalar    *v,tmp;

1238:   v = mat->v; m = A->rmap->n; M = mat->lda; n = A->cmap->n;
1239:   if (reuse == MAT_REUSE_MATRIX && *matout == A) { /* in place transpose */
1240:     if (m != n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can not transpose non-square matrix in place");
1241:     else {
1242:       for (j=0; j<m; j++) {
1243:         for (k=0; k<j; k++) {
1244:           tmp        = v[j + k*M];
1245:           v[j + k*M] = v[k + j*M];
1246:           v[k + j*M] = tmp;
1247:         }
1248:       }
1249:     }
1250:   } else { /* out-of-place transpose */
1251:     Mat          tmat;
1252:     Mat_SeqDense *tmatd;
1253:     PetscScalar  *v2;

1255:     if (reuse == MAT_INITIAL_MATRIX) {
1256:       MatCreate(PetscObjectComm((PetscObject)A),&tmat);
1257:       MatSetSizes(tmat,A->cmap->n,A->rmap->n,A->cmap->n,A->rmap->n);
1258:       MatSetType(tmat,((PetscObject)A)->type_name);
1259:       MatSeqDenseSetPreallocation(tmat,NULL);
1260:     } else {
1261:       tmat = *matout;
1262:     }
1263:     tmatd = (Mat_SeqDense*)tmat->data;
1264:     v = mat->v; v2 = tmatd->v;
1265:     for (j=0; j<n; j++) {
1266:       for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M];
1267:     }
1268:     MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);
1269:     MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);
1270:     *matout = tmat;
1271:   }
1272:   return(0);
1273: }

1277: PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscBool  *flg)
1278: {
1279:   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1280:   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
1281:   PetscInt     i,j;
1282:   PetscScalar  *v1,*v2;

1285:   if (A1->rmap->n != A2->rmap->n) {*flg = PETSC_FALSE; return(0);}
1286:   if (A1->cmap->n != A2->cmap->n) {*flg = PETSC_FALSE; return(0);}
1287:   for (i=0; i<A1->rmap->n; i++) {
1288:     v1 = mat1->v+i; v2 = mat2->v+i;
1289:     for (j=0; j<A1->cmap->n; j++) {
1290:       if (*v1 != *v2) {*flg = PETSC_FALSE; return(0);}
1291:       v1 += mat1->lda; v2 += mat2->lda;
1292:     }
1293:   }
1294:   *flg = PETSC_TRUE;
1295:   return(0);
1296: }

1300: PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v)
1301: {
1302:   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1304:   PetscInt       i,n,len;
1305:   PetscScalar    *x,zero = 0.0;

1308:   VecSet(v,zero);
1309:   VecGetSize(v,&n);
1310:   VecGetArray(v,&x);
1311:   len  = PetscMin(A->rmap->n,A->cmap->n);
1312:   if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
1313:   for (i=0; i<len; i++) {
1314:     x[i] = mat->v[i*mat->lda + i];
1315:   }
1316:   VecRestoreArray(v,&x);
1317:   return(0);
1318: }

1322: PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1323: {
1324:   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
1325:   const PetscScalar *l,*r;
1326:   PetscScalar       x,*v;
1327:   PetscErrorCode    ierr;
1328:   PetscInt          i,j,m = A->rmap->n,n = A->cmap->n;

1331:   if (ll) {
1332:     VecGetSize(ll,&m);
1333:     VecGetArrayRead(ll,&l);
1334:     if (m != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1335:     for (i=0; i<m; i++) {
1336:       x = l[i];
1337:       v = mat->v + i;
1338:       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1339:     }
1340:     VecRestoreArrayRead(ll,&l);
1341:     PetscLogFlops(1.0*n*m);
1342:   }
1343:   if (rr) {
1344:     VecGetSize(rr,&n);
1345:     VecGetArrayRead(rr,&r);
1346:     if (n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1347:     for (i=0; i<n; i++) {
1348:       x = r[i];
1349:       v = mat->v + i*m;
1350:       for (j=0; j<m; j++) (*v++) *= x;
1351:     }
1352:     VecRestoreArrayRead(rr,&r);
1353:     PetscLogFlops(1.0*n*m);
1354:   }
1355:   return(0);
1356: }

1360: PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1361: {
1362:   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1363:   PetscScalar    *v   = mat->v;
1364:   PetscReal      sum  = 0.0;
1365:   PetscInt       lda  =mat->lda,m=A->rmap->n,i,j;

1369:   if (type == NORM_FROBENIUS) {
1370:     if (lda>m) {
1371:       for (j=0; j<A->cmap->n; j++) {
1372:         v = mat->v+j*lda;
1373:         for (i=0; i<m; i++) {
1374:           sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1375:         }
1376:       }
1377:     } else {
1378:       for (i=0; i<A->cmap->n*A->rmap->n; i++) {
1379:         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1380:       }
1381:     }
1382:     *nrm = PetscSqrtReal(sum);
1383:     PetscLogFlops(2.0*A->cmap->n*A->rmap->n);
1384:   } else if (type == NORM_1) {
1385:     *nrm = 0.0;
1386:     for (j=0; j<A->cmap->n; j++) {
1387:       v   = mat->v + j*mat->lda;
1388:       sum = 0.0;
1389:       for (i=0; i<A->rmap->n; i++) {
1390:         sum += PetscAbsScalar(*v);  v++;
1391:       }
1392:       if (sum > *nrm) *nrm = sum;
1393:     }
1394:     PetscLogFlops(1.0*A->cmap->n*A->rmap->n);
1395:   } else if (type == NORM_INFINITY) {
1396:     *nrm = 0.0;
1397:     for (j=0; j<A->rmap->n; j++) {
1398:       v   = mat->v + j;
1399:       sum = 0.0;
1400:       for (i=0; i<A->cmap->n; i++) {
1401:         sum += PetscAbsScalar(*v); v += mat->lda;
1402:       }
1403:       if (sum > *nrm) *nrm = sum;
1404:     }
1405:     PetscLogFlops(1.0*A->cmap->n*A->rmap->n);
1406:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No two norm");
1407:   return(0);
1408: }

1412: PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscBool flg)
1413: {
1414:   Mat_SeqDense   *aij = (Mat_SeqDense*)A->data;

1418:   switch (op) {
1419:   case MAT_ROW_ORIENTED:
1420:     aij->roworiented = flg;
1421:     break;
1422:   case MAT_NEW_NONZERO_LOCATIONS:
1423:   case MAT_NEW_NONZERO_LOCATION_ERR:
1424:   case MAT_NEW_NONZERO_ALLOCATION_ERR:
1425:   case MAT_NEW_DIAGONALS:
1426:   case MAT_KEEP_NONZERO_PATTERN:
1427:   case MAT_IGNORE_OFF_PROC_ENTRIES:
1428:   case MAT_USE_HASH_TABLE:
1429:   case MAT_IGNORE_LOWER_TRIANGULAR:
1430:     PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1431:     break;
1432:   case MAT_SPD:
1433:   case MAT_SYMMETRIC:
1434:   case MAT_STRUCTURALLY_SYMMETRIC:
1435:   case MAT_HERMITIAN:
1436:   case MAT_SYMMETRY_ETERNAL:
1437:     /* These options are handled directly by MatSetOption() */
1438:     break;
1439:   default:
1440:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
1441:   }
1442:   return(0);
1443: }

1447: PetscErrorCode MatZeroEntries_SeqDense(Mat A)
1448: {
1449:   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
1451:   PetscInt       lda=l->lda,m=A->rmap->n,j;

1454:   if (lda>m) {
1455:     for (j=0; j<A->cmap->n; j++) {
1456:       PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));
1457:     }
1458:   } else {
1459:     PetscMemzero(l->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));
1460:   }
1461:   return(0);
1462: }

1466: PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
1467: {
1468:   PetscErrorCode    ierr;
1469:   Mat_SeqDense      *l = (Mat_SeqDense*)A->data;
1470:   PetscInt          m  = l->lda, n = A->cmap->n, i,j;
1471:   PetscScalar       *slot,*bb;
1472:   const PetscScalar *xx;

1475: #if defined(PETSC_USE_DEBUG)
1476:   for (i=0; i<N; i++) {
1477:     if (rows[i] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row requested to be zeroed");
1478:     if (rows[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %D requested to be zeroed greater than or equal number of rows %D",rows[i],A->rmap->n);
1479:   }
1480: #endif

1482:   /* fix right hand side if needed */
1483:   if (x && b) {
1484:     VecGetArrayRead(x,&xx);
1485:     VecGetArray(b,&bb);
1486:     for (i=0; i<N; i++) bb[rows[i]] = diag*xx[rows[i]];
1487:     VecRestoreArrayRead(x,&xx);
1488:     VecRestoreArray(b,&bb);
1489:   }

1491:   for (i=0; i<N; i++) {
1492:     slot = l->v + rows[i];
1493:     for (j=0; j<n; j++) { *slot = 0.0; slot += m;}
1494:   }
1495:   if (diag != 0.0) {
1496:     if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only coded for square matrices");
1497:     for (i=0; i<N; i++) {
1498:       slot  = l->v + (m+1)*rows[i];
1499:       *slot = diag;
1500:     }
1501:   }
1502:   return(0);
1503: }

1507: PetscErrorCode MatDenseGetArray_SeqDense(Mat A,PetscScalar *array[])
1508: {
1509:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;

1512:   if (mat->lda != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot get array for Dense matrices with LDA different from number of rows");
1513:   *array = mat->v;
1514:   return(0);
1515: }

1519: PetscErrorCode MatDenseRestoreArray_SeqDense(Mat A,PetscScalar *array[])
1520: {
1522:   *array = 0; /* user cannot accidently use the array later */
1523:   return(0);
1524: }

1528: /*@C
1529:    MatDenseGetArray - gives access to the array where the data for a SeqDense matrix is stored

1531:    Not Collective

1533:    Input Parameter:
1534: .  mat - a MATSEQDENSE or MATMPIDENSE matrix

1536:    Output Parameter:
1537: .   array - pointer to the data

1539:    Level: intermediate

1541: .seealso: MatDenseRestoreArray()
1542: @*/
1543: PetscErrorCode  MatDenseGetArray(Mat A,PetscScalar **array)
1544: {

1548:   PetscUseMethod(A,"MatDenseGetArray_C",(Mat,PetscScalar**),(A,array));
1549:   return(0);
1550: }

1554: /*@C
1555:    MatDenseRestoreArray - returns access to the array where the data for a dense matrix is stored obtained by MatDenseGetArray()

1557:    Not Collective

1559:    Input Parameters:
1560: .  mat - a MATSEQDENSE or MATMPIDENSE matrix
1561: .  array - pointer to the data

1563:    Level: intermediate

1565: .seealso: MatDenseGetArray()
1566: @*/
1567: PetscErrorCode  MatDenseRestoreArray(Mat A,PetscScalar **array)
1568: {

1572:   PetscUseMethod(A,"MatDenseRestoreArray_C",(Mat,PetscScalar**),(A,array));
1573:   return(0);
1574: }

1578: static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
1579: {
1580:   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1582:   PetscInt       i,j,nrows,ncols;
1583:   const PetscInt *irow,*icol;
1584:   PetscScalar    *av,*bv,*v = mat->v;
1585:   Mat            newmat;

1588:   ISGetIndices(isrow,&irow);
1589:   ISGetIndices(iscol,&icol);
1590:   ISGetLocalSize(isrow,&nrows);
1591:   ISGetLocalSize(iscol,&ncols);

1593:   /* Check submatrixcall */
1594:   if (scall == MAT_REUSE_MATRIX) {
1595:     PetscInt n_cols,n_rows;
1596:     MatGetSize(*B,&n_rows,&n_cols);
1597:     if (n_rows != nrows || n_cols != ncols) {
1598:       /* resize the result matrix to match number of requested rows/columns */
1599:       MatSetSizes(*B,nrows,ncols,nrows,ncols);
1600:     }
1601:     newmat = *B;
1602:   } else {
1603:     /* Create and fill new matrix */
1604:     MatCreate(PetscObjectComm((PetscObject)A),&newmat);
1605:     MatSetSizes(newmat,nrows,ncols,nrows,ncols);
1606:     MatSetType(newmat,((PetscObject)A)->type_name);
1607:     MatSeqDenseSetPreallocation(newmat,NULL);
1608:   }

1610:   /* Now extract the data pointers and do the copy,column at a time */
1611:   bv = ((Mat_SeqDense*)newmat->data)->v;

1613:   for (i=0; i<ncols; i++) {
1614:     av = v + mat->lda*icol[i];
1615:     for (j=0; j<nrows; j++) *bv++ = av[irow[j]];
1616:   }

1618:   /* Assemble the matrices so that the correct flags are set */
1619:   MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
1620:   MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);

1622:   /* Free work space */
1623:   ISRestoreIndices(isrow,&irow);
1624:   ISRestoreIndices(iscol,&icol);
1625:   *B   = newmat;
1626:   return(0);
1627: }

1631: PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1632: {
1634:   PetscInt       i;

1637:   if (scall == MAT_INITIAL_MATRIX) {
1638:     PetscMalloc1(n+1,B);
1639:   }

1641:   for (i=0; i<n; i++) {
1642:     MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);
1643:   }
1644:   return(0);
1645: }

1649: PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode)
1650: {
1652:   return(0);
1653: }

1657: PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode)
1658: {
1660:   return(0);
1661: }

1665: PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
1666: {
1667:   Mat_SeqDense   *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense*)B->data;
1669:   PetscInt       lda1=a->lda,lda2=b->lda, m=A->rmap->n,n=A->cmap->n, j;

1672:   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
1673:   if (A->ops->copy != B->ops->copy) {
1674:     MatCopy_Basic(A,B,str);
1675:     return(0);
1676:   }
1677:   if (m != B->rmap->n || n != B->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1678:   if (lda1>m || lda2>m) {
1679:     for (j=0; j<n; j++) {
1680:       PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));
1681:     }
1682:   } else {
1683:     PetscMemcpy(b->v,a->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));
1684:   }
1685:   return(0);
1686: }

1690: PetscErrorCode MatSetUp_SeqDense(Mat A)
1691: {

1695:    MatSeqDenseSetPreallocation(A,0);
1696:   return(0);
1697: }

1701: static PetscErrorCode MatConjugate_SeqDense(Mat A)
1702: {
1703:   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
1704:   PetscInt     i,nz = A->rmap->n*A->cmap->n;
1705:   PetscScalar  *aa = a->v;

1708:   for (i=0; i<nz; i++) aa[i] = PetscConj(aa[i]);
1709:   return(0);
1710: }

1714: static PetscErrorCode MatRealPart_SeqDense(Mat A)
1715: {
1716:   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
1717:   PetscInt     i,nz = A->rmap->n*A->cmap->n;
1718:   PetscScalar  *aa = a->v;

1721:   for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
1722:   return(0);
1723: }

1727: static PetscErrorCode MatImaginaryPart_SeqDense(Mat A)
1728: {
1729:   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
1730:   PetscInt     i,nz = A->rmap->n*A->cmap->n;
1731:   PetscScalar  *aa = a->v;

1734:   for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
1735:   return(0);
1736: }

1738: /* ----------------------------------------------------------------*/
1741: PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1742: {

1746:   if (scall == MAT_INITIAL_MATRIX) {
1747:     PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);
1748:     MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);
1749:     PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);
1750:   }
1751:   PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);
1752:   MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);
1753:   PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);
1754:   return(0);
1755: }

1759: PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1760: {
1762:   PetscInt       m=A->rmap->n,n=B->cmap->n;
1763:   Mat            Cmat;

1766:   if (A->cmap->n != B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"A->cmap->n %d != B->rmap->n %d\n",A->cmap->n,B->rmap->n);
1767:   MatCreate(PETSC_COMM_SELF,&Cmat);
1768:   MatSetSizes(Cmat,m,n,m,n);
1769:   MatSetType(Cmat,MATSEQDENSE);
1770:   MatSeqDenseSetPreallocation(Cmat,NULL);

1772:   *C = Cmat;
1773:   return(0);
1774: }

1778: PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1779: {
1780:   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1781:   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1782:   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
1783:   PetscBLASInt   m,n,k;
1784:   PetscScalar    _DOne=1.0,_DZero=0.0;

1788:   PetscBLASIntCast(A->rmap->n,&m);
1789:   PetscBLASIntCast(B->cmap->n,&n);
1790:   PetscBLASIntCast(A->cmap->n,&k);
1791:   PetscStackCallBLAS("BLASgemv",BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda));
1792:   return(0);
1793: }

1797: PetscErrorCode MatTransposeMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1798: {

1802:   if (scall == MAT_INITIAL_MATRIX) {
1803:     PetscLogEventBegin(MAT_TransposeMatMultSymbolic,A,B,0,0);
1804:     MatTransposeMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);
1805:     PetscLogEventEnd(MAT_TransposeMatMultSymbolic,A,B,0,0);
1806:   }
1807:   PetscLogEventBegin(MAT_TransposeMatMultNumeric,A,B,0,0);
1808:   MatTransposeMatMultNumeric_SeqDense_SeqDense(A,B,*C);
1809:   PetscLogEventEnd(MAT_TransposeMatMultNumeric,A,B,0,0);
1810:   return(0);
1811: }

1815: PetscErrorCode MatTransposeMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1816: {
1818:   PetscInt       m=A->cmap->n,n=B->cmap->n;
1819:   Mat            Cmat;

1822:   if (A->rmap->n != B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"A->rmap->n %d != B->rmap->n %d\n",A->rmap->n,B->rmap->n);
1823:   MatCreate(PETSC_COMM_SELF,&Cmat);
1824:   MatSetSizes(Cmat,m,n,m,n);
1825:   MatSetType(Cmat,MATSEQDENSE);
1826:   MatSeqDenseSetPreallocation(Cmat,NULL);

1828:   Cmat->assembled = PETSC_TRUE;

1830:   *C = Cmat;
1831:   return(0);
1832: }

1836: PetscErrorCode MatTransposeMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1837: {
1838:   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1839:   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1840:   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
1841:   PetscBLASInt   m,n,k;
1842:   PetscScalar    _DOne=1.0,_DZero=0.0;

1846:   PetscBLASIntCast(A->cmap->n,&m);
1847:   PetscBLASIntCast(B->cmap->n,&n);
1848:   PetscBLASIntCast(A->rmap->n,&k);
1849:   /*
1850:      Note the m and n arguments below are the number rows and columns of A', not A!
1851:   */
1852:   PetscStackCallBLAS("BLASgemv",BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda));
1853:   return(0);
1854: }

1858: PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[])
1859: {
1860:   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1862:   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1863:   PetscScalar    *x;
1864:   MatScalar      *aa = a->v;

1867:   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");

1869:   VecSet(v,0.0);
1870:   VecGetArray(v,&x);
1871:   VecGetLocalSize(v,&p);
1872:   if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1873:   for (i=0; i<m; i++) {
1874:     x[i] = aa[i]; if (idx) idx[i] = 0;
1875:     for (j=1; j<n; j++) {
1876:       if (PetscRealPart(x[i]) < PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
1877:     }
1878:   }
1879:   VecRestoreArray(v,&x);
1880:   return(0);
1881: }

1885: PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[])
1886: {
1887:   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1889:   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1890:   PetscScalar    *x;
1891:   PetscReal      atmp;
1892:   MatScalar      *aa = a->v;

1895:   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");

1897:   VecSet(v,0.0);
1898:   VecGetArray(v,&x);
1899:   VecGetLocalSize(v,&p);
1900:   if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1901:   for (i=0; i<m; i++) {
1902:     x[i] = PetscAbsScalar(aa[i]);
1903:     for (j=1; j<n; j++) {
1904:       atmp = PetscAbsScalar(aa[i+m*j]);
1905:       if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;}
1906:     }
1907:   }
1908:   VecRestoreArray(v,&x);
1909:   return(0);
1910: }

1914: PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[])
1915: {
1916:   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1918:   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1919:   PetscScalar    *x;
1920:   MatScalar      *aa = a->v;

1923:   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");

1925:   VecSet(v,0.0);
1926:   VecGetArray(v,&x);
1927:   VecGetLocalSize(v,&p);
1928:   if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1929:   for (i=0; i<m; i++) {
1930:     x[i] = aa[i]; if (idx) idx[i] = 0;
1931:     for (j=1; j<n; j++) {
1932:       if (PetscRealPart(x[i]) > PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
1933:     }
1934:   }
1935:   VecRestoreArray(v,&x);
1936:   return(0);
1937: }

1941: PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col)
1942: {
1943:   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1945:   PetscScalar    *x;

1948:   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");

1950:   VecGetArray(v,&x);
1951:   PetscMemcpy(x,a->v+col*a->lda,A->rmap->n*sizeof(PetscScalar));
1952:   VecRestoreArray(v,&x);
1953:   return(0);
1954: }


1959: PetscErrorCode MatGetColumnNorms_SeqDense(Mat A,NormType type,PetscReal *norms)
1960: {
1962:   PetscInt       i,j,m,n;
1963:   PetscScalar    *a;

1966:   MatGetSize(A,&m,&n);
1967:   PetscMemzero(norms,n*sizeof(PetscReal));
1968:   MatDenseGetArray(A,&a);
1969:   if (type == NORM_2) {
1970:     for (i=0; i<n; i++) {
1971:       for (j=0; j<m; j++) {
1972:         norms[i] += PetscAbsScalar(a[j]*a[j]);
1973:       }
1974:       a += m;
1975:     }
1976:   } else if (type == NORM_1) {
1977:     for (i=0; i<n; i++) {
1978:       for (j=0; j<m; j++) {
1979:         norms[i] += PetscAbsScalar(a[j]);
1980:       }
1981:       a += m;
1982:     }
1983:   } else if (type == NORM_INFINITY) {
1984:     for (i=0; i<n; i++) {
1985:       for (j=0; j<m; j++) {
1986:         norms[i] = PetscMax(PetscAbsScalar(a[j]),norms[i]);
1987:       }
1988:       a += m;
1989:     }
1990:   } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Unknown NormType");
1991:   MatDenseRestoreArray(A,&a);
1992:   if (type == NORM_2) {
1993:     for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]);
1994:   }
1995:   return(0);
1996: }

2000: static PetscErrorCode  MatSetRandom_SeqDense(Mat x,PetscRandom rctx)
2001: {
2003:   PetscScalar    *a;
2004:   PetscInt       m,n,i;

2007:   MatGetSize(x,&m,&n);
2008:   MatDenseGetArray(x,&a);
2009:   for (i=0; i<m*n; i++) {
2010:     PetscRandomGetValue(rctx,a+i);
2011:   }
2012:   MatDenseRestoreArray(x,&a);
2013:   return(0);
2014: }


2017: /* -------------------------------------------------------------------*/
2018: static struct _MatOps MatOps_Values = { MatSetValues_SeqDense,
2019:                                         MatGetRow_SeqDense,
2020:                                         MatRestoreRow_SeqDense,
2021:                                         MatMult_SeqDense,
2022:                                 /*  4*/ MatMultAdd_SeqDense,
2023:                                         MatMultTranspose_SeqDense,
2024:                                         MatMultTransposeAdd_SeqDense,
2025:                                         0,
2026:                                         0,
2027:                                         0,
2028:                                 /* 10*/ 0,
2029:                                         MatLUFactor_SeqDense,
2030:                                         MatCholeskyFactor_SeqDense,
2031:                                         MatSOR_SeqDense,
2032:                                         MatTranspose_SeqDense,
2033:                                 /* 15*/ MatGetInfo_SeqDense,
2034:                                         MatEqual_SeqDense,
2035:                                         MatGetDiagonal_SeqDense,
2036:                                         MatDiagonalScale_SeqDense,
2037:                                         MatNorm_SeqDense,
2038:                                 /* 20*/ MatAssemblyBegin_SeqDense,
2039:                                         MatAssemblyEnd_SeqDense,
2040:                                         MatSetOption_SeqDense,
2041:                                         MatZeroEntries_SeqDense,
2042:                                 /* 24*/ MatZeroRows_SeqDense,
2043:                                         0,
2044:                                         0,
2045:                                         0,
2046:                                         0,
2047:                                 /* 29*/ MatSetUp_SeqDense,
2048:                                         0,
2049:                                         0,
2050:                                         0,
2051:                                         0,
2052:                                 /* 34*/ MatDuplicate_SeqDense,
2053:                                         0,
2054:                                         0,
2055:                                         0,
2056:                                         0,
2057:                                 /* 39*/ MatAXPY_SeqDense,
2058:                                         MatGetSubMatrices_SeqDense,
2059:                                         0,
2060:                                         MatGetValues_SeqDense,
2061:                                         MatCopy_SeqDense,
2062:                                 /* 44*/ MatGetRowMax_SeqDense,
2063:                                         MatScale_SeqDense,
2064:                                         0,
2065:                                         0,
2066:                                         0,
2067:                                 /* 49*/ MatSetRandom_SeqDense,
2068:                                         0,
2069:                                         0,
2070:                                         0,
2071:                                         0,
2072:                                 /* 54*/ 0,
2073:                                         0,
2074:                                         0,
2075:                                         0,
2076:                                         0,
2077:                                 /* 59*/ 0,
2078:                                         MatDestroy_SeqDense,
2079:                                         MatView_SeqDense,
2080:                                         0,
2081:                                         0,
2082:                                 /* 64*/ 0,
2083:                                         0,
2084:                                         0,
2085:                                         0,
2086:                                         0,
2087:                                 /* 69*/ MatGetRowMaxAbs_SeqDense,
2088:                                         0,
2089:                                         0,
2090:                                         0,
2091:                                         0,
2092:                                 /* 74*/ 0,
2093:                                         0,
2094:                                         0,
2095:                                         0,
2096:                                         0,
2097:                                 /* 79*/ 0,
2098:                                         0,
2099:                                         0,
2100:                                         0,
2101:                                 /* 83*/ MatLoad_SeqDense,
2102:                                         0,
2103:                                         MatIsHermitian_SeqDense,
2104:                                         0,
2105:                                         0,
2106:                                         0,
2107:                                 /* 89*/ MatMatMult_SeqDense_SeqDense,
2108:                                         MatMatMultSymbolic_SeqDense_SeqDense,
2109:                                         MatMatMultNumeric_SeqDense_SeqDense,
2110:                                         0,
2111:                                         0,
2112:                                 /* 94*/ 0,
2113:                                         0,
2114:                                         0,
2115:                                         0,
2116:                                         0,
2117:                                 /* 99*/ 0,
2118:                                         0,
2119:                                         0,
2120:                                         MatConjugate_SeqDense,
2121:                                         0,
2122:                                 /*104*/ 0,
2123:                                         MatRealPart_SeqDense,
2124:                                         MatImaginaryPart_SeqDense,
2125:                                         0,
2126:                                         0,
2127:                                 /*109*/ MatMatSolve_SeqDense,
2128:                                         0,
2129:                                         MatGetRowMin_SeqDense,
2130:                                         MatGetColumnVector_SeqDense,
2131:                                         0,
2132:                                 /*114*/ 0,
2133:                                         0,
2134:                                         0,
2135:                                         0,
2136:                                         0,
2137:                                 /*119*/ 0,
2138:                                         0,
2139:                                         0,
2140:                                         0,
2141:                                         0,
2142:                                 /*124*/ 0,
2143:                                         MatGetColumnNorms_SeqDense,
2144:                                         0,
2145:                                         0,
2146:                                         0,
2147:                                 /*129*/ 0,
2148:                                         MatTransposeMatMult_SeqDense_SeqDense,
2149:                                         MatTransposeMatMultSymbolic_SeqDense_SeqDense,
2150:                                         MatTransposeMatMultNumeric_SeqDense_SeqDense,
2151:                                         0,
2152:                                 /*134*/ 0,
2153:                                         0,
2154:                                         0,
2155:                                         0,
2156:                                         0,
2157:                                 /*139*/ 0,
2158:                                         0,
2159:                                         0
2160: };

2164: /*@C
2165:    MatCreateSeqDense - Creates a sequential dense matrix that
2166:    is stored in column major order (the usual Fortran 77 manner). Many
2167:    of the matrix operations use the BLAS and LAPACK routines.

2169:    Collective on MPI_Comm

2171:    Input Parameters:
2172: +  comm - MPI communicator, set to PETSC_COMM_SELF
2173: .  m - number of rows
2174: .  n - number of columns
2175: -  data - optional location of matrix data in column major order.  Set data=NULL for PETSc
2176:    to control all matrix memory allocation.

2178:    Output Parameter:
2179: .  A - the matrix

2181:    Notes:
2182:    The data input variable is intended primarily for Fortran programmers
2183:    who wish to allocate their own matrix memory space.  Most users should
2184:    set data=NULL.

2186:    Level: intermediate

2188: .keywords: dense, matrix, LAPACK, BLAS

2190: .seealso: MatCreate(), MatCreateDense(), MatSetValues()
2191: @*/
2192: PetscErrorCode  MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A)
2193: {

2197:   MatCreate(comm,A);
2198:   MatSetSizes(*A,m,n,m,n);
2199:   MatSetType(*A,MATSEQDENSE);
2200:   MatSeqDenseSetPreallocation(*A,data);
2201:   return(0);
2202: }

2206: /*@C
2207:    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements

2209:    Collective on MPI_Comm

2211:    Input Parameters:
2212: +  B - the matrix
2213: -  data - the array (or NULL)

2215:    Notes:
2216:    The data input variable is intended primarily for Fortran programmers
2217:    who wish to allocate their own matrix memory space.  Most users should
2218:    need not call this routine.

2220:    Level: intermediate

2222: .keywords: dense, matrix, LAPACK, BLAS

2224: .seealso: MatCreate(), MatCreateDense(), MatSetValues(), MatSeqDenseSetLDA()

2226: @*/
2227: PetscErrorCode  MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
2228: {

2232:   PetscTryMethod(B,"MatSeqDenseSetPreallocation_C",(Mat,PetscScalar[]),(B,data));
2233:   return(0);
2234: }

2238: PetscErrorCode  MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
2239: {
2240:   Mat_SeqDense   *b;

2244:   B->preallocated = PETSC_TRUE;

2246:   PetscLayoutSetUp(B->rmap);
2247:   PetscLayoutSetUp(B->cmap);

2249:   b       = (Mat_SeqDense*)B->data;
2250:   b->Mmax = B->rmap->n;
2251:   b->Nmax = B->cmap->n;
2252:   if (b->lda <= 0 || b->changelda) b->lda = B->rmap->n;

2254:   if (!data) { /* petsc-allocated storage */
2255:     if (!b->user_alloc) { PetscFree(b->v); }
2256:     PetscCalloc1((size_t)b->lda*b->Nmax,&b->v);
2257:     PetscLogObjectMemory((PetscObject)B,b->lda*b->Nmax*sizeof(PetscScalar));

2259:     b->user_alloc = PETSC_FALSE;
2260:   } else { /* user-allocated storage */
2261:     if (!b->user_alloc) { PetscFree(b->v); }
2262:     b->v          = data;
2263:     b->user_alloc = PETSC_TRUE;
2264:   }
2265:   B->assembled = PETSC_TRUE;
2266:   return(0);
2267: }

2269: #if defined(PETSC_HAVE_ELEMENTAL)
2272: PETSC_EXTERN PetscErrorCode MatConvert_SeqDense_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
2273: {
2274:   Mat            mat_elemental;
2276:   PetscScalar    *array,*v_colwise;
2277:   PetscInt       M=A->rmap->N,N=A->cmap->N,i,j,k,*rows,*cols;

2280:   PetscMalloc3(M*N,&v_colwise,M,&rows,N,&cols);
2281:   MatDenseGetArray(A,&array);
2282:   /* convert column-wise array into row-wise v_colwise, see MatSetValues_Elemental() */
2283:   k = 0;
2284:   for (j=0; j<N; j++) {
2285:     cols[j] = j;
2286:     for (i=0; i<M; i++) {
2287:       v_colwise[j*M+i] = array[k++];
2288:     }
2289:   }
2290:   for (i=0; i<M; i++) {
2291:     rows[i] = i;
2292:   }
2293:   MatDenseRestoreArray(A,&array);

2295:   MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);
2296:   MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,M,N);
2297:   MatSetType(mat_elemental,MATELEMENTAL);
2298:   MatSetUp(mat_elemental);

2300:   /* PETSc-Elemental interaface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
2301:   MatSetValues(mat_elemental,M,rows,N,cols,v_colwise,ADD_VALUES);
2302:   MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);
2303:   MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);
2304:   PetscFree3(v_colwise,rows,cols);

2306:   if (reuse == MAT_REUSE_MATRIX) {
2307:     MatHeaderReplace(A,mat_elemental);
2308:   } else {
2309:     *newmat = mat_elemental;
2310:   }
2311:   return(0);
2312: }
2313: #endif

2317: /*@C
2318:   MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array

2320:   Input parameter:
2321: + A - the matrix
2322: - lda - the leading dimension

2324:   Notes:
2325:   This routine is to be used in conjunction with MatSeqDenseSetPreallocation();
2326:   it asserts that the preallocation has a leading dimension (the LDA parameter
2327:   of Blas and Lapack fame) larger than M, the first dimension of the matrix.

2329:   Level: intermediate

2331: .keywords: dense, matrix, LAPACK, BLAS

2333: .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize()

2335: @*/
2336: PetscErrorCode  MatSeqDenseSetLDA(Mat B,PetscInt lda)
2337: {
2338:   Mat_SeqDense *b = (Mat_SeqDense*)B->data;

2341:   if (lda < B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->rmap->n);
2342:   b->lda       = lda;
2343:   b->changelda = PETSC_FALSE;
2344:   b->Mmax      = PetscMax(b->Mmax,lda);
2345:   return(0);
2346: }

2348: /*MC
2349:    MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices.

2351:    Options Database Keys:
2352: . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions()

2354:   Level: beginner

2356: .seealso: MatCreateSeqDense()

2358: M*/

2362: PETSC_EXTERN PetscErrorCode MatCreate_SeqDense(Mat B)
2363: {
2364:   Mat_SeqDense   *b;
2366:   PetscMPIInt    size;

2369:   MPI_Comm_size(PetscObjectComm((PetscObject)B),&size);
2370:   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Comm must be of size 1");

2372:   PetscNewLog(B,&b);
2373:   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
2374:   B->data = (void*)b;

2376:   b->pivots      = 0;
2377:   b->roworiented = PETSC_TRUE;
2378:   b->v           = 0;
2379:   b->changelda   = PETSC_FALSE;

2381:   PetscObjectComposeFunction((PetscObject)B,"MatDenseGetArray_C",MatDenseGetArray_SeqDense);
2382:   PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreArray_C",MatDenseRestoreArray_SeqDense);
2383:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqdense_seqaij_C",MatConvert_SeqDense_SeqAIJ);
2384: #if defined(PETSC_HAVE_ELEMENTAL)
2385:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqdense_elemental_C",MatConvert_SeqDense_Elemental);
2386: #endif
2387:   PetscObjectComposeFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",MatSeqDenseSetPreallocation_SeqDense);
2388:   PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqaij_seqdense_C",MatMatMult_SeqAIJ_SeqDense);
2389:   PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C",MatMatMultSymbolic_SeqAIJ_SeqDense);
2390:   PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C",MatMatMultNumeric_SeqAIJ_SeqDense);
2391:   PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMult_seqaij_seqdense_C",MatTransposeMatMult_SeqAIJ_SeqDense);
2392:   PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultSymbolic_seqaij_seqdense_C",MatTransposeMatMultSymbolic_SeqAIJ_SeqDense);
2393:   PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultNumeric_seqaij_seqdense_C",MatTransposeMatMultNumeric_SeqAIJ_SeqDense);
2394:   PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);
2395:   return(0);
2396: }