xref: /petsc/src/mat/impls/aij/mpi/fdmpiaij.c (revision 74d3cef9fb3a51ffeab1b96bc8a83dd2f44a1f02)
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;
1270e7395fSHong Zhang   PetscInt       k,start,end,l,row,col,colb,**columnsforrow=coloring->columnsforrow,nz;
13fcd7ac73SHong Zhang   PetscScalar    dx,*y,*xx,*w3_array;
14fcd7ac73SHong Zhang   PetscScalar    *vscale_array;
15fcd7ac73SHong Zhang   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin,unorm;
16fcd7ac73SHong Zhang   Vec            w1      = coloring->w1,w2=coloring->w2,w3;
17fcd7ac73SHong Zhang   void           *fctx   = coloring->fctx;
18fcd7ac73SHong Zhang   PetscBool      flg     = PETSC_FALSE;
19*74d3cef9SHong Zhang   PetscInt       ctype   = coloring->ctype,nxloc;
20f5aae955SHong Zhang   Mat_MPIAIJ     *aij = (Mat_MPIAIJ*)J->data;
21f5aae955SHong Zhang   Mat            A=aij->A,B=aij->B;
22f5aae955SHong Zhang   Mat_SeqAIJ     *cspA=(Mat_SeqAIJ*)A->data, *cspB=(Mat_SeqAIJ*)B->data;
23f6af9589SHong Zhang   PetscMPIInt    rank;
24fcd7ac73SHong Zhang 
25fcd7ac73SHong Zhang   PetscFunctionBegin;
26f6af9589SHong Zhang   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)J), &rank);CHKERRQ(ierr);
27fcd7ac73SHong Zhang   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
28fcd7ac73SHong Zhang   ierr = PetscOptionsGetBool(NULL,"-mat_fd_coloring_dont_rezero",&flg,NULL);CHKERRQ(ierr);
29fcd7ac73SHong Zhang   if (flg) {
30fcd7ac73SHong Zhang     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
31fcd7ac73SHong Zhang   } else {
32fcd7ac73SHong Zhang     PetscBool assembled;
33fcd7ac73SHong Zhang     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
34fcd7ac73SHong Zhang     if (assembled) {
35fcd7ac73SHong Zhang       ierr = MatZeroEntries(J);CHKERRQ(ierr);
36fcd7ac73SHong Zhang     }
37fcd7ac73SHong Zhang   }
38fcd7ac73SHong Zhang 
39f6af9589SHong Zhang   if (!coloring->vscale) { //ctype == IS_COLORING_GHOSTED ?
40f6af9589SHong Zhang     printf("[%d] create a non-ghostedcoloring->vscale\n",rank);
41f6af9589SHong Zhang     ierr = VecDuplicate(x1,&coloring->vscale);CHKERRQ(ierr);
42f6af9589SHong Zhang     //SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_NULL,"Null Object: coloring->scale");
43fcd7ac73SHong Zhang   }
44fcd7ac73SHong Zhang 
45fcd7ac73SHong Zhang   /* Set w1 = F(x1) */
46fcd7ac73SHong Zhang   if (!coloring->fset) {
47fcd7ac73SHong Zhang     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
48f6af9589SHong Zhang     ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
49fcd7ac73SHong Zhang     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
50fcd7ac73SHong Zhang   } else {
51fcd7ac73SHong Zhang     coloring->fset = PETSC_FALSE;
52fcd7ac73SHong Zhang   }
53fcd7ac73SHong Zhang 
54fcd7ac73SHong Zhang   if (!coloring->w3) {
55f6af9589SHong Zhang     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
56fcd7ac73SHong Zhang     ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr);
57fcd7ac73SHong Zhang   }
58fcd7ac73SHong Zhang   w3 = coloring->w3;
59fcd7ac73SHong Zhang 
60*74d3cef9SHong Zhang   /* Compute vscale = 1./dx - the local scale factors, including ghost points */
61f6af9589SHong Zhang   ierr = VecGetOwnershipRange(w1,&start,&end);CHKERRQ(ierr); /* =J's row ownershipRange, used by ghosted x! */
62f6af9589SHong Zhang   //printf("[%d] start end %d %d\n",rank,start,end);
63f6af9589SHong Zhang 
64f6af9589SHong Zhang   if (coloring->htype[0] == 'w') {
65f6af9589SHong Zhang     /* tacky test; need to make systematic if we add other approaches to computing h*/
66f6af9589SHong Zhang     ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr);
67*74d3cef9SHong Zhang     dx = 1.0/((1.0 + unorm)*epsilon);
68f6af9589SHong Zhang     ierr = VecSet(coloring->vscale,dx);CHKERRQ(ierr);
6970e7395fSHong Zhang   } else {
70*74d3cef9SHong Zhang     ierr = VecGetLocalSize(x1,&nxloc);CHKERRQ(ierr);
71f6af9589SHong Zhang     ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);
72fcd7ac73SHong Zhang     ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
73*74d3cef9SHong Zhang     for (col=0; col<nxloc; col++) {
74fcd7ac73SHong Zhang       dx = xx[col];
75*74d3cef9SHong Zhang       if (PetscAbsScalar(dx) < umin) {
76*74d3cef9SHong Zhang         if (PetscRealPart(dx) >= 0.0)      dx = umin;
77*74d3cef9SHong Zhang         else if (PetscRealPart(dx) < 0.0 ) dx = -umin;
78*74d3cef9SHong Zhang       }
79fcd7ac73SHong Zhang       dx               *= epsilon;
80*74d3cef9SHong Zhang       vscale_array[col] = 1.0/dx;
81f6af9589SHong Zhang     }
82*74d3cef9SHong Zhang     ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
83fcd7ac73SHong Zhang     ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
8470e7395fSHong Zhang   }
85fcd7ac73SHong Zhang   if (ctype == IS_COLORING_GLOBAL) {
86fcd7ac73SHong Zhang     ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
87fcd7ac73SHong Zhang     ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
88fcd7ac73SHong Zhang   }
89fcd7ac73SHong Zhang 
901b97d346SHong Zhang   /* Loop over each color */
91f5aae955SHong Zhang   nz = 0;
92fcd7ac73SHong Zhang   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
93fcd7ac73SHong Zhang   for (k=0; k<coloring->ncolors; k++) {
94fcd7ac73SHong Zhang     coloring->currentcolor = k;
95fcd7ac73SHong Zhang 
96fcd7ac73SHong Zhang     /*
97fcd7ac73SHong Zhang       Loop over each column associated with color
98f6af9589SHong Zhang       adding the perturbation to the vector w3 = x1 + dx.
99fcd7ac73SHong Zhang     */
100f6af9589SHong Zhang     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
101f6af9589SHong Zhang     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
102f6af9589SHong Zhang     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array - start;
10370e7395fSHong Zhang 
104fcd7ac73SHong Zhang     for (l=0; l<coloring->ncolumns[k]; l++) {
105f6af9589SHong Zhang       col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
10670e7395fSHong Zhang       /* convert global col to local colb */
10770e7395fSHong Zhang       if (col >= start && col < end) {
10870e7395fSHong Zhang         colb = col-start;
10970e7395fSHong Zhang       } else {
11070e7395fSHong Zhang #if defined(PETSC_USE_CTABLE)
11170e7395fSHong Zhang         ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
11270e7395fSHong Zhang         colb--;
11370e7395fSHong Zhang #else
11470e7395fSHong Zhang         colb = aij->colmap[col] - 1; /* local column index */
11570e7395fSHong Zhang #endif
11670e7395fSHong Zhang         colb += end - start;
11770e7395fSHong Zhang       }
118*74d3cef9SHong Zhang       w3_array[col] += 1.0/vscale_array[colb];
119fcd7ac73SHong Zhang     }
120fcd7ac73SHong Zhang     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array + start;
121fcd7ac73SHong Zhang     ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
122fcd7ac73SHong Zhang 
123fcd7ac73SHong Zhang     /*
124fcd7ac73SHong Zhang       Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
125fcd7ac73SHong Zhang                            w2 = F(x1 + dx) - F(x1)
126fcd7ac73SHong Zhang     */
127fcd7ac73SHong Zhang     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
128fcd7ac73SHong Zhang     ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
129fcd7ac73SHong Zhang     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
130fcd7ac73SHong Zhang     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
131fcd7ac73SHong Zhang 
1321b97d346SHong Zhang     /* Loop over rows of vector, putting results into Jacobian matrix */
133fcd7ac73SHong Zhang     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
134fcd7ac73SHong Zhang     for (l=0; l<coloring->nrows[k]; l++) {
135fcd7ac73SHong Zhang       row                        = coloring->rows[k][l];   /* local row index */
136*74d3cef9SHong Zhang       *(coloring->valaddr[nz++]) = y[row]*vscale_array[columnsforrow[k][l]];
137fcd7ac73SHong Zhang     }
138fcd7ac73SHong Zhang     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
139fcd7ac73SHong Zhang   } /* endof for each color */
140f5aae955SHong Zhang   if (nz != cspA->nz + cspB->nz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"nz %d != mat->nz %d\n",nz,cspA->nz+cspB->nz);
141f5aae955SHong Zhang   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
142f5aae955SHong Zhang   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
143fcd7ac73SHong Zhang   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
144fcd7ac73SHong Zhang 
145fcd7ac73SHong Zhang   coloring->currentcolor = -1;
146fcd7ac73SHong Zhang   PetscFunctionReturn(0);
147fcd7ac73SHong Zhang }
148fcd7ac73SHong Zhang 
149fcd7ac73SHong Zhang #undef __FUNCT__
150fcd7ac73SHong Zhang #define __FUNCT__ "MatFDColoringCreate_MPIAIJ_new"
151fcd7ac73SHong Zhang PetscErrorCode MatFDColoringCreate_MPIAIJ_new(Mat mat,ISColoring iscoloring,MatFDColoring c)
152fcd7ac73SHong Zhang {
153fcd7ac73SHong Zhang   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)mat->data;
154fcd7ac73SHong Zhang   PetscErrorCode         ierr;
155fcd7ac73SHong Zhang   PetscMPIInt            size,*ncolsonproc,*disp,nn;
156fcd7ac73SHong Zhang   PetscInt               i,n,nrows,j,k,m,ncols,col;
15772c15787SHong Zhang   const PetscInt         *is,*A_ci,*A_cj,*B_ci,*B_cj,*row = NULL,*ltog=NULL;
158fcd7ac73SHong Zhang   PetscInt               nis = iscoloring->n,nctot,*cols;
159d3825b63SHong Zhang   PetscInt               *rowhit,cstart,cend,colb;
160fcd7ac73SHong Zhang   IS                     *isa;
161fcd7ac73SHong Zhang   ISLocalToGlobalMapping map = mat->cmap->mapping;
162fcd7ac73SHong Zhang   PetscInt               ctype=c->ctype;
163fcd7ac73SHong Zhang   Mat                    A=aij->A,B=aij->B;
164fcd7ac73SHong Zhang   Mat_SeqAIJ             *cspA=(Mat_SeqAIJ*)A->data, *cspB=(Mat_SeqAIJ*)B->data;
165fcd7ac73SHong Zhang   PetscScalar            *A_val=cspA->a,*B_val=cspB->a;
166fcd7ac73SHong Zhang   PetscInt               spidx;
1671b97d346SHong Zhang   PetscInt               *spidxA,*spidxB,nz;
168d3825b63SHong Zhang   PetscScalar            **valaddrhit;
169fcd7ac73SHong Zhang 
170fcd7ac73SHong Zhang   PetscFunctionBegin;
17172c15787SHong Zhang   if (ctype == IS_COLORING_GHOSTED) {
17272c15787SHong 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");
17372c15787SHong Zhang     ierr = ISLocalToGlobalMappingGetIndices(map,&ltog);CHKERRQ(ierr);
17472c15787SHong Zhang   }
175fcd7ac73SHong Zhang 
176d3825b63SHong Zhang   m         = mat->rmap->n;
177fcd7ac73SHong Zhang   cstart    = mat->cmap->rstart;
178fcd7ac73SHong Zhang   cend      = mat->cmap->rend;
179fcd7ac73SHong Zhang   c->M      = mat->rmap->N;         /* set the global rows and columns and local rows */
180fcd7ac73SHong Zhang   c->N      = mat->cmap->N;
181d3825b63SHong Zhang   c->m      = m;
182fcd7ac73SHong Zhang   c->rstart = mat->rmap->rstart;
183fcd7ac73SHong Zhang 
184fcd7ac73SHong Zhang   c->ncolors = nis;
185fcd7ac73SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
186fcd7ac73SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr);
187fcd7ac73SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
188fcd7ac73SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->rows);CHKERRQ(ierr);
189fcd7ac73SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columnsforrow);CHKERRQ(ierr);
190fcd7ac73SHong Zhang   ierr       = PetscLogObjectMemory((PetscObject)c,5*nis*sizeof(PetscInt));CHKERRQ(ierr);
191fcd7ac73SHong Zhang 
192a774a6f1SHong Zhang   nz         = cspA->nz + cspB->nz; /* total nonzero entries of mat */
193a774a6f1SHong Zhang   ierr       = PetscMalloc(nz*sizeof(PetscScalar*),&c->valaddr);CHKERRQ(ierr);
194a774a6f1SHong Zhang   ierr       = PetscLogObjectMemory((PetscObject)c,nz*sizeof(PetscScalar*));CHKERRQ(ierr);
195a774a6f1SHong Zhang 
196d3825b63SHong Zhang   /* Allow access to data structures of local part of matrix
197d3825b63SHong Zhang    - creates aij->colmap which maps global column number to local number in part B */
198fcd7ac73SHong Zhang   if (!aij->colmap) {
199fcd7ac73SHong Zhang     ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
200fcd7ac73SHong Zhang   }
201fcd7ac73SHong Zhang   ierr = MatGetColumnIJ_SeqAIJ_Color(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
202fcd7ac73SHong Zhang   ierr = MatGetColumnIJ_SeqAIJ_Color(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
203fcd7ac73SHong Zhang 
204d3825b63SHong Zhang   ierr = PetscMalloc2(m+1,PetscInt,&rowhit,m+1,PetscScalar*,&valaddrhit);CHKERRQ(ierr);
205fcd7ac73SHong Zhang   nz = 0;
20672c15787SHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
207d3825b63SHong Zhang   for (i=0; i<nis; i++) { /* for each local color */
208fcd7ac73SHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
209fcd7ac73SHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
210fcd7ac73SHong Zhang 
211fcd7ac73SHong Zhang     c->ncolumns[i] = n; /* local number of columns of this color on this process */
212fcd7ac73SHong Zhang     if (n) {
213fcd7ac73SHong Zhang       ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr);
214fcd7ac73SHong Zhang       ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr);
215fcd7ac73SHong Zhang       ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr);
216fcd7ac73SHong Zhang     } else {
217fcd7ac73SHong Zhang       c->columns[i] = 0;
218fcd7ac73SHong Zhang     }
219fcd7ac73SHong Zhang 
220fcd7ac73SHong Zhang     if (ctype == IS_COLORING_GLOBAL) {
221d3825b63SHong Zhang       /* Determine nctot, the total (parallel) number of columns of this color */
222fcd7ac73SHong Zhang       ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr);
223fcd7ac73SHong Zhang       ierr = PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);CHKERRQ(ierr);
224fcd7ac73SHong Zhang 
225d3825b63SHong Zhang       /* ncolsonproc[j]: local ncolumns on proc[j] of this color */
226fcd7ac73SHong Zhang       ierr  = PetscMPIIntCast(n,&nn);CHKERRQ(ierr);
227fcd7ac73SHong Zhang       ierr  = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
228fcd7ac73SHong Zhang       nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j];
229fcd7ac73SHong Zhang       if (!nctot) {
230fcd7ac73SHong Zhang         ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr);
231fcd7ac73SHong Zhang       }
232fcd7ac73SHong Zhang 
233fcd7ac73SHong Zhang       disp[0] = 0;
234fcd7ac73SHong Zhang       for (j=1; j<size; j++) {
235fcd7ac73SHong Zhang         disp[j] = disp[j-1] + ncolsonproc[j-1];
236fcd7ac73SHong Zhang       }
237fcd7ac73SHong Zhang 
238d3825b63SHong Zhang       /* Get cols, the complete list of columns for this color on each process */
239fcd7ac73SHong Zhang       ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
240fcd7ac73SHong Zhang       ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
241fcd7ac73SHong Zhang       ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr);
242fcd7ac73SHong Zhang     } else if (ctype == IS_COLORING_GHOSTED) {
243fcd7ac73SHong Zhang       /* Determine local number of columns of this color on this process, including ghost points */
244fcd7ac73SHong Zhang       nctot = n;
245fcd7ac73SHong Zhang       ierr  = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
246fcd7ac73SHong Zhang       ierr  = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr);
247fcd7ac73SHong Zhang     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type");
248fcd7ac73SHong Zhang 
2491b97d346SHong Zhang     /* Mark all rows affect by these columns */
250d3825b63SHong Zhang     ierr = PetscMemzero(rowhit,m*sizeof(PetscInt));CHKERRQ(ierr);
2511b97d346SHong Zhang     for (j=0; j<nctot; j++) { /* loop over columns*/
252fcd7ac73SHong Zhang       if (ctype == IS_COLORING_GHOSTED) {
253fcd7ac73SHong Zhang         col = ltog[cols[j]];
254fcd7ac73SHong Zhang       } else {
255fcd7ac73SHong Zhang         col = cols[j];
256fcd7ac73SHong Zhang       }
257fcd7ac73SHong Zhang       if (col >= cstart && col < cend) { /* column is in diagonal block of matrix A */
258d3825b63SHong Zhang         row = A_cj + A_ci[col-cstart];
259d3825b63SHong Zhang         nrows = A_ci[col-cstart+1] - A_ci[col-cstart];
260fcd7ac73SHong Zhang         /* loop over columns of A marking them in rowhit */
261d3825b63SHong Zhang         for (k=0; k<nrows; k++) {
262d3825b63SHong Zhang           /* set valaddrhit for part A */
263fcd7ac73SHong Zhang           spidx            = spidxA[A_ci[col-cstart] + k];
264d3825b63SHong Zhang           valaddrhit[*row] = &A_val[spidx];
265a774a6f1SHong Zhang           rowhit[*row++]   = col - cstart + 1; /* local column index */
266fcd7ac73SHong Zhang         }
267fcd7ac73SHong Zhang       } else { /* column is in off-diagonal block of matrix B */
268fcd7ac73SHong Zhang #if defined(PETSC_USE_CTABLE)
269fcd7ac73SHong Zhang         ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
270fcd7ac73SHong Zhang         colb--;
271fcd7ac73SHong Zhang #else
272fcd7ac73SHong Zhang         colb = aij->colmap[col] - 1; /* local column index */
273fcd7ac73SHong Zhang #endif
274fcd7ac73SHong Zhang         if (colb == -1) {
275d3825b63SHong Zhang           nrows = 0;
276fcd7ac73SHong Zhang         } else {
277d3825b63SHong Zhang           row   = B_cj + B_ci[colb];
278d3825b63SHong Zhang           nrows = B_ci[colb+1] - B_ci[colb];
279fcd7ac73SHong Zhang         }
280fcd7ac73SHong Zhang         /* loop over columns of B marking them in rowhit */
281d3825b63SHong Zhang         for (k=0; k<nrows; k++) {
282d3825b63SHong Zhang           /* set valaddrhit for part B */
283fcd7ac73SHong Zhang           spidx            = spidxB[B_ci[colb] + k];
284d3825b63SHong Zhang           valaddrhit[*row] = &B_val[spidx];
28570e7395fSHong Zhang           rowhit[*row++]   = colb + 1 + cend - cstart; /* local column index */
286fcd7ac73SHong Zhang         }
287fcd7ac73SHong Zhang       }
288fcd7ac73SHong Zhang     }
289fcd7ac73SHong Zhang 
290fcd7ac73SHong Zhang     /* count the number of hits */
291fcd7ac73SHong Zhang     nrows = 0;
292d3825b63SHong Zhang     for (j=0; j<m; j++) {
293fcd7ac73SHong Zhang       if (rowhit[j]) nrows++;
294fcd7ac73SHong Zhang     }
295fcd7ac73SHong Zhang     c->nrows[i] = nrows;
296fcd7ac73SHong Zhang     ierr        = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr);
297fcd7ac73SHong Zhang     ierr        = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr);
298fcd7ac73SHong Zhang     ierr        = PetscLogObjectMemory((PetscObject)c,2*(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr);
299fcd7ac73SHong Zhang     nrows       = 0;
300d3825b63SHong Zhang     for (j=0; j<m; j++) {
301fcd7ac73SHong Zhang       if (rowhit[j]) {
302fcd7ac73SHong Zhang         c->rows[i][nrows]          = j;              /* local row index */
303a774a6f1SHong Zhang         c->columnsforrow[i][nrows] = rowhit[j] - 1;  /* local column index */
304d3825b63SHong Zhang         c->valaddr[nz++]           = valaddrhit[j];  /* address of mat value for this entry */
305fcd7ac73SHong Zhang         nrows++;
306fcd7ac73SHong Zhang       }
307fcd7ac73SHong Zhang     }
308fcd7ac73SHong Zhang     ierr = PetscFree(cols);CHKERRQ(ierr);
309fcd7ac73SHong Zhang   }
310fcd7ac73SHong Zhang   if (nz != cspA->nz + cspB->nz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"nz %d != mat->nz %d\n",nz,cspA->nz+cspB->nz);
311fcd7ac73SHong Zhang 
312a774a6f1SHong Zhang   /* create vscale for storing dx */
313fcd7ac73SHong Zhang   if (ctype == IS_COLORING_GLOBAL) {
314d3825b63SHong Zhang     ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),m,PETSC_DETERMINE,B->cmap->n,aij->garray,&c->vscale);CHKERRQ(ierr);
315fcd7ac73SHong Zhang   }
316fcd7ac73SHong Zhang 
317fcd7ac73SHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
318fcd7ac73SHong Zhang 
319d3825b63SHong Zhang   ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr);
320fcd7ac73SHong Zhang   ierr = MatRestoreColumnIJ_SeqAIJ_Color(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
321fcd7ac73SHong Zhang   ierr = MatRestoreColumnIJ_SeqAIJ_Color(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
32272c15787SHong Zhang   if (ctype == IS_COLORING_GHOSTED) {
32372c15787SHong Zhang     ierr = ISLocalToGlobalMappingRestoreIndices(map,&ltog);CHKERRQ(ierr);
32472c15787SHong Zhang   }
325fcd7ac73SHong Zhang 
326fcd7ac73SHong Zhang   mat->ops->fdcoloringapply = MatFDColoringApply_MPIAIJ;
327fcd7ac73SHong Zhang   PetscFunctionReturn(0);
328fcd7ac73SHong Zhang }
329fcd7ac73SHong Zhang 
330fcd7ac73SHong Zhang /*------------------------------------------------------*/
331fcd7ac73SHong Zhang #undef __FUNCT__
3324a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate_MPIAIJ"
333dfbe8321SBarry Smith PetscErrorCode MatFDColoringCreate_MPIAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
334a64fbb32SBarry Smith {
3356eaac0f3SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)mat->data;
3366849ba73SBarry Smith   PetscErrorCode         ierr;
337b1d57f15SBarry Smith   PetscMPIInt            size,*ncolsonproc,*disp,nn;
3381a83f524SJed Brown   PetscInt               i,n,nrows,j,k,m,ncols,col;
339afcb2eb5SJed Brown   const PetscInt         *is,*A_ci,*A_cj,*B_ci,*B_cj,*rows = 0,*ltog;
3401a83f524SJed Brown   PetscInt               nis = iscoloring->n,nctot,*cols;
341f6d58c54SBarry Smith   PetscInt               *rowhit,M,cstart,cend,colb;
342b1d57f15SBarry Smith   PetscInt               *columnsforrow,l;
343b9617806SBarry Smith   IS                     *isa;
344ace3abfcSBarry Smith   PetscBool              done,flg;
345992144d0SBarry Smith   ISLocalToGlobalMapping map   = mat->cmap->mapping;
346afcb2eb5SJed Brown   PetscInt               ctype=c->ctype;
347fcd7ac73SHong Zhang   PetscBool              new_impl=PETSC_FALSE;
348a64fbb32SBarry Smith 
3493a40ed3dSBarry Smith   PetscFunctionBegin;
350fcd7ac73SHong Zhang   ierr = PetscOptionsName("-new","using new impls","",&new_impl);CHKERRQ(ierr);
351fcd7ac73SHong Zhang   if (new_impl){
352fcd7ac73SHong Zhang     ierr =  MatFDColoringCreate_MPIAIJ_new(mat,iscoloring,c);CHKERRQ(ierr);
353fcd7ac73SHong Zhang     PetscFunctionReturn(0);
354fcd7ac73SHong Zhang   }
355ce94432eSBarry 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");
356522c5e43SBarry Smith 
357afcb2eb5SJed Brown   if (map) {ierr = ISLocalToGlobalMappingGetIndices(map,&ltog);CHKERRQ(ierr);}
358afcb2eb5SJed Brown   else     ltog = NULL;
359b9617806SBarry Smith   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
3603acb8795SBarry Smith 
361f6d58c54SBarry Smith   M         = mat->rmap->n;
362f6d58c54SBarry Smith   cstart    = mat->cmap->rstart;
363f6d58c54SBarry Smith   cend      = mat->cmap->rend;
364f6d58c54SBarry Smith   c->M      = mat->rmap->N;         /* set the global rows and columns and local rows */
365f6d58c54SBarry Smith   c->N      = mat->cmap->N;
366f6d58c54SBarry Smith   c->m      = mat->rmap->n;
367f6d58c54SBarry Smith   c->rstart = mat->rmap->rstart;
368005c665bSBarry Smith 
369a64fbb32SBarry Smith   c->ncolors = nis;
370b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
371b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr);
372b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
373b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->rows);CHKERRQ(ierr);
374b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columnsforrow);CHKERRQ(ierr);
3753bb1ff40SBarry Smith   ierr       = PetscLogObjectMemory((PetscObject)c,5*nis*sizeof(PetscInt));CHKERRQ(ierr);
3766eaac0f3SBarry Smith 
3776eaac0f3SBarry Smith   /* Allow access to data structures of local part of matrix */
3786eaac0f3SBarry Smith   if (!aij->colmap) {
379ab9863d7SBarry Smith     ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
3806eaac0f3SBarry Smith   }
3813acb8795SBarry Smith   ierr = MatGetColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr);
3823acb8795SBarry Smith   ierr = MatGetColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr);
3836eaac0f3SBarry Smith 
384b1d57f15SBarry Smith   ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr);
385b1d57f15SBarry Smith   ierr = PetscMalloc((M+1)*sizeof(PetscInt),&columnsforrow);CHKERRQ(ierr);
3866eaac0f3SBarry Smith 
387a64fbb32SBarry Smith   for (i=0; i<nis; i++) {
388b9b97703SBarry Smith     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
389a64fbb32SBarry Smith     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
3902205254eSKarl Rupp 
391fcd7ac73SHong Zhang     c->ncolumns[i] = n; /* local number of columns of this color on this process */
392a64fbb32SBarry Smith     if (n) {
393b1d57f15SBarry Smith       ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr);
3943bb1ff40SBarry Smith       ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr);
395b1d57f15SBarry Smith       ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr);
396a64fbb32SBarry Smith     } else {
397a64fbb32SBarry Smith       c->columns[i] = 0;
398a64fbb32SBarry Smith     }
399a64fbb32SBarry Smith 
4008ee2e534SBarry Smith     if (ctype == IS_COLORING_GLOBAL) {
4016eaac0f3SBarry Smith       /* Determine the total (parallel) number of columns of this color */
402ce94432eSBarry Smith       ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr);
403687f1162SBarry Smith       ierr = PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);CHKERRQ(ierr);
404b8f8c88eSHong Zhang 
4054dc2109aSBarry Smith       ierr  = PetscMPIIntCast(n,&nn);CHKERRQ(ierr);
406ce94432eSBarry Smith       ierr  = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
4072205254eSKarl Rupp       nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j];
4083a7fca6bSBarry Smith       if (!nctot) {
409ae15b995SBarry Smith         ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr);
4103a7fca6bSBarry Smith       }
4116eaac0f3SBarry Smith 
4126eaac0f3SBarry Smith       disp[0] = 0;
4136eaac0f3SBarry Smith       for (j=1; j<size; j++) {
4146eaac0f3SBarry Smith         disp[j] = disp[j-1] + ncolsonproc[j-1];
4156eaac0f3SBarry Smith       }
4166eaac0f3SBarry Smith 
4176eaac0f3SBarry Smith       /* Get complete list of columns for color on each processor */
418b1d57f15SBarry Smith       ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
419ce94432eSBarry Smith       ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
4201d79065fSBarry Smith       ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr);
421b8f8c88eSHong Zhang     } else if (ctype == IS_COLORING_GHOSTED) {
422b8f8c88eSHong Zhang       /* Determine local number of columns of this color on this process, including ghost points */
423b8f8c88eSHong Zhang       nctot = n;
424b8f8c88eSHong Zhang       ierr  = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
425b8f8c88eSHong Zhang       ierr  = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr);
426f23aa3ddSBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type");
4276eaac0f3SBarry Smith 
4286eaac0f3SBarry Smith     /*
4296eaac0f3SBarry Smith        Mark all rows affect by these columns
4306eaac0f3SBarry Smith     */
431b8f8c88eSHong Zhang     /* Temporary option to allow for debugging/testing */
43290d69ab7SBarry Smith     flg  = PETSC_FALSE;
4330298fd71SBarry Smith     ierr = PetscOptionsGetBool(NULL,"-matfdcoloring_slow",&flg,NULL);CHKERRQ(ierr);
434f158e583SBarry Smith     if (!flg) { /*-----------------------------------------------------------------------------*/
435f158e583SBarry Smith       /* crude, fast version */
436b1d57f15SBarry Smith       ierr = PetscMemzero(rowhit,M*sizeof(PetscInt));CHKERRQ(ierr);
437a64fbb32SBarry Smith       /* loop over columns*/
4386eaac0f3SBarry Smith       for (j=0; j<nctot; j++) {
439b8f8c88eSHong Zhang         if (ctype == IS_COLORING_GHOSTED) {
440b8f8c88eSHong Zhang           col = ltog[cols[j]];
441b8f8c88eSHong Zhang         } else {
4426eaac0f3SBarry Smith           col = cols[j];
443b8f8c88eSHong Zhang         }
4446eaac0f3SBarry Smith         if (col >= cstart && col < cend) {
4456eaac0f3SBarry Smith           /* column is in diagonal block of matrix */
4466eaac0f3SBarry Smith           rows = A_cj + A_ci[col-cstart];
4476eaac0f3SBarry Smith           m    = A_ci[col-cstart+1] - A_ci[col-cstart];
4486eaac0f3SBarry Smith         } else {
449aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
450cb9801acSJed Brown           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
451fa46199cSSatish Balay           colb--;
452b3d2dc96SSatish Balay #else
4536eaac0f3SBarry Smith           colb = aij->colmap[col] - 1;
454b3d2dc96SSatish Balay #endif
4556eaac0f3SBarry Smith           if (colb == -1) {
4566eaac0f3SBarry Smith             m = 0;
4576eaac0f3SBarry Smith           } else {
4586eaac0f3SBarry Smith             rows = B_cj + B_ci[colb];
4596eaac0f3SBarry Smith             m    = B_ci[colb+1] - B_ci[colb];
4606eaac0f3SBarry Smith           }
4616eaac0f3SBarry Smith         }
462a64fbb32SBarry Smith         /* loop over columns marking them in rowhit */
463a64fbb32SBarry Smith         for (k=0; k<m; k++) {
464a64fbb32SBarry Smith           rowhit[*rows++] = col + 1;
465a64fbb32SBarry Smith         }
466a64fbb32SBarry Smith       }
4676eaac0f3SBarry Smith 
468a64fbb32SBarry Smith       /* count the number of hits */
469a64fbb32SBarry Smith       nrows = 0;
4706eaac0f3SBarry Smith       for (j=0; j<M; j++) {
471a64fbb32SBarry Smith         if (rowhit[j]) nrows++;
472a64fbb32SBarry Smith       }
473a64fbb32SBarry Smith       c->nrows[i] = nrows;
474b1d57f15SBarry Smith       ierr        = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr);
475b1d57f15SBarry Smith       ierr        = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr);
4763bb1ff40SBarry Smith       ierr        = PetscLogObjectMemory((PetscObject)c,2*(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr);
477a64fbb32SBarry Smith       nrows       = 0;
4786eaac0f3SBarry Smith       for (j=0; j<M; j++) {
479a64fbb32SBarry Smith         if (rowhit[j]) {
480fcd7ac73SHong Zhang           c->rows[i][nrows]          = j;              /* local row index */
481fcd7ac73SHong Zhang           c->columnsforrow[i][nrows] = rowhit[j] - 1;  /* global column index */
482a64fbb32SBarry Smith           nrows++;
483a64fbb32SBarry Smith         }
484a64fbb32SBarry Smith       }
485a64fbb32SBarry Smith     } else { /*-------------------------------------------------------------------------------*/
486f158e583SBarry Smith       /* slow version, using rowhit as a linked list */
487b1d57f15SBarry Smith       PetscInt currentcol,fm,mfm;
4886eaac0f3SBarry Smith       rowhit[M] = M;
489a64fbb32SBarry Smith       nrows     = 0;
490a64fbb32SBarry Smith       /* loop over columns*/
4916eaac0f3SBarry Smith       for (j=0; j<nctot; j++) {
492b8f8c88eSHong Zhang         if (ctype == IS_COLORING_GHOSTED) {
493b8f8c88eSHong Zhang           col = ltog[cols[j]];
494b8f8c88eSHong Zhang         } else {
4956eaac0f3SBarry Smith           col = cols[j];
496b8f8c88eSHong Zhang         }
4976eaac0f3SBarry Smith         if (col >= cstart && col < cend) {
4986eaac0f3SBarry Smith           /* column is in diagonal block of matrix */
4996eaac0f3SBarry Smith           rows = A_cj + A_ci[col-cstart];
5006eaac0f3SBarry Smith           m    = A_ci[col-cstart+1] - A_ci[col-cstart];
5016eaac0f3SBarry Smith         } else {
502aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
5030f5bd95cSBarry Smith           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
504fa46199cSSatish Balay           colb--;
505b3d2dc96SSatish Balay #else
5066eaac0f3SBarry Smith           colb = aij->colmap[col] - 1;
507b3d2dc96SSatish Balay #endif
5086eaac0f3SBarry Smith           if (colb == -1) {
5096eaac0f3SBarry Smith             m = 0;
5106eaac0f3SBarry Smith           } else {
5116eaac0f3SBarry Smith             rows = B_cj + B_ci[colb];
5126eaac0f3SBarry Smith             m    = B_ci[colb+1] - B_ci[colb];
5136eaac0f3SBarry Smith           }
5146eaac0f3SBarry Smith         }
515b8f8c88eSHong Zhang 
516a64fbb32SBarry Smith         /* loop over columns marking them in rowhit */
5176eaac0f3SBarry Smith         fm = M;    /* fm points to first entry in linked list */
518a64fbb32SBarry Smith         for (k=0; k<m; k++) {
519a64fbb32SBarry Smith           currentcol = *rows++;
520a64fbb32SBarry Smith           /* is it already in the list? */
521a64fbb32SBarry Smith           do {
522a64fbb32SBarry Smith             mfm = fm;
523a64fbb32SBarry Smith             fm  = rowhit[fm];
524a64fbb32SBarry Smith           } while (fm < currentcol);
525a64fbb32SBarry Smith           /* not in list so add it */
526a64fbb32SBarry Smith           if (fm != currentcol) {
527a64fbb32SBarry Smith             nrows++;
528a64fbb32SBarry Smith             columnsforrow[currentcol] = col;
529a64fbb32SBarry Smith             /* next three lines insert new entry into linked list */
530a64fbb32SBarry Smith             rowhit[mfm]        = currentcol;
531a64fbb32SBarry Smith             rowhit[currentcol] = fm;
532a64fbb32SBarry Smith             fm                 = currentcol;
533a64fbb32SBarry Smith             /* fm points to present position in list since we know the columns are sorted */
534f23aa3ddSBarry Smith           } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid coloring of matrix detected");
535a64fbb32SBarry Smith         }
536a64fbb32SBarry Smith       }
537a64fbb32SBarry Smith       c->nrows[i] = nrows;
5382205254eSKarl Rupp 
539b1d57f15SBarry Smith       ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr);
540b1d57f15SBarry Smith       ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr);
5413bb1ff40SBarry Smith       ierr = PetscLogObjectMemory((PetscObject)c,(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr);
542a64fbb32SBarry Smith       /* now store the linked list of rows into c->rows[i] */
543a64fbb32SBarry Smith       nrows = 0;
5446eaac0f3SBarry Smith       fm    = rowhit[M];
545a64fbb32SBarry Smith       do {
546a64fbb32SBarry Smith         c->rows[i][nrows]            = fm;
547a64fbb32SBarry Smith         c->columnsforrow[i][nrows++] = columnsforrow[fm];
548a64fbb32SBarry Smith         fm                           = rowhit[fm];
5496eaac0f3SBarry Smith       } while (fm < M);
5506eaac0f3SBarry Smith     } /* ---------------------------------------------------------------------------------------*/
551606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
5526eaac0f3SBarry Smith   }
55330b34957SBarry Smith 
55430b34957SBarry Smith   /* Optimize by adding the vscale, and scaleforrow[][] fields */
55530b34957SBarry Smith   /*
55630b34957SBarry Smith        vscale will contain the "diagonal" on processor scalings followed by the off processor
55730b34957SBarry Smith   */
5588ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
559ce94432eSBarry Smith     ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),aij->A->rmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&c->vscale);CHKERRQ(ierr);
560b1d57f15SBarry Smith     ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr);
56130b34957SBarry Smith     for (k=0; k<c->ncolors; k++) {
562b1d57f15SBarry Smith       ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr);
56330b34957SBarry Smith       for (l=0; l<c->nrows[k]; l++) {
56430b34957SBarry Smith         col = c->columnsforrow[k][l];
56530b34957SBarry Smith         if (col >= cstart && col < cend) {
56630b34957SBarry Smith           /* column is in diagonal block of matrix */
56730b34957SBarry Smith           colb = col - cstart;
56830b34957SBarry Smith         } else {
56930b34957SBarry Smith           /* column  is in "off-processor" part */
57030b34957SBarry Smith #if defined(PETSC_USE_CTABLE)
57130b34957SBarry Smith           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
57230b34957SBarry Smith           colb--;
57330b34957SBarry Smith #else
57430b34957SBarry Smith           colb = aij->colmap[col] - 1;
57530b34957SBarry Smith #endif
57630b34957SBarry Smith           colb += cend - cstart;
57730b34957SBarry Smith         }
57830b34957SBarry Smith         c->vscaleforrow[k][l] = colb;
57930b34957SBarry Smith       }
58030b34957SBarry Smith     }
581b8f8c88eSHong Zhang   } else if (ctype == IS_COLORING_GHOSTED) {
582b8f8c88eSHong Zhang     /* Get gtol mapping */
583afcb2eb5SJed Brown     PetscInt N = mat->cmap->N,nlocal,*gtol;
584b8f8c88eSHong Zhang     ierr = PetscMalloc((N+1)*sizeof(PetscInt),&gtol);CHKERRQ(ierr);
585b8f8c88eSHong Zhang     for (i=0; i<N; i++) gtol[i] = -1;
586afcb2eb5SJed Brown     ierr = ISLocalToGlobalMappingGetSize(map,&nlocal);CHKERRQ(ierr);
587afcb2eb5SJed Brown     for (i=0; i<nlocal; i++) gtol[ltog[i]] = i;
588b8f8c88eSHong Zhang 
589b8f8c88eSHong Zhang     c->vscale = 0; /* will be created in MatFDColoringApply() */
590b8f8c88eSHong Zhang     ierr      = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr);
591b8f8c88eSHong Zhang     for (k=0; k<c->ncolors; k++) {
592b8f8c88eSHong Zhang       ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr);
593b8f8c88eSHong Zhang       for (l=0; l<c->nrows[k]; l++) {
594b8f8c88eSHong Zhang         col = c->columnsforrow[k][l];      /* global column index */
595b8f8c88eSHong Zhang         c->vscaleforrow[k][l] = gtol[col]; /* local column index */
596b8f8c88eSHong Zhang       }
597b8f8c88eSHong Zhang     }
598b8f8c88eSHong Zhang     ierr = PetscFree(gtol);CHKERRQ(ierr);
599b8f8c88eSHong Zhang   }
600b9617806SBarry Smith   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
60130b34957SBarry Smith 
602606d414cSSatish Balay   ierr = PetscFree(rowhit);CHKERRQ(ierr);
603606d414cSSatish Balay   ierr = PetscFree(columnsforrow);CHKERRQ(ierr);
6043acb8795SBarry Smith   ierr = MatRestoreColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr);
6053acb8795SBarry Smith   ierr = MatRestoreColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr);
606afcb2eb5SJed Brown   if (map) {ierr = ISLocalToGlobalMappingRestoreIndices(map,&ltog);CHKERRQ(ierr);}
6073a40ed3dSBarry Smith   PetscFunctionReturn(0);
608a64fbb32SBarry Smith }
609a64fbb32SBarry Smith 
610b9617806SBarry Smith 
611b9617806SBarry Smith 
612b9617806SBarry Smith 
613b9617806SBarry Smith 
614b9617806SBarry Smith 
615