xref: /petsc/src/mat/impls/dense/mpi/mmdense.c (revision 3a40ed3dce77c081171d005ae1a6ff4bb9d13b6f)
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: mmdense.c,v 1.12 1997/07/09 20:53:40 balay 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 #undef __FUNC__
12 #define __FUNC__ "MatSetUpMultiply_MPIDense"
13 int MatSetUpMultiply_MPIDense(Mat mat)
14 {
15   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
16   int          ierr,n;
17   IS           tofrom;
18   Vec          gvec;
19 
20   PetscFunctionBegin;
21   /* Create local vector that is used to scatter into */
22   ierr = VecCreateSeq(PETSC_COMM_SELF,mdn->N,&mdn->lvec); CHKERRQ(ierr);
23 
24   /* Create temporary index set for building scatter gather */
25   ierr = ISCreateStride(PETSC_COMM_SELF,mdn->N,0,1,&tofrom); CHKERRQ(ierr);
26 
27   /* Create temporary global vector to generate scatter context */
28   n    = mdn->cowners[mdn->rank+1] - mdn->cowners[mdn->rank];
29   ierr = VecCreateMPI(mat->comm,n,mdn->N,&gvec); CHKERRQ(ierr);
30 
31   /* Generate the scatter context */
32   ierr = VecScatterCreate(gvec,tofrom,mdn->lvec,tofrom,&mdn->Mvctx); CHKERRQ(ierr);
33   PLogObjectParent(mat,mdn->Mvctx);
34   PLogObjectParent(mat,mdn->lvec);
35   PLogObjectParent(mat,tofrom);
36   PLogObjectParent(mat,gvec);
37 
38   ierr = ISDestroy(tofrom); CHKERRQ(ierr);
39   ierr = VecDestroy(gvec); CHKERRQ(ierr);
40   PetscFunctionReturn(0);
41 }
42 
43 
44 
45