Actual source code: dense.c

petsc-master 2019-05-24
Report Typos and Errors

  2: /*
  3:      Defines the basic matrix operations for sequential dense.
  4: */

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

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

 11: static PetscErrorCode MatSeqDenseSymmetrize_Private(Mat A, PetscBool hermitian)
 12: {
 13:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
 14:   PetscInt      j, k, n = A->rmap->n;

 17:   if (A->rmap->n != A->cmap->n) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot symmetrize a rectangular matrix");
 18:   if (!hermitian) {
 19:     for (k=0;k<n;k++) {
 20:       for (j=k;j<n;j++) {
 21:         mat->v[j*mat->lda + k] = mat->v[k*mat->lda + j];
 22:       }
 23:     }
 24:   } else {
 25:     for (k=0;k<n;k++) {
 26:       for (j=k;j<n;j++) {
 27:         mat->v[j*mat->lda + k] = PetscConj(mat->v[k*mat->lda + j]);
 28:       }
 29:     }
 30:   }
 31:   return(0);
 32: }

 34: PETSC_EXTERN PetscErrorCode MatSeqDenseInvertFactors_Private(Mat A)
 35: {
 36: #if defined(PETSC_MISSING_LAPACK_POTRF)
 38:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable.");
 39: #else
 40:   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
 42:   PetscBLASInt   info,n;

 45:   if (!A->rmap->n || !A->cmap->n) return(0);
 46:   PetscBLASIntCast(A->cmap->n,&n);
 47:   if (A->factortype == MAT_FACTOR_LU) {
 48:     if (!mat->pivots) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Pivots not present");
 49:     if (!mat->fwork) {
 50:       mat->lfwork = n;
 51:       PetscMalloc1(mat->lfwork,&mat->fwork);
 52:       PetscLogObjectMemory((PetscObject)A,mat->lfwork*sizeof(PetscBLASInt));
 53:     }
 54:     PetscFPTrapPush(PETSC_FP_TRAP_OFF);
 55:     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&mat->lfwork,&info));
 56:     PetscFPTrapPop();
 57:     PetscLogFlops((1.0*A->cmap->n*A->cmap->n*A->cmap->n)/3.0); /* TODO CHECK FLOPS */
 58:   } else if (A->factortype == MAT_FACTOR_CHOLESKY) {
 59:     if (A->spd) {
 60:       PetscFPTrapPush(PETSC_FP_TRAP_OFF);
 61:       PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&n,mat->v,&mat->lda,&info));
 62:       PetscFPTrapPop();
 63:       MatSeqDenseSymmetrize_Private(A,PETSC_TRUE);
 64: #if defined(PETSC_USE_COMPLEX)
 65:     } else if (A->hermitian) {
 66:       if (!mat->pivots) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Pivots not present");
 67:       if (!mat->fwork) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Fwork not present");
 68:       PetscFPTrapPush(PETSC_FP_TRAP_OFF);
 69:       PetscStackCallBLAS("LAPACKhetri",LAPACKhetri_("L",&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&info));
 70:       PetscFPTrapPop();
 71:       MatSeqDenseSymmetrize_Private(A,PETSC_TRUE);
 72: #endif
 73:     } else { /* symmetric case */
 74:       if (!mat->pivots) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Pivots not present");
 75:       if (!mat->fwork) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Fwork not present");
 76:       PetscFPTrapPush(PETSC_FP_TRAP_OFF);
 77:       PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&info));
 78:       PetscFPTrapPop();
 79:       MatSeqDenseSymmetrize_Private(A,PETSC_FALSE);
 80:     }
 81:     if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Bad Inversion: zero pivot in row %D",(PetscInt)info-1);
 82:     PetscLogFlops((1.0*A->cmap->n*A->cmap->n*A->cmap->n)/3.0); /* TODO CHECK FLOPS */
 83:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
 84: #endif

 86:   A->ops->solve             = NULL;
 87:   A->ops->matsolve          = NULL;
 88:   A->ops->solvetranspose    = NULL;
 89:   A->ops->matsolvetranspose = NULL;
 90:   A->ops->solveadd          = NULL;
 91:   A->ops->solvetransposeadd = NULL;
 92:   A->factortype             = MAT_FACTOR_NONE;
 93:   PetscFree(A->solvertype);
 94:   return(0);
 95: }

 97: PetscErrorCode MatZeroRowsColumns_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
 98: {
 99:   PetscErrorCode    ierr;
100:   Mat_SeqDense      *l = (Mat_SeqDense*)A->data;
101:   PetscInt          m  = l->lda, n = A->cmap->n,r = A->rmap->n, i,j;
102:   PetscScalar       *slot,*bb;
103:   const PetscScalar *xx;

106: #if defined(PETSC_USE_DEBUG)
107:   for (i=0; i<N; i++) {
108:     if (rows[i] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row requested to be zeroed");
109:     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);
110:     if (rows[i] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Col %D requested to be zeroed greater than or equal number of cols %D",rows[i],A->cmap->n);
111:   }
112: #endif

114:   /* fix right hand side if needed */
115:   if (x && b) {
116:     Vec xt;

118:     if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only coded for square matrices");
119:     VecDuplicate(x,&xt);
120:     VecCopy(x,xt);
121:     VecScale(xt,-1.0);
122:     MatMultAdd(A,xt,b,b);
123:     VecDestroy(&xt);
124:     VecGetArrayRead(x,&xx);
125:     VecGetArray(b,&bb);
126:     for (i=0; i<N; i++) bb[rows[i]] = diag*xx[rows[i]];
127:     VecRestoreArrayRead(x,&xx);
128:     VecRestoreArray(b,&bb);
129:   }

131:   for (i=0; i<N; i++) {
132:     slot = l->v + rows[i]*m;
133:     PetscMemzero(slot,r*sizeof(PetscScalar));
134:   }
135:   for (i=0; i<N; i++) {
136:     slot = l->v + rows[i];
137:     for (j=0; j<n; j++) { *slot = 0.0; slot += m;}
138:   }
139:   if (diag != 0.0) {
140:     if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only coded for square matrices");
141:     for (i=0; i<N; i++) {
142:       slot  = l->v + (m+1)*rows[i];
143:       *slot = diag;
144:     }
145:   }
146:   return(0);
147: }

149: PetscErrorCode MatPtAPNumeric_SeqDense_SeqDense(Mat A,Mat P,Mat C)
150: {
151:   Mat_SeqDense   *c = (Mat_SeqDense*)(C->data);

155:   MatMatMultNumeric_SeqDense_SeqDense(A,P,c->ptapwork);
156:   MatTransposeMatMultNumeric_SeqDense_SeqDense(P,c->ptapwork,C);
157:   return(0);
158: }

160: PetscErrorCode MatPtAPSymbolic_SeqDense_SeqDense(Mat A,Mat P,PetscReal fill,Mat *C)
161: {
162:   Mat_SeqDense   *c;

166:   MatCreateSeqDense(PetscObjectComm((PetscObject)A),P->cmap->N,P->cmap->N,NULL,C);
167:   c = (Mat_SeqDense*)((*C)->data);
168:   MatCreateSeqDense(PetscObjectComm((PetscObject)A),A->rmap->N,P->cmap->N,NULL,&c->ptapwork);
169:   return(0);
170: }

172: PETSC_INTERN PetscErrorCode MatPtAP_SeqDense_SeqDense(Mat A,Mat P,MatReuse reuse,PetscReal fill,Mat *C)
173: {

177:   if (reuse == MAT_INITIAL_MATRIX) {
178:     MatPtAPSymbolic_SeqDense_SeqDense(A,P,fill,C);
179:   }
180:   PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);
181:   (*(*C)->ops->ptapnumeric)(A,P,*C);
182:   PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);
183:   return(0);
184: }

186: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqDense(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
187: {
188:   Mat            B = NULL;
189:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
190:   Mat_SeqDense   *b;
192:   PetscInt       *ai=a->i,*aj=a->j,m=A->rmap->N,n=A->cmap->N,i;
193:   MatScalar      *av=a->a;
194:   PetscBool      isseqdense;

197:   if (reuse == MAT_REUSE_MATRIX) {
198:     PetscObjectTypeCompare((PetscObject)*newmat,MATSEQDENSE,&isseqdense);
199:     if (!isseqdense) SETERRQ1(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_USER,"Cannot reuse matrix of type %s",((PetscObject)(*newmat))->type);
200:   }
201:   if (reuse != MAT_REUSE_MATRIX) {
202:     MatCreate(PetscObjectComm((PetscObject)A),&B);
203:     MatSetSizes(B,m,n,m,n);
204:     MatSetType(B,MATSEQDENSE);
205:     MatSeqDenseSetPreallocation(B,NULL);
206:     b    = (Mat_SeqDense*)(B->data);
207:   } else {
208:     b    = (Mat_SeqDense*)((*newmat)->data);
209:     PetscMemzero(b->v,m*n*sizeof(PetscScalar));
210:   }
211:   for (i=0; i<m; i++) {
212:     PetscInt j;
213:     for (j=0;j<ai[1]-ai[0];j++) {
214:       b->v[*aj*m+i] = *av;
215:       aj++;
216:       av++;
217:     }
218:     ai++;
219:   }

221:   if (reuse == MAT_INPLACE_MATRIX) {
222:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
223:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
224:     MatHeaderReplace(A,&B);
225:   } else {
226:     if (B) *newmat = B;
227:     MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);
228:     MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);
229:   }
230:   return(0);
231: }

233: PETSC_INTERN PetscErrorCode MatConvert_SeqDense_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
234: {
235:   Mat            B;
236:   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
238:   PetscInt       i, j;
239:   PetscInt       *rows, *nnz;
240:   MatScalar      *aa = a->v, *vals;

243:   MatCreate(PetscObjectComm((PetscObject)A),&B);
244:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
245:   MatSetType(B,MATSEQAIJ);
246:   PetscCalloc3(A->rmap->n,&rows,A->rmap->n,&nnz,A->rmap->n,&vals);
247:   for (j=0; j<A->cmap->n; j++) {
248:     for (i=0; i<A->rmap->n; i++) if (aa[i] != 0.0 || i == j) ++nnz[i];
249:     aa += a->lda;
250:   }
251:   MatSeqAIJSetPreallocation(B,PETSC_DETERMINE,nnz);
252:   aa = a->v;
253:   for (j=0; j<A->cmap->n; j++) {
254:     PetscInt numRows = 0;
255:     for (i=0; i<A->rmap->n; i++) if (aa[i] != 0.0 || i == j) {rows[numRows] = i; vals[numRows++] = aa[i];}
256:     MatSetValues(B,numRows,rows,1,&j,vals,INSERT_VALUES);
257:     aa  += a->lda;
258:   }
259:   PetscFree3(rows,nnz,vals);
260:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
261:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

263:   if (reuse == MAT_INPLACE_MATRIX) {
264:     MatHeaderReplace(A,&B);
265:   } else {
266:     *newmat = B;
267:   }
268:   return(0);
269: }

271: static PetscErrorCode MatAXPY_SeqDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str)
272: {
273:   Mat_SeqDense   *x     = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data;
274:   PetscScalar    oalpha = alpha;
275:   PetscInt       j;
276:   PetscBLASInt   N,m,ldax,lday,one = 1;

280:   PetscBLASIntCast(X->rmap->n*X->cmap->n,&N);
281:   PetscBLASIntCast(X->rmap->n,&m);
282:   PetscBLASIntCast(x->lda,&ldax);
283:   PetscBLASIntCast(y->lda,&lday);
284:   if (ldax>m || lday>m) {
285:     for (j=0; j<X->cmap->n; j++) {
286:       PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&m,&oalpha,x->v+j*ldax,&one,y->v+j*lday,&one));
287:     }
288:   } else {
289:     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&N,&oalpha,x->v,&one,y->v,&one));
290:   }
291:   PetscObjectStateIncrease((PetscObject)Y);
292:   PetscLogFlops(PetscMax(2*N-1,0));
293:   return(0);
294: }

296: static PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
297: {
298:   PetscInt N = A->rmap->n*A->cmap->n;

301:   info->block_size        = 1.0;
302:   info->nz_allocated      = (double)N;
303:   info->nz_used           = (double)N;
304:   info->nz_unneeded       = (double)0;
305:   info->assemblies        = (double)A->num_ass;
306:   info->mallocs           = 0;
307:   info->memory            = ((PetscObject)A)->mem;
308:   info->fill_ratio_given  = 0;
309:   info->fill_ratio_needed = 0;
310:   info->factor_mallocs    = 0;
311:   return(0);
312: }

