Actual source code: dense.c

  1: /*$Id: dense.c,v 1.208 2001/09/07 20:09:16 bsmith Exp $*/
  2: /*
  3:      Defines the basic matrix operations for sequential dense.
  4: */

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

 11: int MatAXPY_SeqDense(const PetscScalar *alpha,Mat X,Mat Y,MatStructure str)
 12: {
 13:   Mat_SeqDense *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data;
 14:   int          N = X->m*X->n,m=X->m,ldax=x->lda,lday=y->lda, j,one = 1;

 17:   if (X->m != Y->m || X->n != Y->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
 18:   if (ldax>m || lday>m) {
 19:     for (j=0; j<X->n; j++) {
 20:       BLaxpy_(&m,(PetscScalar*)alpha,x->v+j*ldax,&one,y->v+j*lday,&one);
 21:     }
 22:   } else {
 23:     BLaxpy_(&N,(PetscScalar*)alpha,x->v,&one,y->v,&one);
 24:   }
 25:   PetscLogFlops(2*N-1);
 26:   return(0);
 27: }

 31: int MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
 32: {
 33:   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
 34:   int          i,N = A->m*A->n,count = 0;
 35:   PetscScalar  *v = a->v;

 38:   for (i=0; i<N; i++) {if (*v != 0.0) count++; v++;}

 40:   info->rows_global       = (double)A->m;
 41:   info->columns_global    = (double)A->n;
 42:   info->rows_local        = (double)A->m;
 43:   info->columns_local     = (double)A->n;
 44:   info->block_size        = 1.0;
 45:   info->nz_allocated      = (double)N;
 46:   info->nz_used           = (double)count;
 47:   info->nz_unneeded       = (double)(N-count);
 48:   info->assemblies        = (double)A->num_ass;
 49:   info->mallocs           = 0;
 50:   info->memory            = A->mem;
 51:   info->fill_ratio_given  = 0;
 52:   info->fill_ratio_needed = 0;
 53:   info->factor_mallocs    = 0;

 55:   return(0);
 56: }

 60: int MatScale_SeqDense(const PetscScalar *alpha,Mat A)
 61: {
 62:   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
 63:   int          one = 1,lda = a->lda,j,nz;

 66:   if (lda>A->m) {
 67:     nz = A->m;
 68:     for (j=0; j<A->n; j++) {
 69:       BLscal_(&nz,(PetscScalar*)alpha,a->v+j*lda,&one);
 70:     }
 71:   } else {
 72:     nz = A->m*A->n;
 73:     BLscal_(&nz,(PetscScalar*)alpha,a->v,&one);
 74:   }
 75:   PetscLogFlops(nz);
 76:   return(0);
 77: }
 78: 
 79: /* ---------------------------------------------------------------*/
 80: /* COMMENT: I have chosen to hide row permutation in the pivots,
 81:    rather than put it in the Mat->row slot.*/
 84: int MatLUFactor_SeqDense(Mat A,IS row,IS col,MatFactorInfo *minfo)
 85: {
 86:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
 87:   int          info,ierr;

 90:   if (!mat->pivots) {
 91:     PetscMalloc((A->m+1)*sizeof(int),&mat->pivots);
 92:     PetscLogObjectMemory(A,A->m*sizeof(int));
 93:   }
 94:   A->factor = FACTOR_LU;
 95:   if (!A->m || !A->n) return(0);
 96: #if defined(PETSC_MISSING_LAPACK_GETRF) 
 97:   SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavilable.");
 98: #else
 99:   LAgetrf_(&A->m,&A->n,mat->v,&mat->lda,mat->pivots,&info);
100:   if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization");
101:   if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
102: #endif
103:   PetscLogFlops((2*A->n*A->n*A->n)/3);
104:   return(0);
105: }

109: int MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
110: {
111:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l;
112:   int          lda = mat->lda,j,m,ierr;
113:   Mat          newi;

116:   MatCreateSeqDense(A->comm,A->m,A->n,PETSC_NULL,&newi);
117:   if (cpvalues == MAT_COPY_VALUES) {
118:     l = (Mat_SeqDense*)newi->data;
119:     if (lda>A->m) {
120:       m = A->m;
121:       for (j=0; j<A->n; j++) {
122:         PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));
123:       }
124:     } else {
125:       PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));
126:     }
127:   }
128:   newi->assembled = PETSC_TRUE;
129:   *newmat = newi;
130:   return(0);
131: }

135: int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,MatFactorInfo *info,Mat *fact)
136: {

140:   MatDuplicate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);
141:   return(0);
142: }

146: int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact)
147: {
148:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l = (Mat_SeqDense*)(*fact)->data;
149:   int          lda1=mat->lda,lda2=l->lda, m=A->m,n=A->n, j,ierr;

152:   /* copy the numerical values */
153:   if (lda1>m || lda2>m ) {
154:     for (j=0; j<n; j++) {
155:       PetscMemcpy(l->v+j*lda2,mat->v+j*lda1,m*sizeof(PetscScalar));
156:     }
157:   } else {
158:     PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));
159:   }
160:   (*fact)->factor = 0;
161:   MatLUFactor(*fact,0,0,PETSC_NULL);
162:   return(0);
163: }

167: int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,MatFactorInfo *info,Mat *fact)
168: {

172:   MatConvert(A,MATSAME,fact);
173:   return(0);
174: }

178: int MatCholeskyFactor_SeqDense(Mat A,IS perm,MatFactorInfo *factinfo)
179: {
180:   Mat_SeqDense  *mat = (Mat_SeqDense*)A->data;
181:   int           info,ierr;
182: 
184:   if (mat->pivots) {
185:     PetscFree(mat->pivots);
186:     PetscLogObjectMemory(A,-A->m*sizeof(int));
187:     mat->pivots = 0;
188:   }
189:   if (!A->m || !A->n) return(0);
190: #if defined(PETSC_MISSING_LAPACK_POTRF) 
191:   SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavilable.");
192: #else
193:   LApotrf_("L",&A->n,mat->v,&mat->lda,&info);
194:   if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %d",info-1);
195: #endif
196:   A->factor = FACTOR_CHOLESKY;
197:   PetscLogFlops((A->n*A->n*A->n)/3);
198:   return(0);
199: }

203: int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact)
204: {
206:   MatFactorInfo info;

209:   info.fill = 1.0;
210:   MatCholeskyFactor_SeqDense(*fact,0,&info);
211:   return(0);
212: }

