xref: /petsc/src/mat/impls/baij/mpi/mmbaij.c (revision 70f55243aafb320636e2a54ff30cab5d1e8d3d7b)
18016bdd1SSatish Balay #ifndef lint
2*70f55243SBarry Smith static char vcid[] = "$Id: mmbaij.c,v 1.5 1996/08/04 23:12:57 bsmith Exp bsmith $";
38016bdd1SSatish Balay #endif
48016bdd1SSatish Balay 
58016bdd1SSatish Balay 
68016bdd1SSatish Balay /*
7d9653453SSatish Balay    Support for the parallel BAIJ matrix vector multiply
88016bdd1SSatish Balay */
9*70f55243SBarry Smith #include "src/mat/impls/baij/mpi/mpibaij.h"
108016bdd1SSatish Balay #include "src/vec/vecimpl.h"
118016bdd1SSatish Balay 
12d9653453SSatish Balay int MatSetUpMultiply_MPIBAIJ(Mat mat)
138016bdd1SSatish Balay {
14d9653453SSatish Balay   Mat_MPIBAIJ *baij = (Mat_MPIBAIJ *) mat->data;
15d9653453SSatish Balay   Mat_SeqBAIJ *B = (Mat_SeqBAIJ *) (baij->B->data);
16d9653453SSatish Balay   int        Nbs = baij->Nbs,i,j,*indices,*aj = B->j,ierr,ec = 0,*garray;
1757b952d6SSatish Balay   int        col,bs = baij->bs,*tmp;
188016bdd1SSatish Balay   IS         from,to;
198016bdd1SSatish Balay   Vec        gvec;
208016bdd1SSatish Balay 
218016bdd1SSatish Balay   /* For the first stab we make an array as long as the number of columns */
22d9653453SSatish Balay   /* mark those columns that are in baij->B */
23d9653453SSatish Balay   indices = (int *) PetscMalloc( Nbs*sizeof(int) ); CHKPTRQ(indices);
24d9653453SSatish Balay   PetscMemzero(indices,Nbs*sizeof(int));
25d9653453SSatish Balay   for ( i=0; i<B->mbs; i++ ) {
268016bdd1SSatish Balay     for ( j=0; j<B->ilen[i]; j++ ) {
27d9653453SSatish Balay       if (!indices[aj[B->i[i] + j]]) ec++;
28d9653453SSatish Balay       indices[aj[B->i[i] + j] ] = 1;
298016bdd1SSatish Balay     }
308016bdd1SSatish Balay   }
318016bdd1SSatish Balay 
328016bdd1SSatish Balay   /* form array of columns we need */
338016bdd1SSatish Balay   garray = (int *) PetscMalloc( (ec+1)*sizeof(int) ); CHKPTRQ(garray);
34d9653453SSatish Balay   tmp    = (int *) PetscMalloc( (ec*bs+1)*sizeof(int) ); CHKPTRQ(tmp)
358016bdd1SSatish Balay   ec = 0;
36d9653453SSatish Balay   for ( i=0; i<Nbs; i++ ) {
378016bdd1SSatish Balay     if (indices[i]) garray[ec++] = i;
388016bdd1SSatish Balay   }
398016bdd1SSatish Balay 
408016bdd1SSatish Balay   /* make indices now point into garray */
418016bdd1SSatish Balay   for ( i=0; i<ec; i++ ) {
42d9653453SSatish Balay     indices[garray[i]] = i;
438016bdd1SSatish Balay   }
448016bdd1SSatish Balay 
458016bdd1SSatish Balay   /* compact out the extra columns in B */
46d9653453SSatish Balay   for ( i=0; i<B->mbs; i++ ) {
478016bdd1SSatish Balay     for ( j=0; j<B->ilen[i]; j++ ) {
48d9653453SSatish Balay       aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
498016bdd1SSatish Balay     }
508016bdd1SSatish Balay   }
51d9653453SSatish Balay   B->nbs = ec;
52d9653453SSatish Balay   B->n   = ec*B->bs;
538016bdd1SSatish Balay   PetscFree(indices);
548016bdd1SSatish Balay 
5557b952d6SSatish Balay   for ( i=0,col=0; i<ec; i++ ) {
5657b952d6SSatish Balay     for ( j=0; j<bs; j++,col++) tmp[col] = garray[i]*bs+j;
57d9653453SSatish Balay   }
588016bdd1SSatish Balay   /* create local vector that is used to scatter into */
59d9653453SSatish Balay   ierr = VecCreateSeq(MPI_COMM_SELF,ec*bs,&baij->lvec); CHKERRQ(ierr);
608016bdd1SSatish Balay 
61c16cb8f2SBarry Smith   /* create two temporary index sets for building scatter-gather */
62c16cb8f2SBarry Smith 
63c16cb8f2SBarry Smith   /* ierr = ISCreateSeq(MPI_COMM_SELF,ec*bs,tmp,&from); CHKERRQ(ierr); */
64c16cb8f2SBarry Smith   for ( i=0,col=0; i<ec; i++ ) {
65c16cb8f2SBarry Smith     garray[i] = bs*garray[i];
66c16cb8f2SBarry Smith   }
67c16cb8f2SBarry Smith   ierr = ISCreateBlockSeq(MPI_COMM_SELF,bs,ec,garray,&from);CHKERRQ(ierr);
68c16cb8f2SBarry Smith   for ( i=0,col=0; i<ec; i++ ) {
69c16cb8f2SBarry Smith     garray[i] = garray[i]/bs;
70c16cb8f2SBarry Smith   }
71c16cb8f2SBarry Smith 
72d9653453SSatish Balay   ierr = ISCreateStrideSeq(MPI_COMM_SELF,ec*bs,0,1,&to);CHKERRQ(ierr);
738016bdd1SSatish Balay 
748016bdd1SSatish Balay   /* create temporary global vector to generate scatter context */
758016bdd1SSatish Balay   /* this is inefficient, but otherwise we must do either
768016bdd1SSatish Balay      1) save garray until the first actual scatter when the vector is known or
778016bdd1SSatish Balay      2) have another way of generating a scatter context without a vector.*/
78d9653453SSatish Balay   ierr = VecCreateMPI(mat->comm,baij->n,baij->N,&gvec); CHKERRQ(ierr);
798016bdd1SSatish Balay 
808016bdd1SSatish Balay   /* gnerate the scatter context */
81d9653453SSatish Balay   ierr = VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);CHKERRQ(ierr);
82d9653453SSatish Balay   PLogObjectParent(mat,baij->Mvctx);
83d9653453SSatish Balay   PLogObjectParent(mat,baij->lvec);
848016bdd1SSatish Balay   PLogObjectParent(mat,from);
858016bdd1SSatish Balay   PLogObjectParent(mat,to);
86d9653453SSatish Balay   baij->garray = garray;
878016bdd1SSatish Balay   PLogObjectMemory(mat,(ec+1)*sizeof(int));
888016bdd1SSatish Balay   ierr = ISDestroy(from); CHKERRQ(ierr);
898016bdd1SSatish Balay   ierr = ISDestroy(to); CHKERRQ(ierr);
908016bdd1SSatish Balay   ierr = VecDestroy(gvec);
91d9653453SSatish Balay   PetscFree(tmp);
928016bdd1SSatish Balay   return 0;
938016bdd1SSatish Balay }
948016bdd1SSatish Balay 
958016bdd1SSatish Balay 
968016bdd1SSatish Balay /*
97d9653453SSatish Balay      Takes the local part of an already assembled MPIBAIJ matrix
988016bdd1SSatish Balay    and disassembles it. This is to allow new nonzeros into the matrix
998016bdd1SSatish Balay    that require more communication in the matrix vector multiply.
1008016bdd1SSatish Balay    Thus certain data-structures must be rebuilt.
1018016bdd1SSatish Balay 
1028016bdd1SSatish Balay    Kind of slow! But that's what application programmers get when
1038016bdd1SSatish Balay    they are sloppy.
1048016bdd1SSatish Balay */
105d9653453SSatish Balay int DisAssemble_MPIBAIJ(Mat A)
1068016bdd1SSatish Balay {
107d9653453SSatish Balay   Mat_MPIBAIJ *baij = (Mat_MPIBAIJ *) A->data;
108d9653453SSatish Balay   Mat        B = baij->B,Bnew;
109d9653453SSatish Balay   Mat_SeqBAIJ *Bbaij = (Mat_SeqBAIJ*)B->data;
110d9653453SSatish Balay   int        ierr,i,j,mbs=Bbaij->mbs,n = baij->N,col,*garray=baij->garray;
111d9653453SSatish Balay   int        k,bs=baij->bs,bs2=baij->bs2,*rvals,*nz,ec,m=Bbaij->m;
112d9653453SSatish Balay   Scalar     *a=Bbaij->a;
1138016bdd1SSatish Balay 
1148016bdd1SSatish Balay   /* free stuff related to matrix-vec multiply */
115d9653453SSatish Balay   ierr = VecGetSize(baij->lvec,&ec); /* needed for PLogObjectMemory below */
116d9653453SSatish Balay   ierr = VecDestroy(baij->lvec); CHKERRQ(ierr); baij->lvec = 0;
117d9653453SSatish Balay   ierr = VecScatterDestroy(baij->Mvctx); CHKERRQ(ierr); baij->Mvctx = 0;
118d9653453SSatish Balay   if (baij->colmap) {
119d9653453SSatish Balay     PetscFree(baij->colmap); baij->colmap = 0;
120d9653453SSatish Balay     PLogObjectMemory(A,-Bbaij->nbs*sizeof(int));
1218016bdd1SSatish Balay   }
1228016bdd1SSatish Balay 
1238016bdd1SSatish Balay   /* make sure that B is assembled so we can access its values */
1246d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1256d4a8577SBarry Smith   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1268016bdd1SSatish Balay 
1278016bdd1SSatish Balay   /* invent new B and copy stuff over */
128d9653453SSatish Balay   nz = (int *) PetscMalloc( mbs*sizeof(int) ); CHKPTRQ(nz);
129d9653453SSatish Balay   for ( i=0; i<mbs; i++ ) {
130d9653453SSatish Balay     nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
1318016bdd1SSatish Balay   }
132d9653453SSatish Balay   ierr = MatCreateSeqBAIJ(MPI_COMM_SELF,baij->bs,m,n,0,nz,&Bnew); CHKERRQ(ierr);
1338016bdd1SSatish Balay   PetscFree(nz);
134d9653453SSatish Balay 
135d9653453SSatish Balay   rvals = (int *) PetscMalloc(bs*sizeof(int)); CHKPTRQ(rvals);
136d9653453SSatish Balay   for ( i=0; i<mbs; i++ ) {
137d9653453SSatish Balay     rvals[0] = bs*i;
138d9653453SSatish Balay     for ( j=1; j<bs; j++ ) { rvals[j] = rvals[j-1] + 1; }
139d9653453SSatish Balay     for ( j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++ ) {
14057b952d6SSatish Balay       col = garray[Bbaij->j[j]]*bs;
141d9653453SSatish Balay       for (k=0; k<bs; k++ ) {
142d9653453SSatish Balay         ierr = MatSetValues(Bnew,bs,rvals,1,&col,a+j*bs2,INSERT_VALUES);CHKERRQ(ierr);
143d9653453SSatish Balay         col++;
1448016bdd1SSatish Balay       }
1458016bdd1SSatish Balay     }
146d9653453SSatish Balay   }
147d9653453SSatish Balay   PetscFree(baij->garray); baij->garray = 0;
14857b952d6SSatish Balay   PetscFree(rvals);
1498016bdd1SSatish Balay   PLogObjectMemory(A,-ec*sizeof(int));
1508016bdd1SSatish Balay   ierr = MatDestroy(B); CHKERRQ(ierr);
1518016bdd1SSatish Balay   PLogObjectParent(A,Bnew);
152d9653453SSatish Balay   baij->B = Bnew;
1538016bdd1SSatish Balay   A->was_assembled = PETSC_FALSE;
1548016bdd1SSatish Balay   return 0;
1558016bdd1SSatish Balay }
1568016bdd1SSatish Balay 
1578016bdd1SSatish Balay 
158