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