1cb512458SBarry Smith #ifndef lint 2*48b35521SBarry Smith static char vcid[] = "$Id: mmaij.c,v 1.12 1995/06/30 03:35:15 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; 1644a69424SLois Curfman McInnes Mat_AIJ *B = (Mat_AIJ *) (aij->B->data); 176abc6512SBarry Smith int N = aij->N,i,j,*indices,*aj = B->j; 181eb62cbbSBarry Smith int ierr,ec = 0,*garray; 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); 2578b31e54SBarry Smith PETSCMEMSET(indices,0,N*sizeof(int)); 26d6dfbf8fSBarry Smith for ( i=0; i<B->m; i++ ) { 27d6dfbf8fSBarry Smith for ( j=0; j<B->ilen[i]; j++ ) { 28d6dfbf8fSBarry Smith if (!indices[aj[B->i[i] - 1 + j]-1]) ec++; 29d6dfbf8fSBarry Smith indices[aj[B->i[i] - 1 + j]-1] = 1;} 301eb62cbbSBarry Smith } 318c79f6d3SBarry Smith 321eb62cbbSBarry Smith /* form array of columns we need */ 3378b31e54SBarry 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++ ) { 411eb62cbbSBarry Smith indices[garray[i]] = i+1; 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++ ) { 47d6dfbf8fSBarry Smith aj[B->i[i] - 1 + j] = indices[aj[B->i[i] - 1 + j]-1]; 481eb62cbbSBarry Smith } 49d6dfbf8fSBarry Smith } 50d6dfbf8fSBarry Smith B->n = ec; 5178b31e54SBarry Smith PETSCFREE(indices); 521eb62cbbSBarry Smith 531eb62cbbSBarry Smith /* create local vector that is used to scatter into */ 5478b31e54SBarry Smith ierr = VecCreateSequential(MPI_COMM_SELF,ec,&aij->lvec); CHKERRQ(ierr); 551eb62cbbSBarry Smith 56d6dfbf8fSBarry Smith /* create two temporary Index sets for build scatter gather */ 5778b31e54SBarry Smith ierr = ISCreateSequential(MPI_COMM_SELF,ec,garray,&from); CHKERRQ(ierr); 5878b31e54SBarry Smith ierr = ISCreateStrideSequential(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 661eb62cbbSBarry Smith /* gnerate the scatter context */ 6778b31e54SBarry Smith ierr = VecScatterCtxCreate(gvec,from,aij->lvec,to,&aij->Mvctx); CHKERRQ(ierr); 68a5a9c739SBarry Smith PLogObjectParent(mat,aij->Mvctx); 69a5a9c739SBarry Smith PLogObjectParent(mat,aij->lvec); 709e25ed09SBarry Smith aij->garray = garray; 7178b31e54SBarry Smith ierr = ISDestroy(from); CHKERRQ(ierr); 7278b31e54SBarry Smith ierr = ISDestroy(to); CHKERRQ(ierr); 731eb62cbbSBarry Smith ierr = VecDestroy(gvec); 748c79f6d3SBarry Smith return 0; 758c79f6d3SBarry Smith } 769e25ed09SBarry Smith 779e25ed09SBarry Smith 782493cbb0SBarry Smith /* 792493cbb0SBarry Smith Takes the local part of an already assembled MPIAIJ matrix 802493cbb0SBarry Smith and disassembles it. This is to allow new nonzeros into the matrix 812493cbb0SBarry Smith that require more communication in the matrix vector multiply. 822493cbb0SBarry Smith Thus certain data-structures must be rebuilt. 832493cbb0SBarry Smith 842493cbb0SBarry Smith Kind of slow! But that's what application programmers get when 852493cbb0SBarry Smith they are sloppy. 862493cbb0SBarry Smith */ 872493cbb0SBarry Smith int DisAssemble_MPIAIJ(Mat A) 882493cbb0SBarry Smith { 892493cbb0SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ *) A->data; 902493cbb0SBarry Smith Mat B = aij->B,Bnew; 912493cbb0SBarry Smith Mat_AIJ *Baij = (Mat_AIJ*)B->data; 922493cbb0SBarry Smith int ierr,i,j,m=Baij->m,n = aij->N,col,ct = 0,*garray = aij->garray; 93*48b35521SBarry Smith int *nz; 942493cbb0SBarry Smith Scalar v; 952493cbb0SBarry Smith 96*48b35521SBarry Smith /* 972493cbb0SBarry Smith fprintf(stderr,"Warning: you are adding tricky new nonzeros\n"); 98*48b35521SBarry Smith */ 992493cbb0SBarry Smith 1002493cbb0SBarry Smith /* free stuff related to matrix-vec multiply */ 1012493cbb0SBarry Smith ierr = VecDestroy(aij->lvec); CHKERRQ(ierr); aij->lvec = 0; 1022493cbb0SBarry Smith ierr = VecScatterCtxDestroy(aij->Mvctx); CHKERRQ(ierr); aij->Mvctx = 0; 1032493cbb0SBarry Smith PETSCFREE(aij->colmap); aij->colmap = 0; 1042493cbb0SBarry Smith 1052493cbb0SBarry Smith /* make sure that B is assembled so we can access its values */ 1062493cbb0SBarry Smith ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 1072493cbb0SBarry Smith MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 1082493cbb0SBarry Smith 1092493cbb0SBarry Smith /* invent new B and copy stuff over */ 110*48b35521SBarry Smith nz = (int *) PETSCMALLOC( m*sizeof(int) ); CHKPTRQ(nz); 111*48b35521SBarry Smith for ( i=0; i<m; i++ ) { 112*48b35521SBarry Smith nz[i] = Baij->i[i+1]-Baij->i[i]; 113*48b35521SBarry Smith } 114*48b35521SBarry Smith ierr = MatCreateSequentialAIJ(MPI_COMM_SELF,m,n,0,nz,&Bnew); CHKERRQ(ierr); 115*48b35521SBarry Smith PETSCFREE(nz); 1162493cbb0SBarry Smith for ( i=0; i<m; i++ ) { 1172493cbb0SBarry Smith for ( j=Baij->i[i]-1; j<Baij->i[i+1]-1; j++ ) { 1182493cbb0SBarry Smith col = garray[Baij->j[ct]-1]; 1192493cbb0SBarry Smith v = Baij->a[ct++]; 1202493cbb0SBarry Smith ierr = MatSetValues(Bnew,1,&i,1,&col,&v,INSERTVALUES); CHKERRQ(ierr); 1212493cbb0SBarry Smith } 1222493cbb0SBarry Smith } 1232493cbb0SBarry Smith PETSCFREE(aij->garray); aij->garray = 0; 1242493cbb0SBarry Smith ierr = MatDestroy(B); CHKERRQ(ierr); 1252493cbb0SBarry Smith aij->B = Bnew; 1262493cbb0SBarry Smith aij->assembled = 0; 1272493cbb0SBarry Smith return 0; 1282493cbb0SBarry Smith } 1292493cbb0SBarry Smith 130*48b35521SBarry Smith 131