1a64fbb32SBarry Smith 2c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h> 3a64fbb32SBarry Smith 4ab9863d7SBarry Smith extern PetscErrorCode MatCreateColmap_MPIAIJ_Private(Mat); 5a64fbb32SBarry Smith 64a2ae208SSatish Balay #undef __FUNCT__ 7fcd7ac73SHong Zhang #define __FUNCT__ "MatFDColoringApply_MPIAIJ" 8fcd7ac73SHong Zhang PetscErrorCode MatFDColoringApply_MPIAIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 9fcd7ac73SHong Zhang { 10fcd7ac73SHong Zhang PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void*))coloring->f; 11fcd7ac73SHong Zhang PetscErrorCode ierr; 128bc97078SHong Zhang PetscInt k,cstart,cend,l,row,col,colb,nz; 13*9e917edbSHong Zhang PetscScalar dx=0.0,*y,*xx,*w3_array; 14fcd7ac73SHong Zhang PetscScalar *vscale_array; 15fcd7ac73SHong Zhang PetscReal epsilon=coloring->error_rel,umin=coloring->umin,unorm; 168bc97078SHong Zhang Vec w1=coloring->w1,w2=coloring->w2,w3,vscale=coloring->vscale; 17fcd7ac73SHong Zhang void *fctx=coloring->fctx; 18fcd7ac73SHong Zhang PetscBool flg=PETSC_FALSE; 1974d3cef9SHong Zhang PetscInt ctype=coloring->ctype,nxloc; 20f5aae955SHong Zhang Mat_MPIAIJ *aij=(Mat_MPIAIJ*)J->data; 21573f477fSHong Zhang MatEntry *Jentry=coloring->matentry; 228bc97078SHong Zhang const PetscInt ncolors=coloring->ncolors,*ncolumns=coloring->ncolumns,*nrows=coloring->nrows; 23fcd7ac73SHong Zhang 24fcd7ac73SHong Zhang PetscFunctionBegin; 25fcd7ac73SHong Zhang ierr = MatSetUnfactored(J);CHKERRQ(ierr); 26fcd7ac73SHong Zhang ierr = PetscOptionsGetBool(NULL,"-mat_fd_coloring_dont_rezero",&flg,NULL);CHKERRQ(ierr); 27fcd7ac73SHong Zhang if (flg) { 28fcd7ac73SHong Zhang ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr); 29fcd7ac73SHong Zhang } else { 30fcd7ac73SHong Zhang PetscBool assembled; 31fcd7ac73SHong Zhang ierr = MatAssembled(J,&assembled);CHKERRQ(ierr); 32fcd7ac73SHong Zhang if (assembled) { 33fcd7ac73SHong Zhang ierr = MatZeroEntries(J);CHKERRQ(ierr); 34fcd7ac73SHong Zhang } 35fcd7ac73SHong Zhang } 36fcd7ac73SHong Zhang 37*9e917edbSHong Zhang /* create vscale for storing dx */ 38*9e917edbSHong Zhang if (!vscale) { 39*9e917edbSHong Zhang if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') { 40*9e917edbSHong Zhang ierr = VecCreateGhost(PetscObjectComm((PetscObject)J),J->cmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&vscale);CHKERRQ(ierr); 41*9e917edbSHong Zhang } else if (ctype == IS_COLORING_GHOSTED) { 428bc97078SHong Zhang ierr = VecDuplicate(x1,&vscale);CHKERRQ(ierr); 43*9e917edbSHong Zhang } 448bc97078SHong Zhang coloring->vscale = vscale; 45fcd7ac73SHong Zhang } 46fcd7ac73SHong Zhang 478bc97078SHong Zhang /* (1) Set w1 = F(x1) */ 48fcd7ac73SHong Zhang if (!coloring->fset) { 49fcd7ac73SHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 50f6af9589SHong Zhang ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr); 51fcd7ac73SHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 52fcd7ac73SHong Zhang } else { 53fcd7ac73SHong Zhang coloring->fset = PETSC_FALSE; 54fcd7ac73SHong Zhang } 55fcd7ac73SHong Zhang 568bc97078SHong Zhang /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */ 57f6af9589SHong Zhang if (coloring->htype[0] == 'w') { 58*9e917edbSHong Zhang /* vscale = dx is a constant scalar */ 59f6af9589SHong Zhang ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr); 6074d3cef9SHong Zhang dx = 1.0/((1.0 + unorm)*epsilon); 6170e7395fSHong Zhang } else { 6274d3cef9SHong Zhang ierr = VecGetLocalSize(x1,&nxloc);CHKERRQ(ierr); 63f6af9589SHong Zhang ierr = VecGetArray(x1,&xx);CHKERRQ(ierr); 648bc97078SHong Zhang ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr); 6574d3cef9SHong Zhang for (col=0; col<nxloc; col++) { 66fcd7ac73SHong Zhang dx = xx[col]; 6774d3cef9SHong Zhang if (PetscAbsScalar(dx) < umin) { 6874d3cef9SHong Zhang if (PetscRealPart(dx) >= 0.0) dx = umin; 6974d3cef9SHong Zhang else if (PetscRealPart(dx) < 0.0 ) dx = -umin; 7074d3cef9SHong Zhang } 71fcd7ac73SHong Zhang dx *= epsilon; 7274d3cef9SHong Zhang vscale_array[col] = 1.0/dx; 73f6af9589SHong Zhang } 7474d3cef9SHong Zhang ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 758bc97078SHong Zhang ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr); 7670e7395fSHong Zhang } 77*9e917edbSHong Zhang if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] != 'w') { 788bc97078SHong Zhang ierr = VecGhostUpdateBegin(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 798bc97078SHong Zhang ierr = VecGhostUpdateEnd(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 80fcd7ac73SHong Zhang } 81fcd7ac73SHong Zhang 828bc97078SHong Zhang /* (3) Loop over each color */ 838bc97078SHong Zhang if (!coloring->w3) { 848bc97078SHong Zhang ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 858bc97078SHong Zhang ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr); 868bc97078SHong Zhang } 878bc97078SHong Zhang w3 = coloring->w3; 88fcd7ac73SHong Zhang 898bc97078SHong Zhang ierr = VecGetOwnershipRange(x1,&cstart,&cend);CHKERRQ(ierr); /* used by ghosted vscale */ 90*9e917edbSHong Zhang if (vscale) { 918bc97078SHong Zhang ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr); 92*9e917edbSHong Zhang } 938bc97078SHong Zhang nz = 0; 948bc97078SHong Zhang for (k=0; k<ncolors; k++) { 958bc97078SHong Zhang coloring->currentcolor = k; 96*9e917edbSHong Zhang 97fcd7ac73SHong Zhang /* 988bc97078SHong Zhang (3-1) Loop over each column associated with color 99f6af9589SHong Zhang adding the perturbation to the vector w3 = x1 + dx. 100fcd7ac73SHong Zhang */ 101f6af9589SHong Zhang ierr = VecCopy(x1,w3);CHKERRQ(ierr); 102f6af9589SHong Zhang ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr); 1038bc97078SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array - cstart; /* shift pointer so global index can be used */ 1048bc97078SHong Zhang for (l=0; l<ncolumns[k]; l++) { 105f6af9589SHong Zhang col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 106*9e917edbSHong Zhang if (coloring->htype[0] == 'w') { 107*9e917edbSHong Zhang w3_array[col] += 1.0/dx; 108*9e917edbSHong Zhang continue; 109*9e917edbSHong Zhang } 11070e7395fSHong Zhang /* convert global col to local colb */ 1118bc97078SHong Zhang if (col >= cstart && col < cend) { 1128bc97078SHong Zhang colb = col - cstart; 11370e7395fSHong Zhang } else { 11470e7395fSHong Zhang #if defined(PETSC_USE_CTABLE) 11570e7395fSHong Zhang ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr); 11670e7395fSHong Zhang colb--; 11770e7395fSHong Zhang #else 11870e7395fSHong Zhang colb = aij->colmap[col] - 1; /* local column index */ 11970e7395fSHong Zhang #endif 1208bc97078SHong Zhang colb += cend - cstart; 12170e7395fSHong Zhang } 12274d3cef9SHong Zhang w3_array[col] += 1.0/vscale_array[colb]; 123fcd7ac73SHong Zhang } 1248bc97078SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array + cstart; 125fcd7ac73SHong Zhang ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 126*9e917edbSHong Zhang 127fcd7ac73SHong Zhang /* 1288bc97078SHong Zhang (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 129fcd7ac73SHong Zhang w2 = F(x1 + dx) - F(x1) 130fcd7ac73SHong Zhang */ 131fcd7ac73SHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 132fcd7ac73SHong Zhang ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 133fcd7ac73SHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 134fcd7ac73SHong Zhang ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); 135*9e917edbSHong Zhang 1368bc97078SHong Zhang /* 1378bc97078SHong Zhang (3-3) Loop over rows of vector, putting results into Jacobian matrix 1388bc97078SHong Zhang */ 139fcd7ac73SHong Zhang ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 1408bc97078SHong Zhang for (l=0; l<nrows[k]; l++) { 141573f477fSHong Zhang row = Jentry[nz].row; /* local row index */ 142*9e917edbSHong Zhang if (coloring->htype[0] == 'w') { 143*9e917edbSHong Zhang *(Jentry[nz].valaddr) = y[row]*dx; 144*9e917edbSHong Zhang } else { 145573f477fSHong Zhang *(Jentry[nz].valaddr) = y[row]*vscale_array[Jentry[nz].col]; 146*9e917edbSHong Zhang } 147573f477fSHong Zhang nz++; 148fcd7ac73SHong Zhang } 149fcd7ac73SHong Zhang ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 1508bc97078SHong Zhang } 151f5aae955SHong Zhang ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 152f5aae955SHong Zhang ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 153*9e917edbSHong Zhang if (vscale) { 1548bc97078SHong Zhang ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr); 155*9e917edbSHong Zhang } 156fcd7ac73SHong Zhang 157fcd7ac73SHong Zhang coloring->currentcolor = -1; 158fcd7ac73SHong Zhang PetscFunctionReturn(0); 159fcd7ac73SHong Zhang } 160fcd7ac73SHong Zhang 161fcd7ac73SHong Zhang #undef __FUNCT__ 162fcd7ac73SHong Zhang #define __FUNCT__ "MatFDColoringCreate_MPIAIJ_new" 163fcd7ac73SHong Zhang PetscErrorCode MatFDColoringCreate_MPIAIJ_new(Mat mat,ISColoring iscoloring,MatFDColoring c) 164fcd7ac73SHong Zhang { 165fcd7ac73SHong Zhang Mat_MPIAIJ *aij=(Mat_MPIAIJ*)mat->data; 166fcd7ac73SHong Zhang PetscErrorCode ierr; 167fcd7ac73SHong Zhang PetscMPIInt size,*ncolsonproc,*disp,nn; 168fcd7ac73SHong Zhang PetscInt i,n,nrows,j,k,m,ncols,col; 16972c15787SHong Zhang const PetscInt *is,*A_ci,*A_cj,*B_ci,*B_cj,*row = NULL,*ltog=NULL; 170fcd7ac73SHong Zhang PetscInt nis=iscoloring->n,nctot,*cols; 171d3825b63SHong Zhang PetscInt *rowhit,cstart,cend,colb; 172fcd7ac73SHong Zhang IS *isa; 173fcd7ac73SHong Zhang ISLocalToGlobalMapping map=mat->cmap->mapping; 174fcd7ac73SHong Zhang PetscInt ctype=c->ctype; 175fcd7ac73SHong Zhang Mat A=aij->A,B=aij->B; 1768bc97078SHong Zhang Mat_SeqAIJ *spA=(Mat_SeqAIJ*)A->data,*spB=(Mat_SeqAIJ*)B->data; 1778bc97078SHong Zhang PetscScalar *A_val=spA->a,*B_val=spB->a; 178fcd7ac73SHong Zhang PetscInt spidx; 1791b97d346SHong Zhang PetscInt *spidxA,*spidxB,nz; 180d3825b63SHong Zhang PetscScalar **valaddrhit; 181573f477fSHong Zhang MatEntry *Jentry; 182fcd7ac73SHong Zhang 183fcd7ac73SHong Zhang PetscFunctionBegin; 18472c15787SHong Zhang if (ctype == IS_COLORING_GHOSTED) { 18572c15787SHong Zhang if(!map) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_INCOMP,"When using ghosted differencing matrix must have local to global mapping provided with MatSetLocalToGlobalMapping"); 18672c15787SHong Zhang ierr = ISLocalToGlobalMappingGetIndices(map,<og);CHKERRQ(ierr); 18772c15787SHong Zhang } 188fcd7ac73SHong Zhang 189d3825b63SHong Zhang m = mat->rmap->n; 190fcd7ac73SHong Zhang cstart = mat->cmap->rstart; 191fcd7ac73SHong Zhang cend = mat->cmap->rend; 192fcd7ac73SHong Zhang c->M = mat->rmap->N; /* set the global rows and columns and local rows */ 193fcd7ac73SHong Zhang c->N = mat->cmap->N; 194d3825b63SHong Zhang c->m = m; 195fcd7ac73SHong Zhang c->rstart = mat->rmap->rstart; 196fcd7ac73SHong Zhang 197fcd7ac73SHong Zhang c->ncolors = nis; 198fcd7ac73SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 199fcd7ac73SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr); 200fcd7ac73SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 2018bc97078SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr); 202fcd7ac73SHong Zhang 2038bc97078SHong Zhang nz = spA->nz + spB->nz; /* total nonzero entries of mat */ 204573f477fSHong Zhang ierr = PetscMalloc(nz*sizeof(MatEntry),&Jentry);CHKERRQ(ierr); 2058bc97078SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr); 2068bc97078SHong Zhang c->matentry = Jentry; 207a774a6f1SHong Zhang 208d3825b63SHong Zhang /* Allow access to data structures of local part of matrix 209d3825b63SHong Zhang - creates aij->colmap which maps global column number to local number in part B */ 210fcd7ac73SHong Zhang if (!aij->colmap) { 211fcd7ac73SHong Zhang ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr); 212fcd7ac73SHong Zhang } 213fcd7ac73SHong Zhang ierr = MatGetColumnIJ_SeqAIJ_Color(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); 214fcd7ac73SHong Zhang ierr = MatGetColumnIJ_SeqAIJ_Color(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); 215fcd7ac73SHong Zhang 216d3825b63SHong Zhang ierr = PetscMalloc2(m+1,PetscInt,&rowhit,m+1,PetscScalar*,&valaddrhit);CHKERRQ(ierr); 217fcd7ac73SHong Zhang nz = 0; 21872c15787SHong Zhang ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 219d3825b63SHong Zhang for (i=0; i<nis; i++) { /* for each local color */ 220fcd7ac73SHong Zhang ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 221fcd7ac73SHong Zhang ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 222fcd7ac73SHong Zhang 223fcd7ac73SHong Zhang c->ncolumns[i] = n; /* local number of columns of this color on this process */ 224fcd7ac73SHong Zhang if (n) { 225fcd7ac73SHong Zhang ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr); 226fcd7ac73SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr); 227fcd7ac73SHong Zhang ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr); 228fcd7ac73SHong Zhang } else { 229fcd7ac73SHong Zhang c->columns[i] = 0; 230fcd7ac73SHong Zhang } 231fcd7ac73SHong Zhang 232fcd7ac73SHong Zhang if (ctype == IS_COLORING_GLOBAL) { 233d3825b63SHong Zhang /* Determine nctot, the total (parallel) number of columns of this color */ 234fcd7ac73SHong Zhang ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr); 235fcd7ac73SHong Zhang ierr = PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);CHKERRQ(ierr); 236fcd7ac73SHong Zhang 237d3825b63SHong Zhang /* ncolsonproc[j]: local ncolumns on proc[j] of this color */ 238fcd7ac73SHong Zhang ierr = PetscMPIIntCast(n,&nn);CHKERRQ(ierr); 239fcd7ac73SHong Zhang ierr = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 240fcd7ac73SHong Zhang nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j]; 241fcd7ac73SHong Zhang if (!nctot) { 242fcd7ac73SHong Zhang ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr); 243fcd7ac73SHong Zhang } 244fcd7ac73SHong Zhang 245fcd7ac73SHong Zhang disp[0] = 0; 246fcd7ac73SHong Zhang for (j=1; j<size; j++) { 247fcd7ac73SHong Zhang disp[j] = disp[j-1] + ncolsonproc[j-1]; 248fcd7ac73SHong Zhang } 249fcd7ac73SHong Zhang 250d3825b63SHong Zhang /* Get cols, the complete list of columns for this color on each process */ 251fcd7ac73SHong Zhang ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 252fcd7ac73SHong Zhang ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 253fcd7ac73SHong Zhang ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr); 254fcd7ac73SHong Zhang } else if (ctype == IS_COLORING_GHOSTED) { 255fcd7ac73SHong Zhang /* Determine local number of columns of this color on this process, including ghost points */ 256fcd7ac73SHong Zhang nctot = n; 257fcd7ac73SHong Zhang ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 258fcd7ac73SHong Zhang ierr = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr); 259fcd7ac73SHong Zhang } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type"); 260fcd7ac73SHong Zhang 2611b97d346SHong Zhang /* Mark all rows affect by these columns */ 262d3825b63SHong Zhang ierr = PetscMemzero(rowhit,m*sizeof(PetscInt));CHKERRQ(ierr); 2631b97d346SHong Zhang for (j=0; j<nctot; j++) { /* loop over columns*/ 264fcd7ac73SHong Zhang if (ctype == IS_COLORING_GHOSTED) { 265fcd7ac73SHong Zhang col = ltog[cols[j]]; 266fcd7ac73SHong Zhang } else { 267fcd7ac73SHong Zhang col = cols[j]; 268fcd7ac73SHong Zhang } 269fcd7ac73SHong Zhang if (col >= cstart && col < cend) { /* column is in diagonal block of matrix A */ 270d3825b63SHong Zhang row = A_cj + A_ci[col-cstart]; 271d3825b63SHong Zhang nrows = A_ci[col-cstart+1] - A_ci[col-cstart]; 272fcd7ac73SHong Zhang /* loop over columns of A marking them in rowhit */ 273d3825b63SHong Zhang for (k=0; k<nrows; k++) { 274d3825b63SHong Zhang /* set valaddrhit for part A */ 275fcd7ac73SHong Zhang spidx = spidxA[A_ci[col-cstart] + k]; 276d3825b63SHong Zhang valaddrhit[*row] = &A_val[spidx]; 277a774a6f1SHong Zhang rowhit[*row++] = col - cstart + 1; /* local column index */ 278fcd7ac73SHong Zhang } 279fcd7ac73SHong Zhang } else { /* column is in off-diagonal block of matrix B */ 280fcd7ac73SHong Zhang #if defined(PETSC_USE_CTABLE) 281fcd7ac73SHong Zhang ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr); 282fcd7ac73SHong Zhang colb--; 283fcd7ac73SHong Zhang #else 284fcd7ac73SHong Zhang colb = aij->colmap[col] - 1; /* local column index */ 285fcd7ac73SHong Zhang #endif 286fcd7ac73SHong Zhang if (colb == -1) { 287d3825b63SHong Zhang nrows = 0; 288fcd7ac73SHong Zhang } else { 289d3825b63SHong Zhang row = B_cj + B_ci[colb]; 290d3825b63SHong Zhang nrows = B_ci[colb+1] - B_ci[colb]; 291fcd7ac73SHong Zhang } 292fcd7ac73SHong Zhang /* loop over columns of B marking them in rowhit */ 293d3825b63SHong Zhang for (k=0; k<nrows; k++) { 294d3825b63SHong Zhang /* set valaddrhit for part B */ 295fcd7ac73SHong Zhang spidx = spidxB[B_ci[colb] + k]; 296d3825b63SHong Zhang valaddrhit[*row] = &B_val[spidx]; 29770e7395fSHong Zhang rowhit[*row++] = colb + 1 + cend - cstart; /* local column index */ 298fcd7ac73SHong Zhang } 299fcd7ac73SHong Zhang } 300fcd7ac73SHong Zhang } 301fcd7ac73SHong Zhang 302fcd7ac73SHong Zhang /* count the number of hits */ 303fcd7ac73SHong Zhang nrows = 0; 304d3825b63SHong Zhang for (j=0; j<m; j++) { 305fcd7ac73SHong Zhang if (rowhit[j]) nrows++; 306fcd7ac73SHong Zhang } 307fcd7ac73SHong Zhang c->nrows[i] = nrows; 3088bc97078SHong Zhang 309d3825b63SHong Zhang for (j=0; j<m; j++) { 310fcd7ac73SHong Zhang if (rowhit[j]) { 311573f477fSHong Zhang Jentry[nz].row = j; /* local row index */ 312573f477fSHong Zhang Jentry[nz].col = rowhit[j] - 1; /* local column index */ 313573f477fSHong Zhang Jentry[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */ 314573f477fSHong Zhang nz++; 315fcd7ac73SHong Zhang } 316fcd7ac73SHong Zhang } 317fcd7ac73SHong Zhang ierr = PetscFree(cols);CHKERRQ(ierr); 318fcd7ac73SHong Zhang } 3198bc97078SHong Zhang ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 3208bc97078SHong Zhang if (nz != spA->nz + spB->nz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"nz %d != mat->nz %d\n",nz,spA->nz+spB->nz); 321fcd7ac73SHong Zhang 322d3825b63SHong Zhang ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr); 323fcd7ac73SHong Zhang ierr = MatRestoreColumnIJ_SeqAIJ_Color(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); 324fcd7ac73SHong Zhang ierr = MatRestoreColumnIJ_SeqAIJ_Color(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); 32572c15787SHong Zhang if (ctype == IS_COLORING_GHOSTED) { 32672c15787SHong Zhang ierr = ISLocalToGlobalMappingRestoreIndices(map,<og);CHKERRQ(ierr); 32772c15787SHong Zhang } 328fcd7ac73SHong Zhang 329fcd7ac73SHong Zhang mat->ops->fdcoloringapply = MatFDColoringApply_MPIAIJ; 330fcd7ac73SHong Zhang PetscFunctionReturn(0); 331fcd7ac73SHong Zhang } 332fcd7ac73SHong Zhang 333fcd7ac73SHong Zhang /*------------------------------------------------------*/ 334fcd7ac73SHong Zhang #undef __FUNCT__ 3354a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate_MPIAIJ" 336dfbe8321SBarry Smith PetscErrorCode MatFDColoringCreate_MPIAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c) 337a64fbb32SBarry Smith { 3386eaac0f3SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)mat->data; 3396849ba73SBarry Smith PetscErrorCode ierr; 340b1d57f15SBarry Smith PetscMPIInt size,*ncolsonproc,*disp,nn; 3411a83f524SJed Brown PetscInt i,n,nrows,j,k,m,ncols,col; 342afcb2eb5SJed Brown const PetscInt *is,*A_ci,*A_cj,*B_ci,*B_cj,*rows = 0,*ltog; 3431a83f524SJed Brown PetscInt nis = iscoloring->n,nctot,*cols; 344f6d58c54SBarry Smith PetscInt *rowhit,M,cstart,cend,colb; 345b1d57f15SBarry Smith PetscInt *columnsforrow,l; 346b9617806SBarry Smith IS *isa; 347ace3abfcSBarry Smith PetscBool done,flg; 348992144d0SBarry Smith ISLocalToGlobalMapping map = mat->cmap->mapping; 349afcb2eb5SJed Brown PetscInt ctype=c->ctype; 350fcd7ac73SHong Zhang PetscBool new_impl=PETSC_FALSE; 351a64fbb32SBarry Smith 3523a40ed3dSBarry Smith PetscFunctionBegin; 353fcd7ac73SHong Zhang ierr = PetscOptionsName("-new","using new impls","",&new_impl);CHKERRQ(ierr); 354fcd7ac73SHong Zhang if (new_impl){ 355fcd7ac73SHong Zhang ierr = MatFDColoringCreate_MPIAIJ_new(mat,iscoloring,c);CHKERRQ(ierr); 356fcd7ac73SHong Zhang PetscFunctionReturn(0); 357fcd7ac73SHong Zhang } 358ce94432eSBarry Smith 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"); 359522c5e43SBarry Smith 360afcb2eb5SJed Brown if (map) {ierr = ISLocalToGlobalMappingGetIndices(map,<og);CHKERRQ(ierr);} 361afcb2eb5SJed Brown else ltog = NULL; 362b9617806SBarry Smith ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 3633acb8795SBarry Smith 364f6d58c54SBarry Smith M = mat->rmap->n; 365f6d58c54SBarry Smith cstart = mat->cmap->rstart; 366f6d58c54SBarry Smith cend = mat->cmap->rend; 367f6d58c54SBarry Smith c->M = mat->rmap->N; /* set the global rows and columns and local rows */ 368f6d58c54SBarry Smith c->N = mat->cmap->N; 369f6d58c54SBarry Smith c->m = mat->rmap->n; 370f6d58c54SBarry Smith c->rstart = mat->rmap->rstart; 371005c665bSBarry Smith 372a64fbb32SBarry Smith c->ncolors = nis; 373b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 374b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr); 375b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 376b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->rows);CHKERRQ(ierr); 377b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columnsforrow);CHKERRQ(ierr); 3783bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)c,5*nis*sizeof(PetscInt));CHKERRQ(ierr); 3796eaac0f3SBarry Smith 3806eaac0f3SBarry Smith /* Allow access to data structures of local part of matrix */ 3816eaac0f3SBarry Smith if (!aij->colmap) { 382ab9863d7SBarry Smith ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr); 3836eaac0f3SBarry Smith } 3843acb8795SBarry Smith ierr = MatGetColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr); 3853acb8795SBarry Smith ierr = MatGetColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr); 3866eaac0f3SBarry Smith 387b1d57f15SBarry Smith ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr); 388b1d57f15SBarry Smith ierr = PetscMalloc((M+1)*sizeof(PetscInt),&columnsforrow);CHKERRQ(ierr); 3896eaac0f3SBarry Smith 390a64fbb32SBarry Smith for (i=0; i<nis; i++) { 391b9b97703SBarry Smith ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 392a64fbb32SBarry Smith ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 3932205254eSKarl Rupp 394fcd7ac73SHong Zhang c->ncolumns[i] = n; /* local number of columns of this color on this process */ 395a64fbb32SBarry Smith if (n) { 396b1d57f15SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr); 3973bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr); 398b1d57f15SBarry Smith ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr); 399a64fbb32SBarry Smith } else { 400a64fbb32SBarry Smith c->columns[i] = 0; 401a64fbb32SBarry Smith } 402a64fbb32SBarry Smith 4038ee2e534SBarry Smith if (ctype == IS_COLORING_GLOBAL) { 4046eaac0f3SBarry Smith /* Determine the total (parallel) number of columns of this color */ 405ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr); 406687f1162SBarry Smith ierr = PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);CHKERRQ(ierr); 407b8f8c88eSHong Zhang 4084dc2109aSBarry Smith ierr = PetscMPIIntCast(n,&nn);CHKERRQ(ierr); 409ce94432eSBarry Smith ierr = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 4102205254eSKarl Rupp nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j]; 4113a7fca6bSBarry Smith if (!nctot) { 412ae15b995SBarry Smith ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr); 4133a7fca6bSBarry Smith } 4146eaac0f3SBarry Smith 4156eaac0f3SBarry Smith disp[0] = 0; 4166eaac0f3SBarry Smith for (j=1; j<size; j++) { 4176eaac0f3SBarry Smith disp[j] = disp[j-1] + ncolsonproc[j-1]; 4186eaac0f3SBarry Smith } 4196eaac0f3SBarry Smith 4206eaac0f3SBarry Smith /* Get complete list of columns for color on each processor */ 421b1d57f15SBarry Smith ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 422ce94432eSBarry Smith ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 4231d79065fSBarry Smith ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr); 424b8f8c88eSHong Zhang } else if (ctype == IS_COLORING_GHOSTED) { 425b8f8c88eSHong Zhang /* Determine local number of columns of this color on this process, including ghost points */ 426b8f8c88eSHong Zhang nctot = n; 427b8f8c88eSHong Zhang ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 428b8f8c88eSHong Zhang ierr = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr); 429f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type"); 4306eaac0f3SBarry Smith 4316eaac0f3SBarry Smith /* 4326eaac0f3SBarry Smith Mark all rows affect by these columns 4336eaac0f3SBarry Smith */ 434b8f8c88eSHong Zhang /* Temporary option to allow for debugging/testing */ 43590d69ab7SBarry Smith flg = PETSC_FALSE; 4360298fd71SBarry Smith ierr = PetscOptionsGetBool(NULL,"-matfdcoloring_slow",&flg,NULL);CHKERRQ(ierr); 437f158e583SBarry Smith if (!flg) { /*-----------------------------------------------------------------------------*/ 438f158e583SBarry Smith /* crude, fast version */ 439b1d57f15SBarry Smith ierr = PetscMemzero(rowhit,M*sizeof(PetscInt));CHKERRQ(ierr); 440a64fbb32SBarry Smith /* loop over columns*/ 4416eaac0f3SBarry Smith for (j=0; j<nctot; j++) { 442b8f8c88eSHong Zhang if (ctype == IS_COLORING_GHOSTED) { 443b8f8c88eSHong Zhang col = ltog[cols[j]]; 444b8f8c88eSHong Zhang } else { 4456eaac0f3SBarry Smith col = cols[j]; 446b8f8c88eSHong Zhang } 4476eaac0f3SBarry Smith if (col >= cstart && col < cend) { 4486eaac0f3SBarry Smith /* column is in diagonal block of matrix */ 4496eaac0f3SBarry Smith rows = A_cj + A_ci[col-cstart]; 4506eaac0f3SBarry Smith m = A_ci[col-cstart+1] - A_ci[col-cstart]; 4516eaac0f3SBarry Smith } else { 452aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 453cb9801acSJed Brown ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr); 454fa46199cSSatish Balay colb--; 455b3d2dc96SSatish Balay #else 4566eaac0f3SBarry Smith colb = aij->colmap[col] - 1; 457b3d2dc96SSatish Balay #endif 4586eaac0f3SBarry Smith if (colb == -1) { 4596eaac0f3SBarry Smith m = 0; 4606eaac0f3SBarry Smith } else { 4616eaac0f3SBarry Smith rows = B_cj + B_ci[colb]; 4626eaac0f3SBarry Smith m = B_ci[colb+1] - B_ci[colb]; 4636eaac0f3SBarry Smith } 4646eaac0f3SBarry Smith } 465a64fbb32SBarry Smith /* loop over columns marking them in rowhit */ 466a64fbb32SBarry Smith for (k=0; k<m; k++) { 467a64fbb32SBarry Smith rowhit[*rows++] = col + 1; 468a64fbb32SBarry Smith } 469a64fbb32SBarry Smith } 4706eaac0f3SBarry Smith 471a64fbb32SBarry Smith /* count the number of hits */ 472a64fbb32SBarry Smith nrows = 0; 4736eaac0f3SBarry Smith for (j=0; j<M; j++) { 474a64fbb32SBarry Smith if (rowhit[j]) nrows++; 475a64fbb32SBarry Smith } 476a64fbb32SBarry Smith c->nrows[i] = nrows; 477b1d57f15SBarry Smith ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr); 478b1d57f15SBarry Smith ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr); 4793bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)c,2*(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr); 480a64fbb32SBarry Smith nrows = 0; 4816eaac0f3SBarry Smith for (j=0; j<M; j++) { 482a64fbb32SBarry Smith if (rowhit[j]) { 483fcd7ac73SHong Zhang c->rows[i][nrows] = j; /* local row index */ 484fcd7ac73SHong Zhang c->columnsforrow[i][nrows] = rowhit[j] - 1; /* global column index */ 485a64fbb32SBarry Smith nrows++; 486a64fbb32SBarry Smith } 487a64fbb32SBarry Smith } 488a64fbb32SBarry Smith } else { /*-------------------------------------------------------------------------------*/ 489f158e583SBarry Smith /* slow version, using rowhit as a linked list */ 490b1d57f15SBarry Smith PetscInt currentcol,fm,mfm; 4916eaac0f3SBarry Smith rowhit[M] = M; 492a64fbb32SBarry Smith nrows = 0; 493a64fbb32SBarry Smith /* loop over columns*/ 4946eaac0f3SBarry Smith for (j=0; j<nctot; j++) { 495b8f8c88eSHong Zhang if (ctype == IS_COLORING_GHOSTED) { 496b8f8c88eSHong Zhang col = ltog[cols[j]]; 497b8f8c88eSHong Zhang } else { 4986eaac0f3SBarry Smith col = cols[j]; 499b8f8c88eSHong Zhang } 5006eaac0f3SBarry Smith if (col >= cstart && col < cend) { 5016eaac0f3SBarry Smith /* column is in diagonal block of matrix */ 5026eaac0f3SBarry Smith rows = A_cj + A_ci[col-cstart]; 5036eaac0f3SBarry Smith m = A_ci[col-cstart+1] - A_ci[col-cstart]; 5046eaac0f3SBarry Smith } else { 505aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 5060f5bd95cSBarry Smith ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr); 507fa46199cSSatish Balay colb--; 508b3d2dc96SSatish Balay #else 5096eaac0f3SBarry Smith colb = aij->colmap[col] - 1; 510b3d2dc96SSatish Balay #endif 5116eaac0f3SBarry Smith if (colb == -1) { 5126eaac0f3SBarry Smith m = 0; 5136eaac0f3SBarry Smith } else { 5146eaac0f3SBarry Smith rows = B_cj + B_ci[colb]; 5156eaac0f3SBarry Smith m = B_ci[colb+1] - B_ci[colb]; 5166eaac0f3SBarry Smith } 5176eaac0f3SBarry Smith } 518b8f8c88eSHong Zhang 519a64fbb32SBarry Smith /* loop over columns marking them in rowhit */ 5206eaac0f3SBarry Smith fm = M; /* fm points to first entry in linked list */ 521a64fbb32SBarry Smith for (k=0; k<m; k++) { 522a64fbb32SBarry Smith currentcol = *rows++; 523a64fbb32SBarry Smith /* is it already in the list? */ 524a64fbb32SBarry Smith do { 525a64fbb32SBarry Smith mfm = fm; 526a64fbb32SBarry Smith fm = rowhit[fm]; 527a64fbb32SBarry Smith } while (fm < currentcol); 528a64fbb32SBarry Smith /* not in list so add it */ 529a64fbb32SBarry Smith if (fm != currentcol) { 530a64fbb32SBarry Smith nrows++; 531a64fbb32SBarry Smith columnsforrow[currentcol] = col; 532a64fbb32SBarry Smith /* next three lines insert new entry into linked list */ 533a64fbb32SBarry Smith rowhit[mfm] = currentcol; 534a64fbb32SBarry Smith rowhit[currentcol] = fm; 535a64fbb32SBarry Smith fm = currentcol; 536a64fbb32SBarry Smith /* fm points to present position in list since we know the columns are sorted */ 537f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid coloring of matrix detected"); 538a64fbb32SBarry Smith } 539a64fbb32SBarry Smith } 540a64fbb32SBarry Smith c->nrows[i] = nrows; 5412205254eSKarl Rupp 542b1d57f15SBarry Smith ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr); 543b1d57f15SBarry Smith ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr); 5443bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)c,(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr); 545a64fbb32SBarry Smith /* now store the linked list of rows into c->rows[i] */ 546a64fbb32SBarry Smith nrows = 0; 5476eaac0f3SBarry Smith fm = rowhit[M]; 548a64fbb32SBarry Smith do { 549a64fbb32SBarry Smith c->rows[i][nrows] = fm; 550a64fbb32SBarry Smith c->columnsforrow[i][nrows++] = columnsforrow[fm]; 551a64fbb32SBarry Smith fm = rowhit[fm]; 5526eaac0f3SBarry Smith } while (fm < M); 5536eaac0f3SBarry Smith } /* ---------------------------------------------------------------------------------------*/ 554606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 5556eaac0f3SBarry Smith } 55630b34957SBarry Smith 55730b34957SBarry Smith /* Optimize by adding the vscale, and scaleforrow[][] fields */ 55830b34957SBarry Smith /* 55930b34957SBarry Smith vscale will contain the "diagonal" on processor scalings followed by the off processor 56030b34957SBarry Smith */ 5618ee2e534SBarry Smith if (ctype == IS_COLORING_GLOBAL) { 562ce94432eSBarry Smith ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),aij->A->rmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&c->vscale);CHKERRQ(ierr); 563b1d57f15SBarry Smith ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr); 56430b34957SBarry Smith for (k=0; k<c->ncolors; k++) { 565b1d57f15SBarry Smith ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr); 56630b34957SBarry Smith for (l=0; l<c->nrows[k]; l++) { 56730b34957SBarry Smith col = c->columnsforrow[k][l]; 56830b34957SBarry Smith if (col >= cstart && col < cend) { 56930b34957SBarry Smith /* column is in diagonal block of matrix */ 57030b34957SBarry Smith colb = col - cstart; 57130b34957SBarry Smith } else { 57230b34957SBarry Smith /* column is in "off-processor" part */ 57330b34957SBarry Smith #if defined(PETSC_USE_CTABLE) 57430b34957SBarry Smith ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr); 57530b34957SBarry Smith colb--; 57630b34957SBarry Smith #else 57730b34957SBarry Smith colb = aij->colmap[col] - 1; 57830b34957SBarry Smith #endif 57930b34957SBarry Smith colb += cend - cstart; 58030b34957SBarry Smith } 58130b34957SBarry Smith c->vscaleforrow[k][l] = colb; 58230b34957SBarry Smith } 58330b34957SBarry Smith } 584b8f8c88eSHong Zhang } else if (ctype == IS_COLORING_GHOSTED) { 585b8f8c88eSHong Zhang /* Get gtol mapping */ 586afcb2eb5SJed Brown PetscInt N = mat->cmap->N,nlocal,*gtol; 587b8f8c88eSHong Zhang ierr = PetscMalloc((N+1)*sizeof(PetscInt),>ol);CHKERRQ(ierr); 588b8f8c88eSHong Zhang for (i=0; i<N; i++) gtol[i] = -1; 589afcb2eb5SJed Brown ierr = ISLocalToGlobalMappingGetSize(map,&nlocal);CHKERRQ(ierr); 590afcb2eb5SJed Brown for (i=0; i<nlocal; i++) gtol[ltog[i]] = i; 591b8f8c88eSHong Zhang 592b8f8c88eSHong Zhang c->vscale = 0; /* will be created in MatFDColoringApply() */ 593b8f8c88eSHong Zhang ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr); 594b8f8c88eSHong Zhang for (k=0; k<c->ncolors; k++) { 595b8f8c88eSHong Zhang ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr); 596b8f8c88eSHong Zhang for (l=0; l<c->nrows[k]; l++) { 597b8f8c88eSHong Zhang col = c->columnsforrow[k][l]; /* global column index */ 598b8f8c88eSHong Zhang c->vscaleforrow[k][l] = gtol[col]; /* local column index */ 599b8f8c88eSHong Zhang } 600b8f8c88eSHong Zhang } 601b8f8c88eSHong Zhang ierr = PetscFree(gtol);CHKERRQ(ierr); 602b8f8c88eSHong Zhang } 603b9617806SBarry Smith ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 60430b34957SBarry Smith 605606d414cSSatish Balay ierr = PetscFree(rowhit);CHKERRQ(ierr); 606606d414cSSatish Balay ierr = PetscFree(columnsforrow);CHKERRQ(ierr); 6073acb8795SBarry Smith ierr = MatRestoreColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr); 6083acb8795SBarry Smith ierr = MatRestoreColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr); 609afcb2eb5SJed Brown if (map) {ierr = ISLocalToGlobalMappingRestoreIndices(map,<og);CHKERRQ(ierr);} 6103a40ed3dSBarry Smith PetscFunctionReturn(0); 611a64fbb32SBarry Smith } 612a64fbb32SBarry Smith 613b9617806SBarry Smith 614b9617806SBarry Smith 615b9617806SBarry Smith 616b9617806SBarry Smith 617b9617806SBarry Smith 618