1 #ifdef PETSC_RCS_HEADER 2 static char vcid[] = "$Id: mmaij.c,v 1.48 1999/06/30 23:51:10 balay Exp bsmith $"; 3 #endif 4 5 /* 6 Support for the parallel AIJ matrix vector multiply 7 */ 8 #include "src/mat/impls/aij/mpi/mpiaij.h" 9 #include "src/vec/vecimpl.h" 10 11 #undef __FUNC__ 12 #define __FUNC__ "MatSetUpMultiply_MPIAIJ" 13 int MatSetUpMultiply_MPIAIJ(Mat mat) 14 { 15 Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 16 Mat_SeqAIJ *B = (Mat_SeqAIJ *) (aij->B->data); 17 int N = aij->N,i,j,*indices,*aj = B->j,ierr,ec = 0,*garray; 18 int shift = B->indexshift; 19 IS from,to; 20 Vec gvec; 21 #if defined (PETSC_USE_CTABLE) 22 PetscTable gid1_lid1; 23 PetscTablePosition tpos; 24 int gid, lid; 25 #endif 26 27 PetscFunctionBegin; 28 29 #if defined (PETSC_USE_CTABLE) 30 /* use a table - Mark Adams (this has not been tested with "shift") */ 31 PetscTableCreate(B->m,&gid1_lid1); 32 for ( i=0; i<B->m; i++ ) { 33 for ( j=0; j<B->ilen[i]; j++ ) { 34 int data,gid1 = aj[B->i[i] + shift + j] + 1 + shift; 35 ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr); 36 if (!data) { 37 /* one based table */ 38 ierr = PetscTableAdd(gid1_lid1,gid1,++ec);CHKERRQ(ierr); 39 } 40 } 41 } 42 /* form array of columns we need */ 43 garray = (int *)PetscMalloc((ec+1)*sizeof(int));CHKPTRQ(garray); 44 ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr); 45 while (tpos) { 46 ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr); 47 gid--; lid--; 48 garray[lid] = gid; 49 } 50 ierr = PetscSortInt(ec,garray);CHKERRQ(ierr); /* sort, and rebuild */ 51 /* qsort( garray, ec, sizeof(int), intcomparc ); */ 52 ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr); 53 for ( i=0; i<ec; i++ ) { 54 ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1);CHKERRQ(ierr); 55 } 56 /* compact out the extra columns in B */ 57 for ( i=0; i<B->m; i++ ) { 58 for ( j=0; j<B->ilen[i]; j++ ) { 59 int gid1 = aj[B->i[i] + shift + j] + 1 + shift; 60 ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr); 61 lid --; 62 aj[B->i[i] + shift + j] = lid - shift; 63 } 64 } 65 B->n = aij->B->n = aij->B->N = ec; 66 ierr = PetscTableDelete(gid1_lid1);CHKERRQ(ierr); 67 /* Mark Adams */ 68 #else 69 /* For the first stab we make an array as long as the number of columns */ 70 /* mark those columns that are in aij->B */ 71 indices = (int *) PetscMalloc( (N+1)*sizeof(int) );CHKPTRQ(indices); 72 ierr = PetscMemzero(indices,N*sizeof(int));CHKERRQ(ierr); 73 for ( i=0; i<B->m; i++ ) { 74 for ( j=0; j<B->ilen[i]; j++ ) { 75 if (!indices[aj[B->i[i] +shift + j] + shift]) ec++; 76 indices[aj[B->i[i] + shift + j] + shift] = 1; 77 } 78 } 79 80 /* form array of columns we need */ 81 garray = (int *) PetscMalloc( (ec+1)*sizeof(int) );CHKPTRQ(garray); 82 ec = 0; 83 for ( i=0; i<N; i++ ) { 84 if (indices[i]) garray[ec++] = i; 85 } 86 87 /* make indices now point into garray */ 88 for ( i=0; i<ec; i++ ) { 89 indices[garray[i]] = i-shift; 90 } 91 92 /* compact out the extra columns in B */ 93 for ( i=0; i<B->m; i++ ) { 94 for ( j=0; j<B->ilen[i]; j++ ) { 95 aj[B->i[i] + shift + j] = indices[aj[B->i[i] + shift + j]+shift]; 96 } 97 } 98 B->n = aij->B->n = aij->B->N = ec; 99 ierr = PetscFree(indices);CHKERRQ(ierr); 100 #endif 101 /* create local vector that is used to scatter into */ 102 ierr = VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);CHKERRQ(ierr); 103 104 /* create two temporary Index sets for build scatter gather */ 105 ierr = ISCreateGeneral(PETSC_COMM_SELF,ec,garray,&from);CHKERRQ(ierr); 106 ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr); 107 108 /* create temporary global vector to generate scatter context */ 109 /* this is inefficient, but otherwise we must do either 110 1) save garray until the first actual scatter when the vector is known or 111 2) have another way of generating a scatter context without a vector.*/ 112 ierr = VecCreateMPI(mat->comm,aij->n,aij->N,&gvec);CHKERRQ(ierr); 113 114 /* generate the scatter context */ 115 ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);CHKERRQ(ierr); 116 PLogObjectParent(mat,aij->Mvctx); 117 PLogObjectParent(mat,aij->lvec); 118 PLogObjectParent(mat,from); 119 PLogObjectParent(mat,to); 120 aij->garray = garray; 121 PLogObjectMemory(mat,(ec+1)*sizeof(int)); 122 ierr = ISDestroy(from);CHKERRQ(ierr); 123 ierr = ISDestroy(to);CHKERRQ(ierr); 124 ierr = VecDestroy(gvec);CHKERRQ(ierr); 125 PetscFunctionReturn(0); 126 } 127 128 129 #undef __FUNC__ 130 #define __FUNC__ "DisAssemble_MPIAIJ" 131 /* 132 Takes the local part of an already assembled MPIAIJ matrix 133 and disassembles it. This is to allow new nonzeros into the matrix 134 that require more communication in the matrix vector multiply. 135 Thus certain data-structures must be rebuilt. 136 137 Kind of slow! But that's what application programmers get when 138 they are sloppy. 139 */ 140 int DisAssemble_MPIAIJ(Mat A) 141 { 142 Mat_MPIAIJ *aij = (Mat_MPIAIJ *) A->data; 143 Mat B = aij->B,Bnew; 144 Mat_SeqAIJ *Baij = (Mat_SeqAIJ*)B->data; 145 int ierr,i,j,m = Baij->m,n = aij->N,col,ct = 0,*garray = aij->garray; 146 int *nz,ec,shift = Baij->indexshift; 147 Scalar v; 148 149 PetscFunctionBegin; 150 /* free stuff related to matrix-vec multiply */ 151 ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); /* needed for PLogObjectMemory below */ 152 ierr = VecDestroy(aij->lvec);CHKERRQ(ierr); aij->lvec = 0; 153 ierr = VecScatterDestroy(aij->Mvctx);CHKERRQ(ierr); aij->Mvctx = 0; 154 if (aij->colmap) { 155 #if defined (PETSC_USE_CTABLE) 156 ierr = PetscTableDelete(aij->colmap); CHKERRQ(ierr); 157 aij->colmap = 0; 158 #else 159 ierr = PetscFree(aij->colmap);CHKERRQ(ierr); 160 aij->colmap = 0; 161 PLogObjectMemory(A,-Baij->n*sizeof(int)); 162 #endif 163 } 164 165 /* make sure that B is assembled so we can access its values */ 166 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 167 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 168 169 /* invent new B and copy stuff over */ 170 nz = (int *) PetscMalloc( (m+1)*sizeof(int) );CHKPTRQ(nz); 171 for ( i=0; i<m; i++ ) { 172 nz[i] = Baij->i[i+1] - Baij->i[i]; 173 } 174 ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m,n,0,nz,&Bnew);CHKERRQ(ierr); 175 ierr = PetscFree(nz);CHKERRQ(ierr); 176 for ( i=0; i<m; i++ ) { 177 for ( j=Baij->i[i]+shift; j<Baij->i[i+1]+shift; j++ ) { 178 col = garray[Baij->j[ct]+shift]; 179 v = Baij->a[ct++]; 180 ierr = MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);CHKERRQ(ierr); 181 } 182 } 183 ierr = PetscFree(aij->garray);CHKERRQ(ierr); 184 aij->garray = 0; 185 PLogObjectMemory(A,-ec*sizeof(int)); 186 ierr = MatDestroy(B);CHKERRQ(ierr); 187 PLogObjectParent(A,Bnew); 188 aij->B = Bnew; 189 A->was_assembled = PETSC_FALSE; 190 PetscFunctionReturn(0); 191 } 192 193 194