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