216: int MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
217: {
218:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
219:   int          one = 1,info,ierr;
220:   PetscScalar  *x,*y;
221: 
223:   if (!A->m || !A->n) return(0);
224:   VecGetArray(xx,&x);
225:   VecGetArray(yy,&y);
226:   PetscMemcpy(y,x,A->m*sizeof(PetscScalar));
227:   if (A->factor == FACTOR_LU) {
228: #if defined(PETSC_MISSING_LAPACK_GETRS) 
229:     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable.");
230: #else
231:     LAgetrs_("N",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info);
232:     if (info) SETERRQ(PETSC_ERR_LIB,"MBad solve");
233: #endif
234:   } else if (A->factor == FACTOR_CHOLESKY){
235: #if defined(PETSC_MISSING_LAPACK_POTRS) 
236:     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable.");
237: #else
238:     LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info);
239:     if (info) SETERRQ(PETSC_ERR_LIB,"MBad solve");
240: #endif
241:   }
242:   else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
243:   VecRestoreArray(xx,&x);
244:   VecRestoreArray(yy,&y);
245:   PetscLogFlops(2*A->n*A->n - A->n);
246:   return(0);
247: }

251: int MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
252: {
253:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
254:   int          ierr,one = 1,info;
255:   PetscScalar  *x,*y;
256: 
258:   if (!A->m || !A->n) return(0);
259:   VecGetArray(xx,&x);
260:   VecGetArray(yy,&y);
261:   PetscMemcpy(y,x,A->m*sizeof(PetscScalar));
262:   /* assume if pivots exist then use LU; else Cholesky */
263:   if (mat->pivots) {
264: #if defined(PETSC_MISSING_LAPACK_GETRS) 
265:     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable.");
266: #else
267:     LAgetrs_("T",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info);
268:     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
269: #endif
270:   } else {
271: #if defined(PETSC_MISSING_LAPACK_POTRS) 
272:     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable.");
273: #else
274:     LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info);
275:     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
276: #endif
277:   }
278:   VecRestoreArray(xx,&x);
279:   VecRestoreArray(yy,&y);
280:   PetscLogFlops(2*A->n*A->n - A->n);
281:   return(0);
282: }

286: int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
287: {
288:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
289:   int          ierr,one = 1,info;
290:   PetscScalar  *x,*y,sone = 1.0;
291:   Vec          tmp = 0;
292: 
294:   VecGetArray(xx,&x);
295:   VecGetArray(yy,&y);
296:   if (!A->m || !A->n) return(0);
297:   if (yy == zz) {
298:     VecDuplicate(yy,&tmp);
299:     PetscLogObjectParent(A,tmp);
300:     VecCopy(yy,tmp);
301:   }
302:   PetscMemcpy(y,x,A->m*sizeof(PetscScalar));
303:   /* assume if pivots exist then use LU; else Cholesky */
304:   if (mat->pivots) {
305: #if defined(PETSC_MISSING_LAPACK_GETRS) 
306:     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable.");
307: #else
308:     LAgetrs_("N",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info);
309:     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
310: #endif
311:   } else {
312: #if defined(PETSC_MISSING_LAPACK_POTRS) 
313:     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable.");
314: #else
315:     LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info);
316:     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
317: #endif
318:   }
319:   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
320:   else     {VecAXPY(&sone,zz,yy);}
321:   VecRestoreArray(xx,&x);
322:   VecRestoreArray(yy,&y);
323:   PetscLogFlops(2*A->n*A->n);
324:   return(0);
325: }

329: int MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
330: {
331:   Mat_SeqDense  *mat = (Mat_SeqDense*)A->data;
332:   int           one = 1,info,ierr;
333:   PetscScalar   *x,*y,sone = 1.0;
334:   Vec           tmp;
335: 
337:   if (!A->m || !A->n) return(0);
338:   VecGetArray(xx,&x);
339:   VecGetArray(yy,&y);
340:   if (yy == zz) {
341:     VecDuplicate(yy,&tmp);
342:     PetscLogObjectParent(A,tmp);
343:     VecCopy(yy,tmp);
344:   }
345:   PetscMemcpy(y,x,A->m*sizeof(PetscScalar));
346:   /* assume if pivots exist then use LU; else Cholesky */
347:   if (mat->pivots) {
348: #if defined(PETSC_MISSING_LAPACK_GETRS) 
349:     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable.");
350: #else
351:     LAgetrs_("T",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info);
352:     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
353: #endif
354:   } else {
355: #if defined(PETSC_MISSING_LAPACK_POTRS) 
356:     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable.");
357: #else
358:     LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info);
359:     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
360: #endif
361:   }
362:   if (tmp) {
363:     VecAXPY(&sone,tmp,yy);
364:     VecDestroy(tmp);
365:   } else {
366:     VecAXPY(&sone,zz,yy);
367:   }
368:   VecRestoreArray(xx,&x);
369:   VecRestoreArray(yy,&y);
370:   PetscLogFlops(2*A->n*A->n);
371:   return(0);
372: }
373: /* ------------------------------------------------------------------*/
376: int MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,
377:                           PetscReal shift,int its,int lits,Vec xx)
378: {
379:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
380:   PetscScalar  *x,*b,*v = mat->v,zero = 0.0,xt;
381:   int          ierr,m = A->m,i;
382: #if !defined(PETSC_USE_COMPLEX)
383:   int          o = 1;
384: #endif

387:   if (flag & SOR_ZERO_INITIAL_GUESS) {
388:     /* this is a hack fix, should have another version without the second BLdot */
389:     VecSet(&zero,xx);
390:   }
391:   VecGetArray(xx,&x);
392:   VecGetArray(bb,&b);
393:   its  = its*lits;
394:   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits);
395:   while (its--) {
396:     if (flag & SOR_FORWARD_SWEEP){
397:       for (i=0; i<m; i++) {
398: #if defined(PETSC_USE_COMPLEX)
399:         /* cannot use BLAS dot for complex because compiler/linker is 
400:            not happy about returning a double complex */
401:         int         _i;
402:         PetscScalar sum = b[i];
403:         for (_i=0; _i<m; _i++) {
404:           sum -= PetscConj(v[i+_i*m])*x[_i];
405:         }
406:         xt = sum;
407: #else
408:         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
409: #endif
410:         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
411:       }
412:     }
413:     if (flag & SOR_BACKWARD_SWEEP) {
414:       for (i=m-1; i>=0; i--) {
415: #if defined(PETSC_USE_COMPLEX)
416:         /* cannot use BLAS dot for complex because compiler/linker is 
417:            not happy about returning a double complex */
418:         int         _i;
419:         PetscScalar sum = b[i];
420:         for (_i=0; _i<m; _i++) {
421:           sum -= PetscConj(v[i+_i*m])*x[_i];
422:         }
423:         xt = sum;
424: #else
425:         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
426: #endif
427:         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
428:       }
429:     }
430:   }
431:   VecRestoreArray(bb,&b);
432:   VecRestoreArray(xx,&x);
433:   return(0);
434: }

