Actual source code: dense.c

  1: /*$Id: dense.c,v 1.199 2001/04/10 19:35:15 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"

  9: int MatAXPY_SeqDense(Scalar *alpha,Mat X,Mat Y)
 10: {
 11:   Mat_SeqDense *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data;
 12:   int          N = X->m*X->n,one = 1;

 15:   BLaxpy_(&N,alpha,x->v,&one,y->v,&one);
 16:   PetscLogFlops(2*N-1);
 17:   return(0);
 18: }

 20: int MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
 21: {
 22:   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
 23:   int          i,N = A->m*A->n,count = 0;
 24:   Scalar       *v = a->v;

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

 29:   info->rows_global       = (double)A->m;
 30:   info->columns_global    = (double)A->n;
 31:   info->rows_local        = (double)A->m;
 32:   info->columns_local     = (double)A->n;
 33:   info->block_size        = 1.0;
 34:   info->nz_allocated      = (double)N;
 35:   info->nz_used           = (double)count;
 36:   info->nz_unneeded       = (double)(N-count);
 37:   info->assemblies        = (double)A->num_ass;
 38:   info->mallocs           = 0;
 39:   info->memory            = A->mem;
 40:   info->fill_ratio_given  = 0;
 41:   info->fill_ratio_needed = 0;
 42:   info->factor_mallocs    = 0;

 44:   return(0);
 45: }

 47: int MatScale_SeqDense(Scalar *alpha,Mat A)
 48: {
 49:   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
 50:   int          one = 1,nz;

 53:   nz = A->m*A->n;
 54:   BLscal_(&nz,alpha,a->v,&one);
 55:   PetscLogFlops(nz);
 56:   return(0);
 57: }
 58: 
 59: /* ---------------------------------------------------------------*/
 60: /* COMMENT: I have chosen to hide row permutation in the pivots,
 61:    rather than put it in the Mat->row slot.*/
 62: int MatLUFactor_SeqDense(Mat A,IS row,IS col,MatLUInfo *minfo)
 63: {
 64:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
 65:   int          info,ierr;

 68:   if (!mat->pivots) {
 69:     PetscMalloc((A->m+1)*sizeof(int),&mat->pivots);
 70:     PetscLogObjectMemory(A,A->m*sizeof(int));
 71:   }
 72:   A->factor = FACTOR_LU;
 73:   if (!A->m || !A->n) return(0);
 74: #if defined(PETSC_MISSING_LAPACK_GETRF) 
 75:   SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavilable.");
 76: #else
 77:   LAgetrf_(&A->m,&A->n,mat->v,&A->m,mat->pivots,&info);
 78: #endif
 79:   if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization");
 80:   if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
 81:   PetscLogFlops((2*A->n*A->n*A->n)/3);
 82:   return(0);
 83: }

 85: int MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
 86: {
 87:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l;
 88:   int          ierr;
 89:   Mat          newi;

 92:   MatCreateSeqDense(A->comm,A->m,A->n,PETSC_NULL,&newi);
 93:   l = (Mat_SeqDense*)newi->data;
 94:   if (cpvalues == MAT_COPY_VALUES) {
 95:     PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(Scalar));
 96:   }
 97:   newi->assembled = PETSC_TRUE;
 98:   *newmat = newi;
 99:   return(0);
100: }

102: int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,MatLUInfo *info,Mat *fact)
103: {

107:   MatDuplicate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);
108:   return(0);
109: }

111: int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact)
112: {
113:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l = (Mat_SeqDense*)(*fact)->data;
114:   int          ierr;

117:   /* copy the numerical values */
118:   PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(Scalar));
119:   (*fact)->factor = 0;
120:   MatLUFactor(*fact,0,0,PETSC_NULL);
121:   return(0);
122: }

124: int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,PetscReal f,Mat *fact)
125: {

129:   MatConvert(A,MATSAME,fact);
130:   return(0);
131: }

133: int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact)
134: {

138:   MatCholeskyFactor(*fact,0,1.0);
139:   return(0);
140: }

142: int MatCholeskyFactor_SeqDense(Mat A,IS perm,PetscReal f)
143: {
144:   Mat_SeqDense  *mat = (Mat_SeqDense*)A->data;
145:   int           info,ierr;
146: 
148:   if (mat->pivots) {
149:     PetscFree(mat->pivots);
150:     PetscLogObjectMemory(A,-A->m*sizeof(int));
151:     mat->pivots = 0;
152:   }
153:   if (!A->m || !A->n) return(0);
154: #if defined(PETSC_MISSING_LAPACK_POTRF) 
155:   SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavilable.");
156: #else
157:   LApotrf_("L",&A->n,mat->v,&A->m,&info);
158: #endif
159:   if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %d",info-1);
160:   A->factor = FACTOR_CHOLESKY;
161:   PetscLogFlops((A->n*A->n*A->n)/3);
162:   return(0);
163: }

165: int MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
166: {
167:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
168:   int          one = 1,info,ierr;
169:   Scalar       *x,*y;
170: 
172:   if (!A->m || !A->n) return(0);
173:   VecGetArray(xx,&x);
174:   VecGetArray(yy,&y);
175:   PetscMemcpy(y,x,A->m*sizeof(Scalar));
176:   if (A->factor == FACTOR_LU) {
177: #if defined(PETSC_MISSING_LAPACK_GETRS) 
178:   SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable.");
179: #else
180:     LAgetrs_("N",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info);
181: #endif
182:   } else if (A->factor == FACTOR_CHOLESKY){
183: #if defined(PETSC_MISSING_LAPACK_POTRS) 
184:   SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable.");
185: #else
186:     LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info);
187: #endif
188:   }
189:   else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
190:   if (info) SETERRQ(PETSC_ERR_LIB,"MBad solve");
191:   VecRestoreArray(xx,&x);
192:   VecRestoreArray(yy,&y);
193:   PetscLogFlops(2*A->n*A->n - A->n);
194:   return(0);
195: }

