1 2 #include <../src/mat/impls/aij/mpi/mpiaij.h> 3 4 extern PetscErrorCode CreateColmap_MPIAIJ_Private(Mat); 5 6 #undef __FUNCT__ 7 #define __FUNCT__ "MatFDColoringCreate_MPIAIJ" 8 PetscErrorCode MatFDColoringCreate_MPIAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c) 9 { 10 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)mat->data; 11 PetscErrorCode ierr; 12 PetscMPIInt size,*ncolsonproc,*disp,nn; 13 PetscInt i,n,nrows,j,k,m,*rows = 0,*A_ci,*A_cj,ncols,col; 14 const PetscInt *is; 15 PetscInt nis = iscoloring->n,nctot,*cols,*B_ci,*B_cj; 16 PetscInt *rowhit,M,cstart,cend,colb; 17 PetscInt *columnsforrow,l; 18 IS *isa; 19 PetscBool done,flg; 20 ISLocalToGlobalMapping map = mat->cmap->mapping; 21 PetscInt *ltog = (map ? map->indices : (PetscInt*) PETSC_NULL) ,ctype=c->ctype; 22 23 PetscFunctionBegin; 24 if (!mat->assembled) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled first; MatAssemblyBegin/End();"); 25 if (ctype == IS_COLORING_GHOSTED && !map) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_ARG_INCOMP,"When using ghosted differencing matrix must have local to global mapping provided with MatSetLocalToGlobalMapping"); 26 27 ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 28 29 M = mat->rmap->n; 30 cstart = mat->cmap->rstart; 31 cend = mat->cmap->rend; 32 c->M = mat->rmap->N; /* set the global rows and columns and local rows */ 33 c->N = mat->cmap->N; 34 c->m = mat->rmap->n; 35 c->rstart = mat->rmap->rstart; 36 37 c->ncolors = nis; 38 ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 39 ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr); 40 ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 41 ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->rows);CHKERRQ(ierr); 42 ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columnsforrow);CHKERRQ(ierr); 43 ierr = PetscLogObjectMemory(c,5*nis*sizeof(PetscInt));CHKERRQ(ierr); 44 45 /* Allow access to data structures of local part of matrix */ 46 if (!aij->colmap) { 47 ierr = CreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr); 48 } 49 ierr = MatGetColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr); 50 ierr = MatGetColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr); 51 52 ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr); 53 ierr = PetscMalloc((M+1)*sizeof(PetscInt),&columnsforrow);CHKERRQ(ierr); 54 55 for (i=0; i<nis; i++) { 56 ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 57 ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 58 c->ncolumns[i] = n; 59 if (n) { 60 ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr); 61 ierr = PetscLogObjectMemory(c,n*sizeof(PetscInt));CHKERRQ(ierr); 62 ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr); 63 } else { 64 c->columns[i] = 0; 65 } 66 67 if (ctype == IS_COLORING_GLOBAL){ 68 /* Determine the total (parallel) number of columns of this color */ 69 ierr = MPI_Comm_size(((PetscObject)mat)->comm,&size);CHKERRQ(ierr); 70 ierr = PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);CHKERRQ(ierr); 71 72 nn = PetscMPIIntCast(n); 73 ierr = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,((PetscObject)mat)->comm);CHKERRQ(ierr); 74 nctot = 0; for (j=0; j<size; j++) {nctot += ncolsonproc[j];} 75 if (!nctot) { 76 ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr); 77 } 78 79 disp[0] = 0; 80 for (j=1; j<size; j++) { 81 disp[j] = disp[j-1] + ncolsonproc[j-1]; 82 } 83 84 /* Get complete list of columns for color on each processor */ 85 ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 86 ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,((PetscObject)mat)->comm);CHKERRQ(ierr); 87 ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr); 88 } else if (ctype == IS_COLORING_GHOSTED){ 89 /* Determine local number of columns of this color on this process, including ghost points */ 90 nctot = n; 91 ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 92 ierr = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr); 93 } else { 94 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type"); 95 } 96 97 /* 98 Mark all rows affect by these columns 99 */ 100 /* Temporary option to allow for debugging/testing */ 101 flg = PETSC_FALSE; 102 ierr = PetscOptionsGetBool(PETSC_NULL,"-matfdcoloring_slow",&flg,PETSC_NULL);CHKERRQ(ierr); 103 if (!flg) {/*-----------------------------------------------------------------------------*/ 104 /* crude, fast version */ 105 ierr = PetscMemzero(rowhit,M*sizeof(PetscInt));CHKERRQ(ierr); 106 /* loop over columns*/ 107 for (j=0; j<nctot; j++) { 108 if (ctype == IS_COLORING_GHOSTED) { 109 col = ltog[cols[j]]; 110 } else { 111 col = cols[j]; 112 } 113 if (col >= cstart && col < cend) { 114 /* column is in diagonal block of matrix */ 115 rows = A_cj + A_ci[col-cstart]; 116 m = A_ci[col-cstart+1] - A_ci[col-cstart]; 117 } else { 118 #if defined (PETSC_USE_CTABLE) 119 ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr); 120 colb --; 121 #else 122 colb = aij->colmap[col] - 1; 123 #endif 124 if (colb == -1) { 125 m = 0; 126 } else { 127 rows = B_cj + B_ci[colb]; 128 m = B_ci[colb+1] - B_ci[colb]; 129 } 130 } 131 /* loop over columns marking them in rowhit */ 132 for (k=0; k<m; k++) { 133 rowhit[*rows++] = col + 1; 134 } 135 } 136 137 /* count the number of hits */ 138 nrows = 0; 139 for (j=0; j<M; j++) { 140 if (rowhit[j]) nrows++; 141 } 142 c->nrows[i] = nrows; 143 ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr); 144 ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr); 145 ierr = PetscLogObjectMemory(c,2*(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr); 146 nrows = 0; 147 for (j=0; j<M; j++) { 148 if (rowhit[j]) { 149 c->rows[i][nrows] = j; 150 c->columnsforrow[i][nrows] = rowhit[j] - 1; 151 nrows++; 152 } 153 } 154 } else {/*-------------------------------------------------------------------------------*/ 155 /* slow version, using rowhit as a linked list */ 156 PetscInt currentcol,fm,mfm; 157 rowhit[M] = M; 158 nrows = 0; 159 /* loop over columns*/ 160 for (j=0; j<nctot; j++) { 161 if (ctype == IS_COLORING_GHOSTED) { 162 col = ltog[cols[j]]; 163 } else { 164 col = cols[j]; 165 } 166 if (col >= cstart && col < cend) { 167 /* column is in diagonal block of matrix */ 168 rows = A_cj + A_ci[col-cstart]; 169 m = A_ci[col-cstart+1] - A_ci[col-cstart]; 170 } else { 171 #if defined (PETSC_USE_CTABLE) 172 ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr); 173 colb --; 174 #else 175 colb = aij->colmap[col] - 1; 176 #endif 177 if (colb == -1) { 178 m = 0; 179 } else { 180 rows = B_cj + B_ci[colb]; 181 m = B_ci[colb+1] - B_ci[colb]; 182 } 183 } 184 185 /* loop over columns marking them in rowhit */ 186 fm = M; /* fm points to first entry in linked list */ 187 for (k=0; k<m; k++) { 188 currentcol = *rows++; 189 /* is it already in the list? */ 190 do { 191 mfm = fm; 192 fm = rowhit[fm]; 193 } while (fm < currentcol); 194 /* not in list so add it */ 195 if (fm != currentcol) { 196 nrows++; 197 columnsforrow[currentcol] = col; 198 /* next three lines insert new entry into linked list */ 199 rowhit[mfm] = currentcol; 200 rowhit[currentcol] = fm; 201 fm = currentcol; 202 /* fm points to present position in list since we know the columns are sorted */ 203 } else { 204 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid coloring of matrix detected"); 205 } 206 } 207 } 208 c->nrows[i] = nrows; 209 ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr); 210 ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr); 211 ierr = PetscLogObjectMemory(c,(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr); 212 /* now store the linked list of rows into c->rows[i] */ 213 nrows = 0; 214 fm = rowhit[M]; 215 do { 216 c->rows[i][nrows] = fm; 217 c->columnsforrow[i][nrows++] = columnsforrow[fm]; 218 fm = rowhit[fm]; 219 } while (fm < M); 220 } /* ---------------------------------------------------------------------------------------*/ 221 ierr = PetscFree(cols);CHKERRQ(ierr); 222 } 223 224 /* Optimize by adding the vscale, and scaleforrow[][] fields */ 225 /* 226 vscale will contain the "diagonal" on processor scalings followed by the off processor 227 */ 228 if (ctype == IS_COLORING_GLOBAL) { 229 ierr = VecCreateGhost(((PetscObject)mat)->comm,aij->A->rmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&c->vscale);CHKERRQ(ierr); 230 ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr); 231 for (k=0; k<c->ncolors; k++) { 232 ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr); 233 for (l=0; l<c->nrows[k]; l++) { 234 col = c->columnsforrow[k][l]; 235 if (col >= cstart && col < cend) { 236 /* column is in diagonal block of matrix */ 237 colb = col - cstart; 238 } else { 239 /* column is in "off-processor" part */ 240 #if defined (PETSC_USE_CTABLE) 241 ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr); 242 colb --; 243 #else 244 colb = aij->colmap[col] - 1; 245 #endif 246 colb += cend - cstart; 247 } 248 c->vscaleforrow[k][l] = colb; 249 } 250 } 251 } else if (ctype == IS_COLORING_GHOSTED) { 252 /* Get gtol mapping */ 253 PetscInt N = mat->cmap->N, *gtol; 254 ierr = PetscMalloc((N+1)*sizeof(PetscInt),>ol);CHKERRQ(ierr); 255 for (i=0; i<N; i++) gtol[i] = -1; 256 for (i=0; i<map->n; i++) gtol[ltog[i]] = i; 257 258 c->vscale = 0; /* will be created in MatFDColoringApply() */ 259 ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr); 260 for (k=0; k<c->ncolors; k++) { 261 ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr); 262 for (l=0; l<c->nrows[k]; l++) { 263 col = c->columnsforrow[k][l]; /* global column index */ 264 c->vscaleforrow[k][l] = gtol[col]; /* local column index */ 265 } 266 } 267 ierr = PetscFree(gtol);CHKERRQ(ierr); 268 } 269 ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 270 271 ierr = PetscFree(rowhit);CHKERRQ(ierr); 272 ierr = PetscFree(columnsforrow);CHKERRQ(ierr); 273 ierr = MatRestoreColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr); 274 ierr = MatRestoreColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr); 275 PetscFunctionReturn(0); 276 } 277 278 279 280 281 282 283