1be1d678aSKris Buschelman #define PETSCMAT_DLL 2a64fbb32SBarry Smith 37c4f633dSBarry 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" 9*3acb8795SBarry Smith /* 10*3acb8795SBarry Smith This routine is shared by AIJ and BAIJ matrices, since it operators only on the nonzero structure of the elements or blocks. 11*3acb8795SBarry Smith This is why it has the ugly code with the MatGetBlockSize() 12*3acb8795SBarry Smith */ 13dfbe8321SBarry Smith PetscErrorCode MatFDColoringCreate_MPIAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c) 14a64fbb32SBarry Smith { 156eaac0f3SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)mat->data; 166849ba73SBarry Smith PetscErrorCode ierr; 17b1d57f15SBarry Smith PetscMPIInt size,*ncolsonproc,*disp,nn; 18*3acb8795SBarry Smith PetscInt bs = 1,i,n,nrows,j,k,m,*rows = 0,*A_ci,*A_cj,ncols,col; 195d0c19d7SBarry Smith const PetscInt *is; 20b1d57f15SBarry Smith PetscInt nis = iscoloring->n,nctot,*cols,*B_ci,*B_cj; 21d0f46423SBarry Smith PetscInt *rowhit,M = mat->rmap->n,cstart = mat->cmap->rstart,cend = mat->cmap->rend,colb; 22b1d57f15SBarry Smith PetscInt *columnsforrow,l; 23b9617806SBarry Smith IS *isa; 24f1af5d2fSBarry Smith PetscTruth done,flg; 25b8f8c88eSHong Zhang ISLocalToGlobalMapping map = mat->mapping; 26e60d7cb6SSatish Balay PetscInt *ltog = (map ? map->indices : (PetscInt*) PETSC_NULL) ,ctype=c->ctype; 27*3acb8795SBarry Smith PetscTruth flg1,flg2; 28a64fbb32SBarry Smith 293a40ed3dSBarry Smith PetscFunctionBegin; 30522c5e43SBarry Smith if (!mat->assembled) { 3129bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled first; MatAssemblyBegin/End();"); 32522c5e43SBarry Smith } 337cc688d7SBarry 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"); 34522c5e43SBarry Smith 35b9617806SBarry Smith ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 36*3acb8795SBarry Smith 37*3acb8795SBarry Smith /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */ 38*3acb8795SBarry Smith ierr = PetscTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr); 39*3acb8795SBarry Smith ierr = PetscTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr); 40*3acb8795SBarry Smith if (flg1 || flg2) { 41*3acb8795SBarry Smith ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 42*3acb8795SBarry Smith } 43*3acb8795SBarry Smith 44*3acb8795SBarry Smith c->M = mat->rmap->N/bs; /* set the global rows and columns and local rows */ 45*3acb8795SBarry Smith c->N = mat->cmap->N/bs; 46*3acb8795SBarry Smith c->m = mat->rmap->n/bs; 47*3acb8795SBarry Smith c->rstart = mat->rmap->rstart/bs; 48005c665bSBarry Smith 49a64fbb32SBarry Smith c->ncolors = nis; 50b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 51b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr); 52b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 53b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->rows);CHKERRQ(ierr); 54b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columnsforrow);CHKERRQ(ierr); 5552e6d16bSBarry Smith ierr = PetscLogObjectMemory(c,5*nis*sizeof(PetscInt));CHKERRQ(ierr); 566eaac0f3SBarry Smith 576eaac0f3SBarry Smith /* Allow access to data structures of local part of matrix */ 586eaac0f3SBarry Smith if (!aij->colmap) { 596eaac0f3SBarry Smith ierr = CreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr); 606eaac0f3SBarry Smith } 61*3acb8795SBarry Smith ierr = MatGetColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr); 62*3acb8795SBarry Smith ierr = MatGetColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr); 636eaac0f3SBarry Smith 64b1d57f15SBarry Smith ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr); 65b1d57f15SBarry Smith ierr = PetscMalloc((M+1)*sizeof(PetscInt),&columnsforrow);CHKERRQ(ierr); 666eaac0f3SBarry Smith 67a64fbb32SBarry Smith for (i=0; i<nis; i++) { 68b9b97703SBarry Smith ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 69a64fbb32SBarry Smith ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 70a64fbb32SBarry Smith c->ncolumns[i] = n; 71a64fbb32SBarry Smith if (n) { 72b1d57f15SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr); 7352e6d16bSBarry Smith ierr = PetscLogObjectMemory(c,n*sizeof(PetscInt));CHKERRQ(ierr); 74b1d57f15SBarry Smith ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr); 75a64fbb32SBarry Smith } else { 76a64fbb32SBarry Smith c->columns[i] = 0; 77a64fbb32SBarry Smith } 78a64fbb32SBarry Smith 798ee2e534SBarry Smith if (ctype == IS_COLORING_GLOBAL){ 806eaac0f3SBarry Smith /* Determine the total (parallel) number of columns of this color */ 817adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)mat)->comm,&size);CHKERRQ(ierr); 82b8f8c88eSHong Zhang ierr = PetscMalloc(2*size*sizeof(PetscInt*),&ncolsonproc);CHKERRQ(ierr); 83b8f8c88eSHong Zhang disp = ncolsonproc + size; 84b8f8c88eSHong Zhang 85e44c0bd4SBarry Smith nn = PetscMPIIntCast(n); 867adad957SLisandro Dalcin ierr = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,((PetscObject)mat)->comm);CHKERRQ(ierr); 876eaac0f3SBarry Smith nctot = 0; for (j=0; j<size; j++) {nctot += ncolsonproc[j];} 883a7fca6bSBarry Smith if (!nctot) { 89ae15b995SBarry Smith ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr); 903a7fca6bSBarry Smith } 916eaac0f3SBarry Smith 926eaac0f3SBarry Smith disp[0] = 0; 936eaac0f3SBarry Smith for (j=1; j<size; j++) { 946eaac0f3SBarry Smith disp[j] = disp[j-1] + ncolsonproc[j-1]; 956eaac0f3SBarry Smith } 966eaac0f3SBarry Smith 976eaac0f3SBarry Smith /* Get complete list of columns for color on each processor */ 98b1d57f15SBarry Smith ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 995d0c19d7SBarry Smith ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,((PetscObject)mat)->comm);CHKERRQ(ierr); 100b8f8c88eSHong Zhang ierr = PetscFree(ncolsonproc);CHKERRQ(ierr); 101b8f8c88eSHong Zhang } else if (ctype == IS_COLORING_GHOSTED){ 102b8f8c88eSHong Zhang /* Determine local number of columns of this color on this process, including ghost points */ 103b8f8c88eSHong Zhang nctot = n; 104b8f8c88eSHong Zhang ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 105b8f8c88eSHong Zhang ierr = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr); 106b8f8c88eSHong Zhang } else { 107b8f8c88eSHong Zhang SETERRQ(PETSC_ERR_SUP,"Not provided for this MatFDColoring type"); 108b8f8c88eSHong Zhang } 1096eaac0f3SBarry Smith 1106eaac0f3SBarry Smith /* 1116eaac0f3SBarry Smith Mark all rows affect by these columns 1126eaac0f3SBarry Smith */ 113b8f8c88eSHong Zhang /* Temporary option to allow for debugging/testing */ 11490d69ab7SBarry Smith flg = PETSC_FALSE; 11590d69ab7SBarry Smith ierr = PetscOptionsGetTruth(PETSC_NULL,"-matfdcoloring_slow",&flg,PETSC_NULL);CHKERRQ(ierr); 116f158e583SBarry Smith if (!flg) {/*-----------------------------------------------------------------------------*/ 117f158e583SBarry Smith /* crude, fast version */ 118b1d57f15SBarry Smith ierr = PetscMemzero(rowhit,M*sizeof(PetscInt));CHKERRQ(ierr); 119a64fbb32SBarry Smith /* loop over columns*/ 1206eaac0f3SBarry Smith for (j=0; j<nctot; j++) { 121b8f8c88eSHong Zhang if (ctype == IS_COLORING_GHOSTED) { 122b8f8c88eSHong Zhang col = ltog[cols[j]]; 123b8f8c88eSHong Zhang } else { 1246eaac0f3SBarry Smith col = cols[j]; 125b8f8c88eSHong Zhang } 1266eaac0f3SBarry Smith if (col >= cstart && col < cend) { 1276eaac0f3SBarry Smith /* column is in diagonal block of matrix */ 1286eaac0f3SBarry Smith rows = A_cj + A_ci[col-cstart]; 1296eaac0f3SBarry Smith m = A_ci[col-cstart+1] - A_ci[col-cstart]; 1306eaac0f3SBarry Smith } else { 131aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 1320f5bd95cSBarry Smith ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr) 133fa46199cSSatish Balay colb --; 134b3d2dc96SSatish Balay #else 1356eaac0f3SBarry Smith colb = aij->colmap[col] - 1; 136b3d2dc96SSatish Balay #endif 1376eaac0f3SBarry Smith if (colb == -1) { 1386eaac0f3SBarry Smith m = 0; 1396eaac0f3SBarry Smith } else { 1406eaac0f3SBarry Smith rows = B_cj + B_ci[colb]; 1416eaac0f3SBarry Smith m = B_ci[colb+1] - B_ci[colb]; 1426eaac0f3SBarry Smith } 1436eaac0f3SBarry Smith } 144a64fbb32SBarry Smith /* loop over columns marking them in rowhit */ 145a64fbb32SBarry Smith for (k=0; k<m; k++) { 146a64fbb32SBarry Smith rowhit[*rows++] = col + 1; 147a64fbb32SBarry Smith } 148a64fbb32SBarry Smith } 1496eaac0f3SBarry Smith 150a64fbb32SBarry Smith /* count the number of hits */ 151a64fbb32SBarry Smith nrows = 0; 1526eaac0f3SBarry Smith for (j=0; j<M; j++) { 153a64fbb32SBarry Smith if (rowhit[j]) nrows++; 154a64fbb32SBarry Smith } 155a64fbb32SBarry Smith c->nrows[i] = nrows; 156b1d57f15SBarry Smith ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr); 157b1d57f15SBarry Smith ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr); 15852e6d16bSBarry Smith ierr = PetscLogObjectMemory(c,2*(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr); 159a64fbb32SBarry Smith nrows = 0; 1606eaac0f3SBarry Smith for (j=0; j<M; j++) { 161a64fbb32SBarry Smith if (rowhit[j]) { 162a64fbb32SBarry Smith c->rows[i][nrows] = j; 163a64fbb32SBarry Smith c->columnsforrow[i][nrows] = rowhit[j] - 1; 164a64fbb32SBarry Smith nrows++; 165a64fbb32SBarry Smith } 166a64fbb32SBarry Smith } 167a64fbb32SBarry Smith } else {/*-------------------------------------------------------------------------------*/ 168f158e583SBarry Smith /* slow version, using rowhit as a linked list */ 169b1d57f15SBarry Smith PetscInt currentcol,fm,mfm; 1706eaac0f3SBarry Smith rowhit[M] = M; 171a64fbb32SBarry Smith nrows = 0; 172a64fbb32SBarry Smith /* loop over columns*/ 1736eaac0f3SBarry Smith for (j=0; j<nctot; j++) { 174b8f8c88eSHong Zhang if (ctype == IS_COLORING_GHOSTED) { 175b8f8c88eSHong Zhang col = ltog[cols[j]]; 176b8f8c88eSHong Zhang } else { 1776eaac0f3SBarry Smith col = cols[j]; 178b8f8c88eSHong Zhang } 1796eaac0f3SBarry Smith if (col >= cstart && col < cend) { 1806eaac0f3SBarry Smith /* column is in diagonal block of matrix */ 1816eaac0f3SBarry Smith rows = A_cj + A_ci[col-cstart]; 1826eaac0f3SBarry Smith m = A_ci[col-cstart+1] - A_ci[col-cstart]; 1836eaac0f3SBarry Smith } else { 184aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 1850f5bd95cSBarry Smith ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr); 186fa46199cSSatish Balay colb --; 187b3d2dc96SSatish Balay #else 1886eaac0f3SBarry Smith colb = aij->colmap[col] - 1; 189b3d2dc96SSatish Balay #endif 1906eaac0f3SBarry Smith if (colb == -1) { 1916eaac0f3SBarry Smith m = 0; 1926eaac0f3SBarry Smith } else { 1936eaac0f3SBarry Smith rows = B_cj + B_ci[colb]; 1946eaac0f3SBarry Smith m = B_ci[colb+1] - B_ci[colb]; 1956eaac0f3SBarry Smith } 1966eaac0f3SBarry Smith } 197b8f8c88eSHong Zhang 198a64fbb32SBarry Smith /* loop over columns marking them in rowhit */ 1996eaac0f3SBarry Smith fm = M; /* fm points to first entry in linked list */ 200a64fbb32SBarry Smith for (k=0; k<m; k++) { 201a64fbb32SBarry Smith currentcol = *rows++; 202a64fbb32SBarry Smith /* is it already in the list? */ 203a64fbb32SBarry Smith do { 204a64fbb32SBarry Smith mfm = fm; 205a64fbb32SBarry Smith fm = rowhit[fm]; 206a64fbb32SBarry Smith } while (fm < currentcol); 207a64fbb32SBarry Smith /* not in list so add it */ 208a64fbb32SBarry Smith if (fm != currentcol) { 209a64fbb32SBarry Smith nrows++; 210a64fbb32SBarry Smith columnsforrow[currentcol] = col; 211a64fbb32SBarry Smith /* next three lines insert new entry into linked list */ 212a64fbb32SBarry Smith rowhit[mfm] = currentcol; 213a64fbb32SBarry Smith rowhit[currentcol] = fm; 214a64fbb32SBarry Smith fm = currentcol; 215a64fbb32SBarry Smith /* fm points to present position in list since we know the columns are sorted */ 216a64fbb32SBarry Smith } else { 21729bbc08cSBarry Smith SETERRQ(PETSC_ERR_PLIB,"Invalid coloring of matrix detected"); 218a64fbb32SBarry Smith } 219a64fbb32SBarry Smith } 220a64fbb32SBarry Smith } 221a64fbb32SBarry Smith c->nrows[i] = nrows; 222b1d57f15SBarry Smith ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr); 223b1d57f15SBarry Smith ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr); 22452e6d16bSBarry Smith ierr = PetscLogObjectMemory(c,(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr); 225a64fbb32SBarry Smith /* now store the linked list of rows into c->rows[i] */ 226a64fbb32SBarry Smith nrows = 0; 2276eaac0f3SBarry Smith fm = rowhit[M]; 228a64fbb32SBarry Smith do { 229a64fbb32SBarry Smith c->rows[i][nrows] = fm; 230a64fbb32SBarry Smith c->columnsforrow[i][nrows++] = columnsforrow[fm]; 231a64fbb32SBarry Smith fm = rowhit[fm]; 2326eaac0f3SBarry Smith } while (fm < M); 2336eaac0f3SBarry Smith } /* ---------------------------------------------------------------------------------------*/ 234606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 2356eaac0f3SBarry Smith } 23630b34957SBarry Smith 23730b34957SBarry Smith /* Optimize by adding the vscale, and scaleforrow[][] fields */ 23830b34957SBarry Smith /* 23930b34957SBarry Smith vscale will contain the "diagonal" on processor scalings followed by the off processor 24030b34957SBarry Smith */ 2418ee2e534SBarry Smith if (ctype == IS_COLORING_GLOBAL) { 242d0f46423SBarry Smith ierr = VecCreateGhost(((PetscObject)mat)->comm,aij->A->rmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&c->vscale);CHKERRQ(ierr); 243b1d57f15SBarry Smith ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr); 24430b34957SBarry Smith for (k=0; k<c->ncolors; k++) { 245b1d57f15SBarry Smith ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr); 24630b34957SBarry Smith for (l=0; l<c->nrows[k]; l++) { 24730b34957SBarry Smith col = c->columnsforrow[k][l]; 24830b34957SBarry Smith if (col >= cstart && col < cend) { 24930b34957SBarry Smith /* column is in diagonal block of matrix */ 25030b34957SBarry Smith colb = col - cstart; 25130b34957SBarry Smith } else { 25230b34957SBarry Smith /* column is in "off-processor" part */ 25330b34957SBarry Smith #if defined (PETSC_USE_CTABLE) 25430b34957SBarry Smith ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr); 25530b34957SBarry Smith colb --; 25630b34957SBarry Smith #else 25730b34957SBarry Smith colb = aij->colmap[col] - 1; 25830b34957SBarry Smith #endif 25930b34957SBarry Smith colb += cend - cstart; 26030b34957SBarry Smith } 26130b34957SBarry Smith c->vscaleforrow[k][l] = colb; 26230b34957SBarry Smith } 26330b34957SBarry Smith } 264b8f8c88eSHong Zhang } else if (ctype == IS_COLORING_GHOSTED) { 265b8f8c88eSHong Zhang /* Get gtol mapping */ 266d0f46423SBarry Smith PetscInt N = mat->cmap->N, *gtol; 267b8f8c88eSHong Zhang ierr = PetscMalloc((N+1)*sizeof(PetscInt),>ol);CHKERRQ(ierr); 268b8f8c88eSHong Zhang for (i=0; i<N; i++) gtol[i] = -1; 269b8f8c88eSHong Zhang for (i=0; i<map->n; i++) gtol[ltog[i]] = i; 270b8f8c88eSHong Zhang 271b8f8c88eSHong Zhang c->vscale = 0; /* will be created in MatFDColoringApply() */ 272b8f8c88eSHong Zhang ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr); 273b8f8c88eSHong Zhang for (k=0; k<c->ncolors; k++) { 274b8f8c88eSHong Zhang ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr); 275b8f8c88eSHong Zhang for (l=0; l<c->nrows[k]; l++) { 276b8f8c88eSHong Zhang col = c->columnsforrow[k][l]; /* global column index */ 277b8f8c88eSHong Zhang c->vscaleforrow[k][l] = gtol[col]; /* local column index */ 278b8f8c88eSHong Zhang } 279b8f8c88eSHong Zhang } 280b8f8c88eSHong Zhang ierr = PetscFree(gtol);CHKERRQ(ierr); 281b8f8c88eSHong Zhang } 282b9617806SBarry Smith ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 28330b34957SBarry Smith 284606d414cSSatish Balay ierr = PetscFree(rowhit);CHKERRQ(ierr); 285606d414cSSatish Balay ierr = PetscFree(columnsforrow);CHKERRQ(ierr); 286*3acb8795SBarry Smith ierr = MatRestoreColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr); 287*3acb8795SBarry Smith ierr = MatRestoreColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr); 2883a40ed3dSBarry Smith PetscFunctionReturn(0); 289a64fbb32SBarry Smith } 290a64fbb32SBarry Smith 291b9617806SBarry Smith 292b9617806SBarry Smith 293b9617806SBarry Smith 294b9617806SBarry Smith 295b9617806SBarry Smith 296