197: int MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
198: {
199:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
200:   int          ierr,one = 1,info;
201:   Scalar       *x,*y;
202: 
204:   if (!A->m || !A->n) return(0);
205:   VecGetArray(xx,&x);
206:   VecGetArray(yy,&y);
207:   PetscMemcpy(y,x,A->m*sizeof(Scalar));
208:   /* assume if pivots exist then use LU; else Cholesky */
209:   if (mat->pivots) {
210: #if defined(PETSC_MISSING_LAPACK_GETRS) 
211:   SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable.");
212: #else
213:     LAgetrs_("T",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info);
214: #endif
215:   } else {
216: #if defined(PETSC_MISSING_LAPACK_POTRS) 
217:   SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable.");
218: #else
219:     LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info);
220: #endif
221:   }
222:   if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
223:   VecRestoreArray(xx,&x);
224:   VecRestoreArray(yy,&y);
225:   PetscLogFlops(2*A->n*A->n - A->n);
226:   return(0);
227: }

229: int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
230: {
231:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
232:   int          ierr,one = 1,info;
233:   Scalar       *x,*y,sone = 1.0;
234:   Vec          tmp = 0;
235: 
237:   VecGetArray(xx,&x);
238:   VecGetArray(yy,&y);
239:   if (!A->m || !A->n) return(0);
240:   if (yy == zz) {
241:     VecDuplicate(yy,&tmp);
242:     PetscLogObjectParent(A,tmp);
243:     VecCopy(yy,tmp);
244:   }
245:   PetscMemcpy(y,x,A->m*sizeof(Scalar));
246:   /* assume if pivots exist then use LU; else Cholesky */
247:   if (mat->pivots) {
248: #if defined(PETSC_MISSING_LAPACK_GETRS) 
249:   SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable.");
250: #else
251:     LAgetrs_("N",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info);
252: #endif
253:   } else {
254: #if defined(PETSC_MISSING_LAPACK_POTRS) 
255:   SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable.");
256: #else
257:     LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info);
258: #endif
259:   }
260:   if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
261:   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
262:   else     {VecAXPY(&sone,zz,yy);}
263:   VecRestoreArray(xx,&x);
264:   VecRestoreArray(yy,&y);
265:   PetscLogFlops(2*A->n*A->n);
266:   return(0);
267: }

269: int MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
270: {
271:   Mat_SeqDense  *mat = (Mat_SeqDense*)A->data;
272:   int           one = 1,info,ierr;
273:   Scalar        *x,*y,sone = 1.0;
274:   Vec           tmp;
275: 
277:   if (!A->m || !A->n) return(0);
278:   VecGetArray(xx,&x);
279:   VecGetArray(yy,&y);
280:   if (yy == zz) {
281:     VecDuplicate(yy,&tmp);
282:     PetscLogObjectParent(A,tmp);
283:     VecCopy(yy,tmp);
284:   }
285:   PetscMemcpy(y,x,A->m*sizeof(Scalar));
286:   /* assume if pivots exist then use LU; else Cholesky */
287:   if (mat->pivots) {
288: #if defined(PETSC_MISSING_LAPACK_GETRS) 
289:   SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable.");
290: #else
291:     LAgetrs_("T",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info);
292: #endif
293:   } else {
294: #if defined(PETSC_MISSING_LAPACK_POTRS) 
295:   SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable.");
296: #else
297:     LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info);
298: #endif
299:   }
300:   if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
301:   if (tmp) {
302:     VecAXPY(&sone,tmp,yy);
303:     VecDestroy(tmp);
304:   } else {
305:     VecAXPY(&sone,zz,yy);
306:   }
307:   VecRestoreArray(xx,&x);
308:   VecRestoreArray(yy,&y);
309:   PetscLogFlops(2*A->n*A->n);
310:   return(0);
311: }
312: /* ------------------------------------------------------------------*/
313: int MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,
314:                           PetscReal shift,int its,Vec xx)
315: {
316:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
317:   Scalar       *x,*b,*v = mat->v,zero = 0.0,xt;
318:   int          ierr,m = A->m,i;
319: #if !defined(PETSC_USE_COMPLEX)
320:   int          o = 1;
321: #endif

324:   if (flag & SOR_ZERO_INITIAL_GUESS) {
325:     /* this is a hack fix, should have another version without the second BLdot */
326:     VecSet(&zero,xx);
327:   }
328:   VecGetArray(xx,&x);
329:   VecGetArray(bb,&b);
330:   while (its--) {
331:     if (flag & SOR_FORWARD_SWEEP){
332:       for (i=0; i<m; i++) {
333: #if defined(PETSC_USE_COMPLEX)
334:         /* cannot use BLAS dot for complex because compiler/linker is 
335:            not happy about returning a double complex */
336:         int    _i;
337:         Scalar sum = b[i];
338:         for (_i=0; _i<m; _i++) {
339:           sum -= PetscConj(v[i+_i*m])*x[_i];
340:         }
341:         xt = sum;
342: #else
343:         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
344: #endif
345:         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
346:       }
347:     }
348:     if (flag & SOR_BACKWARD_SWEEP) {
349:       for (i=m-1; i>=0; i--) {
350: #if defined(PETSC_USE_COMPLEX)
351:         /* cannot use BLAS dot for complex because compiler/linker is 
352:            not happy about returning a double complex */
353:         int    _i;
354:         Scalar sum = b[i];
355:         for (_i=0; _i<m; _i++) {
356:           sum -= PetscConj(v[i+_i*m])*x[_i];
357:         }
358:         xt = sum;
359: #else
360:         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
361: #endif
362:         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
363:       }
364:     }
365:   }
366:   VecRestoreArray(bb,&b);
367:   VecRestoreArray(xx,&x);
368:   return(0);
369: }

