Actual source code: mmbdiag.c

  1: /*$Id: mmbdiag.c,v 1.40 2001/03/23 23:22:05 balay Exp $*/

  3: /*
  4:    Support for the MPIBDIAG matrix-vector multiply
  5: */
 6:  #include src/mat/impls/bdiag/mpi/mpibdiag.h
 7:  #include src/vec/vecimpl.h

  9: int MatSetUpMultiply_MPIBDiag(Mat mat)
 10: {
 11:   Mat_MPIBDiag *mbd = (Mat_MPIBDiag*)mat->data;
 12:   Mat_SeqBDiag *lmbd = (Mat_SeqBDiag*)mbd->A->data;
 13:   int          ierr,N = mat->N,*indices,*garray,ec=0;
 14:   int          bs = lmbd->bs,d,i,j,diag;
 15:   IS           to,from;
 16:   Vec          gvec;

 19:   /* We make an array as long as the number of columns */
 20:   /* mark those columns that are in mbd->A */
 21:   PetscMalloc((N+1)*sizeof(int),&indices);
 22:   ierr    = PetscMemzero(indices,N*sizeof(int));

 24:   if (bs == 1) {
 25:     for (d=0; d<lmbd->nd; d++) {
 26:       diag = lmbd->diag[d];
 27:       if (diag > 0) { /* col = loc */
 28:         for (j=0; j<lmbd->bdlen[d]; j++) {
 29:           if (!indices[j]) ec++;
 30:           indices[j] = 1;
 31:         }
 32:       } else { /* col = loc-diag */
 33:         for (j=0; j<lmbd->bdlen[d]; j++) {
 34:           if (!indices[j-diag]) ec++;
 35:           indices[j-diag] = 1;
 36:         }
 37:       }
 38:     }
 39:   } else {
 40:     for (d=0; d<lmbd->nd; d++) {
 41:       diag = lmbd->diag[d];
 42:       if (diag > 0) { /* col = loc */
 43:         for (j=0; j<lmbd->bdlen[d]; j++) {
 44:           if (!indices[bs*j]) ec += bs;
 45:           for (i=0; i<bs; i++) indices[bs*j+i] = 1;
 46:         }
 47:       } else { /* col = loc-diag */
 48:         for (j=0; j<lmbd->bdlen[d]; j++) {
 49:           if (!indices[bs*(j-diag)]) ec += bs;
 50:           for (i=0; i<bs; i++) indices[bs*(j-diag)+i] = 1;
 51:         }
 52:       }
 53:     }
 54:   }

 56:   /* form array of columns we need */
 57:   PetscMalloc((ec+1)*sizeof(int),&garray);
 58:   ec   = 0;
 59:   for (i=0; i<N; i++) {
 60:     if (indices[i]) garray[ec++] = i;
 61:   }
 62:   PetscFree(indices);

 64:   /* create local vector that is used to scatter into */
 65:   VecCreateSeq(PETSC_COMM_SELF,N,&mbd->lvec);

 67:   /* create temporary index set for building scatter-gather */
 68:   ISCreateGeneral(mat->comm,ec,garray,&from);
 69:   ISCreateGeneral(PETSC_COMM_SELF,ec,garray,&to);
 70:   PetscFree(garray);

 72:   /* create temporary global vector to generate scatter context */
 73:   /* this is inefficient, but otherwise we must do either 
 74:      1) save garray until the first actual scatter when the vector is known or
 75:      2) have another way of generating a scatter context without a vector.*/
 76:   /*
 77:      This is not correct for a rectangular matrix mbd->m? 
 78:   */
 79:   VecCreateMPI(mat->comm,mat->m,mat->N,&gvec);

 81:   /* generate the scatter context */
 82:   VecScatterCreate(gvec,from,mbd->lvec,to,&mbd->Mvctx);
 83:   PetscLogObjectParent(mat,mbd->Mvctx);
 84:   PetscLogObjectParent(mat,mbd->lvec);
 85:   PetscLogObjectParent(mat,to);
 86:   PetscLogObjectParent(mat,from);
 87:   PetscLogObjectParent(mat,gvec);

 89:   ISDestroy(to);
 90:   ISDestroy(from);
 91:   VecDestroy(gvec);
 92:   return(0);
 93: }