314: static PetscErrorCode MatScale_SeqDense(Mat A,PetscScalar alpha)
315: {
316:   Mat_SeqDense   *a     = (Mat_SeqDense*)A->data;
317:   PetscScalar    oalpha = alpha;
319:   PetscBLASInt   one = 1,j,nz,lda;

322:   PetscBLASIntCast(a->lda,&lda);
323:   if (lda>A->rmap->n) {
324:     PetscBLASIntCast(A->rmap->n,&nz);
325:     for (j=0; j<A->cmap->n; j++) {
326:       PetscStackCallBLAS("BLASscal",BLASscal_(&nz,&oalpha,a->v+j*lda,&one));
327:     }
328:   } else {
329:     PetscBLASIntCast(A->rmap->n*A->cmap->n,&nz);
330:     PetscStackCallBLAS("BLASscal",BLASscal_(&nz,&oalpha,a->v,&one));
331:   }
332:   PetscLogFlops(nz);
333:   return(0);
334: }

336: static PetscErrorCode MatIsHermitian_SeqDense(Mat A,PetscReal rtol,PetscBool  *fl)
337: {
338:   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
339:   PetscInt     i,j,m = A->rmap->n,N;
340:   PetscScalar  *v = a->v;

343:   *fl = PETSC_FALSE;
344:   if (A->rmap->n != A->cmap->n) return(0);
345:   N = a->lda;

347:   for (i=0; i<m; i++) {
348:     for (j=i+1; j<m; j++) {
349:       if (PetscAbsScalar(v[i+j*N] - PetscConj(v[j+i*N])) > rtol) return(0);
350:     }
351:   }
352:   *fl = PETSC_TRUE;
353:   return(0);
354: }

356: static PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat newi,Mat A,MatDuplicateOption cpvalues)
357: {
358:   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data,*l;
360:   PetscInt       lda = (PetscInt)mat->lda,j,m;

363:   PetscLayoutReference(A->rmap,&newi->rmap);
364:   PetscLayoutReference(A->cmap,&newi->cmap);
365:   MatSeqDenseSetPreallocation(newi,NULL);
366:   if (cpvalues == MAT_COPY_VALUES) {
367:     l = (Mat_SeqDense*)newi->data;
368:     if (lda>A->rmap->n) {
369:       m = A->rmap->n;
370:       for (j=0; j<A->cmap->n; j++) {
371:         PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));
372:       }
373:     } else {
374:       PetscMemcpy(l->v,mat->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));
375:     }
376:   }
377:   newi->assembled = PETSC_TRUE;
378:   return(0);
379: }

381: static PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
382: {

386:   MatCreate(PetscObjectComm((PetscObject)A),newmat);
387:   MatSetSizes(*newmat,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);
388:   MatSetType(*newmat,((PetscObject)A)->type_name);
389:   MatDuplicateNoCreate_SeqDense(*newmat,A,cpvalues);
390:   return(0);
391: }


394: static PetscErrorCode MatLUFactor_SeqDense(Mat,IS,IS,const MatFactorInfo*);

396: static PetscErrorCode MatLUFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy)
397: {
398:   MatFactorInfo  info;

402:   MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);
403:   MatLUFactor_SeqDense(fact,0,0,&info);
404:   return(0);
405: }

407: static PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
408: {
409:   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
410:   PetscErrorCode    ierr;
411:   const PetscScalar *x;
412:   PetscScalar       *y;
413:   PetscBLASInt      one = 1,info,m;

416:   PetscBLASIntCast(A->rmap->n,&m);
417:   VecGetArrayRead(xx,&x);
418:   VecGetArray(yy,&y);
419:   PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));
420:   if (A->factortype == MAT_FACTOR_LU) {
421: #if defined(PETSC_MISSING_LAPACK_GETRS)
422:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
423: #else
424:     PetscFPTrapPush(PETSC_FP_TRAP_OFF);
425:     PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
426:     PetscFPTrapPop();
427:     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve");
428: #endif
429:   } else if (A->factortype == MAT_FACTOR_CHOLESKY) {
430: #if defined(PETSC_MISSING_LAPACK_POTRS)
431:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
432: #else
433:     if (A->spd) {
434:       PetscFPTrapPush(PETSC_FP_TRAP_OFF);
435:       PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info));
436:       PetscFPTrapPop();
437:       if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve");
438: #if defined(PETSC_USE_COMPLEX)
439:     } else if (A->hermitian) {
440:       PetscFPTrapPush(PETSC_FP_TRAP_OFF);
441:       PetscStackCallBLAS("LAPACKhetrs",LAPACKhetrs_("L",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
442:       PetscFPTrapPop();
443:       if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"HETRS Bad solve");
444: #endif
445:     } else { /* symmetric case */
446:       PetscFPTrapPush(PETSC_FP_TRAP_OFF);
447:       PetscStackCallBLAS("LAPACKsytrs",LAPACKsytrs_("L",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
448:       PetscFPTrapPop();
449:       if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SYTRS Bad solve");
450:     }
451: #endif
452:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
453:   VecRestoreArrayRead(xx,&x);
454:   VecRestoreArray(yy,&y);
455:   PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);
456:   return(0);
457: }

459: static PetscErrorCode MatMatSolve_SeqDense(Mat A,Mat B,Mat X)
460: {
461:   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
462:   PetscErrorCode    ierr;
463:   const PetscScalar *b;
464:   PetscScalar       *x;
465:   PetscInt          n;
466:   PetscBLASInt      nrhs,info,m;
467:   PetscBool         flg;

470:   PetscBLASIntCast(A->rmap->n,&m);
471:   PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
472:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
473:   PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
474:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");

476:   MatGetSize(B,NULL,&n);
477:   PetscBLASIntCast(n,&nrhs);
478:   MatDenseGetArrayRead(B,&b);
479:   MatDenseGetArray(X,&x);

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

483:   if (A->factortype == MAT_FACTOR_LU) {
484: #if defined(PETSC_MISSING_LAPACK_GETRS)
485:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
486: #else
487:     PetscFPTrapPush(PETSC_FP_TRAP_OFF);
488:     PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info));
489:     PetscFPTrapPop();
490:     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve");
491: #endif
492:   } else if (A->factortype == MAT_FACTOR_CHOLESKY) {
493: #if defined(PETSC_MISSING_LAPACK_POTRS)
494:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
495: #else
496:     if (A->spd) {
497:       PetscFPTrapPush(PETSC_FP_TRAP_OFF);
498:       PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&m,&nrhs,mat->v,&mat->lda,x,&m,&info));
499:       PetscFPTrapPop();
500:       if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve");
501: #if defined(PETSC_USE_COMPLEX)
502:     } else if (A->hermitian) {
503:       PetscFPTrapPush(PETSC_FP_TRAP_OFF);
504:       PetscStackCallBLAS("LAPACKhetrs",LAPACKhetrs_("L",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info));
505:       PetscFPTrapPop();
506:       if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"HETRS Bad solve");
507: #endif
508:     } else { /* symmetric case */
509:       PetscFPTrapPush(PETSC_FP_TRAP_OFF);
510:       PetscStackCallBLAS("LAPACKsytrs",LAPACKsytrs_("L",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info));
511:       PetscFPTrapPop();
512:       if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SYTRS Bad solve");
513:     }
514: #endif
515:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");

517:   MatDenseRestoreArrayRead(B,&b);
518:   MatDenseRestoreArray(X,&x);
519:   PetscLogFlops(nrhs*(2.0*m*m - m));
520:   return(0);
521: }

523: static PetscErrorCode MatConjugate_SeqDense(Mat);

525: static PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
526: {
527:   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
528:   PetscErrorCode    ierr;
529:   const PetscScalar *x;
530:   PetscScalar       *y;
531:   PetscBLASInt      one = 1,info,m;

534:   PetscBLASIntCast(A->rmap->n,&m);
535:   VecGetArrayRead(xx,&x);
536:   VecGetArray(yy,&y);
537:   PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));
538:   if (A->factortype == MAT_FACTOR_LU) {
539: #if defined(PETSC_MISSING_LAPACK_GETRS)
540:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
541: #else
542:     PetscFPTrapPush(PETSC_FP_TRAP_OFF);
543:     PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
544:     PetscFPTrapPop();
545:     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve");
546: #endif
547:   } else if (A->factortype == MAT_FACTOR_CHOLESKY) {
548: #if defined(PETSC_MISSING_LAPACK_POTRS)
549:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
550: #else
551:     if (A->spd) {
552: #if defined(PETSC_USE_COMPLEX)
553:       MatConjugate_SeqDense(A);
554: #endif
555:       PetscFPTrapPush(PETSC_FP_TRAP_OFF);
556:       PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info));
557:       PetscFPTrapPop();
558: #if defined(PETSC_USE_COMPLEX)
559:       MatConjugate_SeqDense(A);
560: #endif
561:       if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve");
562: #if defined(PETSC_USE_COMPLEX)
563:     } else if (A->hermitian) {
564:       MatConjugate_SeqDense(A);
565:       PetscFPTrapPush(PETSC_FP_TRAP_OFF);
566:       PetscStackCallBLAS("LAPACKhetrs",LAPACKhetrs_("L",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
567:       PetscFPTrapPop();
568:       MatConjugate_SeqDense(A);
569: #endif
570:     } else { /* symmetric case */
571:       PetscFPTrapPush(PETSC_FP_TRAP_OFF);
572:       PetscStackCallBLAS("LAPACKsytrs",LAPACKsytrs_("L",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
573:       PetscFPTrapPop();
574:       if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SYTRS Bad solve");
575:     }
576: #endif
577:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
578:   VecRestoreArrayRead(xx,&x);
579:   VecRestoreArray(yy,&y);
580:   PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);
581:   return(0);
582: }

584: static PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
585: {
586:   PetscErrorCode    ierr;
587:   const PetscScalar *x;
588:   PetscScalar       *y,sone = 1.0;
589:   Vec               tmp = 0;

592:   VecGetArrayRead(xx,&x);
593:   VecGetArray(yy,&y);
594:   if (!A->rmap->n || !A->cmap->n) return(0);
595:   if (yy == zz) {
596:     VecDuplicate(yy,&tmp);
597:     PetscLogObjectParent((PetscObject)A,(PetscObject)tmp);
598:     VecCopy(yy,tmp);
599:   }
600:   MatSolve_SeqDense(A,xx,yy);
601:   if (tmp) {
602:     VecAXPY(yy,sone,tmp);
603:     VecDestroy(&tmp);
604:   } else {
605:     VecAXPY(yy,sone,zz);
606:   }
607:   VecRestoreArrayRead(xx,&x);
608:   VecRestoreArray(yy,&y);
609:   PetscLogFlops(A->cmap->n);
610:   return(0);
611: }

613: static PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
614: {
615:   PetscErrorCode    ierr;
616:   const PetscScalar *x;
617:   PetscScalar       *y,sone = 1.0;
618:   Vec               tmp = NULL;

621:   if (!A->rmap->n || !A->cmap->n) return(0);
622:   VecGetArrayRead(xx,&x);
623:   VecGetArray(yy,&y);
624:   if (yy == zz) {
625:     VecDuplicate(yy,&tmp);
626:     PetscLogObjectParent((PetscObject)A,(PetscObject)tmp);
627:     VecCopy(yy,tmp);
628:   }
629:   MatSolveTranspose_SeqDense(A,xx,yy);
630:   if (tmp) {
631:     VecAXPY(yy,sone,tmp);
632:     VecDestroy(&tmp);
633:   } else {
634:     VecAXPY(yy,sone,zz);
635:   }
636:   VecRestoreArrayRead(xx,&x);
637:   VecRestoreArray(yy,&y);
638:   PetscLogFlops(A->cmap->n);
639:   return(0);
640: }

642: /* ---------------------------------------------------------------*/
643: /* COMMENT: I have chosen to hide row permutation in the pivots,
644:    rather than put it in the Mat->row slot.*/
645: static PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,const MatFactorInfo *minfo)
646: {
647: #if defined(PETSC_MISSING_LAPACK_GETRF)
649:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable.");
650: #else
651:   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
653:   PetscBLASInt   n,m,info;

656:   PetscBLASIntCast(A->cmap->n,&n);
657:   PetscBLASIntCast(A->rmap->n,&m);
658:   if (!mat->pivots) {
659:     PetscMalloc1(A->rmap->n,&mat->pivots);
660:     PetscLogObjectMemory((PetscObject)A,A->rmap->n*sizeof(PetscBLASInt));
661:   }
662:   if (!A->rmap->n || !A->cmap->n) return(0);
663:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
664:   PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info));
665:   PetscFPTrapPop();