371: /* -----------------------------------------------------------------*/
372: int MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
373: {
374:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
375:   Scalar       *v = mat->v,*x,*y;
376:   int          ierr,_One=1;Scalar _DOne=1.0,_DZero=0.0;

379:   if (!A->m || !A->n) return(0);
380:   VecGetArray(xx,&x);
381:   VecGetArray(yy,&y);
382:   LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DZero,y,&_One);
383:   VecRestoreArray(xx,&x);
384:   VecRestoreArray(yy,&y);
385:   PetscLogFlops(2*A->m*A->n - A->n);
386:   return(0);
387: }

389: int MatMult_SeqDense(Mat A,Vec xx,Vec yy)
390: {
391:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
392:   Scalar       *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0;
393:   int          ierr,_One=1;

396:   if (!A->m || !A->n) return(0);
397:   VecGetArray(xx,&x);
398:   VecGetArray(yy,&y);
399:   LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DZero,y,&_One);
400:   VecRestoreArray(xx,&x);
401:   VecRestoreArray(yy,&y);
402:   PetscLogFlops(2*A->m*A->n - A->m);
403:   return(0);
404: }

406: int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
407: {
408:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
409:   Scalar       *v = mat->v,*x,*y,_DOne=1.0;
410:   int          ierr,_One=1;

413:   if (!A->m || !A->n) return(0);
414:   if (zz != yy) {VecCopy(zz,yy);}
415:   VecGetArray(xx,&x);
416:   VecGetArray(yy,&y);
417:   LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DOne,y,&_One);
418:   VecRestoreArray(xx,&x);
419:   VecRestoreArray(yy,&y);
420:   PetscLogFlops(2*A->m*A->n);
421:   return(0);
422: }

424: int MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
425: {
426:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
427:   Scalar       *v = mat->v,*x,*y;
428:   int          ierr,_One=1;
429:   Scalar       _DOne=1.0;

432:   if (!A->m || !A->n) return(0);
433:   if (zz != yy) {VecCopy(zz,yy);}
434:   VecGetArray(xx,&x);
435:   VecGetArray(yy,&y);
436:   LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DOne,y,&_One);
437:   VecRestoreArray(xx,&x);
438:   VecRestoreArray(yy,&y);
439:   PetscLogFlops(2*A->m*A->n);
440:   return(0);
441: }

443: /* -----------------------------------------------------------------*/
444: int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
445: {
446:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
447:   Scalar       *v;
448:   int          i,ierr;
449: 
451:   *ncols = A->n;
452:   if (cols) {
453:     PetscMalloc((A->n+1)*sizeof(int),cols);
454:     for (i=0; i<A->n; i++) (*cols)[i] = i;
455:   }
456:   if (vals) {
457:     PetscMalloc((A->n+1)*sizeof(Scalar),vals);
458:     v    = mat->v + row;
459:     for (i=0; i<A->n; i++) {(*vals)[i] = *v; v += A->m;}
460:   }
461:   return(0);
462: }

464: int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
465: {
468:   if (cols) {PetscFree(*cols);}
469:   if (vals) {PetscFree(*vals); }
470:   return(0);
471: }
472: /* ----------------------------------------------------------------*/
473: int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n,
474:                                     int *indexn,Scalar *v,InsertMode addv)
475: {
476:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
477:   int          i,j;
478: 
480:   if (!mat->roworiented) {
481:     if (addv == INSERT_VALUES) {
482:       for (j=0; j<n; j++) {
483:         if (indexn[j] < 0) {v += m; continue;}
484: #if defined(PETSC_USE_BOPT_g)  
485:         if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
486: #endif
487:         for (i=0; i<m; i++) {
488:           if (indexm[i] < 0) {v++; continue;}
489: #if defined(PETSC_USE_BOPT_g)  
490:           if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
491: #endif
492:           mat->v[indexn[j]*A->m + indexm[i]] = *v++;
493:         }
494:       }
495:     } else {
496:       for (j=0; j<n; j++) {
497:         if (indexn[j] < 0) {v += m; continue;}
498: #if defined(PETSC_USE_BOPT_g)  
499:         if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
500: #endif
501:         for (i=0; i<m; i++) {
502:           if (indexm[i] < 0) {v++; continue;}
503: #if defined(PETSC_USE_BOPT_g)  
504:           if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
505: #endif
506:           mat->v[indexn[j]*A->m + indexm[i]] += *v++;
507:         }
508:       }
509:     }
510:   } else {
511:     if (addv == INSERT_VALUES) {
512:       for (i=0; i<m; i++) {
513:         if (indexm[i] < 0) { v += n; continue;}
514: #if defined(PETSC_USE_BOPT_g)  
515:         if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
516: #endif
517:         for (j=0; j<n; j++) {
518:           if (indexn[j] < 0) { v++; continue;}
519: #if defined(PETSC_USE_BOPT_g)  
520:           if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
521: #endif
522:           mat->v[indexn[j]*A->m + indexm[i]] = *v++;
523:         }
524:       }
525:     } else {
526:       for (i=0; i<m; i++) {
527:         if (indexm[i] < 0) { v += n; continue;}
528: #if defined(PETSC_USE_BOPT_g)  
529:         if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
530: #endif
531:         for (j=0; j<n; j++) {
532:           if (indexn[j] < 0) { v++; continue;}
533: #if defined(PETSC_USE_BOPT_g)  
534:           if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
535: #endif
536:           mat->v[indexn[j]*A->m + indexm[i]] += *v++;
537:         }
538:       }
539:     }
540:   }
541:   return(0);
542: }

544: int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v)
545: {
546:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
547:   int          i,j;
548:   Scalar       *vpt = v;

551:   /* row-oriented output */
552:   for (i=0; i<m; i++) {
553:     for (j=0; j<n; j++) {
554:       *vpt++ = mat->v[indexn[j]*A->m + indexm[i]];
555:     }
556:   }
557:   return(0);
558: }

560: /* -----------------------------------------------------------------*/

562: #include "petscsys.h"

