1be1d678aSKris Buschelman #define PETSCMAT_DLL 2a64fbb32SBarry Smith 37c4f633dSBarry Smith #include "../src/mat/impls/aij/mpi/mpiaij.h" 4a64fbb32SBarry Smith 5*09573ac7SBarry 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; 14f6d58c54SBarry Smith PetscInt i,n,nrows,j,k,m,*rows = 0,*A_ci,*A_cj,ncols,col; 155d0c19d7SBarry Smith const PetscInt *is; 16b1d57f15SBarry Smith PetscInt nis = iscoloring->n,nctot,*cols,*B_ci,*B_cj; 17f6d58c54SBarry Smith PetscInt *rowhit,M,cstart,cend,colb; 18b1d57f15SBarry Smith PetscInt *columnsforrow,l; 19b9617806SBarry Smith IS *isa; 20ace3abfcSBarry Smith PetscBool done,flg; 21784ac674SJed Brown ISLocalToGlobalMapping map = mat->cmapping; 22e60d7cb6SSatish Balay PetscInt *ltog = (map ? map->indices : (PetscInt*) PETSC_NULL) ,ctype=c->ctype; 23a64fbb32SBarry Smith 243a40ed3dSBarry Smith PetscFunctionBegin; 25e7e72b3dSBarry Smith if (!mat->assembled) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled first; MatAssemblyBegin/End();"); 26e7e72b3dSBarry Smith 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"); 27522c5e43SBarry Smith 28b9617806SBarry Smith ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 293acb8795SBarry Smith 30f6d58c54SBarry Smith M = mat->rmap->n; 31f6d58c54SBarry Smith cstart = mat->cmap->rstart; 32f6d58c54SBarry Smith cend = mat->cmap->rend; 33f6d58c54SBarry Smith c->M = mat->rmap->N; /* set the global rows and columns and local rows */ 34f6d58c54SBarry Smith c->N = mat->cmap->N; 35f6d58c54SBarry Smith c->m = mat->rmap->n; 36f6d58c54SBarry Smith c->rstart = mat->rmap->rstart; 37005c665bSBarry Smith 38a64fbb32SBarry Smith c->ncolors = nis; 39b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 40b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr); 41b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 42b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->rows);CHKERRQ(ierr); 43b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columnsforrow);CHKERRQ(ierr); 4452e6d16bSBarry Smith ierr = PetscLogObjectMemory(c,5*nis*sizeof(PetscInt));CHKERRQ(ierr); 456eaac0f3SBarry Smith 466eaac0f3SBarry Smith /* Allow access to data structures of local part of matrix */ 476eaac0f3SBarry Smith if (!aij->colmap) { 486eaac0f3SBarry Smith ierr = CreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr); 496eaac0f3SBarry Smith } 503acb8795SBarry Smith ierr = MatGetColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr); 513acb8795SBarry Smith ierr = MatGetColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr); 526eaac0f3SBarry Smith 53b1d57f15SBarry Smith ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr); 54b1d57f15SBarry Smith ierr = PetscMalloc((M+1)*sizeof(PetscInt),&columnsforrow);CHKERRQ(ierr); 556eaac0f3SBarry Smith 56a64fbb32SBarry Smith for (i=0; i<nis; i++) { 57b9b97703SBarry Smith ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 58a64fbb32SBarry Smith ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 59a64fbb32SBarry Smith c->ncolumns[i] = n; 60a64fbb32SBarry Smith if (n) { 61b1d57f15SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr); 6252e6d16bSBarry Smith ierr = PetscLogObjectMemory(c,n*sizeof(PetscInt));CHKERRQ(ierr); 63b1d57f15SBarry Smith ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr); 64a64fbb32SBarry Smith } else { 65a64fbb32SBarry Smith c->columns[i] = 0; 66a64fbb32SBarry Smith } 67a64fbb32SBarry Smith 688ee2e534SBarry Smith if (ctype == IS_COLORING_GLOBAL){ 696eaac0f3SBarry Smith /* Determine the total (parallel) number of columns of this color */ 707adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)mat)->comm,&size);CHKERRQ(ierr); 71687f1162SBarry Smith ierr = PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);CHKERRQ(ierr); 72b8f8c88eSHong Zhang 73e44c0bd4SBarry Smith nn = PetscMPIIntCast(n); 747adad957SLisandro Dalcin ierr = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,((PetscObject)mat)->comm);CHKERRQ(ierr); 756eaac0f3SBarry Smith nctot = 0; for (j=0; j<size; j++) {nctot += ncolsonproc[j];} 763a7fca6bSBarry Smith if (!nctot) { 77ae15b995SBarry Smith ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr); 783a7fca6bSBarry Smith } 796eaac0f3SBarry Smith 806eaac0f3SBarry Smith disp[0] = 0; 816eaac0f3SBarry Smith for (j=1; j<size; j++) { 826eaac0f3SBarry Smith disp[j] = disp[j-1] + ncolsonproc[j-1]; 836eaac0f3SBarry Smith } 846eaac0f3SBarry Smith 856eaac0f3SBarry Smith /* Get complete list of columns for color on each processor */ 86b1d57f15SBarry Smith ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 875d0c19d7SBarry Smith ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,((PetscObject)mat)->comm);CHKERRQ(ierr); 881d79065fSBarry Smith ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr); 89b8f8c88eSHong Zhang } else if (ctype == IS_COLORING_GHOSTED){ 90b8f8c88eSHong Zhang /* Determine local number of columns of this color on this process, including ghost points */ 91b8f8c88eSHong Zhang nctot = n; 92b8f8c88eSHong Zhang ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 93b8f8c88eSHong Zhang ierr = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr); 94b8f8c88eSHong Zhang } else { 95e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type"); 96b8f8c88eSHong Zhang } 976eaac0f3SBarry Smith 986eaac0f3SBarry Smith /* 996eaac0f3SBarry Smith Mark all rows affect by these columns 1006eaac0f3SBarry Smith */ 101b8f8c88eSHong Zhang /* Temporary option to allow for debugging/testing */ 10290d69ab7SBarry Smith flg = PETSC_FALSE; 103acfcf0e5SJed Brown ierr = PetscOptionsGetBool(PETSC_NULL,"-matfdcoloring_slow",&flg,PETSC_NULL);CHKERRQ(ierr); 104f158e583SBarry Smith if (!flg) {/*-----------------------------------------------------------------------------*/ 105f158e583SBarry Smith /* crude, fast version */ 106b1d57f15SBarry Smith ierr = PetscMemzero(rowhit,M*sizeof(PetscInt));CHKERRQ(ierr); 107a64fbb32SBarry Smith /* loop over columns*/ 1086eaac0f3SBarry Smith for (j=0; j<nctot; j++) { 109b8f8c88eSHong Zhang if (ctype == IS_COLORING_GHOSTED) { 110b8f8c88eSHong Zhang col = ltog[cols[j]]; 111b8f8c88eSHong Zhang } else { 1126eaac0f3SBarry Smith col = cols[j]; 113b8f8c88eSHong Zhang } 1146eaac0f3SBarry Smith if (col >= cstart && col < cend) { 1156eaac0f3SBarry Smith /* column is in diagonal block of matrix */ 1166eaac0f3SBarry Smith rows = A_cj + A_ci[col-cstart]; 1176eaac0f3SBarry Smith m = A_ci[col-cstart+1] - A_ci[col-cstart]; 1186eaac0f3SBarry Smith } else { 119aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 120cb9801acSJed Brown ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr); 121fa46199cSSatish Balay colb --; 122b3d2dc96SSatish Balay #else 1236eaac0f3SBarry Smith colb = aij->colmap[col] - 1; 124b3d2dc96SSatish Balay #endif 1256eaac0f3SBarry Smith if (colb == -1) { 1266eaac0f3SBarry Smith m = 0; 1276eaac0f3SBarry Smith } else { 1286eaac0f3SBarry Smith rows = B_cj + B_ci[colb]; 1296eaac0f3SBarry Smith m = B_ci[colb+1] - B_ci[colb]; 1306eaac0f3SBarry Smith } 1316eaac0f3SBarry Smith } 132a64fbb32SBarry Smith /* loop over columns marking them in rowhit */ 133a64fbb32SBarry Smith for (k=0; k<m; k++) { 134a64fbb32SBarry Smith rowhit[*rows++] = col + 1; 135a64fbb32SBarry Smith } 136a64fbb32SBarry Smith } 1376eaac0f3SBarry Smith 138a64fbb32SBarry Smith /* count the number of hits */ 139a64fbb32SBarry Smith nrows = 0; 1406eaac0f3SBarry Smith for (j=0; j<M; j++) { 141a64fbb32SBarry Smith if (rowhit[j]) nrows++; 142a64fbb32SBarry Smith } 143a64fbb32SBarry Smith c->nrows[i] = nrows; 144b1d57f15SBarry Smith ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr); 145b1d57f15SBarry Smith ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr); 14652e6d16bSBarry Smith ierr = PetscLogObjectMemory(c,2*(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr); 147a64fbb32SBarry Smith nrows = 0; 1486eaac0f3SBarry Smith for (j=0; j<M; j++) { 149a64fbb32SBarry Smith if (rowhit[j]) { 150a64fbb32SBarry Smith c->rows[i][nrows] = j; 151a64fbb32SBarry Smith c->columnsforrow[i][nrows] = rowhit[j] - 1; 152a64fbb32SBarry Smith nrows++; 153a64fbb32SBarry Smith } 154a64fbb32SBarry Smith } 155a64fbb32SBarry Smith } else {/*-------------------------------------------------------------------------------*/ 156f158e583SBarry Smith /* slow version, using rowhit as a linked list */ 157b1d57f15SBarry Smith PetscInt currentcol,fm,mfm; 1586eaac0f3SBarry Smith rowhit[M] = M; 159a64fbb32SBarry Smith nrows = 0; 160a64fbb32SBarry Smith /* loop over columns*/ 1616eaac0f3SBarry Smith for (j=0; j<nctot; j++) { 162b8f8c88eSHong Zhang if (ctype == IS_COLORING_GHOSTED) { 163b8f8c88eSHong Zhang col = ltog[cols[j]]; 164b8f8c88eSHong Zhang } else { 1656eaac0f3SBarry Smith col = cols[j]; 166b8f8c88eSHong Zhang } 1676eaac0f3SBarry Smith if (col >= cstart && col < cend) { 1686eaac0f3SBarry Smith /* column is in diagonal block of matrix */ 1696eaac0f3SBarry Smith rows = A_cj + A_ci[col-cstart]; 1706eaac0f3SBarry Smith m = A_ci[col-cstart+1] - A_ci[col-cstart]; 1716eaac0f3SBarry Smith } else { 172aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 1730f5bd95cSBarry Smith ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr); 174fa46199cSSatish Balay colb --; 175b3d2dc96SSatish Balay #else 1766eaac0f3SBarry Smith colb = aij->colmap[col] - 1; 177b3d2dc96SSatish Balay #endif 1786eaac0f3SBarry Smith if (colb == -1) { 1796eaac0f3SBarry Smith m = 0; 1806eaac0f3SBarry Smith } else { 1816eaac0f3SBarry Smith rows = B_cj + B_ci[colb]; 1826eaac0f3SBarry Smith m = B_ci[colb+1] - B_ci[colb]; 1836eaac0f3SBarry Smith } 1846eaac0f3SBarry Smith } 185b8f8c88eSHong Zhang 186a64fbb32SBarry Smith /* loop over columns marking them in rowhit */ 1876eaac0f3SBarry Smith fm = M; /* fm points to first entry in linked list */ 188a64fbb32SBarry Smith for (k=0; k<m; k++) { 189a64fbb32SBarry Smith currentcol = *rows++; 190a64fbb32SBarry Smith /* is it already in the list? */ 191a64fbb32SBarry Smith do { 192a64fbb32SBarry Smith mfm = fm; 193a64fbb32SBarry Smith fm = rowhit[fm]; 194a64fbb32SBarry Smith } while (fm < currentcol); 195a64fbb32SBarry Smith /* not in list so add it */ 196a64fbb32SBarry Smith if (fm != currentcol) { 197a64fbb32SBarry Smith nrows++; 198a64fbb32SBarry Smith columnsforrow[currentcol] = col; 199a64fbb32SBarry Smith /* next three lines insert new entry into linked list */ 200a64fbb32SBarry Smith rowhit[mfm] = currentcol; 201a64fbb32SBarry Smith rowhit[currentcol] = fm; 202a64fbb32SBarry Smith fm = currentcol; 203a64fbb32SBarry Smith /* fm points to present position in list since we know the columns are sorted */ 204a64fbb32SBarry Smith } else { 205e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid coloring of matrix detected"); 206a64fbb32SBarry Smith } 207a64fbb32SBarry Smith } 208a64fbb32SBarry Smith } 209a64fbb32SBarry Smith c->nrows[i] = nrows; 210b1d57f15SBarry Smith ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr); 211b1d57f15SBarry Smith ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr); 21252e6d16bSBarry Smith ierr = PetscLogObjectMemory(c,(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr); 213a64fbb32SBarry Smith /* now store the linked list of rows into c->rows[i] */ 214a64fbb32SBarry Smith nrows = 0; 2156eaac0f3SBarry Smith fm = rowhit[M]; 216a64fbb32SBarry Smith do { 217a64fbb32SBarry Smith c->rows[i][nrows] = fm; 218a64fbb32SBarry Smith c->columnsforrow[i][nrows++] = columnsforrow[fm]; 219a64fbb32SBarry Smith fm = rowhit[fm]; 2206eaac0f3SBarry Smith } while (fm < M); 2216eaac0f3SBarry Smith } /* ---------------------------------------------------------------------------------------*/ 222606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 2236eaac0f3SBarry Smith } 22430b34957SBarry Smith 22530b34957SBarry Smith /* Optimize by adding the vscale, and scaleforrow[][] fields */ 22630b34957SBarry Smith /* 22730b34957SBarry Smith vscale will contain the "diagonal" on processor scalings followed by the off processor 22830b34957SBarry Smith */ 2298ee2e534SBarry Smith if (ctype == IS_COLORING_GLOBAL) { 230d0f46423SBarry Smith ierr = VecCreateGhost(((PetscObject)mat)->comm,aij->A->rmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&c->vscale);CHKERRQ(ierr); 231b1d57f15SBarry Smith ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr); 23230b34957SBarry Smith for (k=0; k<c->ncolors; k++) { 233b1d57f15SBarry Smith ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr); 23430b34957SBarry Smith for (l=0; l<c->nrows[k]; l++) { 23530b34957SBarry Smith col = c->columnsforrow[k][l]; 23630b34957SBarry Smith if (col >= cstart && col < cend) { 23730b34957SBarry Smith /* column is in diagonal block of matrix */ 23830b34957SBarry Smith colb = col - cstart; 23930b34957SBarry Smith } else { 24030b34957SBarry Smith /* column is in "off-processor" part */ 24130b34957SBarry Smith #if defined (PETSC_USE_CTABLE) 24230b34957SBarry Smith ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr); 24330b34957SBarry Smith colb --; 24430b34957SBarry Smith #else 24530b34957SBarry Smith colb = aij->colmap[col] - 1; 24630b34957SBarry Smith #endif 24730b34957SBarry Smith colb += cend - cstart; 24830b34957SBarry Smith } 24930b34957SBarry Smith c->vscaleforrow[k][l] = colb; 25030b34957SBarry Smith } 25130b34957SBarry Smith } 252b8f8c88eSHong Zhang } else if (ctype == IS_COLORING_GHOSTED) { 253b8f8c88eSHong Zhang /* Get gtol mapping */ 254d0f46423SBarry Smith PetscInt N = mat->cmap->N, *gtol; 255b8f8c88eSHong Zhang ierr = PetscMalloc((N+1)*sizeof(PetscInt),>ol);CHKERRQ(ierr); 256b8f8c88eSHong Zhang for (i=0; i<N; i++) gtol[i] = -1; 257b8f8c88eSHong Zhang for (i=0; i<map->n; i++) gtol[ltog[i]] = i; 258b8f8c88eSHong Zhang 259b8f8c88eSHong Zhang c->vscale = 0; /* will be created in MatFDColoringApply() */ 260b8f8c88eSHong Zhang ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr); 261b8f8c88eSHong Zhang for (k=0; k<c->ncolors; k++) { 262b8f8c88eSHong Zhang ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr); 263b8f8c88eSHong Zhang for (l=0; l<c->nrows[k]; l++) { 264b8f8c88eSHong Zhang col = c->columnsforrow[k][l]; /* global column index */ 265b8f8c88eSHong Zhang c->vscaleforrow[k][l] = gtol[col]; /* local column index */ 266b8f8c88eSHong Zhang } 267b8f8c88eSHong Zhang } 268b8f8c88eSHong Zhang ierr = PetscFree(gtol);CHKERRQ(ierr); 269b8f8c88eSHong Zhang } 270b9617806SBarry Smith ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 27130b34957SBarry Smith 272606d414cSSatish Balay ierr = PetscFree(rowhit);CHKERRQ(ierr); 273606d414cSSatish Balay ierr = PetscFree(columnsforrow);CHKERRQ(ierr); 2743acb8795SBarry Smith ierr = MatRestoreColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr); 2753acb8795SBarry Smith ierr = MatRestoreColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr); 2763a40ed3dSBarry Smith PetscFunctionReturn(0); 277a64fbb32SBarry Smith } 278a64fbb32SBarry Smith 279b9617806SBarry Smith 280b9617806SBarry Smith 281b9617806SBarry Smith 282b9617806SBarry Smith 283b9617806SBarry Smith 284