xref: /petsc/src/mat/impls/dense/mpi/mmdense.c (revision 08480c60afa5ef1d2e4e27b9ebdf48b02c6a2186)
1 #ifndef lint
2 static char vcid[] = "$Id: mmdense.c,v 1.3 1995/10/23 22:02:53 curfman Exp curfman $";
3 #endif
4 
5 /*
6    Support for the parallel dense matrix vector multiply
7 */
8 #include "mpidense.h"
9 #include "vec/vecimpl.h"
10 
11 int MatSetUpMultiply_MPIDense(Mat mat)
12 {
13   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
14   int          ierr;
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 = ISCreateStrideSeq(MPI_COMM_SELF,mdn->N,0,1,&tofrom); CHKERRQ(ierr);
23 
24   /* Create temporary global vector to generate scatter context */
25   ierr = VecCreateMPI(mat->comm,PETSC_DECIDE,mdn->N,&gvec); CHKERRQ(ierr);
26 
27   /* Generate the scatter context */
28   ierr = VecScatterCreate(gvec,tofrom,mdn->lvec,tofrom,&mdn->Mvctx); CHKERRQ(ierr);
29   PLogObjectParent(mat,mdn->Mvctx);
30   PLogObjectParent(mat,mdn->lvec);
31   PLogObjectParent(mat,tofrom);
32   PLogObjectParent(mat,gvec);
33 
34   ierr = ISDestroy(tofrom); CHKERRQ(ierr);
35   ierr = VecDestroy(gvec); CHKERRQ(ierr);
36   return 0;
37 }
38 
39 
40 
41