564: EXTERN_C_BEGIN
565: int MatLoad_SeqDense(PetscViewer viewer,MatType type,Mat *A)
566: {
567:   Mat_SeqDense *a;
568:   Mat          B;
569:   int          *scols,i,j,nz,ierr,fd,header[4],size;
570:   int          *rowlengths = 0,M,N,*cols;
571:   Scalar       *vals,*svals,*v,*w;
572:   MPI_Comm     comm = ((PetscObject)viewer)->comm;

575:   MPI_Comm_size(comm,&size);
576:   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
577:   PetscViewerBinaryGetDescriptor(viewer,&fd);
578:   PetscBinaryRead(fd,header,4,PETSC_INT);
579:   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
580:   M = header[1]; N = header[2]; nz = header[3];

582:   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
583:     MatCreateSeqDense(comm,M,N,PETSC_NULL,A);
584:     B    = *A;
585:     a    = (Mat_SeqDense*)B->data;
586:     v    = a->v;
587:     /* Allocate some temp space to read in the values and then flip them
588:        from row major to column major */
589:     PetscMalloc((M*N+1)*sizeof(Scalar),&w);
590:     /* read in nonzero values */
591:     PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);
592:     /* now flip the values and store them in the matrix*/
593:     for (j=0; j<N; j++) {
594:       for (i=0; i<M; i++) {
595:         *v++ =w[i*N+j];
596:       }
597:     }
598:     PetscFree(w);
599:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
600:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
601:   } else {
602:     /* read row lengths */
603:     PetscMalloc((M+1)*sizeof(int),&rowlengths);
604:     PetscBinaryRead(fd,rowlengths,M,PETSC_INT);

606:     /* create our matrix */
607:     MatCreateSeqDense(comm,M,N,PETSC_NULL,A);
608:     B = *A;
609:     a = (Mat_SeqDense*)B->data;
610:     v = a->v;

612:     /* read column indices and nonzeros */
613:     PetscMalloc((nz+1)*sizeof(int),&scols);
614:     cols = scols;
615:     PetscBinaryRead(fd,cols,nz,PETSC_INT);
616:     PetscMalloc((nz+1)*sizeof(Scalar),&svals);
617:     vals = svals;
618:     PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);

620:     /* insert into matrix */
621:     for (i=0; i<M; i++) {
622:       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
623:       svals += rowlengths[i]; scols += rowlengths[i];
624:     }
625:     PetscFree(vals);
626:     PetscFree(cols);
627:     PetscFree(rowlengths);

629:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
630:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
631:   }
632:   return(0);
633: }
634: EXTERN_C_END

636: #include "petscsys.h"

638: static int MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
639: {
640:   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
641:   int               ierr,i,j;
642:   char              *name;
643:   Scalar            *v;
644:   PetscViewerFormat format;

647:   PetscObjectGetName((PetscObject)A,&name);
648:   PetscViewerGetFormat(viewer,&format);
649:   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_LONG) {
650:     return(0);  /* do nothing for now */
651:   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
652:     PetscViewerASCIIUseTabs(viewer,PETSC_NO);
653:     for (i=0; i<A->m; i++) {
654:       v = a->v + i;
655:       PetscViewerASCIIPrintf(viewer,"row %d:",i);
656:       for (j=0; j<A->n; j++) {
657: #if defined(PETSC_USE_COMPLEX)
658:         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
659:           PetscViewerASCIIPrintf(viewer," %d %g + %g i",j,PetscRealPart(*v),PetscImaginaryPart(*v));
660:         } else if (PetscRealPart(*v)) {
661:           PetscViewerASCIIPrintf(viewer," %d %g ",j,PetscRealPart(*v));
662:         }
663: #else
664:         if (*v) {
665:           PetscViewerASCIIPrintf(viewer," %d %g ",j,*v);
666:         }
667: #endif
668:         v += A->m;
669:       }
670:       PetscViewerASCIIPrintf(viewer,"n");
671:     }
672:     PetscViewerASCIIUseTabs(viewer,PETSC_YES);
673:   } else {
674:     PetscViewerASCIIUseTabs(viewer,PETSC_NO);
675: #if defined(PETSC_USE_COMPLEX)
676:     PetscTruth allreal = PETSC_TRUE;
677:     /* determine if matrix has all real values */
678:     v = a->v;
679:     for (i=0; i<A->m*A->n; i++) {
680:       if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;}
681:     }
682: #endif
683:     if (format == PETSC_VIEWER_ASCII_MATLAB) {
684:       PetscObjectGetName((PetscObject)viewer,&name);
685:       PetscViewerASCIIPrintf(viewer,"%% Size = %d %d n",A->m,A->n);
686:       PetscViewerASCIIPrintf(viewer,"%s = zeros(%d,%d);n",name,A->m,A->n);
687:       PetscViewerASCIIPrintf(viewer,"%s = [n",name);
688:     }

690:     for (i=0; i<A->m; i++) {
691:       v = a->v + i;
692:       for (j=0; j<A->n; j++) {
693: #if defined(PETSC_USE_COMPLEX)
694:         if (allreal) {
695:           PetscViewerASCIIPrintf(viewer,"%6.4e ",PetscRealPart(*v)); v += A->m;
696:         } else {
697:           PetscViewerASCIIPrintf(viewer,"%6.4e + %6.4e i ",PetscRealPart(*v),PetscImaginaryPart(*v)); v += A->m;
698:         }
699: #else
700:         PetscViewerASCIIPrintf(viewer,"%6.4e ",*v); v += A->m;
701: #endif
702:       }
703:       PetscViewerASCIIPrintf(viewer,"n");
704:     }
705:     if (format == PETSC_VIEWER_ASCII_MATLAB) {
706:       PetscViewerASCIIPrintf(viewer,"];n");
707:     }
708:     PetscViewerASCIIUseTabs(viewer,PETSC_YES);
709:   }
710:   PetscViewerFlush(viewer);
711:   return(0);
712: }

