xref: /petsc/src/mat/impls/dense/mpi/mmdense.c (revision f5eb4b8140e51d60f0c4bc35d8fe3d9a5fa43604)
10500aeceSBarry Smith #ifndef lint
2*f5eb4b81SSatish Balay static char vcid[] = "$Id: mmdense.c,v 1.5 1995/11/03 02:59:00 bsmith Exp balay $";
30500aeceSBarry Smith #endif
40500aeceSBarry Smith 
50500aeceSBarry Smith /*
60500aeceSBarry Smith    Support for the parallel dense matrix vector multiply
70500aeceSBarry Smith */
80500aeceSBarry Smith #include "mpidense.h"
9*f5eb4b81SSatish Balay #include "src/vec/vecimpl.h"
100500aeceSBarry Smith 
110500aeceSBarry Smith int MatSetUpMultiply_MPIDense(Mat mat)
120500aeceSBarry Smith {
130500aeceSBarry Smith   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
14ed3cc1f0SBarry Smith   int          ierr,n;
15d9f96c7cSLois Curfman McInnes   IS           tofrom;
160500aeceSBarry Smith   Vec          gvec;
170500aeceSBarry Smith 
180500aeceSBarry Smith   /* Create local vector that is used to scatter into */
190500aeceSBarry Smith   ierr = VecCreateSeq(MPI_COMM_SELF,mdn->N,&mdn->lvec); CHKERRQ(ierr);
200500aeceSBarry Smith 
210500aeceSBarry Smith   /* Create temporary index set for building scatter gather */
22d9f96c7cSLois Curfman McInnes   ierr = ISCreateStrideSeq(MPI_COMM_SELF,mdn->N,0,1,&tofrom); CHKERRQ(ierr);
230500aeceSBarry Smith 
240500aeceSBarry Smith   /* Create temporary global vector to generate scatter context */
25ed3cc1f0SBarry Smith   n    = mdn->cowners[mdn->rank+1] - mdn->cowners[mdn->rank];
26ed3cc1f0SBarry Smith   ierr = VecCreateMPI(mat->comm,n,mdn->N,&gvec); CHKERRQ(ierr);
27e19c7942SLois Curfman McInnes 
280500aeceSBarry Smith   /* Generate the scatter context */
29d9f96c7cSLois Curfman McInnes   ierr = VecScatterCreate(gvec,tofrom,mdn->lvec,tofrom,&mdn->Mvctx); CHKERRQ(ierr);
300500aeceSBarry Smith   PLogObjectParent(mat,mdn->Mvctx);
310500aeceSBarry Smith   PLogObjectParent(mat,mdn->lvec);
32d9f96c7cSLois Curfman McInnes   PLogObjectParent(mat,tofrom);
33d9f96c7cSLois Curfman McInnes   PLogObjectParent(mat,gvec);
34e19c7942SLois Curfman McInnes 
35d9f96c7cSLois Curfman McInnes   ierr = ISDestroy(tofrom); CHKERRQ(ierr);
36d9f96c7cSLois Curfman McInnes   ierr = VecDestroy(gvec); CHKERRQ(ierr);
370500aeceSBarry Smith   return 0;
380500aeceSBarry Smith }
390500aeceSBarry Smith 
400500aeceSBarry Smith 
410500aeceSBarry Smith 
42