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