714: static int MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
715: {
716:   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
717:   int               ict,j,n = A->n,m = A->m,i,fd,*col_lens,ierr,nz = m*n;
718:   Scalar            *v,*anonz,*vals;
719:   PetscViewerFormat format;
720: 
722:   PetscViewerBinaryGetDescriptor(viewer,&fd);

724:   PetscViewerGetFormat(viewer,&format);
725:   if (format == PETSC_VIEWER_BINARY_NATIVE) {
726:     /* store the matrix as a dense matrix */
727:     PetscMalloc(4*sizeof(int),&col_lens);
728:     col_lens[0] = MAT_COOKIE;
729:     col_lens[1] = m;
730:     col_lens[2] = n;
731:     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
732:     PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1);
733:     PetscFree(col_lens);

735:     /* write out matrix, by rows */
736:     PetscMalloc((m*n+1)*sizeof(Scalar),&vals);
737:     v    = a->v;
738:     for (i=0; i<m; i++) {
739:       for (j=0; j<n; j++) {
740:         vals[i + j*m] = *v++;
741:       }
742:     }
743:     PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0);
744:     PetscFree(vals);
745:   } else {
746:     PetscMalloc((4+nz)*sizeof(int),&col_lens);
747:     col_lens[0] = MAT_COOKIE;
748:     col_lens[1] = m;
749:     col_lens[2] = n;
750:     col_lens[3] = nz;

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

756:     /* Possibly should write in smaller increments, not whole matrix at once? */
757:     /* store column indices (zero start index) */
758:     ict = 0;
759:     for (i=0; i<m; i++) {
760:       for (j=0; j<n; j++) col_lens[ict++] = j;
761:     }
762:     PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0);
763:     PetscFree(col_lens);

765:     /* store nonzero values */
766:     PetscMalloc((nz+1)*sizeof(Scalar),&anonz);
767:     ict  = 0;
768:     for (i=0; i<m; i++) {
769:       v = a->v + i;
770:       for (j=0; j<n; j++) {
771:         anonz[ict++] = *v; v += A->m;
772:       }
773:     }
774:     PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0);
775:     PetscFree(anonz);
776:   }
777:   return(0);
778: }

780: int MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
781: {
782:   Mat               A = (Mat) Aa;
783:   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
784:   int               m = A->m,n = A->n,color,i,j,ierr;
785:   Scalar            *v = a->v;
786:   PetscViewer       viewer;
787:   PetscDraw         popup;
788:   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0;
789:   PetscViewerFormat format;


793:   PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);
794:   PetscViewerGetFormat(viewer,&format);
795:   PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);

797:   /* Loop over matrix elements drawing boxes */
798:   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
799:     /* Blue for negative and Red for positive */
800:     color = PETSC_DRAW_BLUE;
801:     for(j = 0; j < n; j++) {
802:       x_l = j;
803:       x_r = x_l + 1.0;
804:       for(i = 0; i < m; i++) {
805:         y_l = m - i - 1.0;
806:         y_r = y_l + 1.0;
807: #if defined(PETSC_USE_COMPLEX)
808:         if (PetscRealPart(v[j*m+i]) >  0.) {
809:           color = PETSC_DRAW_RED;
810:         } else if (PetscRealPart(v[j*m+i]) <  0.) {
811:           color = PETSC_DRAW_BLUE;
812:         } else {
813:           continue;
814:         }
815: #else
816:         if (v[j*m+i] >  0.) {
817:           color = PETSC_DRAW_RED;
818:         } else if (v[j*m+i] <  0.) {
819:           color = PETSC_DRAW_BLUE;
820:         } else {
821:           continue;
822:         }
823: #endif
824:         PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
825:       }
826:     }
827:   } else {
828:     /* use contour shading to indicate magnitude of values */
829:     /* first determine max of all nonzero values */
830:     for(i = 0; i < m*n; i++) {
831:       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
832:     }
833:     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
834:     ierr  = PetscDrawGetPopup(draw,&popup);
835:     if (popup) {ierr  = PetscDrawScalePopup(popup,0.0,maxv);}
836:     for(j = 0; j < n; j++) {
837:       x_l = j;
838:       x_r = x_l + 1.0;
839:       for(i = 0; i < m; i++) {
840:         y_l   = m - i - 1.0;
841:         y_r   = y_l + 1.0;
842:         color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i]));
843:         ierr  = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
844:       }
845:     }
846:   }
847:   return(0);
848: }

850: int MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
851: {
852:   PetscDraw  draw;
853:   PetscTruth isnull;
854:   PetscReal  xr,yr,xl,yl,h,w;
855:   int        ierr;

858:   PetscViewerDrawGetDraw(viewer,0,&draw);
859:   PetscDrawIsNull(draw,&isnull);
860:   if (isnull == PETSC_TRUE) return(0);

862:   PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
863:   xr  = A->n; yr = A->m; h = yr/10.0; w = xr/10.0;
864:   xr += w;    yr += h;  xl = -w;     yl = -h;
865:   PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
866:   PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);
867:   PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);
868:   return(0);
869: }

871: int MatView_SeqDense(Mat A,PetscViewer viewer)
872: {
873:   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
874:   int          ierr;
875:   PetscTruth   issocket,isascii,isbinary,isdraw;

878:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);
879:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
880:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
881:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);

883:   if (issocket) {
884:     PetscViewerSocketPutScalar(viewer,A->m,A->n,a->v);
885:   } else if (isascii) {
886:     MatView_SeqDense_ASCII(A,viewer);
887:   } else if (isbinary) {
888:     MatView_SeqDense_Binary(A,viewer);
889:   } else if (isdraw) {
890:     MatView_SeqDense_Draw(A,viewer);
891:   } else {
892:     SETERRQ1(1,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name);
893:   }
894:   return(0);
895: }