667:   if (info<0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to LU factorization");
668:   if (info>0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");

670:   A->ops->solve             = MatSolve_SeqDense;
671:   A->ops->matsolve          = MatMatSolve_SeqDense;
672:   A->ops->solvetranspose    = MatSolveTranspose_SeqDense;
673:   A->ops->solveadd          = MatSolveAdd_SeqDense;
674:   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense;
675:   A->factortype             = MAT_FACTOR_LU;

677:   PetscFree(A->solvertype);
678:   PetscStrallocpy(MATSOLVERPETSC,&A->solvertype);

680:   PetscLogFlops((2.0*A->cmap->n*A->cmap->n*A->cmap->n)/3);
681: #endif
682:   return(0);
683: }

685: /* Cholesky as L*L^T or L*D*L^T and the symmetric/hermitian complex variants */
686: static PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,const MatFactorInfo *factinfo)
687: {
688: #if defined(PETSC_MISSING_LAPACK_POTRF)
690:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable.");
691: #else
692:   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
694:   PetscBLASInt   info,n;

697:   PetscBLASIntCast(A->cmap->n,&n);
698:   if (!A->rmap->n || !A->cmap->n) return(0);
699:   if (A->spd) {
700:     PetscFPTrapPush(PETSC_FP_TRAP_OFF);
701:     PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info));
702:     PetscFPTrapPop();
703: #if defined(PETSC_USE_COMPLEX)
704:   } else if (A->hermitian) {
705:     if (!mat->pivots) {
706:       PetscMalloc1(A->rmap->n,&mat->pivots);
707:       PetscLogObjectMemory((PetscObject)A,A->rmap->n*sizeof(PetscBLASInt));
708:     }
709:     if (!mat->fwork) {
710:       PetscScalar dummy;

712:       mat->lfwork = -1;
713:       PetscFPTrapPush(PETSC_FP_TRAP_OFF);
714:       PetscStackCallBLAS("LAPACKhetrf",LAPACKhetrf_("L",&n,mat->v,&mat->lda,mat->pivots,&dummy,&mat->lfwork,&info));
715:       PetscFPTrapPop();
716:       mat->lfwork = (PetscInt)PetscRealPart(dummy);
717:       PetscMalloc1(mat->lfwork,&mat->fwork);
718:       PetscLogObjectMemory((PetscObject)A,mat->lfwork*sizeof(PetscBLASInt));
719:     }
720:     PetscFPTrapPush(PETSC_FP_TRAP_OFF);
721:     PetscStackCallBLAS("LAPACKhetrf",LAPACKhetrf_("L",&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&mat->lfwork,&info));
722:     PetscFPTrapPop();
723: #endif
724:   } else { /* symmetric case */
725:     if (!mat->pivots) {
726:       PetscMalloc1(A->rmap->n,&mat->pivots);
727:       PetscLogObjectMemory((PetscObject)A,A->rmap->n*sizeof(PetscBLASInt));
728:     }
729:     if (!mat->fwork) {
730:       PetscScalar dummy;

732:       mat->lfwork = -1;
733:       PetscFPTrapPush(PETSC_FP_TRAP_OFF);
734:       PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&n,mat->v,&mat->lda,mat->pivots,&dummy,&mat->lfwork,&info));
735:       PetscFPTrapPop();
736:       mat->lfwork = (PetscInt)PetscRealPart(dummy);
737:       PetscMalloc1(mat->lfwork,&mat->fwork);
738:       PetscLogObjectMemory((PetscObject)A,mat->lfwork*sizeof(PetscBLASInt));
739:     }
740:     PetscFPTrapPush(PETSC_FP_TRAP_OFF);
741:     PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&mat->lfwork,&info));
742:     PetscFPTrapPop();
743:   }
744:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1);

746:   A->ops->solve             = MatSolve_SeqDense;
747:   A->ops->matsolve          = MatMatSolve_SeqDense;
748:   A->ops->solvetranspose    = MatSolveTranspose_SeqDense;
749:   A->ops->solveadd          = MatSolveAdd_SeqDense;
750:   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense;
751:   A->factortype             = MAT_FACTOR_CHOLESKY;

753:   PetscFree(A->solvertype);
754:   PetscStrallocpy(MATSOLVERPETSC,&A->solvertype);

756:   PetscLogFlops((1.0*A->cmap->n*A->cmap->n*A->cmap->n)/3.0);
757: #endif
758:   return(0);
759: }


762: PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy)
763: {
765:   MatFactorInfo  info;

768:   info.fill = 1.0;

770:   MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);
771:   MatCholeskyFactor_SeqDense(fact,0,&info);
772:   return(0);
773: }

775: static PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,const MatFactorInfo *info)
776: {
778:   fact->assembled                  = PETSC_TRUE;
779:   fact->preallocated               = PETSC_TRUE;
780:   fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqDense;
781:   fact->ops->solve                 = MatSolve_SeqDense;
782:   fact->ops->matsolve              = MatMatSolve_SeqDense;
783:   fact->ops->solvetranspose        = MatSolveTranspose_SeqDense;
784:   fact->ops->solveadd              = MatSolveAdd_SeqDense;
785:   fact->ops->solvetransposeadd     = MatSolveTransposeAdd_SeqDense;
786:   return(0);
787: }

789: static PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info)
790: {
792:   fact->preallocated           = PETSC_TRUE;
793:   fact->assembled              = PETSC_TRUE;
794:   fact->ops->lufactornumeric   = MatLUFactorNumeric_SeqDense;
795:   fact->ops->solve             = MatSolve_SeqDense;
796:   fact->ops->matsolve          = MatMatSolve_SeqDense;
797:   fact->ops->solvetranspose    = MatSolveTranspose_SeqDense;
798:   fact->ops->solveadd          = MatSolveAdd_SeqDense;
799:   fact->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense;
800:   return(0);
801: }

803: PETSC_INTERN PetscErrorCode MatGetFactor_seqdense_petsc(Mat A,MatFactorType ftype,Mat *fact)
804: {

808:   MatCreate(PetscObjectComm((PetscObject)A),fact);
809:   MatSetSizes(*fact,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);
810:   MatSetType(*fact,((PetscObject)A)->type_name);
811:   if (ftype == MAT_FACTOR_LU) {
812:     (*fact)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqDense;
813:   } else {
814:     (*fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqDense;
815:   }
816:   (*fact)->factortype = ftype;

818:   PetscFree((*fact)->solvertype);
819:   PetscStrallocpy(MATSOLVERPETSC,&(*fact)->solvertype);
820:   return(0);
821: }

823: /* ------------------------------------------------------------------*/
824: static PetscErrorCode MatSOR_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx)
825: {
826:   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
827:   PetscScalar       *x,*v = mat->v,zero = 0.0,xt;
828:   const PetscScalar *b;
829:   PetscErrorCode    ierr;
830:   PetscInt          m = A->rmap->n,i;
831:   PetscBLASInt      o = 1,bm;

834:   if (shift == -1) shift = 0.0; /* negative shift indicates do not error on zero diagonal; this code never zeros on zero diagonal */
835:   PetscBLASIntCast(m,&bm);
836:   if (flag & SOR_ZERO_INITIAL_GUESS) {
837:     /* this is a hack fix, should have another version without the second BLASdotu */
838:     VecSet(xx,zero);
839:   }
840:   VecGetArray(xx,&x);
841:   VecGetArrayRead(bb,&b);
842:   its  = its*lits;
843:   if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
844:   while (its--) {
845:     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
846:       for (i=0; i<m; i++) {
847:         PetscStackCallBLAS("BLASdotu",xt   = b[i] - BLASdotu_(&bm,v+i,&bm,x,&o));
848:         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
849:       }
850:     }
851:     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
852:       for (i=m-1; i>=0; i--) {
853:         PetscStackCallBLAS("BLASdotu",xt   = b[i] - BLASdotu_(&bm,v+i,&bm,x,&o));
854:         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
855:       }
856:     }
857:   }
858:   VecRestoreArrayRead(bb,&b);
859:   VecRestoreArray(xx,&x);
860:   return(0);
861: }

863: /* -----------------------------------------------------------------*/
864: PETSC_INTERN PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
865: {
866:   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
867:   const PetscScalar *v   = mat->v,*x;
868:   PetscScalar       *y;
869:   PetscErrorCode    ierr;
870:   PetscBLASInt      m, n,_One=1;
871:   PetscScalar       _DOne=1.0,_DZero=0.0;

874:   PetscBLASIntCast(A->rmap->n,&m);
875:   PetscBLASIntCast(A->cmap->n,&n);
876:   VecGetArrayRead(xx,&x);
877:   VecGetArray(yy,&y);
878:   if (!A->rmap->n || !A->cmap->n) {
879:     PetscBLASInt i;
880:     for (i=0; i<n; i++) y[i] = 0.0;
881:   } else {
882:     PetscStackCallBLAS("BLASgemv",BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One));
883:     PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->cmap->n);
884:   }
885:   VecRestoreArrayRead(xx,&x);
886:   VecRestoreArray(yy,&y);
887:   return(0);
888: }

890: PETSC_INTERN PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy)
891: {
892:   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
893:   PetscScalar       *y,_DOne=1.0,_DZero=0.0;
894:   PetscErrorCode    ierr;
895:   PetscBLASInt      m, n, _One=1;
896:   const PetscScalar *v = mat->v,*x;

899:   PetscBLASIntCast(A->rmap->n,&m);
900:   PetscBLASIntCast(A->cmap->n,&n);
901:   VecGetArrayRead(xx,&x);
902:   VecGetArray(yy,&y);
903:   if (!A->rmap->n || !A->cmap->n) {
904:     PetscBLASInt i;
905:     for (i=0; i<m; i++) y[i] = 0.0;
906:   } else {
907:     PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One));
908:     PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->rmap->n);
909:   }
910:   VecRestoreArrayRead(xx,&x);
911:   VecRestoreArray(yy,&y);
912:   return(0);
913: }

915: PETSC_INTERN PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
916: {
917:   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
918:   const PetscScalar *v = mat->v,*x;
919:   PetscScalar       *y,_DOne=1.0;
920:   PetscErrorCode    ierr;
921:   PetscBLASInt      m, n, _One=1;

924:   PetscBLASIntCast(A->rmap->n,&m);
925:   PetscBLASIntCast(A->cmap->n,&n);
926:   if (!A->rmap->n || !A->cmap->n) return(0);
927:   if (zz != yy) {VecCopy(zz,yy);}
928:   VecGetArrayRead(xx,&x);
929:   VecGetArray(yy,&y);
930:   PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One));
931:   VecRestoreArrayRead(xx,&x);
932:   VecRestoreArray(yy,&y);
933:   PetscLogFlops(2.0*A->rmap->n*A->cmap->n);
934:   return(0);
935: }

937: static PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
938: {
939:   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
940:   const PetscScalar *v = mat->v,*x;
941:   PetscScalar       *y;
942:   PetscErrorCode    ierr;
943:   PetscBLASInt      m, n, _One=1;
944:   PetscScalar       _DOne=1.0;

947:   PetscBLASIntCast(A->rmap->n,&m);
948:   PetscBLASIntCast(A->cmap->n,&n);
949:   if (!A->rmap->n || !A->cmap->n) return(0);
950:   if (zz != yy) {VecCopy(zz,yy);}
951:   VecGetArrayRead(xx,&x);
952:   VecGetArray(yy,&y);
953:   PetscStackCallBLAS("BLASgemv",BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One));
954:   VecRestoreArrayRead(xx,&x);
955:   VecRestoreArray(yy,&y);
956:   PetscLogFlops(2.0*A->rmap->n*A->cmap->n);
957:   return(0);
958: }

960: /* -----------------------------------------------------------------*/
961: static PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
962: {
963:   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
964:   PetscScalar    *v;
966:   PetscInt       i;

969:   *ncols = A->cmap->n;
970:   if (cols) {
971:     PetscMalloc1(A->cmap->n+1,cols);
972:     for (i=0; i<A->cmap->n; i++) (*cols)[i] = i;
973:   }
974:   if (vals) {
975:     PetscMalloc1(A->cmap->n+1,vals);
976:     v    = mat->v + row;
977:     for (i=0; i<A->cmap->n; i++) {(*vals)[i] = *v; v += mat->lda;}
978:   }
979:   return(0);
980: }

982: static PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
983: {

987:   if (cols) {PetscFree(*cols);}
988:   if (vals) {PetscFree(*vals); }
989:   return(0);
990: }
991: /* ----------------------------------------------------------------*/
992: static PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv)
993: {
994:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
995:   PetscInt     i,j,idx=0;

998:   if (!mat->roworiented) {
999:     if (addv == INSERT_VALUES) {
1000:       for (j=0; j<n; j++) {
1001:         if (indexn[j] < 0) {idx += m; continue;}
1002: #if defined(PETSC_USE_DEBUG)
1003:         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);
1004: #endif
1005:         for (i=0; i<m; i++) {
1006:           if (indexm[i] < 0) {idx++; continue;}
1007: #if defined(PETSC_USE_DEBUG)
1008:           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);
1009: #endif
1010:           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
1011:         }
1012:       }
1013:     } else {
1014:       for (j=0; j<n; j++) {
1015:         if (indexn[j] < 0) {idx += m; continue;}
1016: #if defined(PETSC_USE_DEBUG)
1017:         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);
1018: #endif
1019:         for (i=0; i<m; i++) {
1020:           if (indexm[i] < 0) {idx++; continue;}
1021: #if defined(PETSC_USE_DEBUG)
1022:           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);
1023: #endif
1024:           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
1025:         }
1026:       }
1027:     }
1028:   } else {
1029:     if (addv == INSERT_VALUES) {
1030:       for (i=0; i<m; i++) {
1031:         if (indexm[i] < 0) { idx += n; continue;}
1032: #if defined(PETSC_USE_DEBUG)
1033:         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);
1034: #endif
1035:         for (j=0; j<n; j++) {
1036:           if (indexn[j] < 0) { idx++; continue;}
1037: #if defined(PETSC_USE_DEBUG)
1038:           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);
1039: #endif
1040:           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
1041:         }
1042:       }
1043:     } else {
1044:       for (i=0; i<m; i++) {
1045:         if (indexm[i] < 0) { idx += n; continue;}
1046: #if defined(PETSC_USE_DEBUG)
1047:         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);
1048: #endif
1049:         for (j=0; j<n; j++) {
1050:           if (indexn[j] < 0) { idx++; continue;}
1051: #if defined(PETSC_USE_DEBUG)
1052:           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);
1053: #endif
1054:           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
1055:         }
1056:       }
1057:     }
1058:   }
1059:   return(0);
1060: }