436: /* -----------------------------------------------------------------*/
439: int MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
440: {
441:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
442:   PetscScalar  *v = mat->v,*x,*y;
443:   int          ierr,_One=1;
444:   PetscScalar  _DOne=1.0,_DZero=0.0;

447:   if (!A->m || !A->n) return(0);
448:   VecGetArray(xx,&x);
449:   VecGetArray(yy,&y);
450:   LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One);
451:   VecRestoreArray(xx,&x);
452:   VecRestoreArray(yy,&y);
453:   PetscLogFlops(2*A->m*A->n - A->n);
454:   return(0);
455: }

459: int MatMult_SeqDense(Mat A,Vec xx,Vec yy)
460: {
461:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
462:   PetscScalar  *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0;
463:   int          ierr,_One=1;

466:   if (!A->m || !A->n) return(0);
467:   VecGetArray(xx,&x);
468:   VecGetArray(yy,&y);
469:   LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One);
470:   VecRestoreArray(xx,&x);
471:   VecRestoreArray(yy,&y);
472:   PetscLogFlops(2*A->m*A->n - A->m);
473:   return(0);
474: }

478: int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
479: {
480:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
481:   PetscScalar  *v = mat->v,*x,*y,_DOne=1.0;
482:   int          ierr,_One=1;

485:   if (!A->m || !A->n) return(0);
486:   if (zz != yy) {VecCopy(zz,yy);}
487:   VecGetArray(xx,&x);
488:   VecGetArray(yy,&y);
489:   LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
490:   VecRestoreArray(xx,&x);
491:   VecRestoreArray(yy,&y);
492:   PetscLogFlops(2*A->m*A->n);
493:   return(0);
494: }

498: int MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
499: {
500:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
501:   PetscScalar  *v = mat->v,*x,*y;
502:   int          ierr,_One=1;
503:   PetscScalar  _DOne=1.0;

506:   if (!A->m || !A->n) return(0);
507:   if (zz != yy) {VecCopy(zz,yy);}
508:   VecGetArray(xx,&x);
509:   VecGetArray(yy,&y);
510:   LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
511:   VecRestoreArray(xx,&x);
512:   VecRestoreArray(yy,&y);
513:   PetscLogFlops(2*A->m*A->n);
514:   return(0);
515: }

517: /* -----------------------------------------------------------------*/
520: int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,PetscScalar **vals)
521: {
522:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
523:   PetscScalar  *v;
524:   int          i,ierr;
525: 
527:   *ncols = A->n;
528:   if (cols) {
529:     PetscMalloc((A->n+1)*sizeof(int),cols);
530:     for (i=0; i<A->n; i++) (*cols)[i] = i;
531:   }
532:   if (vals) {
533:     PetscMalloc((A->n+1)*sizeof(PetscScalar),vals);
534:     v    = mat->v + row;
535:     for (i=0; i<A->n; i++) {(*vals)[i] = *v; v += mat->lda;}
536:   }
537:   return(0);
538: }

542: int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,PetscScalar **vals)
543: {
546:   if (cols) {PetscFree(*cols);}
547:   if (vals) {PetscFree(*vals); }
548:   return(0);
549: }
550: /* ----------------------------------------------------------------*/
553: int MatSetValues_SeqDense(Mat A,int m,const int indexm[],int n,const int indexn[],const PetscScalar v[],InsertMode addv)
554: {
555:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
556:   int          i,j,idx=0;
557: 
559:   if (!mat->roworiented) {
560:     if (addv == INSERT_VALUES) {
561:       for (j=0; j<n; j++) {
562:         if (indexn[j] < 0) {idx += m; continue;}
563: #if defined(PETSC_USE_BOPT_g)  
564:         if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",indexn[j],A->n-1);
565: #endif
566:         for (i=0; i<m; i++) {
567:           if (indexm[i] < 0) {idx++; continue;}
568: #if defined(PETSC_USE_BOPT_g)  
569:           if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",indexm[i],A->m-1);
570: #endif
571:           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
572:         }
573:       }
574:     } else {
575:       for (j=0; j<n; j++) {
576:         if (indexn[j] < 0) {idx += m; continue;}
577: #if defined(PETSC_USE_BOPT_g)  
578:         if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",indexn[j],A->n-1);
579: #endif
580:         for (i=0; i<m; i++) {
581:           if (indexm[i] < 0) {idx++; continue;}
582: #if defined(PETSC_USE_BOPT_g)  
583:           if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",indexm[i],A->m-1);
584: #endif
585:           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
586:         }
587:       }
588:     }
589:   } else {
590:     if (addv == INSERT_VALUES) {
591:       for (i=0; i<m; i++) {
592:         if (indexm[i] < 0) { idx += n; continue;}
593: #if defined(PETSC_USE_BOPT_g)  
594:         if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",indexm[j],A->m-1);
595: #endif
596:         for (j=0; j<n; j++) {
597:           if (indexn[j] < 0) { idx++; continue;}
598: #if defined(PETSC_USE_BOPT_g)  
599:           if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",indexn[i],A->n-1);
600: #endif
601:           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
602:         }
603:       }
604:     } else {
605:       for (i=0; i<m; i++) {
606:         if (indexm[i] < 0) { idx += n; continue;}
607: #if defined(PETSC_USE_BOPT_g)  
608:         if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",indexm[j],A->m-1);
609: #endif
610:         for (j=0; j<n; j++) {
611:           if (indexn[j] < 0) { idx++; continue;}
612: #if defined(PETSC_USE_BOPT_g)  
613:           if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",indexn[i],A->n-1);
614: #endif
615:           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
616:         }
617:       }
618:     }
619:   }
620:   return(0);
621: }