897: int MatDestroy_SeqDense(Mat mat)
898: {
899:   Mat_SeqDense *l = (Mat_SeqDense*)mat->data;
900:   int          ierr;

903: #if defined(PETSC_USE_LOG)
904:   PetscLogObjectState((PetscObject)mat,"Rows %d Cols %d",mat->m,mat->n);
905: #endif
906:   if (l->pivots) {PetscFree(l->pivots);}
907:   if (!l->user_alloc) {PetscFree(l->v);}
908:   PetscFree(l);
909:   return(0);
910: }

912: int MatTranspose_SeqDense(Mat A,Mat *matout)
913: {
914:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
915:   int          k,j,m,n,ierr;
916:   Scalar       *v,tmp;

919:   v = mat->v; m = A->m; n = A->n;
920:   if (!matout) { /* in place transpose */
921:     if (m != n) { /* malloc temp to hold transpose */
922:       Scalar *w;
923:       PetscMalloc((m+1)*(n+1)*sizeof(Scalar),&w);
924:       for (j=0; j<m; j++) {
925:         for (k=0; k<n; k++) {
926:           w[k + j*n] = v[j + k*m];
927:         }
928:       }
929:       PetscMemcpy(v,w,m*n*sizeof(Scalar));
930:       PetscFree(w);
931:     } else {
932:       for (j=0; j<m; j++) {
933:         for (k=0; k<j; k++) {
934:           tmp = v[j + k*n];
935:           v[j + k*n] = v[k + j*n];
936:           v[k + j*n] = tmp;
937:         }
938:       }
939:     }
940:   } else { /* out-of-place transpose */
941:     Mat          tmat;
942:     Mat_SeqDense *tmatd;
943:     Scalar       *v2;
944:     ierr  = MatCreateSeqDense(A->comm,A->n,A->m,PETSC_NULL,&tmat);
945:     tmatd = (Mat_SeqDense*)tmat->data;
946:     v = mat->v; v2 = tmatd->v;
947:     for (j=0; j<n; j++) {
948:       for (k=0; k<m; k++) v2[j + k*n] = v[k + j*m];
949:     }
950:     MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);
951:     MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);
952:     *matout = tmat;
953:   }
954:   return(0);
955: }

957: int MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg)
958: {
959:   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
960:   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
961:   int          ierr,i;
962:   Scalar       *v1 = mat1->v,*v2 = mat2->v;
963:   PetscTruth   flag;

966:   PetscTypeCompare((PetscObject)A2,MATSEQDENSE,&flag);
967:   if (!flag) SETERRQ(PETSC_ERR_SUP,"Matrices should be of type  MATSEQDENSE");
968:   if (A1->m != A2->m) {*flg = PETSC_FALSE; return(0);}
969:   if (A1->n != A2->n) {*flg = PETSC_FALSE; return(0);}
970:   for (i=0; i<A1->m*A1->n; i++) {
971:     if (*v1 != *v2) {*flg = PETSC_FALSE; return(0);}
972:     v1++; v2++;
973:   }
974:   *flg = PETSC_TRUE;
975:   return(0);
976: }

978: int MatGetDiagonal_SeqDense(Mat A,Vec v)
979: {
980:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
981:   int          ierr,i,n,len;
982:   Scalar       *x,zero = 0.0;

985:   VecSet(&zero,v);
986:   VecGetSize(v,&n);
987:   VecGetArray(v,&x);
988:   len = PetscMin(A->m,A->n);
989:   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
990:   for (i=0; i<len; i++) {
991:     x[i] = mat->v[i*A->m + i];
992:   }
993:   VecRestoreArray(v,&x);
994:   return(0);
995: }

997: int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
998: {
999:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1000:   Scalar       *l,*r,x,*v;
1001:   int          ierr,i,j,m = A->m,n = A->n;

1004:   if (ll) {
1005:     VecGetSize(ll,&m);
1006:     VecGetArray(ll,&l);
1007:     if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1008:     for (i=0; i<m; i++) {
1009:       x = l[i];
1010:       v = mat->v + i;
1011:       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1012:     }
1013:     VecRestoreArray(ll,&l);
1014:     PetscLogFlops(n*m);
1015:   }
1016:   if (rr) {
1017:     VecGetSize(rr,&n);
1018:     VecGetArray(rr,&r);
1019:     if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1020:     for (i=0; i<n; i++) {
1021:       x = r[i];
1022:       v = mat->v + i*m;
1023:       for (j=0; j<m; j++) { (*v++) *= x;}
1024:     }
1025:     VecRestoreArray(rr,&r);
1026:     PetscLogFlops(n*m);
1027:   }
1028:   return(0);
1029: }

1031: int MatNorm_SeqDense(Mat A,NormType type,PetscReal *norm)
1032: {
1033:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1034:   Scalar       *v = mat->v;
1035:   PetscReal    sum = 0.0;
1036:   int          i,j;

1039:   if (type == NORM_FROBENIUS) {
1040:     for (i=0; i<A->n*A->m; i++) {
1041: #if defined(PETSC_USE_COMPLEX)
1042:       sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1043: #else
1044:       sum += (*v)*(*v); v++;
1045: #endif
1046:     }
1047:     *norm = sqrt(sum);
1048:     PetscLogFlops(2*A->n*A->m);
1049:   } else if (type == NORM_1) {
1050:     *norm = 0.0;
1051:     for (j=0; j<A->n; j++) {
1052:       sum = 0.0;
1053:       for (i=0; i<A->m; i++) {
1054:         sum += PetscAbsScalar(*v);  v++;
1055:       }
1056:       if (sum > *norm) *norm = sum;
1057:     }
1058:     PetscLogFlops(A->n*A->m);
1059:   } else if (type == NORM_INFINITY) {
1060:     *norm = 0.0;
1061:     for (j=0; j<A->m; j++) {
1062:       v = mat->v + j;
1063:       sum = 0.0;
1064:       for (i=0; i<A->n; i++) {
1065:         sum += PetscAbsScalar(*v); v += A->m;
1066:       }
1067:       if (sum > *norm) *norm = sum;
1068:     }
1069:     PetscLogFlops(A->n*A->m);
1070:   } else {
1071:     SETERRQ(PETSC_ERR_SUP,"No two norm");
1072:   }
1073:   return(0);
1074: }