1062: static PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[])
1063: {
1064:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1065:   PetscInt     i,j;

1068:   /* row-oriented output */
1069:   for (i=0; i<m; i++) {
1070:     if (indexm[i] < 0) {v += n;continue;}
1071:     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);
1072:     for (j=0; j<n; j++) {
1073:       if (indexn[j] < 0) {v++; continue;}
1074:       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);
1075:       *v++ = mat->v[indexn[j]*mat->lda + indexm[i]];
1076:     }
1077:   }
1078:   return(0);
1079: }

1081: /* -----------------------------------------------------------------*/

1083: static PetscErrorCode MatLoad_SeqDense(Mat newmat,PetscViewer viewer)
1084: {
1085:   Mat_SeqDense   *a;
1087:   PetscInt       *scols,i,j,nz,header[4];
1088:   int            fd;
1089:   PetscMPIInt    size;
1090:   PetscInt       *rowlengths = 0,M,N,*cols,grows,gcols;
1091:   PetscScalar    *vals,*svals,*v,*w;
1092:   MPI_Comm       comm;
1093:   PetscBool      isbinary;

1096:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
1097:   if (!isbinary) SETERRQ2(PetscObjectComm((PetscObject)newmat),PETSC_ERR_SUP,"Viewer type %s not yet supported for reading %s matrices",((PetscObject)viewer)->type_name,((PetscObject)newmat)->type_name);

1099:   /* force binary viewer to load .info file if it has not yet done so */
1100:   PetscViewerSetUp(viewer);
1101:   PetscObjectGetComm((PetscObject)viewer,&comm);
1102:   MPI_Comm_size(comm,&size);
1103:   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"view must have one processor");
1104:   PetscViewerBinaryGetDescriptor(viewer,&fd);
1105:   PetscBinaryRead(fd,header,4,NULL,PETSC_INT);
1106:   if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
1107:   M = header[1]; N = header[2]; nz = header[3];

1109:   /* set global size if not set already*/
1110:   if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) {
1111:     MatSetSizes(newmat,M,N,M,N);
1112:   } else {
1113:     /* if sizes and type are already set, check if the vector global sizes are correct */
1114:     MatGetSize(newmat,&grows,&gcols);
1115:     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);
1116:   }
1117:   a = (Mat_SeqDense*)newmat->data;
1118:   if (!a->user_alloc) {
1119:     MatSeqDenseSetPreallocation(newmat,NULL);
1120:   }

1122:   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
1123:     a = (Mat_SeqDense*)newmat->data;
1124:     v = a->v;
1125:     /* Allocate some temp space to read in the values and then flip them
1126:        from row major to column major */
1127:     PetscMalloc1(M*N > 0 ? M*N : 1,&w);
1128:     /* read in nonzero values */
1129:     PetscBinaryRead(fd,w,M*N,NULL,PETSC_SCALAR);
1130:     /* now flip the values and store them in the matrix*/
1131:     for (j=0; j<N; j++) {
1132:       for (i=0; i<M; i++) {
1133:         *v++ =w[i*N+j];
1134:       }
1135:     }
1136:     PetscFree(w);
1137:   } else {
1138:     /* read row lengths */
1139:     PetscMalloc1(M+1,&rowlengths);
1140:     PetscBinaryRead(fd,rowlengths,M,NULL,PETSC_INT);

1142:     a = (Mat_SeqDense*)newmat->data;
1143:     v = a->v;

1145:     /* read column indices and nonzeros */
1146:     PetscMalloc1(nz+1,&scols);
1147:     cols = scols;
1148:     PetscBinaryRead(fd,cols,nz,NULL,PETSC_INT);
1149:     PetscMalloc1(nz+1,&svals);
1150:     vals = svals;
1151:     PetscBinaryRead(fd,vals,nz,NULL,PETSC_SCALAR);

1153:     /* insert into matrix */
1154:     for (i=0; i<M; i++) {
1155:       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
1156:       svals += rowlengths[i]; scols += rowlengths[i];
1157:     }
1158:     PetscFree(vals);
1159:     PetscFree(cols);
1160:     PetscFree(rowlengths);
1161:   }
1162:   MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
1163:   MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);
1164:   return(0);
1165: }

1167: static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
1168: {
1169:   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
1170:   PetscErrorCode    ierr;
1171:   PetscInt          i,j;
1172:   const char        *name;
1173:   PetscScalar       *v;
1174:   PetscViewerFormat format;
1175: #if defined(PETSC_USE_COMPLEX)
1176:   PetscBool         allreal = PETSC_TRUE;
1177: #endif

1180:   PetscViewerGetFormat(viewer,&format);
1181:   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1182:     return(0);  /* do nothing for now */
1183:   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
1184:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
1185:     for (i=0; i<A->rmap->n; i++) {
1186:       v    = a->v + i;
1187:       PetscViewerASCIIPrintf(viewer,"row %D:",i);
1188:       for (j=0; j<A->cmap->n; j++) {
1189: #if defined(PETSC_USE_COMPLEX)
1190:         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
1191:           PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",j,(double)PetscRealPart(*v),(double)PetscImaginaryPart(*v));
1192:         } else if (PetscRealPart(*v)) {
1193:           PetscViewerASCIIPrintf(viewer," (%D, %g) ",j,(double)PetscRealPart(*v));
1194:         }
1195: #else
1196:         if (*v) {
1197:           PetscViewerASCIIPrintf(viewer," (%D, %g) ",j,(double)*v);
1198:         }
1199: #endif
1200:         v += a->lda;
1201:       }
1202:       PetscViewerASCIIPrintf(viewer,"\n");
1203:     }
1204:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
1205:   } else {
1206:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
1207: #if defined(PETSC_USE_COMPLEX)
1208:     /* determine if matrix has all real values */
1209:     v = a->v;
1210:     for (i=0; i<A->rmap->n*A->cmap->n; i++) {
1211:       if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break;}
1212:     }
1213: #endif
1214:     if (format == PETSC_VIEWER_ASCII_MATLAB) {
1215:       PetscObjectGetName((PetscObject)A,&name);
1216:       PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap->n,A->cmap->n);
1217:       PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap->n,A->cmap->n);
1218:       PetscViewerASCIIPrintf(viewer,"%s = [\n",name);
1219:     }

1221:     for (i=0; i<A->rmap->n; i++) {
1222:       v = a->v + i;
1223:       for (j=0; j<A->cmap->n; j++) {
1224: #if defined(PETSC_USE_COMPLEX)
1225:         if (allreal) {
1226:           PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)PetscRealPart(*v));
1227:         } else {
1228:           PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16ei ",(double)PetscRealPart(*v),(double)PetscImaginaryPart(*v));
1229:         }
1230: #else
1231:         PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)*v);
1232: #endif
1233:         v += a->lda;
1234:       }
1235:       PetscViewerASCIIPrintf(viewer,"\n");
1236:     }
1237:     if (format == PETSC_VIEWER_ASCII_MATLAB) {
1238:       PetscViewerASCIIPrintf(viewer,"];\n");
1239:     }
1240:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
1241:   }
1242:   PetscViewerFlush(viewer);
1243:   return(0);
1244: }

1246: static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
1247: {
1248:   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
1249:   PetscErrorCode    ierr;
1250:   int               fd;
1251:   PetscInt          ict,j,n = A->cmap->n,m = A->rmap->n,i,*col_lens,nz = m*n;
1252:   PetscScalar       *v,*anonz,*vals;
1253:   PetscViewerFormat format;

1256:   PetscViewerBinaryGetDescriptor(viewer,&fd);

1258:   PetscViewerGetFormat(viewer,&format);
1259:   if (format == PETSC_VIEWER_NATIVE) {
1260:     /* store the matrix as a dense matrix */
1261:     PetscMalloc1(4,&col_lens);

1263:     col_lens[0] = MAT_FILE_CLASSID;
1264:     col_lens[1] = m;
1265:     col_lens[2] = n;
1266:     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;

1268:     PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);
1269:     PetscFree(col_lens);

1271:     /* write out matrix, by rows */
1272:     PetscMalloc1(m*n+1,&vals);
1273:     v    = a->v;
1274:     for (j=0; j<n; j++) {
1275:       for (i=0; i<m; i++) {
1276:         vals[j + i*n] = *v++;
1277:       }
1278:     }
1279:     PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);
1280:     PetscFree(vals);
1281:   } else {
1282:     PetscMalloc1(4+nz,&col_lens);

1284:     col_lens[0] = MAT_FILE_CLASSID;
1285:     col_lens[1] = m;
1286:     col_lens[2] = n;
1287:     col_lens[3] = nz;

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

1293:     /* Possibly should write in smaller increments, not whole matrix at once? */
1294:     /* store column indices (zero start index) */
1295:     ict = 0;
1296:     for (i=0; i<m; i++) {
1297:       for (j=0; j<n; j++) col_lens[ict++] = j;
1298:     }
1299:     PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);
1300:     PetscFree(col_lens);

1302:     /* store nonzero values */
1303:     PetscMalloc1(nz+1,&anonz);
1304:     ict  = 0;
1305:     for (i=0; i<m; i++) {
1306:       v = a->v + i;
1307:       for (j=0; j<n; j++) {
1308:         anonz[ict++] = *v; v += a->lda;
1309:       }
1310:     }
1311:     PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);
1312:     PetscFree(anonz);
1313:   }
1314:   return(0);
1315: }

1317:  #include <petscdraw.h>
1318: static PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
1319: {
1320:   Mat               A  = (Mat) Aa;
1321:   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
1322:   PetscErrorCode    ierr;
1323:   PetscInt          m  = A->rmap->n,n = A->cmap->n,i,j;
1324:   int               color = PETSC_DRAW_WHITE;
1325:   PetscScalar       *v = a->v;
1326:   PetscViewer       viewer;
1327:   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r;
1328:   PetscViewerFormat format;

1331:   PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);
1332:   PetscViewerGetFormat(viewer,&format);
1333:   PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);

1335:   /* Loop over matrix elements drawing boxes */

1337:   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
1338:     PetscDrawCollectiveBegin(draw);
1339:     /* Blue for negative and Red for positive */
1340:     for (j = 0; j < n; j++) {
1341:       x_l = j; x_r = x_l + 1.0;
1342:       for (i = 0; i < m; i++) {
1343:         y_l = m - i - 1.0;
1344:         y_r = y_l + 1.0;
1345:         if (PetscRealPart(v[j*m+i]) >  0.) {
1346:           color = PETSC_DRAW_RED;
1347:         } else if (PetscRealPart(v[j*m+i]) <  0.) {
1348:           color = PETSC_DRAW_BLUE;
1349:         } else {
1350:           continue;
1351:         }
1352:         PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
1353:       }
1354:     }
1355:     PetscDrawCollectiveEnd(draw);
1356:   } else {
1357:     /* use contour shading to indicate magnitude of values */
1358:     /* first determine max of all nonzero values */
1359:     PetscReal minv = 0.0, maxv = 0.0;
1360:     PetscDraw popup;

1362:     for (i=0; i < m*n; i++) {
1363:       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
1364:     }
1365:     if (minv >= maxv) maxv = minv + PETSC_SMALL;
1366:     PetscDrawGetPopup(draw,&popup);
1367:     PetscDrawScalePopup(popup,minv,maxv);

1369:     PetscDrawCollectiveBegin(draw);
1370:     for (j=0; j<n; j++) {
1371:       x_l = j;
1372:       x_r = x_l + 1.0;
1373:       for (i=0; i<m; i++) {
1374:         y_l = m - i - 1.0;
1375:         y_r = y_l + 1.0;
1376:         color = PetscDrawRealToColor(PetscAbsScalar(v[j*m+i]),minv,maxv);
1377:         PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
1378:       }
1379:     }
1380:     PetscDrawCollectiveEnd(draw);
1381:   }
1382:   return(0);
1383: }