625: int MatGetValues_SeqDense(Mat A,int m,const int indexm[],int n,const int indexn[],PetscScalar v[])
626: {
627:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
628:   int          i,j;
629:   PetscScalar  *vpt = v;

632:   /* row-oriented output */
633:   for (i=0; i<m; i++) {
634:     for (j=0; j<n; j++) {
635:       *vpt++ = mat->v[indexn[j]*mat->lda + indexm[i]];
636:     }
637:   }
638:   return(0);
639: }

641: /* -----------------------------------------------------------------*/

643:  #include petscsys.h

647: int MatLoad_SeqDense(PetscViewer viewer,MatType type,Mat *A)
648: {
649:   Mat_SeqDense *a;
650:   Mat          B;
651:   int          *scols,i,j,nz,ierr,fd,header[4],size;
652:   int          *rowlengths = 0,M,N,*cols;
653:   PetscScalar  *vals,*svals,*v,*w;
654:   MPI_Comm     comm = ((PetscObject)viewer)->comm;

657:   MPI_Comm_size(comm,&size);
658:   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
659:   PetscViewerBinaryGetDescriptor(viewer,&fd);
660:   PetscBinaryRead(fd,header,4,PETSC_INT);
661:   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
662:   M = header[1]; N = header[2]; nz = header[3];

664:   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
665:     MatCreateSeqDense(comm,M,N,PETSC_NULL,A);
666:     B    = *A;
667:     a    = (Mat_SeqDense*)B->data;
668:     v    = a->v;
669:     /* Allocate some temp space to read in the values and then flip them
670:        from row major to column major */
671:     PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);
672:     /* read in nonzero values */
673:     PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);
674:     /* now flip the values and store them in the matrix*/
675:     for (j=0; j<N; j++) {
676:       for (i=0; i<M; i++) {
677:         *v++ =w[i*N+j];
678:       }
679:     }
680:     PetscFree(w);
681:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
682:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
683:   } else {
684:     /* read row lengths */
685:     PetscMalloc((M+1)*sizeof(int),&rowlengths);
686:     PetscBinaryRead(fd,rowlengths,M,PETSC_INT);

688:     /* create our matrix */
689:     MatCreateSeqDense(comm,M,N,PETSC_NULL,A);
690:     B = *A;
691:     a = (Mat_SeqDense*)B->data;
692:     v = a->v;

694:     /* read column indices and nonzeros */
695:     PetscMalloc((nz+1)*sizeof(int),&scols);
696:     cols = scols;
697:     PetscBinaryRead(fd,cols,nz,PETSC_INT);
698:     PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);
699:     vals = svals;
700:     PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);

702:     /* insert into matrix */
703:     for (i=0; i<M; i++) {
704:       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
705:       svals += rowlengths[i]; scols += rowlengths[i];
706:     }
707:     PetscFree(vals);
708:     PetscFree(cols);
709:     PetscFree(rowlengths);

711:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
712:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
713:   }
714:   return(0);
715: }

717:  #include petscsys.h

721: static int MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
722: {
723:   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
724:   int               ierr,i,j;
725:   char              *name;
726:   PetscScalar       *v;
727:   PetscViewerFormat format;

730:   PetscObjectGetName((PetscObject)A,&name);
731:   PetscViewerGetFormat(viewer,&format);
732:   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
733:     return(0);  /* do nothing for now */
734:   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
735:     PetscViewerASCIIUseTabs(viewer,PETSC_NO);
736:     for (i=0; i<A->m; i++) {
737:       v = a->v + i;
738:       PetscViewerASCIIPrintf(viewer,"row %d:",i);
739:       for (j=0; j<A->n; j++) {
740: #if defined(PETSC_USE_COMPLEX)
741:         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
742:           PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));
743:         } else if (PetscRealPart(*v)) {
744:           PetscViewerASCIIPrintf(viewer," (%d, %g) ",j,PetscRealPart(*v));
745:         }
746: #else
747:         if (*v) {
748:           PetscViewerASCIIPrintf(viewer," (%d, %g) ",j,*v);
749:         }
750: #endif
751:         v += a->lda;
752:       }
753:       PetscViewerASCIIPrintf(viewer,"\n");
754:     }
755:     PetscViewerASCIIUseTabs(viewer,PETSC_YES);
756:   } else {
757:     PetscViewerASCIIUseTabs(viewer,PETSC_NO);
758: #if defined(PETSC_USE_COMPLEX)
759:     PetscTruth allreal = PETSC_TRUE;
760:     /* determine if matrix has all real values */
761:     v = a->v;
762:     for (i=0; i<A->m*A->n; i++) {
763:         if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;}
764:     }
765: #endif
766:     if (format == PETSC_VIEWER_ASCII_MATLAB) {
767:       PetscObjectGetName((PetscObject)A,&name);
768:       PetscViewerASCIIPrintf(viewer,"%% Size = %d %d \n",A->m,A->n);
769:       PetscViewerASCIIPrintf(viewer,"%s = zeros(%d,%d);\n",name,A->m,A->n);
770:       PetscViewerASCIIPrintf(viewer,"%s = [\n",name);
771:     }

773:     for (i=0; i<A->m; i++) {
774:       v = a->v + i;
775:       for (j=0; j<A->n; j++) {
776: #if defined(PETSC_USE_COMPLEX)
777:         if (allreal) {
778:           PetscViewerASCIIPrintf(viewer,"%6.4e ",PetscRealPart(*v));
779:         } else {
780:           PetscViewerASCIIPrintf(viewer,"%6.4e + %6.4e i ",PetscRealPart(*v),PetscImaginaryPart(*v));
781:         }
782: #else
783:         PetscViewerASCIIPrintf(viewer,"%6.4e ",*v);
784: #endif
785:         v += a->lda;
786:       }
787:       PetscViewerASCIIPrintf(viewer,"\n");
788:     }
789:     if (format == PETSC_VIEWER_ASCII_MATLAB) {
790:       PetscViewerASCIIPrintf(viewer,"];\n");
791:     }
792:     PetscViewerASCIIUseTabs(viewer,PETSC_YES);
793:   }
794:   PetscViewerFlush(viewer);
795:   return(0);
796: }