1076: int MatSetOption_SeqDense(Mat A,MatOption op)
1077: {
1078:   Mat_SeqDense *aij = (Mat_SeqDense*)A->data;
1079: 
1081:   if (op == MAT_ROW_ORIENTED)            aij->roworiented = PETSC_TRUE;
1082:   else if (op == MAT_COLUMN_ORIENTED)    aij->roworiented = PETSC_FALSE;
1083:   else if (op == MAT_ROWS_SORTED ||
1084:            op == MAT_ROWS_UNSORTED ||
1085:            op == MAT_COLUMNS_SORTED ||
1086:            op == MAT_COLUMNS_UNSORTED ||
1087:            op == MAT_NO_NEW_NONZERO_LOCATIONS ||
1088:            op == MAT_YES_NEW_NONZERO_LOCATIONS ||
1089:            op == MAT_NEW_NONZERO_LOCATION_ERR ||
1090:            op == MAT_NO_NEW_DIAGONALS ||
1091:            op == MAT_YES_NEW_DIAGONALS ||
1092:            op == MAT_IGNORE_OFF_PROC_ENTRIES ||
1093:            op == MAT_USE_HASH_TABLE) {
1094:     PetscLogInfo(A,"MatSetOption_SeqDense:Option ignoredn");
1095:   } else {
1096:     SETERRQ(PETSC_ERR_SUP,"unknown option");
1097:   }
1098:   return(0);
1099: }

1101: int MatZeroEntries_SeqDense(Mat A)
1102: {
1103:   Mat_SeqDense *l = (Mat_SeqDense*)A->data;
1104:   int          ierr;

1107:   PetscMemzero(l->v,A->m*A->n*sizeof(Scalar));
1108:   return(0);
1109: }

1111: int MatGetBlockSize_SeqDense(Mat A,int *bs)
1112: {
1114:   *bs = 1;
1115:   return(0);
1116: }

1118: int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
1119: {
1120:   Mat_SeqDense *l = (Mat_SeqDense*)A->data;
1121:   int          n = A->n,i,j,ierr,N,*rows;
1122:   Scalar       *slot;

1125:   ISGetLocalSize(is,&N);
1126:   ISGetIndices(is,&rows);
1127:   for (i=0; i<N; i++) {
1128:     slot = l->v + rows[i];
1129:     for (j=0; j<n; j++) { *slot = 0.0; slot += n;}
1130:   }
1131:   if (diag) {
1132:     for (i=0; i<N; i++) {
1133:       slot = l->v + (n+1)*rows[i];
1134:       *slot = *diag;
1135:     }
1136:   }
1137:   ISRestoreIndices(is,&rows);
1138:   return(0);
1139: }

1141: int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
1142: {
1144:   if (m) *m = 0;
1145:   if (n) *n = A->m;
1146:   return(0);
1147: }

1149: int MatGetArray_SeqDense(Mat A,Scalar **array)
1150: {
1151:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;

1154:   *array = mat->v;
1155:   return(0);
1156: }

1158: int MatRestoreArray_SeqDense(Mat A,Scalar **array)
1159: {
1161:   *array = 0; /* user cannot accidently use the array later */
1162:   return(0);
1163: }

1165: static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
1166: {
1167:   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1168:   int          i,j,ierr,m = A->m,*irow,*icol,nrows,ncols;
1169:   Scalar       *av,*bv,*v = mat->v;
1170:   Mat          newmat;

1173:   ISGetIndices(isrow,&irow);
1174:   ISGetIndices(iscol,&icol);
1175:   ISGetLocalSize(isrow,&nrows);
1176:   ISGetLocalSize(iscol,&ncols);
1177: 
1178:   /* Check submatrixcall */
1179:   if (scall == MAT_REUSE_MATRIX) {
1180:     int n_cols,n_rows;
1181:     MatGetSize(*B,&n_rows,&n_cols);
1182:     if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
1183:     newmat = *B;
1184:   } else {
1185:     /* Create and fill new matrix */
1186:     MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat);
1187:   }

1189:   /* Now extract the data pointers and do the copy,column at a time */
1190:   bv = ((Mat_SeqDense*)newmat->data)->v;
1191: 
1192:   for (i=0; i<ncols; i++) {
1193:     av = v + m*icol[i];
1194:     for (j=0; j<nrows; j++) {
1195:       *bv++ = av[irow[j]];
1196:     }
1197:   }

1199:   /* Assemble the matrices so that the correct flags are set */
1200:   MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
1201:   MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);

1203:   /* Free work space */
1204:   ISRestoreIndices(isrow,&irow);
1205:   ISRestoreIndices(iscol,&icol);
1206:   *B = newmat;
1207:   return(0);
1208: }

1210: int MatGetSubMatrices_SeqDense(Mat A,int n,IS *irow,IS *icol,MatReuse scall,Mat **B)
1211: {
1212:   int ierr,i;

1215:   if (scall == MAT_INITIAL_MATRIX) {
1216:     PetscMalloc((n+1)*sizeof(Mat),B);
1217:   }

1219:   for (i=0; i<n; i++) {
1220:     MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);
1221:   }
1222:   return(0);
1223: }

