xref: /petsc/src/mat/impls/aij/mpi/fdmpiaij.c (revision 0d1c53f1c6b58c59919d2989d565a8bf51fcaee8)
1a64fbb32SBarry Smith 
2c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h>
3*0d1c53f1SHong Zhang #include <../src/mat/impls/baij/mpi/mpibaij.h>
4*0d1c53f1SHong Zhang 
5*0d1c53f1SHong Zhang #undef __FUNCT__
6*0d1c53f1SHong Zhang #define __FUNCT__ "MatFDColoringApply_BAIJ_new"
7*0d1c53f1SHong Zhang PetscErrorCode  MatFDColoringApply_BAIJ_new(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
8*0d1c53f1SHong Zhang {
9*0d1c53f1SHong Zhang   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void*))coloring->f;
10*0d1c53f1SHong Zhang   PetscErrorCode ierr;
11*0d1c53f1SHong Zhang   PetscInt       k,cstart,cend,l,row,col,nz,spidx,i,j;
12*0d1c53f1SHong Zhang   PetscScalar    dx=0.0,*y,*xx,*w3_array,*dy_i,*dy=coloring->dy;
13*0d1c53f1SHong Zhang   PetscScalar    *vscale_array;
14*0d1c53f1SHong Zhang   PetscReal      epsilon=coloring->error_rel,umin=coloring->umin,unorm;
15*0d1c53f1SHong Zhang   Vec            w1=coloring->w1,w2=coloring->w2,w3,vscale=coloring->vscale;
16*0d1c53f1SHong Zhang   void           *fctx=coloring->fctx;
17*0d1c53f1SHong Zhang   PetscBool      flg=PETSC_FALSE;
18*0d1c53f1SHong Zhang   PetscInt       ctype=coloring->ctype,nxloc,nrows_k;
19*0d1c53f1SHong Zhang   Mat_MPIBAIJ    *aij=(Mat_MPIBAIJ*)J->data;
20*0d1c53f1SHong Zhang   PetscScalar    *valaddr;
21*0d1c53f1SHong Zhang   MatEntry       *Jentry=coloring->matentry;
22*0d1c53f1SHong Zhang   const PetscInt ncolors=coloring->ncolors,*ncolumns=coloring->ncolumns,*nrows=coloring->nrows;
23*0d1c53f1SHong Zhang   PetscInt       bs=J->rmap->bs;
24*0d1c53f1SHong Zhang 
25*0d1c53f1SHong Zhang   PetscFunctionBegin;
26*0d1c53f1SHong Zhang   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
27*0d1c53f1SHong Zhang   ierr = PetscOptionsGetBool(NULL,"-mat_fd_coloring_dont_rezero",&flg,NULL);CHKERRQ(ierr);
28*0d1c53f1SHong Zhang   if (flg) {
29*0d1c53f1SHong Zhang     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
30*0d1c53f1SHong Zhang   } else {
31*0d1c53f1SHong Zhang     PetscBool assembled;
32*0d1c53f1SHong Zhang     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
33*0d1c53f1SHong Zhang     if (assembled) {
34*0d1c53f1SHong Zhang       ierr = MatZeroEntries(J);CHKERRQ(ierr);
35*0d1c53f1SHong Zhang     }
36*0d1c53f1SHong Zhang   }
37*0d1c53f1SHong Zhang 
38*0d1c53f1SHong Zhang   /* create vscale for storing dx */
39*0d1c53f1SHong Zhang   if (!vscale) {
40*0d1c53f1SHong Zhang     if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') {
41*0d1c53f1SHong Zhang       ierr = VecCreateGhost(PetscObjectComm((PetscObject)J),J->cmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&vscale);CHKERRQ(ierr);
42*0d1c53f1SHong Zhang 
43*0d1c53f1SHong Zhang     } else if (ctype == IS_COLORING_GHOSTED) {
44*0d1c53f1SHong Zhang       ierr = VecDuplicate(x1,&vscale);CHKERRQ(ierr);
45*0d1c53f1SHong Zhang     }
46*0d1c53f1SHong Zhang     coloring->vscale = vscale;
47*0d1c53f1SHong Zhang   }
48*0d1c53f1SHong Zhang 
49*0d1c53f1SHong Zhang   /* (1) Set w1 = F(x1) */
50*0d1c53f1SHong Zhang   if (!coloring->fset) {
51*0d1c53f1SHong Zhang     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
52*0d1c53f1SHong Zhang     ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
53*0d1c53f1SHong Zhang     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
54*0d1c53f1SHong Zhang   } else {
55*0d1c53f1SHong Zhang     coloring->fset = PETSC_FALSE;
56*0d1c53f1SHong Zhang   }
57*0d1c53f1SHong Zhang 
58*0d1c53f1SHong Zhang   /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */
59*0d1c53f1SHong Zhang   ierr = VecGetLocalSize(x1,&nxloc);CHKERRQ(ierr);
60*0d1c53f1SHong Zhang   //PetscMPIInt rank;
61*0d1c53f1SHong Zhang   //ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
62*0d1c53f1SHong Zhang   //printf("[%d] nxloc %d\n",rank,nxloc);
63*0d1c53f1SHong Zhang   if (coloring->htype[0] == 'w') {
64*0d1c53f1SHong Zhang     /* vscale = dx is a constant scalar */
65*0d1c53f1SHong Zhang     ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr);
66*0d1c53f1SHong Zhang     dx = 1.0/(PetscSqrtReal(1.0 + unorm)*epsilon);
67*0d1c53f1SHong Zhang   } else {
68*0d1c53f1SHong Zhang     ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);
69*0d1c53f1SHong Zhang     ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr);
70*0d1c53f1SHong Zhang     for (col=0; col<nxloc; col++) {
71*0d1c53f1SHong Zhang       dx = xx[col];
72*0d1c53f1SHong Zhang       if (PetscAbsScalar(dx) < umin) {
73*0d1c53f1SHong Zhang         if (PetscRealPart(dx) >= 0.0)      dx = umin;
74*0d1c53f1SHong Zhang         else if (PetscRealPart(dx) < 0.0 ) dx = -umin;
75*0d1c53f1SHong Zhang       }
76*0d1c53f1SHong Zhang       dx               *= epsilon;
77*0d1c53f1SHong Zhang       vscale_array[col] = 1.0/dx;
78*0d1c53f1SHong Zhang     }
79*0d1c53f1SHong Zhang     ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
80*0d1c53f1SHong Zhang     ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr);
81*0d1c53f1SHong Zhang   }
82*0d1c53f1SHong Zhang   if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] != 'w') {
83*0d1c53f1SHong Zhang     ierr = VecGhostUpdateBegin(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
84*0d1c53f1SHong Zhang     ierr = VecGhostUpdateEnd(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
85*0d1c53f1SHong Zhang   }
86*0d1c53f1SHong Zhang 
87*0d1c53f1SHong Zhang   /* (3) Loop over each color */
88*0d1c53f1SHong Zhang   if (!coloring->w3) {
89*0d1c53f1SHong Zhang     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
90*0d1c53f1SHong Zhang     ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr);
91*0d1c53f1SHong Zhang   }
92*0d1c53f1SHong Zhang   w3 = coloring->w3;
93*0d1c53f1SHong Zhang 
94*0d1c53f1SHong Zhang   ierr = VecGetOwnershipRange(x1,&cstart,&cend);CHKERRQ(ierr); /* used by ghosted vscale */
95*0d1c53f1SHong Zhang   if (vscale) {
96*0d1c53f1SHong Zhang     ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr);
97*0d1c53f1SHong Zhang   }
98*0d1c53f1SHong Zhang   nz   = 0;
99*0d1c53f1SHong Zhang   for (k=0; k<ncolors; k++) {
100*0d1c53f1SHong Zhang     coloring->currentcolor = k;
101*0d1c53f1SHong Zhang 
102*0d1c53f1SHong Zhang     /*
103*0d1c53f1SHong Zhang       (3-1) Loop over each column associated with color
104*0d1c53f1SHong Zhang       adding the perturbation to the vector w3 = x1 + dx.
105*0d1c53f1SHong Zhang     */
106*0d1c53f1SHong Zhang     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
107*0d1c53f1SHong Zhang     dy_i = dy;
108*0d1c53f1SHong Zhang     for (i=0; i<bs; i++) {     /* Loop over a block of columns */
109*0d1c53f1SHong Zhang       //------------------------------------------------
110*0d1c53f1SHong Zhang       ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
111*0d1c53f1SHong Zhang       if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */
112*0d1c53f1SHong Zhang       if (coloring->htype[0] == 'w') {
113*0d1c53f1SHong Zhang         for (l=0; l<ncolumns[k]; l++) {
114*0d1c53f1SHong Zhang           col            = i + bs*coloring->columns[k][l];  /* local column (in global index!) of the matrix we are probing for */
115*0d1c53f1SHong Zhang           w3_array[col] += 1.0/dx;
116*0d1c53f1SHong Zhang 
117*0d1c53f1SHong Zhang           if (i) {
118*0d1c53f1SHong Zhang             w3_array[col-1] -= 1.0/dx; /* resume original w3[col-1] */
119*0d1c53f1SHong Zhang           }
120*0d1c53f1SHong Zhang 
121*0d1c53f1SHong Zhang         }
122*0d1c53f1SHong Zhang       } else { /* htype == 'ds' */
123*0d1c53f1SHong Zhang         vscale_array -= cstart; /* shift pointer so global index can be used */
124*0d1c53f1SHong Zhang         for (l=0; l<ncolumns[k]; l++) {
125*0d1c53f1SHong Zhang           col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
126*0d1c53f1SHong Zhang           w3_array[col] += 1.0/vscale_array[col];
127*0d1c53f1SHong Zhang 
128*0d1c53f1SHong Zhang           if (i) {
129*0d1c53f1SHong Zhang             w3_array[col-1] -=  1.0/vscale_array[col-1]; /* resume original w3[col-1] */
130*0d1c53f1SHong Zhang           }
131*0d1c53f1SHong Zhang 
132*0d1c53f1SHong Zhang         }
133*0d1c53f1SHong Zhang         vscale_array += cstart;
134*0d1c53f1SHong Zhang       }
135*0d1c53f1SHong Zhang       if (ctype == IS_COLORING_GLOBAL) w3_array += cstart;
136*0d1c53f1SHong Zhang       ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
137*0d1c53f1SHong Zhang 
138*0d1c53f1SHong Zhang       /*
139*0d1c53f1SHong Zhang        (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
140*0d1c53f1SHong Zhang                            w2 = F(x1 + dx) - F(x1)
141*0d1c53f1SHong Zhang        */
142*0d1c53f1SHong Zhang       ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
143*0d1c53f1SHong Zhang       ierr = VecPlaceArray(w2,dy_i);CHKERRQ(ierr); /* place w2 to the array dy_i */
144*0d1c53f1SHong Zhang       ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
145*0d1c53f1SHong Zhang       ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
146*0d1c53f1SHong Zhang       ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
147*0d1c53f1SHong Zhang       //---------------------------------------------
148*0d1c53f1SHong Zhang       ierr = VecResetArray(w2);CHKERRQ(ierr);
149*0d1c53f1SHong Zhang       dy_i += nxloc; /* points to dy+i*nxloc */
150*0d1c53f1SHong Zhang     }
151*0d1c53f1SHong Zhang 
152*0d1c53f1SHong Zhang     /*
153*0d1c53f1SHong Zhang      (3-3) Loop over rows of vector, putting results into Jacobian matrix
154*0d1c53f1SHong Zhang     */
155*0d1c53f1SHong Zhang     nrows_k = nrows[k];
156*0d1c53f1SHong Zhang     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
157*0d1c53f1SHong Zhang     if (coloring->htype[0] == 'w') {
158*0d1c53f1SHong Zhang       for (l=0; l<nrows_k; l++) {
159*0d1c53f1SHong Zhang         row     = bs*Jentry[nz].row;   /* local row index */
160*0d1c53f1SHong Zhang         valaddr = Jentry[nz].valaddr;
161*0d1c53f1SHong Zhang         nz++;
162*0d1c53f1SHong Zhang         spidx = 0;
163*0d1c53f1SHong Zhang         dy_i  = dy;
164*0d1c53f1SHong Zhang         for (i=0; i<bs; i++) {   /* column of the block */
165*0d1c53f1SHong Zhang           for (j=0; j<bs; j++) { /* row of the block */
166*0d1c53f1SHong Zhang             valaddr[spidx++] = dy_i[row+j]*dx;
167*0d1c53f1SHong Zhang           }
168*0d1c53f1SHong Zhang           dy_i += nxloc; /* points to dy+i*N */
169*0d1c53f1SHong Zhang         }
170*0d1c53f1SHong Zhang       }
171*0d1c53f1SHong Zhang     } else { /* htype == 'ds' */
172*0d1c53f1SHong Zhang       for (l=0; l<nrows_k; l++) {
173*0d1c53f1SHong Zhang         row                   = bs*Jentry[nz].row;   /* local row index */
174*0d1c53f1SHong Zhang         col                   = bs*Jentry[nz].col;   /* local column index */
175*0d1c53f1SHong Zhang         //*(Jentry[nz].valaddr) = y[row]*vscale_array[col];
176*0d1c53f1SHong Zhang         nz++;
177*0d1c53f1SHong Zhang       }
178*0d1c53f1SHong Zhang     }
179*0d1c53f1SHong Zhang     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
180*0d1c53f1SHong Zhang   }
181*0d1c53f1SHong Zhang   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
182*0d1c53f1SHong Zhang   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
183*0d1c53f1SHong Zhang   if (vscale) {
184*0d1c53f1SHong Zhang     ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr);
185*0d1c53f1SHong Zhang   }
186*0d1c53f1SHong Zhang 
187*0d1c53f1SHong Zhang   coloring->currentcolor = -1;
188*0d1c53f1SHong Zhang   PetscFunctionReturn(0);
189*0d1c53f1SHong Zhang }
190a64fbb32SBarry Smith 
191ab9863d7SBarry Smith extern PetscErrorCode MatCreateColmap_MPIAIJ_Private(Mat);
1924a2ae208SSatish Balay #undef __FUNCT__
193c53567a0SHong Zhang #define __FUNCT__ "MatFDColoringApply_AIJ_new"
194c53567a0SHong Zhang PetscErrorCode  MatFDColoringApply_AIJ_new(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
195fcd7ac73SHong Zhang {
196fcd7ac73SHong Zhang   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void*))coloring->f;
197fcd7ac73SHong Zhang   PetscErrorCode ierr;
198a2f2d239SHong Zhang   PetscInt       k,cstart,cend,l,row,col,nz;
1999e917edbSHong Zhang   PetscScalar    dx=0.0,*y,*xx,*w3_array;
200fcd7ac73SHong Zhang   PetscScalar    *vscale_array;
201fcd7ac73SHong Zhang   PetscReal      epsilon=coloring->error_rel,umin=coloring->umin,unorm;
2028bc97078SHong Zhang   Vec            w1=coloring->w1,w2=coloring->w2,w3,vscale=coloring->vscale;
203fcd7ac73SHong Zhang   void           *fctx=coloring->fctx;
204fcd7ac73SHong Zhang   PetscBool      flg=PETSC_FALSE;
205c53567a0SHong Zhang   PetscInt       ctype=coloring->ctype,nxloc,nrows_k;
206f5aae955SHong Zhang   Mat_MPIAIJ     *aij=(Mat_MPIAIJ*)J->data;
207573f477fSHong Zhang   MatEntry       *Jentry=coloring->matentry;
2088bc97078SHong Zhang   const PetscInt ncolors=coloring->ncolors,*ncolumns=coloring->ncolumns,*nrows=coloring->nrows;
209fcd7ac73SHong Zhang 
210fcd7ac73SHong Zhang   PetscFunctionBegin;
211fcd7ac73SHong Zhang   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
212fcd7ac73SHong Zhang   ierr = PetscOptionsGetBool(NULL,"-mat_fd_coloring_dont_rezero",&flg,NULL);CHKERRQ(ierr);
213fcd7ac73SHong Zhang   if (flg) {
214fcd7ac73SHong Zhang     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
215fcd7ac73SHong Zhang   } else {
216fcd7ac73SHong Zhang     PetscBool assembled;
217fcd7ac73SHong Zhang     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
218fcd7ac73SHong Zhang     if (assembled) {
219fcd7ac73SHong Zhang       ierr = MatZeroEntries(J);CHKERRQ(ierr);
220fcd7ac73SHong Zhang     }
221fcd7ac73SHong Zhang   }
222fcd7ac73SHong Zhang 
2239e917edbSHong Zhang   /* create vscale for storing dx */
2249e917edbSHong Zhang   if (!vscale) {
2259e917edbSHong Zhang     if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') {
2269e917edbSHong Zhang       ierr = VecCreateGhost(PetscObjectComm((PetscObject)J),J->cmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&vscale);CHKERRQ(ierr);
2279e917edbSHong Zhang     } else if (ctype == IS_COLORING_GHOSTED) {
2288bc97078SHong Zhang       ierr = VecDuplicate(x1,&vscale);CHKERRQ(ierr);
2299e917edbSHong Zhang     }
2308bc97078SHong Zhang     coloring->vscale = vscale;
231fcd7ac73SHong Zhang   }
232fcd7ac73SHong Zhang 
2338bc97078SHong Zhang   /* (1) Set w1 = F(x1) */
234fcd7ac73SHong Zhang   if (!coloring->fset) {
235fcd7ac73SHong Zhang     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
236f6af9589SHong Zhang     ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
237fcd7ac73SHong Zhang     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
238fcd7ac73SHong Zhang   } else {
239fcd7ac73SHong Zhang     coloring->fset = PETSC_FALSE;
240fcd7ac73SHong Zhang   }
241fcd7ac73SHong Zhang 
2428bc97078SHong Zhang   /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */
243f6af9589SHong Zhang   if (coloring->htype[0] == 'w') {
2449e917edbSHong Zhang     /* vscale = dx is a constant scalar */
245f6af9589SHong Zhang     ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr);
246c53567a0SHong Zhang     dx = 1.0/(PetscSqrtReal(1.0 + unorm)*epsilon);
24770e7395fSHong Zhang   } else {
24874d3cef9SHong Zhang     ierr = VecGetLocalSize(x1,&nxloc);CHKERRQ(ierr);
249f6af9589SHong Zhang     ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);
2508bc97078SHong Zhang     ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr);
25174d3cef9SHong Zhang     for (col=0; col<nxloc; col++) {
252fcd7ac73SHong Zhang       dx = xx[col];
25374d3cef9SHong Zhang       if (PetscAbsScalar(dx) < umin) {
25474d3cef9SHong Zhang         if (PetscRealPart(dx) >= 0.0)      dx = umin;
25574d3cef9SHong Zhang         else if (PetscRealPart(dx) < 0.0 ) dx = -umin;
25674d3cef9SHong Zhang       }
257fcd7ac73SHong Zhang       dx               *= epsilon;
25874d3cef9SHong Zhang       vscale_array[col] = 1.0/dx;
259f6af9589SHong Zhang     }
26074d3cef9SHong Zhang     ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
2618bc97078SHong Zhang     ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr);
26270e7395fSHong Zhang   }
2639e917edbSHong Zhang   if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] != 'w') {
2648bc97078SHong Zhang     ierr = VecGhostUpdateBegin(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2658bc97078SHong Zhang     ierr = VecGhostUpdateEnd(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
266fcd7ac73SHong Zhang   }
267fcd7ac73SHong Zhang 
2688bc97078SHong Zhang   /* (3) Loop over each color */
2698bc97078SHong Zhang   if (!coloring->w3) {
2708bc97078SHong Zhang     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
2718bc97078SHong Zhang     ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr);
2728bc97078SHong Zhang   }
2738bc97078SHong Zhang   w3 = coloring->w3;
274fcd7ac73SHong Zhang 
2758bc97078SHong Zhang   ierr = VecGetOwnershipRange(x1,&cstart,&cend);CHKERRQ(ierr); /* used by ghosted vscale */
2769e917edbSHong Zhang   if (vscale) {
2778bc97078SHong Zhang     ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr);
2789e917edbSHong Zhang   }
2798bc97078SHong Zhang   nz   = 0;
2808bc97078SHong Zhang   for (k=0; k<ncolors; k++) {
2818bc97078SHong Zhang     coloring->currentcolor = k;
2829e917edbSHong Zhang 
283fcd7ac73SHong Zhang     /*
2848bc97078SHong Zhang       (3-1) Loop over each column associated with color
285f6af9589SHong Zhang       adding the perturbation to the vector w3 = x1 + dx.
286fcd7ac73SHong Zhang     */
287f6af9589SHong Zhang     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
288f6af9589SHong Zhang     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
289a2f2d239SHong Zhang     if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */
290b039d6c4SHong Zhang     if (coloring->htype[0] == 'w') {
2918bc97078SHong Zhang       for (l=0; l<ncolumns[k]; l++) {
292f6af9589SHong Zhang         col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
2939e917edbSHong Zhang         w3_array[col] += 1.0/dx;
2949e917edbSHong Zhang       }
295b039d6c4SHong Zhang     } else { /* htype == 'ds' */
296a2f2d239SHong Zhang       vscale_array -= cstart; /* shift pointer so global index can be used */
297b039d6c4SHong Zhang       for (l=0; l<ncolumns[k]; l++) {
298b039d6c4SHong Zhang         col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
299a2f2d239SHong Zhang         w3_array[col] += 1.0/vscale_array[col];
30070e7395fSHong Zhang       }
301a2f2d239SHong Zhang       vscale_array += cstart;
302fcd7ac73SHong Zhang     }
303a2f2d239SHong Zhang     if (ctype == IS_COLORING_GLOBAL) w3_array += cstart;
304fcd7ac73SHong Zhang     ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
3059e917edbSHong Zhang 
306fcd7ac73SHong Zhang     /*
3078bc97078SHong Zhang       (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
308fcd7ac73SHong Zhang                            w2 = F(x1 + dx) - F(x1)
309fcd7ac73SHong Zhang     */
310fcd7ac73SHong Zhang     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
311fcd7ac73SHong Zhang     ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
312fcd7ac73SHong Zhang     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
313fcd7ac73SHong Zhang     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
3149e917edbSHong Zhang 
3158bc97078SHong Zhang     /*
3168bc97078SHong Zhang      (3-3) Loop over rows of vector, putting results into Jacobian matrix
3178bc97078SHong Zhang     */
318c53567a0SHong Zhang     nrows_k = nrows[k];
319fcd7ac73SHong Zhang     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
320b039d6c4SHong Zhang     if (coloring->htype[0] == 'w') {
321c53567a0SHong Zhang       for (l=0; l<nrows_k; l++) {
322573f477fSHong Zhang         row                     = Jentry[nz].row;   /* local row index */
323b039d6c4SHong Zhang         *(Jentry[nz++].valaddr) = y[row]*dx;
3249e917edbSHong Zhang       }
325b039d6c4SHong Zhang     } else { /* htype == 'ds' */
326c53567a0SHong Zhang       for (l=0; l<nrows_k; l++) {
327b039d6c4SHong Zhang         row                   = Jentry[nz].row;   /* local row index */
328b039d6c4SHong Zhang         *(Jentry[nz].valaddr) = y[row]*vscale_array[Jentry[nz].col];
329573f477fSHong Zhang         nz++;
330fcd7ac73SHong Zhang       }
331b039d6c4SHong Zhang     }
332fcd7ac73SHong Zhang     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
3338bc97078SHong Zhang   }
334f5aae955SHong Zhang   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
335f5aae955SHong Zhang   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3369e917edbSHong Zhang   if (vscale) {
3378bc97078SHong Zhang     ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr);
3389e917edbSHong Zhang   }
339fcd7ac73SHong Zhang 
340fcd7ac73SHong Zhang   coloring->currentcolor = -1;
341fcd7ac73SHong Zhang   PetscFunctionReturn(0);
342fcd7ac73SHong Zhang }
343fcd7ac73SHong Zhang 
344fcd7ac73SHong Zhang #undef __FUNCT__
345fcd7ac73SHong Zhang #define __FUNCT__ "MatFDColoringCreate_MPIAIJ_new"
346fcd7ac73SHong Zhang PetscErrorCode MatFDColoringCreate_MPIAIJ_new(Mat mat,ISColoring iscoloring,MatFDColoring c)
347fcd7ac73SHong Zhang {
348*0d1c53f1SHong Zhang   //Mat_MPIAIJ             *aij=(Mat_MPIAIJ*)mat->data;
349fcd7ac73SHong Zhang   PetscErrorCode         ierr;
350fcd7ac73SHong Zhang   PetscMPIInt            size,*ncolsonproc,*disp,nn;
3514b2e90caSHong Zhang   PetscInt               i,n,nrows,nrows_i,j,k,m,ncols,col;
35272c15787SHong Zhang   const PetscInt         *is,*A_ci,*A_cj,*B_ci,*B_cj,*row = NULL,*ltog=NULL;
353fcd7ac73SHong Zhang   PetscInt               nis=iscoloring->n,nctot,*cols;
354d3825b63SHong Zhang   PetscInt               *rowhit,cstart,cend,colb;
355fcd7ac73SHong Zhang   IS                     *isa;
356fcd7ac73SHong Zhang   ISLocalToGlobalMapping map=mat->cmap->mapping;
357fcd7ac73SHong Zhang   PetscInt               ctype=c->ctype;
358*0d1c53f1SHong Zhang   Mat                    A,B;
359*0d1c53f1SHong Zhang   //Mat                    A=aij->A,B=aij->B;
360*0d1c53f1SHong Zhang   //Mat_SeqAIJ             *spA=(Mat_SeqAIJ*)A->data,*spB=(Mat_SeqAIJ*)B->data;
361*0d1c53f1SHong Zhang   //PetscScalar            *A_val=spA->a,*B_val=spB->a;
362*0d1c53f1SHong Zhang   PetscScalar            *A_val,*B_val;
363fcd7ac73SHong Zhang   PetscInt               spidx;
364*0d1c53f1SHong Zhang   PetscInt               *spidxA,*spidxB,nz,bs,bs2;
365d3825b63SHong Zhang   PetscScalar            **valaddrhit;
366573f477fSHong Zhang   MatEntry               *Jentry;
367*0d1c53f1SHong Zhang   PetscBool              isBAIJ;
368*0d1c53f1SHong Zhang #if defined(PETSC_USE_CTABLE)
369*0d1c53f1SHong Zhang   PetscTable             colmap=NULL;
370*0d1c53f1SHong Zhang #else
371*0d1c53f1SHong Zhang   PetscInt               *colmap=NULL;     /* local col number of off-diag col */
372*0d1c53f1SHong Zhang #endif
373fcd7ac73SHong Zhang 
374fcd7ac73SHong Zhang   PetscFunctionBegin;
37572c15787SHong Zhang   if (ctype == IS_COLORING_GHOSTED) {
37672c15787SHong 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");
37772c15787SHong Zhang     ierr = ISLocalToGlobalMappingGetIndices(map,&ltog);CHKERRQ(ierr);
37872c15787SHong Zhang   }
379fcd7ac73SHong Zhang 
380*0d1c53f1SHong Zhang   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
381*0d1c53f1SHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr);
382*0d1c53f1SHong Zhang   if (!isBAIJ) {
383*0d1c53f1SHong Zhang     bs = 1; /* only bs=1 is supported for non MPIBAIJ matrix */
384*0d1c53f1SHong Zhang     Mat_MPIAIJ             *aij=(Mat_MPIAIJ*)mat->data;
385*0d1c53f1SHong Zhang     A=aij->A; B=aij->B;
386*0d1c53f1SHong Zhang     Mat_SeqAIJ             *spA=(Mat_SeqAIJ*)A->data,*spB=(Mat_SeqAIJ*)B->data;
387*0d1c53f1SHong Zhang     A_val=spA->a; B_val=spB->a;
388*0d1c53f1SHong Zhang     nz         = spA->nz + spB->nz; /* total nonzero entries of mat */
389*0d1c53f1SHong Zhang     if (!aij->colmap) {
390*0d1c53f1SHong Zhang       ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
391*0d1c53f1SHong Zhang       colmap = aij->colmap;
392*0d1c53f1SHong Zhang     }
393*0d1c53f1SHong Zhang     ierr = MatGetColumnIJ_SeqAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
394*0d1c53f1SHong Zhang     ierr = MatGetColumnIJ_SeqAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
395*0d1c53f1SHong Zhang   } else {
396*0d1c53f1SHong Zhang     Mat_MPIBAIJ             *aij=(Mat_MPIBAIJ*)mat->data;
397*0d1c53f1SHong Zhang     A=aij->A; B=aij->B;
398*0d1c53f1SHong Zhang     Mat_SeqBAIJ             *spA=(Mat_SeqBAIJ*)A->data,*spB=(Mat_SeqBAIJ*)B->data;
399*0d1c53f1SHong Zhang     A_val=spA->a; B_val=spB->a;
400*0d1c53f1SHong Zhang     nz         = spA->nz + spB->nz; /* total nonzero entries of mat */
401*0d1c53f1SHong Zhang     if (!aij->colmap) {
402*0d1c53f1SHong Zhang       ierr = MatCreateColmap_MPIBAIJ_Private(mat);CHKERRQ(ierr);
403*0d1c53f1SHong Zhang       colmap = aij->colmap;
404*0d1c53f1SHong Zhang     }
405*0d1c53f1SHong Zhang     ierr = MatGetColumnIJ_SeqBAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
406*0d1c53f1SHong Zhang     ierr = MatGetColumnIJ_SeqBAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
407*0d1c53f1SHong Zhang   }
408*0d1c53f1SHong Zhang 
409*0d1c53f1SHong Zhang   m         = mat->rmap->n/bs;
410*0d1c53f1SHong Zhang   cstart    = mat->cmap->rstart/bs;
411*0d1c53f1SHong Zhang   cend      = mat->cmap->rend/bs;
412*0d1c53f1SHong Zhang   c->M      = mat->rmap->N/bs;         /* set the global rows and columns and local rows */
413*0d1c53f1SHong Zhang   c->N      = mat->cmap->N/bs;
414d3825b63SHong Zhang   c->m      = m;
415*0d1c53f1SHong Zhang   c->rstart = mat->rmap->rstart/bs;
416fcd7ac73SHong Zhang 
417fcd7ac73SHong Zhang   c->ncolors = nis;
418fcd7ac73SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
419fcd7ac73SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr);
420fcd7ac73SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
4218bc97078SHong Zhang   ierr       = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr);
422fcd7ac73SHong Zhang 
423*0d1c53f1SHong Zhang   //nz         = spA->nz + spB->nz; /* total nonzero entries of mat */
424573f477fSHong Zhang   ierr       = PetscMalloc(nz*sizeof(MatEntry),&Jentry);CHKERRQ(ierr);
4258bc97078SHong Zhang   ierr       = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr);
4268bc97078SHong Zhang   c->matentry = Jentry;
427a774a6f1SHong Zhang 
428d3825b63SHong Zhang   /* Allow access to data structures of local part of matrix
429d3825b63SHong Zhang    - creates aij->colmap which maps global column number to local number in part B */
430*0d1c53f1SHong Zhang   //if (!aij->colmap) {
431*0d1c53f1SHong Zhang   //  ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
432*0d1c53f1SHong Zhang   //}
433*0d1c53f1SHong Zhang   //ierr = MatGetColumnIJ_SeqAIJ_Color(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
434*0d1c53f1SHong Zhang   //ierr = MatGetColumnIJ_SeqAIJ_Color(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
435fcd7ac73SHong Zhang 
436d3825b63SHong Zhang   ierr = PetscMalloc2(m+1,PetscInt,&rowhit,m+1,PetscScalar*,&valaddrhit);CHKERRQ(ierr);
437fcd7ac73SHong Zhang   nz = 0;
43872c15787SHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
439d3825b63SHong Zhang   for (i=0; i<nis; i++) { /* for each local color */
440fcd7ac73SHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
441fcd7ac73SHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
442fcd7ac73SHong Zhang 
443fcd7ac73SHong Zhang     c->ncolumns[i] = n; /* local number of columns of this color on this process */
444fcd7ac73SHong Zhang     if (n) {
445fcd7ac73SHong Zhang       ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr);
446fcd7ac73SHong Zhang       ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr);
447fcd7ac73SHong Zhang       ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr);
448a2f2d239SHong Zhang       /* convert global column indices to local ones! */
449a2f2d239SHong Zhang 
450fcd7ac73SHong Zhang     } else {
451fcd7ac73SHong Zhang       c->columns[i] = 0;
452fcd7ac73SHong Zhang     }
453fcd7ac73SHong Zhang 
454fcd7ac73SHong Zhang     if (ctype == IS_COLORING_GLOBAL) {
455d3825b63SHong Zhang       /* Determine nctot, the total (parallel) number of columns of this color */
456fcd7ac73SHong Zhang       ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr);
457fcd7ac73SHong Zhang       ierr = PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);CHKERRQ(ierr);
458fcd7ac73SHong Zhang 
459d3825b63SHong Zhang       /* ncolsonproc[j]: local ncolumns on proc[j] of this color */
460fcd7ac73SHong Zhang       ierr  = PetscMPIIntCast(n,&nn);CHKERRQ(ierr);
461fcd7ac73SHong Zhang       ierr  = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
462fcd7ac73SHong Zhang       nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j];
463fcd7ac73SHong Zhang       if (!nctot) {
464fcd7ac73SHong Zhang         ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr);
465fcd7ac73SHong Zhang       }
466fcd7ac73SHong Zhang 
467fcd7ac73SHong Zhang       disp[0] = 0;
468fcd7ac73SHong Zhang       for (j=1; j<size; j++) {
469fcd7ac73SHong Zhang         disp[j] = disp[j-1] + ncolsonproc[j-1];
470fcd7ac73SHong Zhang       }
471fcd7ac73SHong Zhang 
472d3825b63SHong Zhang       /* Get cols, the complete list of columns for this color on each process */
473fcd7ac73SHong Zhang       ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
474fcd7ac73SHong Zhang       ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
475fcd7ac73SHong Zhang       ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr);
476fcd7ac73SHong Zhang     } else if (ctype == IS_COLORING_GHOSTED) {
477fcd7ac73SHong Zhang       /* Determine local number of columns of this color on this process, including ghost points */
478fcd7ac73SHong Zhang       nctot = n;
479fcd7ac73SHong Zhang       ierr  = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
480fcd7ac73SHong Zhang       ierr  = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr);
481fcd7ac73SHong Zhang     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type");
482fcd7ac73SHong Zhang 
4831b97d346SHong Zhang     /* Mark all rows affect by these columns */
484d3825b63SHong Zhang     ierr = PetscMemzero(rowhit,m*sizeof(PetscInt));CHKERRQ(ierr);
485*0d1c53f1SHong Zhang     bs2     = bs*bs;
4864b2e90caSHong Zhang     nrows_i = 0;
4871b97d346SHong Zhang     for (j=0; j<nctot; j++) { /* loop over columns*/
488fcd7ac73SHong Zhang       if (ctype == IS_COLORING_GHOSTED) {
489fcd7ac73SHong Zhang         col = ltog[cols[j]];
490fcd7ac73SHong Zhang       } else {
491fcd7ac73SHong Zhang         col = cols[j];
492fcd7ac73SHong Zhang       }
493fcd7ac73SHong Zhang       if (col >= cstart && col < cend) { /* column is in diagonal block of matrix A */
494d3825b63SHong Zhang         row      = A_cj + A_ci[col-cstart];
495d3825b63SHong Zhang         nrows    = A_ci[col-cstart+1] - A_ci[col-cstart];
4964b2e90caSHong Zhang         nrows_i += nrows;
497fcd7ac73SHong Zhang         /* loop over columns of A marking them in rowhit */
498d3825b63SHong Zhang         for (k=0; k<nrows; k++) {
499d3825b63SHong Zhang           /* set valaddrhit for part A */
500fcd7ac73SHong Zhang           spidx            = spidxA[A_ci[col-cstart] + k];
501*0d1c53f1SHong Zhang           valaddrhit[*row] = &A_val[bs2*spidx];
502a774a6f1SHong Zhang           rowhit[*row++]   = col - cstart + 1; /* local column index */
503fcd7ac73SHong Zhang         }
504fcd7ac73SHong Zhang       } else { /* column is in off-diagonal block of matrix B */
505fcd7ac73SHong Zhang #if defined(PETSC_USE_CTABLE)
506*0d1c53f1SHong Zhang         //ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
507*0d1c53f1SHong Zhang         ierr = PetscTableFind(colmap,col+1,&colb);CHKERRQ(ierr);
508fcd7ac73SHong Zhang         colb--;
509fcd7ac73SHong Zhang #else
510*0d1c53f1SHong Zhang         //colb = aij->colmap[col] - 1; /* local column index */
511*0d1c53f1SHong Zhang         colb = colmap[col] - 1; /* local column index */
512fcd7ac73SHong Zhang #endif
513fcd7ac73SHong Zhang         if (colb == -1) {
514d3825b63SHong Zhang           nrows = 0;
515fcd7ac73SHong Zhang         } else {
516*0d1c53f1SHong Zhang           colb  = colb/bs;
517d3825b63SHong Zhang           row   = B_cj + B_ci[colb];
518d3825b63SHong Zhang           nrows = B_ci[colb+1] - B_ci[colb];
519fcd7ac73SHong Zhang         }
5204b2e90caSHong Zhang         nrows_i += nrows;
521fcd7ac73SHong Zhang         /* loop over columns of B marking them in rowhit */
522d3825b63SHong Zhang         for (k=0; k<nrows; k++) {
523d3825b63SHong Zhang           /* set valaddrhit for part B */
524fcd7ac73SHong Zhang           spidx            = spidxB[B_ci[colb] + k];
525*0d1c53f1SHong Zhang           valaddrhit[*row] = &B_val[bs2*spidx];
52670e7395fSHong Zhang           rowhit[*row++]   = colb + 1 + cend - cstart; /* local column index */
527fcd7ac73SHong Zhang         }
528fcd7ac73SHong Zhang       }
529fcd7ac73SHong Zhang     }
5304b2e90caSHong Zhang     c->nrows[i] = nrows_i;
5318bc97078SHong Zhang 
532d3825b63SHong Zhang     for (j=0; j<m; j++) {
533fcd7ac73SHong Zhang       if (rowhit[j]) {
534573f477fSHong Zhang         Jentry[nz].row     = j;              /* local row index */
535573f477fSHong Zhang         Jentry[nz].col     = rowhit[j] - 1;  /* local column index */
536573f477fSHong Zhang         Jentry[nz].valaddr = valaddrhit[j];  /* address of mat value for this entry */
537573f477fSHong Zhang         nz++;
538fcd7ac73SHong Zhang       }
539fcd7ac73SHong Zhang     }
540fcd7ac73SHong Zhang     ierr = PetscFree(cols);CHKERRQ(ierr);
541fcd7ac73SHong Zhang   }
5428bc97078SHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
543*0d1c53f1SHong 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);
544fcd7ac73SHong Zhang 
545d3825b63SHong Zhang   ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr);
546*0d1c53f1SHong Zhang   if (isBAIJ) {
547*0d1c53f1SHong Zhang     ierr = MatRestoreColumnIJ_SeqBAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
548*0d1c53f1SHong Zhang     ierr = MatRestoreColumnIJ_SeqBAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
549*0d1c53f1SHong Zhang     ierr = PetscMalloc(bs*mat->cmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr);
550*0d1c53f1SHong Zhang     mat->ops->fdcoloringapply = MatFDColoringApply_BAIJ_new;
551*0d1c53f1SHong Zhang   } else {
552*0d1c53f1SHong Zhang     ierr = MatRestoreColumnIJ_SeqAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
553*0d1c53f1SHong Zhang     ierr = MatRestoreColumnIJ_SeqAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
554*0d1c53f1SHong Zhang     mat->ops->fdcoloringapply = MatFDColoringApply_AIJ_new;
555*0d1c53f1SHong Zhang   }
556*0d1c53f1SHong Zhang 
55772c15787SHong Zhang   if (ctype == IS_COLORING_GHOSTED) {
55872c15787SHong Zhang     ierr = ISLocalToGlobalMappingRestoreIndices(map,&ltog);CHKERRQ(ierr);
55972c15787SHong Zhang   }
560fcd7ac73SHong Zhang   PetscFunctionReturn(0);
561fcd7ac73SHong Zhang }
562fcd7ac73SHong Zhang 
563fcd7ac73SHong Zhang /*------------------------------------------------------*/
564fcd7ac73SHong Zhang #undef __FUNCT__
5654a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate_MPIAIJ"
566dfbe8321SBarry Smith PetscErrorCode MatFDColoringCreate_MPIAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
567a64fbb32SBarry Smith {
5686eaac0f3SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)mat->data;
5696849ba73SBarry Smith   PetscErrorCode         ierr;
570b1d57f15SBarry Smith   PetscMPIInt            size,*ncolsonproc,*disp,nn;
5711a83f524SJed Brown   PetscInt               i,n,nrows,j,k,m,ncols,col;
572afcb2eb5SJed Brown   const PetscInt         *is,*A_ci,*A_cj,*B_ci,*B_cj,*rows = 0,*ltog;
5731a83f524SJed Brown   PetscInt               nis = iscoloring->n,nctot,*cols;
574f6d58c54SBarry Smith   PetscInt               *rowhit,M,cstart,cend,colb;
575b1d57f15SBarry Smith   PetscInt               *columnsforrow,l;
576b9617806SBarry Smith   IS                     *isa;
577ace3abfcSBarry Smith   PetscBool              done,flg;
578992144d0SBarry Smith   ISLocalToGlobalMapping map   = mat->cmap->mapping;
579afcb2eb5SJed Brown   PetscInt               ctype=c->ctype;
580fcd7ac73SHong Zhang   PetscBool              new_impl=PETSC_FALSE;
581a64fbb32SBarry Smith 
5823a40ed3dSBarry Smith   PetscFunctionBegin;
583fcd7ac73SHong Zhang   ierr = PetscOptionsName("-new","using new impls","",&new_impl);CHKERRQ(ierr);
584fcd7ac73SHong Zhang   if (new_impl){
585fcd7ac73SHong Zhang     ierr =  MatFDColoringCreate_MPIAIJ_new(mat,iscoloring,c);CHKERRQ(ierr);
586fcd7ac73SHong Zhang     PetscFunctionReturn(0);
587fcd7ac73SHong Zhang   }
588ce94432eSBarry 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");
589522c5e43SBarry Smith 
590afcb2eb5SJed Brown   if (map) {ierr = ISLocalToGlobalMappingGetIndices(map,&ltog);CHKERRQ(ierr);}
591afcb2eb5SJed Brown   else     ltog = NULL;
592b9617806SBarry Smith   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
5933acb8795SBarry Smith 
594f6d58c54SBarry Smith   M         = mat->rmap->n;
595f6d58c54SBarry Smith   cstart    = mat->cmap->rstart;
596f6d58c54SBarry Smith   cend      = mat->cmap->rend;
597f6d58c54SBarry Smith   c->M      = mat->rmap->N;         /* set the global rows and columns and local rows */
598f6d58c54SBarry Smith   c->N      = mat->cmap->N;
599f6d58c54SBarry Smith   c->m      = mat->rmap->n;
600f6d58c54SBarry Smith   c->rstart = mat->rmap->rstart;
601005c665bSBarry Smith 
602a64fbb32SBarry Smith   c->ncolors = nis;
603b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
604b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr);
605b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
606b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->rows);CHKERRQ(ierr);
607b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columnsforrow);CHKERRQ(ierr);
6083bb1ff40SBarry Smith   ierr       = PetscLogObjectMemory((PetscObject)c,5*nis*sizeof(PetscInt));CHKERRQ(ierr);
6096eaac0f3SBarry Smith 
6106eaac0f3SBarry Smith   /* Allow access to data structures of local part of matrix */
6116eaac0f3SBarry Smith   if (!aij->colmap) {
612ab9863d7SBarry Smith     ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
6136eaac0f3SBarry Smith   }
6143acb8795SBarry Smith   ierr = MatGetColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr);
6153acb8795SBarry Smith   ierr = MatGetColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr);
6166eaac0f3SBarry Smith 
617b1d57f15SBarry Smith   ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr);
618b1d57f15SBarry Smith   ierr = PetscMalloc((M+1)*sizeof(PetscInt),&columnsforrow);CHKERRQ(ierr);
6196eaac0f3SBarry Smith 
620a64fbb32SBarry Smith   for (i=0; i<nis; i++) {
621b9b97703SBarry Smith     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
622a64fbb32SBarry Smith     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
6232205254eSKarl Rupp 
624fcd7ac73SHong Zhang     c->ncolumns[i] = n; /* local number of columns of this color on this process */
625a64fbb32SBarry Smith     if (n) {
626b1d57f15SBarry Smith       ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr);
6273bb1ff40SBarry Smith       ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr);
628b1d57f15SBarry Smith       ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr);
629a64fbb32SBarry Smith     } else {
630a64fbb32SBarry Smith       c->columns[i] = 0;
631a64fbb32SBarry Smith     }
632a64fbb32SBarry Smith 
6338ee2e534SBarry Smith     if (ctype == IS_COLORING_GLOBAL) {
6346eaac0f3SBarry Smith       /* Determine the total (parallel) number of columns of this color */
635ce94432eSBarry Smith       ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr);
636687f1162SBarry Smith       ierr = PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);CHKERRQ(ierr);
637b8f8c88eSHong Zhang 
6384dc2109aSBarry Smith       ierr  = PetscMPIIntCast(n,&nn);CHKERRQ(ierr);
639ce94432eSBarry Smith       ierr  = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
6402205254eSKarl Rupp       nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j];
6413a7fca6bSBarry Smith       if (!nctot) {
642ae15b995SBarry Smith         ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr);
6433a7fca6bSBarry Smith       }
6446eaac0f3SBarry Smith 
6456eaac0f3SBarry Smith       disp[0] = 0;
6466eaac0f3SBarry Smith       for (j=1; j<size; j++) {
6476eaac0f3SBarry Smith         disp[j] = disp[j-1] + ncolsonproc[j-1];
6486eaac0f3SBarry Smith       }
6496eaac0f3SBarry Smith 
6506eaac0f3SBarry Smith       /* Get complete list of columns for color on each processor */
651b1d57f15SBarry Smith       ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
652ce94432eSBarry Smith       ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
6531d79065fSBarry Smith       ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr);
654b8f8c88eSHong Zhang     } else if (ctype == IS_COLORING_GHOSTED) {
655b8f8c88eSHong Zhang       /* Determine local number of columns of this color on this process, including ghost points */
656b8f8c88eSHong Zhang       nctot = n;
657b8f8c88eSHong Zhang       ierr  = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
658b8f8c88eSHong Zhang       ierr  = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr);
659f23aa3ddSBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type");
6606eaac0f3SBarry Smith 
6616eaac0f3SBarry Smith     /*
6626eaac0f3SBarry Smith        Mark all rows affect by these columns
6636eaac0f3SBarry Smith     */
664b8f8c88eSHong Zhang     /* Temporary option to allow for debugging/testing */
66590d69ab7SBarry Smith     flg  = PETSC_FALSE;
6660298fd71SBarry Smith     ierr = PetscOptionsGetBool(NULL,"-matfdcoloring_slow",&flg,NULL);CHKERRQ(ierr);
667f158e583SBarry Smith     if (!flg) { /*-----------------------------------------------------------------------------*/
668f158e583SBarry Smith       /* crude, fast version */
669b1d57f15SBarry Smith       ierr = PetscMemzero(rowhit,M*sizeof(PetscInt));CHKERRQ(ierr);
670a64fbb32SBarry Smith       /* loop over columns*/
6716eaac0f3SBarry Smith       for (j=0; j<nctot; j++) {
672b8f8c88eSHong Zhang         if (ctype == IS_COLORING_GHOSTED) {
673b8f8c88eSHong Zhang           col = ltog[cols[j]];
674b8f8c88eSHong Zhang         } else {
6756eaac0f3SBarry Smith           col = cols[j];
676b8f8c88eSHong Zhang         }
6776eaac0f3SBarry Smith         if (col >= cstart && col < cend) {
6786eaac0f3SBarry Smith           /* column is in diagonal block of matrix */
6796eaac0f3SBarry Smith           rows = A_cj + A_ci[col-cstart];
6806eaac0f3SBarry Smith           m    = A_ci[col-cstart+1] - A_ci[col-cstart];
6816eaac0f3SBarry Smith         } else {
682aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
683cb9801acSJed Brown           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
684fa46199cSSatish Balay           colb--;
685b3d2dc96SSatish Balay #else
6866eaac0f3SBarry Smith           colb = aij->colmap[col] - 1;
687b3d2dc96SSatish Balay #endif
6886eaac0f3SBarry Smith           if (colb == -1) {
6896eaac0f3SBarry Smith             m = 0;
6906eaac0f3SBarry Smith           } else {
6916eaac0f3SBarry Smith             rows = B_cj + B_ci[colb];
6926eaac0f3SBarry Smith             m    = B_ci[colb+1] - B_ci[colb];
6936eaac0f3SBarry Smith           }
6946eaac0f3SBarry Smith         }
695a64fbb32SBarry Smith         /* loop over columns marking them in rowhit */
696a64fbb32SBarry Smith         for (k=0; k<m; k++) {
697a64fbb32SBarry Smith           rowhit[*rows++] = col + 1;
698a64fbb32SBarry Smith         }
699a64fbb32SBarry Smith       }
7006eaac0f3SBarry Smith 
701a64fbb32SBarry Smith       /* count the number of hits */
702a64fbb32SBarry Smith       nrows = 0;
7036eaac0f3SBarry Smith       for (j=0; j<M; j++) {
704a64fbb32SBarry Smith         if (rowhit[j]) nrows++;
705a64fbb32SBarry Smith       }
706a64fbb32SBarry Smith       c->nrows[i] = nrows;
707b1d57f15SBarry Smith       ierr        = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr);
708b1d57f15SBarry Smith       ierr        = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr);
7093bb1ff40SBarry Smith       ierr        = PetscLogObjectMemory((PetscObject)c,2*(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr);
710a64fbb32SBarry Smith       nrows       = 0;
7116eaac0f3SBarry Smith       for (j=0; j<M; j++) {
712a64fbb32SBarry Smith         if (rowhit[j]) {
713fcd7ac73SHong Zhang           c->rows[i][nrows]          = j;              /* local row index */
714fcd7ac73SHong Zhang           c->columnsforrow[i][nrows] = rowhit[j] - 1;  /* global column index */
715a64fbb32SBarry Smith           nrows++;
716a64fbb32SBarry Smith         }
717a64fbb32SBarry Smith       }
718a64fbb32SBarry Smith     } else { /*-------------------------------------------------------------------------------*/
719f158e583SBarry Smith       /* slow version, using rowhit as a linked list */
720b1d57f15SBarry Smith       PetscInt currentcol,fm,mfm;
7216eaac0f3SBarry Smith       rowhit[M] = M;
722a64fbb32SBarry Smith       nrows     = 0;
723a64fbb32SBarry Smith       /* loop over columns*/
7246eaac0f3SBarry Smith       for (j=0; j<nctot; j++) {
725b8f8c88eSHong Zhang         if (ctype == IS_COLORING_GHOSTED) {
726b8f8c88eSHong Zhang           col = ltog[cols[j]];
727b8f8c88eSHong Zhang         } else {
7286eaac0f3SBarry Smith           col = cols[j];
729b8f8c88eSHong Zhang         }
7306eaac0f3SBarry Smith         if (col >= cstart && col < cend) {
7316eaac0f3SBarry Smith           /* column is in diagonal block of matrix */
7326eaac0f3SBarry Smith           rows = A_cj + A_ci[col-cstart];
7336eaac0f3SBarry Smith           m    = A_ci[col-cstart+1] - A_ci[col-cstart];
7346eaac0f3SBarry Smith         } else {
735aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
7360f5bd95cSBarry Smith           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
737fa46199cSSatish Balay           colb--;
738b3d2dc96SSatish Balay #else
7396eaac0f3SBarry Smith           colb = aij->colmap[col] - 1;
740b3d2dc96SSatish Balay #endif
7416eaac0f3SBarry Smith           if (colb == -1) {
7426eaac0f3SBarry Smith             m = 0;
7436eaac0f3SBarry Smith           } else {
7446eaac0f3SBarry Smith             rows = B_cj + B_ci[colb];
7456eaac0f3SBarry Smith             m    = B_ci[colb+1] - B_ci[colb];
7466eaac0f3SBarry Smith           }
7476eaac0f3SBarry Smith         }
748b8f8c88eSHong Zhang 
749a64fbb32SBarry Smith         /* loop over columns marking them in rowhit */
7506eaac0f3SBarry Smith         fm = M;    /* fm points to first entry in linked list */
751a64fbb32SBarry Smith         for (k=0; k<m; k++) {
752a64fbb32SBarry Smith           currentcol = *rows++;
753a64fbb32SBarry Smith           /* is it already in the list? */
754a64fbb32SBarry Smith           do {
755a64fbb32SBarry Smith             mfm = fm;
756a64fbb32SBarry Smith             fm  = rowhit[fm];
757a64fbb32SBarry Smith           } while (fm < currentcol);
758a64fbb32SBarry Smith           /* not in list so add it */
759a64fbb32SBarry Smith           if (fm != currentcol) {
760a64fbb32SBarry Smith             nrows++;
761a64fbb32SBarry Smith             columnsforrow[currentcol] = col;
762a64fbb32SBarry Smith             /* next three lines insert new entry into linked list */
763a64fbb32SBarry Smith             rowhit[mfm]        = currentcol;
764a64fbb32SBarry Smith             rowhit[currentcol] = fm;
765a64fbb32SBarry Smith             fm                 = currentcol;
766a64fbb32SBarry Smith             /* fm points to present position in list since we know the columns are sorted */
767f23aa3ddSBarry Smith           } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid coloring of matrix detected");
768a64fbb32SBarry Smith         }
769a64fbb32SBarry Smith       }
770a64fbb32SBarry Smith       c->nrows[i] = nrows;
7712205254eSKarl Rupp 
772b1d57f15SBarry Smith       ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr);
773b1d57f15SBarry Smith       ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr);
7743bb1ff40SBarry Smith       ierr = PetscLogObjectMemory((PetscObject)c,(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr);
775a64fbb32SBarry Smith       /* now store the linked list of rows into c->rows[i] */
776a64fbb32SBarry Smith       nrows = 0;
7776eaac0f3SBarry Smith       fm    = rowhit[M];
778a64fbb32SBarry Smith       do {
779a64fbb32SBarry Smith         c->rows[i][nrows]            = fm;
780a64fbb32SBarry Smith         c->columnsforrow[i][nrows++] = columnsforrow[fm];
781a64fbb32SBarry Smith         fm                           = rowhit[fm];
7826eaac0f3SBarry Smith       } while (fm < M);
7836eaac0f3SBarry Smith     } /* ---------------------------------------------------------------------------------------*/
784606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
7856eaac0f3SBarry Smith   }
78630b34957SBarry Smith 
78730b34957SBarry Smith   /* Optimize by adding the vscale, and scaleforrow[][] fields */
78830b34957SBarry Smith   /*
78930b34957SBarry Smith        vscale will contain the "diagonal" on processor scalings followed by the off processor
79030b34957SBarry Smith   */
7918ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
792ce94432eSBarry Smith     ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),aij->A->rmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&c->vscale);CHKERRQ(ierr);
793b1d57f15SBarry Smith     ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr);
79430b34957SBarry Smith     for (k=0; k<c->ncolors; k++) {
795b1d57f15SBarry Smith       ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr);
79630b34957SBarry Smith       for (l=0; l<c->nrows[k]; l++) {
79730b34957SBarry Smith         col = c->columnsforrow[k][l];
79830b34957SBarry Smith         if (col >= cstart && col < cend) {
79930b34957SBarry Smith           /* column is in diagonal block of matrix */
80030b34957SBarry Smith           colb = col - cstart;
80130b34957SBarry Smith         } else {
80230b34957SBarry Smith           /* column  is in "off-processor" part */
80330b34957SBarry Smith #if defined(PETSC_USE_CTABLE)
80430b34957SBarry Smith           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
80530b34957SBarry Smith           colb--;
80630b34957SBarry Smith #else
80730b34957SBarry Smith           colb = aij->colmap[col] - 1;
80830b34957SBarry Smith #endif
80930b34957SBarry Smith           colb += cend - cstart;
81030b34957SBarry Smith         }
81130b34957SBarry Smith         c->vscaleforrow[k][l] = colb;
81230b34957SBarry Smith       }
81330b34957SBarry Smith     }
814b8f8c88eSHong Zhang   } else if (ctype == IS_COLORING_GHOSTED) {
815b8f8c88eSHong Zhang     /* Get gtol mapping */
816afcb2eb5SJed Brown     PetscInt N = mat->cmap->N,nlocal,*gtol;
817b8f8c88eSHong Zhang     ierr = PetscMalloc((N+1)*sizeof(PetscInt),&gtol);CHKERRQ(ierr);
818b8f8c88eSHong Zhang     for (i=0; i<N; i++) gtol[i] = -1;
819afcb2eb5SJed Brown     ierr = ISLocalToGlobalMappingGetSize(map,&nlocal);CHKERRQ(ierr);
820afcb2eb5SJed Brown     for (i=0; i<nlocal; i++) gtol[ltog[i]] = i;
821b8f8c88eSHong Zhang 
822b8f8c88eSHong Zhang     c->vscale = 0; /* will be created in MatFDColoringApply() */
823b8f8c88eSHong Zhang     ierr      = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr);
824b8f8c88eSHong Zhang     for (k=0; k<c->ncolors; k++) {
825b8f8c88eSHong Zhang       ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr);
826b8f8c88eSHong Zhang       for (l=0; l<c->nrows[k]; l++) {
827b8f8c88eSHong Zhang         col = c->columnsforrow[k][l];      /* global column index */
828b8f8c88eSHong Zhang         c->vscaleforrow[k][l] = gtol[col]; /* local column index */
829b8f8c88eSHong Zhang       }
830b8f8c88eSHong Zhang     }
831b8f8c88eSHong Zhang     ierr = PetscFree(gtol);CHKERRQ(ierr);
832b8f8c88eSHong Zhang   }
833b9617806SBarry Smith   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
83430b34957SBarry Smith 
835606d414cSSatish Balay   ierr = PetscFree(rowhit);CHKERRQ(ierr);
836606d414cSSatish Balay   ierr = PetscFree(columnsforrow);CHKERRQ(ierr);
8373acb8795SBarry Smith   ierr = MatRestoreColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr);
8383acb8795SBarry Smith   ierr = MatRestoreColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr);
839afcb2eb5SJed Brown   if (map) {ierr = ISLocalToGlobalMappingRestoreIndices(map,&ltog);CHKERRQ(ierr);}
8403a40ed3dSBarry Smith   PetscFunctionReturn(0);
841a64fbb32SBarry Smith }
842a64fbb32SBarry Smith 
843b9617806SBarry Smith 
844b9617806SBarry Smith 
845b9617806SBarry Smith 
846b9617806SBarry Smith 
847b9617806SBarry Smith 
848