xref: /petsc/src/mat/impls/dense/mpi/mmdense.c (revision 6a67db7e1512b5d9ff9ba8f94e48d05f0a9eae3d)
1 #ifndef lint
2 static char vcid[] = "$Id: mmdense.c,v 1.7 1996/08/08 14:42:42 bsmith Exp bsmith $";
3 #endif
4 
5 /*
6    Support for the parallel dense matrix vector multiply
7 */
8 #include "src/mat/impls/dense/mpi/mpidense.h"
9 #include "src/vec/vecimpl.h"
10 
11 int MatSetUpMultiply_MPIDense(Mat mat)
12 {
13   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
14   int          ierr,n;
15   IS           tofrom;
16   Vec          gvec;
17 
18   /* Create local vector that is used to scatter into */
19   ierr = VecCreateSeq(MPI_COMM_SELF,mdn->N,&mdn->lvec); CHKERRQ(ierr);
20 
21   /* Create temporary index set for building scatter gather */
22   ierr = ISCreateStride(MPI_COMM_SELF,mdn->N,0,1,&tofrom); CHKERRQ(ierr);
23 
24   /* Create temporary global vector to generate scatter context */
25   n    = mdn->cowners[mdn->rank+1] - mdn->cowners[mdn->rank];
26   ierr = VecCreateMPI(mat->comm,n,mdn->N,&gvec); CHKERRQ(ierr);
27 
28   /* Generate the scatter context */
29   ierr = VecScatterCreate(gvec,tofrom,mdn->lvec,tofrom,&mdn->Mvctx); CHKERRQ(ierr);
30   PLogObjectParent(mat,mdn->Mvctx);
31   PLogObjectParent(mat,mdn->lvec);
32   PLogObjectParent(mat,tofrom);
33   PLogObjectParent(mat,gvec);
34 
35   ierr = ISDestroy(tofrom); CHKERRQ(ierr);
36   ierr = VecDestroy(gvec); CHKERRQ(ierr);
37   return 0;
38 }
39 
40 
41 
42