1cb512458SBarry Smith #ifndef lint 2*639f9d9dSBarry Smith static char vcid[] = "$Id: mmaij.c,v 1.31 1996/09/19 20:21:32 bsmith Exp bsmith $"; 3cb512458SBarry Smith #endif 48c79f6d3SBarry Smith 5d6dfbf8fSBarry Smith 68c79f6d3SBarry Smith /* 78c79f6d3SBarry Smith Support for the parallel AIJ matrix vector multiply 88c79f6d3SBarry Smith */ 970f55243SBarry Smith #include "src/mat/impls/aij/mpi/mpiaij.h" 10f5eb4b81SSatish Balay #include "src/vec/vecimpl.h" 118c79f6d3SBarry Smith 1244a69424SLois Curfman McInnes int MatSetUpMultiply_MPIAIJ(Mat mat) 138c79f6d3SBarry Smith { 1444a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 15ec8511deSBarry Smith Mat_SeqAIJ *B = (Mat_SeqAIJ *) (aij->B->data); 16416022c9SBarry Smith int N = aij->N,i,j,*indices,*aj = B->j,ierr,ec = 0,*garray; 17416022c9SBarry Smith int shift = B->indexshift; 181eb62cbbSBarry Smith IS from,to; 191eb62cbbSBarry Smith Vec gvec; 208c79f6d3SBarry Smith 218c79f6d3SBarry Smith /* For the first stab we make an array as long as the number of columns */ 221eb62cbbSBarry Smith /* mark those columns that are in aij->B */ 230452661fSBarry Smith indices = (int *) PetscMalloc( N*sizeof(int) ); CHKPTRQ(indices); 24cddf8d76SBarry Smith PetscMemzero(indices,N*sizeof(int)); 25d6dfbf8fSBarry Smith for ( i=0; i<B->m; i++ ) { 26d6dfbf8fSBarry Smith for ( j=0; j<B->ilen[i]; j++ ) { 27dbb450caSBarry Smith if (!indices[aj[B->i[i] +shift + j] + shift]) ec++; 28416022c9SBarry Smith indices[aj[B->i[i] + shift + j] + shift] = 1; 29416022c9SBarry Smith } 301eb62cbbSBarry Smith } 318c79f6d3SBarry Smith 321eb62cbbSBarry Smith /* form array of columns we need */ 330452661fSBarry Smith garray = (int *) PetscMalloc( (ec+1)*sizeof(int) ); CHKPTRQ(garray); 341eb62cbbSBarry Smith ec = 0; 351eb62cbbSBarry Smith for ( i=0; i<N; i++ ) { 361eb62cbbSBarry Smith if (indices[i]) garray[ec++] = i; 371eb62cbbSBarry Smith } 381eb62cbbSBarry Smith 391eb62cbbSBarry Smith /* make indices now point into garray */ 401eb62cbbSBarry Smith for ( i=0; i<ec; i++ ) { 41dbb450caSBarry Smith indices[garray[i]] = i-shift; 421eb62cbbSBarry Smith } 431eb62cbbSBarry Smith 441eb62cbbSBarry Smith /* compact out the extra columns in B */ 45d6dfbf8fSBarry Smith for ( i=0; i<B->m; i++ ) { 46d6dfbf8fSBarry Smith for ( j=0; j<B->ilen[i]; j++ ) { 47dbb450caSBarry Smith aj[B->i[i] + shift + j] = indices[aj[B->i[i] + shift + j]+shift]; 481eb62cbbSBarry Smith } 49d6dfbf8fSBarry Smith } 50*639f9d9dSBarry Smith B->n = aij->B->n = aij->B->N = ec; 510452661fSBarry Smith PetscFree(indices); 521eb62cbbSBarry Smith 531eb62cbbSBarry Smith /* create local vector that is used to scatter into */ 54fafbff53SBarry Smith ierr = VecCreateSeq(MPI_COMM_SELF,ec,&aij->lvec); CHKERRQ(ierr); 551eb62cbbSBarry Smith 56d6dfbf8fSBarry Smith /* create two temporary Index sets for build scatter gather */ 57537820f0SBarry Smith ierr = ISCreateGeneral(MPI_COMM_SELF,ec,garray,&from); CHKERRQ(ierr); 58537820f0SBarry Smith ierr = ISCreateStride(MPI_COMM_SELF,ec,0,1,&to); CHKERRQ(ierr); 591eb62cbbSBarry Smith 601eb62cbbSBarry Smith /* create temporary global vector to generate scatter context */ 611eb62cbbSBarry Smith /* this is inefficient, but otherwise we must do either 621eb62cbbSBarry Smith 1) save garray until the first actual scatter when the vector is known or 631eb62cbbSBarry Smith 2) have another way of generating a scatter context without a vector.*/ 6478b31e54SBarry Smith ierr = VecCreateMPI(mat->comm,aij->n,aij->N,&gvec); CHKERRQ(ierr); 651eb62cbbSBarry Smith 662d336d48SLois Curfman McInnes /* generate the scatter context */ 6708480c60SBarry Smith ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx); CHKERRQ(ierr); 68a5a9c739SBarry Smith PLogObjectParent(mat,aij->Mvctx); 69a5a9c739SBarry Smith PLogObjectParent(mat,aij->lvec); 70464493b3SBarry Smith PLogObjectParent(mat,from); 71464493b3SBarry Smith PLogObjectParent(mat,to); 729e25ed09SBarry Smith aij->garray = garray; 73464493b3SBarry Smith PLogObjectMemory(mat,(ec+1)*sizeof(int)); 7478b31e54SBarry Smith ierr = ISDestroy(from); CHKERRQ(ierr); 7578b31e54SBarry Smith ierr = ISDestroy(to); CHKERRQ(ierr); 761eb62cbbSBarry Smith ierr = VecDestroy(gvec); 778c79f6d3SBarry Smith return 0; 788c79f6d3SBarry Smith } 799e25ed09SBarry Smith 809e25ed09SBarry Smith 812493cbb0SBarry Smith /* 822493cbb0SBarry Smith Takes the local part of an already assembled MPIAIJ matrix 832493cbb0SBarry Smith and disassembles it. This is to allow new nonzeros into the matrix 842493cbb0SBarry Smith that require more communication in the matrix vector multiply. 852493cbb0SBarry Smith Thus certain data-structures must be rebuilt. 862493cbb0SBarry Smith 872493cbb0SBarry Smith Kind of slow! But that's what application programmers get when 882493cbb0SBarry Smith they are sloppy. 892493cbb0SBarry Smith */ 902493cbb0SBarry Smith int DisAssemble_MPIAIJ(Mat A) 912493cbb0SBarry Smith { 922493cbb0SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ *) A->data; 932493cbb0SBarry Smith Mat B = aij->B,Bnew; 94ec8511deSBarry Smith Mat_SeqAIJ *Baij = (Mat_SeqAIJ*)B->data; 952493cbb0SBarry Smith int ierr,i,j,m = Baij->m,n = aij->N,col,ct = 0,*garray = aij->garray; 96416022c9SBarry Smith int *nz,ec,shift = Baij->indexshift; 972493cbb0SBarry Smith Scalar v; 982493cbb0SBarry Smith 992493cbb0SBarry Smith /* free stuff related to matrix-vec multiply */ 100464493b3SBarry Smith ierr = VecGetSize(aij->lvec,&ec); /* needed for PLogObjectMemory below */ 1012493cbb0SBarry Smith ierr = VecDestroy(aij->lvec); CHKERRQ(ierr); aij->lvec = 0; 10208480c60SBarry Smith ierr = VecScatterDestroy(aij->Mvctx); CHKERRQ(ierr); aij->Mvctx = 0; 103464493b3SBarry Smith if (aij->colmap) { 1040452661fSBarry Smith PetscFree(aij->colmap); aij->colmap = 0; 105464493b3SBarry Smith PLogObjectMemory(A,-Baij->n*sizeof(int)); 106464493b3SBarry Smith } 1072493cbb0SBarry Smith 1082493cbb0SBarry Smith /* make sure that B is assembled so we can access its values */ 1096d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 110fe2f2677SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 1112493cbb0SBarry Smith 1122493cbb0SBarry Smith /* invent new B and copy stuff over */ 113*639f9d9dSBarry Smith nz = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(nz); 11448b35521SBarry Smith for ( i=0; i<m; i++ ) { 11548b35521SBarry Smith nz[i] = Baij->i[i+1] - Baij->i[i]; 11648b35521SBarry Smith } 117fafbff53SBarry Smith ierr = MatCreateSeqAIJ(MPI_COMM_SELF,m,n,0,nz,&Bnew); CHKERRQ(ierr); 1180452661fSBarry Smith PetscFree(nz); 1192493cbb0SBarry Smith for ( i=0; i<m; i++ ) { 120dbb450caSBarry Smith for ( j=Baij->i[i]+shift; j<Baij->i[i+1]+shift; j++ ) { 121dbb450caSBarry Smith col = garray[Baij->j[ct]+shift]; 1222493cbb0SBarry Smith v = Baij->a[ct++]; 123dbb450caSBarry Smith ierr = MatSetValues(Bnew,1,&i,1,&col,&v,INSERT_VALUES); CHKERRQ(ierr); 1242493cbb0SBarry Smith } 1252493cbb0SBarry Smith } 1260452661fSBarry Smith PetscFree(aij->garray); aij->garray = 0; 127464493b3SBarry Smith PLogObjectMemory(A,-ec*sizeof(int)); 1282493cbb0SBarry Smith ierr = MatDestroy(B); CHKERRQ(ierr); 129464493b3SBarry Smith PLogObjectParent(A,Bnew); 1302493cbb0SBarry Smith aij->B = Bnew; 131227d817aSBarry Smith A->was_assembled = PETSC_FALSE; 1322493cbb0SBarry Smith return 0; 1332493cbb0SBarry Smith } 1342493cbb0SBarry Smith 13548b35521SBarry Smith 136