800: static int MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
801: {
802:   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
803:   int               ict,j,n = A->n,m = A->m,i,fd,*col_lens,ierr,nz = m*n;
804:   PetscScalar       *v,*anonz,*vals;
805:   PetscViewerFormat format;
806: 
808:   PetscViewerBinaryGetDescriptor(viewer,&fd);

810:   PetscViewerGetFormat(viewer,&format);
811:   if (format == PETSC_VIEWER_BINARY_NATIVE) {
812:     /* store the matrix as a dense matrix */
813:     PetscMalloc(4*sizeof(int),&col_lens);
814:     col_lens[0] = MAT_FILE_COOKIE;
815:     col_lens[1] = m;
816:     col_lens[2] = n;
817:     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
818:     PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1);
819:     PetscFree(col_lens);

821:     /* write out matrix, by rows */
822:     PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);
823:     v    = a->v;
824:     for (i=0; i<m; i++) {
825:       for (j=0; j<n; j++) {
826:         vals[i + j*m] = *v++;
827:       }
828:     }
829:     PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0);
830:     PetscFree(vals);
831:   } else {
832:     PetscMalloc((4+nz)*sizeof(int),&col_lens);
833:     col_lens[0] = MAT_FILE_COOKIE;
834:     col_lens[1] = m;
835:     col_lens[2] = n;
836:     col_lens[3] = nz;

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

842:     /* Possibly should write in smaller increments, not whole matrix at once? */
843:     /* store column indices (zero start index) */
844:     ict = 0;
845:     for (i=0; i<m; i++) {
846:       for (j=0; j<n; j++) col_lens[ict++] = j;
847:     }
848:     PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0);
849:     PetscFree(col_lens);

851:     /* store nonzero values */
852:     PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);
853:     ict  = 0;
854:     for (i=0; i<m; i++) {
855:       v = a->v + i;
856:       for (j=0; j<n; j++) {
857:         anonz[ict++] = *v; v += a->lda;
858:       }
859:     }
860:     PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0);
861:     PetscFree(anonz);
862:   }
863:   return(0);
864: }

868: int MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
869: {
870:   Mat               A = (Mat) Aa;
871:   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
872:   int               m = A->m,n = A->n,color,i,j,ierr;
873:   PetscScalar       *v = a->v;
874:   PetscViewer       viewer;
875:   PetscDraw         popup;
876:   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0;
877:   PetscViewerFormat format;


881:   PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);
882:   PetscViewerGetFormat(viewer,&format);
883:   PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);

885:   /* Loop over matrix elements drawing boxes */
886:   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
887:     /* Blue for negative and Red for positive */
888:     color = PETSC_DRAW_BLUE;
889:     for(j = 0; j < n; j++) {
890:       x_l = j;
891:       x_r = x_l + 1.0;
892:       for(i = 0; i < m; i++) {
893:         y_l = m - i - 1.0;
894:         y_r = y_l + 1.0;
895: #if defined(PETSC_USE_COMPLEX)
896:         if (PetscRealPart(v[j*m+i]) >  0.) {
897:           color = PETSC_DRAW_RED;
898:         } else if (PetscRealPart(v[j*m+i]) <  0.) {
899:           color = PETSC_DRAW_BLUE;
900:         } else {
901:           continue;
902:         }
903: #else
904:         if (v[j*m+i] >  0.) {
905:           color = PETSC_DRAW_RED;
906:         } else if (v[j*m+i] <  0.) {
907:           color = PETSC_DRAW_BLUE;
908:         } else {
909:           continue;
910:         }
911: #endif
912:         PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
913:       }
914:     }
915:   } else {
916:     /* use contour shading to indicate magnitude of values */
917:     /* first determine max of all nonzero values */
918:     for(i = 0; i < m*n; i++) {
919:       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
920:     }
921:     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
922:     PetscDrawGetPopup(draw,&popup);
923:     if (popup) {PetscDrawScalePopup(popup,0.0,maxv);}
924:     for(j = 0; j < n; j++) {
925:       x_l = j;
926:       x_r = x_l + 1.0;
927:       for(i = 0; i < m; i++) {
928:         y_l   = m - i - 1.0;
929:         y_r   = y_l + 1.0;
930:         color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i]));
931:         PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
932:       }
933:     }
934:   }
935:   return(0);
936: }

940: int MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
941: {
942:   PetscDraw  draw;
943:   PetscTruth isnull;
944:   PetscReal  xr,yr,xl,yl,h,w;
945:   int        ierr;

948:   PetscViewerDrawGetDraw(viewer,0,&draw);
949:   PetscDrawIsNull(draw,&isnull);
950:   if (isnull == PETSC_TRUE) return(0);

952:   PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
953:   xr  = A->n; yr = A->m; h = yr/10.0; w = xr/10.0;
954:   xr += w;    yr += h;  xl = -w;     yl = -h;
955:   PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
956:   PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);
957:   PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);
958:   return(0);
959: }

963: int MatView_SeqDense(Mat A,PetscViewer viewer)
964: {
965:   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
966:   int          ierr;
967:   PetscTruth   issocket,isascii,isbinary,isdraw;

970:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);
971:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
972:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
973:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);

975:   if (issocket) {
976:     if (a->lda>A->m) SETERRQ(1,"Case can not handle LDA");
977:     PetscViewerSocketPutScalar(viewer,A->m,A->n,a->v);
978:   } else if (isascii) {
979:     MatView_SeqDense_ASCII(A,viewer);
980:   } else if (isbinary) {
981:     MatView_SeqDense_Binary(A,viewer);
982:   } else if (isdraw) {
983:     MatView_SeqDense_Draw(A,viewer);
984:   } else {
985:     SETERRQ1(1,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name);
986:   }
987:   return(0);
988: }

992: int MatDestroy_SeqDense(Mat mat)
993: {
994:   Mat_SeqDense *l = (Mat_SeqDense*)mat->data;
995:   int          ierr;

998: #if defined(PETSC_USE_LOG)
999:   PetscLogObjectState((PetscObject)mat,"Rows %d Cols %d",mat->m,mat->n);
1000: #endif
1001:   if (l->pivots) {PetscFree(l->pivots);}
1002:   if (!l->user_alloc) {PetscFree(l->v);}
1003:   PetscFree(l);
1004:   return(0);
1005: }

