xref: /petsc/src/mat/impls/aij/mpi/fdmpiaij.c (revision 4b2e90ca22e5d236f72fef601ed2fa3dbdcf4eae)
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;
12a2f2d239SHong Zhang   PetscInt       k,cstart,cend,l,row,col,nz;
139e917edbSHong 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 
379e917edbSHong Zhang   /* create vscale for storing dx */
389e917edbSHong Zhang   if (!vscale) {
399e917edbSHong Zhang     if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') {
409e917edbSHong Zhang       ierr = VecCreateGhost(PetscObjectComm((PetscObject)J),J->cmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&vscale);CHKERRQ(ierr);
419e917edbSHong Zhang     } else if (ctype == IS_COLORING_GHOSTED) {
428bc97078SHong Zhang       ierr = VecDuplicate(x1,&vscale);CHKERRQ(ierr);
439e917edbSHong 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') {
589e917edbSHong 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   }
779e917edbSHong 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 */
909e917edbSHong Zhang   if (vscale) {
918bc97078SHong Zhang     ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr);
929e917edbSHong Zhang   }
938bc97078SHong Zhang   nz   = 0;
948bc97078SHong Zhang   for (k=0; k<ncolors; k++) {
958bc97078SHong Zhang     coloring->currentcolor = k;
969e917edbSHong 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);
103a2f2d239SHong Zhang     if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */
104b039d6c4SHong Zhang     if (coloring->htype[0] == 'w') {
1058bc97078SHong Zhang       for (l=0; l<ncolumns[k]; l++) {
106f6af9589SHong Zhang         col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
1079e917edbSHong Zhang         w3_array[col] += 1.0/dx;
1089e917edbSHong Zhang       }
109b039d6c4SHong Zhang     } else { /* htype == 'ds' */
110a2f2d239SHong Zhang       vscale_array -= cstart; /* shift pointer so global index can be used */
111b039d6c4SHong Zhang       for (l=0; l<ncolumns[k]; l++) {
112b039d6c4SHong Zhang         col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
113a2f2d239SHong Zhang         w3_array[col] += 1.0/vscale_array[col];
11470e7395fSHong Zhang       }
115a2f2d239SHong Zhang       vscale_array += cstart;
116fcd7ac73SHong Zhang     }
117a2f2d239SHong Zhang     if (ctype == IS_COLORING_GLOBAL) w3_array += cstart;
118fcd7ac73SHong Zhang     ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
1199e917edbSHong Zhang 
120fcd7ac73SHong Zhang     /*
1218bc97078SHong Zhang       (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
122fcd7ac73SHong Zhang                            w2 = F(x1 + dx) - F(x1)
123fcd7ac73SHong Zhang     */
124fcd7ac73SHong Zhang     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
125fcd7ac73SHong Zhang     ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
126fcd7ac73SHong Zhang     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
127fcd7ac73SHong Zhang     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
1289e917edbSHong Zhang 
1298bc97078SHong Zhang     /*
1308bc97078SHong Zhang      (3-3) Loop over rows of vector, putting results into Jacobian matrix
1318bc97078SHong Zhang     */
132fcd7ac73SHong Zhang     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
133b039d6c4SHong Zhang     if (coloring->htype[0] == 'w') {
1348bc97078SHong Zhang       for (l=0; l<nrows[k]; l++) {
135573f477fSHong Zhang         row                     = Jentry[nz].row;   /* local row index */
136b039d6c4SHong Zhang         *(Jentry[nz++].valaddr) = y[row]*dx;
1379e917edbSHong Zhang       }
138b039d6c4SHong Zhang     } else { /* htype == 'ds' */
139b039d6c4SHong Zhang       for (l=0; l<nrows[k]; l++) {
140b039d6c4SHong Zhang         row                   = Jentry[nz].row;   /* local row index */
141b039d6c4SHong Zhang         *(Jentry[nz].valaddr) = y[row]*vscale_array[Jentry[nz].col];
142573f477fSHong Zhang         nz++;
143fcd7ac73SHong Zhang       }
144b039d6c4SHong Zhang     }
145fcd7ac73SHong Zhang     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
1468bc97078SHong Zhang   }
147f5aae955SHong Zhang   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
148f5aae955SHong Zhang   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1499e917edbSHong Zhang   if (vscale) {
1508bc97078SHong Zhang     ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr);
1519e917edbSHong Zhang   }
152fcd7ac73SHong Zhang 
153fcd7ac73SHong Zhang   coloring->currentcolor = -1;
154fcd7ac73SHong Zhang   PetscFunctionReturn(0);
155fcd7ac73SHong Zhang }
156fcd7ac73SHong Zhang 
157fcd7ac73SHong Zhang #undef __FUNCT__
158fcd7ac73SHong Zhang #define __FUNCT__ "MatFDColoringCreate_MPIAIJ_new"
159fcd7ac73SHong Zhang PetscErrorCode MatFDColoringCreate_MPIAIJ_new(Mat mat,ISColoring iscoloring,MatFDColoring c)
160fcd7ac73SHong Zhang {
161fcd7ac73SHong Zhang   Mat_MPIAIJ             *aij=(Mat_MPIAIJ*)mat->data;
162fcd7ac73SHong Zhang   PetscErrorCode         ierr;
163fcd7ac73SHong Zhang   PetscMPIInt            size,*ncolsonproc,*disp,nn;
164*4b2e90caSHong Zhang   PetscInt               i,n,nrows,nrows_i,j,k,m,ncols,col;
16572c15787SHong Zhang   const PetscInt         *is,*A_ci,*A_cj,*B_ci,*B_cj,*row = NULL,*ltog=NULL;
166fcd7ac73SHong Zhang   PetscInt               nis=iscoloring->n,nctot,*cols;
167d3825b63SHong Zhang   PetscInt               *rowhit,cstart,cend,colb;
168fcd7ac73SHong Zhang   IS                     *isa;
169fcd7ac73SHong Zhang   ISLocalToGlobalMapping map=mat->cmap->mapping;
170fcd7ac73SHong Zhang   PetscInt               ctype=c->ctype;
171fcd7ac73SHong Zhang   Mat                    A=aij->A,B=aij->B;
1728bc97078SHong Zhang   Mat_SeqAIJ             *spA=(Mat_SeqAIJ*)A->data,*spB=(Mat_SeqAIJ*)B->data;
1738bc97078SHong Zhang   PetscScalar            *A_val=spA->a,*B_val=spB->a;
174fcd7ac73SHong Zhang   PetscInt               spidx;
1751b97d346SHong Zhang   PetscInt               *spidxA,*spidxB,nz;
176d3825b63SHong Zhang   PetscScalar            **valaddrhit;
177573f477fSHong Zhang   MatEntry               *Jentry;
178fcd7ac73SHong Zhang 
179fcd7ac73SHong Zhang   PetscFunctionBegin;
18072c15787SHong Zhang   if (ctype == IS_COLORING_GHOSTED) {
18172c15787SHong 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");
18272c15787SHong Zhang     ierr = ISLocalToGlobalMappingGetIndices(map,&ltog);CHKERRQ(ierr);
18372c15787SHong Zhang   }
184fcd7ac73SHong Zhang 
185d3825b63SHong Zhang   m         = mat->rmap->n;
186fcd7ac73SHong Zhang   cstart    = mat->cmap->rstart;
187fcd7ac73SHong Zhang   cend      = mat->cmap->rend;
188fcd7ac73SHong Zhang   c->M      = mat->rmap->N;         /* set the global rows and columns and local rows */
189fcd7ac73SHong Zhang   c->N      = mat->cmap->N;
190d3825b63SHong Zhang   c->m      = m;
191fcd7ac73SHong Zhang   c->rstart = mat->rmap->rstart;
192fcd7ac73SHong Zhang 
193fcd7ac73SHong Zhang   c->ncolors = nis;
194fcd7ac73SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
195fcd7ac73SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr);
196fcd7ac73SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
1978bc97078SHong Zhang   ierr       = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr);
198fcd7ac73SHong Zhang 
1998bc97078SHong Zhang   nz         = spA->nz + spB->nz; /* total nonzero entries of mat */
200573f477fSHong Zhang   ierr       = PetscMalloc(nz*sizeof(MatEntry),&Jentry);CHKERRQ(ierr);
2018bc97078SHong Zhang   ierr       = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr);
2028bc97078SHong Zhang   c->matentry = Jentry;
203a774a6f1SHong Zhang 
204d3825b63SHong Zhang   /* Allow access to data structures of local part of matrix
205d3825b63SHong Zhang    - creates aij->colmap which maps global column number to local number in part B */
206fcd7ac73SHong Zhang   if (!aij->colmap) {
207fcd7ac73SHong Zhang     ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
208fcd7ac73SHong Zhang   }
209fcd7ac73SHong Zhang   ierr = MatGetColumnIJ_SeqAIJ_Color(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
210fcd7ac73SHong Zhang   ierr = MatGetColumnIJ_SeqAIJ_Color(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
211fcd7ac73SHong Zhang 
212d3825b63SHong Zhang   ierr = PetscMalloc2(m+1,PetscInt,&rowhit,m+1,PetscScalar*,&valaddrhit);CHKERRQ(ierr);
213fcd7ac73SHong Zhang   nz = 0;
21472c15787SHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
215d3825b63SHong Zhang   for (i=0; i<nis; i++) { /* for each local color */
216fcd7ac73SHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
217fcd7ac73SHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
218fcd7ac73SHong Zhang 
219fcd7ac73SHong Zhang     c->ncolumns[i] = n; /* local number of columns of this color on this process */
220fcd7ac73SHong Zhang     if (n) {
221fcd7ac73SHong Zhang       ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr);
222fcd7ac73SHong Zhang       ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr);
223fcd7ac73SHong Zhang       ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr);
224a2f2d239SHong Zhang       /* convert global column indices to local ones! */
225a2f2d239SHong Zhang 
226fcd7ac73SHong Zhang     } else {
227fcd7ac73SHong Zhang       c->columns[i] = 0;
228fcd7ac73SHong Zhang     }
229fcd7ac73SHong Zhang 
230fcd7ac73SHong Zhang     if (ctype == IS_COLORING_GLOBAL) {
231d3825b63SHong Zhang       /* Determine nctot, the total (parallel) number of columns of this color */
232fcd7ac73SHong Zhang       ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr);
233fcd7ac73SHong Zhang       ierr = PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);CHKERRQ(ierr);
234fcd7ac73SHong Zhang 
235d3825b63SHong Zhang       /* ncolsonproc[j]: local ncolumns on proc[j] of this color */
236fcd7ac73SHong Zhang       ierr  = PetscMPIIntCast(n,&nn);CHKERRQ(ierr);
237fcd7ac73SHong Zhang       ierr  = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
238fcd7ac73SHong Zhang       nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j];
239fcd7ac73SHong Zhang       if (!nctot) {
240fcd7ac73SHong Zhang         ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr);
241fcd7ac73SHong Zhang       }
242fcd7ac73SHong Zhang 
243fcd7ac73SHong Zhang       disp[0] = 0;
244fcd7ac73SHong Zhang       for (j=1; j<size; j++) {
245fcd7ac73SHong Zhang         disp[j] = disp[j-1] + ncolsonproc[j-1];
246fcd7ac73SHong Zhang       }
247fcd7ac73SHong Zhang 
248d3825b63SHong Zhang       /* Get cols, the complete list of columns for this color on each process */
249fcd7ac73SHong Zhang       ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
250fcd7ac73SHong Zhang       ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
251fcd7ac73SHong Zhang       ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr);
252fcd7ac73SHong Zhang     } else if (ctype == IS_COLORING_GHOSTED) {
253fcd7ac73SHong Zhang       /* Determine local number of columns of this color on this process, including ghost points */
254fcd7ac73SHong Zhang       nctot = n;
255fcd7ac73SHong Zhang       ierr  = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
256fcd7ac73SHong Zhang       ierr  = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr);
257fcd7ac73SHong Zhang     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type");
258fcd7ac73SHong Zhang 
2591b97d346SHong Zhang     /* Mark all rows affect by these columns */
260d3825b63SHong Zhang     ierr = PetscMemzero(rowhit,m*sizeof(PetscInt));CHKERRQ(ierr);
261*4b2e90caSHong Zhang     nrows_i = 0;
2621b97d346SHong Zhang     for (j=0; j<nctot; j++) { /* loop over columns*/
263fcd7ac73SHong Zhang       if (ctype == IS_COLORING_GHOSTED) {
264fcd7ac73SHong Zhang         col = ltog[cols[j]];
265fcd7ac73SHong Zhang       } else {
266fcd7ac73SHong Zhang         col = cols[j];
267fcd7ac73SHong Zhang       }
268fcd7ac73SHong Zhang       if (col >= cstart && col < cend) { /* column is in diagonal block of matrix A */
269d3825b63SHong Zhang         row      = A_cj + A_ci[col-cstart];
270d3825b63SHong Zhang         nrows    = A_ci[col-cstart+1] - A_ci[col-cstart];
271*4b2e90caSHong Zhang         nrows_i += nrows;
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         }
292*4b2e90caSHong Zhang         nrows_i += nrows;
293fcd7ac73SHong Zhang         /* loop over columns of B marking them in rowhit */
294d3825b63SHong Zhang         for (k=0; k<nrows; k++) {
295d3825b63SHong Zhang           /* set valaddrhit for part B */
296fcd7ac73SHong Zhang           spidx            = spidxB[B_ci[colb] + k];
297d3825b63SHong Zhang           valaddrhit[*row] = &B_val[spidx];
29870e7395fSHong Zhang           rowhit[*row++]   = colb + 1 + cend - cstart; /* local column index */
299fcd7ac73SHong Zhang         }
300fcd7ac73SHong Zhang       }
301fcd7ac73SHong Zhang     }
302*4b2e90caSHong Zhang     c->nrows[i] = nrows_i;
3038bc97078SHong Zhang 
304d3825b63SHong Zhang     for (j=0; j<m; j++) {
305fcd7ac73SHong Zhang       if (rowhit[j]) {
306573f477fSHong Zhang         Jentry[nz].row     = j;              /* local row index */
307573f477fSHong Zhang         Jentry[nz].col     = rowhit[j] - 1;  /* local column index */
308573f477fSHong Zhang         Jentry[nz].valaddr = valaddrhit[j];  /* address of mat value for this entry */
309573f477fSHong Zhang         nz++;
310fcd7ac73SHong Zhang       }
311fcd7ac73SHong Zhang     }
312fcd7ac73SHong Zhang     ierr = PetscFree(cols);CHKERRQ(ierr);
313fcd7ac73SHong Zhang   }
3148bc97078SHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
3158bc97078SHong 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);
316fcd7ac73SHong Zhang 
317d3825b63SHong Zhang   ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr);
318fcd7ac73SHong Zhang   ierr = MatRestoreColumnIJ_SeqAIJ_Color(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
319fcd7ac73SHong Zhang   ierr = MatRestoreColumnIJ_SeqAIJ_Color(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
32072c15787SHong Zhang   if (ctype == IS_COLORING_GHOSTED) {
32172c15787SHong Zhang     ierr = ISLocalToGlobalMappingRestoreIndices(map,&ltog);CHKERRQ(ierr);
32272c15787SHong Zhang   }
323fcd7ac73SHong Zhang 
324fcd7ac73SHong Zhang   mat->ops->fdcoloringapply = MatFDColoringApply_MPIAIJ;
325fcd7ac73SHong Zhang   PetscFunctionReturn(0);
326fcd7ac73SHong Zhang }
327fcd7ac73SHong Zhang 
328fcd7ac73SHong Zhang /*------------------------------------------------------*/
329fcd7ac73SHong Zhang #undef __FUNCT__
3304a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate_MPIAIJ"
331dfbe8321SBarry Smith PetscErrorCode MatFDColoringCreate_MPIAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
332a64fbb32SBarry Smith {
3336eaac0f3SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)mat->data;
3346849ba73SBarry Smith   PetscErrorCode         ierr;
335b1d57f15SBarry Smith   PetscMPIInt            size,*ncolsonproc,*disp,nn;
3361a83f524SJed Brown   PetscInt               i,n,nrows,j,k,m,ncols,col;
337afcb2eb5SJed Brown   const PetscInt         *is,*A_ci,*A_cj,*B_ci,*B_cj,*rows = 0,*ltog;
3381a83f524SJed Brown   PetscInt               nis = iscoloring->n,nctot,*cols;
339f6d58c54SBarry Smith   PetscInt               *rowhit,M,cstart,cend,colb;
340b1d57f15SBarry Smith   PetscInt               *columnsforrow,l;
341b9617806SBarry Smith   IS                     *isa;
342ace3abfcSBarry Smith   PetscBool              done,flg;
343992144d0SBarry Smith   ISLocalToGlobalMapping map   = mat->cmap->mapping;
344afcb2eb5SJed Brown   PetscInt               ctype=c->ctype;
345fcd7ac73SHong Zhang   PetscBool              new_impl=PETSC_FALSE;
346a64fbb32SBarry Smith 
3473a40ed3dSBarry Smith   PetscFunctionBegin;
348fcd7ac73SHong Zhang   ierr = PetscOptionsName("-new","using new impls","",&new_impl);CHKERRQ(ierr);
349fcd7ac73SHong Zhang   if (new_impl){
350fcd7ac73SHong Zhang     ierr =  MatFDColoringCreate_MPIAIJ_new(mat,iscoloring,c);CHKERRQ(ierr);
351fcd7ac73SHong Zhang     PetscFunctionReturn(0);
352fcd7ac73SHong Zhang   }
353ce94432eSBarry 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");
354522c5e43SBarry Smith 
355afcb2eb5SJed Brown   if (map) {ierr = ISLocalToGlobalMappingGetIndices(map,&ltog);CHKERRQ(ierr);}
356afcb2eb5SJed Brown   else     ltog = NULL;
357b9617806SBarry Smith   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
3583acb8795SBarry Smith 
359f6d58c54SBarry Smith   M         = mat->rmap->n;
360f6d58c54SBarry Smith   cstart    = mat->cmap->rstart;
361f6d58c54SBarry Smith   cend      = mat->cmap->rend;
362f6d58c54SBarry Smith   c->M      = mat->rmap->N;         /* set the global rows and columns and local rows */
363f6d58c54SBarry Smith   c->N      = mat->cmap->N;
364f6d58c54SBarry Smith   c->m      = mat->rmap->n;
365f6d58c54SBarry Smith   c->rstart = mat->rmap->rstart;
366005c665bSBarry Smith 
367a64fbb32SBarry Smith   c->ncolors = nis;
368b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
369b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr);
370b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
371b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->rows);CHKERRQ(ierr);
372b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columnsforrow);CHKERRQ(ierr);
3733bb1ff40SBarry Smith   ierr       = PetscLogObjectMemory((PetscObject)c,5*nis*sizeof(PetscInt));CHKERRQ(ierr);
3746eaac0f3SBarry Smith 
3756eaac0f3SBarry Smith   /* Allow access to data structures of local part of matrix */
3766eaac0f3SBarry Smith   if (!aij->colmap) {
377ab9863d7SBarry Smith     ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
3786eaac0f3SBarry Smith   }
3793acb8795SBarry Smith   ierr = MatGetColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr);
3803acb8795SBarry Smith   ierr = MatGetColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr);
3816eaac0f3SBarry Smith 
382b1d57f15SBarry Smith   ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr);
383b1d57f15SBarry Smith   ierr = PetscMalloc((M+1)*sizeof(PetscInt),&columnsforrow);CHKERRQ(ierr);
3846eaac0f3SBarry Smith 
385a64fbb32SBarry Smith   for (i=0; i<nis; i++) {
386b9b97703SBarry Smith     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
387a64fbb32SBarry Smith     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
3882205254eSKarl Rupp 
389fcd7ac73SHong Zhang     c->ncolumns[i] = n; /* local number of columns of this color on this process */
390a64fbb32SBarry Smith     if (n) {
391b1d57f15SBarry Smith       ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr);
3923bb1ff40SBarry Smith       ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr);
393b1d57f15SBarry Smith       ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr);
394a64fbb32SBarry Smith     } else {
395a64fbb32SBarry Smith       c->columns[i] = 0;
396a64fbb32SBarry Smith     }
397a64fbb32SBarry Smith 
3988ee2e534SBarry Smith     if (ctype == IS_COLORING_GLOBAL) {
3996eaac0f3SBarry Smith       /* Determine the total (parallel) number of columns of this color */
400ce94432eSBarry Smith       ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr);
401687f1162SBarry Smith       ierr = PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);CHKERRQ(ierr);
402b8f8c88eSHong Zhang 
4034dc2109aSBarry Smith       ierr  = PetscMPIIntCast(n,&nn);CHKERRQ(ierr);
404ce94432eSBarry Smith       ierr  = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
4052205254eSKarl Rupp       nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j];
4063a7fca6bSBarry Smith       if (!nctot) {
407ae15b995SBarry Smith         ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr);
4083a7fca6bSBarry Smith       }
4096eaac0f3SBarry Smith 
4106eaac0f3SBarry Smith       disp[0] = 0;
4116eaac0f3SBarry Smith       for (j=1; j<size; j++) {
4126eaac0f3SBarry Smith         disp[j] = disp[j-1] + ncolsonproc[j-1];
4136eaac0f3SBarry Smith       }
4146eaac0f3SBarry Smith 
4156eaac0f3SBarry Smith       /* Get complete list of columns for color on each processor */
416b1d57f15SBarry Smith       ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
417ce94432eSBarry Smith       ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
4181d79065fSBarry Smith       ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr);
419b8f8c88eSHong Zhang     } else if (ctype == IS_COLORING_GHOSTED) {
420b8f8c88eSHong Zhang       /* Determine local number of columns of this color on this process, including ghost points */
421b8f8c88eSHong Zhang       nctot = n;
422b8f8c88eSHong Zhang       ierr  = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
423b8f8c88eSHong Zhang       ierr  = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr);
424f23aa3ddSBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type");
4256eaac0f3SBarry Smith 
4266eaac0f3SBarry Smith     /*
4276eaac0f3SBarry Smith        Mark all rows affect by these columns
4286eaac0f3SBarry Smith     */
429b8f8c88eSHong Zhang     /* Temporary option to allow for debugging/testing */
43090d69ab7SBarry Smith     flg  = PETSC_FALSE;
4310298fd71SBarry Smith     ierr = PetscOptionsGetBool(NULL,"-matfdcoloring_slow",&flg,NULL);CHKERRQ(ierr);
432f158e583SBarry Smith     if (!flg) { /*-----------------------------------------------------------------------------*/
433f158e583SBarry Smith       /* crude, fast version */
434b1d57f15SBarry Smith       ierr = PetscMemzero(rowhit,M*sizeof(PetscInt));CHKERRQ(ierr);
435a64fbb32SBarry Smith       /* loop over columns*/
4366eaac0f3SBarry Smith       for (j=0; j<nctot; j++) {
437b8f8c88eSHong Zhang         if (ctype == IS_COLORING_GHOSTED) {
438b8f8c88eSHong Zhang           col = ltog[cols[j]];
439b8f8c88eSHong Zhang         } else {
4406eaac0f3SBarry Smith           col = cols[j];
441b8f8c88eSHong Zhang         }
4426eaac0f3SBarry Smith         if (col >= cstart && col < cend) {
4436eaac0f3SBarry Smith           /* column is in diagonal block of matrix */
4446eaac0f3SBarry Smith           rows = A_cj + A_ci[col-cstart];
4456eaac0f3SBarry Smith           m    = A_ci[col-cstart+1] - A_ci[col-cstart];
4466eaac0f3SBarry Smith         } else {
447aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
448cb9801acSJed Brown           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
449fa46199cSSatish Balay           colb--;
450b3d2dc96SSatish Balay #else
4516eaac0f3SBarry Smith           colb = aij->colmap[col] - 1;
452b3d2dc96SSatish Balay #endif
4536eaac0f3SBarry Smith           if (colb == -1) {
4546eaac0f3SBarry Smith             m = 0;
4556eaac0f3SBarry Smith           } else {
4566eaac0f3SBarry Smith             rows = B_cj + B_ci[colb];
4576eaac0f3SBarry Smith             m    = B_ci[colb+1] - B_ci[colb];
4586eaac0f3SBarry Smith           }
4596eaac0f3SBarry Smith         }
460a64fbb32SBarry Smith         /* loop over columns marking them in rowhit */
461a64fbb32SBarry Smith         for (k=0; k<m; k++) {
462a64fbb32SBarry Smith           rowhit[*rows++] = col + 1;
463a64fbb32SBarry Smith         }
464a64fbb32SBarry Smith       }
4656eaac0f3SBarry Smith 
466a64fbb32SBarry Smith       /* count the number of hits */
467a64fbb32SBarry Smith       nrows = 0;
4686eaac0f3SBarry Smith       for (j=0; j<M; j++) {
469a64fbb32SBarry Smith         if (rowhit[j]) nrows++;
470a64fbb32SBarry Smith       }
471a64fbb32SBarry Smith       c->nrows[i] = nrows;
472b1d57f15SBarry Smith       ierr        = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr);
473b1d57f15SBarry Smith       ierr        = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr);
4743bb1ff40SBarry Smith       ierr        = PetscLogObjectMemory((PetscObject)c,2*(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr);
475a64fbb32SBarry Smith       nrows       = 0;
4766eaac0f3SBarry Smith       for (j=0; j<M; j++) {
477a64fbb32SBarry Smith         if (rowhit[j]) {
478fcd7ac73SHong Zhang           c->rows[i][nrows]          = j;              /* local row index */
479fcd7ac73SHong Zhang           c->columnsforrow[i][nrows] = rowhit[j] - 1;  /* global column index */
480a64fbb32SBarry Smith           nrows++;
481a64fbb32SBarry Smith         }
482a64fbb32SBarry Smith       }
483a64fbb32SBarry Smith     } else { /*-------------------------------------------------------------------------------*/
484f158e583SBarry Smith       /* slow version, using rowhit as a linked list */
485b1d57f15SBarry Smith       PetscInt currentcol,fm,mfm;
4866eaac0f3SBarry Smith       rowhit[M] = M;
487a64fbb32SBarry Smith       nrows     = 0;
488a64fbb32SBarry Smith       /* loop over columns*/
4896eaac0f3SBarry Smith       for (j=0; j<nctot; j++) {
490b8f8c88eSHong Zhang         if (ctype == IS_COLORING_GHOSTED) {
491b8f8c88eSHong Zhang           col = ltog[cols[j]];
492b8f8c88eSHong Zhang         } else {
4936eaac0f3SBarry Smith           col = cols[j];
494b8f8c88eSHong Zhang         }
4956eaac0f3SBarry Smith         if (col >= cstart && col < cend) {
4966eaac0f3SBarry Smith           /* column is in diagonal block of matrix */
4976eaac0f3SBarry Smith           rows = A_cj + A_ci[col-cstart];
4986eaac0f3SBarry Smith           m    = A_ci[col-cstart+1] - A_ci[col-cstart];
4996eaac0f3SBarry Smith         } else {
500aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
5010f5bd95cSBarry Smith           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
502fa46199cSSatish Balay           colb--;
503b3d2dc96SSatish Balay #else
5046eaac0f3SBarry Smith           colb = aij->colmap[col] - 1;
505b3d2dc96SSatish Balay #endif
5066eaac0f3SBarry Smith           if (colb == -1) {
5076eaac0f3SBarry Smith             m = 0;
5086eaac0f3SBarry Smith           } else {
5096eaac0f3SBarry Smith             rows = B_cj + B_ci[colb];
5106eaac0f3SBarry Smith             m    = B_ci[colb+1] - B_ci[colb];
5116eaac0f3SBarry Smith           }
5126eaac0f3SBarry Smith         }
513b8f8c88eSHong Zhang 
514a64fbb32SBarry Smith         /* loop over columns marking them in rowhit */
5156eaac0f3SBarry Smith         fm = M;    /* fm points to first entry in linked list */
516a64fbb32SBarry Smith         for (k=0; k<m; k++) {
517a64fbb32SBarry Smith           currentcol = *rows++;
518a64fbb32SBarry Smith           /* is it already in the list? */
519a64fbb32SBarry Smith           do {
520a64fbb32SBarry Smith             mfm = fm;
521a64fbb32SBarry Smith             fm  = rowhit[fm];
522a64fbb32SBarry Smith           } while (fm < currentcol);
523a64fbb32SBarry Smith           /* not in list so add it */
524a64fbb32SBarry Smith           if (fm != currentcol) {
525a64fbb32SBarry Smith             nrows++;
526a64fbb32SBarry Smith             columnsforrow[currentcol] = col;
527a64fbb32SBarry Smith             /* next three lines insert new entry into linked list */
528a64fbb32SBarry Smith             rowhit[mfm]        = currentcol;
529a64fbb32SBarry Smith             rowhit[currentcol] = fm;
530a64fbb32SBarry Smith             fm                 = currentcol;
531a64fbb32SBarry Smith             /* fm points to present position in list since we know the columns are sorted */
532f23aa3ddSBarry Smith           } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid coloring of matrix detected");
533a64fbb32SBarry Smith         }
534a64fbb32SBarry Smith       }
535a64fbb32SBarry Smith       c->nrows[i] = nrows;
5362205254eSKarl Rupp 
537b1d57f15SBarry Smith       ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr);
538b1d57f15SBarry Smith       ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr);
5393bb1ff40SBarry Smith       ierr = PetscLogObjectMemory((PetscObject)c,(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr);
540a64fbb32SBarry Smith       /* now store the linked list of rows into c->rows[i] */
541a64fbb32SBarry Smith       nrows = 0;
5426eaac0f3SBarry Smith       fm    = rowhit[M];
543a64fbb32SBarry Smith       do {
544a64fbb32SBarry Smith         c->rows[i][nrows]            = fm;
545a64fbb32SBarry Smith         c->columnsforrow[i][nrows++] = columnsforrow[fm];
546a64fbb32SBarry Smith         fm                           = rowhit[fm];
5476eaac0f3SBarry Smith       } while (fm < M);
5486eaac0f3SBarry Smith     } /* ---------------------------------------------------------------------------------------*/
549606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
5506eaac0f3SBarry Smith   }
55130b34957SBarry Smith 
55230b34957SBarry Smith   /* Optimize by adding the vscale, and scaleforrow[][] fields */
55330b34957SBarry Smith   /*
55430b34957SBarry Smith        vscale will contain the "diagonal" on processor scalings followed by the off processor
55530b34957SBarry Smith   */
5568ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
557ce94432eSBarry Smith     ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),aij->A->rmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&c->vscale);CHKERRQ(ierr);
558b1d57f15SBarry Smith     ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr);
55930b34957SBarry Smith     for (k=0; k<c->ncolors; k++) {
560b1d57f15SBarry Smith       ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr);
56130b34957SBarry Smith       for (l=0; l<c->nrows[k]; l++) {
56230b34957SBarry Smith         col = c->columnsforrow[k][l];
56330b34957SBarry Smith         if (col >= cstart && col < cend) {
56430b34957SBarry Smith           /* column is in diagonal block of matrix */
56530b34957SBarry Smith           colb = col - cstart;
56630b34957SBarry Smith         } else {
56730b34957SBarry Smith           /* column  is in "off-processor" part */
56830b34957SBarry Smith #if defined(PETSC_USE_CTABLE)
56930b34957SBarry Smith           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
57030b34957SBarry Smith           colb--;
57130b34957SBarry Smith #else
57230b34957SBarry Smith           colb = aij->colmap[col] - 1;
57330b34957SBarry Smith #endif
57430b34957SBarry Smith           colb += cend - cstart;
57530b34957SBarry Smith         }
57630b34957SBarry Smith         c->vscaleforrow[k][l] = colb;
57730b34957SBarry Smith       }
57830b34957SBarry Smith     }
579b8f8c88eSHong Zhang   } else if (ctype == IS_COLORING_GHOSTED) {
580b8f8c88eSHong Zhang     /* Get gtol mapping */
581afcb2eb5SJed Brown     PetscInt N = mat->cmap->N,nlocal,*gtol;
582b8f8c88eSHong Zhang     ierr = PetscMalloc((N+1)*sizeof(PetscInt),&gtol);CHKERRQ(ierr);
583b8f8c88eSHong Zhang     for (i=0; i<N; i++) gtol[i] = -1;
584afcb2eb5SJed Brown     ierr = ISLocalToGlobalMappingGetSize(map,&nlocal);CHKERRQ(ierr);
585afcb2eb5SJed Brown     for (i=0; i<nlocal; i++) gtol[ltog[i]] = i;
586b8f8c88eSHong Zhang 
587b8f8c88eSHong Zhang     c->vscale = 0; /* will be created in MatFDColoringApply() */
588b8f8c88eSHong Zhang     ierr      = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr);
589b8f8c88eSHong Zhang     for (k=0; k<c->ncolors; k++) {
590b8f8c88eSHong Zhang       ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr);
591b8f8c88eSHong Zhang       for (l=0; l<c->nrows[k]; l++) {
592b8f8c88eSHong Zhang         col = c->columnsforrow[k][l];      /* global column index */
593b8f8c88eSHong Zhang         c->vscaleforrow[k][l] = gtol[col]; /* local column index */
594b8f8c88eSHong Zhang       }
595b8f8c88eSHong Zhang     }
596b8f8c88eSHong Zhang     ierr = PetscFree(gtol);CHKERRQ(ierr);
597b8f8c88eSHong Zhang   }
598b9617806SBarry Smith   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
59930b34957SBarry Smith 
600606d414cSSatish Balay   ierr = PetscFree(rowhit);CHKERRQ(ierr);
601606d414cSSatish Balay   ierr = PetscFree(columnsforrow);CHKERRQ(ierr);
6023acb8795SBarry Smith   ierr = MatRestoreColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr);
6033acb8795SBarry Smith   ierr = MatRestoreColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr);
604afcb2eb5SJed Brown   if (map) {ierr = ISLocalToGlobalMappingRestoreIndices(map,&ltog);CHKERRQ(ierr);}
6053a40ed3dSBarry Smith   PetscFunctionReturn(0);
606a64fbb32SBarry Smith }
607a64fbb32SBarry Smith 
608b9617806SBarry Smith 
609b9617806SBarry Smith 
610b9617806SBarry Smith 
611b9617806SBarry Smith 
612b9617806SBarry Smith 
613