1385: static PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
1386: {
1387:   PetscDraw      draw;
1388:   PetscBool      isnull;
1389:   PetscReal      xr,yr,xl,yl,h,w;

1393:   PetscViewerDrawGetDraw(viewer,0,&draw);
1394:   PetscDrawIsNull(draw,&isnull);
1395:   if (isnull) return(0);

1397:   xr   = A->cmap->n; yr = A->rmap->n; h = yr/10.0; w = xr/10.0;
1398:   xr  += w;          yr += h;        xl = -w;     yl = -h;
1399:   PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
1400:   PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
1401:   PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);
1402:   PetscObjectCompose((PetscObject)A,"Zoomviewer",NULL);
1403:   PetscDrawSave(draw);
1404:   return(0);
1405: }

1407: PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer)
1408: {
1410:   PetscBool      iascii,isbinary,isdraw;

1413:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1414:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
1415:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);

1417:   if (iascii) {
1418:     MatView_SeqDense_ASCII(A,viewer);
1419:   } else if (isbinary) {
1420:     MatView_SeqDense_Binary(A,viewer);
1421:   } else if (isdraw) {
1422:     MatView_SeqDense_Draw(A,viewer);
1423:   }
1424:   return(0);
1425: }

1427: static PetscErrorCode MatDensePlaceArray_SeqDense(Mat A,const PetscScalar array[])
1428: {
1429:   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;

1432:   a->unplacedarray       = a->v;
1433:   a->unplaced_user_alloc = a->user_alloc;
1434:   a->v                   = (PetscScalar*) array;
1435:   return(0);
1436: }

1438: static PetscErrorCode MatDenseResetArray_SeqDense(Mat A)
1439: {
1440:   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;

1443:   a->v             = a->unplacedarray;
1444:   a->user_alloc    = a->unplaced_user_alloc;
1445:   a->unplacedarray = NULL;
1446:   return(0);
1447: }

1449: static PetscErrorCode MatDestroy_SeqDense(Mat mat)
1450: {
1451:   Mat_SeqDense   *l = (Mat_SeqDense*)mat->data;

1455: #if defined(PETSC_USE_LOG)
1456:   PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap->n,mat->cmap->n);
1457: #endif
1458:   PetscFree(l->pivots);
1459:   PetscFree(l->fwork);
1460:   MatDestroy(&l->ptapwork);
1461:   if (!l->user_alloc) {PetscFree(l->v);}
1462:   PetscFree(mat->data);

1464:   PetscObjectChangeTypeName((PetscObject)mat,0);
1465:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetLDA_C",NULL);
1466:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArray_C",NULL);
1467:   PetscObjectComposeFunction((PetscObject)mat,"MatDensePlaceArray_C",NULL);
1468:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseResetArray_C",NULL);
1469:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArray_C",NULL);
1470:   PetscObjectComposeFunction((PetscObject)mat,"MatConvert_seqdense_seqaij_C",NULL);
1471: #if defined(PETSC_HAVE_ELEMENTAL)
1472:   PetscObjectComposeFunction((PetscObject)mat,"MatConvert_seqdense_elemental_C",NULL);
1473: #endif
1474:   PetscObjectComposeFunction((PetscObject)mat,"MatSeqDenseSetPreallocation_C",NULL);
1475:   PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_seqaij_seqdense_C",NULL);
1476:   PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C",NULL);
1477:   PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C",NULL);
1478:   PetscObjectComposeFunction((PetscObject)mat,"MatPtAP_seqaij_seqdense_C",NULL);
1479:   PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMult_seqaij_seqdense_C",NULL);
1480:   PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultSymbolic_seqaij_seqdense_C",NULL);
1481:   PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultNumeric_seqaij_seqdense_C",NULL);
1482:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumn_C",NULL);
1483:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumn_C",NULL);
1484:   return(0);
1485: }

1487: static PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout)
1488: {
1489:   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1491:   PetscInt       k,j,m,n,M;
1492:   PetscScalar    *v,tmp;

1495:   v = mat->v; m = A->rmap->n; M = mat->lda; n = A->cmap->n;
1496:   if (reuse == MAT_INPLACE_MATRIX && m == n) { /* in place transpose */
1497:     for (j=0; j<m; j++) {
1498:       for (k=0; k<j; k++) {
1499:         tmp        = v[j + k*M];
1500:         v[j + k*M] = v[k + j*M];
1501:         v[k + j*M] = tmp;
1502:       }
1503:     }
1504:   } else { /* out-of-place transpose */
1505:     Mat          tmat;
1506:     Mat_SeqDense *tmatd;
1507:     PetscScalar  *v2;
1508:     PetscInt     M2;

1510:     if (reuse != MAT_REUSE_MATRIX) {
1511:       MatCreate(PetscObjectComm((PetscObject)A),&tmat);
1512:       MatSetSizes(tmat,A->cmap->n,A->rmap->n,A->cmap->n,A->rmap->n);
1513:       MatSetType(tmat,((PetscObject)A)->type_name);
1514:       MatSeqDenseSetPreallocation(tmat,NULL);
1515:     } else {
1516:       tmat = *matout;
1517:     }
1518:     tmatd = (Mat_SeqDense*)tmat->data;
1519:     v = mat->v; v2 = tmatd->v; M2 = tmatd->lda;
1520:     for (j=0; j<n; j++) {
1521:       for (k=0; k<m; k++) v2[j + k*M2] = v[k + j*M];
1522:     }
1523:     MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);
1524:     MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);
1525:     if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX) *matout = tmat;
1526:     else {
1527:       MatHeaderMerge(A,&tmat);
1528:     }
1529:   }
1530:   return(0);
1531: }

1533: static PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscBool  *flg)
1534: {
1535:   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1536:   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
1537:   PetscInt     i,j;
1538:   PetscScalar  *v1,*v2;

1541:   if (A1->rmap->n != A2->rmap->n) {*flg = PETSC_FALSE; return(0);}
1542:   if (A1->cmap->n != A2->cmap->n) {*flg = PETSC_FALSE; return(0);}
1543:   for (i=0; i<A1->rmap->n; i++) {
1544:     v1 = mat1->v+i; v2 = mat2->v+i;
1545:     for (j=0; j<A1->cmap->n; j++) {
1546:       if (*v1 != *v2) {*flg = PETSC_FALSE; return(0);}
1547:       v1 += mat1->lda; v2 += mat2->lda;
1548:     }
1549:   }
1550:   *flg = PETSC_TRUE;
1551:   return(0);
1552: }

1554: static PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v)
1555: {
1556:   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1558:   PetscInt       i,n,len;
1559:   PetscScalar    *x,zero = 0.0;

1562:   VecSet(v,zero);
1563:   VecGetSize(v,&n);
1564:   VecGetArray(v,&x);
1565:   len  = PetscMin(A->rmap->n,A->cmap->n);
1566:   if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
1567:   for (i=0; i<len; i++) {
1568:     x[i] = mat->v[i*mat->lda + i];
1569:   }
1570:   VecRestoreArray(v,&x);
1571:   return(0);
1572: }

1574: static PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1575: {
1576:   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
1577:   const PetscScalar *l,*r;
1578:   PetscScalar       x,*v;
1579:   PetscErrorCode    ierr;
1580:   PetscInt          i,j,m = A->rmap->n,n = A->cmap->n;

1583:   if (ll) {
1584:     VecGetSize(ll,&m);
1585:     VecGetArrayRead(ll,&l);
1586:     if (m != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1587:     for (i=0; i<m; i++) {
1588:       x = l[i];
1589:       v = mat->v + i;
1590:       for (j=0; j<n; j++) { (*v) *= x; v+= mat->lda;}
1591:     }
1592:     VecRestoreArrayRead(ll,&l);
1593:     PetscLogFlops(1.0*n*m);
1594:   }
1595:   if (rr) {
1596:     VecGetSize(rr,&n);
1597:     VecGetArrayRead(rr,&r);
1598:     if (n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1599:     for (i=0; i<n; i++) {
1600:       x = r[i];
1601:       v = mat->v + i*mat->lda;
1602:       for (j=0; j<m; j++) (*v++) *= x;
1603:     }
1604:     VecRestoreArrayRead(rr,&r);
1605:     PetscLogFlops(1.0*n*m);
1606:   }
1607:   return(0);
1608: }

1610: static PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1611: {
1612:   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1613:   PetscScalar    *v   = mat->v;
1614:   PetscReal      sum  = 0.0;
1615:   PetscInt       lda  =mat->lda,m=A->rmap->n,i,j;

1619:   if (type == NORM_FROBENIUS) {
1620:     if (lda>m) {
1621:       for (j=0; j<A->cmap->n; j++) {
1622:         v = mat->v+j*lda;
1623:         for (i=0; i<m; i++) {
1624:           sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1625:         }
1626:       }
1627:     } else {
1628: #if defined(PETSC_USE_REAL___FP16)
1629:       PetscBLASInt one = 1,cnt = A->cmap->n*A->rmap->n;
1630:       *nrm = BLASnrm2_(&cnt,v,&one);
1631:     }
1632: #else
1633:       for (i=0; i<A->cmap->n*A->rmap->n; i++) {
1634:         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1635:       }
1636:     }
1637:     *nrm = PetscSqrtReal(sum);
1638: #endif
1639:     PetscLogFlops(2.0*A->cmap->n*A->rmap->n);
1640:   } else if (type == NORM_1) {
1641:     *nrm = 0.0;
1642:     for (j=0; j<A->cmap->n; j++) {
1643:       v   = mat->v + j*mat->lda;
1644:       sum = 0.0;
1645:       for (i=0; i<A->rmap->n; i++) {
1646:         sum += PetscAbsScalar(*v);  v++;
1647:       }
1648:       if (sum > *nrm) *nrm = sum;
1649:     }
1650:     PetscLogFlops(1.0*A->cmap->n*A->rmap->n);
1651:   } else if (type == NORM_INFINITY) {
1652:     *nrm = 0.0;
1653:     for (j=0; j<A->rmap->n; j++) {
1654:       v   = mat->v + j;
1655:       sum = 0.0;
1656:       for (i=0; i<A->cmap->n; i++) {
1657:         sum += PetscAbsScalar(*v); v += mat->lda;
1658:       }
1659:       if (sum > *nrm) *nrm = sum;
1660:     }
1661:     PetscLogFlops(1.0*A->cmap->n*A->rmap->n);
1662:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No two norm");
1663:   return(0);
1664: }

1666: static PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscBool flg)
1667: {
1668:   Mat_SeqDense   *aij = (Mat_SeqDense*)A->data;

1672:   switch (op) {
1673:   case MAT_ROW_ORIENTED:
1674:     aij->roworiented = flg;
1675:     break;
1676:   case MAT_NEW_NONZERO_LOCATIONS:
1677:   case MAT_NEW_NONZERO_LOCATION_ERR:
1678:   case MAT_NEW_NONZERO_ALLOCATION_ERR:
1679:   case MAT_NEW_DIAGONALS:
1680:   case MAT_KEEP_NONZERO_PATTERN:
1681:   case MAT_IGNORE_OFF_PROC_ENTRIES:
1682:   case MAT_USE_HASH_TABLE:
1683:   case MAT_IGNORE_ZERO_ENTRIES:
1684:   case MAT_IGNORE_LOWER_TRIANGULAR:
1685:     PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1686:     break;
1687:   case MAT_SPD:
1688:   case MAT_SYMMETRIC:
1689:   case MAT_STRUCTURALLY_SYMMETRIC:
1690:   case MAT_HERMITIAN:
1691:   case MAT_SYMMETRY_ETERNAL:
1692:     /* These options are handled directly by MatSetOption() */
1693:     break;
1694:   default:
1695:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
1696:   }
1697:   return(0);
1698: }

1700: static PetscErrorCode MatZeroEntries_SeqDense(Mat A)
1701: {
1702:   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
1704:   PetscInt       lda=l->lda,m=A->rmap->n,j;

1707:   if (lda>m) {
1708:     for (j=0; j<A->cmap->n; j++) {
1709:       PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));
1710:     }
1711:   } else {
1712:     PetscMemzero(l->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));
1713:   }
1714:   return(0);
1715: }

1717: static PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
1718: {
1719:   PetscErrorCode    ierr;
1720:   Mat_SeqDense      *l = (Mat_SeqDense*)A->data;
1721:   PetscInt          m  = l->lda, n = A->cmap->n, i,j;
1722:   PetscScalar       *slot,*bb;
1723:   const PetscScalar *xx;

1726: #if defined(PETSC_USE_DEBUG)
1727:   for (i=0; i<N; i++) {
1728:     if (rows[i] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row requested to be zeroed");
1729:     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);
1730:   }
1731: #endif

1733:   /* fix right hand side if needed */
1734:   if (x && b) {
1735:     VecGetArrayRead(x,&xx);
1736:     VecGetArray(b,&bb);
1737:     for (i=0; i<N; i++) bb[rows[i]] = diag*xx[rows[i]];
1738:     VecRestoreArrayRead(x,&xx);
1739:     VecRestoreArray(b,&bb);
1740:   }

1742:   for (i=0; i<N; i++) {
1743:     slot = l->v + rows[i];
1744:     for (j=0; j<n; j++) { *slot = 0.0; slot += m;}
1745:   }
1746:   if (diag != 0.0) {
1747:     if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only coded for square matrices");
1748:     for (i=0; i<N; i++) {
1749:       slot  = l->v + (m+1)*rows[i];
1750:       *slot = diag;
1751:     }
1752:   }
1753:   return(0);
1754: }