1009: int MatTranspose_SeqDense(Mat A,Mat *matout)
1010: {
1011:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1012:   int          k,j,m,n,M,ierr;
1013:   PetscScalar  *v,tmp;

1016:   v = mat->v; m = A->m; M = mat->lda; n = A->n;
1017:   if (!matout) { /* in place transpose */
1018:     if (m != n) {
1019:       SETERRQ(1,"Can not transpose non-square matrix in place");
1020:     } else {
1021:       for (j=0; j<m; j++) {
1022:         for (k=0; k<j; k++) {
1023:           tmp = v[j + k*M];
1024:           v[j + k*M] = v[k + j*M];
1025:           v[k + j*M] = tmp;
1026:         }
1027:       }
1028:     }
1029:   } else { /* out-of-place transpose */
1030:     Mat          tmat;
1031:     Mat_SeqDense *tmatd;
1032:     PetscScalar  *v2;

1034:     MatCreateSeqDense(A->comm,A->n,A->m,PETSC_NULL,&tmat);
1035:     tmatd = (Mat_SeqDense*)tmat->data;
1036:     v = mat->v; v2 = tmatd->v;
1037:     for (j=0; j<n; j++) {
1038:       for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M];
1039:     }
1040:     MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);
1041:     MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);
1042:     *matout = tmat;
1043:   }
1044:   return(0);
1045: }

1049: int MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg)
1050: {
1051:   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1052:   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
1053:   int          i,j;
1054:   PetscScalar  *v1 = mat1->v,*v2 = mat2->v;

1057:   if (A1->m != A2->m) {*flg = PETSC_FALSE; return(0);}
1058:   if (A1->n != A2->n) {*flg = PETSC_FALSE; return(0);}
1059:   for (i=0; i<A1->m; i++) {
1060:     v1 = mat1->v+i; v2 = mat2->v+i;
1061:     for (j=0; j<A1->n; j++) {
1062:       if (*v1 != *v2) {*flg = PETSC_FALSE; return(0);}
1063:       v1 += mat1->lda; v2 += mat2->lda;
1064:     }
1065:   }
1066:   *flg = PETSC_TRUE;
1067:   return(0);
1068: }

1072: int MatGetDiagonal_SeqDense(Mat A,Vec v)
1073: {
1074:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1075:   int          ierr,i,n,len;
1076:   PetscScalar  *x,zero = 0.0;

1079:   VecSet(&zero,v);
1080:   VecGetSize(v,&n);
1081:   VecGetArray(v,&x);
1082:   len = PetscMin(A->m,A->n);
1083:   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
1084:   for (i=0; i<len; i++) {
1085:     x[i] = mat->v[i*mat->lda + i];
1086:   }
1087:   VecRestoreArray(v,&x);
1088:   return(0);
1089: }

1093: int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1094: {
1095:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1096:   PetscScalar  *l,*r,x,*v;
1097:   int          ierr,i,j,m = A->m,n = A->n;

1100:   if (ll) {
1101:     VecGetSize(ll,&m);
1102:     VecGetArray(ll,&l);
1103:     if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1104:     for (i=0; i<m; i++) {
1105:       x = l[i];
1106:       v = mat->v + i;
1107:       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1108:     }
1109:     VecRestoreArray(ll,&l);
1110:     PetscLogFlops(n*m);
1111:   }
1112:   if (rr) {
1113:     VecGetSize(rr,&n);
1114:     VecGetArray(rr,&r);
1115:     if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1116:     for (i=0; i<n; i++) {
1117:       x = r[i];
1118:       v = mat->v + i*m;
1119:       for (j=0; j<m; j++) { (*v++) *= x;}
1120:     }
1121:     VecRestoreArray(rr,&r);
1122:     PetscLogFlops(n*m);
1123:   }
1124:   return(0);
1125: }

1129: int MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1130: {
1131:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1132:   PetscScalar  *v = mat->v;
1133:   PetscReal    sum = 0.0;
1134:   int          lda=mat->lda,m=A->m,i,j;

1137:   if (type == NORM_FROBENIUS) {
1138:     if (lda>m) {
1139:       for (j=0; j<A->n; j++) {
1140:         v = mat->v+j*lda;
1141:         for (i=0; i<m; i++) {
1142: #if defined(PETSC_USE_COMPLEX)
1143:           sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1144: #else
1145:           sum += (*v)*(*v); v++;
1146: #endif
1147:         }
1148:       }
1149:     } else {
1150:       for (i=0; i<A->n*A->m; i++) {
1151: #if defined(PETSC_USE_COMPLEX)
1152:         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1153: #else
1154:         sum += (*v)*(*v); v++;
1155: #endif
1156:       }
1157:     }
1158:     *nrm = sqrt(sum);
1159:     PetscLogFlops(2*A->n*A->m);
1160:   } else if (type == NORM_1) {
1161:     *nrm = 0.0;
1162:     for (j=0; j<A->n; j++) {
1163:       v = mat->v + j*mat->lda;
1164:       sum = 0.0;
1165:       for (i=0; i<A->m; i++) {
1166:         sum += PetscAbsScalar(*v);  v++;
1167:       }
1168:       if (sum > *nrm) *nrm = sum;
1169:     }
1170:     PetscLogFlops(A->n*A->m);
1171:   } else if (type == NORM_INFINITY) {
1172:     *nrm = 0.0;
1173:     for (j=0; j<A->m; j++) {
1174:       v = mat->v + j;
1175:       sum = 0.0;
1176:       for (i=0; i<A->n; i++) {
1177:         sum += PetscAbsScalar(*v); v += mat->lda;
1178:       }
1179:       if (sum > *nrm) *nrm = sum;
1180:     }
1181:     PetscLogFlops(A->n*A->m);
1182:   } else {
1183:     SETERRQ(PETSC_ERR_SUP,"No two norm");
1184:   }
1185:   return(0);
1186: }

1190: int MatSetOption_SeqDense(Mat A,MatOption op)
1191: {
1192:   Mat_SeqDense *aij = (Mat_SeqDense*)A->data;
1193: 
1195:   switch (op) {
1196:   case MAT_ROW_ORIENTED:
1197:     aij->roworiented = PETSC_TRUE;
1198:     break;
1199:   case MAT_COLUMN_ORIENTED:
1200:     aij->roworiented = PETSC_FALSE;
1201:     break;
1202:   case MAT_ROWS_SORTED:
1203:   case MAT_ROWS_UNSORTED:
1204:   case MAT_COLUMNS_SORTED:
1205:   case MAT_COLUMNS_UNSORTED:
1206:   case MAT_NO_NEW_NONZERO_LOCATIONS:
1207:   case MAT_YES_NEW_NONZERO_LOCATIONS:
1208:   case MAT_NEW_NONZERO_LOCATION_ERR:
1209:   case MAT_NO_NEW_DIAGONALS:
1210:   case MAT_YES_NEW_DIAGONALS:
1211:   case MAT_IGNORE_OFF_PROC_ENTRIES:
1212:   case MAT_USE_HASH_TABLE:
1213:     PetscLogInfo(A,"MatSetOption_SeqDense:Option ignored\n");
1214:     break;
1215:   default:
1216:     SETERRQ(PETSC_ERR_SUP,"unknown option");
1217:   }
1218:   return(0);
1219: }

