1*8016bdd1SSatish Balay #ifndef lint 2*8016bdd1SSatish Balay static char vcid[] = "$Id: mmaij.c,v 1.26 1996/04/26 00:03:57 balay Exp $"; 3*8016bdd1SSatish Balay #endif 4*8016bdd1SSatish Balay 5*8016bdd1SSatish Balay 6*8016bdd1SSatish Balay /* 7*8016bdd1SSatish Balay Support for the parallel AIJ matrix vector multiply 8*8016bdd1SSatish Balay */ 9*8016bdd1SSatish Balay #include "mpiaij.h" 10*8016bdd1SSatish Balay #include "src/vec/vecimpl.h" 11*8016bdd1SSatish Balay #include "../seq/aij.h" 12*8016bdd1SSatish Balay 13*8016bdd1SSatish Balay int MatSetUpMultiply_MPIAIJ(Mat mat) 14*8016bdd1SSatish Balay { 15*8016bdd1SSatish Balay Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 16*8016bdd1SSatish Balay Mat_SeqAIJ *B = (Mat_SeqAIJ *) (aij->B->data); 17*8016bdd1SSatish Balay int N = aij->N,i,j,*indices,*aj = B->j,ierr,ec = 0,*garray; 18*8016bdd1SSatish Balay int shift = B->indexshift; 19*8016bdd1SSatish Balay IS from,to; 20*8016bdd1SSatish Balay Vec gvec; 21*8016bdd1SSatish Balay 22*8016bdd1SSatish Balay /* For the first stab we make an array as long as the number of columns */ 23*8016bdd1SSatish Balay /* mark those columns that are in aij->B */ 24*8016bdd1SSatish Balay indices = (int *) PetscMalloc( N*sizeof(int) ); CHKPTRQ(indices); 25*8016bdd1SSatish Balay PetscMemzero(indices,N*sizeof(int)); 26*8016bdd1SSatish Balay for ( i=0; i<B->m; i++ ) { 27*8016bdd1SSatish Balay for ( j=0; j<B->ilen[i]; j++ ) { 28*8016bdd1SSatish Balay if (!indices[aj[B->i[i] +shift + j] + shift]) ec++; 29*8016bdd1SSatish Balay indices[aj[B->i[i] + shift + j] + shift] = 1; 30*8016bdd1SSatish Balay } 31*8016bdd1SSatish Balay } 32*8016bdd1SSatish Balay 33*8016bdd1SSatish Balay /* form array of columns we need */ 34*8016bdd1SSatish Balay garray = (int *) PetscMalloc( (ec+1)*sizeof(int) ); CHKPTRQ(garray); 35*8016bdd1SSatish Balay ec = 0; 36*8016bdd1SSatish Balay for ( i=0; i<N; i++ ) { 37*8016bdd1SSatish Balay if (indices[i]) garray[ec++] = i; 38*8016bdd1SSatish Balay } 39*8016bdd1SSatish Balay 40*8016bdd1SSatish Balay /* make indices now point into garray */ 41*8016bdd1SSatish Balay for ( i=0; i<ec; i++ ) { 42*8016bdd1SSatish Balay indices[garray[i]] = i-shift; 43*8016bdd1SSatish Balay } 44*8016bdd1SSatish Balay 45*8016bdd1SSatish Balay /* compact out the extra columns in B */ 46*8016bdd1SSatish Balay for ( i=0; i<B->m; i++ ) { 47*8016bdd1SSatish Balay for ( j=0; j<B->ilen[i]; j++ ) { 48*8016bdd1SSatish Balay aj[B->i[i] + shift + j] = indices[aj[B->i[i] + shift + j]+shift]; 49*8016bdd1SSatish Balay } 50*8016bdd1SSatish Balay } 51*8016bdd1SSatish Balay B->n = ec; 52*8016bdd1SSatish Balay PetscFree(indices); 53*8016bdd1SSatish Balay 54*8016bdd1SSatish Balay /* create local vector that is used to scatter into */ 55*8016bdd1SSatish Balay ierr = VecCreateSeq(MPI_COMM_SELF,ec,&aij->lvec); CHKERRQ(ierr); 56*8016bdd1SSatish Balay 57*8016bdd1SSatish Balay /* create two temporary Index sets for build scatter gather */ 58*8016bdd1SSatish Balay ierr = ISCreateSeq(MPI_COMM_SELF,ec,garray,&from); CHKERRQ(ierr); 59*8016bdd1SSatish Balay ierr = ISCreateStrideSeq(MPI_COMM_SELF,ec,0,1,&to); CHKERRQ(ierr); 60*8016bdd1SSatish Balay 61*8016bdd1SSatish Balay /* create temporary global vector to generate scatter context */ 62*8016bdd1SSatish Balay /* this is inefficient, but otherwise we must do either 63*8016bdd1SSatish Balay 1) save garray until the first actual scatter when the vector is known or 64*8016bdd1SSatish Balay 2) have another way of generating a scatter context without a vector.*/ 65*8016bdd1SSatish Balay ierr = VecCreateMPI(mat->comm,aij->n,aij->N,&gvec); CHKERRQ(ierr); 66*8016bdd1SSatish Balay 67*8016bdd1SSatish Balay /* gnerate the scatter context */ 68*8016bdd1SSatish Balay ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx); CHKERRQ(ierr); 69*8016bdd1SSatish Balay PLogObjectParent(mat,aij->Mvctx); 70*8016bdd1SSatish Balay PLogObjectParent(mat,aij->lvec); 71*8016bdd1SSatish Balay PLogObjectParent(mat,from); 72*8016bdd1SSatish Balay PLogObjectParent(mat,to); 73*8016bdd1SSatish Balay aij->garray = garray; 74*8016bdd1SSatish Balay PLogObjectMemory(mat,(ec+1)*sizeof(int)); 75*8016bdd1SSatish Balay ierr = ISDestroy(from); CHKERRQ(ierr); 76*8016bdd1SSatish Balay ierr = ISDestroy(to); CHKERRQ(ierr); 77*8016bdd1SSatish Balay ierr = VecDestroy(gvec); 78*8016bdd1SSatish Balay return 0; 79*8016bdd1SSatish Balay } 80*8016bdd1SSatish Balay 81*8016bdd1SSatish Balay 82*8016bdd1SSatish Balay /* 83*8016bdd1SSatish Balay Takes the local part of an already assembled MPIAIJ matrix 84*8016bdd1SSatish Balay and disassembles it. This is to allow new nonzeros into the matrix 85*8016bdd1SSatish Balay that require more communication in the matrix vector multiply. 86*8016bdd1SSatish Balay Thus certain data-structures must be rebuilt. 87*8016bdd1SSatish Balay 88*8016bdd1SSatish Balay Kind of slow! But that's what application programmers get when 89*8016bdd1SSatish Balay they are sloppy. 90*8016bdd1SSatish Balay */ 91*8016bdd1SSatish Balay int DisAssemble_MPIAIJ(Mat A) 92*8016bdd1SSatish Balay { 93*8016bdd1SSatish Balay Mat_MPIAIJ *aij = (Mat_MPIAIJ *) A->data; 94*8016bdd1SSatish Balay Mat B = aij->B,Bnew; 95*8016bdd1SSatish Balay Mat_SeqAIJ *Baij = (Mat_SeqAIJ*)B->data; 96*8016bdd1SSatish Balay int ierr,i,j,m=Baij->m,n = aij->N,col,ct = 0,*garray = aij->garray; 97*8016bdd1SSatish Balay int *nz,ec,shift = Baij->indexshift; 98*8016bdd1SSatish Balay Scalar v; 99*8016bdd1SSatish Balay 100*8016bdd1SSatish Balay /* free stuff related to matrix-vec multiply */ 101*8016bdd1SSatish Balay ierr = VecGetSize(aij->lvec,&ec); /* needed for PLogObjectMemory below */ 102*8016bdd1SSatish Balay ierr = VecDestroy(aij->lvec); CHKERRQ(ierr); aij->lvec = 0; 103*8016bdd1SSatish Balay ierr = VecScatterDestroy(aij->Mvctx); CHKERRQ(ierr); aij->Mvctx = 0; 104*8016bdd1SSatish Balay if (aij->colmap) { 105*8016bdd1SSatish Balay PetscFree(aij->colmap); aij->colmap = 0; 106*8016bdd1SSatish Balay PLogObjectMemory(A,-Baij->n*sizeof(int)); 107*8016bdd1SSatish Balay } 108*8016bdd1SSatish Balay 109*8016bdd1SSatish Balay /* make sure that B is assembled so we can access its values */ 110*8016bdd1SSatish Balay ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 111*8016bdd1SSatish Balay MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 112*8016bdd1SSatish Balay 113*8016bdd1SSatish Balay /* invent new B and copy stuff over */ 114*8016bdd1SSatish Balay nz = (int *) PetscMalloc( m*sizeof(int) ); CHKPTRQ(nz); 115*8016bdd1SSatish Balay for ( i=0; i<m; i++ ) { 116*8016bdd1SSatish Balay nz[i] = Baij->i[i+1]-Baij->i[i]; 117*8016bdd1SSatish Balay } 118*8016bdd1SSatish Balay ierr = MatCreateSeqAIJ(MPI_COMM_SELF,m,n,0,nz,&Bnew); CHKERRQ(ierr); 119*8016bdd1SSatish Balay PetscFree(nz); 120*8016bdd1SSatish Balay for ( i=0; i<m; i++ ) { 121*8016bdd1SSatish Balay for ( j=Baij->i[i]+shift; j<Baij->i[i+1]+shift; j++ ) { 122*8016bdd1SSatish Balay col = garray[Baij->j[ct]+shift]; 123*8016bdd1SSatish Balay v = Baij->a[ct++]; 124*8016bdd1SSatish Balay ierr = MatSetValues(Bnew,1,&i,1,&col,&v,INSERT_VALUES); CHKERRQ(ierr); 125*8016bdd1SSatish Balay } 126*8016bdd1SSatish Balay } 127*8016bdd1SSatish Balay PetscFree(aij->garray); aij->garray = 0; 128*8016bdd1SSatish Balay PLogObjectMemory(A,-ec*sizeof(int)); 129*8016bdd1SSatish Balay ierr = MatDestroy(B); CHKERRQ(ierr); 130*8016bdd1SSatish Balay PLogObjectParent(A,Bnew); 131*8016bdd1SSatish Balay aij->B = Bnew; 132*8016bdd1SSatish Balay A->was_assembled = PETSC_FALSE; 133*8016bdd1SSatish Balay return 0; 134*8016bdd1SSatish Balay } 135*8016bdd1SSatish Balay 136*8016bdd1SSatish Balay 137