1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*3a40ed3dSBarry Smith static char vcid[] = "$Id: mmaij.c,v 1.37 1997/07/09 20:54:04 balay 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 125615d1e5SSatish Balay #undef __FUNC__ 135615d1e5SSatish Balay #define __FUNC__ "MatSetUpMultiply_MPIAIJ" 1444a69424SLois Curfman McInnes int MatSetUpMultiply_MPIAIJ(Mat mat) 158c79f6d3SBarry Smith { 1644a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 17ec8511deSBarry Smith Mat_SeqAIJ *B = (Mat_SeqAIJ *) (aij->B->data); 18416022c9SBarry Smith int N = aij->N,i,j,*indices,*aj = B->j,ierr,ec = 0,*garray; 19416022c9SBarry Smith int shift = B->indexshift; 201eb62cbbSBarry Smith IS from,to; 211eb62cbbSBarry Smith Vec gvec; 228c79f6d3SBarry Smith 23*3a40ed3dSBarry Smith PetscFunctionBegin; 248c79f6d3SBarry Smith /* For the first stab we make an array as long as the number of columns */ 251eb62cbbSBarry Smith /* mark those columns that are in aij->B */ 26ac03b0baSBarry Smith indices = (int *) PetscMalloc( (N+1)*sizeof(int) ); CHKPTRQ(indices); 27cddf8d76SBarry Smith PetscMemzero(indices,N*sizeof(int)); 28d6dfbf8fSBarry Smith for ( i=0; i<B->m; i++ ) { 29d6dfbf8fSBarry Smith for ( j=0; j<B->ilen[i]; j++ ) { 30dbb450caSBarry Smith if (!indices[aj[B->i[i] +shift + j] + shift]) ec++; 31416022c9SBarry Smith indices[aj[B->i[i] + shift + j] + shift] = 1; 32416022c9SBarry Smith } 331eb62cbbSBarry Smith } 348c79f6d3SBarry Smith 351eb62cbbSBarry Smith /* form array of columns we need */ 360452661fSBarry Smith garray = (int *) PetscMalloc( (ec+1)*sizeof(int) ); CHKPTRQ(garray); 371eb62cbbSBarry Smith ec = 0; 381eb62cbbSBarry Smith for ( i=0; i<N; i++ ) { 391eb62cbbSBarry Smith if (indices[i]) garray[ec++] = i; 401eb62cbbSBarry Smith } 411eb62cbbSBarry Smith 421eb62cbbSBarry Smith /* make indices now point into garray */ 431eb62cbbSBarry Smith for ( i=0; i<ec; i++ ) { 44dbb450caSBarry Smith indices[garray[i]] = i-shift; 451eb62cbbSBarry Smith } 461eb62cbbSBarry Smith 471eb62cbbSBarry Smith /* compact out the extra columns in B */ 48d6dfbf8fSBarry Smith for ( i=0; i<B->m; i++ ) { 49d6dfbf8fSBarry Smith for ( j=0; j<B->ilen[i]; j++ ) { 50dbb450caSBarry Smith aj[B->i[i] + shift + j] = indices[aj[B->i[i] + shift + j]+shift]; 511eb62cbbSBarry Smith } 52d6dfbf8fSBarry Smith } 53639f9d9dSBarry Smith B->n = aij->B->n = aij->B->N = ec; 540452661fSBarry Smith PetscFree(indices); 551eb62cbbSBarry Smith 561eb62cbbSBarry Smith /* create local vector that is used to scatter into */ 57029af93fSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec); CHKERRQ(ierr); 581eb62cbbSBarry Smith 59d6dfbf8fSBarry Smith /* create two temporary Index sets for build scatter gather */ 60029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,ec,garray,&from); CHKERRQ(ierr); 61029af93fSBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to); CHKERRQ(ierr); 621eb62cbbSBarry Smith 631eb62cbbSBarry Smith /* create temporary global vector to generate scatter context */ 641eb62cbbSBarry Smith /* this is inefficient, but otherwise we must do either 651eb62cbbSBarry Smith 1) save garray until the first actual scatter when the vector is known or 661eb62cbbSBarry Smith 2) have another way of generating a scatter context without a vector.*/ 6778b31e54SBarry Smith ierr = VecCreateMPI(mat->comm,aij->n,aij->N,&gvec); CHKERRQ(ierr); 681eb62cbbSBarry Smith 692d336d48SLois Curfman McInnes /* generate the scatter context */ 7008480c60SBarry Smith ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx); CHKERRQ(ierr); 71a5a9c739SBarry Smith PLogObjectParent(mat,aij->Mvctx); 72a5a9c739SBarry Smith PLogObjectParent(mat,aij->lvec); 73464493b3SBarry Smith PLogObjectParent(mat,from); 74464493b3SBarry Smith PLogObjectParent(mat,to); 759e25ed09SBarry Smith aij->garray = garray; 76464493b3SBarry Smith PLogObjectMemory(mat,(ec+1)*sizeof(int)); 7778b31e54SBarry Smith ierr = ISDestroy(from); CHKERRQ(ierr); 7878b31e54SBarry Smith ierr = ISDestroy(to); CHKERRQ(ierr); 791eb62cbbSBarry Smith ierr = VecDestroy(gvec); 80*3a40ed3dSBarry Smith PetscFunctionReturn(0); 818c79f6d3SBarry Smith } 829e25ed09SBarry Smith 839e25ed09SBarry Smith 845615d1e5SSatish Balay #undef __FUNC__ 855615d1e5SSatish Balay #define __FUNC__ "DisAssemble_MPIAIJ" 862493cbb0SBarry Smith /* 872493cbb0SBarry Smith Takes the local part of an already assembled MPIAIJ matrix 882493cbb0SBarry Smith and disassembles it. This is to allow new nonzeros into the matrix 892493cbb0SBarry Smith that require more communication in the matrix vector multiply. 902493cbb0SBarry Smith Thus certain data-structures must be rebuilt. 912493cbb0SBarry Smith 922493cbb0SBarry Smith Kind of slow! But that's what application programmers get when 932493cbb0SBarry Smith they are sloppy. 942493cbb0SBarry Smith */ 952493cbb0SBarry Smith int DisAssemble_MPIAIJ(Mat A) 962493cbb0SBarry Smith { 972493cbb0SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ *) A->data; 982493cbb0SBarry Smith Mat B = aij->B,Bnew; 99ec8511deSBarry Smith Mat_SeqAIJ *Baij = (Mat_SeqAIJ*)B->data; 1002493cbb0SBarry Smith int ierr,i,j,m = Baij->m,n = aij->N,col,ct = 0,*garray = aij->garray; 101416022c9SBarry Smith int *nz,ec,shift = Baij->indexshift; 1022493cbb0SBarry Smith Scalar v; 1032493cbb0SBarry Smith 104*3a40ed3dSBarry Smith PetscFunctionBegin; 1052493cbb0SBarry Smith /* free stuff related to matrix-vec multiply */ 106464493b3SBarry Smith ierr = VecGetSize(aij->lvec,&ec); /* needed for PLogObjectMemory below */ 1072493cbb0SBarry Smith ierr = VecDestroy(aij->lvec); CHKERRQ(ierr); aij->lvec = 0; 10808480c60SBarry Smith ierr = VecScatterDestroy(aij->Mvctx); CHKERRQ(ierr); aij->Mvctx = 0; 109464493b3SBarry Smith if (aij->colmap) { 1100452661fSBarry Smith PetscFree(aij->colmap); aij->colmap = 0; 111464493b3SBarry Smith PLogObjectMemory(A,-Baij->n*sizeof(int)); 112464493b3SBarry Smith } 1132493cbb0SBarry Smith 1142493cbb0SBarry Smith /* make sure that B is assembled so we can access its values */ 1156d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 116fe2f2677SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 1172493cbb0SBarry Smith 1182493cbb0SBarry Smith /* invent new B and copy stuff over */ 119639f9d9dSBarry Smith nz = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(nz); 12048b35521SBarry Smith for ( i=0; i<m; i++ ) { 12148b35521SBarry Smith nz[i] = Baij->i[i+1] - Baij->i[i]; 12248b35521SBarry Smith } 123029af93fSBarry Smith ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m,n,0,nz,&Bnew); CHKERRQ(ierr); 1240452661fSBarry Smith PetscFree(nz); 1252493cbb0SBarry Smith for ( i=0; i<m; i++ ) { 126dbb450caSBarry Smith for ( j=Baij->i[i]+shift; j<Baij->i[i+1]+shift; j++ ) { 127dbb450caSBarry Smith col = garray[Baij->j[ct]+shift]; 1282493cbb0SBarry Smith v = Baij->a[ct++]; 129dbb450caSBarry Smith ierr = MatSetValues(Bnew,1,&i,1,&col,&v,INSERT_VALUES); CHKERRQ(ierr); 1302493cbb0SBarry Smith } 1312493cbb0SBarry Smith } 1320452661fSBarry Smith PetscFree(aij->garray); aij->garray = 0; 133464493b3SBarry Smith PLogObjectMemory(A,-ec*sizeof(int)); 1342493cbb0SBarry Smith ierr = MatDestroy(B); CHKERRQ(ierr); 135464493b3SBarry Smith PLogObjectParent(A,Bnew); 1362493cbb0SBarry Smith aij->B = Bnew; 137227d817aSBarry Smith A->was_assembled = PETSC_FALSE; 138*3a40ed3dSBarry Smith PetscFunctionReturn(0); 1392493cbb0SBarry Smith } 1402493cbb0SBarry Smith 14148b35521SBarry Smith 142