xref: /petsc/src/mat/impls/aij/mpi/fdmpiaij.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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;
25b470e4b4SRichard Tran Mills   ierr = VecBindToCPU(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);
64ce911d59SRichard Tran Mills     /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */
65b470e4b4SRichard Tran Mills     ierr = VecBindToCPU(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;
158b470e4b4SRichard Tran Mills   ierr = VecBindToCPU(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;
179b294de21SRichard Tran Mills   PetscBool         alreadyboundtocpu;
180fcd7ac73SHong Zhang 
181fcd7ac73SHong Zhang   PetscFunctionBegin;
182b294de21SRichard Tran Mills   ierr = VecBoundToCPU(x1,&alreadyboundtocpu);CHKERRQ(ierr);
183b470e4b4SRichard Tran Mills   ierr = VecBindToCPU(x1,PETSC_TRUE);CHKERRQ(ierr);
184*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse((ctype == IS_COLORING_LOCAL) && (J->ops->fdcoloringapply == MatFDColoringApply_AIJ),PetscObjectComm((PetscObject)J),PETSC_ERR_SUP,"Must call MatColoringUseDM() with IS_COLORING_LOCAL");
1858bc97078SHong Zhang   /* (1) Set w1 = F(x1) */
186fcd7ac73SHong Zhang   if (!coloring->fset) {
187fcd7ac73SHong Zhang     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
188f6af9589SHong Zhang     ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
189fcd7ac73SHong Zhang     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
190fcd7ac73SHong Zhang   } else {
191fcd7ac73SHong Zhang     coloring->fset = PETSC_FALSE;
192fcd7ac73SHong Zhang   }
193fcd7ac73SHong Zhang 
1948bc97078SHong Zhang   /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */
195f6af9589SHong Zhang   if (coloring->htype[0] == 'w') {
196684f2004SHong Zhang     /* vscale = 1./dx is a constant scalar */
197f6af9589SHong Zhang     ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr);
198c53567a0SHong Zhang     dx = 1.0/(PetscSqrtReal(1.0 + unorm)*epsilon);
19970e7395fSHong Zhang   } else {
20074d3cef9SHong Zhang     ierr = VecGetLocalSize(x1,&nxloc);CHKERRQ(ierr);
2015edff71fSBarry Smith     ierr = VecGetArrayRead(x1,&xx);CHKERRQ(ierr);
2028bc97078SHong Zhang     ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr);
20374d3cef9SHong Zhang     for (col=0; col<nxloc; col++) {
204fcd7ac73SHong Zhang       dx = xx[col];
20574d3cef9SHong Zhang       if (PetscAbsScalar(dx) < umin) {
20674d3cef9SHong Zhang         if (PetscRealPart(dx) >= 0.0)      dx = umin;
20774d3cef9SHong Zhang         else if (PetscRealPart(dx) < 0.0) dx = -umin;
20874d3cef9SHong Zhang       }
209fcd7ac73SHong Zhang       dx               *= epsilon;
21074d3cef9SHong Zhang       vscale_array[col] = 1.0/dx;
211f6af9589SHong Zhang     }
2125edff71fSBarry Smith     ierr = VecRestoreArrayRead(x1,&xx);CHKERRQ(ierr);
2138bc97078SHong Zhang     ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr);
21470e7395fSHong Zhang   }
215684f2004SHong Zhang   if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') {
2168bc97078SHong Zhang     ierr = VecGhostUpdateBegin(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2178bc97078SHong Zhang     ierr = VecGhostUpdateEnd(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
218fcd7ac73SHong Zhang   }
219fcd7ac73SHong Zhang 
2208bc97078SHong Zhang   /* (3) Loop over each color */
2218bc97078SHong Zhang   if (!coloring->w3) {
2228bc97078SHong Zhang     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
2238bc97078SHong Zhang     ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr);
2248bc97078SHong Zhang   }
2258bc97078SHong Zhang   w3 = coloring->w3;
226fcd7ac73SHong Zhang 
2278bc97078SHong Zhang   ierr = VecGetOwnershipRange(x1,&cstart,&cend);CHKERRQ(ierr); /* used by ghosted vscale */
2289e917edbSHong Zhang   if (vscale) {
2298bc97078SHong Zhang     ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr);
2309e917edbSHong Zhang   }
2318bc97078SHong Zhang   nz = 0;
232e2c857f8SHong Zhang 
233b3e1f37bSHong Zhang   if (coloring->bcols > 1) { /* use blocked insertion of Jentry */
234b3e1f37bSHong Zhang     PetscInt    i,m=J->rmap->n,nbcols,bcols=coloring->bcols;
235e2c857f8SHong Zhang     PetscScalar *dy=coloring->dy,*dy_k;
236e2c857f8SHong Zhang 
237e2c857f8SHong Zhang     nbcols = 0;
238e2c857f8SHong Zhang     for (k=0; k<ncolors; k+=bcols) {
239e2c857f8SHong Zhang 
240e2c857f8SHong Zhang       /*
241e2c857f8SHong Zhang        (3-1) Loop over each column associated with color
242e2c857f8SHong Zhang        adding the perturbation to the vector w3 = x1 + dx.
243e2c857f8SHong Zhang        */
244e2c857f8SHong Zhang 
245e2c857f8SHong Zhang       dy_k = dy;
246e2c857f8SHong Zhang       if (k + bcols > ncolors) bcols = ncolors - k;
247e2c857f8SHong Zhang       for (i=0; i<bcols; i++) {
248ceef8a49SBarry Smith         coloring->currentcolor = k+i;
249e2c857f8SHong Zhang 
250e2c857f8SHong Zhang         ierr = VecCopy(x1,w3);CHKERRQ(ierr);
251e2c857f8SHong Zhang         ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
252e2c857f8SHong Zhang         if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */
253e2c857f8SHong Zhang         if (coloring->htype[0] == 'w') {
254e2c857f8SHong Zhang           for (l=0; l<ncolumns[k+i]; l++) {
255e2c857f8SHong Zhang             col = coloring->columns[k+i][l]; /* local column (in global index!) of the matrix we are probing for */
256e2c857f8SHong Zhang             w3_array[col] += 1.0/dx;
257e2c857f8SHong Zhang           }
258e2c857f8SHong Zhang         } else { /* htype == 'ds' */
259e2c857f8SHong Zhang           vscale_array -= cstart; /* shift pointer so global index can be used */
260e2c857f8SHong Zhang           for (l=0; l<ncolumns[k+i]; l++) {
261e2c857f8SHong Zhang             col = coloring->columns[k+i][l]; /* local column (in global index!) of the matrix we are probing for */
262e2c857f8SHong Zhang             w3_array[col] += 1.0/vscale_array[col];
263e2c857f8SHong Zhang           }
264e2c857f8SHong Zhang           vscale_array += cstart;
265e2c857f8SHong Zhang         }
266e2c857f8SHong Zhang         if (ctype == IS_COLORING_GLOBAL) w3_array += cstart;
267e2c857f8SHong Zhang         ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
268e2c857f8SHong Zhang 
269e2c857f8SHong Zhang         /*
270e2c857f8SHong Zhang          (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
271e2c857f8SHong Zhang                            w2 = F(x1 + dx) - F(x1)
272e2c857f8SHong Zhang          */
273e2c857f8SHong Zhang         ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
274e2c857f8SHong Zhang         ierr = VecPlaceArray(w2,dy_k);CHKERRQ(ierr); /* place w2 to the array dy_i */
275e2c857f8SHong Zhang         ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
276e2c857f8SHong Zhang         ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
277e2c857f8SHong Zhang         ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
278e2c857f8SHong Zhang         ierr = VecResetArray(w2);CHKERRQ(ierr);
279e2c857f8SHong Zhang         dy_k += m; /* points to dy+i*nxloc */
280e2c857f8SHong Zhang       }
281e2c857f8SHong Zhang 
282e2c857f8SHong Zhang       /*
283e2c857f8SHong Zhang        (3-3) Loop over block rows of vector, putting results into Jacobian matrix
284e2c857f8SHong Zhang        */
285d880da65SHong Zhang       nrows_k = nrows[nbcols++];
286d880da65SHong Zhang 
287e2c857f8SHong Zhang       if (coloring->htype[0] == 'w') {
288e2c857f8SHong Zhang         for (l=0; l<nrows_k; l++) {
2890df34763SHong Zhang           row  = Jentry2[nz].row;   /* local row index */
2905059b303SJunchao Zhang           /* The 'useless' ifdef is due to a bug in NVIDIA nvc 21.11, which triggers a segfault on this line. We write it in
2915059b303SJunchao Zhang              another way, and it seems work. See https://lists.mcs.anl.gov/pipermail/petsc-users/2021-December/045158.html
2925059b303SJunchao Zhang            */
2935059b303SJunchao Zhang          #if defined(PETSC_USE_COMPLEX)
2945059b303SJunchao Zhang           PetscScalar *tmp = Jentry2[nz].valaddr;
2955059b303SJunchao Zhang           *tmp = dy[row]*dx;
2965059b303SJunchao Zhang          #else
2975059b303SJunchao Zhang           *(Jentry2[nz].valaddr) = dy[row]*dx;
2985059b303SJunchao Zhang          #endif
2995059b303SJunchao Zhang           nz++;
300e2c857f8SHong Zhang         }
301e2c857f8SHong Zhang       } else { /* htype == 'ds' */
302e2c857f8SHong Zhang         for (l=0; l<nrows_k; l++) {
303e2c857f8SHong Zhang           row = Jentry[nz].row;   /* local row index */
3045059b303SJunchao Zhang          #if defined(PETSC_USE_COMPLEX) /* See https://lists.mcs.anl.gov/pipermail/petsc-users/2021-December/045158.html */
3055059b303SJunchao Zhang           PetscScalar *tmp = Jentry[nz].valaddr;
3065059b303SJunchao Zhang           *tmp = dy[row]*vscale_array[Jentry[nz].col];
3075059b303SJunchao Zhang          #else
308d880da65SHong Zhang           *(Jentry[nz].valaddr) = dy[row]*vscale_array[Jentry[nz].col];
3095059b303SJunchao Zhang          #endif
310e2c857f8SHong Zhang           nz++;
311e2c857f8SHong Zhang         }
312e2c857f8SHong Zhang       }
313e2c857f8SHong Zhang     }
314b3e1f37bSHong Zhang   } else { /* bcols == 1 */
3158bc97078SHong Zhang     for (k=0; k<ncolors; k++) {
3168bc97078SHong Zhang       coloring->currentcolor = k;
3179e917edbSHong Zhang 
318fcd7ac73SHong Zhang       /*
3198bc97078SHong Zhang        (3-1) Loop over each column associated with color
320f6af9589SHong Zhang        adding the perturbation to the vector w3 = x1 + dx.
321fcd7ac73SHong Zhang        */
322f6af9589SHong Zhang       ierr = VecCopy(x1,w3);CHKERRQ(ierr);
323f6af9589SHong Zhang       ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
324a2f2d239SHong Zhang       if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */
325b039d6c4SHong Zhang       if (coloring->htype[0] == 'w') {
3268bc97078SHong Zhang         for (l=0; l<ncolumns[k]; l++) {
327f6af9589SHong Zhang           col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
3289e917edbSHong Zhang           w3_array[col] += 1.0/dx;
3299e917edbSHong Zhang         }
330b039d6c4SHong Zhang       } else { /* htype == 'ds' */
331a2f2d239SHong Zhang         vscale_array -= cstart; /* shift pointer so global index can be used */
332b039d6c4SHong Zhang         for (l=0; l<ncolumns[k]; l++) {
333b039d6c4SHong Zhang           col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
334a2f2d239SHong Zhang           w3_array[col] += 1.0/vscale_array[col];
33570e7395fSHong Zhang         }
336a2f2d239SHong Zhang         vscale_array += cstart;
337fcd7ac73SHong Zhang       }
338a2f2d239SHong Zhang       if (ctype == IS_COLORING_GLOBAL) w3_array += cstart;
339fcd7ac73SHong Zhang       ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
3409e917edbSHong Zhang 
341fcd7ac73SHong Zhang       /*
3428bc97078SHong Zhang        (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
343fcd7ac73SHong Zhang                            w2 = F(x1 + dx) - F(x1)
344fcd7ac73SHong Zhang        */
345fcd7ac73SHong Zhang       ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
346fcd7ac73SHong Zhang       ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
347fcd7ac73SHong Zhang       ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
348fcd7ac73SHong Zhang       ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
3499e917edbSHong Zhang 
3508bc97078SHong Zhang       /*
3518bc97078SHong Zhang        (3-3) Loop over rows of vector, putting results into Jacobian matrix
3528bc97078SHong Zhang        */
353c53567a0SHong Zhang       nrows_k = nrows[k];
354fcd7ac73SHong Zhang       ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
355b039d6c4SHong Zhang       if (coloring->htype[0] == 'w') {
356c53567a0SHong Zhang         for (l=0; l<nrows_k; l++) {
3570df34763SHong Zhang           row  = Jentry2[nz].row;   /* local row index */
3585059b303SJunchao Zhang          #if defined(PETSC_USE_COMPLEX) /* See https://lists.mcs.anl.gov/pipermail/petsc-users/2021-December/045158.html */
3595059b303SJunchao Zhang           PetscScalar *tmp = Jentry2[nz].valaddr;
3605059b303SJunchao Zhang           *tmp = y[row]*dx;
3615059b303SJunchao Zhang          #else
3625059b303SJunchao Zhang           *(Jentry2[nz].valaddr) = y[row]*dx;
3635059b303SJunchao Zhang          #endif
3645059b303SJunchao Zhang           nz++;
3659e917edbSHong Zhang         }
366b039d6c4SHong Zhang       } else { /* htype == 'ds' */
367c53567a0SHong Zhang         for (l=0; l<nrows_k; l++) {
368b039d6c4SHong Zhang           row  = Jentry[nz].row;   /* local row index */
3695059b303SJunchao Zhang          #if defined(PETSC_USE_COMPLEX) /* See https://lists.mcs.anl.gov/pipermail/petsc-users/2021-December/045158.html */
3705059b303SJunchao Zhang           PetscScalar *tmp = Jentry[nz].valaddr;
3715059b303SJunchao Zhang           *tmp = y[row]*vscale_array[Jentry[nz].col];
3725059b303SJunchao Zhang          #else
373b039d6c4SHong Zhang           *(Jentry[nz].valaddr) = y[row]*vscale_array[Jentry[nz].col];
3745059b303SJunchao Zhang          #endif
375573f477fSHong Zhang           nz++;
376fcd7ac73SHong Zhang         }
377b039d6c4SHong Zhang       }
378fcd7ac73SHong Zhang       ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
3798bc97078SHong Zhang     }
380b3e1f37bSHong Zhang   }
381b3e1f37bSHong Zhang 
382e2cf4d64SStefano Zampini #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
383c70f7ee4SJunchao Zhang   if (J->offloadmask != PETSC_OFFLOAD_UNALLOCATED) J->offloadmask = PETSC_OFFLOAD_CPU;
384e2cf4d64SStefano Zampini #endif
385f5aae955SHong Zhang   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
386f5aae955SHong Zhang   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3879e917edbSHong Zhang   if (vscale) {
3888bc97078SHong Zhang     ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr);
3899e917edbSHong Zhang   }
390fcd7ac73SHong Zhang   coloring->currentcolor = -1;
391b294de21SRichard Tran Mills   if (!alreadyboundtocpu) {ierr = VecBindToCPU(x1,PETSC_FALSE);CHKERRQ(ierr);}
392fcd7ac73SHong Zhang   PetscFunctionReturn(0);
393fcd7ac73SHong Zhang }
394fcd7ac73SHong Zhang 
395f86b9fbaSHong Zhang PetscErrorCode MatFDColoringSetUp_MPIXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
396fcd7ac73SHong Zhang {
397fcd7ac73SHong Zhang   PetscErrorCode         ierr;
398fcd7ac73SHong Zhang   PetscMPIInt            size,*ncolsonproc,*disp,nn;
399e3225e9dSHong Zhang   PetscInt               i,n,nrows,nrows_i,j,k,m,ncols,col,*rowhit,cstart,cend,colb;
40072c15787SHong Zhang   const PetscInt         *is,*A_ci,*A_cj,*B_ci,*B_cj,*row=NULL,*ltog=NULL;
401071fcb05SBarry Smith   PetscInt               nis=iscoloring->n,nctot,*cols,tmp = 0;
402fcd7ac73SHong Zhang   ISLocalToGlobalMapping map=mat->cmap->mapping;
403a8971b87SHong Zhang   PetscInt               ctype=c->ctype,*spidxA,*spidxB,nz,bs,bs2,spidx;
4040d1c53f1SHong Zhang   Mat                    A,B;
405e3225e9dSHong Zhang   PetscScalar            *A_val,*B_val,**valaddrhit;
406573f477fSHong Zhang   MatEntry               *Jentry;
4070df34763SHong Zhang   MatEntry2              *Jentry2;
408d4002b98SHong Zhang   PetscBool              isBAIJ,isSELL;
409531e53bdSHong Zhang   PetscInt               bcols=c->bcols;
4100d1c53f1SHong Zhang #if defined(PETSC_USE_CTABLE)
4110d1c53f1SHong Zhang   PetscTable             colmap=NULL;
4120d1c53f1SHong Zhang #else
4130d1c53f1SHong Zhang   PetscInt               *colmap=NULL;     /* local col number of off-diag col */
4140d1c53f1SHong Zhang #endif
415fcd7ac73SHong Zhang 
416fcd7ac73SHong Zhang   PetscFunctionBegin;
4175bdb020cSBarry Smith   if (ctype == IS_COLORING_LOCAL) {
418*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!map,PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_INCOMP,"When using ghosted differencing matrix must have local to global mapping provided with MatSetLocalToGlobalMapping");
41972c15787SHong Zhang     ierr = ISLocalToGlobalMappingGetIndices(map,&ltog);CHKERRQ(ierr);
42072c15787SHong Zhang   }
421fcd7ac73SHong Zhang 
4220d1c53f1SHong Zhang   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
423b9e7e5c1SBarry Smith   ierr = PetscObjectBaseTypeCompare((PetscObject)mat,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr);
424d4002b98SHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPISELL,&isSELL);CHKERRQ(ierr);
425e3225e9dSHong Zhang   if (isBAIJ) {
426195687c5SHong Zhang     Mat_MPIBAIJ *baij=(Mat_MPIBAIJ*)mat->data;
4276c145f34SHong Zhang     Mat_SeqBAIJ *spA,*spB;
428195687c5SHong Zhang     A = baij->A;  spA = (Mat_SeqBAIJ*)A->data; A_val = spA->a;
429195687c5SHong Zhang     B = baij->B;  spB = (Mat_SeqBAIJ*)B->data; B_val = spB->a;
4300d1c53f1SHong Zhang     nz = spA->nz + spB->nz; /* total nonzero entries of mat */
431195687c5SHong Zhang     if (!baij->colmap) {
4320d1c53f1SHong Zhang       ierr = MatCreateColmap_MPIBAIJ_Private(mat);CHKERRQ(ierr);
4330d1c53f1SHong Zhang     }
434097e26fcSJed Brown     colmap = baij->colmap;
4350d1c53f1SHong Zhang     ierr = MatGetColumnIJ_SeqBAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
4360d1c53f1SHong Zhang     ierr = MatGetColumnIJ_SeqBAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
437195687c5SHong Zhang 
438195687c5SHong Zhang     if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') {  /* create vscale for storing dx */
439195687c5SHong Zhang       PetscInt    *garray;
440785e854fSJed Brown       ierr = PetscMalloc1(B->cmap->n,&garray);CHKERRQ(ierr);
441195687c5SHong Zhang       for (i=0; i<baij->B->cmap->n/bs; i++) {
442195687c5SHong Zhang         for (j=0; j<bs; j++) {
443195687c5SHong Zhang           garray[i*bs+j] = bs*baij->garray[i]+j;
444195687c5SHong Zhang         }
445195687c5SHong Zhang       }
446195687c5SHong Zhang       ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->cmap->n,PETSC_DETERMINE,B->cmap->n,garray,&c->vscale);CHKERRQ(ierr);
447b470e4b4SRichard Tran Mills       ierr = VecBindToCPU(c->vscale,PETSC_TRUE);CHKERRQ(ierr);
448195687c5SHong Zhang       ierr = PetscFree(garray);CHKERRQ(ierr);
449195687c5SHong Zhang     }
450d4002b98SHong Zhang   } else if (isSELL) {
451d4002b98SHong Zhang     Mat_MPISELL *sell=(Mat_MPISELL*)mat->data;
452d4002b98SHong Zhang     Mat_SeqSELL *spA,*spB;
453d4002b98SHong Zhang     A = sell->A;  spA = (Mat_SeqSELL*)A->data; A_val = spA->val;
454d4002b98SHong Zhang     B = sell->B;  spB = (Mat_SeqSELL*)B->data; B_val = spB->val;
455f4b713efSHong Zhang     nz = spA->nz + spB->nz; /* total nonzero entries of mat */
456d4002b98SHong Zhang     if (!sell->colmap) {
457f4b713efSHong Zhang       /* Allow access to data structures of local part of matrix
458f4b713efSHong Zhang        - creates aij->colmap which maps global column number to local number in part B */
459d4002b98SHong Zhang       ierr = MatCreateColmap_MPISELL_Private(mat);CHKERRQ(ierr);
460f4b713efSHong Zhang     }
461d4002b98SHong Zhang     colmap = sell->colmap;
462d4002b98SHong Zhang     ierr = MatGetColumnIJ_SeqSELL_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
463d4002b98SHong Zhang     ierr = MatGetColumnIJ_SeqSELL_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
464f4b713efSHong Zhang 
465f4b713efSHong Zhang     bs = 1; /* only bs=1 is supported for non MPIBAIJ matrix */
466f4b713efSHong Zhang 
467f4b713efSHong Zhang     if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') { /* create vscale for storing dx */
468d4002b98SHong Zhang       ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->cmap->n,PETSC_DETERMINE,B->cmap->n,sell->garray,&c->vscale);CHKERRQ(ierr);
469b470e4b4SRichard Tran Mills       ierr = VecBindToCPU(c->vscale,PETSC_TRUE);CHKERRQ(ierr);
470f4b713efSHong Zhang     }
471e3225e9dSHong Zhang   } else {
472e3225e9dSHong Zhang     Mat_MPIAIJ *aij=(Mat_MPIAIJ*)mat->data;
473e3225e9dSHong Zhang     Mat_SeqAIJ *spA,*spB;
474e3225e9dSHong Zhang     A = aij->A;  spA = (Mat_SeqAIJ*)A->data; A_val = spA->a;
475e3225e9dSHong Zhang     B = aij->B;  spB = (Mat_SeqAIJ*)B->data; B_val = spB->a;
476e3225e9dSHong Zhang     nz = spA->nz + spB->nz; /* total nonzero entries of mat */
477e3225e9dSHong Zhang     if (!aij->colmap) {
478e3225e9dSHong Zhang       /* Allow access to data structures of local part of matrix
479e3225e9dSHong Zhang        - creates aij->colmap which maps global column number to local number in part B */
480e3225e9dSHong Zhang       ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
481e3225e9dSHong Zhang     }
482097e26fcSJed Brown     colmap = aij->colmap;
483e3225e9dSHong Zhang     ierr = MatGetColumnIJ_SeqAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
484e3225e9dSHong Zhang     ierr = MatGetColumnIJ_SeqAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
485e3225e9dSHong Zhang 
486e3225e9dSHong Zhang     bs = 1; /* only bs=1 is supported for non MPIBAIJ matrix */
487195687c5SHong Zhang 
488195687c5SHong Zhang     if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') { /* create vscale for storing dx */
489195687c5SHong Zhang       ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->cmap->n,PETSC_DETERMINE,B->cmap->n,aij->garray,&c->vscale);CHKERRQ(ierr);
490b470e4b4SRichard Tran Mills       ierr = VecBindToCPU(c->vscale,PETSC_TRUE);CHKERRQ(ierr);
491195687c5SHong Zhang     }
4920d1c53f1SHong Zhang   }
4930d1c53f1SHong Zhang 
4940d1c53f1SHong Zhang   m      = mat->rmap->n/bs;
4950d1c53f1SHong Zhang   cstart = mat->cmap->rstart/bs;
4960d1c53f1SHong Zhang   cend   = mat->cmap->rend/bs;
497fcd7ac73SHong Zhang 
498071fcb05SBarry Smith   ierr = PetscMalloc2(nis,&c->ncolumns,nis,&c->columns);CHKERRQ(ierr);
499071fcb05SBarry Smith   ierr = PetscMalloc1(nis,&c->nrows);CHKERRQ(ierr);
5008bc97078SHong Zhang   ierr = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr);
501fcd7ac73SHong Zhang 
5020df34763SHong Zhang   if (c->htype[0] == 'd') {
503785e854fSJed Brown     ierr        = PetscMalloc1(nz,&Jentry);CHKERRQ(ierr);
5048bc97078SHong Zhang     ierr        = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr);
5058bc97078SHong Zhang     c->matentry = Jentry;
5060df34763SHong Zhang   } else if (c->htype[0] == 'w') {
507785e854fSJed Brown     ierr         = PetscMalloc1(nz,&Jentry2);CHKERRQ(ierr);
5080df34763SHong Zhang     ierr         = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry2));CHKERRQ(ierr);
5090df34763SHong Zhang     c->matentry2 = Jentry2;
5100df34763SHong Zhang   } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"htype is not supported");
511a774a6f1SHong Zhang 
512dcca6d9dSJed Brown   ierr = PetscMalloc2(m+1,&rowhit,m+1,&valaddrhit);CHKERRQ(ierr);
513fcd7ac73SHong Zhang   nz   = 0;
514071fcb05SBarry Smith   ierr = ISColoringGetIS(iscoloring,PETSC_OWN_POINTER, PETSC_IGNORE,&c->isa);CHKERRQ(ierr);
515071fcb05SBarry Smith 
516071fcb05SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
517ffc4695bSBarry Smith     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRMPI(ierr);
518071fcb05SBarry Smith     ierr = PetscMalloc2(size,&ncolsonproc,size,&disp);CHKERRQ(ierr);
519071fcb05SBarry Smith   }
520071fcb05SBarry Smith 
521d3825b63SHong Zhang   for (i=0; i<nis; i++) { /* for each local color */
522071fcb05SBarry Smith     ierr = ISGetLocalSize(c->isa[i],&n);CHKERRQ(ierr);
523071fcb05SBarry Smith     ierr = ISGetIndices(c->isa[i],&is);CHKERRQ(ierr);
524fcd7ac73SHong Zhang 
525fcd7ac73SHong Zhang     c->ncolumns[i] = n; /* local number of columns of this color on this process */
526071fcb05SBarry Smith     c->columns[i]  = (PetscInt*)is;
527fcd7ac73SHong Zhang 
528fcd7ac73SHong Zhang     if (ctype == IS_COLORING_GLOBAL) {
529d3825b63SHong Zhang       /* Determine nctot, the total (parallel) number of columns of this color */
530d3825b63SHong Zhang       /* ncolsonproc[j]: local ncolumns on proc[j] of this color */
531fcd7ac73SHong Zhang       ierr  = PetscMPIIntCast(n,&nn);CHKERRQ(ierr);
53255b25c41SPierre Jolivet       ierr  = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRMPI(ierr);
533fcd7ac73SHong Zhang       nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j];
534fcd7ac73SHong Zhang       if (!nctot) {
535fcd7ac73SHong Zhang         ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr);
536fcd7ac73SHong Zhang       }
537fcd7ac73SHong Zhang 
538fcd7ac73SHong Zhang       disp[0] = 0;
539fcd7ac73SHong Zhang       for (j=1; j<size; j++) {
540fcd7ac73SHong Zhang         disp[j] = disp[j-1] + ncolsonproc[j-1];
541fcd7ac73SHong Zhang       }
542fcd7ac73SHong Zhang 
543d3825b63SHong Zhang       /* Get cols, the complete list of columns for this color on each process */
544854ce69bSBarry Smith       ierr = PetscMalloc1(nctot+1,&cols);CHKERRQ(ierr);
545ffc4695bSBarry Smith       ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRMPI(ierr);
5465bdb020cSBarry Smith     } else if (ctype == IS_COLORING_LOCAL) {
547fcd7ac73SHong Zhang       /* Determine local number of columns of this color on this process, including ghost points */
548fcd7ac73SHong Zhang       nctot = n;
549071fcb05SBarry Smith       cols  = (PetscInt*)is;
550fcd7ac73SHong Zhang     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type");
551fcd7ac73SHong Zhang 
5521b97d346SHong Zhang     /* Mark all rows affect by these columns */
553580bdb30SBarry Smith     ierr    = PetscArrayzero(rowhit,m);CHKERRQ(ierr);
5540d1c53f1SHong Zhang     bs2     = bs*bs;
5554b2e90caSHong Zhang     nrows_i = 0;
5561b97d346SHong Zhang     for (j=0; j<nctot; j++) { /* loop over columns*/
5575bdb020cSBarry Smith       if (ctype == IS_COLORING_LOCAL) {
558fcd7ac73SHong Zhang         col = ltog[cols[j]];
559fcd7ac73SHong Zhang       } else {
560fcd7ac73SHong Zhang         col = cols[j];
561fcd7ac73SHong Zhang       }
562e3225e9dSHong Zhang       if (col >= cstart && col < cend) { /* column is in A, diagonal block of mat */
563071fcb05SBarry Smith         tmp      = A_ci[col-cstart];
564071fcb05SBarry Smith         row      = A_cj + tmp;
565071fcb05SBarry Smith         nrows    = A_ci[col-cstart+1] - tmp;
5664b2e90caSHong Zhang         nrows_i += nrows;
567071fcb05SBarry Smith 
568fcd7ac73SHong Zhang         /* loop over columns of A marking them in rowhit */
569d3825b63SHong Zhang         for (k=0; k<nrows; k++) {
570d3825b63SHong Zhang           /* set valaddrhit for part A */
571071fcb05SBarry Smith           spidx            = bs2*spidxA[tmp + k];
572e3225e9dSHong Zhang           valaddrhit[*row] = &A_val[spidx];
573a774a6f1SHong Zhang           rowhit[*row++]   = col - cstart + 1; /* local column index */
574fcd7ac73SHong Zhang         }
575e3225e9dSHong Zhang       } else { /* column is in B, off-diagonal block of mat */
576fcd7ac73SHong Zhang #if defined(PETSC_USE_CTABLE)
5770d1c53f1SHong Zhang         ierr = PetscTableFind(colmap,col+1,&colb);CHKERRQ(ierr);
578fcd7ac73SHong Zhang         colb--;
579fcd7ac73SHong Zhang #else
5800d1c53f1SHong Zhang         colb = colmap[col] - 1; /* local column index */
581fcd7ac73SHong Zhang #endif
582fcd7ac73SHong Zhang         if (colb == -1) {
583d3825b63SHong Zhang           nrows = 0;
584fcd7ac73SHong Zhang         } else {
5850d1c53f1SHong Zhang           colb  = colb/bs;
586071fcb05SBarry Smith           tmp   = B_ci[colb];
587071fcb05SBarry Smith           row   = B_cj + tmp;
588071fcb05SBarry Smith           nrows = B_ci[colb+1] - tmp;
589fcd7ac73SHong Zhang         }
5904b2e90caSHong Zhang         nrows_i += nrows;
591fcd7ac73SHong Zhang         /* loop over columns of B marking them in rowhit */
592d3825b63SHong Zhang         for (k=0; k<nrows; k++) {
593d3825b63SHong Zhang           /* set valaddrhit for part B */
594071fcb05SBarry Smith           spidx            = bs2*spidxB[tmp + k];
595e3225e9dSHong Zhang           valaddrhit[*row] = &B_val[spidx];
59670e7395fSHong Zhang           rowhit[*row++]   = colb + 1 + cend - cstart; /* local column index */
597fcd7ac73SHong Zhang         }
598fcd7ac73SHong Zhang       }
599fcd7ac73SHong Zhang     }
6004b2e90caSHong Zhang     c->nrows[i] = nrows_i;
6018bc97078SHong Zhang 
6020df34763SHong Zhang     if (c->htype[0] == 'd') {
603d3825b63SHong Zhang       for (j=0; j<m; j++) {
604fcd7ac73SHong Zhang         if (rowhit[j]) {
605573f477fSHong Zhang           Jentry[nz].row     = j;              /* local row index */
606573f477fSHong Zhang           Jentry[nz].col     = rowhit[j] - 1;  /* local column index */
607573f477fSHong Zhang           Jentry[nz].valaddr = valaddrhit[j];  /* address of mat value for this entry */
608573f477fSHong Zhang           nz++;
609fcd7ac73SHong Zhang         }
610fcd7ac73SHong Zhang       }
6110df34763SHong Zhang     } else { /* c->htype == 'wp' */
6120df34763SHong Zhang       for (j=0; j<m; j++) {
6130df34763SHong Zhang         if (rowhit[j]) {
6140df34763SHong Zhang           Jentry2[nz].row     = j;              /* local row index */
6150df34763SHong Zhang           Jentry2[nz].valaddr = valaddrhit[j];  /* address of mat value for this entry */
6160df34763SHong Zhang           nz++;
6170df34763SHong Zhang         }
6180df34763SHong Zhang       }
6190df34763SHong Zhang     }
620071fcb05SBarry Smith     if (ctype == IS_COLORING_GLOBAL) {
621fcd7ac73SHong Zhang       ierr = PetscFree(cols);CHKERRQ(ierr);
622fcd7ac73SHong Zhang     }
623071fcb05SBarry Smith   }
624071fcb05SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
625071fcb05SBarry Smith     ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr);
626071fcb05SBarry Smith   }
627fcd7ac73SHong Zhang 
628a8971b87SHong Zhang   if (bcols > 1) { /* reorder Jentry for faster MatFDColoringApply() */
629a8971b87SHong Zhang     ierr = MatFDColoringSetUpBlocked_AIJ_Private(mat,c,nz);CHKERRQ(ierr);
630531e53bdSHong Zhang   }
631531e53bdSHong Zhang 
6320d1c53f1SHong Zhang   if (isBAIJ) {
6330d1c53f1SHong Zhang     ierr = MatRestoreColumnIJ_SeqBAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
6340d1c53f1SHong Zhang     ierr = MatRestoreColumnIJ_SeqBAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
635785e854fSJed Brown     ierr = PetscMalloc1(bs*mat->rmap->n,&c->dy);CHKERRQ(ierr);
636d4002b98SHong Zhang   } else if (isSELL) {
637d4002b98SHong Zhang     ierr = MatRestoreColumnIJ_SeqSELL_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
638d4002b98SHong Zhang     ierr = MatRestoreColumnIJ_SeqSELL_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
6390d1c53f1SHong Zhang   } else {
6400d1c53f1SHong Zhang     ierr = MatRestoreColumnIJ_SeqAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
6410d1c53f1SHong Zhang     ierr = MatRestoreColumnIJ_SeqAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
6420d1c53f1SHong Zhang   }
6430d1c53f1SHong Zhang 
644071fcb05SBarry Smith   ierr = ISColoringRestoreIS(iscoloring,PETSC_OWN_POINTER,&c->isa);CHKERRQ(ierr);
645a8971b87SHong Zhang   ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr);
646a8971b87SHong Zhang 
6475bdb020cSBarry Smith   if (ctype == IS_COLORING_LOCAL) {
64872c15787SHong Zhang     ierr = ISLocalToGlobalMappingRestoreIndices(map,&ltog);CHKERRQ(ierr);
64972c15787SHong Zhang   }
6507d3de750SJacob Faibussowitsch   ierr = PetscInfo(c,"ncolors %" PetscInt_FMT ", brows %" PetscInt_FMT " and bcols %" PetscInt_FMT " are used.\n",c->ncolors,c->brows,c->bcols);CHKERRQ(ierr);
651fcd7ac73SHong Zhang   PetscFunctionReturn(0);
652fcd7ac73SHong Zhang }
65393dfae19SHong Zhang 
65493dfae19SHong Zhang PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
65593dfae19SHong Zhang {
656531e53bdSHong Zhang   PetscErrorCode ierr;
657a8971b87SHong Zhang   PetscInt       bs,nis=iscoloring->n,m=mat->rmap->n;
658d4002b98SHong Zhang   PetscBool      isBAIJ,isSELL;
659531e53bdSHong Zhang 
66093dfae19SHong Zhang   PetscFunctionBegin;
661531e53bdSHong Zhang   /* set default brows and bcols for speedup inserting the dense matrix into sparse Jacobian;
662531e53bdSHong Zhang    bcols is chosen s.t. dy-array takes 50% of memory space as mat */
663531e53bdSHong Zhang   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
664b9e7e5c1SBarry Smith   ierr = PetscObjectBaseTypeCompare((PetscObject)mat,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr);
665d4002b98SHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPISELL,&isSELL);CHKERRQ(ierr);
666b8dfa574SBarry Smith   if (isBAIJ || m == 0) {
667a8971b87SHong Zhang     c->brows = m;
668a8971b87SHong Zhang     c->bcols = 1;
669d4002b98SHong Zhang   } else if (isSELL) {
670f4b713efSHong Zhang     /* bcols is chosen s.t. dy-array takes 50% of local memory space as mat */
671d4002b98SHong Zhang     Mat_MPISELL *sell=(Mat_MPISELL*)mat->data;
672d4002b98SHong Zhang     Mat_SeqSELL *spA,*spB;
673f4b713efSHong Zhang     Mat        A,B;
674f4b713efSHong Zhang     PetscInt   nz,brows,bcols;
675f4b713efSHong Zhang     PetscReal  mem;
676f4b713efSHong Zhang 
677d4002b98SHong Zhang     bs    = 1; /* only bs=1 is supported for MPISELL matrix */
678f4b713efSHong Zhang 
679d4002b98SHong Zhang     A = sell->A;  spA = (Mat_SeqSELL*)A->data;
680d4002b98SHong Zhang     B = sell->B;  spB = (Mat_SeqSELL*)B->data;
681f4b713efSHong Zhang     nz = spA->nz + spB->nz; /* total local nonzero entries of mat */
682f4b713efSHong Zhang     mem = nz*(sizeof(PetscScalar) + sizeof(PetscInt)) + 3*m*sizeof(PetscInt);
683f4b713efSHong Zhang     bcols = (PetscInt)(0.5*mem /(m*sizeof(PetscScalar)));
684f4b713efSHong Zhang     brows = 1000/bcols;
685f4b713efSHong Zhang     if (bcols > nis) bcols = nis;
686f4b713efSHong Zhang     if (brows == 0 || brows > m) brows = m;
687f4b713efSHong Zhang     c->brows = brows;
688f4b713efSHong Zhang     c->bcols = bcols;
689531e53bdSHong Zhang   } else { /* mpiaij matrix */
690531e53bdSHong Zhang     /* bcols is chosen s.t. dy-array takes 50% of local memory space as mat */
691531e53bdSHong Zhang     Mat_MPIAIJ *aij=(Mat_MPIAIJ*)mat->data;
692531e53bdSHong Zhang     Mat_SeqAIJ *spA,*spB;
693531e53bdSHong Zhang     Mat        A,B;
694a8971b87SHong Zhang     PetscInt   nz,brows,bcols;
695531e53bdSHong Zhang     PetscReal  mem;
696531e53bdSHong Zhang 
697531e53bdSHong Zhang     bs    = 1; /* only bs=1 is supported for MPIAIJ matrix */
698531e53bdSHong Zhang 
699531e53bdSHong Zhang     A = aij->A;  spA = (Mat_SeqAIJ*)A->data;
700531e53bdSHong Zhang     B = aij->B;  spB = (Mat_SeqAIJ*)B->data;
701531e53bdSHong Zhang     nz = spA->nz + spB->nz; /* total local nonzero entries of mat */
702531e53bdSHong Zhang     mem = nz*(sizeof(PetscScalar) + sizeof(PetscInt)) + 3*m*sizeof(PetscInt);
703531e53bdSHong Zhang     bcols = (PetscInt)(0.5*mem /(m*sizeof(PetscScalar)));
704531e53bdSHong Zhang     brows = 1000/bcols;
705531e53bdSHong Zhang     if (bcols > nis) bcols = nis;
706531e53bdSHong Zhang     if (brows == 0 || brows > m) brows = m;
707531e53bdSHong Zhang     c->brows = brows;
708531e53bdSHong Zhang     c->bcols = bcols;
709531e53bdSHong Zhang   }
710a8971b87SHong Zhang 
711a8971b87SHong Zhang   c->M       = mat->rmap->N/bs;         /* set the global rows and columns and local rows */
712a8971b87SHong Zhang   c->N       = mat->cmap->N/bs;
713a8971b87SHong Zhang   c->m       = mat->rmap->n/bs;
714a8971b87SHong Zhang   c->rstart  = mat->rmap->rstart/bs;
715a8971b87SHong Zhang   c->ncolors = nis;
71693dfae19SHong Zhang   PetscFunctionReturn(0);
71793dfae19SHong Zhang }
718bdaf1daeSBarry Smith 
719bdaf1daeSBarry Smith /*@C
720bdaf1daeSBarry Smith 
721bdaf1daeSBarry Smith     MatFDColoringSetValues - takes a matrix in compressed color format and enters the matrix into a PETSc Mat
722bdaf1daeSBarry Smith 
723bdaf1daeSBarry Smith    Collective on J
724bdaf1daeSBarry Smith 
725bdaf1daeSBarry Smith    Input Parameters:
726bdaf1daeSBarry Smith +    J - the sparse matrix
727bdaf1daeSBarry Smith .    coloring - created with MatFDColoringCreate() and a local coloring
728bdaf1daeSBarry Smith -    y - column major storage of matrix values with one color of values per column, the number of rows of y should match
729bdaf1daeSBarry Smith          the number of local rows of J and the number of columns is the number of colors.
730bdaf1daeSBarry Smith 
731bdaf1daeSBarry Smith    Level: intermediate
732bdaf1daeSBarry Smith 
7336a9b8d82SBarry Smith    Notes: the matrix in compressed color format may come from an Automatic Differentiation code
734bdaf1daeSBarry Smith 
735bdaf1daeSBarry Smith    The code will be slightly faster if MatFDColoringSetBlockSize(coloring,PETSC_DEFAULT,nc); is called immediately after creating the coloring
736bdaf1daeSBarry Smith 
737bdaf1daeSBarry Smith .seealso: MatFDColoringCreate(), ISColoring, ISColoringCreate(), ISColoringSetType(), IS_COLORING_LOCAL, MatFDColoringSetBlockSize()
738bdaf1daeSBarry Smith 
739bdaf1daeSBarry Smith @*/
740bdaf1daeSBarry Smith PetscErrorCode  MatFDColoringSetValues(Mat J,MatFDColoring coloring,const PetscScalar *y)
741bdaf1daeSBarry Smith {
742bdaf1daeSBarry Smith   PetscErrorCode    ierr;
743bdaf1daeSBarry Smith   MatEntry2         *Jentry2;
744bdaf1daeSBarry Smith   PetscInt          row,i,nrows_k,l,ncolors,nz = 0,bcols,nbcols = 0;
745bdaf1daeSBarry Smith   const PetscInt    *nrows;
746bdaf1daeSBarry Smith   PetscBool         eq;
747bdaf1daeSBarry Smith 
748bdaf1daeSBarry Smith   PetscFunctionBegin;
749bdaf1daeSBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
750bdaf1daeSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
751bdaf1daeSBarry Smith   ierr = PetscObjectCompareId((PetscObject)J,coloring->matid,&eq);CHKERRQ(ierr);
752*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!eq,PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONG,"Matrix used with MatFDColoringSetValues() must be that used with MatFDColoringCreate()");
753bdaf1daeSBarry Smith   Jentry2 = coloring->matentry2;
754bdaf1daeSBarry Smith   nrows   = coloring->nrows;
755bdaf1daeSBarry Smith   ncolors = coloring->ncolors;
756bdaf1daeSBarry Smith   bcols   = coloring->bcols;
757bdaf1daeSBarry Smith 
758bdaf1daeSBarry Smith   for (i=0; i<ncolors; i+=bcols) {
759bdaf1daeSBarry Smith     nrows_k = nrows[nbcols++];
760bdaf1daeSBarry Smith     for (l=0; l<nrows_k; l++) {
761bdaf1daeSBarry Smith       row                      = Jentry2[nz].row;   /* local row index */
762bdaf1daeSBarry Smith       *(Jentry2[nz++].valaddr) = y[row];
763bdaf1daeSBarry Smith     }
764bdaf1daeSBarry Smith     y += bcols*coloring->m;
765bdaf1daeSBarry Smith   }
766bdaf1daeSBarry Smith   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
767bdaf1daeSBarry Smith   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
768bdaf1daeSBarry Smith   PetscFunctionReturn(0);
769bdaf1daeSBarry Smith }
770