1225: int MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
1226: {
1227:   Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
1228:   int          ierr;
1229:   PetscTruth   flag;

1232:   PetscTypeCompare((PetscObject)B,MATSEQDENSE,&flag);
1233:   if (!flag) {
1234:     MatCopy_Basic(A,B,str);
1235:     return(0);
1236:   }
1237:   if (A->m != B->m || A->n != B->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1238:   PetscMemcpy(b->v,a->v,A->m*A->n*sizeof(Scalar));
1239:   return(0);
1240: }

1242: int MatSetUpPreallocation_SeqDense(Mat A)
1243: {
1244:   int        ierr;

1247:    MatSeqDenseSetPreallocation(A,0);
1248:   return(0);
1249: }

1251: /* -------------------------------------------------------------------*/
1252: static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1253:        MatGetRow_SeqDense,
1254:        MatRestoreRow_SeqDense,
1255:        MatMult_SeqDense,
1256:        MatMultAdd_SeqDense,
1257:        MatMultTranspose_SeqDense,
1258:        MatMultTransposeAdd_SeqDense,
1259:        MatSolve_SeqDense,
1260:        MatSolveAdd_SeqDense,
1261:        MatSolveTranspose_SeqDense,
1262:        MatSolveTransposeAdd_SeqDense,
1263:        MatLUFactor_SeqDense,
1264:        MatCholeskyFactor_SeqDense,
1265:        MatRelax_SeqDense,
1266:        MatTranspose_SeqDense,
1267:        MatGetInfo_SeqDense,
1268:        MatEqual_SeqDense,
1269:        MatGetDiagonal_SeqDense,
1270:        MatDiagonalScale_SeqDense,
1271:        MatNorm_SeqDense,
1272:        0,
1273:        0,
1274:        0,
1275:        MatSetOption_SeqDense,
1276:        MatZeroEntries_SeqDense,
1277:        MatZeroRows_SeqDense,
1278:        MatLUFactorSymbolic_SeqDense,
1279:        MatLUFactorNumeric_SeqDense,
1280:        MatCholeskyFactorSymbolic_SeqDense,
1281:        MatCholeskyFactorNumeric_SeqDense,
1282:        MatSetUpPreallocation_SeqDense,
1283:        0,
1284:        MatGetOwnershipRange_SeqDense,
1285:        0,
1286:        0,
1287:        MatGetArray_SeqDense,
1288:        MatRestoreArray_SeqDense,
1289:        MatDuplicate_SeqDense,
1290:        0,
1291:        0,
1292:        0,
1293:        0,
1294:        MatAXPY_SeqDense,
1295:        MatGetSubMatrices_SeqDense,
1296:        0,
1297:        MatGetValues_SeqDense,
1298:        MatCopy_SeqDense,
1299:        0,
1300:        MatScale_SeqDense,
1301:        0,
1302:        0,
1303:        0,
1304:        MatGetBlockSize_SeqDense,
1305:        0,
1306:        0,
1307:        0,
1308:        0,
1309:        0,
1310:        0,
1311:        0,
1312:        0,
1313:        0,
1314:        0,
1315:        MatDestroy_SeqDense,
1316:        MatView_SeqDense,
1317:        MatGetMaps_Petsc};

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

1324:    Collective on MPI_Comm

1326:    Input Parameters:
1327: +  comm - MPI communicator, set to PETSC_COMM_SELF
1328: .  m - number of rows
1329: .  n - number of columns
1330: -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1331:    to control all matrix memory allocation.

1333:    Output Parameter:
1334: .  A - the matrix

1336:    Notes:
1337:    The data input variable is intended primarily for Fortran programmers
1338:    who wish to allocate their own matrix memory space.  Most users should
1339:    set data=PETSC_NULL.

1341:    Level: intermediate

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

1345: .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
1346: @*/
1347: int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A)
1348: {

1352:   MatCreate(comm,m,n,m,n,A);
1353:   MatSetType(*A,MATSEQDENSE);
1354:   MatSeqDenseSetPreallocation(*A,data);
1355:   return(0);
1356: }

1358: /*@C
1359:    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements

1361:    Collective on MPI_Comm

1363:    Input Parameters:
1364: +  A - the matrix
1365: -  data - the array (or PETSC_NULL)

1367:    Notes:
1368:    The data input variable is intended primarily for Fortran programmers
1369:    who wish to allocate their own matrix memory space.  Most users should
1370:    set data=PETSC_NULL.

1372:    Level: intermediate

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

1376: .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
1377: @*/
1378: int MatSeqDenseSetPreallocation(Mat B,Scalar *data)
1379: {
1380:   Mat_SeqDense *b;
1381:   int          ierr;
1382:   PetscTruth   flg2;

1385:   PetscTypeCompare((PetscObject)B,MATSEQDENSE,&flg2);
1386:   if (!flg2) return(0);
1387:   B->preallocated = PETSC_TRUE;
1388:   b               = (Mat_SeqDense*)B->data;
1389:   if (!data) {
1390:     ierr          = PetscMalloc((B->m*B->n+1)*sizeof(Scalar),&b->v);
1391:     ierr          = PetscMemzero(b->v,B->m*B->n*sizeof(Scalar));
1392:     b->user_alloc = PETSC_FALSE;
1393:     PetscLogObjectMemory(B,B->n*B->m*sizeof(Scalar));
1394:   } else { /* user-allocated storage */
1395:     b->v          = data;
1396:     b->user_alloc = PETSC_TRUE;
1397:   }
1398:   return(0);
1399: }

1401: EXTERN_C_BEGIN
1402: int MatCreate_SeqDense(Mat B)
1403: {
1404:   Mat_SeqDense *b;
1405:   int          ierr,size;

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

1411:   B->m = B->M = PetscMax(B->m,B->M);
1412:   B->n = B->N = PetscMax(B->n,B->N);

1414:   ierr            = PetscNew(Mat_SeqDense,&b);
1415:   ierr            = PetscMemzero(b,sizeof(Mat_SeqDense));
1416:   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
1417:   B->factor       = 0;
1418:   B->mapping      = 0;
1419:   PetscLogObjectMemory(B,sizeof(struct _p_Mat));
1420:   B->data         = (void*)b;

1422:   MapCreateMPI(B->comm,B->m,B->m,&B->rmap);
1423:   MapCreateMPI(B->comm,B->n,B->n,&B->cmap);

1425:   b->pivots       = 0;
1426:   b->roworiented  = PETSC_TRUE;
1427:   b->v            = 0;

1429:   return(0);
1430: }
1431: EXTERN_C_END