1 /*$Id: mmbaij.c,v 1.34 2000/07/10 03:39:45 bsmith Exp bsmith $*/ 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 EXTERN int MatSetValues_SeqBAIJ(Mat,int,int*,int,int*,Scalar*,InsertMode); 9 10 #undef __FUNC__ 11 #define __FUNC__ /*<a name=""></a>*/"MatSetUpMultiply_MPIBAIJ" 12 int MatSetUpMultiply_MPIBAIJ(Mat mat) 13 { 14 Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)mat->data; 15 Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)(baij->B->data); 16 int Nbs = baij->Nbs,i,j,*indices,*aj = B->j,ierr,ec = 0,*garray; 17 int bs = baij->bs,*stmp; 18 IS from,to; 19 Vec gvec; 20 #if defined (PETSC_USE_CTABLE) 21 PetscTable gid1_lid1; 22 PetscTablePosition tpos; 23 int gid,lid; 24 #endif 25 26 PetscFunctionBegin; 27 28 #if defined (PETSC_USE_CTABLE) 29 /* use a table - Mark Adams */ 30 PetscTableCreate(B->mbs,&gid1_lid1); 31 for (i=0; i<B->mbs; i++) { 32 for (j=0; j<B->ilen[i]; j++) { 33 int data,gid1 = aj[B->i[i]+j] + 1; 34 ierr = PetscTableFind(gid1_lid1,gid1,&data) ;CHKERRQ(ierr); 35 if (!data) { 36 /* one based table */ 37 ierr = PetscTableAdd(gid1_lid1,gid1,++ec);CHKERRQ(ierr); 38 } 39 } 40 } 41 /* form array of columns we need */ 42 garray = (int *)PetscMalloc((ec+1)*sizeof(int));CHKPTRQ(garray); 43 ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr); 44 while (tpos) { 45 ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr); 46 gid--; lid--; 47 garray[lid] = gid; 48 } 49 ierr = PetscSortInt(ec,garray);CHKERRQ(ierr); 50 /* qsort(garray, ec, sizeof(int), intcomparcarc); */ 51 ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr); 52 for (i=0; i<ec; i++) { 53 ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1);CHKERRQ(ierr); 54 } 55 /* compact out the extra columns in B */ 56 for (i=0; i<B->mbs; i++) { 57 for (j=0; j<B->ilen[i]; j++) { 58 int gid1 = aj[B->i[i] + j] + 1; 59 ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr); 60 lid --; 61 aj[B->i[i]+j] = lid; 62 } 63 } 64 B->nbs = ec; 65 B->n = ec*B->bs; 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 baij->B */ 71 indices = (int*)PetscMalloc((Nbs+1)*sizeof(int));CHKPTRQ(indices); 72 ierr = PetscMemzero(indices,Nbs*sizeof(int));CHKERRQ(ierr); 73 for (i=0; i<B->mbs; i++) { 74 for (j=0; j<B->ilen[i]; j++) { 75 if (!indices[aj[B->i[i] + j]]) ec++; 76 indices[aj[B->i[i] + j]] = 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<Nbs; i++) { 84 if (indices[i]) { 85 garray[ec++] = i; 86 } 87 } 88 89 /* make indices now point into garray */ 90 for (i=0; i<ec; i++) { 91 indices[garray[i]] = i; 92 } 93 94 /* compact out the extra columns in B */ 95 for (i=0; i<B->mbs; i++) { 96 for (j=0; j<B->ilen[i]; j++) { 97 aj[B->i[i] + j] = indices[aj[B->i[i] + j]]; 98 } 99 } 100 B->nbs = ec; 101 B->n = ec*B->bs; 102 ierr = PetscFree(indices);CHKERRQ(ierr); 103 #endif 104 105 /* create local vector that is used to scatter into */ 106 ierr = VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);CHKERRQ(ierr); 107 108 /* create two temporary index sets for building scatter-gather */ 109 for (i=0; i<ec; i++) { 110 garray[i] = bs*garray[i]; 111 } 112 ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,&from);CHKERRQ(ierr); 113 for (i=0; i<ec; i++) { 114 garray[i] = garray[i]/bs; 115 } 116 117 stmp = (int*)PetscMalloc((ec+1)*sizeof(int));CHKPTRQ(stmp); 118 for (i=0; i<ec; i++) { stmp[i] = bs*i; } 119 ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,&to);CHKERRQ(ierr); 120 ierr = PetscFree(stmp);CHKERRQ(ierr); 121 122 /* create temporary global vector to generate scatter context */ 123 /* this is inefficient, but otherwise we must do either 124 1) save garray until the first actual scatter when the vector is known or 125 2) have another way of generating a scatter context without a vector.*/ 126 ierr = VecCreateMPI(mat->comm,baij->n,baij->N,&gvec);CHKERRQ(ierr); 127 128 /* gnerate the scatter context */ 129 ierr = VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);CHKERRQ(ierr); 130 131 /* 132 Post the receives for the first matrix vector product. We sync-chronize after 133 this on the chance that the user immediately calls MatMult() after assemblying 134 the matrix. 135 */ 136 ierr = VecScatterPostRecvs(gvec,baij->lvec,INSERT_VALUES,SCATTER_FORWARD,baij->Mvctx);CHKERRQ(ierr); 137 ierr = MPI_Barrier(mat->comm);CHKERRQ(ierr); 138 139 PLogObjectParent(mat,baij->Mvctx); 140 PLogObjectParent(mat,baij->lvec); 141 PLogObjectParent(mat,from); 142 PLogObjectParent(mat,to); 143 baij->garray = garray; 144 PLogObjectMemory(mat,(ec+1)*sizeof(int)); 145 ierr = ISDestroy(from);CHKERRQ(ierr); 146 ierr = ISDestroy(to);CHKERRQ(ierr); 147 ierr = VecDestroy(gvec);CHKERRQ(ierr); 148 PetscFunctionReturn(0); 149 } 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 __FUNC__ 162 #define __FUNC__ /*<a name=""></a>*/"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 = baij->N,col,*garray=baij->garray; 169 int k,bs=baij->bs,bs2=baij->bs2,*rvals,*nz,ec,m=Bbaij->m; 170 MatScalar *a = Bbaij->a; 171 #if defined(PETSC_USE_MAT_SINGLE) 172 Scalar *atmp = (Scalar*)PetscMalloc(baij->bs*sizeof(Scalar)); 173 int l; 174 #else 175 Scalar *atmp; 176 #endif 177 178 PetscFunctionBegin; 179 /* free stuff related to matrix-vec multiply */ 180 ierr = VecGetSize(baij->lvec,&ec);CHKERRQ(ierr); /* needed for PLogObjectMemory below */ 181 ierr = VecDestroy(baij->lvec);CHKERRQ(ierr); baij->lvec = 0; 182 ierr = VecScatterDestroy(baij->Mvctx);CHKERRQ(ierr); baij->Mvctx = 0; 183 if (baij->colmap) { 184 #if defined (PETSC_USE_CTABLE) 185 ierr = PetscTableDelete(baij->colmap); baij->colmap = 0;CHKERRQ(ierr); 186 #else 187 ierr = PetscFree(baij->colmap);CHKERRQ(ierr); 188 baij->colmap = 0; 189 PLogObjectMemory(A,-Bbaij->nbs*sizeof(int)); 190 #endif 191 } 192 193 /* make sure that B is assembled so we can access its values */ 194 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 195 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 196 197 /* invent new B and copy stuff over */ 198 nz = (int*)PetscMalloc(mbs*sizeof(int));CHKPTRQ(nz); 199 for (i=0; i<mbs; i++) { 200 nz[i] = Bbaij->i[i+1]-Bbaij->i[i]; 201 } 202 ierr = MatCreateSeqBAIJ(PETSC_COMM_SELF,baij->bs,m,n,0,nz,&Bnew);CHKERRQ(ierr); 203 ierr = PetscFree(nz);CHKERRQ(ierr); 204 205 rvals = (int*)PetscMalloc(bs*sizeof(int));CHKPTRQ(rvals); 206 for (i=0; i<mbs; i++) { 207 rvals[0] = bs*i; 208 for (j=1; j<bs; j++) { rvals[j] = rvals[j-1] + 1; } 209 for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) { 210 col = garray[Bbaij->j[j]]*bs; 211 for (k=0; k<bs; k++) { 212 #if defined(PETSC_USE_MAT_SINGLE) 213 for (l=0; l<bs; l++) atmp[l] = a[j*bs2+l]; 214 #else 215 atmp = a+j*bs2; 216 #endif 217 ierr = MatSetValues_SeqBAIJ(Bnew,bs,rvals,1,&col,atmp,B->insertmode);CHKERRQ(ierr); 218 col++; 219 } 220 } 221 } 222 #if defined(PETSC_USE_MAT_SINGLE) 223 ierr = PetscFree(atmp);CHKERRQ(ierr); 224 #endif 225 ierr = PetscFree(baij->garray);CHKERRQ(ierr); 226 baij->garray = 0; 227 ierr = PetscFree(rvals);CHKERRQ(ierr); 228 PLogObjectMemory(A,-ec*sizeof(int)); 229 ierr = MatDestroy(B);CHKERRQ(ierr); 230 PLogObjectParent(A,Bnew); 231 baij->B = Bnew; 232 A->was_assembled = PETSC_FALSE; 233 PetscFunctionReturn(0); 234 } 235 236 237 238