1be1d678aSKris Buschelman #define PETSCMAT_DLL 2a64fbb32SBarry Smith 36eaac0f3SBarry Smith #include "src/mat/impls/aij/mpi/mpiaij.h" 4a64fbb32SBarry Smith 5dfbe8321SBarry Smith EXTERN PetscErrorCode CreateColmap_MPIAIJ_Private(Mat); 6a64fbb32SBarry Smith 74a2ae208SSatish Balay #undef __FUNCT__ 84a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate_MPIAIJ" 9dfbe8321SBarry Smith PetscErrorCode MatFDColoringCreate_MPIAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c) 10a64fbb32SBarry Smith { 116eaac0f3SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)mat->data; 126849ba73SBarry Smith PetscErrorCode ierr; 13b1d57f15SBarry Smith PetscMPIInt size,*ncolsonproc,*disp,nn; 14*5d0c19d7SBarry Smith PetscInt i,n,nrows,j,k,m,*rows = 0,*A_ci,*A_cj,ncols,col; 15*5d0c19d7SBarry Smith const PetscInt *is; 16b1d57f15SBarry Smith PetscInt nis = iscoloring->n,nctot,*cols,*B_ci,*B_cj; 17d0f46423SBarry Smith PetscInt *rowhit,M = mat->rmap->n,cstart = mat->cmap->rstart,cend = mat->cmap->rend,colb; 18b1d57f15SBarry Smith PetscInt *columnsforrow,l; 19b9617806SBarry Smith IS *isa; 20f1af5d2fSBarry Smith PetscTruth done,flg; 21b8f8c88eSHong Zhang ISLocalToGlobalMapping map = mat->mapping; 22e60d7cb6SSatish Balay PetscInt *ltog = (map ? map->indices : (PetscInt*) PETSC_NULL) ,ctype=c->ctype; 23a64fbb32SBarry Smith 243a40ed3dSBarry Smith PetscFunctionBegin; 25522c5e43SBarry Smith if (!mat->assembled) { 2629bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled first; MatAssemblyBegin/End();"); 27522c5e43SBarry Smith } 287cc688d7SBarry Smith if (ctype == IS_COLORING_GHOSTED && !map) SETERRQ(PETSC_ERR_ARG_INCOMP,"When using ghosted differencing matrix must have local to global mapping provided with MatSetLocalToGlobalMapping"); 29522c5e43SBarry Smith 30b9617806SBarry Smith ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 31d0f46423SBarry Smith c->M = mat->rmap->N; /* set the global rows and columns and local rows */ 32d0f46423SBarry Smith c->N = mat->cmap->N; 33d0f46423SBarry Smith c->m = mat->rmap->n; 34d0f46423SBarry Smith c->rstart = mat->rmap->rstart; 35005c665bSBarry Smith 36a64fbb32SBarry Smith c->ncolors = nis; 37b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 38b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr); 39b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 40b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->rows);CHKERRQ(ierr); 41b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columnsforrow);CHKERRQ(ierr); 4252e6d16bSBarry Smith ierr = PetscLogObjectMemory(c,5*nis*sizeof(PetscInt));CHKERRQ(ierr); 436eaac0f3SBarry Smith 446eaac0f3SBarry Smith /* Allow access to data structures of local part of matrix */ 456eaac0f3SBarry Smith if (!aij->colmap) { 466eaac0f3SBarry Smith ierr = CreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr); 476eaac0f3SBarry Smith } 4843a90d84SBarry Smith /* 4943a90d84SBarry Smith Calls the _SeqAIJ() version of these routines to make sure it does not 5043a90d84SBarry Smith get the reduced (by inodes) version of I and J 5143a90d84SBarry Smith */ 528f7157efSSatish Balay ierr = MatGetColumnIJ_SeqAIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr); 538f7157efSSatish Balay ierr = MatGetColumnIJ_SeqAIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr); 546eaac0f3SBarry Smith 55b1d57f15SBarry Smith ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr); 56b1d57f15SBarry Smith ierr = PetscMalloc((M+1)*sizeof(PetscInt),&columnsforrow);CHKERRQ(ierr); 576eaac0f3SBarry Smith 58a64fbb32SBarry Smith for (i=0; i<nis; i++) { 59b9b97703SBarry Smith ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 60a64fbb32SBarry Smith ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 61a64fbb32SBarry Smith c->ncolumns[i] = n; 62a64fbb32SBarry Smith if (n) { 63b1d57f15SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr); 6452e6d16bSBarry Smith ierr = PetscLogObjectMemory(c,n*sizeof(PetscInt));CHKERRQ(ierr); 65b1d57f15SBarry Smith ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr); 66a64fbb32SBarry Smith } else { 67a64fbb32SBarry Smith c->columns[i] = 0; 68a64fbb32SBarry Smith } 69a64fbb32SBarry Smith 708ee2e534SBarry Smith if (ctype == IS_COLORING_GLOBAL){ 716eaac0f3SBarry Smith /* Determine the total (parallel) number of columns of this color */ 727adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)mat)->comm,&size);CHKERRQ(ierr); 73b8f8c88eSHong Zhang ierr = PetscMalloc(2*size*sizeof(PetscInt*),&ncolsonproc);CHKERRQ(ierr); 74b8f8c88eSHong Zhang disp = ncolsonproc + size; 75b8f8c88eSHong Zhang 76e44c0bd4SBarry Smith nn = PetscMPIIntCast(n); 777adad957SLisandro Dalcin ierr = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,((PetscObject)mat)->comm);CHKERRQ(ierr); 786eaac0f3SBarry Smith nctot = 0; for (j=0; j<size; j++) {nctot += ncolsonproc[j];} 793a7fca6bSBarry Smith if (!nctot) { 80ae15b995SBarry Smith ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr); 813a7fca6bSBarry Smith } 826eaac0f3SBarry Smith 836eaac0f3SBarry Smith disp[0] = 0; 846eaac0f3SBarry Smith for (j=1; j<size; j++) { 856eaac0f3SBarry Smith disp[j] = disp[j-1] + ncolsonproc[j-1]; 866eaac0f3SBarry Smith } 876eaac0f3SBarry Smith 886eaac0f3SBarry Smith /* Get complete list of columns for color on each processor */ 89b1d57f15SBarry Smith ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 90*5d0c19d7SBarry Smith ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,((PetscObject)mat)->comm);CHKERRQ(ierr); 91b8f8c88eSHong Zhang ierr = PetscFree(ncolsonproc);CHKERRQ(ierr); 92b8f8c88eSHong Zhang } else if (ctype == IS_COLORING_GHOSTED){ 93b8f8c88eSHong Zhang /* Determine local number of columns of this color on this process, including ghost points */ 94b8f8c88eSHong Zhang nctot = n; 95b8f8c88eSHong Zhang ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 96b8f8c88eSHong Zhang ierr = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr); 97b8f8c88eSHong Zhang } else { 98b8f8c88eSHong Zhang SETERRQ(PETSC_ERR_SUP,"Not provided for this MatFDColoring type"); 99b8f8c88eSHong Zhang } 1006eaac0f3SBarry Smith 1016eaac0f3SBarry Smith /* 1026eaac0f3SBarry Smith Mark all rows affect by these columns 1036eaac0f3SBarry Smith */ 104b8f8c88eSHong Zhang /* Temporary option to allow for debugging/testing */ 105b8f8c88eSHong Zhang ierr = PetscOptionsHasName(PETSC_NULL,"-matfdcoloring_slow",&flg);CHKERRQ(ierr); 106f158e583SBarry Smith if (!flg) {/*-----------------------------------------------------------------------------*/ 107f158e583SBarry Smith /* crude, fast version */ 108b1d57f15SBarry Smith ierr = PetscMemzero(rowhit,M*sizeof(PetscInt));CHKERRQ(ierr); 109a64fbb32SBarry Smith /* loop over columns*/ 1106eaac0f3SBarry Smith for (j=0; j<nctot; j++) { 111b8f8c88eSHong Zhang if (ctype == IS_COLORING_GHOSTED) { 112b8f8c88eSHong Zhang col = ltog[cols[j]]; 113b8f8c88eSHong Zhang } else { 1146eaac0f3SBarry Smith col = cols[j]; 115b8f8c88eSHong Zhang } 1166eaac0f3SBarry Smith if (col >= cstart && col < cend) { 1176eaac0f3SBarry Smith /* column is in diagonal block of matrix */ 1186eaac0f3SBarry Smith rows = A_cj + A_ci[col-cstart]; 1196eaac0f3SBarry Smith m = A_ci[col-cstart+1] - A_ci[col-cstart]; 1206eaac0f3SBarry Smith } else { 121aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 1220f5bd95cSBarry Smith ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr) 123fa46199cSSatish Balay colb --; 124b3d2dc96SSatish Balay #else 1256eaac0f3SBarry Smith colb = aij->colmap[col] - 1; 126b3d2dc96SSatish Balay #endif 1276eaac0f3SBarry Smith if (colb == -1) { 1286eaac0f3SBarry Smith m = 0; 1296eaac0f3SBarry Smith } else { 1306eaac0f3SBarry Smith rows = B_cj + B_ci[colb]; 1316eaac0f3SBarry Smith m = B_ci[colb+1] - B_ci[colb]; 1326eaac0f3SBarry Smith } 1336eaac0f3SBarry Smith } 134a64fbb32SBarry Smith /* loop over columns marking them in rowhit */ 135a64fbb32SBarry Smith for (k=0; k<m; k++) { 136a64fbb32SBarry Smith rowhit[*rows++] = col + 1; 137a64fbb32SBarry Smith } 138a64fbb32SBarry Smith } 1396eaac0f3SBarry Smith 140a64fbb32SBarry Smith /* count the number of hits */ 141a64fbb32SBarry Smith nrows = 0; 1426eaac0f3SBarry Smith for (j=0; j<M; j++) { 143a64fbb32SBarry Smith if (rowhit[j]) nrows++; 144a64fbb32SBarry Smith } 145a64fbb32SBarry Smith c->nrows[i] = nrows; 146b1d57f15SBarry Smith ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr); 147b1d57f15SBarry Smith ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr); 14852e6d16bSBarry Smith ierr = PetscLogObjectMemory(c,2*(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr); 149a64fbb32SBarry Smith nrows = 0; 1506eaac0f3SBarry Smith for (j=0; j<M; j++) { 151a64fbb32SBarry Smith if (rowhit[j]) { 152a64fbb32SBarry Smith c->rows[i][nrows] = j; 153a64fbb32SBarry Smith c->columnsforrow[i][nrows] = rowhit[j] - 1; 154a64fbb32SBarry Smith nrows++; 155a64fbb32SBarry Smith } 156a64fbb32SBarry Smith } 157a64fbb32SBarry Smith } else {/*-------------------------------------------------------------------------------*/ 158f158e583SBarry Smith /* slow version, using rowhit as a linked list */ 159b1d57f15SBarry Smith PetscInt currentcol,fm,mfm; 1606eaac0f3SBarry Smith rowhit[M] = M; 161a64fbb32SBarry Smith nrows = 0; 162a64fbb32SBarry Smith /* loop over columns*/ 1636eaac0f3SBarry Smith for (j=0; j<nctot; j++) { 164b8f8c88eSHong Zhang if (ctype == IS_COLORING_GHOSTED) { 165b8f8c88eSHong Zhang col = ltog[cols[j]]; 166b8f8c88eSHong Zhang } else { 1676eaac0f3SBarry Smith col = cols[j]; 168b8f8c88eSHong Zhang } 1696eaac0f3SBarry Smith if (col >= cstart && col < cend) { 1706eaac0f3SBarry Smith /* column is in diagonal block of matrix */ 1716eaac0f3SBarry Smith rows = A_cj + A_ci[col-cstart]; 1726eaac0f3SBarry Smith m = A_ci[col-cstart+1] - A_ci[col-cstart]; 1736eaac0f3SBarry Smith } else { 174aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 1750f5bd95cSBarry Smith ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr); 176fa46199cSSatish Balay colb --; 177b3d2dc96SSatish Balay #else 1786eaac0f3SBarry Smith colb = aij->colmap[col] - 1; 179b3d2dc96SSatish Balay #endif 1806eaac0f3SBarry Smith if (colb == -1) { 1816eaac0f3SBarry Smith m = 0; 1826eaac0f3SBarry Smith } else { 1836eaac0f3SBarry Smith rows = B_cj + B_ci[colb]; 1846eaac0f3SBarry Smith m = B_ci[colb+1] - B_ci[colb]; 1856eaac0f3SBarry Smith } 1866eaac0f3SBarry Smith } 187b8f8c88eSHong Zhang 188a64fbb32SBarry Smith /* loop over columns marking them in rowhit */ 1896eaac0f3SBarry Smith fm = M; /* fm points to first entry in linked list */ 190a64fbb32SBarry Smith for (k=0; k<m; k++) { 191a64fbb32SBarry Smith currentcol = *rows++; 192a64fbb32SBarry Smith /* is it already in the list? */ 193a64fbb32SBarry Smith do { 194a64fbb32SBarry Smith mfm = fm; 195a64fbb32SBarry Smith fm = rowhit[fm]; 196a64fbb32SBarry Smith } while (fm < currentcol); 197a64fbb32SBarry Smith /* not in list so add it */ 198a64fbb32SBarry Smith if (fm != currentcol) { 199a64fbb32SBarry Smith nrows++; 200a64fbb32SBarry Smith columnsforrow[currentcol] = col; 201a64fbb32SBarry Smith /* next three lines insert new entry into linked list */ 202a64fbb32SBarry Smith rowhit[mfm] = currentcol; 203a64fbb32SBarry Smith rowhit[currentcol] = fm; 204a64fbb32SBarry Smith fm = currentcol; 205a64fbb32SBarry Smith /* fm points to present position in list since we know the columns are sorted */ 206a64fbb32SBarry Smith } else { 20729bbc08cSBarry Smith SETERRQ(PETSC_ERR_PLIB,"Invalid coloring of matrix detected"); 208a64fbb32SBarry Smith } 209a64fbb32SBarry Smith } 210a64fbb32SBarry Smith } 211a64fbb32SBarry Smith c->nrows[i] = nrows; 212b1d57f15SBarry Smith ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr); 213b1d57f15SBarry Smith ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr); 21452e6d16bSBarry Smith ierr = PetscLogObjectMemory(c,(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr); 215a64fbb32SBarry Smith /* now store the linked list of rows into c->rows[i] */ 216a64fbb32SBarry Smith nrows = 0; 2176eaac0f3SBarry Smith fm = rowhit[M]; 218a64fbb32SBarry Smith do { 219a64fbb32SBarry Smith c->rows[i][nrows] = fm; 220a64fbb32SBarry Smith c->columnsforrow[i][nrows++] = columnsforrow[fm]; 221a64fbb32SBarry Smith fm = rowhit[fm]; 2226eaac0f3SBarry Smith } while (fm < M); 2236eaac0f3SBarry Smith } /* ---------------------------------------------------------------------------------------*/ 224606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 2256eaac0f3SBarry Smith } 22630b34957SBarry Smith 22730b34957SBarry Smith /* Optimize by adding the vscale, and scaleforrow[][] fields */ 22830b34957SBarry Smith /* 22930b34957SBarry Smith vscale will contain the "diagonal" on processor scalings followed by the off processor 23030b34957SBarry Smith */ 2318ee2e534SBarry Smith if (ctype == IS_COLORING_GLOBAL) { 232d0f46423SBarry Smith ierr = VecCreateGhost(((PetscObject)mat)->comm,aij->A->rmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&c->vscale);CHKERRQ(ierr); 233b1d57f15SBarry Smith ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr); 23430b34957SBarry Smith for (k=0; k<c->ncolors; k++) { 235b1d57f15SBarry Smith ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr); 23630b34957SBarry Smith for (l=0; l<c->nrows[k]; l++) { 23730b34957SBarry Smith col = c->columnsforrow[k][l]; 23830b34957SBarry Smith if (col >= cstart && col < cend) { 23930b34957SBarry Smith /* column is in diagonal block of matrix */ 24030b34957SBarry Smith colb = col - cstart; 24130b34957SBarry Smith } else { 24230b34957SBarry Smith /* column is in "off-processor" part */ 24330b34957SBarry Smith #if defined (PETSC_USE_CTABLE) 24430b34957SBarry Smith ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr); 24530b34957SBarry Smith colb --; 24630b34957SBarry Smith #else 24730b34957SBarry Smith colb = aij->colmap[col] - 1; 24830b34957SBarry Smith #endif 24930b34957SBarry Smith colb += cend - cstart; 25030b34957SBarry Smith } 25130b34957SBarry Smith c->vscaleforrow[k][l] = colb; 25230b34957SBarry Smith } 25330b34957SBarry Smith } 254b8f8c88eSHong Zhang } else if (ctype == IS_COLORING_GHOSTED) { 255b8f8c88eSHong Zhang /* Get gtol mapping */ 256d0f46423SBarry Smith PetscInt N = mat->cmap->N, *gtol; 257b8f8c88eSHong Zhang ierr = PetscMalloc((N+1)*sizeof(PetscInt),>ol);CHKERRQ(ierr); 258b8f8c88eSHong Zhang for (i=0; i<N; i++) gtol[i] = -1; 259b8f8c88eSHong Zhang for (i=0; i<map->n; i++) gtol[ltog[i]] = i; 260b8f8c88eSHong Zhang 261b8f8c88eSHong Zhang c->vscale = 0; /* will be created in MatFDColoringApply() */ 262b8f8c88eSHong Zhang ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr); 263b8f8c88eSHong Zhang for (k=0; k<c->ncolors; k++) { 264b8f8c88eSHong Zhang ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr); 265b8f8c88eSHong Zhang for (l=0; l<c->nrows[k]; l++) { 266b8f8c88eSHong Zhang col = c->columnsforrow[k][l]; /* global column index */ 267b8f8c88eSHong Zhang c->vscaleforrow[k][l] = gtol[col]; /* local column index */ 268b8f8c88eSHong Zhang } 269b8f8c88eSHong Zhang } 270b8f8c88eSHong Zhang ierr = PetscFree(gtol);CHKERRQ(ierr); 271b8f8c88eSHong Zhang } 272b9617806SBarry Smith ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 27330b34957SBarry Smith 274606d414cSSatish Balay ierr = PetscFree(rowhit);CHKERRQ(ierr); 275606d414cSSatish Balay ierr = PetscFree(columnsforrow);CHKERRQ(ierr); 2768f7157efSSatish Balay ierr = MatRestoreColumnIJ_SeqAIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr); 2778f7157efSSatish Balay ierr = MatRestoreColumnIJ_SeqAIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr); 2783a40ed3dSBarry Smith PetscFunctionReturn(0); 279a64fbb32SBarry Smith } 280a64fbb32SBarry Smith 281b9617806SBarry Smith 282b9617806SBarry Smith 283b9617806SBarry Smith 284b9617806SBarry Smith 285b9617806SBarry Smith 286