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 #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 #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 #if defined(PETSC_USE_CTABLE) 28 /* use a table */ 29 ierr = PetscTableCreate(aij->B->rmap->n,mat->cmap->N+1,&gid1_lid1);CHKERRQ(ierr); 30 for (i=0; i<aij->B->rmap->n; i++) { 31 for (j=0; j<B->ilen[i]; j++) { 32 PetscInt data,gid1 = aj[B->i[i] + j] + 1; 33 ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr); 34 if (!data) { 35 /* one based table */ 36 ierr = PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);CHKERRQ(ierr); 37 } 38 } 39 } 40 /* form array of columns we need */ 41 ierr = PetscMalloc1((ec+1),&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,INSERT_VALUES);CHKERRQ(ierr); 53 } 54 /* compact out the extra columns in B */ 55 for (i=0; i<aij->B->rmap->n; i++) { 56 for (j=0; j<B->ilen[i]; j++) { 57 PetscInt gid1 = aj[B->i[i] + j] + 1; 58 ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr); 59 lid--; 60 aj[B->i[i] + j] = lid; 61 } 62 } 63 aij->B->cmap->n = aij->B->cmap->N = ec; 64 65 ierr = PetscLayoutSetUp((aij->B->cmap));CHKERRQ(ierr); 66 ierr = PetscTableDestroy(&gid1_lid1);CHKERRQ(ierr); 67 #else 68 /* Make an array as long as the number of columns */ 69 /* mark those columns that are in aij->B */ 70 ierr = PetscCalloc1(N+1,&indices);CHKERRQ(ierr); 71 for (i=0; i<aij->B->rmap->n; i++) { 72 for (j=0; j<B->ilen[i]; j++) { 73 if (!indices[aj[B->i[i] + j]]) ec++; 74 indices[aj[B->i[i] + j]] = 1; 75 } 76 } 77 78 /* form array of columns we need */ 79 ierr = PetscMalloc1((ec+1),&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; 88 } 89 90 /* compact out the extra columns in B */ 91 for (i=0; i<aij->B->rmap->n; i++) { 92 for (j=0; j<B->ilen[i]; j++) { 93 aj[B->i[i] + j] = indices[aj[B->i[i] + j]]; 94 } 95 } 96 aij->B->cmap->n = aij->B->cmap->N = ec; 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 #undef __FUNCT__ 131 #define __FUNCT__ "MatDisAssemble_MPIAIJ" 132 /* 133 Takes the local part of an already assembled MPIAIJ matrix 134 and disassembles it. This is to allow new nonzeros into the matrix 135 that require more communication in the matrix vector multiply. 136 Thus certain data-structures must be rebuilt. 137 138 Kind of slow! But that's what application programmers get when 139 they are sloppy. 140 */ 141 PetscErrorCode MatDisAssemble_MPIAIJ(Mat A) 142 { 143 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 144 Mat B = aij->B,Bnew; 145 Mat_SeqAIJ *Baij = (Mat_SeqAIJ*)B->data; 146 PetscErrorCode ierr; 147 PetscInt i,j,m = B->rmap->n,n = A->cmap->N,col,ct = 0,*garray = aij->garray,*nz,ec; 148 PetscScalar v; 149 150 PetscFunctionBegin; 151 /* free stuff related to matrix-vec multiply */ 152 ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */ 153 ierr = VecDestroy(&aij->lvec);CHKERRQ(ierr); 154 ierr = VecScatterDestroy(&aij->Mvctx);CHKERRQ(ierr); 155 if (aij->colmap) { 156 #if defined(PETSC_USE_CTABLE) 157 ierr = PetscTableDestroy(&aij->colmap);CHKERRQ(ierr); 158 #else 159 ierr = PetscFree(aij->colmap);CHKERRQ(ierr); 160 ierr = PetscLogObjectMemory((PetscObject)A,-aij->B->cmap->n*sizeof(PetscInt));CHKERRQ(ierr); 161 #endif 162 } 163 164 /* make sure that B is assembled so we can access its values */ 165 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 166 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 167 168 /* invent new B and copy stuff over */ 169 ierr = PetscMalloc1((m+1),&nz);CHKERRQ(ierr); 170 for (i=0; i<m; i++) { 171 nz[i] = Baij->i[i+1] - Baij->i[i]; 172 } 173 ierr = MatCreate(PETSC_COMM_SELF,&Bnew);CHKERRQ(ierr); 174 ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr); 175 ierr = MatSetBlockSizes(Bnew,A->rmap->bs,A->cmap->bs);CHKERRQ(ierr); 176 ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr); 177 ierr = MatSeqAIJSetPreallocation(Bnew,0,nz);CHKERRQ(ierr); 178 179 ((Mat_SeqAIJ*)Bnew->data)->nonew = Baij->nonew; /* Inherit insertion error options. */ 180 181 ierr = PetscFree(nz);CHKERRQ(ierr); 182 for (i=0; i<m; i++) { 183 for (j=Baij->i[i]; j<Baij->i[i+1]; j++) { 184 col = garray[Baij->j[ct]]; 185 v = Baij->a[ct++]; 186 ierr = MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);CHKERRQ(ierr); 187 } 188 } 189 ierr = PetscFree(aij->garray);CHKERRQ(ierr); 190 ierr = PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));CHKERRQ(ierr); 191 ierr = MatDestroy(&B);CHKERRQ(ierr); 192 ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);CHKERRQ(ierr); 193 194 aij->B = Bnew; 195 A->was_assembled = PETSC_FALSE; 196 PetscFunctionReturn(0); 197 } 198 199 /* ugly stuff added for Glenn someday we should fix this up */ 200 201 static PetscInt *auglyrmapd = 0,*auglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */ 202 static Vec auglydd = 0,auglyoo = 0; /* work vectors used to scale the two parts of the local matrix */ 203 204 205 #undef __FUNCT__ 206 #define __FUNCT__ "MatMPIAIJDiagonalScaleLocalSetUp" 207 PetscErrorCode MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale) 208 { 209 Mat_MPIAIJ *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */ 210 PetscErrorCode ierr; 211 PetscInt i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices; 212 PetscInt *r_rmapd,*r_rmapo; 213 214 PetscFunctionBegin; 215 ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr); 216 ierr = MatGetSize(ina->A,NULL,&n);CHKERRQ(ierr); 217 ierr = PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapd);CHKERRQ(ierr); 218 nt = 0; 219 for (i=0; i<inA->rmap->mapping->n; i++) { 220 if (inA->rmap->mapping->indices[i] >= cstart && inA->rmap->mapping->indices[i] < cend) { 221 nt++; 222 r_rmapd[i] = inA->rmap->mapping->indices[i] + 1; 223 } 224 } 225 if (nt != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D n %D",nt,n); 226 ierr = PetscMalloc1((n+1),&auglyrmapd);CHKERRQ(ierr); 227 for (i=0; i<inA->rmap->mapping->n; i++) { 228 if (r_rmapd[i]) { 229 auglyrmapd[(r_rmapd[i]-1)-cstart] = i; 230 } 231 } 232 ierr = PetscFree(r_rmapd);CHKERRQ(ierr); 233 ierr = VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);CHKERRQ(ierr); 234 235 ierr = PetscCalloc1(inA->cmap->N+1,&lindices);CHKERRQ(ierr); 236 for (i=0; i<ina->B->cmap->n; i++) { 237 lindices[garray[i]] = i+1; 238 } 239 no = inA->rmap->mapping->n - nt; 240 ierr = PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapo);CHKERRQ(ierr); 241 nt = 0; 242 for (i=0; i<inA->rmap->mapping->n; i++) { 243 if (lindices[inA->rmap->mapping->indices[i]]) { 244 nt++; 245 r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]]; 246 } 247 } 248 if (nt > no) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n); 249 ierr = PetscFree(lindices);CHKERRQ(ierr); 250 ierr = PetscMalloc1((nt+1),&auglyrmapo);CHKERRQ(ierr); 251 for (i=0; i<inA->rmap->mapping->n; i++) { 252 if (r_rmapo[i]) { 253 auglyrmapo[(r_rmapo[i]-1)] = i; 254 } 255 } 256 ierr = PetscFree(r_rmapo);CHKERRQ(ierr); 257 ierr = VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);CHKERRQ(ierr); 258 PetscFunctionReturn(0); 259 } 260 261 #undef __FUNCT__ 262 #define __FUNCT__ "MatMPIAIJDiagonalScaleLocal" 263 PetscErrorCode MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale) 264 { 265 /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */ 266 PetscErrorCode ierr; 267 268 PetscFunctionBegin; 269 ierr = PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));CHKERRQ(ierr); 270 PetscFunctionReturn(0); 271 } 272 273 #undef __FUNCT__ 274 #define __FUNCT__ "MatDiagonalScaleLocal_MPIAIJ" 275 PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat A,Vec scale) 276 { 277 Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */ 278 PetscErrorCode ierr; 279 PetscInt n,i; 280 PetscScalar *d,*o,*s; 281 282 PetscFunctionBegin; 283 if (!auglyrmapd) { 284 ierr = MatMPIAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr); 285 } 286 287 ierr = VecGetArray(scale,&s);CHKERRQ(ierr); 288 289 ierr = VecGetLocalSize(auglydd,&n);CHKERRQ(ierr); 290 ierr = VecGetArray(auglydd,&d);CHKERRQ(ierr); 291 for (i=0; i<n; i++) { 292 d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */ 293 } 294 ierr = VecRestoreArray(auglydd,&d);CHKERRQ(ierr); 295 /* column scale "diagonal" portion of local matrix */ 296 ierr = MatDiagonalScale(a->A,NULL,auglydd);CHKERRQ(ierr); 297 298 ierr = VecGetLocalSize(auglyoo,&n);CHKERRQ(ierr); 299 ierr = VecGetArray(auglyoo,&o);CHKERRQ(ierr); 300 for (i=0; i<n; i++) { 301 o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */ 302 } 303 ierr = VecRestoreArray(scale,&s);CHKERRQ(ierr); 304 ierr = VecRestoreArray(auglyoo,&o);CHKERRQ(ierr); 305 /* column scale "off-diagonal" portion of local matrix */ 306 ierr = MatDiagonalScale(a->B,NULL,auglyoo);CHKERRQ(ierr); 307 PetscFunctionReturn(0); 308 } 309 310 311