1cb512458SBarry Smith #ifndef lint 2*08480c60SBarry Smith static char vcid[] = "$Id: mmaij.c,v 1.20 1995/09/30 19:28:49 bsmith Exp bsmith $"; 3cb512458SBarry Smith #endif 48c79f6d3SBarry Smith 5d6dfbf8fSBarry Smith 68c79f6d3SBarry Smith /* 78c79f6d3SBarry Smith Support for the parallel AIJ matrix vector multiply 88c79f6d3SBarry Smith */ 98c79f6d3SBarry Smith #include "mpiaij.h" 108c79f6d3SBarry Smith #include "vec/vecimpl.h" 118c79f6d3SBarry Smith #include "../seq/aij.h" 128c79f6d3SBarry Smith 1344a69424SLois Curfman McInnes int MatSetUpMultiply_MPIAIJ(Mat mat) 148c79f6d3SBarry Smith { 1544a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 16ec8511deSBarry Smith Mat_SeqAIJ *B = (Mat_SeqAIJ *) (aij->B->data); 17416022c9SBarry Smith int N = aij->N,i,j,*indices,*aj = B->j,ierr,ec = 0,*garray; 18416022c9SBarry Smith int shift = B->indexshift; 191eb62cbbSBarry Smith IS from,to; 201eb62cbbSBarry Smith Vec gvec; 218c79f6d3SBarry Smith 228c79f6d3SBarry Smith /* For the first stab we make an array as long as the number of columns */ 231eb62cbbSBarry Smith /* mark those columns that are in aij->B */ 2478b31e54SBarry Smith indices = (int *) PETSCMALLOC( N*sizeof(int) ); CHKPTRQ(indices); 25416022c9SBarry Smith PetscZero(indices,N*sizeof(int)); 26d6dfbf8fSBarry Smith for ( i=0; i<B->m; i++ ) { 27d6dfbf8fSBarry Smith for ( j=0; j<B->ilen[i]; j++ ) { 28dbb450caSBarry Smith if (!indices[aj[B->i[i] +shift + j] + shift]) ec++; 29416022c9SBarry Smith indices[aj[B->i[i] + shift + j] + shift] = 1; 30416022c9SBarry Smith } 311eb62cbbSBarry Smith } 328c79f6d3SBarry Smith 331eb62cbbSBarry Smith /* form array of columns we need */ 3478b31e54SBarry Smith garray = (int *) PETSCMALLOC( (ec+1)*sizeof(int) ); CHKPTRQ(garray); 351eb62cbbSBarry Smith ec = 0; 361eb62cbbSBarry Smith for ( i=0; i<N; i++ ) { 371eb62cbbSBarry Smith if (indices[i]) garray[ec++] = i; 381eb62cbbSBarry Smith } 391eb62cbbSBarry Smith 401eb62cbbSBarry Smith /* make indices now point into garray */ 411eb62cbbSBarry Smith for ( i=0; i<ec; i++ ) { 42dbb450caSBarry Smith indices[garray[i]] = i-shift; 431eb62cbbSBarry Smith } 441eb62cbbSBarry Smith 451eb62cbbSBarry Smith /* compact out the extra columns in B */ 46d6dfbf8fSBarry Smith for ( i=0; i<B->m; i++ ) { 47d6dfbf8fSBarry Smith for ( j=0; j<B->ilen[i]; j++ ) { 48dbb450caSBarry Smith aj[B->i[i] + shift + j] = indices[aj[B->i[i] + shift + j]+shift]; 491eb62cbbSBarry Smith } 50d6dfbf8fSBarry Smith } 51d6dfbf8fSBarry Smith B->n = ec; 5278b31e54SBarry Smith PETSCFREE(indices); 531eb62cbbSBarry Smith 541eb62cbbSBarry Smith /* create local vector that is used to scatter into */ 55fafbff53SBarry Smith ierr = VecCreateSeq(MPI_COMM_SELF,ec,&aij->lvec); CHKERRQ(ierr); 561eb62cbbSBarry Smith 57d6dfbf8fSBarry Smith /* create two temporary Index sets for build scatter gather */ 58fafbff53SBarry Smith ierr = ISCreateSeq(MPI_COMM_SELF,ec,garray,&from); CHKERRQ(ierr); 59fafbff53SBarry Smith ierr = ISCreateStrideSeq(MPI_COMM_SELF,ec,0,1,&to); CHKERRQ(ierr); 601eb62cbbSBarry Smith 611eb62cbbSBarry Smith /* create temporary global vector to generate scatter context */ 621eb62cbbSBarry Smith /* this is inefficient, but otherwise we must do either 631eb62cbbSBarry Smith 1) save garray until the first actual scatter when the vector is known or 641eb62cbbSBarry Smith 2) have another way of generating a scatter context without a vector.*/ 6578b31e54SBarry Smith ierr = VecCreateMPI(mat->comm,aij->n,aij->N,&gvec); CHKERRQ(ierr); 661eb62cbbSBarry Smith 671eb62cbbSBarry Smith /* gnerate the scatter context */ 68*08480c60SBarry Smith ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx); CHKERRQ(ierr); 69a5a9c739SBarry Smith PLogObjectParent(mat,aij->Mvctx); 70a5a9c739SBarry Smith PLogObjectParent(mat,aij->lvec); 71464493b3SBarry Smith PLogObjectParent(mat,from); 72464493b3SBarry Smith PLogObjectParent(mat,to); 739e25ed09SBarry Smith aij->garray = garray; 74464493b3SBarry Smith PLogObjectMemory(mat,(ec+1)*sizeof(int)); 7578b31e54SBarry Smith ierr = ISDestroy(from); CHKERRQ(ierr); 7678b31e54SBarry Smith ierr = ISDestroy(to); CHKERRQ(ierr); 771eb62cbbSBarry Smith ierr = VecDestroy(gvec); 788c79f6d3SBarry Smith return 0; 798c79f6d3SBarry Smith } 809e25ed09SBarry Smith 819e25ed09SBarry Smith 822493cbb0SBarry Smith /* 832493cbb0SBarry Smith Takes the local part of an already assembled MPIAIJ matrix 842493cbb0SBarry Smith and disassembles it. This is to allow new nonzeros into the matrix 852493cbb0SBarry Smith that require more communication in the matrix vector multiply. 862493cbb0SBarry Smith Thus certain data-structures must be rebuilt. 872493cbb0SBarry Smith 882493cbb0SBarry Smith Kind of slow! But that's what application programmers get when 892493cbb0SBarry Smith they are sloppy. 902493cbb0SBarry Smith */ 912493cbb0SBarry Smith int DisAssemble_MPIAIJ(Mat A) 922493cbb0SBarry Smith { 932493cbb0SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ *) A->data; 942493cbb0SBarry Smith Mat B = aij->B,Bnew; 95ec8511deSBarry Smith Mat_SeqAIJ *Baij = (Mat_SeqAIJ*)B->data; 962493cbb0SBarry Smith int ierr,i,j,m=Baij->m,n = aij->N,col,ct = 0,*garray = aij->garray; 97416022c9SBarry Smith int *nz,ec,shift = Baij->indexshift; 982493cbb0SBarry Smith Scalar v; 992493cbb0SBarry Smith 1002493cbb0SBarry Smith /* free stuff related to matrix-vec multiply */ 101464493b3SBarry Smith ierr = VecGetSize(aij->lvec,&ec); /* needed for PLogObjectMemory below */ 1022493cbb0SBarry Smith ierr = VecDestroy(aij->lvec); CHKERRQ(ierr); aij->lvec = 0; 103*08480c60SBarry Smith ierr = VecScatterDestroy(aij->Mvctx); CHKERRQ(ierr); aij->Mvctx = 0; 104464493b3SBarry Smith if (aij->colmap) { 105464493b3SBarry Smith PETSCFREE(aij->colmap); aij->colmap = 0; 106464493b3SBarry Smith PLogObjectMemory(A,-Baij->n*sizeof(int)); 107464493b3SBarry Smith } 1082493cbb0SBarry Smith 1092493cbb0SBarry Smith /* make sure that B is assembled so we can access its values */ 1102493cbb0SBarry Smith ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 1112493cbb0SBarry Smith MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 1122493cbb0SBarry Smith 1132493cbb0SBarry Smith /* invent new B and copy stuff over */ 11448b35521SBarry Smith nz = (int *) PETSCMALLOC( m*sizeof(int) ); CHKPTRQ(nz); 11548b35521SBarry Smith for ( i=0; i<m; i++ ) { 11648b35521SBarry Smith nz[i] = Baij->i[i+1]-Baij->i[i]; 11748b35521SBarry Smith } 118fafbff53SBarry Smith ierr = MatCreateSeqAIJ(MPI_COMM_SELF,m,n,0,nz,&Bnew); CHKERRQ(ierr); 11948b35521SBarry Smith PETSCFREE(nz); 1202493cbb0SBarry Smith for ( i=0; i<m; i++ ) { 121dbb450caSBarry Smith for ( j=Baij->i[i]+shift; j<Baij->i[i+1]+shift; j++ ) { 122dbb450caSBarry Smith col = garray[Baij->j[ct]+shift]; 1232493cbb0SBarry Smith v = Baij->a[ct++]; 124dbb450caSBarry Smith ierr = MatSetValues(Bnew,1,&i,1,&col,&v,INSERT_VALUES); CHKERRQ(ierr); 1252493cbb0SBarry Smith } 1262493cbb0SBarry Smith } 1272493cbb0SBarry Smith PETSCFREE(aij->garray); aij->garray = 0; 128464493b3SBarry Smith PLogObjectMemory(A,-ec*sizeof(int)); 1292493cbb0SBarry Smith ierr = MatDestroy(B); CHKERRQ(ierr); 130464493b3SBarry Smith PLogObjectParent(A,Bnew); 1312493cbb0SBarry Smith aij->B = Bnew; 1322493cbb0SBarry Smith aij->assembled = 0; 1332493cbb0SBarry Smith return 0; 1342493cbb0SBarry Smith } 1352493cbb0SBarry Smith 13648b35521SBarry Smith 137