1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*606d414cSSatish Balay static char vcid[] = "$Id: mmaij.c,v 1.47 1999/06/08 22:55:49 balay Exp balay $"; 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; 22aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 232066d6f7SSatish Balay Table gid1_lid1; 242066d6f7SSatish Balay CTablePos tpos; 252066d6f7SSatish Balay int gid, lid; 262066d6f7SSatish Balay #endif 272066d6f7SSatish Balay 283a40ed3dSBarry Smith PetscFunctionBegin; 292066d6f7SSatish Balay 30aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 312066d6f7SSatish Balay /* use a table - Mark Adams (this has not been tested with "shift") */ 32fa46199cSSatish Balay TableCreate(B->m,&gid1_lid1); 332066d6f7SSatish Balay for ( i=0; i<B->m; i++ ) { 342066d6f7SSatish Balay for ( j=0; j<B->ilen[i]; j++ ) { 3508c73f0fSSatish Balay int data,gid1 = aj[B->i[i] + shift + j] + 1 + shift; 36fa46199cSSatish Balay ierr = TableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr); 37fa46199cSSatish Balay if (!data) { 382066d6f7SSatish Balay /* one based table */ 392066d6f7SSatish Balay ierr = TableAdd(gid1_lid1,gid1,++ec);CHKERRQ(ierr); 402066d6f7SSatish Balay } 412066d6f7SSatish Balay } 422066d6f7SSatish Balay } 432066d6f7SSatish Balay /* form array of columns we need */ 442066d6f7SSatish Balay garray = (int *)PetscMalloc((ec+1)*sizeof(int));CHKPTRQ(garray); 452066d6f7SSatish Balay ierr = TableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr); 462066d6f7SSatish Balay while (tpos) { 472066d6f7SSatish Balay ierr = TableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr); 482066d6f7SSatish Balay gid--; lid--; 492066d6f7SSatish Balay garray[lid] = gid; 502066d6f7SSatish Balay } 510064e2bbSSatish Balay ierr = PetscSortInt(ec,garray);CHKERRQ(ierr); /* sort, and rebuild */ 520064e2bbSSatish Balay /* qsort( garray, ec, sizeof(int), intcomparc ); */ 532066d6f7SSatish Balay TableRemoveAll(gid1_lid1); 542066d6f7SSatish Balay for ( i=0; i<ec; i++ ) { 552066d6f7SSatish Balay ierr = TableAdd(gid1_lid1,garray[i]+1,i+1);CHKERRQ(ierr); 562066d6f7SSatish Balay } 572066d6f7SSatish Balay /* compact out the extra columns in B */ 582066d6f7SSatish Balay for ( i=0; i<B->m; i++ ) { 592066d6f7SSatish Balay for ( j=0; j<B->ilen[i]; j++ ) { 6008c73f0fSSatish Balay int gid1 = aj[B->i[i] + shift + j] + 1 + shift; 61fa46199cSSatish Balay ierr = TableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr); 62fa46199cSSatish Balay lid --; 6308c73f0fSSatish Balay aj[B->i[i] + shift + j] = lid - shift; 642066d6f7SSatish Balay } 652066d6f7SSatish Balay } 662066d6f7SSatish Balay B->n = aij->B->n = aij->B->N = ec; 672066d6f7SSatish Balay TableDelete(gid1_lid1); 682066d6f7SSatish Balay /* Mark Adams */ 692066d6f7SSatish Balay #else 708c79f6d3SBarry Smith /* For the first stab we make an array as long as the number of columns */ 711eb62cbbSBarry Smith /* mark those columns that are in aij->B */ 72ac03b0baSBarry Smith indices = (int *) PetscMalloc( (N+1)*sizeof(int) );CHKPTRQ(indices); 73549d3d68SSatish Balay ierr = PetscMemzero(indices,N*sizeof(int));CHKERRQ(ierr); 74d6dfbf8fSBarry Smith for ( i=0; i<B->m; i++ ) { 75d6dfbf8fSBarry Smith for ( j=0; j<B->ilen[i]; j++ ) { 76dbb450caSBarry Smith if (!indices[aj[B->i[i] +shift + j] + shift]) ec++; 77416022c9SBarry Smith indices[aj[B->i[i] + shift + j] + shift] = 1; 78416022c9SBarry Smith } 791eb62cbbSBarry Smith } 808c79f6d3SBarry Smith 811eb62cbbSBarry Smith /* form array of columns we need */ 820452661fSBarry Smith garray = (int *) PetscMalloc( (ec+1)*sizeof(int) );CHKPTRQ(garray); 831eb62cbbSBarry Smith ec = 0; 841eb62cbbSBarry Smith for ( i=0; i<N; i++ ) { 851eb62cbbSBarry Smith if (indices[i]) garray[ec++] = i; 861eb62cbbSBarry Smith } 871eb62cbbSBarry Smith 881eb62cbbSBarry Smith /* make indices now point into garray */ 891eb62cbbSBarry Smith for ( i=0; i<ec; i++ ) { 90dbb450caSBarry Smith indices[garray[i]] = i-shift; 911eb62cbbSBarry Smith } 921eb62cbbSBarry Smith 931eb62cbbSBarry Smith /* compact out the extra columns in B */ 94d6dfbf8fSBarry Smith for ( i=0; i<B->m; i++ ) { 95d6dfbf8fSBarry Smith for ( j=0; j<B->ilen[i]; j++ ) { 96dbb450caSBarry Smith aj[B->i[i] + shift + j] = indices[aj[B->i[i] + shift + j]+shift]; 971eb62cbbSBarry Smith } 98d6dfbf8fSBarry Smith } 99639f9d9dSBarry Smith B->n = aij->B->n = aij->B->N = ec; 100*606d414cSSatish Balay ierr = PetscFree(indices);CHKERRQ(ierr); 1012066d6f7SSatish Balay #endif 1021eb62cbbSBarry Smith /* create local vector that is used to scatter into */ 103029af93fSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);CHKERRQ(ierr); 1041eb62cbbSBarry Smith 105d6dfbf8fSBarry Smith /* create two temporary Index sets for build scatter gather */ 106029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,ec,garray,&from);CHKERRQ(ierr); 107029af93fSBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr); 1081eb62cbbSBarry Smith 1091eb62cbbSBarry Smith /* create temporary global vector to generate scatter context */ 1101eb62cbbSBarry Smith /* this is inefficient, but otherwise we must do either 1111eb62cbbSBarry Smith 1) save garray until the first actual scatter when the vector is known or 1121eb62cbbSBarry Smith 2) have another way of generating a scatter context without a vector.*/ 11378b31e54SBarry Smith ierr = VecCreateMPI(mat->comm,aij->n,aij->N,&gvec);CHKERRQ(ierr); 1141eb62cbbSBarry Smith 1152d336d48SLois Curfman McInnes /* generate the scatter context */ 11608480c60SBarry Smith ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);CHKERRQ(ierr); 117a5a9c739SBarry Smith PLogObjectParent(mat,aij->Mvctx); 118a5a9c739SBarry Smith PLogObjectParent(mat,aij->lvec); 119464493b3SBarry Smith PLogObjectParent(mat,from); 120464493b3SBarry Smith PLogObjectParent(mat,to); 1219e25ed09SBarry Smith aij->garray = garray; 122464493b3SBarry Smith PLogObjectMemory(mat,(ec+1)*sizeof(int)); 12378b31e54SBarry Smith ierr = ISDestroy(from);CHKERRQ(ierr); 12478b31e54SBarry Smith ierr = ISDestroy(to);CHKERRQ(ierr); 125888f2ed8SSatish Balay ierr = VecDestroy(gvec);CHKERRQ(ierr); 1263a40ed3dSBarry Smith PetscFunctionReturn(0); 1278c79f6d3SBarry Smith } 1289e25ed09SBarry Smith 1299e25ed09SBarry Smith 1305615d1e5SSatish Balay #undef __FUNC__ 1315615d1e5SSatish Balay #define __FUNC__ "DisAssemble_MPIAIJ" 1322493cbb0SBarry Smith /* 1332493cbb0SBarry Smith Takes the local part of an already assembled MPIAIJ matrix 1342493cbb0SBarry Smith and disassembles it. This is to allow new nonzeros into the matrix 1352493cbb0SBarry Smith that require more communication in the matrix vector multiply. 1362493cbb0SBarry Smith Thus certain data-structures must be rebuilt. 1372493cbb0SBarry Smith 1382493cbb0SBarry Smith Kind of slow! But that's what application programmers get when 1392493cbb0SBarry Smith they are sloppy. 1402493cbb0SBarry Smith */ 1412493cbb0SBarry Smith int DisAssemble_MPIAIJ(Mat A) 1422493cbb0SBarry Smith { 1432493cbb0SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ *) A->data; 1442493cbb0SBarry Smith Mat B = aij->B,Bnew; 145ec8511deSBarry Smith Mat_SeqAIJ *Baij = (Mat_SeqAIJ*)B->data; 1462493cbb0SBarry Smith int ierr,i,j,m = Baij->m,n = aij->N,col,ct = 0,*garray = aij->garray; 147416022c9SBarry Smith int *nz,ec,shift = Baij->indexshift; 1482493cbb0SBarry Smith Scalar v; 1492493cbb0SBarry Smith 1503a40ed3dSBarry Smith PetscFunctionBegin; 1512493cbb0SBarry Smith /* free stuff related to matrix-vec multiply */ 152888f2ed8SSatish Balay ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); /* needed for PLogObjectMemory below */ 1532493cbb0SBarry Smith ierr = VecDestroy(aij->lvec);CHKERRQ(ierr); aij->lvec = 0; 15408480c60SBarry Smith ierr = VecScatterDestroy(aij->Mvctx);CHKERRQ(ierr); aij->Mvctx = 0; 155464493b3SBarry Smith if (aij->colmap) { 156aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 1572066d6f7SSatish Balay TableDelete(aij->colmap); aij->colmap = 0; 1582066d6f7SSatish Balay #else 159*606d414cSSatish Balay ierr = PetscFree(aij->colmap);CHKERRQ(ierr); 160*606d414cSSatish Balay aij->colmap = 0; 161464493b3SBarry Smith PLogObjectMemory(A,-Baij->n*sizeof(int)); 1622066d6f7SSatish Balay #endif 163464493b3SBarry Smith } 1642493cbb0SBarry Smith 1652493cbb0SBarry Smith /* make sure that B is assembled so we can access its values */ 1666d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 167fe2f2677SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1682493cbb0SBarry Smith 1692493cbb0SBarry Smith /* invent new B and copy stuff over */ 170639f9d9dSBarry Smith nz = (int *) PetscMalloc( (m+1)*sizeof(int) );CHKPTRQ(nz); 17148b35521SBarry Smith for ( i=0; i<m; i++ ) { 17248b35521SBarry Smith nz[i] = Baij->i[i+1] - Baij->i[i]; 17348b35521SBarry Smith } 174029af93fSBarry Smith ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m,n,0,nz,&Bnew);CHKERRQ(ierr); 175*606d414cSSatish Balay ierr = PetscFree(nz);CHKERRQ(ierr); 1762493cbb0SBarry Smith for ( i=0; i<m; i++ ) { 177dbb450caSBarry Smith for ( j=Baij->i[i]+shift; j<Baij->i[i+1]+shift; j++ ) { 178dbb450caSBarry Smith col = garray[Baij->j[ct]+shift]; 1792493cbb0SBarry Smith v = Baij->a[ct++]; 18083271157SBarry Smith ierr = MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);CHKERRQ(ierr); 1812493cbb0SBarry Smith } 1822493cbb0SBarry Smith } 183*606d414cSSatish Balay ierr = PetscFree(aij->garray);CHKERRQ(ierr); 184*606d414cSSatish Balay aij->garray = 0; 185464493b3SBarry Smith PLogObjectMemory(A,-ec*sizeof(int)); 1862493cbb0SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 187464493b3SBarry Smith PLogObjectParent(A,Bnew); 1882493cbb0SBarry Smith aij->B = Bnew; 189227d817aSBarry Smith A->was_assembled = PETSC_FALSE; 1903a40ed3dSBarry Smith PetscFunctionReturn(0); 1912493cbb0SBarry Smith } 1922493cbb0SBarry Smith 19348b35521SBarry Smith 194