1223: int MatZeroEntries_SeqDense(Mat A)
1224: {
1225:   Mat_SeqDense *l = (Mat_SeqDense*)A->data;
1226:   int          lda=l->lda,m=A->m,j, ierr;

1229:   if (lda>m) {
1230:     for (j=0; j<A->n; j++) {
1231:       PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));
1232:     }
1233:   } else {
1234:     PetscMemzero(l->v,A->m*A->n*sizeof(PetscScalar));
1235:   }
1236:   return(0);
1237: }

1241: int MatGetBlockSize_SeqDense(Mat A,int *bs)
1242: {
1244:   *bs = 1;
1245:   return(0);
1246: }

1250: int MatZeroRows_SeqDense(Mat A,IS is,const PetscScalar *diag)
1251: {
1252:   Mat_SeqDense *l = (Mat_SeqDense*)A->data;
1253:   int          n = A->n,i,j,ierr,N,*rows;
1254:   PetscScalar  *slot;

1257:   ISGetLocalSize(is,&N);
1258:   ISGetIndices(is,&rows);
1259:   for (i=0; i<N; i++) {
1260:     slot = l->v + rows[i];
1261:     for (j=0; j<n; j++) { *slot = 0.0; slot += n;}
1262:   }
1263:   if (diag) {
1264:     for (i=0; i<N; i++) {
1265:       slot = l->v + (n+1)*rows[i];
1266:       *slot = *diag;
1267:     }
1268:   }
1269:   ISRestoreIndices(is,&rows);
1270:   return(0);
1271: }

1275: int MatGetArray_SeqDense(Mat A,PetscScalar *array[])
1276: {
1277:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;

1280:   *array = mat->v;
1281:   return(0);
1282: }

1286: int MatRestoreArray_SeqDense(Mat A,PetscScalar *array[])
1287: {
1289:   *array = 0; /* user cannot accidently use the array later */
1290:   return(0);
1291: }

1295: static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
1296: {
1297:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1298:   int          i,j,ierr,m = A->m,*irow,*icol,nrows,ncols;
1299:   PetscScalar  *av,*bv,*v = mat->v;
1300:   Mat          newmat;

1303:   ISGetIndices(isrow,&irow);
1304:   ISGetIndices(iscol,&icol);
1305:   ISGetLocalSize(isrow,&nrows);
1306:   ISGetLocalSize(iscol,&ncols);
1307: 
1308:   /* Check submatrixcall */
1309:   if (scall == MAT_REUSE_MATRIX) {
1310:     int n_cols,n_rows;
1311:     MatGetSize(*B,&n_rows,&n_cols);
1312:     if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
1313:     newmat = *B;
1314:   } else {
1315:     /* Create and fill new matrix */
1316:     MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat);
1317:   }

1319:   /* Now extract the data pointers and do the copy,column at a time */
1320:   bv = ((Mat_SeqDense*)newmat->data)->v;
1321: 
1322:   for (i=0; i<ncols; i++) {
1323:     av = v + m*icol[i];
1324:     for (j=0; j<nrows; j++) {
1325:       *bv++ = av[irow[j]];
1326:     }
1327:   }

1329:   /* Assemble the matrices so that the correct flags are set */
1330:   MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
1331:   MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);

1333:   /* Free work space */
1334:   ISRestoreIndices(isrow,&irow);
1335:   ISRestoreIndices(iscol,&icol);
1336:   *B = newmat;
1337:   return(0);
1338: }

1342: int MatGetSubMatrices_SeqDense(Mat A,int n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1343: {
1344:   int ierr,i;

1347:   if (scall == MAT_INITIAL_MATRIX) {
1348:     PetscMalloc((n+1)*sizeof(Mat),B);
1349:   }

1351:   for (i=0; i<n; i++) {
1352:     MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);
1353:   }
1354:   return(0);
1355: }

1359: int MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
1360: {
1361:   Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
1362:   int          lda1=a->lda,lda2=b->lda, m=A->m,n=A->n, j,ierr;

1365:   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
1366:   if (A->ops->copy != B->ops->copy) {
1367:     MatCopy_Basic(A,B,str);
1368:     return(0);
1369:   }
1370:   if (m != B->m || n != B->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1371:   if (lda1>m || lda2>m) {
1372:     for (j=0; j<n; j++) {
1373:       PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));
1374:     }
1375:   } else {
1376:     PetscMemcpy(b->v,a->v,A->m*A->n*sizeof(PetscScalar));
1377:   }
1378:   return(0);
1379: }

1383: int MatSetUpPreallocation_SeqDense(Mat A)
1384: {
1385:   int        ierr;

1388:    MatSeqDenseSetPreallocation(A,0);
1389:   return(0);
1390: }

