Actual source code: mmbdiag.c

  1: #define PETSCMAT_DLL

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

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

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

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

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

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

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

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

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

 91:   ISDestroy(to);
 92:   ISDestroy(from);
 93:   VecDestroy(gvec);
 94:   return(0);
 95: }