1 /* 2 Support for the parallel SELL matrix vector multiply 3 */ 4 #include <../src/mat/impls/sell/mpi/mpisell.h> 5 #include <petsc/private/isimpl.h> /* needed because accesses data structure of ISLocalToGlobalMapping directly */ 6 7 /* 8 Takes the local part of an already assembled MPISELL matrix 9 and disassembles it. This is to allow new nonzeros into the matrix 10 that require more communication in the matrix vector multiply. 11 Thus certain data-structures must be rebuilt. 12 13 Kind of slow! But that's what application programmers get when 14 they are sloppy. 15 */ 16 PetscErrorCode MatDisAssemble_MPISELL(Mat A) 17 { 18 Mat_MPISELL *sell=(Mat_MPISELL*)A->data; 19 Mat B=sell->B,Bnew; 20 Mat_SeqSELL *Bsell=(Mat_SeqSELL*)B->data; 21 PetscInt i,j,totalslices,N=A->cmap->N,ec,row; 22 PetscBool isnonzero; 23 PetscErrorCode ierr; 24 25 PetscFunctionBegin; 26 /* free stuff related to matrix-vec multiply */ 27 ierr = VecGetSize(sell->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */ 28 ierr = VecDestroy(&sell->lvec);CHKERRQ(ierr); 29 ierr = VecScatterDestroy(&sell->Mvctx);CHKERRQ(ierr); 30 if (sell->colmap) { 31 #if defined(PETSC_USE_CTABLE) 32 ierr = PetscTableDestroy(&sell->colmap);CHKERRQ(ierr); 33 #else 34 ierr = PetscFree(sell->colmap);CHKERRQ(ierr); 35 ierr = PetscLogObjectMemory((PetscObject)A,-sell->B->cmap->n*sizeof(PetscInt));CHKERRQ(ierr); 36 #endif 37 } 38 39 /* make sure that B is assembled so we can access its values */ 40 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 41 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 42 43 /* invent new B and copy stuff over */ 44 ierr = MatCreate(PETSC_COMM_SELF,&Bnew);CHKERRQ(ierr); 45 ierr = MatSetSizes(Bnew,B->rmap->n,N,B->rmap->n,N);CHKERRQ(ierr); 46 ierr = MatSetBlockSizesFromMats(Bnew,A,A);CHKERRQ(ierr); 47 ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr); 48 ierr = MatSeqSELLSetPreallocation(Bnew,0,Bsell->rlen);CHKERRQ(ierr); 49 if (Bsell->nonew >= 0) { /* Inherit insertion error options (if positive). */ 50 ((Mat_SeqSELL*)Bnew->data)->nonew = Bsell->nonew; 51 } 52 53 /* 54 Ensure that B's nonzerostate is monotonically increasing. 55 Or should this follow the MatSetValues() loop to preserve B's nonzerstate across a MatDisAssemble() call? 56 */ 57 Bnew->nonzerostate = B->nonzerostate; 58 59 totalslices = B->rmap->n/8+((B->rmap->n & 0x07)?1:0); /* floor(n/8) */ 60 for (i=0; i<totalslices; i++) { /* loop over slices */ 61 for (j=Bsell->sliidx[i],row=0; j<Bsell->sliidx[i+1]; j++,row=((row+1)&0x07)) { 62 isnonzero = (PetscBool)((j-Bsell->sliidx[i])/8 < Bsell->rlen[8*i+row]); 63 if (isnonzero) { 64 ierr = MatSetValue(Bnew,8*i+row,sell->garray[Bsell->colidx[j]],Bsell->val[j],B->insertmode);CHKERRQ(ierr); 65 } 66 } 67 } 68 69 ierr = PetscFree(sell->garray);CHKERRQ(ierr); 70 ierr = PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));CHKERRQ(ierr); 71 ierr = MatDestroy(&B);CHKERRQ(ierr); 72 ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);CHKERRQ(ierr); 73 74 sell->B = Bnew; 75 A->was_assembled = PETSC_FALSE; 76 A->assembled = PETSC_FALSE; 77 PetscFunctionReturn(0); 78 } 79 80 PetscErrorCode MatSetUpMultiply_MPISELL(Mat mat) 81 { 82 Mat_MPISELL *sell=(Mat_MPISELL*)mat->data; 83 Mat_SeqSELL *B=(Mat_SeqSELL*)(sell->B->data); 84 PetscErrorCode ierr; 85 PetscInt i,j,*bcolidx=B->colidx,ec=0,*garray,totalslices; 86 IS from,to; 87 Vec gvec; 88 PetscBool isnonzero; 89 #if defined(PETSC_USE_CTABLE) 90 PetscTable gid1_lid1; 91 PetscTablePosition tpos; 92 PetscInt gid,lid; 93 #else 94 PetscInt N = mat->cmap->N,*indices; 95 #endif 96 97 PetscFunctionBegin; 98 totalslices = sell->B->rmap->n/8+((sell->B->rmap->n & 0x07)?1:0); /* floor(n/8) */ 99 100 /* ec counts the number of columns that contain nonzeros */ 101 #if defined(PETSC_USE_CTABLE) 102 /* use a table */ 103 ierr = PetscTableCreate(sell->B->rmap->n,mat->cmap->N+1,&gid1_lid1);CHKERRQ(ierr); 104 for (i=0; i<totalslices; i++) { /* loop over slices */ 105 for (j=B->sliidx[i]; j<B->sliidx[i+1]; j++) { 106 isnonzero = (PetscBool)((j-B->sliidx[i])/8 < B->rlen[(i<<3)+(j&0x07)]); 107 if (isnonzero) { /* check the mask bit */ 108 PetscInt data,gid1 = bcolidx[j] + 1; 109 ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr); 110 if (!data) { 111 /* one based table */ 112 ierr = PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);CHKERRQ(ierr); 113 } 114 } 115 } 116 } 117 118 /* form array of columns we need */ 119 ierr = PetscMalloc1(ec+1,&garray);CHKERRQ(ierr); 120 ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr); 121 while (tpos) { 122 ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr); 123 gid--; 124 lid--; 125 garray[lid] = gid; 126 } 127 ierr = PetscSortInt(ec,garray);CHKERRQ(ierr); /* sort, and rebuild */ 128 ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr); 129 for (i=0; i<ec; i++) { 130 ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);CHKERRQ(ierr); 131 } 132 133 /* compact out the extra columns in B */ 134 for (i=0; i<totalslices; i++) { /* loop over slices */ 135 for (j=B->sliidx[i]; j<B->sliidx[i+1]; j++) { 136 isnonzero = (PetscBool)((j-B->sliidx[i])/8 < B->rlen[(i<<3)+(j&0x07)]); 137 if (isnonzero) { 138 PetscInt gid1 = bcolidx[j] + 1; 139 ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr); 140 lid--; 141 bcolidx[j] = lid; 142 } 143 } 144 } 145 sell->B->cmap->n = sell->B->cmap->N = ec; 146 sell->B->cmap->bs = 1; 147 148 ierr = PetscLayoutSetUp((sell->B->cmap));CHKERRQ(ierr); 149 ierr = PetscTableDestroy(&gid1_lid1);CHKERRQ(ierr); 150 #else 151 /* Make an array as long as the number of columns */ 152 ierr = PetscCalloc1(N+1,&indices);CHKERRQ(ierr); 153 /* mark those columns that are in sell->B */ 154 for (i=0; i<totalslices; i++) { /* loop over slices */ 155 for (j=B->sliidx[i]; j<B->sliidx[i+1]; j++) { 156 isnonzero = (PetscBool)((j-B->sliidx[i])/8 < B->rlen[(i<<3)+(j&0x07)]); 157 if (isnonzero) { 158 if (!indices[bcolidx[j]]) ec++; 159 indices[bcolidx[j]] = 1; 160 } 161 } 162 } 163 164 /* form array of columns we need */ 165 ierr = PetscMalloc1(ec+1,&garray);CHKERRQ(ierr); 166 ec = 0; 167 for (i=0; i<N; i++) { 168 if (indices[i]) garray[ec++] = i; 169 } 170 171 /* make indices now point into garray */ 172 for (i=0; i<ec; i++) { 173 indices[garray[i]] = i; 174 } 175 176 /* compact out the extra columns in B */ 177 for (i=0; i<totalslices; i++) { /* loop over slices */ 178 for (j=B->sliidx[i]; j<B->sliidx[i+1]; j++) { 179 isnonzero = (PetscBool)((j-B->sliidx[i])/8 < B->rlen[(i<<3)+(j&0x07)]); 180 if (isnonzero) bcolidx[j] = indices[bcolidx[j]]; 181 } 182 } 183 sell->B->cmap->n = sell->B->cmap->N = ec; /* number of columns that are not all zeros */ 184 sell->B->cmap->bs = 1; 185 186 ierr = PetscLayoutSetUp((sell->B->cmap));CHKERRQ(ierr); 187 ierr = PetscFree(indices);CHKERRQ(ierr); 188 #endif 189 /* create local vector that is used to scatter into */ 190 ierr = VecCreateSeq(PETSC_COMM_SELF,ec,&sell->lvec);CHKERRQ(ierr); 191 /* create two temporary Index sets for build scatter gather */ 192 ierr = ISCreateGeneral(PETSC_COMM_SELF,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 193 ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr); 194 195 /* create temporary global vector to generate scatter context */ 196 /* This does not allocate the array's memory so is efficient */ 197 ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);CHKERRQ(ierr); 198 199 /* generate the scatter context */ 200 ierr = VecScatterCreate(gvec,from,sell->lvec,to,&sell->Mvctx);CHKERRQ(ierr); 201 ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sell->Mvctx);CHKERRQ(ierr); 202 ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sell->lvec);CHKERRQ(ierr); 203 ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr); 204 ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr); 205 206 sell->garray = garray; 207 208 ierr = PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr); 209 ierr = ISDestroy(&from);CHKERRQ(ierr); 210 ierr = ISDestroy(&to);CHKERRQ(ierr); 211 ierr = VecDestroy(&gvec);CHKERRQ(ierr); 212 PetscFunctionReturn(0); 213 } 214 215 /* ugly stuff added for Glenn someday we should fix this up */ 216 static PetscInt *auglyrmapd = 0,*auglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */ 217 static Vec auglydd = 0,auglyoo = 0; /* work vectors used to scale the two parts of the local matrix */ 218 219 PetscErrorCode MatMPISELLDiagonalScaleLocalSetUp(Mat inA,Vec scale) 220 { 221 Mat_MPISELL *ina=(Mat_MPISELL*)inA->data; /*access private part of matrix */ 222 PetscErrorCode ierr; 223 PetscInt i,n,nt,cstart,cend,no,*garray=ina->garray,*lindices; 224 PetscInt *r_rmapd,*r_rmapo; 225 226 PetscFunctionBegin; 227 ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr); 228 ierr = MatGetSize(ina->A,NULL,&n);CHKERRQ(ierr); 229 ierr = PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapd);CHKERRQ(ierr); 230 nt = 0; 231 for (i=0; i<inA->rmap->mapping->n; i++) { 232 if (inA->rmap->mapping->indices[i] >= cstart && inA->rmap->mapping->indices[i] < cend) { 233 nt++; 234 r_rmapd[i] = inA->rmap->mapping->indices[i] + 1; 235 } 236 } 237 if (nt != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D n %D",nt,n); 238 ierr = PetscMalloc1(n+1,&auglyrmapd);CHKERRQ(ierr); 239 for (i=0; i<inA->rmap->mapping->n; i++) { 240 if (r_rmapd[i]) { 241 auglyrmapd[(r_rmapd[i]-1)-cstart] = i; 242 } 243 } 244 ierr = PetscFree(r_rmapd);CHKERRQ(ierr); 245 ierr = VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);CHKERRQ(ierr); 246 ierr = PetscCalloc1(inA->cmap->N+1,&lindices);CHKERRQ(ierr); 247 for (i=0; i<ina->B->cmap->n; i++) { 248 lindices[garray[i]] = i+1; 249 } 250 no = inA->rmap->mapping->n - nt; 251 ierr = PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapo);CHKERRQ(ierr); 252 nt = 0; 253 for (i=0; i<inA->rmap->mapping->n; i++) { 254 if (lindices[inA->rmap->mapping->indices[i]]) { 255 nt++; 256 r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]]; 257 } 258 } 259 if (nt > no) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n); 260 ierr = PetscFree(lindices);CHKERRQ(ierr); 261 ierr = PetscMalloc1(nt+1,&auglyrmapo);CHKERRQ(ierr); 262 for (i=0; i<inA->rmap->mapping->n; i++) { 263 if (r_rmapo[i]) { 264 auglyrmapo[(r_rmapo[i]-1)] = i; 265 } 266 } 267 ierr = PetscFree(r_rmapo);CHKERRQ(ierr); 268 ierr = VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);CHKERRQ(ierr); 269 PetscFunctionReturn(0); 270 } 271 272 PetscErrorCode MatDiagonalScaleLocal_MPISELL(Mat A,Vec scale) 273 { 274 Mat_MPISELL *a=(Mat_MPISELL*)A->data; /*access private part of matrix */ 275 PetscErrorCode ierr; 276 PetscInt n,i; 277 PetscScalar *d,*o; 278 const PetscScalar *s; 279 280 PetscFunctionBegin; 281 if (!auglyrmapd) { 282 ierr = MatMPISELLDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr); 283 } 284 ierr = VecGetArrayRead(scale,&s);CHKERRQ(ierr); 285 ierr = VecGetLocalSize(auglydd,&n);CHKERRQ(ierr); 286 ierr = VecGetArray(auglydd,&d);CHKERRQ(ierr); 287 for (i=0; i<n; i++) { 288 d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */ 289 } 290 ierr = VecRestoreArray(auglydd,&d);CHKERRQ(ierr); 291 /* column scale "diagonal" portion of local matrix */ 292 ierr = MatDiagonalScale(a->A,NULL,auglydd);CHKERRQ(ierr); 293 ierr = VecGetLocalSize(auglyoo,&n);CHKERRQ(ierr); 294 ierr = VecGetArray(auglyoo,&o);CHKERRQ(ierr); 295 for (i=0; i<n; i++) { 296 o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */ 297 } 298 ierr = VecRestoreArrayRead(scale,&s);CHKERRQ(ierr); 299 ierr = VecRestoreArray(auglyoo,&o);CHKERRQ(ierr); 300 /* column scale "off-diagonal" portion of local matrix */ 301 ierr = MatDiagonalScale(a->B,NULL,auglyoo);CHKERRQ(ierr); 302 PetscFunctionReturn(0); 303 } 304