1 /*$Id: mmbaij.c,v 1.46 2001/09/25 00:31:36 balay Exp $*/ 2 3 /* 4 Support for the parallel BAIJ matrix vector multiply 5 */ 6 #include "src/mat/impls/baij/mpi/mpibaij.h" 7 #include "src/vec/vecimpl.h" 8 9 EXTERN int MatSetValuesBlocked_SeqBAIJ(Mat,int,int*,int,int*,PetscScalar*,InsertMode); 10 11 #undef __FUNCT__ 12 #define __FUNCT__ "MatSetUpMultiply_MPIBAIJ" 13 int MatSetUpMultiply_MPIBAIJ(Mat mat) 14 { 15 Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)mat->data; 16 Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)(baij->B->data); 17 int Nbs = baij->Nbs,i,j,*indices,*aj = B->j,ierr,ec = 0,*garray; 18 int bs = baij->bs,*stmp; 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 */ 31 ierr = PetscTableCreate(B->mbs,&gid1_lid1);CHKERRQ(ierr); 32 for (i=0; i<B->mbs; i++) { 33 for (j=0; j<B->ilen[i]; j++) { 34 int data,gid1 = aj[B->i[i]+j] + 1; 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 ierr = PetscMalloc((ec+1)*sizeof(int),&garray);CHKERRQ(ierr); 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); 51 /* qsort(garray, ec, sizeof(int), intcomparcarc); */ 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->mbs; i++) { 58 for (j=0; j<B->ilen[i]; j++) { 59 int gid1 = aj[B->i[i] + j] + 1; 60 ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr); 61 lid --; 62 aj[B->i[i]+j] = lid; 63 } 64 } 65 B->nbs = ec; 66 baij->B->n = ec*B->bs; 67 ierr = PetscTableDelete(gid1_lid1);CHKERRQ(ierr); 68 /* Mark Adams */ 69 #else 70 /* Make an array as long as the number of columns */ 71 /* mark those columns that are in baij->B */ 72 ierr = PetscMalloc((Nbs+1)*sizeof(int),&indices);CHKERRQ(ierr); 73 ierr = PetscMemzero(indices,Nbs*sizeof(int));CHKERRQ(ierr); 74 for (i=0; i<B->mbs; i++) { 75 for (j=0; j<B->ilen[i]; j++) { 76 if (!indices[aj[B->i[i] + j]]) ec++; 77 indices[aj[B->i[i] + j]] = 1; 78 } 79 } 80 81 /* form array of columns we need */ 82 ierr = PetscMalloc((ec+1)*sizeof(int),&garray);CHKERRQ(ierr); 83 ec = 0; 84 for (i=0; i<Nbs; i++) { 85 if (indices[i]) { 86 garray[ec++] = i; 87 } 88 } 89 90 /* make indices now point into garray */ 91 for (i=0; i<ec; i++) { 92 indices[garray[i]] = i; 93 } 94 95 /* compact out the extra columns in B */ 96 for (i=0; i<B->mbs; i++) { 97 for (j=0; j<B->ilen[i]; j++) { 98 aj[B->i[i] + j] = indices[aj[B->i[i] + j]]; 99 } 100 } 101 B->nbs = ec; 102 baij->B->n = ec*B->bs; 103 ierr = PetscFree(indices);CHKERRQ(ierr); 104 #endif 105 106 /* create local vector that is used to scatter into */ 107 ierr = VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);CHKERRQ(ierr); 108 109 /* create two temporary index sets for building scatter-gather */ 110 for (i=0; i<ec; i++) { 111 garray[i] = bs*garray[i]; 112 } 113 ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,&from);CHKERRQ(ierr); 114 for (i=0; i<ec; i++) { 115 garray[i] = garray[i]/bs; 116 } 117 118 ierr = PetscMalloc((ec+1)*sizeof(int),&stmp);CHKERRQ(ierr); 119 for (i=0; i<ec; i++) { stmp[i] = bs*i; } 120 ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,&to);CHKERRQ(ierr); 121 ierr = PetscFree(stmp);CHKERRQ(ierr); 122 123 /* create temporary global vector to generate scatter context */ 124 /* this is inefficient, but otherwise we must do either 125 1) save garray until the first actual scatter when the vector is known or 126 2) have another way of generating a scatter context without a vector.*/ 127 ierr = VecCreateMPI(mat->comm,mat->n,mat->N,&gvec);CHKERRQ(ierr); 128 129 /* gnerate the scatter context */ 130 ierr = VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);CHKERRQ(ierr); 131 132 /* 133 Post the receives for the first matrix vector product. We sync-chronize after 134 this on the chance that the user immediately calls MatMult() after assemblying 135 the matrix. 136 */ 137 ierr = VecScatterPostRecvs(gvec,baij->lvec,INSERT_VALUES,SCATTER_FORWARD,baij->Mvctx);CHKERRQ(ierr); 138 ierr = MPI_Barrier(mat->comm);CHKERRQ(ierr); 139 140 PetscLogObjectParent(mat,baij->Mvctx); 141 PetscLogObjectParent(mat,baij->lvec); 142 PetscLogObjectParent(mat,from); 143 PetscLogObjectParent(mat,to); 144 baij->garray = garray; 145 PetscLogObjectMemory(mat,(ec+1)*sizeof(int)); 146 ierr = ISDestroy(from);CHKERRQ(ierr); 147 ierr = ISDestroy(to);CHKERRQ(ierr); 148 ierr = VecDestroy(gvec);CHKERRQ(ierr); 149 PetscFunctionReturn(0); 150 } 151 152 /* 153 Takes the local part of an already assembled MPIBAIJ matrix 154 and disassembles it. This is to allow new nonzeros into the matrix 155 that require more communication in the matrix vector multiply. 156 Thus certain data-structures must be rebuilt. 157 158 Kind of slow! But that's what application programmers get when 159 they are sloppy. 160 */ 161 #undef __FUNCT__ 162 #define __FUNCT__ "DisAssemble_MPIBAIJ" 163 int DisAssemble_MPIBAIJ(Mat A) 164 { 165 Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)A->data; 166 Mat B = baij->B,Bnew; 167 Mat_SeqBAIJ *Bbaij = (Mat_SeqBAIJ*)B->data; 168 int ierr,i,j,mbs=Bbaij->mbs,n = A->N,col,*garray=baij->garray; 169 int bs2=baij->bs2,*nz,ec,m = A->m; 170 MatScalar *a = Bbaij->a; 171 PetscScalar *atmp; 172 173 PetscFunctionBegin; 174 /* free stuff related to matrix-vec multiply */ 175 ierr = VecGetSize(baij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */ 176 ierr = VecDestroy(baij->lvec);CHKERRQ(ierr); baij->lvec = 0; 177 ierr = VecScatterDestroy(baij->Mvctx);CHKERRQ(ierr); baij->Mvctx = 0; 178 if (baij->colmap) { 179 #if defined (PETSC_USE_CTABLE) 180 ierr = PetscTableDelete(baij->colmap); baij->colmap = 0;CHKERRQ(ierr); 181 #else 182 ierr = PetscFree(baij->colmap);CHKERRQ(ierr); 183 baij->colmap = 0; 184 PetscLogObjectMemory(A,-Bbaij->nbs*sizeof(int)); 185 #endif 186 } 187 188 /* make sure that B is assembled so we can access its values */ 189 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 190 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 191 192 /* invent new B and copy stuff over */ 193 ierr = PetscMalloc(mbs*sizeof(int),&nz);CHKERRQ(ierr); 194 for (i=0; i<mbs; i++) { 195 nz[i] = Bbaij->i[i+1]-Bbaij->i[i]; 196 } 197 ierr = MatCreateSeqBAIJ(PETSC_COMM_SELF,baij->bs,m,n,0,nz,&Bnew);CHKERRQ(ierr); 198 ierr = MatSetOption(Bnew,MAT_COLUMN_ORIENTED);CHKERRQ(ierr); 199 200 #if defined(PETSC_USE_MAT_SINGLE) 201 ierr = PetscMalloc(bs2*sizeof(PetscScalar),&atmp);CHKERRQ(ierr); 202 #endif 203 for (i=0; i<mbs; i++) { 204 for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) { 205 col = garray[Bbaij->j[j]]; 206 #if defined(PETSC_USE_MAT_SINGLE) 207 for (k=0; k<bs2; k++) atmp[k] = a[j*bs2+k]; 208 #else 209 atmp = a + j*bs2; 210 #endif 211 ierr = MatSetValuesBlocked_SeqBAIJ(Bnew,1,&i,1,&col,atmp,B->insertmode);CHKERRQ(ierr); 212 } 213 } 214 ierr = MatSetOption(Bnew,MAT_ROW_ORIENTED);CHKERRQ(ierr); 215 216 #if defined(PETSC_USE_MAT_SINGLE) 217 ierr = PetscFree(atmp);CHKERRQ(ierr); 218 #endif 219 220 ierr = PetscFree(nz);CHKERRQ(ierr); 221 ierr = PetscFree(baij->garray);CHKERRQ(ierr); 222 baij->garray = 0; 223 PetscLogObjectMemory(A,-ec*sizeof(int)); 224 ierr = MatDestroy(B);CHKERRQ(ierr); 225 PetscLogObjectParent(A,Bnew); 226 baij->B = Bnew; 227 A->was_assembled = PETSC_FALSE; 228 PetscFunctionReturn(0); 229 } 230 231 232 233