1756: static PetscErrorCode MatDenseGetLDA_SeqDense(Mat A,PetscInt *lda)
1757: {
1758:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;

1761:   *lda = mat->lda;
1762:   return(0);
1763: }

1765: static PetscErrorCode MatDenseGetArray_SeqDense(Mat A,PetscScalar *array[])
1766: {
1767:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;

1770:   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");
1771:   *array = mat->v;
1772:   return(0);
1773: }

1775: static PetscErrorCode MatDenseRestoreArray_SeqDense(Mat A,PetscScalar *array[])
1776: {
1778:   return(0);
1779: }

1781: /*@C
1782:    MatDenseGetLDA - gets the leading dimension of the array returned from MatDenseGetArray()

1784:    Logically Collective on Mat

1786:    Input Parameter:
1787: .  mat - a MATSEQDENSE or MATMPIDENSE matrix

1789:    Output Parameter:
1790: .   lda - the leading dimension

1792:    Level: intermediate

1794: .seealso: MatDenseGetArray(), MatDenseRestoreArray(), MatDenseGetArrayRead(), MatDenseRestoreArrayRead(), MatSeqDenseSetLDA()
1795: @*/
1796: PetscErrorCode  MatDenseGetLDA(Mat A,PetscInt *lda)
1797: {

1801:   PetscUseMethod(A,"MatDenseGetLDA_C",(Mat,PetscInt*),(A,lda));
1802:   return(0);
1803: }

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

1808:    Logically Collective on Mat

1810:    Input Parameter:
1811: .  mat - a MATSEQDENSE or MATMPIDENSE matrix

1813:    Output Parameter:
1814: .   array - pointer to the data

1816:    Level: intermediate

1818: .seealso: MatDenseRestoreArray(), MatDenseGetArrayRead(), MatDenseRestoreArrayRead()
1819: @*/
1820: PetscErrorCode  MatDenseGetArray(Mat A,PetscScalar **array)
1821: {

1825:   PetscUseMethod(A,"MatDenseGetArray_C",(Mat,PetscScalar**),(A,array));
1826:   return(0);
1827: }

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

1832:    Logically Collective on Mat

1834:    Input Parameters:
1835: .  mat - a MATSEQDENSE or MATMPIDENSE matrix
1836: .  array - pointer to the data

1838:    Level: intermediate

1840: .seealso: MatDenseGetArray(), MatDenseGetArrayRead(), MatDenseRestoreArrayRead()
1841: @*/
1842: PetscErrorCode  MatDenseRestoreArray(Mat A,PetscScalar **array)
1843: {

1847:   PetscUseMethod(A,"MatDenseRestoreArray_C",(Mat,PetscScalar**),(A,array));
1848:   if (array) *array = NULL;
1849:   PetscObjectStateIncrease((PetscObject)A);
1850:   return(0);
1851: }

1853: /*@C
1854:    MatDenseGetArrayRead - gives access to the array where the data for a SeqDense matrix is stored

1856:    Not Collective

1858:    Input Parameter:
1859: .  mat - a MATSEQDENSE or MATMPIDENSE matrix

1861:    Output Parameter:
1862: .   array - pointer to the data

1864:    Level: intermediate

1866: .seealso: MatDenseRestoreArray(), MatDenseGetArray(), MatDenseRestoreArrayRead()
1867: @*/
1868: PetscErrorCode  MatDenseGetArrayRead(Mat A,const PetscScalar **array)
1869: {

1873:   PetscUseMethod(A,"MatDenseGetArrayRead_C",(Mat,const PetscScalar**),(A,array));
1874:   return(0);
1875: }

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

1880:    Not Collective

1882:    Input Parameters:
1883: .  mat - a MATSEQDENSE or MATMPIDENSE matrix
1884: .  array - pointer to the data

1886:    Level: intermediate

1888: .seealso: MatDenseGetArray(), MatDenseGetArrayRead(), MatDenseRestoreArray()
1889: @*/
1890: PetscErrorCode  MatDenseRestoreArrayRead(Mat A,const PetscScalar **array)
1891: {

1895:   PetscUseMethod(A,"MatDenseRestoreArrayRead_C",(Mat,const PetscScalar**),(A,array));
1896:   if (array) *array = NULL;
1897:   return(0);
1898: }

1900: static PetscErrorCode MatCreateSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
1901: {
1902:   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1904:   PetscInt       i,j,nrows,ncols;
1905:   const PetscInt *irow,*icol;
1906:   PetscScalar    *av,*bv,*v = mat->v;
1907:   Mat            newmat;

1910:   ISGetIndices(isrow,&irow);
1911:   ISGetIndices(iscol,&icol);
1912:   ISGetLocalSize(isrow,&nrows);
1913:   ISGetLocalSize(iscol,&ncols);

1915:   /* Check submatrixcall */
1916:   if (scall == MAT_REUSE_MATRIX) {
1917:     PetscInt n_cols,n_rows;
1918:     MatGetSize(*B,&n_rows,&n_cols);
1919:     if (n_rows != nrows || n_cols != ncols) {
1920:       /* resize the result matrix to match number of requested rows/columns */
1921:       MatSetSizes(*B,nrows,ncols,nrows,ncols);
1922:     }
1923:     newmat = *B;
1924:   } else {
1925:     /* Create and fill new matrix */
1926:     MatCreate(PetscObjectComm((PetscObject)A),&newmat);
1927:     MatSetSizes(newmat,nrows,ncols,nrows,ncols);
1928:     MatSetType(newmat,((PetscObject)A)->type_name);
1929:     MatSeqDenseSetPreallocation(newmat,NULL);
1930:   }

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

1935:   for (i=0; i<ncols; i++) {
1936:     av = v + mat->lda*icol[i];
1937:     for (j=0; j<nrows; j++) *bv++ = av[irow[j]];
1938:   }

1940:   /* Assemble the matrices so that the correct flags are set */
1941:   MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
1942:   MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);

1944:   /* Free work space */
1945:   ISRestoreIndices(isrow,&irow);
1946:   ISRestoreIndices(iscol,&icol);
1947:   *B   = newmat;
1948:   return(0);
1949: }

1951: static PetscErrorCode MatCreateSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1952: {
1954:   PetscInt       i;

1957:   if (scall == MAT_INITIAL_MATRIX) {
1958:     PetscCalloc1(n+1,B);
1959:   }

1961:   for (i=0; i<n; i++) {
1962:     MatCreateSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);
1963:   }
1964:   return(0);
1965: }

1967: static PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode)
1968: {
1970:   return(0);
1971: }

1973: static PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode)
1974: {
1976:   return(0);
1977: }

1979: static PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
1980: {
1981:   Mat_SeqDense   *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense*)B->data;
1983:   PetscInt       lda1=a->lda,lda2=b->lda, m=A->rmap->n,n=A->cmap->n, j;

1986:   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
1987:   if (A->ops->copy != B->ops->copy) {
1988:     MatCopy_Basic(A,B,str);
1989:     return(0);
1990:   }
1991:   if (m != B->rmap->n || n != B->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1992:   if (lda1>m || lda2>m) {
1993:     for (j=0; j<n; j++) {
1994:       PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));
1995:     }
1996:   } else {
1997:     PetscMemcpy(b->v,a->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));
1998:   }
1999:   PetscObjectStateIncrease((PetscObject)B);
2000:   return(0);
2001: }

2003: static PetscErrorCode MatSetUp_SeqDense(Mat A)
2004: {

2008:    MatSeqDenseSetPreallocation(A,0);
2009:   return(0);
2010: }

2012: static PetscErrorCode MatConjugate_SeqDense(Mat A)
2013: {
2014:   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2015:   PetscInt     i,nz = A->rmap->n*A->cmap->n;
2016:   PetscScalar  *aa = a->v;

2019:   for (i=0; i<nz; i++) aa[i] = PetscConj(aa[i]);
2020:   return(0);
2021: }

2023: static PetscErrorCode MatRealPart_SeqDense(Mat A)
2024: {
2025:   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2026:   PetscInt     i,nz = A->rmap->n*A->cmap->n;
2027:   PetscScalar  *aa = a->v;

2030:   for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
2031:   return(0);
2032: }

2034: static PetscErrorCode MatImaginaryPart_SeqDense(Mat A)
2035: {
2036:   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2037:   PetscInt     i,nz = A->rmap->n*A->cmap->n;
2038:   PetscScalar  *aa = a->v;

2041:   for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
2042:   return(0);
2043: }

2045: /* ----------------------------------------------------------------*/
2046: PETSC_INTERN PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
2047: {

2051:   if (scall == MAT_INITIAL_MATRIX) {
2052:     PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);
2053:     MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);
2054:     PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);
2055:   }
2056:   PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);
2057:   MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);
2058:   PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);
2059:   return(0);
2060: }

2062: PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
2063: {
2065:   PetscInt       m=A->rmap->n,n=B->cmap->n;
2066:   Mat            Cmat;

2069:   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);
2070:   MatCreate(PETSC_COMM_SELF,&Cmat);
2071:   MatSetSizes(Cmat,m,n,m,n);
2072:   MatSetType(Cmat,MATSEQDENSE);
2073:   MatSeqDenseSetPreallocation(Cmat,NULL);

2075:   *C = Cmat;
2076:   return(0);
2077: }

2079: PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
2080: {
2081:   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
2082:   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
2083:   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
2084:   PetscBLASInt   m,n,k;
2085:   PetscScalar    _DOne=1.0,_DZero=0.0;
2087:   PetscBool      flg;

2090:   PetscObjectTypeCompare((PetscObject)B,MATSEQDENSE,&flg);
2091:   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Second matrix must be dense");

2093:   /* Handle case where where user provided the final C matrix rather than calling MatMatMult() with MAT_INITIAL_MATRIX*/
2094:   PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);
2095:   if (flg) {
2096:     C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqDense;
2097:     (*C->ops->matmultnumeric)(A,B,C);
2098:     return(0);
2099:   }

2101:   PetscObjectTypeCompare((PetscObject)A,MATSEQDENSE,&flg);
2102:   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"First matrix must be dense");
2103:   PetscBLASIntCast(C->rmap->n,&m);
2104:   PetscBLASIntCast(C->cmap->n,&n);
2105:   PetscBLASIntCast(A->cmap->n,&k);
2106:   if (!m || !n || !k) return(0);
2107:   PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda));
2108:   return(0);
2109: }

2111: PetscErrorCode MatMatTransposeMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
2112: {

2116:   if (scall == MAT_INITIAL_MATRIX) {
2117:     PetscLogEventBegin(MAT_MatTransposeMultSymbolic,A,B,0,0);
2118:     MatMatTransposeMultSymbolic_SeqDense_SeqDense(A,B,fill,C);
2119:     PetscLogEventEnd(MAT_MatTransposeMultSymbolic,A,B,0,0);
2120:   }
2121:   PetscLogEventBegin(MAT_MatTransposeMultNumeric,A,B,0,0);
2122:   MatMatTransposeMultNumeric_SeqDense_SeqDense(A,B,*C);
2123:   PetscLogEventEnd(MAT_MatTransposeMultNumeric,A,B,0,0);
2124:   return(0);
2125: }

2127: PetscErrorCode MatMatTransposeMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
2128: {
2130:   PetscInt       m=A->rmap->n,n=B->rmap->n;
2131:   Mat            Cmat;

2134:   if (A->cmap->n != B->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"A->cmap->n %d != B->cmap->n %d\n",A->cmap->n,B->cmap->n);
2135:   MatCreate(PETSC_COMM_SELF,&Cmat);
2136:   MatSetSizes(Cmat,m,n,m,n);
2137:   MatSetType(Cmat,MATSEQDENSE);
2138:   MatSeqDenseSetPreallocation(Cmat,NULL);

2140:   Cmat->assembled = PETSC_TRUE;

2142:   *C = Cmat;
2143:   return(0);
2144: }

2146: PetscErrorCode MatMatTransposeMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
2147: {
2148:   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
2149:   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
2150:   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
2151:   PetscBLASInt   m,n,k;
2152:   PetscScalar    _DOne=1.0,_DZero=0.0;

2156:   PetscBLASIntCast(C->rmap->n,&m);
2157:   PetscBLASIntCast(C->cmap->n,&n);
2158:   PetscBLASIntCast(A->cmap->n,&k);
2159:   if (!m || !n || !k) return(0);
2160:   PetscStackCallBLAS("BLASgemm",BLASgemm_("N","T",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda));
2161:   return(0);
2162: }

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

