1 /*$Id: mmaij.c,v 1.59 2001/08/07 03:02:49 balay Exp $*/ 2 3 /* 4 Support for the parallel AIJ matrix vector multiply 5 */ 6 #include "src/mat/impls/aij/mpi/mpiaij.h" 7 #include "src/vec/vecimpl.h" 8 9 #undef __FUNCT__ 10 #define __FUNCT__ "MatSetUpMultiply_MPIAIJ" 11 int MatSetUpMultiply_MPIAIJ(Mat mat) 12 { 13 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)mat->data; 14 Mat_SeqAIJ *B = (Mat_SeqAIJ*)(aij->B->data); 15 int N = mat->N,i,j,*indices,*aj = B->j,ierr,ec = 0,*garray; 16 int shift = B->indexshift; 17 IS from,to; 18 Vec gvec; 19 #if defined (PETSC_USE_CTABLE) 20 PetscTable gid1_lid1; 21 PetscTablePosition tpos; 22 int gid,lid; 23 #endif 24 25 PetscFunctionBegin; 26 27 #if defined (PETSC_USE_CTABLE) 28 /* use a table - Mark Adams (this has not been tested with "shift") */ 29 ierr = PetscTableCreate(aij->B->m,&gid1_lid1);CHKERRQ(ierr); 30 for (i=0; i<aij->B->m; i++) { 31 for (j=0; j<B->ilen[i]; j++) { 32 int data,gid1 = aj[B->i[i] + shift + j] + 1 + shift; 33 ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr); 34 if (!data) { 35 /* one based table */ 36 ierr = PetscTableAdd(gid1_lid1,gid1,++ec);CHKERRQ(ierr); 37 } 38 } 39 } 40 /* form array of columns we need */ 41 ierr = PetscMalloc((ec+1)*sizeof(int),&garray);CHKERRQ(ierr); 42 ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr); 43 while (tpos) { 44 ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr); 45 gid--; 46 lid--; 47 garray[lid] = gid; 48 } 49 ierr = PetscSortInt(ec,garray);CHKERRQ(ierr); /* sort, and rebuild */ 50 ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr); 51 for (i=0; i<ec; i++) { 52 ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1);CHKERRQ(ierr); 53 } 54 /* compact out the extra columns in B */ 55 for (i=0; i<aij->B->m; i++) { 56 for (j=0; j<B->ilen[i]; j++) { 57 int gid1 = aj[B->i[i] + shift + j] + 1 + shift; 58 ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr); 59 lid --; 60 aj[B->i[i] + shift + j] = lid - shift; 61 } 62 } 63 aij->B->n = aij->B->N = ec; 64 ierr = PetscTableDelete(gid1_lid1);CHKERRQ(ierr); 65 /* Mark Adams */ 66 #else 67 /* For the first stab we make an array as long as the number of columns */ 68 /* mark those columns that are in aij->B */ 69 ierr = PetscMalloc((N+1)*sizeof(int),&indices);CHKERRQ(ierr); 70 ierr = PetscMemzero(indices,N*sizeof(int));CHKERRQ(ierr); 71 for (i=0; i<aij->B->m; i++) { 72 for (j=0; j<B->ilen[i]; j++) { 73 if (!indices[aj[B->i[i] +shift + j] + shift]) ec++; 74 indices[aj[B->i[i] + shift + j] + shift] = 1; 75 } 76 } 77 78 /* form array of columns we need */ 79 ierr = PetscMalloc((ec+1)*sizeof(int),&garray);CHKERRQ(ierr); 80 ec = 0; 81 for (i=0; i<N; i++) { 82 if (indices[i]) garray[ec++] = i; 83 } 84 85 /* make indices now point into garray */ 86 for (i=0; i<ec; i++) { 87 indices[garray[i]] = i-shift; 88 } 89 90 /* compact out the extra columns in B */ 91 for (i=0; i<aij->B->m; i++) { 92 for (j=0; j<B->ilen[i]; j++) { 93 aj[B->i[i] + shift + j] = indices[aj[B->i[i] + shift + j]+shift]; 94 } 95 } 96 aij->B->n = aij->B->N = ec; 97 ierr = PetscFree(indices);CHKERRQ(ierr); 98 #endif 99 /* create local vector that is used to scatter into */ 100 ierr = VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);CHKERRQ(ierr); 101 102 /* create two temporary Index sets for build scatter gather */ 103 ierr = ISCreateGeneral(mat->comm,ec,garray,&from);CHKERRQ(ierr); 104 ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr); 105 106 /* create temporary global vector to generate scatter context */ 107 /* this is inefficient, but otherwise we must do either 108 1) save garray until the first actual scatter when the vector is known or 109 2) have another way of generating a scatter context without a vector.*/ 110 ierr = VecCreateMPI(mat->comm,mat->n,mat->N,&gvec);CHKERRQ(ierr); 111 112 /* generate the scatter context */ 113 ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);CHKERRQ(ierr); 114 PetscLogObjectParent(mat,aij->Mvctx); 115 PetscLogObjectParent(mat,aij->lvec); 116 PetscLogObjectParent(mat,from); 117 PetscLogObjectParent(mat,to); 118 aij->garray = garray; 119 PetscLogObjectMemory(mat,(ec+1)*sizeof(int)); 120 ierr = ISDestroy(from);CHKERRQ(ierr); 121 ierr = ISDestroy(to);CHKERRQ(ierr); 122 ierr = VecDestroy(gvec);CHKERRQ(ierr); 123 PetscFunctionReturn(0); 124 } 125 126 127 #undef __FUNCT__ 128 #define __FUNCT__ "DisAssemble_MPIAIJ" 129 /* 130 Takes the local part of an already assembled MPIAIJ matrix 131 and disassembles it. This is to allow new nonzeros into the matrix 132 that require more communication in the matrix vector multiply. 133 Thus certain data-structures must be rebuilt. 134 135 Kind of slow! But that's what application programmers get when 136 they are sloppy. 137 */ 138 int DisAssemble_MPIAIJ(Mat A) 139 { 140 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 141 Mat B = aij->B,Bnew; 142 Mat_SeqAIJ *Baij = (Mat_SeqAIJ*)B->data; 143 int ierr,i,j,m = B->m,n = A->N,col,ct = 0,*garray = aij->garray; 144 int *nz,ec,shift = Baij->indexshift; 145 PetscScalar v; 146 147 PetscFunctionBegin; 148 /* free stuff related to matrix-vec multiply */ 149 ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */ 150 ierr = VecDestroy(aij->lvec);CHKERRQ(ierr); aij->lvec = 0; 151 ierr = VecScatterDestroy(aij->Mvctx);CHKERRQ(ierr); aij->Mvctx = 0; 152 if (aij->colmap) { 153 #if defined (PETSC_USE_CTABLE) 154 ierr = PetscTableDelete(aij->colmap);CHKERRQ(ierr); 155 aij->colmap = 0; 156 #else 157 ierr = PetscFree(aij->colmap);CHKERRQ(ierr); 158 aij->colmap = 0; 159 PetscLogObjectMemory(A,-aij->B->n*sizeof(int)); 160 #endif 161 } 162 163 /* make sure that B is assembled so we can access its values */ 164 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 165 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 166 167 /* invent new B and copy stuff over */ 168 ierr = PetscMalloc((m+1)*sizeof(int),&nz);CHKERRQ(ierr); 169 for (i=0; i<m; i++) { 170 nz[i] = Baij->i[i+1] - Baij->i[i]; 171 } 172 ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m,n,0,nz,&Bnew);CHKERRQ(ierr); 173 ierr = PetscFree(nz);CHKERRQ(ierr); 174 for (i=0; i<m; i++) { 175 for (j=Baij->i[i]+shift; j<Baij->i[i+1]+shift; j++) { 176 col = garray[Baij->j[ct]+shift]; 177 v = Baij->a[ct++]; 178 ierr = MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);CHKERRQ(ierr); 179 } 180 } 181 ierr = PetscFree(aij->garray);CHKERRQ(ierr); 182 aij->garray = 0; 183 PetscLogObjectMemory(A,-ec*sizeof(int)); 184 ierr = MatDestroy(B);CHKERRQ(ierr); 185 PetscLogObjectParent(A,Bnew); 186 aij->B = Bnew; 187 A->was_assembled = PETSC_FALSE; 188 PetscFunctionReturn(0); 189 } 190 191 /* ugly stuff added for Glenn someday we should fix this up */ 192 193 static int *auglyrmapd = 0,*auglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" 194 parts of the local matrix */ 195 static Vec auglydd = 0,auglyoo = 0; /* work vectors used to scale the two parts of the local matrix */ 196 197 198 #undef __FUNCT__ 199 #define __FUNCT__ "MatMPIAIJDiagonalScaleLocalSetUp" 200 int MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale) 201 { 202 Mat_MPIAIJ *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */ 203 Mat_SeqAIJ *A = (Mat_SeqAIJ*)ina->A->data; 204 Mat_SeqAIJ *B = (Mat_SeqAIJ*)ina->B->data; 205 int ierr,i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices; 206 int *r_rmapd,*r_rmapo; 207 208 PetscFunctionBegin; 209 ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr); 210 ierr = MatGetSize(ina->A,PETSC_NULL,&n);CHKERRQ(ierr); 211 ierr = PetscMalloc((inA->mapping->n+1)*sizeof(int),&r_rmapd);CHKERRQ(ierr); 212 ierr = PetscMemzero(r_rmapd,inA->mapping->n*sizeof(int));CHKERRQ(ierr); 213 nt = 0; 214 for (i=0; i<inA->mapping->n; i++) { 215 if (inA->mapping->indices[i] >= cstart && inA->mapping->indices[i] < cend) { 216 nt++; 217 r_rmapd[i] = inA->mapping->indices[i] + 1; 218 } 219 } 220 if (nt != n) SETERRQ2(1,"Hmm nt %d n %d",nt,n); 221 ierr = PetscMalloc((n+1)*sizeof(int),&auglyrmapd);CHKERRQ(ierr); 222 for (i=0; i<inA->mapping->n; i++) { 223 if (r_rmapd[i]){ 224 auglyrmapd[(r_rmapd[i]-1)-cstart] = i; 225 } 226 } 227 ierr = PetscFree(r_rmapd);CHKERRQ(ierr); 228 ierr = VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);CHKERRQ(ierr); 229 230 ierr = PetscMalloc((inA->N+1)*sizeof(int),&lindices);CHKERRQ(ierr); 231 ierr = PetscMemzero(lindices,inA->N*sizeof(int));CHKERRQ(ierr); 232 for (i=0; i<ina->B->n; i++) { 233 lindices[garray[i]] = i+1; 234 } 235 no = inA->mapping->n - nt; 236 ierr = PetscMalloc((inA->mapping->n+1)*sizeof(int),&r_rmapo);CHKERRQ(ierr); 237 ierr = PetscMemzero(r_rmapo,inA->mapping->n*sizeof(int));CHKERRQ(ierr); 238 nt = 0; 239 for (i=0; i<inA->mapping->n; i++) { 240 if (lindices[inA->mapping->indices[i]]) { 241 nt++; 242 r_rmapo[i] = lindices[inA->mapping->indices[i]]; 243 } 244 } 245 if (nt > no) SETERRQ2(1,"Hmm nt %d no %d",nt,n); 246 ierr = PetscFree(lindices);CHKERRQ(ierr); 247 ierr = PetscMalloc((nt+1)*sizeof(int),&auglyrmapo);CHKERRQ(ierr); 248 for (i=0; i<inA->mapping->n; i++) { 249 if (r_rmapo[i]){ 250 auglyrmapo[(r_rmapo[i]-1)] = i; 251 } 252 } 253 ierr = PetscFree(r_rmapo);CHKERRQ(ierr); 254 ierr = VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);CHKERRQ(ierr); 255 256 PetscFunctionReturn(0); 257 } 258 259 #undef __FUNCT__ 260 #define __FUNCT__ "MatMPIAIJDiagonalScaleLocal" 261 int MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale) 262 { 263 Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */ 264 int ierr,n,i; 265 PetscScalar *d,*o,*s; 266 267 PetscFunctionBegin; 268 if (!auglyrmapd) { 269 ierr = MatMPIAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr); 270 } 271 272 ierr = VecGetArray(scale,&s);CHKERRQ(ierr); 273 274 ierr = VecGetLocalSize(auglydd,&n);CHKERRQ(ierr); 275 ierr = VecGetArray(auglydd,&d);CHKERRQ(ierr); 276 for (i=0; i<n; i++) { 277 d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */ 278 } 279 ierr = VecRestoreArray(auglydd,&d);CHKERRQ(ierr); 280 /* column scale "diagonal" portion of local matrix */ 281 ierr = MatDiagonalScale(a->A,PETSC_NULL,auglydd);CHKERRQ(ierr); 282 283 ierr = VecGetLocalSize(auglyoo,&n);CHKERRQ(ierr); 284 ierr = VecGetArray(auglyoo,&o);CHKERRQ(ierr); 285 for (i=0; i<n; i++) { 286 o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */ 287 } 288 ierr = VecRestoreArray(scale,&s);CHKERRQ(ierr); 289 ierr = VecRestoreArray(auglyoo,&o);CHKERRQ(ierr); 290 /* column scale "off-diagonal" portion of local matrix */ 291 ierr = MatDiagonalScale(a->B,PETSC_NULL,auglyoo);CHKERRQ(ierr); 292 293 PetscFunctionReturn(0); 294 } 295 296 297 298 299