xref: /petsc/src/mat/impls/aij/mpi/fdmpiaij.c (revision b3e1f37b5be90d2fa371e91c89aba21b3af857c2)
1a64fbb32SBarry Smith 
2c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h>
30d1c53f1SHong Zhang #include <../src/mat/impls/baij/mpi/mpibaij.h>
40d1c53f1SHong Zhang 
50d1c53f1SHong Zhang #undef __FUNCT__
6040ebd07SHong Zhang #define __FUNCT__ "MatFDColoringApply_BAIJ"
7040ebd07SHong Zhang PetscErrorCode  MatFDColoringApply_BAIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
80d1c53f1SHong Zhang {
90d1c53f1SHong Zhang   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void*))coloring->f;
100d1c53f1SHong Zhang   PetscErrorCode ierr;
110d1c53f1SHong Zhang   PetscInt       k,cstart,cend,l,row,col,nz,spidx,i,j;
12e2c857f8SHong Zhang   PetscScalar    dx=0.0,*xx,*w3_array,*dy_i,*dy=coloring->dy;
130d1c53f1SHong Zhang   PetscScalar    *vscale_array;
140d1c53f1SHong Zhang   PetscReal      epsilon=coloring->error_rel,umin=coloring->umin,unorm;
150d1c53f1SHong Zhang   Vec            w1=coloring->w1,w2=coloring->w2,w3,vscale=coloring->vscale;
160d1c53f1SHong Zhang   void           *fctx=coloring->fctx;
170d1c53f1SHong Zhang   PetscInt       ctype=coloring->ctype,nxloc,nrows_k;
180d1c53f1SHong Zhang   PetscScalar    *valaddr;
190d1c53f1SHong Zhang   MatEntry       *Jentry=coloring->matentry;
200d1c53f1SHong Zhang   const PetscInt ncolors=coloring->ncolors,*ncolumns=coloring->ncolumns,*nrows=coloring->nrows;
210d1c53f1SHong Zhang   PetscInt       bs=J->rmap->bs;
220d1c53f1SHong Zhang 
230d1c53f1SHong Zhang   PetscFunctionBegin;
240d1c53f1SHong Zhang   /* create vscale for storing dx */
250d1c53f1SHong Zhang   if (!vscale) {
260d1c53f1SHong Zhang     if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') {
276c145f34SHong Zhang       /* contain the "diagonal" on processor scalings followed by the off processor* - garray must be non-bloced! */
286c145f34SHong Zhang       Mat_MPIBAIJ    *baij=(Mat_MPIBAIJ*)J->data;
296c145f34SHong Zhang       PetscInt       *garray;
306c145f34SHong Zhang       ierr = PetscMalloc(baij->B->cmap->n*sizeof(PetscInt),&garray);CHKERRQ(ierr);
316c145f34SHong Zhang       for (i=0; i<baij->B->cmap->n/bs; i++) {
326c145f34SHong Zhang         for (j=0; j<bs; j++) {
336c145f34SHong Zhang           garray[i*bs+j] = bs*baij->garray[i]+j;
346c145f34SHong Zhang         }
356c145f34SHong Zhang       }
366c145f34SHong Zhang       ierr = VecCreateGhost(PetscObjectComm((PetscObject)J),J->cmap->n,PETSC_DETERMINE,baij->B->cmap->n,garray,&vscale);CHKERRQ(ierr);
376c145f34SHong Zhang       ierr = PetscFree(garray);CHKERRQ(ierr);
380d1c53f1SHong Zhang     } else if (ctype == IS_COLORING_GHOSTED) {
390d1c53f1SHong Zhang       ierr = VecDuplicate(x1,&vscale);CHKERRQ(ierr);
400d1c53f1SHong Zhang     }
410d1c53f1SHong Zhang     coloring->vscale = vscale;
420d1c53f1SHong Zhang   }
430d1c53f1SHong Zhang 
440d1c53f1SHong Zhang   /* (1) Set w1 = F(x1) */
450d1c53f1SHong Zhang   if (!coloring->fset) {
460d1c53f1SHong Zhang     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
470d1c53f1SHong Zhang     ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
480d1c53f1SHong Zhang     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
490d1c53f1SHong Zhang   } else {
500d1c53f1SHong Zhang     coloring->fset = PETSC_FALSE;
510d1c53f1SHong Zhang   }
520d1c53f1SHong Zhang 
530d1c53f1SHong Zhang   /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */
540d1c53f1SHong Zhang   ierr = VecGetLocalSize(x1,&nxloc);CHKERRQ(ierr);
550d1c53f1SHong Zhang   if (coloring->htype[0] == 'w') {
560d1c53f1SHong Zhang     /* vscale = dx is a constant scalar */
570d1c53f1SHong Zhang     ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr);
580d1c53f1SHong Zhang     dx = 1.0/(PetscSqrtReal(1.0 + unorm)*epsilon);
590d1c53f1SHong Zhang   } else {
600d1c53f1SHong Zhang     ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);
610d1c53f1SHong Zhang     ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr);
620d1c53f1SHong Zhang     for (col=0; col<nxloc; col++) {
630d1c53f1SHong Zhang       dx = xx[col];
640d1c53f1SHong Zhang       if (PetscAbsScalar(dx) < umin) {
650d1c53f1SHong Zhang         if (PetscRealPart(dx) >= 0.0)      dx = umin;
660d1c53f1SHong Zhang         else if (PetscRealPart(dx) < 0.0 ) dx = -umin;
670d1c53f1SHong Zhang       }
680d1c53f1SHong Zhang       dx               *= epsilon;
690d1c53f1SHong Zhang       vscale_array[col] = 1.0/dx;
700d1c53f1SHong Zhang     }
710d1c53f1SHong Zhang     ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
720d1c53f1SHong Zhang     ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr);
730d1c53f1SHong Zhang   }
74684f2004SHong Zhang   if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') {
750d1c53f1SHong Zhang     ierr = VecGhostUpdateBegin(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
760d1c53f1SHong Zhang     ierr = VecGhostUpdateEnd(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
770d1c53f1SHong Zhang   }
780d1c53f1SHong Zhang 
790d1c53f1SHong Zhang   /* (3) Loop over each color */
800d1c53f1SHong Zhang   if (!coloring->w3) {
810d1c53f1SHong Zhang     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
820d1c53f1SHong Zhang     ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr);
830d1c53f1SHong Zhang   }
840d1c53f1SHong Zhang   w3 = coloring->w3;
850d1c53f1SHong Zhang 
860d1c53f1SHong Zhang   ierr = VecGetOwnershipRange(x1,&cstart,&cend);CHKERRQ(ierr); /* used by ghosted vscale */
870d1c53f1SHong Zhang   if (vscale) {
880d1c53f1SHong Zhang     ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr);
890d1c53f1SHong Zhang   }
900d1c53f1SHong Zhang   nz   = 0;
910d1c53f1SHong Zhang   for (k=0; k<ncolors; k++) {
920d1c53f1SHong Zhang     coloring->currentcolor = k;
930d1c53f1SHong Zhang 
940d1c53f1SHong Zhang     /*
950d1c53f1SHong Zhang       (3-1) Loop over each column associated with color
960d1c53f1SHong Zhang       adding the perturbation to the vector w3 = x1 + dx.
970d1c53f1SHong Zhang     */
980d1c53f1SHong Zhang     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
990d1c53f1SHong Zhang     dy_i = dy;
1000d1c53f1SHong Zhang     for (i=0; i<bs; i++) {     /* Loop over a block of columns */
1010d1c53f1SHong Zhang       ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
1020d1c53f1SHong Zhang       if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */
1030d1c53f1SHong Zhang       if (coloring->htype[0] == 'w') {
1040d1c53f1SHong Zhang         for (l=0; l<ncolumns[k]; l++) {
1050d1c53f1SHong Zhang           col            = i + bs*coloring->columns[k][l];  /* local column (in global index!) of the matrix we are probing for */
1060d1c53f1SHong Zhang           w3_array[col] += 1.0/dx;
107684f2004SHong Zhang           if (i) w3_array[col-1] -= 1.0/dx; /* resume original w3[col-1] */
1080d1c53f1SHong Zhang         }
1090d1c53f1SHong Zhang       } else { /* htype == 'ds' */
1100d1c53f1SHong Zhang         vscale_array -= cstart; /* shift pointer so global index can be used */
1110d1c53f1SHong Zhang         for (l=0; l<ncolumns[k]; l++) {
112f8c2866eSHong Zhang           col = i + bs*coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
1130d1c53f1SHong Zhang           w3_array[col] += 1.0/vscale_array[col];
114684f2004SHong Zhang           if (i) w3_array[col-1] -=  1.0/vscale_array[col-1]; /* resume original w3[col-1] */
1150d1c53f1SHong Zhang         }
1160d1c53f1SHong Zhang         vscale_array += cstart;
1170d1c53f1SHong Zhang       }
1180d1c53f1SHong Zhang       if (ctype == IS_COLORING_GLOBAL) w3_array += cstart;
1190d1c53f1SHong Zhang       ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
1200d1c53f1SHong Zhang 
1210d1c53f1SHong Zhang       /*
1220d1c53f1SHong Zhang        (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
1230d1c53f1SHong Zhang                            w2 = F(x1 + dx) - F(x1)
1240d1c53f1SHong Zhang        */
1250d1c53f1SHong Zhang       ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
1260d1c53f1SHong Zhang       ierr = VecPlaceArray(w2,dy_i);CHKERRQ(ierr); /* place w2 to the array dy_i */
1270d1c53f1SHong Zhang       ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
1280d1c53f1SHong Zhang       ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
1290d1c53f1SHong Zhang       ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
1300d1c53f1SHong Zhang       ierr = VecResetArray(w2);CHKERRQ(ierr);
1310d1c53f1SHong Zhang       dy_i += nxloc; /* points to dy+i*nxloc */
1320d1c53f1SHong Zhang     }
1330d1c53f1SHong Zhang 
1340d1c53f1SHong Zhang     /*
1350d1c53f1SHong Zhang      (3-3) Loop over rows of vector, putting results into Jacobian matrix
1360d1c53f1SHong Zhang     */
1370d1c53f1SHong Zhang     nrows_k = nrows[k];
1380d1c53f1SHong Zhang     if (coloring->htype[0] == 'w') {
1390d1c53f1SHong Zhang       for (l=0; l<nrows_k; l++) {
1400d1c53f1SHong Zhang         row     = bs*Jentry[nz].row;   /* local row index */
141684f2004SHong Zhang         valaddr = Jentry[nz++].valaddr;
1420d1c53f1SHong Zhang         spidx   = 0;
1430d1c53f1SHong Zhang         dy_i    = dy;
1440d1c53f1SHong Zhang         for (i=0; i<bs; i++) {   /* column of the block */
1450d1c53f1SHong Zhang           for (j=0; j<bs; j++) { /* row of the block */
1460d1c53f1SHong Zhang             valaddr[spidx++] = dy_i[row+j]*dx;
1470d1c53f1SHong Zhang           }
148f8c2866eSHong Zhang           dy_i += nxloc; /* points to dy+i*nxloc */
1490d1c53f1SHong Zhang         }
1500d1c53f1SHong Zhang       }
1510d1c53f1SHong Zhang     } else { /* htype == 'ds' */
1520d1c53f1SHong Zhang       for (l=0; l<nrows_k; l++) {
1530d1c53f1SHong Zhang         row     = bs*Jentry[nz].row;   /* local row index */
1540d1c53f1SHong Zhang         col     = bs*Jentry[nz].col;   /* local column index */
155684f2004SHong Zhang         valaddr = Jentry[nz++].valaddr;
156f8c2866eSHong Zhang         spidx   = 0;
157f8c2866eSHong Zhang         dy_i    = dy;
158f8c2866eSHong Zhang         for (i=0; i<bs; i++) {   /* column of the block */
159f8c2866eSHong Zhang           for (j=0; j<bs; j++) { /* row of the block */
160f8c2866eSHong Zhang             valaddr[spidx++] = dy_i[row+j]*vscale_array[col+i];
161f8c2866eSHong Zhang           }
162f8c2866eSHong Zhang           dy_i += nxloc; /* points to dy+i*nxloc */
163f8c2866eSHong Zhang         }
1640d1c53f1SHong Zhang       }
1650d1c53f1SHong Zhang     }
1660d1c53f1SHong Zhang   }
1670d1c53f1SHong Zhang   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1680d1c53f1SHong Zhang   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1690d1c53f1SHong Zhang   if (vscale) {
1700d1c53f1SHong Zhang     ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr);
1710d1c53f1SHong Zhang   }
1720d1c53f1SHong Zhang 
1730d1c53f1SHong Zhang   coloring->currentcolor = -1;
1740d1c53f1SHong Zhang   PetscFunctionReturn(0);
1750d1c53f1SHong Zhang }
176a64fbb32SBarry Smith 
1774a2ae208SSatish Balay #undef __FUNCT__
178040ebd07SHong Zhang #define __FUNCT__ "MatFDColoringApply_AIJ"
179040ebd07SHong Zhang PetscErrorCode  MatFDColoringApply_AIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
180fcd7ac73SHong Zhang {
181fcd7ac73SHong Zhang   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void*))coloring->f;
182fcd7ac73SHong Zhang   PetscErrorCode ierr;
183a2f2d239SHong Zhang   PetscInt       k,cstart,cend,l,row,col,nz;
1849e917edbSHong Zhang   PetscScalar    dx=0.0,*y,*xx,*w3_array;
185fcd7ac73SHong Zhang   PetscScalar    *vscale_array;
186fcd7ac73SHong Zhang   PetscReal      epsilon=coloring->error_rel,umin=coloring->umin,unorm;
1878bc97078SHong Zhang   Vec            w1=coloring->w1,w2=coloring->w2,w3,vscale=coloring->vscale;
188fcd7ac73SHong Zhang   void           *fctx=coloring->fctx;
189c53567a0SHong Zhang   PetscInt       ctype=coloring->ctype,nxloc,nrows_k;
190573f477fSHong Zhang   MatEntry       *Jentry=coloring->matentry;
1918bc97078SHong Zhang   const PetscInt ncolors=coloring->ncolors,*ncolumns=coloring->ncolumns,*nrows=coloring->nrows;
192fcd7ac73SHong Zhang 
193fcd7ac73SHong Zhang   PetscFunctionBegin;
1949e917edbSHong Zhang   /* create vscale for storing dx */
1959e917edbSHong Zhang   if (!vscale) {
1969e917edbSHong Zhang     if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') {
1976c145f34SHong Zhang       Mat_MPIAIJ     *aij=(Mat_MPIAIJ*)J->data;
1989e917edbSHong Zhang       ierr = VecCreateGhost(PetscObjectComm((PetscObject)J),J->cmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&vscale);CHKERRQ(ierr);
1999e917edbSHong Zhang     } else if (ctype == IS_COLORING_GHOSTED) {
2008bc97078SHong Zhang       ierr = VecDuplicate(x1,&vscale);CHKERRQ(ierr);
2019e917edbSHong Zhang     }
2028bc97078SHong Zhang     coloring->vscale = vscale;
203fcd7ac73SHong Zhang   }
204fcd7ac73SHong Zhang 
2058bc97078SHong Zhang   /* (1) Set w1 = F(x1) */
206fcd7ac73SHong Zhang   if (!coloring->fset) {
207fcd7ac73SHong Zhang     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
208f6af9589SHong Zhang     ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
209fcd7ac73SHong Zhang     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
210fcd7ac73SHong Zhang   } else {
211fcd7ac73SHong Zhang     coloring->fset = PETSC_FALSE;
212fcd7ac73SHong Zhang   }
213fcd7ac73SHong Zhang 
2148bc97078SHong Zhang   /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */
215f6af9589SHong Zhang   if (coloring->htype[0] == 'w') {
216684f2004SHong Zhang     /* vscale = 1./dx is a constant scalar */
217f6af9589SHong Zhang     ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr);
218c53567a0SHong Zhang     dx = 1.0/(PetscSqrtReal(1.0 + unorm)*epsilon);
21970e7395fSHong Zhang   } else {
22074d3cef9SHong Zhang     ierr = VecGetLocalSize(x1,&nxloc);CHKERRQ(ierr);
221f6af9589SHong Zhang     ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);
2228bc97078SHong Zhang     ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr);
22374d3cef9SHong Zhang     for (col=0; col<nxloc; col++) {
224fcd7ac73SHong Zhang       dx = xx[col];
22574d3cef9SHong Zhang       if (PetscAbsScalar(dx) < umin) {
22674d3cef9SHong Zhang         if (PetscRealPart(dx) >= 0.0)      dx = umin;
22774d3cef9SHong Zhang         else if (PetscRealPart(dx) < 0.0 ) dx = -umin;
22874d3cef9SHong Zhang       }
229fcd7ac73SHong Zhang       dx               *= epsilon;
23074d3cef9SHong Zhang       vscale_array[col] = 1.0/dx;
231f6af9589SHong Zhang     }
23274d3cef9SHong Zhang     ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
2338bc97078SHong Zhang     ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr);
23470e7395fSHong Zhang   }
235684f2004SHong Zhang   if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') {
2368bc97078SHong Zhang     ierr = VecGhostUpdateBegin(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2378bc97078SHong Zhang     ierr = VecGhostUpdateEnd(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
238fcd7ac73SHong Zhang   }
239fcd7ac73SHong Zhang 
2408bc97078SHong Zhang   /* (3) Loop over each color */
2418bc97078SHong Zhang   if (!coloring->w3) {
2428bc97078SHong Zhang     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
2438bc97078SHong Zhang     ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr);
2448bc97078SHong Zhang   }
2458bc97078SHong Zhang   w3 = coloring->w3;
246fcd7ac73SHong Zhang 
2478bc97078SHong Zhang   ierr = VecGetOwnershipRange(x1,&cstart,&cend);CHKERRQ(ierr); /* used by ghosted vscale */
2489e917edbSHong Zhang   if (vscale) {
2498bc97078SHong Zhang     ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr);
2509e917edbSHong Zhang   }
2518bc97078SHong Zhang   nz   = 0;
252e2c857f8SHong Zhang 
253*b3e1f37bSHong Zhang   if (coloring->bcols > 1) { /* use blocked insertion of Jentry */
254*b3e1f37bSHong Zhang     PetscInt    i,m=J->rmap->n,nbcols,bcols=coloring->bcols;
255e2c857f8SHong Zhang     PetscScalar *dy=coloring->dy,*dy_k;
256e2c857f8SHong Zhang 
257e2c857f8SHong Zhang     nbcols = 0;
258e2c857f8SHong Zhang     for (k=0; k<ncolors; k+=bcols) {
259e2c857f8SHong Zhang       coloring->currentcolor = k;
260e2c857f8SHong Zhang 
261e2c857f8SHong Zhang       /*
262e2c857f8SHong Zhang        (3-1) Loop over each column associated with color
263e2c857f8SHong Zhang        adding the perturbation to the vector w3 = x1 + dx.
264e2c857f8SHong Zhang        */
265e2c857f8SHong Zhang 
266e2c857f8SHong Zhang       dy_k = dy;
267e2c857f8SHong Zhang       if (k + bcols > ncolors) bcols = ncolors - k;
268e2c857f8SHong Zhang       for (i=0; i<bcols; i++) {
269e2c857f8SHong Zhang 
270e2c857f8SHong Zhang         ierr = VecCopy(x1,w3);CHKERRQ(ierr);
271e2c857f8SHong Zhang         ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
272e2c857f8SHong Zhang         if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */
273e2c857f8SHong Zhang         if (coloring->htype[0] == 'w') {
274e2c857f8SHong Zhang           for (l=0; l<ncolumns[k+i]; l++) {
275e2c857f8SHong Zhang             col = coloring->columns[k+i][l]; /* local column (in global index!) of the matrix we are probing for */
276e2c857f8SHong Zhang             w3_array[col] += 1.0/dx;
277e2c857f8SHong Zhang           }
278e2c857f8SHong Zhang         } else { /* htype == 'ds' */
279e2c857f8SHong Zhang           vscale_array -= cstart; /* shift pointer so global index can be used */
280e2c857f8SHong Zhang           for (l=0; l<ncolumns[k+i]; l++) {
281e2c857f8SHong Zhang             col = coloring->columns[k+i][l]; /* local column (in global index!) of the matrix we are probing for */
282e2c857f8SHong Zhang             w3_array[col] += 1.0/vscale_array[col];
283e2c857f8SHong Zhang           }
284e2c857f8SHong Zhang           vscale_array += cstart;
285e2c857f8SHong Zhang         }
286e2c857f8SHong Zhang         if (ctype == IS_COLORING_GLOBAL) w3_array += cstart;
287e2c857f8SHong Zhang         ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
288e2c857f8SHong Zhang 
289e2c857f8SHong Zhang         /*
290e2c857f8SHong Zhang          (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
291e2c857f8SHong Zhang                            w2 = F(x1 + dx) - F(x1)
292e2c857f8SHong Zhang          */
293e2c857f8SHong Zhang         ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
294e2c857f8SHong Zhang         ierr = VecPlaceArray(w2,dy_k);CHKERRQ(ierr); /* place w2 to the array dy_i */
295e2c857f8SHong Zhang         ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
296e2c857f8SHong Zhang         ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
297e2c857f8SHong Zhang         ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
298e2c857f8SHong Zhang         ierr = VecResetArray(w2);CHKERRQ(ierr);
299e2c857f8SHong Zhang         dy_k += m; /* points to dy+i*nxloc */
300e2c857f8SHong Zhang 
301e2c857f8SHong Zhang       }
302e2c857f8SHong Zhang 
303e2c857f8SHong Zhang       /*
304e2c857f8SHong Zhang        (3-3) Loop over block rows of vector, putting results into Jacobian matrix
305e2c857f8SHong Zhang        */
306d880da65SHong Zhang       nrows_k = nrows[nbcols++];
307e2c857f8SHong Zhang       ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
308d880da65SHong Zhang 
309e2c857f8SHong Zhang       if (coloring->htype[0] == 'w') {
310e2c857f8SHong Zhang         for (l=0; l<nrows_k; l++) {
311d880da65SHong Zhang           row                     = Jentry[nz].row;   /* local row index */
312d880da65SHong Zhang           *(Jentry[nz++].valaddr) = dy[row]*dx;
313e2c857f8SHong Zhang         }
314e2c857f8SHong Zhang       } else { /* htype == 'ds' */
315e2c857f8SHong Zhang         for (l=0; l<nrows_k; l++) {
316e2c857f8SHong Zhang           row                   = Jentry[nz].row;   /* local row index */
317d880da65SHong Zhang           *(Jentry[nz].valaddr) = dy[row]*vscale_array[Jentry[nz].col];
318e2c857f8SHong Zhang           nz++;
319e2c857f8SHong Zhang         }
320e2c857f8SHong Zhang       }
321e2c857f8SHong Zhang       ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
322e2c857f8SHong Zhang     }
323*b3e1f37bSHong Zhang   } else { /* bcols == 1 */
3248bc97078SHong Zhang     for (k=0; k<ncolors; k++) {
3258bc97078SHong Zhang       coloring->currentcolor = k;
3269e917edbSHong Zhang 
327fcd7ac73SHong Zhang       /*
3288bc97078SHong Zhang        (3-1) Loop over each column associated with color
329f6af9589SHong Zhang        adding the perturbation to the vector w3 = x1 + dx.
330fcd7ac73SHong Zhang        */
331f6af9589SHong Zhang       ierr = VecCopy(x1,w3);CHKERRQ(ierr);
332f6af9589SHong Zhang       ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
333a2f2d239SHong Zhang       if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */
334b039d6c4SHong Zhang       if (coloring->htype[0] == 'w') {
3358bc97078SHong Zhang         for (l=0; l<ncolumns[k]; l++) {
336f6af9589SHong Zhang           col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
3379e917edbSHong Zhang           w3_array[col] += 1.0/dx;
3389e917edbSHong Zhang         }
339b039d6c4SHong Zhang       } else { /* htype == 'ds' */
340a2f2d239SHong Zhang         vscale_array -= cstart; /* shift pointer so global index can be used */
341b039d6c4SHong Zhang         for (l=0; l<ncolumns[k]; l++) {
342b039d6c4SHong Zhang           col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
343a2f2d239SHong Zhang           w3_array[col] += 1.0/vscale_array[col];
34470e7395fSHong Zhang         }
345a2f2d239SHong Zhang         vscale_array += cstart;
346fcd7ac73SHong Zhang       }
347a2f2d239SHong Zhang       if (ctype == IS_COLORING_GLOBAL) w3_array += cstart;
348fcd7ac73SHong Zhang       ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
3499e917edbSHong Zhang 
350fcd7ac73SHong Zhang       /*
3518bc97078SHong Zhang        (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
352fcd7ac73SHong Zhang                            w2 = F(x1 + dx) - F(x1)
353fcd7ac73SHong Zhang        */
354fcd7ac73SHong Zhang       ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
355fcd7ac73SHong Zhang       ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
356fcd7ac73SHong Zhang       ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
357fcd7ac73SHong Zhang       ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
3589e917edbSHong Zhang 
3598bc97078SHong Zhang       /*
3608bc97078SHong Zhang        (3-3) Loop over rows of vector, putting results into Jacobian matrix
3618bc97078SHong Zhang        */
362c53567a0SHong Zhang       nrows_k = nrows[k];
363fcd7ac73SHong Zhang       ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
364b039d6c4SHong Zhang       if (coloring->htype[0] == 'w') {
365c53567a0SHong Zhang         for (l=0; l<nrows_k; l++) {
366573f477fSHong Zhang           row                     = Jentry[nz].row;   /* local row index */
367b039d6c4SHong Zhang           *(Jentry[nz++].valaddr) = y[row]*dx;
3689e917edbSHong Zhang         }
369b039d6c4SHong Zhang       } else { /* htype == 'ds' */
370c53567a0SHong Zhang         for (l=0; l<nrows_k; l++) {
371b039d6c4SHong Zhang           row                   = Jentry[nz].row;   /* local row index */
372b039d6c4SHong Zhang           *(Jentry[nz].valaddr) = y[row]*vscale_array[Jentry[nz].col];
373573f477fSHong Zhang           nz++;
374fcd7ac73SHong Zhang         }
375b039d6c4SHong Zhang       }
376fcd7ac73SHong Zhang       ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
3778bc97078SHong Zhang     }
378*b3e1f37bSHong Zhang   }
379*b3e1f37bSHong Zhang 
380f5aae955SHong Zhang   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
381f5aae955SHong Zhang   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3829e917edbSHong Zhang   if (vscale) {
3838bc97078SHong Zhang     ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr);
3849e917edbSHong Zhang   }
385fcd7ac73SHong Zhang   coloring->currentcolor = -1;
386fcd7ac73SHong Zhang   PetscFunctionReturn(0);
387fcd7ac73SHong Zhang }
388fcd7ac73SHong Zhang 
389fcd7ac73SHong Zhang #undef __FUNCT__
390f86b9fbaSHong Zhang #define __FUNCT__ "MatFDColoringSetUp_MPIXAIJ"
391f86b9fbaSHong Zhang PetscErrorCode MatFDColoringSetUp_MPIXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
392fcd7ac73SHong Zhang {
393fcd7ac73SHong Zhang   PetscErrorCode         ierr;
394fcd7ac73SHong Zhang   PetscMPIInt            size,*ncolsonproc,*disp,nn;
395e3225e9dSHong Zhang   PetscInt               i,n,nrows,nrows_i,j,k,m,ncols,col,*rowhit,cstart,cend,colb;
39672c15787SHong Zhang   const PetscInt         *is,*A_ci,*A_cj,*B_ci,*B_cj,*row=NULL,*ltog=NULL;
397fcd7ac73SHong Zhang   PetscInt               nis=iscoloring->n,nctot,*cols;
398fcd7ac73SHong Zhang   IS                     *isa;
399fcd7ac73SHong Zhang   ISLocalToGlobalMapping map=mat->cmap->mapping;
400a8971b87SHong Zhang   PetscInt               ctype=c->ctype,*spidxA,*spidxB,nz,bs,bs2,spidx;
4010d1c53f1SHong Zhang   Mat                    A,B;
402e3225e9dSHong Zhang   PetscScalar            *A_val,*B_val,**valaddrhit;
403573f477fSHong Zhang   MatEntry               *Jentry;
4040d1c53f1SHong Zhang   PetscBool              isBAIJ;
405531e53bdSHong Zhang   PetscInt               bcols=c->bcols;
4060d1c53f1SHong Zhang #if defined(PETSC_USE_CTABLE)
4070d1c53f1SHong Zhang   PetscTable             colmap=NULL;
4080d1c53f1SHong Zhang #else
4090d1c53f1SHong Zhang   PetscInt               *colmap=NULL;     /* local col number of off-diag col */
4100d1c53f1SHong Zhang #endif
411fcd7ac73SHong Zhang 
412fcd7ac73SHong Zhang   PetscFunctionBegin;
41372c15787SHong Zhang   if (ctype == IS_COLORING_GHOSTED) {
41472c15787SHong 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");
41572c15787SHong Zhang     ierr = ISLocalToGlobalMappingGetIndices(map,&ltog);CHKERRQ(ierr);
41672c15787SHong Zhang   }
417fcd7ac73SHong Zhang 
4180d1c53f1SHong Zhang   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
4190d1c53f1SHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr);
420e3225e9dSHong Zhang   if (isBAIJ) {
4210d1c53f1SHong Zhang     Mat_MPIBAIJ *aij=(Mat_MPIBAIJ*)mat->data;
4226c145f34SHong Zhang     Mat_SeqBAIJ *spA,*spB;
4236c145f34SHong Zhang     A = aij->A;  spA = (Mat_SeqBAIJ*)A->data; A_val = spA->a;
4246c145f34SHong Zhang     B = aij->B;  spB = (Mat_SeqBAIJ*)B->data; B_val = spB->a;
4250d1c53f1SHong Zhang     nz = spA->nz + spB->nz; /* total nonzero entries of mat */
4260d1c53f1SHong Zhang     if (!aij->colmap) {
4270d1c53f1SHong Zhang       ierr = MatCreateColmap_MPIBAIJ_Private(mat);CHKERRQ(ierr);
4280d1c53f1SHong Zhang       colmap = aij->colmap;
4290d1c53f1SHong Zhang     }
4300d1c53f1SHong Zhang     ierr = MatGetColumnIJ_SeqBAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
4310d1c53f1SHong Zhang     ierr = MatGetColumnIJ_SeqBAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
432e3225e9dSHong Zhang   } else {
433e3225e9dSHong Zhang     Mat_MPIAIJ *aij=(Mat_MPIAIJ*)mat->data;
434e3225e9dSHong Zhang     Mat_SeqAIJ *spA,*spB;
435e3225e9dSHong Zhang     A = aij->A;  spA = (Mat_SeqAIJ*)A->data; A_val = spA->a;
436e3225e9dSHong Zhang     B = aij->B;  spB = (Mat_SeqAIJ*)B->data; B_val = spB->a;
437e3225e9dSHong Zhang     nz = spA->nz + spB->nz; /* total nonzero entries of mat */
438e3225e9dSHong Zhang     if (!aij->colmap) {
439e3225e9dSHong Zhang       /* Allow access to data structures of local part of matrix
440e3225e9dSHong Zhang        - creates aij->colmap which maps global column number to local number in part B */
441e3225e9dSHong Zhang       ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
442e3225e9dSHong Zhang       colmap = aij->colmap;
443e3225e9dSHong Zhang     }
444e3225e9dSHong Zhang     ierr = MatGetColumnIJ_SeqAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
445e3225e9dSHong Zhang     ierr = MatGetColumnIJ_SeqAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
446e3225e9dSHong Zhang 
447e3225e9dSHong Zhang     bs = 1; /* only bs=1 is supported for non MPIBAIJ matrix */
4480d1c53f1SHong Zhang   }
4490d1c53f1SHong Zhang 
4500d1c53f1SHong Zhang   m         = mat->rmap->n/bs;
4510d1c53f1SHong Zhang   cstart    = mat->cmap->rstart/bs;
4520d1c53f1SHong Zhang   cend      = mat->cmap->rend/bs;
453fcd7ac73SHong Zhang 
454fcd7ac73SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
455fcd7ac73SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr);
456fcd7ac73SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
4578bc97078SHong Zhang   ierr       = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr);
458fcd7ac73SHong Zhang 
459573f477fSHong Zhang   ierr       = PetscMalloc(nz*sizeof(MatEntry),&Jentry);CHKERRQ(ierr);
4608bc97078SHong Zhang   ierr       = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr);
4618bc97078SHong Zhang   c->matentry = Jentry;
462a774a6f1SHong Zhang 
463a8971b87SHong Zhang   ierr = PetscMalloc2(m+1,PetscInt,&rowhit,m+1,PetscScalar*,&valaddrhit);CHKERRQ(ierr);
464fcd7ac73SHong Zhang   nz = 0;
46572c15787SHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
466d3825b63SHong Zhang   for (i=0; i<nis; i++) { /* for each local color */
467fcd7ac73SHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
468fcd7ac73SHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
469fcd7ac73SHong Zhang 
470fcd7ac73SHong Zhang     c->ncolumns[i] = n; /* local number of columns of this color on this process */
471fcd7ac73SHong Zhang     if (n) {
472fcd7ac73SHong Zhang       ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr);
473fcd7ac73SHong Zhang       ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr);
474fcd7ac73SHong Zhang       ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr);
475fcd7ac73SHong Zhang     } else {
476fcd7ac73SHong Zhang       c->columns[i] = 0;
477fcd7ac73SHong Zhang     }
478fcd7ac73SHong Zhang 
479fcd7ac73SHong Zhang     if (ctype == IS_COLORING_GLOBAL) {
480d3825b63SHong Zhang       /* Determine nctot, the total (parallel) number of columns of this color */
481fcd7ac73SHong Zhang       ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr);
482fcd7ac73SHong Zhang       ierr = PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);CHKERRQ(ierr);
483fcd7ac73SHong Zhang 
484d3825b63SHong Zhang       /* ncolsonproc[j]: local ncolumns on proc[j] of this color */
485fcd7ac73SHong Zhang       ierr  = PetscMPIIntCast(n,&nn);CHKERRQ(ierr);
486fcd7ac73SHong Zhang       ierr  = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
487fcd7ac73SHong Zhang       nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j];
488fcd7ac73SHong Zhang       if (!nctot) {
489fcd7ac73SHong Zhang         ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr);
490fcd7ac73SHong Zhang       }
491fcd7ac73SHong Zhang 
492fcd7ac73SHong Zhang       disp[0] = 0;
493fcd7ac73SHong Zhang       for (j=1; j<size; j++) {
494fcd7ac73SHong Zhang         disp[j] = disp[j-1] + ncolsonproc[j-1];
495fcd7ac73SHong Zhang       }
496fcd7ac73SHong Zhang 
497d3825b63SHong Zhang       /* Get cols, the complete list of columns for this color on each process */
498fcd7ac73SHong Zhang       ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
499fcd7ac73SHong Zhang       ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
500fcd7ac73SHong Zhang       ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr);
501fcd7ac73SHong Zhang     } else if (ctype == IS_COLORING_GHOSTED) {
502fcd7ac73SHong Zhang       /* Determine local number of columns of this color on this process, including ghost points */
503fcd7ac73SHong Zhang       nctot = n;
504fcd7ac73SHong Zhang       ierr  = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
505fcd7ac73SHong Zhang       ierr  = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr);
506fcd7ac73SHong Zhang     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type");
507fcd7ac73SHong Zhang 
5081b97d346SHong Zhang     /* Mark all rows affect by these columns */
509d3825b63SHong Zhang     ierr    = PetscMemzero(rowhit,m*sizeof(PetscInt));CHKERRQ(ierr);
5100d1c53f1SHong Zhang     bs2     = bs*bs;
5114b2e90caSHong Zhang     nrows_i = 0;
5121b97d346SHong Zhang     for (j=0; j<nctot; j++) { /* loop over columns*/
513fcd7ac73SHong Zhang       if (ctype == IS_COLORING_GHOSTED) {
514fcd7ac73SHong Zhang         col = ltog[cols[j]];
515fcd7ac73SHong Zhang       } else {
516fcd7ac73SHong Zhang         col = cols[j];
517fcd7ac73SHong Zhang       }
518e3225e9dSHong Zhang       if (col >= cstart && col < cend) { /* column is in A, diagonal block of mat */
519d3825b63SHong Zhang         row      = A_cj + A_ci[col-cstart];
520d3825b63SHong Zhang         nrows    = A_ci[col-cstart+1] - A_ci[col-cstart];
5214b2e90caSHong Zhang         nrows_i += nrows;
522fcd7ac73SHong Zhang         /* loop over columns of A marking them in rowhit */
523d3825b63SHong Zhang         for (k=0; k<nrows; k++) {
524d3825b63SHong Zhang           /* set valaddrhit for part A */
525e3225e9dSHong Zhang           spidx            = bs2*spidxA[A_ci[col-cstart] + k];
526e3225e9dSHong Zhang           valaddrhit[*row] = &A_val[spidx];
527a774a6f1SHong Zhang           rowhit[*row++]   = col - cstart + 1; /* local column index */
528fcd7ac73SHong Zhang         }
529e3225e9dSHong Zhang       } else { /* column is in B, off-diagonal block of mat */
530fcd7ac73SHong Zhang #if defined(PETSC_USE_CTABLE)
5310d1c53f1SHong Zhang         ierr = PetscTableFind(colmap,col+1,&colb);CHKERRQ(ierr);
532fcd7ac73SHong Zhang         colb--;
533fcd7ac73SHong Zhang #else
5340d1c53f1SHong Zhang         colb = colmap[col] - 1; /* local column index */
535fcd7ac73SHong Zhang #endif
536fcd7ac73SHong Zhang         if (colb == -1) {
537d3825b63SHong Zhang           nrows = 0;
538fcd7ac73SHong Zhang         } else {
5390d1c53f1SHong Zhang           colb  = colb/bs;
540d3825b63SHong Zhang           row   = B_cj + B_ci[colb];
541d3825b63SHong Zhang           nrows = B_ci[colb+1] - B_ci[colb];
542fcd7ac73SHong Zhang         }
5434b2e90caSHong Zhang         nrows_i += nrows;
544fcd7ac73SHong Zhang         /* loop over columns of B marking them in rowhit */
545d3825b63SHong Zhang         for (k=0; k<nrows; k++) {
546d3825b63SHong Zhang           /* set valaddrhit for part B */
547e3225e9dSHong Zhang           spidx            = bs2*spidxB[B_ci[colb] + k];
548e3225e9dSHong Zhang           valaddrhit[*row] = &B_val[spidx];
54970e7395fSHong Zhang           rowhit[*row++]   = colb + 1 + cend - cstart; /* local column index */
550fcd7ac73SHong Zhang         }
551fcd7ac73SHong Zhang       }
552fcd7ac73SHong Zhang     }
5534b2e90caSHong Zhang     c->nrows[i] = nrows_i;
5548bc97078SHong Zhang 
555d3825b63SHong Zhang     for (j=0; j<m; j++) {
556fcd7ac73SHong Zhang       if (rowhit[j]) {
557573f477fSHong Zhang         Jentry[nz].row     = j;              /* local row index */
558573f477fSHong Zhang         Jentry[nz].col     = rowhit[j] - 1;  /* local column index */
559573f477fSHong Zhang         Jentry[nz].valaddr = valaddrhit[j];  /* address of mat value for this entry */
560573f477fSHong Zhang         nz++;
561fcd7ac73SHong Zhang       }
562fcd7ac73SHong Zhang     }
563fcd7ac73SHong Zhang     ierr = PetscFree(cols);CHKERRQ(ierr);
564fcd7ac73SHong Zhang   }
565fcd7ac73SHong Zhang 
566a8971b87SHong Zhang   if (bcols > 1) { /* reorder Jentry for faster MatFDColoringApply() */
567a8971b87SHong Zhang     ierr = MatFDColoringSetUpBlocked_AIJ_Private(mat,c,nz);CHKERRQ(ierr);
568531e53bdSHong Zhang   }
569531e53bdSHong Zhang 
5700d1c53f1SHong Zhang   if (isBAIJ) {
5710d1c53f1SHong Zhang     ierr = MatRestoreColumnIJ_SeqBAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
5720d1c53f1SHong Zhang     ierr = MatRestoreColumnIJ_SeqBAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
573e3225e9dSHong Zhang     ierr = PetscMalloc(bs*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr);
5740d1c53f1SHong Zhang   } else {
5750d1c53f1SHong Zhang     ierr = MatRestoreColumnIJ_SeqAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
5760d1c53f1SHong Zhang     ierr = MatRestoreColumnIJ_SeqAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
5770d1c53f1SHong Zhang   }
5780d1c53f1SHong Zhang 
579a8971b87SHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
580a8971b87SHong Zhang   ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr);
581a8971b87SHong Zhang 
58272c15787SHong Zhang   if (ctype == IS_COLORING_GHOSTED) {
58372c15787SHong Zhang     ierr = ISLocalToGlobalMappingRestoreIndices(map,&ltog);CHKERRQ(ierr);
58472c15787SHong Zhang   }
585*b3e1f37bSHong Zhang #if defined(PETSC_USE_INFO)
586*b3e1f37bSHong Zhang   ierr = PetscInfo3(c,"ncolors %D, brows %D and bcols %D are used.\n",c->ncolors,c->brows,c->bcols);CHKERRQ(ierr);
587*b3e1f37bSHong Zhang #endif
588fcd7ac73SHong Zhang   PetscFunctionReturn(0);
589fcd7ac73SHong Zhang }
59093dfae19SHong Zhang 
59193dfae19SHong Zhang #undef __FUNCT__
59293dfae19SHong Zhang #define __FUNCT__ "MatFDColoringCreate_MPIXAIJ"
59393dfae19SHong Zhang PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
59493dfae19SHong Zhang {
595531e53bdSHong Zhang   PetscErrorCode ierr;
596a8971b87SHong Zhang   PetscInt       bs,nis=iscoloring->n,m=mat->rmap->n;
597531e53bdSHong Zhang   PetscBool      isBAIJ;
598531e53bdSHong Zhang 
59993dfae19SHong Zhang   PetscFunctionBegin;
600531e53bdSHong Zhang   /* set default brows and bcols for speedup inserting the dense matrix into sparse Jacobian;
601531e53bdSHong Zhang    bcols is chosen s.t. dy-array takes 50% of memory space as mat */
602531e53bdSHong Zhang   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
603531e53bdSHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr);
604a8971b87SHong Zhang   if (isBAIJ) {
605a8971b87SHong Zhang     c->brows = m;
606a8971b87SHong Zhang     c->bcols = 1;
607531e53bdSHong Zhang   } else { /* mpiaij matrix */
608531e53bdSHong Zhang     /* bcols is chosen s.t. dy-array takes 50% of local memory space as mat */
609531e53bdSHong Zhang     Mat_MPIAIJ *aij=(Mat_MPIAIJ*)mat->data;
610531e53bdSHong Zhang     Mat_SeqAIJ *spA,*spB;
611531e53bdSHong Zhang     Mat        A,B;
612a8971b87SHong Zhang     PetscInt   nz,brows,bcols;
613531e53bdSHong Zhang     PetscReal  mem;
614531e53bdSHong Zhang 
615531e53bdSHong Zhang     bs    = 1; /* only bs=1 is supported for MPIAIJ matrix */
616531e53bdSHong Zhang 
617531e53bdSHong Zhang     A = aij->A;  spA = (Mat_SeqAIJ*)A->data;
618531e53bdSHong Zhang     B = aij->B;  spB = (Mat_SeqAIJ*)B->data;
619531e53bdSHong Zhang     nz = spA->nz + spB->nz; /* total local nonzero entries of mat */
620531e53bdSHong Zhang     mem = nz*(sizeof(PetscScalar) + sizeof(PetscInt)) + 3*m*sizeof(PetscInt);
621531e53bdSHong Zhang     bcols = (PetscInt)(0.5*mem /(m*sizeof(PetscScalar)));
622531e53bdSHong Zhang     brows = 1000/bcols;
623531e53bdSHong Zhang     if (bcols > nis) bcols = nis;
624531e53bdSHong Zhang     if (brows == 0 || brows > m) brows = m;
625531e53bdSHong Zhang     c->brows = brows;
626531e53bdSHong Zhang     c->bcols = bcols;
627531e53bdSHong Zhang   }
628a8971b87SHong Zhang 
629a8971b87SHong Zhang   c->M       = mat->rmap->N/bs;         /* set the global rows and columns and local rows */
630a8971b87SHong Zhang   c->N       = mat->cmap->N/bs;
631a8971b87SHong Zhang   c->m       = mat->rmap->n/bs;
632a8971b87SHong Zhang   c->rstart  = mat->rmap->rstart/bs;
633a8971b87SHong Zhang   c->ncolors = nis;
63493dfae19SHong Zhang   PetscFunctionReturn(0);
63593dfae19SHong Zhang }
636