2169:   if (scall == MAT_INITIAL_MATRIX) {
2170:     PetscLogEventBegin(MAT_TransposeMatMultSymbolic,A,B,0,0);
2171:     MatTransposeMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);
2172:     PetscLogEventEnd(MAT_TransposeMatMultSymbolic,A,B,0,0);
2173:   }
2174:   PetscLogEventBegin(MAT_TransposeMatMultNumeric,A,B,0,0);
2175:   MatTransposeMatMultNumeric_SeqDense_SeqDense(A,B,*C);
2176:   PetscLogEventEnd(MAT_TransposeMatMultNumeric,A,B,0,0);
2177:   return(0);
2178: }

2180: PetscErrorCode MatTransposeMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
2181: {
2183:   PetscInt       m=A->cmap->n,n=B->cmap->n;
2184:   Mat            Cmat;

2187:   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);
2188:   MatCreate(PETSC_COMM_SELF,&Cmat);
2189:   MatSetSizes(Cmat,m,n,m,n);
2190:   MatSetType(Cmat,MATSEQDENSE);
2191:   MatSeqDenseSetPreallocation(Cmat,NULL);

2193:   Cmat->assembled = PETSC_TRUE;

2195:   *C = Cmat;
2196:   return(0);
2197: }

2199: PetscErrorCode MatTransposeMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
2200: {
2201:   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
2202:   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
2203:   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
2204:   PetscBLASInt   m,n,k;
2205:   PetscScalar    _DOne=1.0,_DZero=0.0;

2209:   PetscBLASIntCast(C->rmap->n,&m);
2210:   PetscBLASIntCast(C->cmap->n,&n);
2211:   PetscBLASIntCast(A->rmap->n,&k);
2212:   if (!m || !n || !k) return(0);
2213:   PetscStackCallBLAS("BLASgemm",BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda));
2214:   return(0);
2215: }

2217: static PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[])
2218: {
2219:   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
2221:   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
2222:   PetscScalar    *x;
2223:   MatScalar      *aa = a->v;

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

2228:   VecSet(v,0.0);
2229:   VecGetArray(v,&x);
2230:   VecGetLocalSize(v,&p);
2231:   if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2232:   for (i=0; i<m; i++) {
2233:     x[i] = aa[i]; if (idx) idx[i] = 0;
2234:     for (j=1; j<n; j++) {
2235:       if (PetscRealPart(x[i]) < PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
2236:     }
2237:   }
2238:   VecRestoreArray(v,&x);
2239:   return(0);
2240: }

2242: static PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[])
2243: {
2244:   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
2246:   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
2247:   PetscScalar    *x;
2248:   PetscReal      atmp;
2249:   MatScalar      *aa = a->v;

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

2254:   VecSet(v,0.0);
2255:   VecGetArray(v,&x);
2256:   VecGetLocalSize(v,&p);
2257:   if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2258:   for (i=0; i<m; i++) {
2259:     x[i] = PetscAbsScalar(aa[i]);
2260:     for (j=1; j<n; j++) {
2261:       atmp = PetscAbsScalar(aa[i+m*j]);
2262:       if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;}
2263:     }
2264:   }
2265:   VecRestoreArray(v,&x);
2266:   return(0);
2267: }

2269: static PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[])
2270: {
2271:   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
2273:   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
2274:   PetscScalar    *x;
2275:   MatScalar      *aa = a->v;

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

2280:   VecSet(v,0.0);
2281:   VecGetArray(v,&x);
2282:   VecGetLocalSize(v,&p);
2283:   if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2284:   for (i=0; i<m; i++) {
2285:     x[i] = aa[i]; if (idx) idx[i] = 0;
2286:     for (j=1; j<n; j++) {
2287:       if (PetscRealPart(x[i]) > PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
2288:     }
2289:   }
2290:   VecRestoreArray(v,&x);
2291:   return(0);
2292: }

2294: static PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col)
2295: {
2296:   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
2298:   PetscScalar    *x;

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

2303:   VecGetArray(v,&x);
2304:   PetscMemcpy(x,a->v+col*a->lda,A->rmap->n*sizeof(PetscScalar));
2305:   VecRestoreArray(v,&x);
2306:   return(0);
2307: }

2309: PetscErrorCode MatGetColumnNorms_SeqDense(Mat A,NormType type,PetscReal *norms)
2310: {
2311:   PetscErrorCode    ierr;
2312:   PetscInt          i,j,m,n;
2313:   const PetscScalar *a;

2316:   MatGetSize(A,&m,&n);
2317:   PetscMemzero(norms,n*sizeof(PetscReal));
2318:   MatDenseGetArrayRead(A,&a);
2319:   if (type == NORM_2) {
2320:     for (i=0; i<n; i++) {
2321:       for (j=0; j<m; j++) {
2322:         norms[i] += PetscAbsScalar(a[j]*a[j]);
2323:       }
2324:       a += m;
2325:     }
2326:   } else if (type == NORM_1) {
2327:     for (i=0; i<n; i++) {
2328:       for (j=0; j<m; j++) {
2329:         norms[i] += PetscAbsScalar(a[j]);
2330:       }
2331:       a += m;
2332:     }
2333:   } else if (type == NORM_INFINITY) {
2334:     for (i=0; i<n; i++) {
2335:       for (j=0; j<m; j++) {
2336:         norms[i] = PetscMax(PetscAbsScalar(a[j]),norms[i]);
2337:       }
2338:       a += m;
2339:     }
2340:   } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Unknown NormType");
2341:   MatDenseRestoreArrayRead(A,&a);
2342:   if (type == NORM_2) {
2343:     for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]);
2344:   }
2345:   return(0);
2346: }

2348: static PetscErrorCode  MatSetRandom_SeqDense(Mat x,PetscRandom rctx)
2349: {
2351:   PetscScalar    *a;
2352:   PetscInt       m,n,i;

2355:   MatGetSize(x,&m,&n);
2356:   MatDenseGetArray(x,&a);
2357:   for (i=0; i<m*n; i++) {
2358:     PetscRandomGetValue(rctx,a+i);
2359:   }
2360:   MatDenseRestoreArray(x,&a);
2361:   return(0);
2362: }

2364: static PetscErrorCode MatMissingDiagonal_SeqDense(Mat A,PetscBool  *missing,PetscInt *d)
2365: {
2367:   *missing = PETSC_FALSE;
2368:   return(0);
2369: }

2371: static PetscErrorCode MatDenseGetColumn_SeqDense(Mat A,PetscInt col,PetscScalar **vals)
2372: {
2373:   Mat_SeqDense *a = (Mat_SeqDense*)A->data;

2376:   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2377:   *vals = a->v+col*a->lda;
2378:   return(0);
2379: }

2381: static PetscErrorCode MatDenseRestoreColumn_SeqDense(Mat A,PetscScalar **vals)
2382: {
2384:   *vals = 0; /* user cannot accidently use the array later */
2385:   return(0);
2386: }

2388: /* -------------------------------------------------------------------*/
2389: static struct _MatOps MatOps_Values = { MatSetValues_SeqDense,
2390:                                         MatGetRow_SeqDense,
2391:                                         MatRestoreRow_SeqDense,
2392:                                         MatMult_SeqDense,
2393:                                 /*  4*/ MatMultAdd_SeqDense,
2394:                                         MatMultTranspose_SeqDense,
2395:                                         MatMultTransposeAdd_SeqDense,
2396:                                         0,
2397:                                         0,
2398:                                         0,
2399:                                 /* 10*/ 0,
2400:                                         MatLUFactor_SeqDense,
2401:                                         MatCholeskyFactor_SeqDense,
2402:                                         MatSOR_SeqDense,
2403:                                         MatTranspose_SeqDense,
2404:                                 /* 15*/ MatGetInfo_SeqDense,
2405:                                         MatEqual_SeqDense,
2406:                                         MatGetDiagonal_SeqDense,
2407:                                         MatDiagonalScale_SeqDense,
2408:                                         MatNorm_SeqDense,
2409:                                 /* 20*/ MatAssemblyBegin_SeqDense,
2410:                                         MatAssemblyEnd_SeqDense,
2411:                                         MatSetOption_SeqDense,
2412:                                         MatZeroEntries_SeqDense,
2413:                                 /* 24*/ MatZeroRows_SeqDense,
2414:                                         0,
2415:                                         0,
2416:                                         0,
2417:                                         0,
2418:                                 /* 29*/ MatSetUp_SeqDense,
2419:                                         0,
2420:                                         0,
2421:                                         0,
2422:                                         0,
2423:                                 /* 34*/ MatDuplicate_SeqDense,
2424:                                         0,
2425:                                         0,
2426:                                         0,
2427:                                         0,
2428:                                 /* 39*/ MatAXPY_SeqDense,
2429:                                         MatCreateSubMatrices_SeqDense,
2430:                                         0,
2431:                                         MatGetValues_SeqDense,
2432:                                         MatCopy_SeqDense,
2433:                                 /* 44*/ MatGetRowMax_SeqDense,
2434:                                         MatScale_SeqDense,
2435:                                         MatShift_Basic,
2436:                                         0,
2437:                                         MatZeroRowsColumns_SeqDense,
2438:                                 /* 49*/ MatSetRandom_SeqDense,
2439:                                         0,
2440:                                         0,
2441:                                         0,
2442:                                         0,
2443:                                 /* 54*/ 0,
2444:                                         0,
2445:                                         0,
2446:                                         0,
2447:                                         0,
2448:                                 /* 59*/ 0,
2449:                                         MatDestroy_SeqDense,
2450:                                         MatView_SeqDense,
2451:                                         0,
2452:                                         0,
2453:                                 /* 64*/ 0,
2454:                                         0,
2455:                                         0,
2456:                                         0,
2457:                                         0,
2458:                                 /* 69*/ MatGetRowMaxAbs_SeqDense,
2459:                                         0,
2460:                                         0,
2461:                                         0,
2462:                                         0,
2463:                                 /* 74*/ 0,
2464:                                         0,
2465:                                         0,
2466:                                         0,
2467:                                         0,
2468:                                 /* 79*/ 0,
2469:                                         0,
2470:                                         0,
2471:                                         0,
2472:                                 /* 83*/ MatLoad_SeqDense,
2473:                                         0,
2474:                                         MatIsHermitian_SeqDense,
2475:                                         0,
2476:                                         0,
2477:                                         0,
2478:                                 /* 89*/ MatMatMult_SeqDense_SeqDense,
2479:                                         MatMatMultSymbolic_SeqDense_SeqDense,
2480:                                         MatMatMultNumeric_SeqDense_SeqDense,
2481:                                         MatPtAP_SeqDense_SeqDense,
2482:                                         MatPtAPSymbolic_SeqDense_SeqDense,
2483:                                 /* 94*/ MatPtAPNumeric_SeqDense_SeqDense,
2484:                                         MatMatTransposeMult_SeqDense_SeqDense,
2485:                                         MatMatTransposeMultSymbolic_SeqDense_SeqDense,
2486:                                         MatMatTransposeMultNumeric_SeqDense_SeqDense,
2487:                                         0,
2488:                                 /* 99*/ 0,
2489:                                         0,
2490:                                         0,
2491:                                         MatConjugate_SeqDense,
2492:                                         0,
2493:                                 /*104*/ 0,
2494:                                         MatRealPart_SeqDense,
2495:                                         MatImaginaryPart_SeqDense,
2496:                                         0,
2497:                                         0,
2498:                                 /*109*/ 0,
2499:                                         0,
2500:                                         MatGetRowMin_SeqDense,
2501:                                         MatGetColumnVector_SeqDense,
2502:                                         MatMissingDiagonal_SeqDense,
2503:                                 /*114*/ 0,
2504:                                         0,
2505:                                         0,
2506:                                         0,
2507:                                         0,
2508:                                 /*119*/ 0,
2509:                                         0,
2510:                                         0,
2511:                                         0,
2512:                                         0,
2513:                                 /*124*/ 0,
2514:                                         MatGetColumnNorms_SeqDense,
2515:                                         0,
2516:                                         0,
2517:                                         0,
2518:                                 /*129*/ 0,
2519:                                         MatTransposeMatMult_SeqDense_SeqDense,
2520:                                         MatTransposeMatMultSymbolic_SeqDense_SeqDense,
2521:                                         MatTransposeMatMultNumeric_SeqDense_SeqDense,
2522:                                         0,
2523:                                 /*134*/ 0,
2524:                                         0,
2525:                                         0,
2526:                                         0,
2527:                                         0,
2528:                                 /*139*/ 0,
2529:                                         0,
2530:                                         0,
2531:                                         0,
2532:                                         0,
2533:                                 /*144*/ MatCreateMPIMatConcatenateSeqMat_SeqDense
2534: };

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

2541:    Collective on MPI_Comm