1392: /* -------------------------------------------------------------------*/
1393: static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1394:        MatGetRow_SeqDense,
1395:        MatRestoreRow_SeqDense,
1396:        MatMult_SeqDense,
1397: /* 4*/ MatMultAdd_SeqDense,
1398:        MatMultTranspose_SeqDense,
1399:        MatMultTransposeAdd_SeqDense,
1400:        MatSolve_SeqDense,
1401:        MatSolveAdd_SeqDense,
1402:        MatSolveTranspose_SeqDense,
1403: /*10*/ MatSolveTransposeAdd_SeqDense,
1404:        MatLUFactor_SeqDense,
1405:        MatCholeskyFactor_SeqDense,
1406:        MatRelax_SeqDense,
1407:        MatTranspose_SeqDense,
1408: /*15*/ MatGetInfo_SeqDense,
1409:        MatEqual_SeqDense,
1410:        MatGetDiagonal_SeqDense,
1411:        MatDiagonalScale_SeqDense,
1412:        MatNorm_SeqDense,
1413: /*20*/ 0,
1414:        0,
1415:        0,
1416:        MatSetOption_SeqDense,
1417:        MatZeroEntries_SeqDense,
1418: /*25*/ MatZeroRows_SeqDense,
1419:        MatLUFactorSymbolic_SeqDense,
1420:        MatLUFactorNumeric_SeqDense,
1421:        MatCholeskyFactorSymbolic_SeqDense,
1422:        MatCholeskyFactorNumeric_SeqDense,
1423: /*30*/ MatSetUpPreallocation_SeqDense,
1424:        0,
1425:        0,
1426:        MatGetArray_SeqDense,
1427:        MatRestoreArray_SeqDense,
1428: /*35*/ MatDuplicate_SeqDense,
1429:        0,
1430:        0,
1431:        0,
1432:        0,
1433: /*40*/ MatAXPY_SeqDense,
1434:        MatGetSubMatrices_SeqDense,
1435:        0,
1436:        MatGetValues_SeqDense,
1437:        MatCopy_SeqDense,
1438: /*45*/ 0,
1439:        MatScale_SeqDense,
1440:        0,
1441:        0,
1442:        0,
1443: /*50*/ MatGetBlockSize_SeqDense,
1444:        0,
1445:        0,
1446:        0,
1447:        0,
1448: /*55*/ 0,
1449:        0,
1450:        0,
1451:        0,
1452:        0,
1453: /*60*/ 0,
1454:        MatDestroy_SeqDense,
1455:        MatView_SeqDense,
1456:        MatGetPetscMaps_Petsc,
1457:        0,
1458: /*65*/ 0,
1459:        0,
1460:        0,
1461:        0,
1462:        0,
1463: /*70*/ 0,
1464:        0,
1465:        0,
1466:        0,
1467:        0,
1468: /*75*/ 0,
1469:        0,
1470:        0,
1471:        0,
1472:        0,
1473: /*80*/ 0,
1474:        0,
1475:        0,
1476:        0,
1477:        0,
1478: /*85*/ MatLoad_SeqDense};

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

1487:    Collective on MPI_Comm

1489:    Input Parameters:
1490: +  comm - MPI communicator, set to PETSC_COMM_SELF
1491: .  m - number of rows
1492: .  n - number of columns
1493: -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1494:    to control all matrix memory allocation.

1496:    Output Parameter:
1497: .  A - the matrix

1499:    Notes:
1500:    The data input variable is intended primarily for Fortran programmers
1501:    who wish to allocate their own matrix memory space.  Most users should
1502:    set data=PETSC_NULL.

1504:    Level: intermediate

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

1508: .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
1509: @*/
1510: int MatCreateSeqDense(MPI_Comm comm,int m,int n,PetscScalar *data,Mat *A)
1511: {

1515:   MatCreate(comm,m,n,m,n,A);
1516:   MatSetType(*A,MATSEQDENSE);
1517:   MatSeqDenseSetPreallocation(*A,data);
1518:   return(0);
1519: }

1523: /*@C
1524:    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements

1526:    Collective on MPI_Comm

1528:    Input Parameters:
1529: +  A - the matrix
1530: -  data - the array (or PETSC_NULL)

1532:    Notes:
1533:    The data input variable is intended primarily for Fortran programmers
1534:    who wish to allocate their own matrix memory space.  Most users should
1535:    set data=PETSC_NULL.

1537:    Level: intermediate

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

1541: .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
1542: @*/
1543: int MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
1544: {
1545:   int ierr,(*f)(Mat,PetscScalar[]);

1548:   PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);
1549:   if (f) {
1550:     (*f)(B,data);
1551:   }
1552:   return(0);
1553: }

1555: EXTERN_C_BEGIN
1558: int MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
1559: {
1560:   Mat_SeqDense *b;
1561:   int          ierr;

1564:   B->preallocated = PETSC_TRUE;
1565:   b               = (Mat_SeqDense*)B->data;
1566:   if (!data) {
1567:     PetscMalloc((B->m*B->n+1)*sizeof(PetscScalar),&b->v);
1568:     PetscMemzero(b->v,B->m*B->n*sizeof(PetscScalar));
1569:     b->user_alloc = PETSC_FALSE;
1570:     PetscLogObjectMemory(B,B->n*B->m*sizeof(PetscScalar));
1571:   } else { /* user-allocated storage */
1572:     b->v          = data;
1573:     b->user_alloc = PETSC_TRUE;
1574:   }
1575:   return(0);
1576: }
1577: EXTERN_C_END

1581: /*@C
1582:   MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array

1584:   Input parameter:
1585: + A - the matrix
1586: - lda - the leading dimension

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

1593:   Level: intermediate

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

1597: .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation()
1598: @*/
1599: int MatSeqDenseSetLDA(Mat B,int lda)
1600: {
1601:   Mat_SeqDense *b = (Mat_SeqDense*)B->data;
1603:   if (lda<B->m) SETERRQ(1,"LDA must be at least matrix i dimension");
1604:   b->lda = lda;
1605:   return(0);
1606: }

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

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

1614:   Level: beginner

1616: .seealso: MatCreateSeqDense
1617: M*/

1619: EXTERN_C_BEGIN
1622: int MatCreate_SeqDense(Mat B)
1623: {
1624:   Mat_SeqDense *b;
1625:   int          ierr,size;

1628:   MPI_Comm_size(B->comm,&size);
1629:   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");

1631:   B->m = B->M = PetscMax(B->m,B->M);
1632:   B->n = B->N = PetscMax(B->n,B->N);

1634:   PetscNew(Mat_SeqDense,&b);
1635:   PetscMemzero(b,sizeof(Mat_SeqDense));
1636:   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
1637:   B->factor       = 0;
1638:   B->mapping      = 0;
1639:   PetscLogObjectMemory(B,sizeof(struct _p_Mat));
1640:   B->data         = (void*)b;

1642:   PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);
1643:   PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);

1645:   b->pivots       = 0;
1646:   b->roworiented  = PETSC_TRUE;
1647:   b->v            = 0;
1648:   b->lda          = B->m;

1650:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C",
1651:                                     "MatSeqDenseSetPreallocation_SeqDense",
1652:                                      MatSeqDenseSetPreallocation_SeqDense);
1653:   return(0);
1654: }
1655: EXTERN_C_END