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