2543:    Input Parameters:
2544: +  comm - MPI communicator, set to PETSC_COMM_SELF
2545: .  m - number of rows
2546: .  n - number of columns
2547: -  data - optional location of matrix data in column major order.  Set data=NULL for PETSc
2548:    to control all matrix memory allocation.

2550:    Output Parameter:
2551: .  A - the matrix

2553:    Notes:
2554:    The data input variable is intended primarily for Fortran programmers
2555:    who wish to allocate their own matrix memory space.  Most users should
2556:    set data=NULL.

2558:    Level: intermediate

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

2562: .seealso: MatCreate(), MatCreateDense(), MatSetValues()
2563: @*/
2564: PetscErrorCode  MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A)
2565: {

2569:   MatCreate(comm,A);
2570:   MatSetSizes(*A,m,n,m,n);
2571:   MatSetType(*A,MATSEQDENSE);
2572:   MatSeqDenseSetPreallocation(*A,data);
2573:   return(0);
2574: }

2576: /*@C
2577:    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements

2579:    Collective on MPI_Comm

2581:    Input Parameters:
2582: +  B - the matrix
2583: -  data - the array (or NULL)

2585:    Notes:
2586:    The data input variable is intended primarily for Fortran programmers
2587:    who wish to allocate their own matrix memory space.  Most users should
2588:    need not call this routine.

2590:    Level: intermediate

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

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

2596: @*/
2597: PetscErrorCode  MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
2598: {

2602:   PetscTryMethod(B,"MatSeqDenseSetPreallocation_C",(Mat,PetscScalar[]),(B,data));
2603:   return(0);
2604: }

2606: PetscErrorCode  MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
2607: {
2608:   Mat_SeqDense   *b;

2612:   B->preallocated = PETSC_TRUE;

2614:   PetscLayoutSetUp(B->rmap);
2615:   PetscLayoutSetUp(B->cmap);

2617:   b       = (Mat_SeqDense*)B->data;
2618:   b->Mmax = B->rmap->n;
2619:   b->Nmax = B->cmap->n;
2620:   if (b->lda <= 0 || b->changelda) b->lda = B->rmap->n;

2622:   PetscIntMultError(b->lda,b->Nmax,NULL);
2623:   if (!data) { /* petsc-allocated storage */
2624:     if (!b->user_alloc) { PetscFree(b->v); }
2625:     PetscCalloc1((size_t)b->lda*b->Nmax,&b->v);
2626:     PetscLogObjectMemory((PetscObject)B,b->lda*b->Nmax*sizeof(PetscScalar));

2628:     b->user_alloc = PETSC_FALSE;
2629:   } else { /* user-allocated storage */
2630:     if (!b->user_alloc) { PetscFree(b->v); }
2631:     b->v          = data;
2632:     b->user_alloc = PETSC_TRUE;
2633:   }
2634:   B->assembled = PETSC_TRUE;
2635:   return(0);
2636: }

2638: #if defined(PETSC_HAVE_ELEMENTAL)
2639: PETSC_INTERN PetscErrorCode MatConvert_SeqDense_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
2640: {
2641:   Mat               mat_elemental;
2642:   PetscErrorCode    ierr;
2643:   const PetscScalar *array;
2644:   PetscScalar       *v_colwise;
2645:   PetscInt          M=A->rmap->N,N=A->cmap->N,i,j,k,*rows,*cols;

2648:   PetscMalloc3(M*N,&v_colwise,M,&rows,N,&cols);
2649:   MatDenseGetArrayRead(A,&array);
2650:   /* convert column-wise array into row-wise v_colwise, see MatSetValues_Elemental() */
2651:   k = 0;
2652:   for (j=0; j<N; j++) {
2653:     cols[j] = j;
2654:     for (i=0; i<M; i++) {
2655:       v_colwise[j*M+i] = array[k++];
2656:     }
2657:   }
2658:   for (i=0; i<M; i++) {
2659:     rows[i] = i;
2660:   }
2661:   MatDenseRestoreArrayRead(A,&array);

2663:   MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);
2664:   MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,M,N);
2665:   MatSetType(mat_elemental,MATELEMENTAL);
2666:   MatSetUp(mat_elemental);

2668:   /* PETSc-Elemental interaface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
2669:   MatSetValues(mat_elemental,M,rows,N,cols,v_colwise,ADD_VALUES);
2670:   MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);
2671:   MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);
2672:   PetscFree3(v_colwise,rows,cols);

2674:   if (reuse == MAT_INPLACE_MATRIX) {
2675:     MatHeaderReplace(A,&mat_elemental);
2676:   } else {
2677:     *newmat = mat_elemental;
2678:   }
2679:   return(0);
2680: }
2681: #endif

2683: /*@C
2684:   MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array

2686:   Input parameter:
2687: + A - the matrix
2688: - lda - the leading dimension

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

2695:   Level: intermediate

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

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

2701: @*/
2702: PetscErrorCode  MatSeqDenseSetLDA(Mat B,PetscInt lda)
2703: {
2704:   Mat_SeqDense *b = (Mat_SeqDense*)B->data;

2707:   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);
2708:   b->lda       = lda;
2709:   b->changelda = PETSC_FALSE;
2710:   b->Mmax      = PetscMax(b->Mmax,lda);
2711:   return(0);
2712: }

2714: PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqDense(MPI_Comm comm,Mat inmat,PetscInt n,MatReuse scall,Mat *outmat)
2715: {
2717:   PetscMPIInt    size;

2720:   MPI_Comm_size(comm,&size);
2721:   if (size == 1) {
2722:     if (scall == MAT_INITIAL_MATRIX) {
2723:       MatDuplicate(inmat,MAT_COPY_VALUES,outmat);
2724:     } else {
2725:       MatCopy(inmat,*outmat,SAME_NONZERO_PATTERN);
2726:     }
2727:   } else {
2728:     MatCreateMPIMatConcatenateSeqMat_MPIDense(comm,inmat,n,scall,outmat);
2729:   }
2730:   return(0);
2731: }

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

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

2739:   Level: beginner

2741: .seealso: MatCreateSeqDense()

2743: M*/

2745: PETSC_EXTERN PetscErrorCode MatCreate_SeqDense(Mat B)
2746: {
2747:   Mat_SeqDense   *b;
2749:   PetscMPIInt    size;

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

2755:   PetscNewLog(B,&b);
2756:   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
2757:   B->data = (void*)b;

2759:   b->roworiented = PETSC_TRUE;

2761:   PetscObjectComposeFunction((PetscObject)B,"MatDenseGetLDA_C",MatDenseGetLDA_SeqDense);
2762:   PetscObjectComposeFunction((PetscObject)B,"MatDenseGetArray_C",MatDenseGetArray_SeqDense);
2763:   PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreArray_C",MatDenseRestoreArray_SeqDense);
2764:   PetscObjectComposeFunction((PetscObject)B,"MatDensePlaceArray_C",MatDensePlaceArray_SeqDense);
2765:   PetscObjectComposeFunction((PetscObject)B,"MatDenseResetArray_C",MatDenseResetArray_SeqDense);
2766:   PetscObjectComposeFunction((PetscObject)B,"MatDenseGetArrayRead_C",MatDenseGetArray_SeqDense);
2767:   PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreArrayRead_C",MatDenseRestoreArray_SeqDense);
2768:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqdense_seqaij_C",MatConvert_SeqDense_SeqAIJ);
2769: #if defined(PETSC_HAVE_ELEMENTAL)
2770:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqdense_elemental_C",MatConvert_SeqDense_Elemental);
2771: #endif
2772:   PetscObjectComposeFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",MatSeqDenseSetPreallocation_SeqDense);
2773:   PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqaij_seqdense_C",MatMatMult_SeqAIJ_SeqDense);
2774:   PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C",MatMatMultSymbolic_SeqAIJ_SeqDense);
2775:   PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C",MatMatMultNumeric_SeqAIJ_SeqDense);
2776:   PetscObjectComposeFunction((PetscObject)B,"MatPtAP_seqaij_seqdense_C",MatPtAP_SeqDense_SeqDense);
2777:   PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqaijperm_seqdense_C",MatMatMult_SeqAIJ_SeqDense);
2778:   PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqaijperm_seqdense_C",MatMatMultSymbolic_SeqAIJ_SeqDense);
2779:   PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqaijperm_seqdense_C",MatMatMultNumeric_SeqAIJ_SeqDense);
2780:   PetscObjectComposeFunction((PetscObject)B,"MatPtAP_seqaijperm_seqdense_C",MatPtAP_SeqDense_SeqDense);
2781:   PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqaijsell_seqdense_C",MatMatMult_SeqAIJ_SeqDense);
2782:   PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqaijsell_seqdense_C",MatMatMultSymbolic_SeqAIJ_SeqDense);
2783:   PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqaijsell_seqdense_C",MatMatMultNumeric_SeqAIJ_SeqDense);
2784:   PetscObjectComposeFunction((PetscObject)B,"MatPtAP_seqaijsell_seqdense_C",MatPtAP_SeqDense_SeqDense);
2785:   PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqaijmkl_seqdense_C",MatMatMult_SeqAIJ_SeqDense);
2786:   PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqaijmkl_seqdense_C",MatMatMultSymbolic_SeqAIJ_SeqDense);
2787:   PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqaijmkl_seqdense_C",MatMatMultNumeric_SeqAIJ_SeqDense);
2788:   PetscObjectComposeFunction((PetscObject)B,"MatPtAP_seqaijmkl_seqdense_C",MatPtAP_SeqDense_SeqDense);

2790:   PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMult_seqaij_seqdense_C",MatTransposeMatMult_SeqAIJ_SeqDense);
2791:   PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultSymbolic_seqaij_seqdense_C",MatTransposeMatMultSymbolic_SeqAIJ_SeqDense);
2792:   PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultNumeric_seqaij_seqdense_C",MatTransposeMatMultNumeric_SeqAIJ_SeqDense);
2793:   PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMult_seqaijperm_seqdense_C",MatTransposeMatMult_SeqAIJ_SeqDense);
2794:   PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultSymbolic_seqaijperm_seqdense_C",MatTransposeMatMultSymbolic_SeqAIJ_SeqDense);
2795:   PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultNumeric_seqaijperm_seqdense_C",MatTransposeMatMultNumeric_SeqAIJ_SeqDense);
2796:   PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMult_seqaijsell_seqdense_C",MatTransposeMatMult_SeqAIJ_SeqDense);
2797:   PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultSymbolic_seqaijsell_seqdense_C",MatTransposeMatMultSymbolic_SeqAIJ_SeqDense);
2798:   PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultNumeric_seqaijsell_seqdense_C",MatTransposeMatMultNumeric_SeqAIJ_SeqDense);

2800:   PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMult_seqaijmkl_seqdense_C",MatTransposeMatMult_SeqAIJ_SeqDense);
2801:   PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultSymbolic_seqaijmkl_seqdense_C",MatTransposeMatMultSymbolic_SeqAIJ_SeqDense);
2802:   PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultNumeric_seqaijmkl_seqdense_C",MatTransposeMatMultNumeric_SeqAIJ_SeqDense);
2803:   PetscObjectComposeFunction((PetscObject)B,"MatDenseGetColumn_C",MatDenseGetColumn_SeqDense);
2804:   PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreColumn_C",MatDenseRestoreColumn_SeqDense);
2805:   PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);
2806:   return(0);
2807: }

2809: /*@C
2810:    MatDenseGetColumn - gives access to a column of a dense matrix. This is only the local part of the column. You MUST call MatDenseRestoreColumn() to avoid memory bleeding.

2812:    Not Collective

2814:    Input Parameter:
2815: +  mat - a MATSEQDENSE or MATMPIDENSE matrix
2816: -  col - column index

2818:    Output Parameter:
2819: .  vals - pointer to the data

2821:    Level: intermediate

2823: .seealso: MatDenseRestoreColumn()
2824: @*/
2825: PetscErrorCode MatDenseGetColumn(Mat A,PetscInt col,PetscScalar **vals)
2826: {

2830:   PetscUseMethod(A,"MatDenseGetColumn_C",(Mat,PetscInt,PetscScalar**),(A,col,vals));
2831:   return(0);
2832: }

2834: /*@C
2835:    MatDenseRestoreColumn - returns access to a column of a dense matrix which is returned by MatDenseGetColumn().

2837:    Not Collective

2839:    Input Parameter:
2840: .  mat - a MATSEQDENSE or MATMPIDENSE matrix

2842:    Output Parameter:
2843: .  vals - pointer to the data

2845:    Level: intermediate

2847: .seealso: MatDenseGetColumn()
2848: @*/
2849: PetscErrorCode MatDenseRestoreColumn(Mat A,PetscScalar **vals)
2850: {

2854:   PetscUseMethod(A,"MatDenseRestoreColumn_C",(Mat,PetscScalar**),(A,vals));
2855:   return(0);
2856: }