xref: /petsc/src/mat/impls/dense/mpi/mmdense.c (revision 3a40ed3dce77c081171d005ae1a6ff4bb9d13b6f)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*3a40ed3dSBarry Smith static char vcid[] = "$Id: mmdense.c,v 1.12 1997/07/09 20:53:40 balay Exp bsmith $";
30500aeceSBarry Smith #endif
40500aeceSBarry Smith 
50500aeceSBarry Smith /*
60500aeceSBarry Smith    Support for the parallel dense matrix vector multiply
70500aeceSBarry Smith */
870f55243SBarry Smith #include "src/mat/impls/dense/mpi/mpidense.h"
9f5eb4b81SSatish Balay #include "src/vec/vecimpl.h"
100500aeceSBarry Smith 
115615d1e5SSatish Balay #undef __FUNC__
125615d1e5SSatish Balay #define __FUNC__ "MatSetUpMultiply_MPIDense"
130500aeceSBarry Smith int MatSetUpMultiply_MPIDense(Mat mat)
140500aeceSBarry Smith {
150500aeceSBarry Smith   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
16ed3cc1f0SBarry Smith   int          ierr,n;
17d9f96c7cSLois Curfman McInnes   IS           tofrom;
180500aeceSBarry Smith   Vec          gvec;
190500aeceSBarry Smith 
20*3a40ed3dSBarry Smith   PetscFunctionBegin;
210500aeceSBarry Smith   /* Create local vector that is used to scatter into */
22029af93fSBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,mdn->N,&mdn->lvec); CHKERRQ(ierr);
230500aeceSBarry Smith 
240500aeceSBarry Smith   /* Create temporary index set for building scatter gather */
25029af93fSBarry Smith   ierr = ISCreateStride(PETSC_COMM_SELF,mdn->N,0,1,&tofrom); CHKERRQ(ierr);
260500aeceSBarry Smith 
270500aeceSBarry Smith   /* Create temporary global vector to generate scatter context */
28ed3cc1f0SBarry Smith   n    = mdn->cowners[mdn->rank+1] - mdn->cowners[mdn->rank];
29ed3cc1f0SBarry Smith   ierr = VecCreateMPI(mat->comm,n,mdn->N,&gvec); CHKERRQ(ierr);
30e19c7942SLois Curfman McInnes 
310500aeceSBarry Smith   /* Generate the scatter context */
32d9f96c7cSLois Curfman McInnes   ierr = VecScatterCreate(gvec,tofrom,mdn->lvec,tofrom,&mdn->Mvctx); CHKERRQ(ierr);
330500aeceSBarry Smith   PLogObjectParent(mat,mdn->Mvctx);
340500aeceSBarry Smith   PLogObjectParent(mat,mdn->lvec);
35d9f96c7cSLois Curfman McInnes   PLogObjectParent(mat,tofrom);
36d9f96c7cSLois Curfman McInnes   PLogObjectParent(mat,gvec);
37e19c7942SLois Curfman McInnes 
38d9f96c7cSLois Curfman McInnes   ierr = ISDestroy(tofrom); CHKERRQ(ierr);
39d9f96c7cSLois Curfman McInnes   ierr = VecDestroy(gvec); CHKERRQ(ierr);
40*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
410500aeceSBarry Smith }
420500aeceSBarry Smith 
430500aeceSBarry Smith 
440500aeceSBarry Smith 
45