1 2 /* 3 Support for the parallel AIJ matrix vector multiply 4 */ 5 #include <../src/mat/impls/aij/mpi/mpiaij.h> 6 7 #undef __FUNCT__ 8 #define __FUNCT__ "MatSetUpMultiply_MPIAIJ" 9 PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat mat) 10 { 11 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)mat->data; 12 Mat_SeqAIJ *B = (Mat_SeqAIJ*)(aij->B->data); 13 PetscErrorCode ierr; 14 PetscInt i,j,*aj = B->j,ec = 0,*garray; 15 IS from,to; 16 Vec gvec; 17 PetscBool useblockis; 18 #if defined (PETSC_USE_CTABLE) 19 PetscTable gid1_lid1; 20 PetscTablePosition tpos; 21 PetscInt gid,lid; 22 #else 23 PetscInt N = mat->cmap->N,*indices; 24 #endif 25 26 PetscFunctionBegin; 27 28 #if defined (PETSC_USE_CTABLE) 29 /* use a table - Mark Adams */ 30 ierr = PetscTableCreate(aij->B->rmap->n,&gid1_lid1);CHKERRQ(ierr); 31 for (i=0; i<aij->B->rmap->n; i++) { 32 for (j=0; j<B->ilen[i]; j++) { 33 PetscInt 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 ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&garray);CHKERRQ(ierr); 43 ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr); 44 while (tpos) { 45 ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr); 46 gid--; 47 lid--; 48 garray[lid] = gid; 49 } 50 ierr = PetscSortInt(ec,garray);CHKERRQ(ierr); /* sort, and rebuild */ 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<aij->B->rmap->n; i++) { 57 for (j=0; j<B->ilen[i]; j++) { 58 PetscInt 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 aij->B->cmap->n = aij->B->cmap->N = ec; 65 ierr = PetscLayoutSetUp((aij->B->cmap));CHKERRQ(ierr); 66 ierr = PetscTableDestroy(&gid1_lid1);CHKERRQ(ierr); 67 /* Mark Adams */ 68 #else 69 /* Make an array as long as the number of columns */ 70 /* mark those columns that are in aij->B */ 71 ierr = PetscMalloc((N+1)*sizeof(PetscInt),&indices);CHKERRQ(ierr); 72 ierr = PetscMemzero(indices,N*sizeof(PetscInt));CHKERRQ(ierr); 73 for (i=0; i<aij->B->rmap->n; 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 ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&garray);CHKERRQ(ierr); 82 ec = 0; 83 for (i=0; i<N; i++) { 84 if (indices[i]) garray[ec++] = i; 85 } 86 87 /* make indices now point into garray */ 88 for (i=0; i<ec; i++) { 89 indices[garray[i]] = i; 90 } 91 92 /* compact out the extra columns in B */ 93 for (i=0; i<aij->B->rmap->n; i++) { 94 for (j=0; j<B->ilen[i]; j++) { 95 aj[B->i[i] + j] = indices[aj[B->i[i] + j]]; 96 } 97 } 98 aij->B->cmap->n = aij->B->cmap->N = ec; 99 ierr = PetscLayoutSetUp((aij->B->cmap));CHKERRQ(ierr); 100 ierr = PetscFree(indices);CHKERRQ(ierr); 101 #endif 102 /* create local vector that is used to scatter into */ 103 ierr = VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);CHKERRQ(ierr); 104 105 /* create two temporary Index sets for build scatter gather */ 106 /* check for the special case where blocks are communicated for faster VecScatterXXX */ 107 useblockis = PETSC_TRUE; 108 if (mat->rmap->bs > 1) { 109 PetscInt bs = mat->rmap->bs,ibs,ga; 110 if (!(ec % bs)) { 111 for (i=0; i<ec/bs; i++) { 112 if ((ga = garray[ibs = i*bs]) % bs) { 113 useblockis = PETSC_FALSE; 114 break; 115 } 116 for (j=1; j<bs; j++) { 117 if (garray[ibs+j] != ga+j) { 118 useblockis = PETSC_FALSE; 119 break; 120 } 121 } 122 if (!useblockis) break; 123 } 124 } 125 } 126 if (useblockis) { 127 PetscInt *ga,bs = mat->rmap->bs,iec = ec/bs; 128 ierr = PetscInfo(mat,"Using block index set to define scatter\n"); 129 ierr = PetscMalloc((ec/mat->rmap->bs)*sizeof(PetscInt),&ga);CHKERRQ(ierr); 130 for (i=0; i<iec; i++) ga[i] = garray[i*bs]/bs; 131 ierr = ISCreateBlock(((PetscObject)mat)->comm,bs,iec,ga,PETSC_OWN_POINTER,&from);CHKERRQ(ierr); 132 } else { 133 ierr = ISCreateGeneral(((PetscObject)mat)->comm,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 134 } 135 ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr); 136 137 /* create temporary global vector to generate scatter context */ 138 /* This does not allocate the array's memory so is efficient */ 139 ierr = VecCreateMPIWithArray(((PetscObject)mat)->comm,mat->cmap->n,mat->cmap->N,PETSC_NULL,&gvec);CHKERRQ(ierr); 140 141 /* generate the scatter context */ 142 ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);CHKERRQ(ierr); 143 ierr = PetscLogObjectParent(mat,aij->Mvctx);CHKERRQ(ierr); 144 ierr = PetscLogObjectParent(mat,aij->lvec);CHKERRQ(ierr); 145 ierr = PetscLogObjectParent(mat,from);CHKERRQ(ierr); 146 ierr = PetscLogObjectParent(mat,to);CHKERRQ(ierr); 147 aij->garray = garray; 148 ierr = PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr); 149 ierr = ISDestroy(&from);CHKERRQ(ierr); 150 ierr = ISDestroy(&to);CHKERRQ(ierr); 151 ierr = VecDestroy(&gvec);CHKERRQ(ierr); 152 PetscFunctionReturn(0); 153 } 154 155 156 #undef __FUNCT__ 157 #define __FUNCT__ "DisAssemble_MPIAIJ" 158 /* 159 Takes the local part of an already assembled MPIAIJ matrix 160 and disassembles it. This is to allow new nonzeros into the matrix 161 that require more communication in the matrix vector multiply. 162 Thus certain data-structures must be rebuilt. 163 164 Kind of slow! But that's what application programmers get when 165 they are sloppy. 166 */ 167 PetscErrorCode DisAssemble_MPIAIJ(Mat A) 168 { 169 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 170 Mat B = aij->B,Bnew; 171 Mat_SeqAIJ *Baij = (Mat_SeqAIJ*)B->data; 172 PetscErrorCode ierr; 173 PetscInt i,j,m = B->rmap->n,n = A->cmap->N,col,ct = 0,*garray = aij->garray,*nz,ec; 174 PetscScalar v; 175 176 PetscFunctionBegin; 177 /* free stuff related to matrix-vec multiply */ 178 ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */ 179 ierr = VecDestroy(&aij->lvec);CHKERRQ(ierr); aij->lvec = 0; 180 ierr = VecScatterDestroy(&aij->Mvctx);CHKERRQ(ierr); aij->Mvctx = 0; 181 if (aij->colmap) { 182 #if defined (PETSC_USE_CTABLE) 183 ierr = PetscTableDestroy(&aij->colmap);CHKERRQ(ierr); 184 #else 185 ierr = PetscFree(aij->colmap);CHKERRQ(ierr); 186 ierr = PetscLogObjectMemory(A,-aij->B->cmap->n*sizeof(PetscInt));CHKERRQ(ierr); 187 #endif 188 } 189 190 /* make sure that B is assembled so we can access its values */ 191 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 192 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 193 194 /* invent new B and copy stuff over */ 195 ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nz);CHKERRQ(ierr); 196 for (i=0; i<m; i++) { 197 nz[i] = Baij->i[i+1] - Baij->i[i]; 198 } 199 ierr = MatCreate(PETSC_COMM_SELF,&Bnew);CHKERRQ(ierr); 200 ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr); 201 ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr); 202 ierr = MatSeqAIJSetPreallocation(Bnew,0,nz);CHKERRQ(ierr); 203 ierr = PetscFree(nz);CHKERRQ(ierr); 204 for (i=0; i<m; i++) { 205 for (j=Baij->i[i]; j<Baij->i[i+1]; j++) { 206 col = garray[Baij->j[ct]]; 207 v = Baij->a[ct++]; 208 ierr = MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);CHKERRQ(ierr); 209 } 210 } 211 ierr = PetscFree(aij->garray);CHKERRQ(ierr); 212 ierr = PetscLogObjectMemory(A,-ec*sizeof(PetscInt));CHKERRQ(ierr); 213 ierr = MatDestroy(&B);CHKERRQ(ierr); 214 ierr = PetscLogObjectParent(A,Bnew);CHKERRQ(ierr); 215 aij->B = Bnew; 216 A->was_assembled = PETSC_FALSE; 217 PetscFunctionReturn(0); 218 } 219 220 /* ugly stuff added for Glenn someday we should fix this up */ 221 222 static PetscInt *auglyrmapd = 0,*auglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" 223 parts of the local matrix */ 224 static Vec auglydd = 0,auglyoo = 0; /* work vectors used to scale the two parts of the local matrix */ 225 226 227 #undef __FUNCT__ 228 #define __FUNCT__ "MatMPIAIJDiagonalScaleLocalSetUp" 229 PetscErrorCode MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale) 230 { 231 Mat_MPIAIJ *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */ 232 PetscErrorCode ierr; 233 PetscInt i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices; 234 PetscInt *r_rmapd,*r_rmapo; 235 236 PetscFunctionBegin; 237 ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr); 238 ierr = MatGetSize(ina->A,PETSC_NULL,&n);CHKERRQ(ierr); 239 ierr = PetscMalloc((inA->rmap->mapping->n+1)*sizeof(PetscInt),&r_rmapd);CHKERRQ(ierr); 240 ierr = PetscMemzero(r_rmapd,inA->rmap->mapping->n*sizeof(PetscInt));CHKERRQ(ierr); 241 nt = 0; 242 for (i=0; i<inA->rmap->mapping->n; i++) { 243 if (inA->rmap->mapping->indices[i] >= cstart && inA->rmap->mapping->indices[i] < cend) { 244 nt++; 245 r_rmapd[i] = inA->rmap->mapping->indices[i] + 1; 246 } 247 } 248 if (nt != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D n %D",nt,n); 249 ierr = PetscMalloc((n+1)*sizeof(PetscInt),&auglyrmapd);CHKERRQ(ierr); 250 for (i=0; i<inA->rmap->mapping->n; i++) { 251 if (r_rmapd[i]){ 252 auglyrmapd[(r_rmapd[i]-1)-cstart] = i; 253 } 254 } 255 ierr = PetscFree(r_rmapd);CHKERRQ(ierr); 256 ierr = VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);CHKERRQ(ierr); 257 258 ierr = PetscMalloc((inA->cmap->N+1)*sizeof(PetscInt),&lindices);CHKERRQ(ierr); 259 ierr = PetscMemzero(lindices,inA->cmap->N*sizeof(PetscInt));CHKERRQ(ierr); 260 for (i=0; i<ina->B->cmap->n; i++) { 261 lindices[garray[i]] = i+1; 262 } 263 no = inA->rmap->mapping->n - nt; 264 ierr = PetscMalloc((inA->rmap->mapping->n+1)*sizeof(PetscInt),&r_rmapo);CHKERRQ(ierr); 265 ierr = PetscMemzero(r_rmapo,inA->rmap->mapping->n*sizeof(PetscInt));CHKERRQ(ierr); 266 nt = 0; 267 for (i=0; i<inA->rmap->mapping->n; i++) { 268 if (lindices[inA->rmap->mapping->indices[i]]) { 269 nt++; 270 r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]]; 271 } 272 } 273 if (nt > no) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n); 274 ierr = PetscFree(lindices);CHKERRQ(ierr); 275 ierr = PetscMalloc((nt+1)*sizeof(PetscInt),&auglyrmapo);CHKERRQ(ierr); 276 for (i=0; i<inA->rmap->mapping->n; i++) { 277 if (r_rmapo[i]){ 278 auglyrmapo[(r_rmapo[i]-1)] = i; 279 } 280 } 281 ierr = PetscFree(r_rmapo);CHKERRQ(ierr); 282 ierr = VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);CHKERRQ(ierr); 283 284 PetscFunctionReturn(0); 285 } 286 287 #undef __FUNCT__ 288 #define __FUNCT__ "MatMPIAIJDiagonalScaleLocal" 289 PetscErrorCode MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale) 290 { 291 /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */ 292 PetscErrorCode ierr; 293 294 PetscFunctionBegin; 295 ierr = PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));CHKERRQ(ierr); 296 PetscFunctionReturn(0); 297 } 298 299 EXTERN_C_BEGIN 300 #undef __FUNCT__ 301 #define __FUNCT__ "MatDiagonalScaleLocal_MPIAIJ" 302 PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat A,Vec scale) 303 { 304 Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */ 305 PetscErrorCode ierr; 306 PetscInt n,i; 307 PetscScalar *d,*o,*s; 308 309 PetscFunctionBegin; 310 if (!auglyrmapd) { 311 ierr = MatMPIAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr); 312 } 313 314 ierr = VecGetArray(scale,&s);CHKERRQ(ierr); 315 316 ierr = VecGetLocalSize(auglydd,&n);CHKERRQ(ierr); 317 ierr = VecGetArray(auglydd,&d);CHKERRQ(ierr); 318 for (i=0; i<n; i++) { 319 d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */ 320 } 321 ierr = VecRestoreArray(auglydd,&d);CHKERRQ(ierr); 322 /* column scale "diagonal" portion of local matrix */ 323 ierr = MatDiagonalScale(a->A,PETSC_NULL,auglydd);CHKERRQ(ierr); 324 325 ierr = VecGetLocalSize(auglyoo,&n);CHKERRQ(ierr); 326 ierr = VecGetArray(auglyoo,&o);CHKERRQ(ierr); 327 for (i=0; i<n; i++) { 328 o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */ 329 } 330 ierr = VecRestoreArray(scale,&s);CHKERRQ(ierr); 331 ierr = VecRestoreArray(auglyoo,&o);CHKERRQ(ierr); 332 /* column scale "off-diagonal" portion of local matrix */ 333 ierr = MatDiagonalScale(a->B,PETSC_NULL,auglyoo);CHKERRQ(ierr); 334 335 PetscFunctionReturn(0); 336 } 337 EXTERN_C_END 338 339 340