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 #if defined (PETSC_USE_CTABLE) 18 PetscTable gid1_lid1; 19 PetscTablePosition tpos; 20 PetscInt gid,lid; 21 #else 22 PetscInt N = mat->cmap->N,*indices; 23 #endif 24 25 PetscFunctionBegin; 26 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 = PetscMalloc((ec+1)*sizeof(PetscInt),&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 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 = PetscMalloc((N+1)*sizeof(PetscInt),&indices);CHKERRQ(ierr); 70 ierr = PetscMemzero(indices,N*sizeof(PetscInt));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 = PetscMalloc((ec+1)*sizeof(PetscInt),&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 ierr = PetscLayoutSetUp((aij->B->cmap));CHKERRQ(ierr); 98 ierr = PetscFree(indices);CHKERRQ(ierr); 99 #endif 100 /* create local vector that is used to scatter into */ 101 ierr = VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);CHKERRQ(ierr); 102 103 /* create two temporary Index sets for build scatter gather */ 104 ierr = ISCreateGeneral(((PetscObject)mat)->comm,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 105 106 ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr); 107 108 /* create temporary global vector to generate scatter context */ 109 /* This does not allocate the array's memory so is efficient */ 110 ierr = VecCreateMPIWithArray(((PetscObject)mat)->comm,1,mat->cmap->n,mat->cmap->N,PETSC_NULL,&gvec);CHKERRQ(ierr); 111 112 /* generate the scatter context */ 113 ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);CHKERRQ(ierr); 114 ierr = PetscLogObjectParent(mat,aij->Mvctx);CHKERRQ(ierr); 115 ierr = PetscLogObjectParent(mat,aij->lvec);CHKERRQ(ierr); 116 ierr = PetscLogObjectParent(mat,from);CHKERRQ(ierr); 117 ierr = PetscLogObjectParent(mat,to);CHKERRQ(ierr); 118 aij->garray = garray; 119 ierr = PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr); 120 ierr = ISDestroy(&from);CHKERRQ(ierr); 121 ierr = ISDestroy(&to);CHKERRQ(ierr); 122 ierr = VecDestroy(&gvec);CHKERRQ(ierr); 123 PetscFunctionReturn(0); 124 } 125 126 127 #undef __FUNCT__ 128 #define __FUNCT__ "MatDisAssemble_MPIAIJ" 129 /* 130 Takes the local part of an already assembled MPIAIJ matrix 131 and disassembles it. This is to allow new nonzeros into the matrix 132 that require more communication in the matrix vector multiply. 133 Thus certain data-structures must be rebuilt. 134 135 Kind of slow! But that's what application programmers get when 136 they are sloppy. 137 */ 138 PetscErrorCode MatDisAssemble_MPIAIJ(Mat A) 139 { 140 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 141 Mat B = aij->B,Bnew; 142 Mat_SeqAIJ *Baij = (Mat_SeqAIJ*)B->data; 143 PetscErrorCode ierr; 144 PetscInt i,j,m = B->rmap->n,n = A->cmap->N,col,ct = 0,*garray = aij->garray,*nz,ec; 145 PetscScalar v; 146 147 PetscFunctionBegin; 148 /* free stuff related to matrix-vec multiply */ 149 ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */ 150 ierr = VecDestroy(&aij->lvec);CHKERRQ(ierr); 151 ierr = VecScatterDestroy(&aij->Mvctx);CHKERRQ(ierr); 152 if (aij->colmap) { 153 #if defined (PETSC_USE_CTABLE) 154 ierr = PetscTableDestroy(&aij->colmap);CHKERRQ(ierr); 155 #else 156 ierr = PetscFree(aij->colmap);CHKERRQ(ierr); 157 ierr = PetscLogObjectMemory(A,-aij->B->cmap->n*sizeof(PetscInt));CHKERRQ(ierr); 158 #endif 159 } 160 161 /* make sure that B is assembled so we can access its values */ 162 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 163 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 164 165 /* invent new B and copy stuff over */ 166 ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nz);CHKERRQ(ierr); 167 for (i=0; i<m; i++) { 168 nz[i] = Baij->i[i+1] - Baij->i[i]; 169 } 170 ierr = MatCreate(PETSC_COMM_SELF,&Bnew);CHKERRQ(ierr); 171 ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr); 172 ierr = MatSetBlockSizes(Bnew,A->rmap->bs,A->cmap->bs);CHKERRQ(ierr); 173 ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr); 174 ierr = MatSeqAIJSetPreallocation(Bnew,0,nz);CHKERRQ(ierr); 175 ((Mat_SeqAIJ*)Bnew->data)->nonew = Baij->nonew; /* Inherit insertion error options. */ 176 ierr = PetscFree(nz);CHKERRQ(ierr); 177 for (i=0; i<m; i++) { 178 for (j=Baij->i[i]; j<Baij->i[i+1]; j++) { 179 col = garray[Baij->j[ct]]; 180 v = Baij->a[ct++]; 181 ierr = MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);CHKERRQ(ierr); 182 } 183 } 184 ierr = PetscFree(aij->garray);CHKERRQ(ierr); 185 ierr = PetscLogObjectMemory(A,-ec*sizeof(PetscInt));CHKERRQ(ierr); 186 ierr = MatDestroy(&B);CHKERRQ(ierr); 187 ierr = PetscLogObjectParent(A,Bnew);CHKERRQ(ierr); 188 aij->B = Bnew; 189 A->was_assembled = PETSC_FALSE; 190 PetscFunctionReturn(0); 191 } 192 193 /* ugly stuff added for Glenn someday we should fix this up */ 194 195 static PetscInt *auglyrmapd = 0,*auglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" 196 parts of the local matrix */ 197 static Vec auglydd = 0,auglyoo = 0; /* work vectors used to scale the two parts of the local matrix */ 198 199 200 #undef __FUNCT__ 201 #define __FUNCT__ "MatMPIAIJDiagonalScaleLocalSetUp" 202 PetscErrorCode MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale) 203 { 204 Mat_MPIAIJ *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */ 205 PetscErrorCode ierr; 206 PetscInt i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices; 207 PetscInt *r_rmapd,*r_rmapo; 208 209 PetscFunctionBegin; 210 ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr); 211 ierr = MatGetSize(ina->A,PETSC_NULL,&n);CHKERRQ(ierr); 212 ierr = PetscMalloc((inA->rmap->mapping->n+1)*sizeof(PetscInt),&r_rmapd);CHKERRQ(ierr); 213 ierr = PetscMemzero(r_rmapd,inA->rmap->mapping->n*sizeof(PetscInt));CHKERRQ(ierr); 214 nt = 0; 215 for (i=0; i<inA->rmap->mapping->n; i++) { 216 if (inA->rmap->mapping->indices[i] >= cstart && inA->rmap->mapping->indices[i] < cend) { 217 nt++; 218 r_rmapd[i] = inA->rmap->mapping->indices[i] + 1; 219 } 220 } 221 if (nt != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D n %D",nt,n); 222 ierr = PetscMalloc((n+1)*sizeof(PetscInt),&auglyrmapd);CHKERRQ(ierr); 223 for (i=0; i<inA->rmap->mapping->n; i++) { 224 if (r_rmapd[i]){ 225 auglyrmapd[(r_rmapd[i]-1)-cstart] = i; 226 } 227 } 228 ierr = PetscFree(r_rmapd);CHKERRQ(ierr); 229 ierr = VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);CHKERRQ(ierr); 230 231 ierr = PetscMalloc((inA->cmap->N+1)*sizeof(PetscInt),&lindices);CHKERRQ(ierr); 232 ierr = PetscMemzero(lindices,inA->cmap->N*sizeof(PetscInt));CHKERRQ(ierr); 233 for (i=0; i<ina->B->cmap->n; i++) { 234 lindices[garray[i]] = i+1; 235 } 236 no = inA->rmap->mapping->n - nt; 237 ierr = PetscMalloc((inA->rmap->mapping->n+1)*sizeof(PetscInt),&r_rmapo);CHKERRQ(ierr); 238 ierr = PetscMemzero(r_rmapo,inA->rmap->mapping->n*sizeof(PetscInt));CHKERRQ(ierr); 239 nt = 0; 240 for (i=0; i<inA->rmap->mapping->n; i++) { 241 if (lindices[inA->rmap->mapping->indices[i]]) { 242 nt++; 243 r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]]; 244 } 245 } 246 if (nt > no) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n); 247 ierr = PetscFree(lindices);CHKERRQ(ierr); 248 ierr = PetscMalloc((nt+1)*sizeof(PetscInt),&auglyrmapo);CHKERRQ(ierr); 249 for (i=0; i<inA->rmap->mapping->n; i++) { 250 if (r_rmapo[i]){ 251 auglyrmapo[(r_rmapo[i]-1)] = i; 252 } 253 } 254 ierr = PetscFree(r_rmapo);CHKERRQ(ierr); 255 ierr = VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);CHKERRQ(ierr); 256 257 PetscFunctionReturn(0); 258 } 259 260 #undef __FUNCT__ 261 #define __FUNCT__ "MatMPIAIJDiagonalScaleLocal" 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 EXTERN_C_BEGIN 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,PETSC_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,PETSC_NULL,auglyoo);CHKERRQ(ierr); 307 308 PetscFunctionReturn(0); 309 } 310 EXTERN_C_END 311 312 313