xref: /petsc/src/mat/impls/aij/mpi/fdmpiaij.c (revision e7e920445bafd8671dd3ea579405cd582ea97ef2)
1d4002b98SHong Zhang #include <../src/mat/impls/sell/mpi/mpisell.h>
2c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h>
30d1c53f1SHong Zhang #include <../src/mat/impls/baij/mpi/mpibaij.h>
4af0996ceSBarry Smith #include <petsc/private/isimpl.h>
50d1c53f1SHong Zhang 
6d1e9a80fSBarry Smith PetscErrorCode MatFDColoringApply_BAIJ(Mat J,MatFDColoring coloring,Vec x1,void *sctx)
70d1c53f1SHong Zhang {
80d1c53f1SHong Zhang   PetscErrorCode    (*f)(void*,Vec,Vec,void*)=(PetscErrorCode (*)(void*,Vec,Vec,void*))coloring->f;
90d1c53f1SHong Zhang   PetscErrorCode    ierr;
100d1c53f1SHong Zhang   PetscInt          k,cstart,cend,l,row,col,nz,spidx,i,j;
115edff71fSBarry Smith   PetscScalar       dx=0.0,*w3_array,*dy_i,*dy=coloring->dy;
120d1c53f1SHong Zhang   PetscScalar       *vscale_array;
135edff71fSBarry Smith   const PetscScalar *xx;
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;
200df34763SHong Zhang   MatEntry2         *Jentry2=coloring->matentry2;
210d1c53f1SHong Zhang   const PetscInt    ncolors=coloring->ncolors,*ncolumns=coloring->ncolumns,*nrows=coloring->nrows;
220d1c53f1SHong Zhang   PetscInt          bs=J->rmap->bs;
230d1c53f1SHong Zhang 
240d1c53f1SHong Zhang   PetscFunctionBegin;
25*e7e92044SBarry Smith   ierr = VecPinToCPU(x1,PETSC_TRUE);CHKERRQ(ierr);
260d1c53f1SHong Zhang   /* (1) Set w1 = F(x1) */
270d1c53f1SHong Zhang   if (!coloring->fset) {
28217044c2SLisandro Dalcin     ierr = PetscLogEventBegin(MAT_FDColoringFunction,coloring,0,0,0);CHKERRQ(ierr);
290d1c53f1SHong Zhang     ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
30217044c2SLisandro Dalcin     ierr = PetscLogEventEnd(MAT_FDColoringFunction,coloring,0,0,0);CHKERRQ(ierr);
310d1c53f1SHong Zhang   } else {
320d1c53f1SHong Zhang     coloring->fset = PETSC_FALSE;
330d1c53f1SHong Zhang   }
340d1c53f1SHong Zhang 
350d1c53f1SHong Zhang   /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */
360d1c53f1SHong Zhang   ierr = VecGetLocalSize(x1,&nxloc);CHKERRQ(ierr);
370d1c53f1SHong Zhang   if (coloring->htype[0] == 'w') {
380d1c53f1SHong Zhang     /* vscale = dx is a constant scalar */
390d1c53f1SHong Zhang     ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr);
400d1c53f1SHong Zhang     dx = 1.0/(PetscSqrtReal(1.0 + unorm)*epsilon);
410d1c53f1SHong Zhang   } else {
425edff71fSBarry Smith     ierr = VecGetArrayRead(x1,&xx);CHKERRQ(ierr);
430d1c53f1SHong Zhang     ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr);
440d1c53f1SHong Zhang     for (col=0; col<nxloc; col++) {
450d1c53f1SHong Zhang       dx = xx[col];
460d1c53f1SHong Zhang       if (PetscAbsScalar(dx) < umin) {
470d1c53f1SHong Zhang         if (PetscRealPart(dx) >= 0.0)      dx = umin;
480d1c53f1SHong Zhang         else if (PetscRealPart(dx) < 0.0 ) dx = -umin;
490d1c53f1SHong Zhang       }
500d1c53f1SHong Zhang       dx               *= epsilon;
510d1c53f1SHong Zhang       vscale_array[col] = 1.0/dx;
520d1c53f1SHong Zhang     }
535edff71fSBarry Smith     ierr = VecRestoreArrayRead(x1,&xx);CHKERRQ(ierr);
540d1c53f1SHong Zhang     ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr);
550d1c53f1SHong Zhang   }
56684f2004SHong Zhang   if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') {
570d1c53f1SHong Zhang     ierr = VecGhostUpdateBegin(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
580d1c53f1SHong Zhang     ierr = VecGhostUpdateEnd(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
590d1c53f1SHong Zhang   }
600d1c53f1SHong Zhang 
610d1c53f1SHong Zhang   /* (3) Loop over each color */
620d1c53f1SHong Zhang   if (!coloring->w3) {
630d1c53f1SHong Zhang     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
64*e7e92044SBarry Smith     /* Vec is used instensively in particular piece of scalar CPU code; won't benifit from bouncing back and forth to the GPU */
65*e7e92044SBarry Smith     ierr = VecPinToCPU(coloring->w3,PETSC_TRUE);CHKERRQ(ierr);
660d1c53f1SHong Zhang     ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr);
670d1c53f1SHong Zhang   }
680d1c53f1SHong Zhang   w3 = coloring->w3;
690d1c53f1SHong Zhang 
700d1c53f1SHong Zhang   ierr = VecGetOwnershipRange(x1,&cstart,&cend);CHKERRQ(ierr); /* used by ghosted vscale */
710d1c53f1SHong Zhang   if (vscale) {
720d1c53f1SHong Zhang     ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr);
730d1c53f1SHong Zhang   }
740d1c53f1SHong Zhang   nz   = 0;
750d1c53f1SHong Zhang   for (k=0; k<ncolors; k++) {
760d1c53f1SHong Zhang     coloring->currentcolor = k;
770d1c53f1SHong Zhang 
780d1c53f1SHong Zhang     /*
790d1c53f1SHong Zhang       (3-1) Loop over each column associated with color
800d1c53f1SHong Zhang       adding the perturbation to the vector w3 = x1 + dx.
810d1c53f1SHong Zhang     */
820d1c53f1SHong Zhang     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
830d1c53f1SHong Zhang     dy_i = dy;
840d1c53f1SHong Zhang     for (i=0; i<bs; i++) {     /* Loop over a block of columns */
850d1c53f1SHong Zhang       ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
860d1c53f1SHong Zhang       if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */
870d1c53f1SHong Zhang       if (coloring->htype[0] == 'w') {
880d1c53f1SHong Zhang         for (l=0; l<ncolumns[k]; l++) {
890d1c53f1SHong Zhang           col            = i + bs*coloring->columns[k][l];  /* local column (in global index!) of the matrix we are probing for */
900d1c53f1SHong Zhang           w3_array[col] += 1.0/dx;
91684f2004SHong Zhang           if (i) w3_array[col-1] -= 1.0/dx; /* resume original w3[col-1] */
920d1c53f1SHong Zhang         }
930d1c53f1SHong Zhang       } else { /* htype == 'ds' */
940d1c53f1SHong Zhang         vscale_array -= cstart; /* shift pointer so global index can be used */
950d1c53f1SHong Zhang         for (l=0; l<ncolumns[k]; l++) {
96f8c2866eSHong Zhang           col = i + bs*coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
970d1c53f1SHong Zhang           w3_array[col] += 1.0/vscale_array[col];
98684f2004SHong Zhang           if (i) w3_array[col-1] -=  1.0/vscale_array[col-1]; /* resume original w3[col-1] */
990d1c53f1SHong Zhang         }
1000d1c53f1SHong Zhang         vscale_array += cstart;
1010d1c53f1SHong Zhang       }
1020d1c53f1SHong Zhang       if (ctype == IS_COLORING_GLOBAL) w3_array += cstart;
1030d1c53f1SHong Zhang       ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
1040d1c53f1SHong Zhang 
1050d1c53f1SHong Zhang       /*
1060d1c53f1SHong Zhang        (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
1070d1c53f1SHong Zhang                            w2 = F(x1 + dx) - F(x1)
1080d1c53f1SHong Zhang        */
1090d1c53f1SHong Zhang       ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
1100d1c53f1SHong Zhang       ierr = VecPlaceArray(w2,dy_i);CHKERRQ(ierr); /* place w2 to the array dy_i */
1110d1c53f1SHong Zhang       ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
1120d1c53f1SHong Zhang       ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
1130d1c53f1SHong Zhang       ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
1140d1c53f1SHong Zhang       ierr = VecResetArray(w2);CHKERRQ(ierr);
1150d1c53f1SHong Zhang       dy_i += nxloc; /* points to dy+i*nxloc */
1160d1c53f1SHong Zhang     }
1170d1c53f1SHong Zhang 
1180d1c53f1SHong Zhang     /*
1190d1c53f1SHong Zhang      (3-3) Loop over rows of vector, putting results into Jacobian matrix
1200d1c53f1SHong Zhang     */
1210d1c53f1SHong Zhang     nrows_k = nrows[k];
1220d1c53f1SHong Zhang     if (coloring->htype[0] == 'w') {
1230d1c53f1SHong Zhang       for (l=0; l<nrows_k; l++) {
1240df34763SHong Zhang         row     = bs*Jentry2[nz].row;   /* local row index */
1250df34763SHong Zhang         valaddr = Jentry2[nz++].valaddr;
1260d1c53f1SHong Zhang         spidx   = 0;
1270d1c53f1SHong Zhang         dy_i    = dy;
1280d1c53f1SHong Zhang         for (i=0; i<bs; i++) {   /* column of the block */
1290d1c53f1SHong Zhang           for (j=0; j<bs; j++) { /* row of the block */
1300d1c53f1SHong Zhang             valaddr[spidx++] = dy_i[row+j]*dx;
1310d1c53f1SHong Zhang           }
132f8c2866eSHong Zhang           dy_i += nxloc; /* points to dy+i*nxloc */
1330d1c53f1SHong Zhang         }
1340d1c53f1SHong Zhang       }
1350d1c53f1SHong Zhang     } else { /* htype == 'ds' */
1360d1c53f1SHong Zhang       for (l=0; l<nrows_k; l++) {
1370d1c53f1SHong Zhang         row     = bs*Jentry[nz].row;   /* local row index */
1380d1c53f1SHong Zhang         col     = bs*Jentry[nz].col;   /* local column index */
139684f2004SHong Zhang         valaddr = Jentry[nz++].valaddr;
140f8c2866eSHong Zhang         spidx   = 0;
141f8c2866eSHong Zhang         dy_i    = dy;
142f8c2866eSHong Zhang         for (i=0; i<bs; i++) {   /* column of the block */
143f8c2866eSHong Zhang           for (j=0; j<bs; j++) { /* row of the block */
144f8c2866eSHong Zhang             valaddr[spidx++] = dy_i[row+j]*vscale_array[col+i];
145f8c2866eSHong Zhang           }
146f8c2866eSHong Zhang           dy_i += nxloc; /* points to dy+i*nxloc */
147f8c2866eSHong Zhang         }
1480d1c53f1SHong Zhang       }
1490d1c53f1SHong Zhang     }
1500d1c53f1SHong Zhang   }
1510d1c53f1SHong Zhang   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1520d1c53f1SHong Zhang   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1530d1c53f1SHong Zhang   if (vscale) {
1540d1c53f1SHong Zhang     ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr);
1550d1c53f1SHong Zhang   }
1560d1c53f1SHong Zhang 
1570d1c53f1SHong Zhang   coloring->currentcolor = -1;
158*e7e92044SBarry Smith   ierr = VecPinToCPU(x1,PETSC_FALSE);CHKERRQ(ierr);
1590d1c53f1SHong Zhang   PetscFunctionReturn(0);
1600d1c53f1SHong Zhang }
161a64fbb32SBarry Smith 
162531c7667SBarry Smith /* this is declared PETSC_EXTERN because it is used by MatFDColoringUseDM() which is in the DM library */
163d1e9a80fSBarry Smith PetscErrorCode  MatFDColoringApply_AIJ(Mat J,MatFDColoring coloring,Vec x1,void *sctx)
164fcd7ac73SHong Zhang {
165fcd7ac73SHong Zhang   PetscErrorCode    (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void*))coloring->f;
166fcd7ac73SHong Zhang   PetscErrorCode    ierr;
167a2f2d239SHong Zhang   PetscInt          k,cstart,cend,l,row,col,nz;
1685edff71fSBarry Smith   PetscScalar       dx=0.0,*y,*w3_array;
1695edff71fSBarry Smith   const PetscScalar *xx;
170fcd7ac73SHong Zhang   PetscScalar       *vscale_array;
171fcd7ac73SHong Zhang   PetscReal         epsilon=coloring->error_rel,umin=coloring->umin,unorm;
1728bc97078SHong Zhang   Vec               w1=coloring->w1,w2=coloring->w2,w3,vscale=coloring->vscale;
173fcd7ac73SHong Zhang   void              *fctx=coloring->fctx;
17491e7fa0fSBarry Smith   ISColoringType    ctype=coloring->ctype;
17591e7fa0fSBarry Smith   PetscInt          nxloc,nrows_k;
176573f477fSHong Zhang   MatEntry          *Jentry=coloring->matentry;
1770df34763SHong Zhang   MatEntry2         *Jentry2=coloring->matentry2;
1788bc97078SHong Zhang   const PetscInt    ncolors=coloring->ncolors,*ncolumns=coloring->ncolumns,*nrows=coloring->nrows;
179fcd7ac73SHong Zhang 
180fcd7ac73SHong Zhang   PetscFunctionBegin;
181*e7e92044SBarry Smith   ierr = VecPinToCPU(x1,PETSC_TRUE);CHKERRQ(ierr);
182531c7667SBarry Smith   if ((ctype == IS_COLORING_LOCAL) && (J->ops->fdcoloringapply == MatFDColoringApply_AIJ)) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_SUP,"Must call MatColoringUseDM() with IS_COLORING_LOCAL");
1838bc97078SHong Zhang   /* (1) Set w1 = F(x1) */
184fcd7ac73SHong Zhang   if (!coloring->fset) {
185fcd7ac73SHong Zhang     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
186f6af9589SHong Zhang     ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
187fcd7ac73SHong Zhang     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
188fcd7ac73SHong Zhang   } else {
189fcd7ac73SHong Zhang     coloring->fset = PETSC_FALSE;
190fcd7ac73SHong Zhang   }
191fcd7ac73SHong Zhang 
1928bc97078SHong Zhang   /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */
193f6af9589SHong Zhang   if (coloring->htype[0] == 'w') {
194684f2004SHong Zhang     /* vscale = 1./dx is a constant scalar */
195f6af9589SHong Zhang     ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr);
196c53567a0SHong Zhang     dx = 1.0/(PetscSqrtReal(1.0 + unorm)*epsilon);
19770e7395fSHong Zhang   } else {
19874d3cef9SHong Zhang     ierr = VecGetLocalSize(x1,&nxloc);CHKERRQ(ierr);
1995edff71fSBarry Smith     ierr = VecGetArrayRead(x1,&xx);CHKERRQ(ierr);
2008bc97078SHong Zhang     ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr);
20174d3cef9SHong Zhang     for (col=0; col<nxloc; col++) {
202fcd7ac73SHong Zhang       dx = xx[col];
20374d3cef9SHong Zhang       if (PetscAbsScalar(dx) < umin) {
20474d3cef9SHong Zhang         if (PetscRealPart(dx) >= 0.0)      dx = umin;
20574d3cef9SHong Zhang         else if (PetscRealPart(dx) < 0.0 ) dx = -umin;
20674d3cef9SHong Zhang       }
207fcd7ac73SHong Zhang       dx               *= epsilon;
20874d3cef9SHong Zhang       vscale_array[col] = 1.0/dx;
209f6af9589SHong Zhang     }
2105edff71fSBarry Smith     ierr = VecRestoreArrayRead(x1,&xx);CHKERRQ(ierr);
2118bc97078SHong Zhang     ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr);
21270e7395fSHong Zhang   }
213684f2004SHong Zhang   if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') {
2148bc97078SHong Zhang     ierr = VecGhostUpdateBegin(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2158bc97078SHong Zhang     ierr = VecGhostUpdateEnd(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
216fcd7ac73SHong Zhang   }
217fcd7ac73SHong Zhang 
2188bc97078SHong Zhang   /* (3) Loop over each color */
2198bc97078SHong Zhang   if (!coloring->w3) {
2208bc97078SHong Zhang     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
2218bc97078SHong Zhang     ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr);
2228bc97078SHong Zhang   }
2238bc97078SHong Zhang   w3 = coloring->w3;
224fcd7ac73SHong Zhang 
2258bc97078SHong Zhang   ierr = VecGetOwnershipRange(x1,&cstart,&cend);CHKERRQ(ierr); /* used by ghosted vscale */
2269e917edbSHong Zhang   if (vscale) {
2278bc97078SHong Zhang     ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr);
2289e917edbSHong Zhang   }
2298bc97078SHong Zhang   nz   = 0;
230e2c857f8SHong Zhang 
231b3e1f37bSHong Zhang   if (coloring->bcols > 1) { /* use blocked insertion of Jentry */
232b3e1f37bSHong Zhang     PetscInt    i,m=J->rmap->n,nbcols,bcols=coloring->bcols;
233e2c857f8SHong Zhang     PetscScalar *dy=coloring->dy,*dy_k;
234e2c857f8SHong Zhang 
235e2c857f8SHong Zhang     nbcols = 0;
236e2c857f8SHong Zhang     for (k=0; k<ncolors; k+=bcols) {
237e2c857f8SHong Zhang 
238e2c857f8SHong Zhang       /*
239e2c857f8SHong Zhang        (3-1) Loop over each column associated with color
240e2c857f8SHong Zhang        adding the perturbation to the vector w3 = x1 + dx.
241e2c857f8SHong Zhang        */
242e2c857f8SHong Zhang 
243e2c857f8SHong Zhang       dy_k = dy;
244e2c857f8SHong Zhang       if (k + bcols > ncolors) bcols = ncolors - k;
245e2c857f8SHong Zhang       for (i=0; i<bcols; i++) {
246ceef8a49SBarry Smith         coloring->currentcolor = k+i;
247e2c857f8SHong Zhang 
248e2c857f8SHong Zhang         ierr = VecCopy(x1,w3);CHKERRQ(ierr);
249e2c857f8SHong Zhang         ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
250e2c857f8SHong Zhang         if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */
251e2c857f8SHong Zhang         if (coloring->htype[0] == 'w') {
252e2c857f8SHong Zhang           for (l=0; l<ncolumns[k+i]; l++) {
253e2c857f8SHong Zhang             col = coloring->columns[k+i][l]; /* local column (in global index!) of the matrix we are probing for */
254e2c857f8SHong Zhang             w3_array[col] += 1.0/dx;
255e2c857f8SHong Zhang           }
256e2c857f8SHong Zhang         } else { /* htype == 'ds' */
257e2c857f8SHong Zhang           vscale_array -= cstart; /* shift pointer so global index can be used */
258e2c857f8SHong Zhang           for (l=0; l<ncolumns[k+i]; l++) {
259e2c857f8SHong Zhang             col = coloring->columns[k+i][l]; /* local column (in global index!) of the matrix we are probing for */
260e2c857f8SHong Zhang             w3_array[col] += 1.0/vscale_array[col];
261e2c857f8SHong Zhang           }
262e2c857f8SHong Zhang           vscale_array += cstart;
263e2c857f8SHong Zhang         }
264e2c857f8SHong Zhang         if (ctype == IS_COLORING_GLOBAL) w3_array += cstart;
265e2c857f8SHong Zhang         ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
266e2c857f8SHong Zhang 
267e2c857f8SHong Zhang         /*
268e2c857f8SHong Zhang          (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
269e2c857f8SHong Zhang                            w2 = F(x1 + dx) - F(x1)
270e2c857f8SHong Zhang          */
271e2c857f8SHong Zhang         ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
272e2c857f8SHong Zhang         ierr = VecPlaceArray(w2,dy_k);CHKERRQ(ierr); /* place w2 to the array dy_i */
273e2c857f8SHong Zhang         ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
274e2c857f8SHong Zhang         ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
275e2c857f8SHong Zhang         ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
276e2c857f8SHong Zhang         ierr = VecResetArray(w2);CHKERRQ(ierr);
277e2c857f8SHong Zhang         dy_k += m; /* points to dy+i*nxloc */
278e2c857f8SHong Zhang       }
279e2c857f8SHong Zhang 
280e2c857f8SHong Zhang       /*
281e2c857f8SHong Zhang        (3-3) Loop over block rows of vector, putting results into Jacobian matrix
282e2c857f8SHong Zhang        */
283d880da65SHong Zhang       nrows_k = nrows[nbcols++];
284e2c857f8SHong Zhang       ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
285d880da65SHong Zhang 
286e2c857f8SHong Zhang       if (coloring->htype[0] == 'w') {
287e2c857f8SHong Zhang         for (l=0; l<nrows_k; l++) {
2880df34763SHong Zhang           row                      = Jentry2[nz].row;   /* local row index */
2890df34763SHong Zhang           *(Jentry2[nz++].valaddr) = dy[row]*dx;
290e2c857f8SHong Zhang         }
291e2c857f8SHong Zhang       } else { /* htype == 'ds' */
292e2c857f8SHong Zhang         for (l=0; l<nrows_k; l++) {
293e2c857f8SHong Zhang           row                   = Jentry[nz].row;   /* local row index */
294d880da65SHong Zhang           *(Jentry[nz].valaddr) = dy[row]*vscale_array[Jentry[nz].col];
295e2c857f8SHong Zhang           nz++;
296e2c857f8SHong Zhang         }
297e2c857f8SHong Zhang       }
298e2c857f8SHong Zhang       ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
299e2c857f8SHong Zhang     }
300b3e1f37bSHong Zhang   } else { /* bcols == 1 */
3018bc97078SHong Zhang     for (k=0; k<ncolors; k++) {
3028bc97078SHong Zhang       coloring->currentcolor = k;
3039e917edbSHong Zhang 
304fcd7ac73SHong Zhang       /*
3058bc97078SHong Zhang        (3-1) Loop over each column associated with color
306f6af9589SHong Zhang        adding the perturbation to the vector w3 = x1 + dx.
307fcd7ac73SHong Zhang        */
308f6af9589SHong Zhang       ierr = VecCopy(x1,w3);CHKERRQ(ierr);
309f6af9589SHong Zhang       ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
310a2f2d239SHong Zhang       if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */
311b039d6c4SHong Zhang       if (coloring->htype[0] == 'w') {
3128bc97078SHong Zhang         for (l=0; l<ncolumns[k]; l++) {
313f6af9589SHong Zhang           col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
3149e917edbSHong Zhang           w3_array[col] += 1.0/dx;
3159e917edbSHong Zhang         }
316b039d6c4SHong Zhang       } else { /* htype == 'ds' */
317a2f2d239SHong Zhang         vscale_array -= cstart; /* shift pointer so global index can be used */
318b039d6c4SHong Zhang         for (l=0; l<ncolumns[k]; l++) {
319b039d6c4SHong Zhang           col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
320a2f2d239SHong Zhang           w3_array[col] += 1.0/vscale_array[col];
32170e7395fSHong Zhang         }
322a2f2d239SHong Zhang         vscale_array += cstart;
323fcd7ac73SHong Zhang       }
324a2f2d239SHong Zhang       if (ctype == IS_COLORING_GLOBAL) w3_array += cstart;
325fcd7ac73SHong Zhang       ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
3269e917edbSHong Zhang 
327fcd7ac73SHong Zhang       /*
3288bc97078SHong Zhang        (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
329fcd7ac73SHong Zhang                            w2 = F(x1 + dx) - F(x1)
330fcd7ac73SHong Zhang        */
331fcd7ac73SHong Zhang       ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
332fcd7ac73SHong Zhang       ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
333fcd7ac73SHong Zhang       ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
334fcd7ac73SHong Zhang       ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
3359e917edbSHong Zhang 
3368bc97078SHong Zhang       /*
3378bc97078SHong Zhang        (3-3) Loop over rows of vector, putting results into Jacobian matrix
3388bc97078SHong Zhang        */
339c53567a0SHong Zhang       nrows_k = nrows[k];
340fcd7ac73SHong Zhang       ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
341b039d6c4SHong Zhang       if (coloring->htype[0] == 'w') {
342c53567a0SHong Zhang         for (l=0; l<nrows_k; l++) {
3430df34763SHong Zhang           row                      = Jentry2[nz].row;   /* local row index */
3440df34763SHong Zhang           *(Jentry2[nz++].valaddr) = y[row]*dx;
3459e917edbSHong Zhang         }
346b039d6c4SHong Zhang       } else { /* htype == 'ds' */
347c53567a0SHong Zhang         for (l=0; l<nrows_k; l++) {
348b039d6c4SHong Zhang           row                   = Jentry[nz].row;   /* local row index */
349b039d6c4SHong Zhang           *(Jentry[nz].valaddr) = y[row]*vscale_array[Jentry[nz].col];
350573f477fSHong Zhang           nz++;
351fcd7ac73SHong Zhang         }
352b039d6c4SHong Zhang       }
353fcd7ac73SHong Zhang       ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
3548bc97078SHong Zhang     }
355b3e1f37bSHong Zhang   }
356b3e1f37bSHong Zhang 
357f5aae955SHong Zhang   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
358f5aae955SHong Zhang   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3599e917edbSHong Zhang   if (vscale) {
3608bc97078SHong Zhang     ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr);
3619e917edbSHong Zhang   }
362fcd7ac73SHong Zhang   coloring->currentcolor = -1;
363*e7e92044SBarry Smith   ierr = VecPinToCPU(x1,PETSC_FALSE);CHKERRQ(ierr);
364fcd7ac73SHong Zhang   PetscFunctionReturn(0);
365fcd7ac73SHong Zhang }
366fcd7ac73SHong Zhang 
367f86b9fbaSHong Zhang PetscErrorCode MatFDColoringSetUp_MPIXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
368fcd7ac73SHong Zhang {
369fcd7ac73SHong Zhang   PetscErrorCode         ierr;
370fcd7ac73SHong Zhang   PetscMPIInt            size,*ncolsonproc,*disp,nn;
371e3225e9dSHong Zhang   PetscInt               i,n,nrows,nrows_i,j,k,m,ncols,col,*rowhit,cstart,cend,colb;
37272c15787SHong Zhang   const PetscInt         *is,*A_ci,*A_cj,*B_ci,*B_cj,*row=NULL,*ltog=NULL;
373fcd7ac73SHong Zhang   PetscInt               nis=iscoloring->n,nctot,*cols;
374fcd7ac73SHong Zhang   IS                     *isa;
375fcd7ac73SHong Zhang   ISLocalToGlobalMapping map=mat->cmap->mapping;
376a8971b87SHong Zhang   PetscInt               ctype=c->ctype,*spidxA,*spidxB,nz,bs,bs2,spidx;
3770d1c53f1SHong Zhang   Mat                    A,B;
378e3225e9dSHong Zhang   PetscScalar            *A_val,*B_val,**valaddrhit;
379573f477fSHong Zhang   MatEntry               *Jentry;
3800df34763SHong Zhang   MatEntry2              *Jentry2;
381d4002b98SHong Zhang   PetscBool              isBAIJ,isSELL;
382531e53bdSHong Zhang   PetscInt               bcols=c->bcols;
3830d1c53f1SHong Zhang #if defined(PETSC_USE_CTABLE)
3840d1c53f1SHong Zhang   PetscTable             colmap=NULL;
3850d1c53f1SHong Zhang #else
3860d1c53f1SHong Zhang   PetscInt               *colmap=NULL;     /* local col number of off-diag col */
3870d1c53f1SHong Zhang #endif
388fcd7ac73SHong Zhang 
389fcd7ac73SHong Zhang   PetscFunctionBegin;
3905bdb020cSBarry Smith   if (ctype == IS_COLORING_LOCAL) {
39172c15787SHong 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");
39272c15787SHong Zhang     ierr = ISLocalToGlobalMappingGetIndices(map,&ltog);CHKERRQ(ierr);
39372c15787SHong Zhang   }
394fcd7ac73SHong Zhang 
3950d1c53f1SHong Zhang   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
3960d1c53f1SHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr);
397d4002b98SHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPISELL,&isSELL);CHKERRQ(ierr);
398e3225e9dSHong Zhang   if (isBAIJ) {
399195687c5SHong Zhang     Mat_MPIBAIJ *baij=(Mat_MPIBAIJ*)mat->data;
4006c145f34SHong Zhang     Mat_SeqBAIJ *spA,*spB;
401195687c5SHong Zhang     A = baij->A;  spA = (Mat_SeqBAIJ*)A->data; A_val = spA->a;
402195687c5SHong Zhang     B = baij->B;  spB = (Mat_SeqBAIJ*)B->data; B_val = spB->a;
4030d1c53f1SHong Zhang     nz = spA->nz + spB->nz; /* total nonzero entries of mat */
404195687c5SHong Zhang     if (!baij->colmap) {
4050d1c53f1SHong Zhang       ierr = MatCreateColmap_MPIBAIJ_Private(mat);CHKERRQ(ierr);
4060d1c53f1SHong Zhang     }
407097e26fcSJed Brown     colmap = baij->colmap;
4080d1c53f1SHong Zhang     ierr = MatGetColumnIJ_SeqBAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
4090d1c53f1SHong Zhang     ierr = MatGetColumnIJ_SeqBAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
410195687c5SHong Zhang 
411195687c5SHong Zhang     if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') {  /* create vscale for storing dx */
412195687c5SHong Zhang       PetscInt    *garray;
413785e854fSJed Brown       ierr = PetscMalloc1(B->cmap->n,&garray);CHKERRQ(ierr);
414195687c5SHong Zhang       for (i=0; i<baij->B->cmap->n/bs; i++) {
415195687c5SHong Zhang         for (j=0; j<bs; j++) {
416195687c5SHong Zhang           garray[i*bs+j] = bs*baij->garray[i]+j;
417195687c5SHong Zhang         }
418195687c5SHong Zhang       }
419195687c5SHong Zhang       ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->cmap->n,PETSC_DETERMINE,B->cmap->n,garray,&c->vscale);CHKERRQ(ierr);
420*e7e92044SBarry Smith       ierr = VecPinToCPU(c->vscale,PETSC_TRUE);CHKERRQ(ierr);
421195687c5SHong Zhang       ierr = PetscFree(garray);CHKERRQ(ierr);
422195687c5SHong Zhang     }
423d4002b98SHong Zhang   } else if (isSELL) {
424d4002b98SHong Zhang     Mat_MPISELL *sell=(Mat_MPISELL*)mat->data;
425d4002b98SHong Zhang     Mat_SeqSELL *spA,*spB;
426d4002b98SHong Zhang     A = sell->A;  spA = (Mat_SeqSELL*)A->data; A_val = spA->val;
427d4002b98SHong Zhang     B = sell->B;  spB = (Mat_SeqSELL*)B->data; B_val = spB->val;
428f4b713efSHong Zhang     nz = spA->nz + spB->nz; /* total nonzero entries of mat */
429d4002b98SHong Zhang     if (!sell->colmap) {
430f4b713efSHong Zhang       /* Allow access to data structures of local part of matrix
431f4b713efSHong Zhang        - creates aij->colmap which maps global column number to local number in part B */
432d4002b98SHong Zhang       ierr = MatCreateColmap_MPISELL_Private(mat);CHKERRQ(ierr);
433f4b713efSHong Zhang     }
434d4002b98SHong Zhang     colmap = sell->colmap;
435d4002b98SHong Zhang     ierr = MatGetColumnIJ_SeqSELL_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
436d4002b98SHong Zhang     ierr = MatGetColumnIJ_SeqSELL_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
437f4b713efSHong Zhang 
438f4b713efSHong Zhang     bs = 1; /* only bs=1 is supported for non MPIBAIJ matrix */
439f4b713efSHong Zhang 
440f4b713efSHong Zhang     if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') { /* create vscale for storing dx */
441d4002b98SHong Zhang       ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->cmap->n,PETSC_DETERMINE,B->cmap->n,sell->garray,&c->vscale);CHKERRQ(ierr);
442*e7e92044SBarry Smith       ierr = VecPinToCPU(c->vscale,PETSC_TRUE);CHKERRQ(ierr);
443f4b713efSHong Zhang     }
444e3225e9dSHong Zhang   } else {
445e3225e9dSHong Zhang     Mat_MPIAIJ *aij=(Mat_MPIAIJ*)mat->data;
446e3225e9dSHong Zhang     Mat_SeqAIJ *spA,*spB;
447e3225e9dSHong Zhang     A = aij->A;  spA = (Mat_SeqAIJ*)A->data; A_val = spA->a;
448e3225e9dSHong Zhang     B = aij->B;  spB = (Mat_SeqAIJ*)B->data; B_val = spB->a;
449e3225e9dSHong Zhang     nz = spA->nz + spB->nz; /* total nonzero entries of mat */
450e3225e9dSHong Zhang     if (!aij->colmap) {
451e3225e9dSHong Zhang       /* Allow access to data structures of local part of matrix
452e3225e9dSHong Zhang        - creates aij->colmap which maps global column number to local number in part B */
453e3225e9dSHong Zhang       ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
454e3225e9dSHong Zhang     }
455097e26fcSJed Brown     colmap = aij->colmap;
456e3225e9dSHong Zhang     ierr = MatGetColumnIJ_SeqAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
457e3225e9dSHong Zhang     ierr = MatGetColumnIJ_SeqAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
458e3225e9dSHong Zhang 
459e3225e9dSHong Zhang     bs = 1; /* only bs=1 is supported for non MPIBAIJ matrix */
460195687c5SHong Zhang 
461195687c5SHong Zhang     if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') { /* create vscale for storing dx */
462195687c5SHong Zhang       ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->cmap->n,PETSC_DETERMINE,B->cmap->n,aij->garray,&c->vscale);CHKERRQ(ierr);
463*e7e92044SBarry Smith       ierr = VecPinToCPU(c->vscale,PETSC_TRUE);CHKERRQ(ierr);
464195687c5SHong Zhang     }
4650d1c53f1SHong Zhang   }
4660d1c53f1SHong Zhang 
4670d1c53f1SHong Zhang   m         = mat->rmap->n/bs;
4680d1c53f1SHong Zhang   cstart    = mat->cmap->rstart/bs;
4690d1c53f1SHong Zhang   cend      = mat->cmap->rend/bs;
470fcd7ac73SHong Zhang 
471785e854fSJed Brown   ierr       = PetscMalloc1(nis,&c->ncolumns);CHKERRQ(ierr);
472785e854fSJed Brown   ierr       = PetscMalloc1(nis,&c->columns);CHKERRQ(ierr);
4735bdb020cSBarry Smith   ierr       = PetscCalloc1(nis,&c->nrows);CHKERRQ(ierr);
4748bc97078SHong Zhang   ierr       = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr);
475fcd7ac73SHong Zhang 
4760df34763SHong Zhang   if (c->htype[0] == 'd') {
477785e854fSJed Brown     ierr       = PetscMalloc1(nz,&Jentry);CHKERRQ(ierr);
4788bc97078SHong Zhang     ierr       = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr);
4798bc97078SHong Zhang     c->matentry = Jentry;
4800df34763SHong Zhang   } else if (c->htype[0] == 'w') {
481785e854fSJed Brown     ierr       = PetscMalloc1(nz,&Jentry2);CHKERRQ(ierr);
4820df34763SHong Zhang     ierr       = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry2));CHKERRQ(ierr);
4830df34763SHong Zhang     c->matentry2 = Jentry2;
4840df34763SHong Zhang   } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"htype is not supported");
485a774a6f1SHong Zhang 
486dcca6d9dSJed Brown   ierr = PetscMalloc2(m+1,&rowhit,m+1,&valaddrhit);CHKERRQ(ierr);
487fcd7ac73SHong Zhang   nz = 0;
48872c15787SHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
489d3825b63SHong Zhang   for (i=0; i<nis; i++) { /* for each local color */
490fcd7ac73SHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
491fcd7ac73SHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
492fcd7ac73SHong Zhang 
493fcd7ac73SHong Zhang     c->ncolumns[i] = n; /* local number of columns of this color on this process */
494fcd7ac73SHong Zhang     if (n) {
495785e854fSJed Brown       ierr = PetscMalloc1(n,&c->columns[i]);CHKERRQ(ierr);
496fcd7ac73SHong Zhang       ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr);
497fcd7ac73SHong Zhang       ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr);
498fcd7ac73SHong Zhang     } else {
499fcd7ac73SHong Zhang       c->columns[i] = 0;
500fcd7ac73SHong Zhang     }
501fcd7ac73SHong Zhang 
502fcd7ac73SHong Zhang     if (ctype == IS_COLORING_GLOBAL) {
503d3825b63SHong Zhang       /* Determine nctot, the total (parallel) number of columns of this color */
504fcd7ac73SHong Zhang       ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr);
505dcca6d9dSJed Brown       ierr = PetscMalloc2(size,&ncolsonproc,size,&disp);CHKERRQ(ierr);
506fcd7ac73SHong Zhang 
507d3825b63SHong Zhang       /* ncolsonproc[j]: local ncolumns on proc[j] of this color */
508fcd7ac73SHong Zhang       ierr  = PetscMPIIntCast(n,&nn);CHKERRQ(ierr);
509fcd7ac73SHong Zhang       ierr  = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
510fcd7ac73SHong Zhang       nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j];
511fcd7ac73SHong Zhang       if (!nctot) {
512fcd7ac73SHong Zhang         ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr);
513fcd7ac73SHong Zhang       }
514fcd7ac73SHong Zhang 
515fcd7ac73SHong Zhang       disp[0] = 0;
516fcd7ac73SHong Zhang       for (j=1; j<size; j++) {
517fcd7ac73SHong Zhang         disp[j] = disp[j-1] + ncolsonproc[j-1];
518fcd7ac73SHong Zhang       }
519fcd7ac73SHong Zhang 
520d3825b63SHong Zhang       /* Get cols, the complete list of columns for this color on each process */
521854ce69bSBarry Smith       ierr = PetscMalloc1(nctot+1,&cols);CHKERRQ(ierr);
522fcd7ac73SHong Zhang       ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
523fcd7ac73SHong Zhang       ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr);
5245bdb020cSBarry Smith     } else if (ctype == IS_COLORING_LOCAL) {
525fcd7ac73SHong Zhang       /* Determine local number of columns of this color on this process, including ghost points */
526fcd7ac73SHong Zhang       nctot = n;
527854ce69bSBarry Smith       ierr  = PetscMalloc1(nctot+1,&cols);CHKERRQ(ierr);
528fcd7ac73SHong Zhang       ierr  = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr);
529fcd7ac73SHong Zhang     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type");
530fcd7ac73SHong Zhang 
5311b97d346SHong Zhang     /* Mark all rows affect by these columns */
532d3825b63SHong Zhang     ierr    = PetscMemzero(rowhit,m*sizeof(PetscInt));CHKERRQ(ierr);
5330d1c53f1SHong Zhang     bs2     = bs*bs;
5344b2e90caSHong Zhang     nrows_i = 0;
5351b97d346SHong Zhang     for (j=0; j<nctot; j++) { /* loop over columns*/
5365bdb020cSBarry Smith       if (ctype == IS_COLORING_LOCAL) {
537fcd7ac73SHong Zhang         col = ltog[cols[j]];
538fcd7ac73SHong Zhang       } else {
539fcd7ac73SHong Zhang         col = cols[j];
540fcd7ac73SHong Zhang       }
541e3225e9dSHong Zhang       if (col >= cstart && col < cend) { /* column is in A, diagonal block of mat */
542d3825b63SHong Zhang         row      = A_cj + A_ci[col-cstart];
543d3825b63SHong Zhang         nrows    = A_ci[col-cstart+1] - A_ci[col-cstart];
5444b2e90caSHong Zhang         nrows_i += nrows;
545fcd7ac73SHong Zhang         /* loop over columns of A marking them in rowhit */
546d3825b63SHong Zhang         for (k=0; k<nrows; k++) {
547d3825b63SHong Zhang           /* set valaddrhit for part A */
548e3225e9dSHong Zhang           spidx            = bs2*spidxA[A_ci[col-cstart] + k];
549e3225e9dSHong Zhang           valaddrhit[*row] = &A_val[spidx];
550a774a6f1SHong Zhang           rowhit[*row++]   = col - cstart + 1; /* local column index */
551fcd7ac73SHong Zhang         }
552e3225e9dSHong Zhang       } else { /* column is in B, off-diagonal block of mat */
553fcd7ac73SHong Zhang #if defined(PETSC_USE_CTABLE)
5540d1c53f1SHong Zhang         ierr = PetscTableFind(colmap,col+1,&colb);CHKERRQ(ierr);
555fcd7ac73SHong Zhang         colb--;
556fcd7ac73SHong Zhang #else
5570d1c53f1SHong Zhang         colb = colmap[col] - 1; /* local column index */
558fcd7ac73SHong Zhang #endif
559fcd7ac73SHong Zhang         if (colb == -1) {
560d3825b63SHong Zhang           nrows = 0;
561fcd7ac73SHong Zhang         } else {
5620d1c53f1SHong Zhang           colb  = colb/bs;
563d3825b63SHong Zhang           row   = B_cj + B_ci[colb];
564d3825b63SHong Zhang           nrows = B_ci[colb+1] - B_ci[colb];
565fcd7ac73SHong Zhang         }
5664b2e90caSHong Zhang         nrows_i += nrows;
567fcd7ac73SHong Zhang         /* loop over columns of B marking them in rowhit */
568d3825b63SHong Zhang         for (k=0; k<nrows; k++) {
569d3825b63SHong Zhang           /* set valaddrhit for part B */
570e3225e9dSHong Zhang           spidx            = bs2*spidxB[B_ci[colb] + k];
571e3225e9dSHong Zhang           valaddrhit[*row] = &B_val[spidx];
57270e7395fSHong Zhang           rowhit[*row++]   = colb + 1 + cend - cstart; /* local column index */
573fcd7ac73SHong Zhang         }
574fcd7ac73SHong Zhang       }
575fcd7ac73SHong Zhang     }
5764b2e90caSHong Zhang     c->nrows[i] = nrows_i;
5778bc97078SHong Zhang 
5780df34763SHong Zhang     if (c->htype[0] == 'd') {
579d3825b63SHong Zhang       for (j=0; j<m; j++) {
580fcd7ac73SHong Zhang         if (rowhit[j]) {
581573f477fSHong Zhang           Jentry[nz].row     = j;              /* local row index */
582573f477fSHong Zhang           Jentry[nz].col     = rowhit[j] - 1;  /* local column index */
583573f477fSHong Zhang           Jentry[nz].valaddr = valaddrhit[j];  /* address of mat value for this entry */
584573f477fSHong Zhang           nz++;
585fcd7ac73SHong Zhang         }
586fcd7ac73SHong Zhang       }
5870df34763SHong Zhang     } else { /* c->htype == 'wp' */
5880df34763SHong Zhang       for (j=0; j<m; j++) {
5890df34763SHong Zhang         if (rowhit[j]) {
5900df34763SHong Zhang           Jentry2[nz].row     = j;              /* local row index */
5910df34763SHong Zhang           Jentry2[nz].valaddr = valaddrhit[j];  /* address of mat value for this entry */
5920df34763SHong Zhang           nz++;
5930df34763SHong Zhang         }
5940df34763SHong Zhang       }
5950df34763SHong Zhang     }
596fcd7ac73SHong Zhang     ierr = PetscFree(cols);CHKERRQ(ierr);
597fcd7ac73SHong Zhang   }
598fcd7ac73SHong Zhang 
599a8971b87SHong Zhang   if (bcols > 1) { /* reorder Jentry for faster MatFDColoringApply() */
600a8971b87SHong Zhang     ierr = MatFDColoringSetUpBlocked_AIJ_Private(mat,c,nz);CHKERRQ(ierr);
601531e53bdSHong Zhang   }
602531e53bdSHong Zhang 
6030d1c53f1SHong Zhang   if (isBAIJ) {
6040d1c53f1SHong Zhang     ierr = MatRestoreColumnIJ_SeqBAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
6050d1c53f1SHong Zhang     ierr = MatRestoreColumnIJ_SeqBAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
606785e854fSJed Brown     ierr = PetscMalloc1(bs*mat->rmap->n,&c->dy);CHKERRQ(ierr);
607d4002b98SHong Zhang   } else if (isSELL) {
608d4002b98SHong Zhang     ierr = MatRestoreColumnIJ_SeqSELL_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
609d4002b98SHong Zhang     ierr = MatRestoreColumnIJ_SeqSELL_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
6100d1c53f1SHong Zhang   }else {
6110d1c53f1SHong Zhang     ierr = MatRestoreColumnIJ_SeqAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
6120d1c53f1SHong Zhang     ierr = MatRestoreColumnIJ_SeqAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
6130d1c53f1SHong Zhang   }
6140d1c53f1SHong Zhang 
615a8971b87SHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
616a8971b87SHong Zhang   ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr);
617a8971b87SHong Zhang 
6185bdb020cSBarry Smith   if (ctype == IS_COLORING_LOCAL) {
61972c15787SHong Zhang     ierr = ISLocalToGlobalMappingRestoreIndices(map,&ltog);CHKERRQ(ierr);
62072c15787SHong Zhang   }
621b3e1f37bSHong Zhang   ierr = PetscInfo3(c,"ncolors %D, brows %D and bcols %D are used.\n",c->ncolors,c->brows,c->bcols);CHKERRQ(ierr);
622fcd7ac73SHong Zhang   PetscFunctionReturn(0);
623fcd7ac73SHong Zhang }
62493dfae19SHong Zhang 
62593dfae19SHong Zhang PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
62693dfae19SHong Zhang {
627531e53bdSHong Zhang   PetscErrorCode ierr;
628a8971b87SHong Zhang   PetscInt       bs,nis=iscoloring->n,m=mat->rmap->n;
629d4002b98SHong Zhang   PetscBool      isBAIJ,isSELL;
630531e53bdSHong Zhang 
63193dfae19SHong Zhang   PetscFunctionBegin;
632531e53bdSHong Zhang   /* set default brows and bcols for speedup inserting the dense matrix into sparse Jacobian;
633531e53bdSHong Zhang    bcols is chosen s.t. dy-array takes 50% of memory space as mat */
634531e53bdSHong Zhang   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
635531e53bdSHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr);
636d4002b98SHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPISELL,&isSELL);CHKERRQ(ierr);
637b8dfa574SBarry Smith   if (isBAIJ || m == 0) {
638a8971b87SHong Zhang     c->brows = m;
639a8971b87SHong Zhang     c->bcols = 1;
640d4002b98SHong Zhang   } else if (isSELL) {
641f4b713efSHong Zhang     /* bcols is chosen s.t. dy-array takes 50% of local memory space as mat */
642d4002b98SHong Zhang     Mat_MPISELL *sell=(Mat_MPISELL*)mat->data;
643d4002b98SHong Zhang     Mat_SeqSELL *spA,*spB;
644f4b713efSHong Zhang     Mat        A,B;
645f4b713efSHong Zhang     PetscInt   nz,brows,bcols;
646f4b713efSHong Zhang     PetscReal  mem;
647f4b713efSHong Zhang 
648d4002b98SHong Zhang     bs    = 1; /* only bs=1 is supported for MPISELL matrix */
649f4b713efSHong Zhang 
650d4002b98SHong Zhang     A = sell->A;  spA = (Mat_SeqSELL*)A->data;
651d4002b98SHong Zhang     B = sell->B;  spB = (Mat_SeqSELL*)B->data;
652f4b713efSHong Zhang     nz = spA->nz + spB->nz; /* total local nonzero entries of mat */
653f4b713efSHong Zhang     mem = nz*(sizeof(PetscScalar) + sizeof(PetscInt)) + 3*m*sizeof(PetscInt);
654f4b713efSHong Zhang     bcols = (PetscInt)(0.5*mem /(m*sizeof(PetscScalar)));
655f4b713efSHong Zhang     brows = 1000/bcols;
656f4b713efSHong Zhang     if (bcols > nis) bcols = nis;
657f4b713efSHong Zhang     if (brows == 0 || brows > m) brows = m;
658f4b713efSHong Zhang     c->brows = brows;
659f4b713efSHong Zhang     c->bcols = bcols;
660531e53bdSHong Zhang   } else { /* mpiaij matrix */
661531e53bdSHong Zhang     /* bcols is chosen s.t. dy-array takes 50% of local memory space as mat */
662531e53bdSHong Zhang     Mat_MPIAIJ *aij=(Mat_MPIAIJ*)mat->data;
663531e53bdSHong Zhang     Mat_SeqAIJ *spA,*spB;
664531e53bdSHong Zhang     Mat        A,B;
665a8971b87SHong Zhang     PetscInt   nz,brows,bcols;
666531e53bdSHong Zhang     PetscReal  mem;
667531e53bdSHong Zhang 
668531e53bdSHong Zhang     bs    = 1; /* only bs=1 is supported for MPIAIJ matrix */
669531e53bdSHong Zhang 
670531e53bdSHong Zhang     A = aij->A;  spA = (Mat_SeqAIJ*)A->data;
671531e53bdSHong Zhang     B = aij->B;  spB = (Mat_SeqAIJ*)B->data;
672531e53bdSHong Zhang     nz = spA->nz + spB->nz; /* total local nonzero entries of mat */
673531e53bdSHong Zhang     mem = nz*(sizeof(PetscScalar) + sizeof(PetscInt)) + 3*m*sizeof(PetscInt);
674531e53bdSHong Zhang     bcols = (PetscInt)(0.5*mem /(m*sizeof(PetscScalar)));
675531e53bdSHong Zhang     brows = 1000/bcols;
676531e53bdSHong Zhang     if (bcols > nis) bcols = nis;
677531e53bdSHong Zhang     if (brows == 0 || brows > m) brows = m;
678531e53bdSHong Zhang     c->brows = brows;
679531e53bdSHong Zhang     c->bcols = bcols;
680531e53bdSHong Zhang   }
681a8971b87SHong Zhang 
682a8971b87SHong Zhang   c->M       = mat->rmap->N/bs;         /* set the global rows and columns and local rows */
683a8971b87SHong Zhang   c->N       = mat->cmap->N/bs;
684a8971b87SHong Zhang   c->m       = mat->rmap->n/bs;
685a8971b87SHong Zhang   c->rstart  = mat->rmap->rstart/bs;
686a8971b87SHong Zhang   c->ncolors = nis;
68793dfae19SHong Zhang   PetscFunctionReturn(0);
68893dfae19SHong Zhang }
689