xref: /petsc/src/mat/impls/aij/mpi/fdmpiaij.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
69371c9d4SSatish Balay PetscErrorCode MatFDColoringApply_BAIJ(Mat J, MatFDColoring coloring, Vec x1, void *sctx) {
70d1c53f1SHong Zhang   PetscErrorCode (*f)(void *, Vec, Vec, void *) = (PetscErrorCode(*)(void *, Vec, Vec, void *))coloring->f;
80d1c53f1SHong Zhang   PetscInt           k, cstart, cend, l, row, col, nz, spidx, i, j;
95edff71fSBarry Smith   PetscScalar        dx = 0.0, *w3_array, *dy_i, *dy = coloring->dy;
100d1c53f1SHong Zhang   PetscScalar       *vscale_array;
115edff71fSBarry Smith   const PetscScalar *xx;
120d1c53f1SHong Zhang   PetscReal          epsilon = coloring->error_rel, umin = coloring->umin, unorm;
130d1c53f1SHong Zhang   Vec                w1 = coloring->w1, w2 = coloring->w2, w3, vscale = coloring->vscale;
140d1c53f1SHong Zhang   void              *fctx  = coloring->fctx;
150d1c53f1SHong Zhang   PetscInt           ctype = coloring->ctype, nxloc, nrows_k;
160d1c53f1SHong Zhang   PetscScalar       *valaddr;
170d1c53f1SHong Zhang   MatEntry          *Jentry  = coloring->matentry;
180df34763SHong Zhang   MatEntry2         *Jentry2 = coloring->matentry2;
190d1c53f1SHong Zhang   const PetscInt     ncolors = coloring->ncolors, *ncolumns = coloring->ncolumns, *nrows = coloring->nrows;
200d1c53f1SHong Zhang   PetscInt           bs = J->rmap->bs;
210d1c53f1SHong Zhang 
220d1c53f1SHong Zhang   PetscFunctionBegin;
239566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(x1, PETSC_TRUE));
240d1c53f1SHong Zhang   /* (1) Set w1 = F(x1) */
250d1c53f1SHong Zhang   if (!coloring->fset) {
269566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(MAT_FDColoringFunction, coloring, 0, 0, 0));
279566063dSJacob Faibussowitsch     PetscCall((*f)(sctx, x1, w1, fctx));
289566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(MAT_FDColoringFunction, coloring, 0, 0, 0));
290d1c53f1SHong Zhang   } else {
300d1c53f1SHong Zhang     coloring->fset = PETSC_FALSE;
310d1c53f1SHong Zhang   }
320d1c53f1SHong Zhang 
330d1c53f1SHong Zhang   /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */
349566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(x1, &nxloc));
350d1c53f1SHong Zhang   if (coloring->htype[0] == 'w') {
360d1c53f1SHong Zhang     /* vscale = dx is a constant scalar */
379566063dSJacob Faibussowitsch     PetscCall(VecNorm(x1, NORM_2, &unorm));
380d1c53f1SHong Zhang     dx = 1.0 / (PetscSqrtReal(1.0 + unorm) * epsilon);
390d1c53f1SHong Zhang   } else {
409566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(x1, &xx));
419566063dSJacob Faibussowitsch     PetscCall(VecGetArray(vscale, &vscale_array));
420d1c53f1SHong Zhang     for (col = 0; col < nxloc; col++) {
430d1c53f1SHong Zhang       dx = xx[col];
440d1c53f1SHong Zhang       if (PetscAbsScalar(dx) < umin) {
450d1c53f1SHong Zhang         if (PetscRealPart(dx) >= 0.0) dx = umin;
460d1c53f1SHong Zhang         else if (PetscRealPart(dx) < 0.0) dx = -umin;
470d1c53f1SHong Zhang       }
480d1c53f1SHong Zhang       dx *= epsilon;
490d1c53f1SHong Zhang       vscale_array[col] = 1.0 / dx;
500d1c53f1SHong Zhang     }
519566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(x1, &xx));
529566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(vscale, &vscale_array));
530d1c53f1SHong Zhang   }
54684f2004SHong Zhang   if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') {
559566063dSJacob Faibussowitsch     PetscCall(VecGhostUpdateBegin(vscale, INSERT_VALUES, SCATTER_FORWARD));
569566063dSJacob Faibussowitsch     PetscCall(VecGhostUpdateEnd(vscale, INSERT_VALUES, SCATTER_FORWARD));
570d1c53f1SHong Zhang   }
580d1c53f1SHong Zhang 
590d1c53f1SHong Zhang   /* (3) Loop over each color */
600d1c53f1SHong Zhang   if (!coloring->w3) {
619566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x1, &coloring->w3));
62ce911d59SRichard Tran Mills     /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */
639566063dSJacob Faibussowitsch     PetscCall(VecBindToCPU(coloring->w3, PETSC_TRUE));
649566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectParent((PetscObject)coloring, (PetscObject)coloring->w3));
650d1c53f1SHong Zhang   }
660d1c53f1SHong Zhang   w3 = coloring->w3;
670d1c53f1SHong Zhang 
689566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(x1, &cstart, &cend)); /* used by ghosted vscale */
69*48a46eb9SPierre Jolivet   if (vscale) PetscCall(VecGetArray(vscale, &vscale_array));
700d1c53f1SHong Zhang   nz = 0;
710d1c53f1SHong Zhang   for (k = 0; k < ncolors; k++) {
720d1c53f1SHong Zhang     coloring->currentcolor = k;
730d1c53f1SHong Zhang 
740d1c53f1SHong Zhang     /*
750d1c53f1SHong Zhang       (3-1) Loop over each column associated with color
760d1c53f1SHong Zhang       adding the perturbation to the vector w3 = x1 + dx.
770d1c53f1SHong Zhang     */
789566063dSJacob Faibussowitsch     PetscCall(VecCopy(x1, w3));
790d1c53f1SHong Zhang     dy_i = dy;
800d1c53f1SHong Zhang     for (i = 0; i < bs; i++) { /* Loop over a block of columns */
819566063dSJacob Faibussowitsch       PetscCall(VecGetArray(w3, &w3_array));
820d1c53f1SHong Zhang       if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */
830d1c53f1SHong Zhang       if (coloring->htype[0] == 'w') {
840d1c53f1SHong Zhang         for (l = 0; l < ncolumns[k]; l++) {
850d1c53f1SHong Zhang           col = i + bs * coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
860d1c53f1SHong Zhang           w3_array[col] += 1.0 / dx;
87684f2004SHong Zhang           if (i) w3_array[col - 1] -= 1.0 / dx; /* resume original w3[col-1] */
880d1c53f1SHong Zhang         }
890d1c53f1SHong Zhang       } else {                  /* htype == 'ds' */
900d1c53f1SHong Zhang         vscale_array -= cstart; /* shift pointer so global index can be used */
910d1c53f1SHong Zhang         for (l = 0; l < ncolumns[k]; l++) {
92f8c2866eSHong Zhang           col = i + bs * coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
930d1c53f1SHong Zhang           w3_array[col] += 1.0 / vscale_array[col];
94684f2004SHong Zhang           if (i) w3_array[col - 1] -= 1.0 / vscale_array[col - 1]; /* resume original w3[col-1] */
950d1c53f1SHong Zhang         }
960d1c53f1SHong Zhang         vscale_array += cstart;
970d1c53f1SHong Zhang       }
980d1c53f1SHong Zhang       if (ctype == IS_COLORING_GLOBAL) w3_array += cstart;
999566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(w3, &w3_array));
1000d1c53f1SHong Zhang 
1010d1c53f1SHong Zhang       /*
1020d1c53f1SHong Zhang        (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
1030d1c53f1SHong Zhang                            w2 = F(x1 + dx) - F(x1)
1040d1c53f1SHong Zhang        */
1059566063dSJacob Faibussowitsch       PetscCall(PetscLogEventBegin(MAT_FDColoringFunction, 0, 0, 0, 0));
1069566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(w2, dy_i)); /* place w2 to the array dy_i */
1079566063dSJacob Faibussowitsch       PetscCall((*f)(sctx, w3, w2, fctx));
1089566063dSJacob Faibussowitsch       PetscCall(PetscLogEventEnd(MAT_FDColoringFunction, 0, 0, 0, 0));
1099566063dSJacob Faibussowitsch       PetscCall(VecAXPY(w2, -1.0, w1));
1109566063dSJacob Faibussowitsch       PetscCall(VecResetArray(w2));
1110d1c53f1SHong Zhang       dy_i += nxloc; /* points to dy+i*nxloc */
1120d1c53f1SHong Zhang     }
1130d1c53f1SHong Zhang 
1140d1c53f1SHong Zhang     /*
1150d1c53f1SHong Zhang      (3-3) Loop over rows of vector, putting results into Jacobian matrix
1160d1c53f1SHong Zhang     */
1170d1c53f1SHong Zhang     nrows_k = nrows[k];
1180d1c53f1SHong Zhang     if (coloring->htype[0] == 'w') {
1190d1c53f1SHong Zhang       for (l = 0; l < nrows_k; l++) {
1200df34763SHong Zhang         row     = bs * Jentry2[nz].row; /* local row index */
1210df34763SHong Zhang         valaddr = Jentry2[nz++].valaddr;
1220d1c53f1SHong Zhang         spidx   = 0;
1230d1c53f1SHong Zhang         dy_i    = dy;
1240d1c53f1SHong Zhang         for (i = 0; i < bs; i++) {   /* column of the block */
1250d1c53f1SHong Zhang           for (j = 0; j < bs; j++) { /* row of the block */
1260d1c53f1SHong Zhang             valaddr[spidx++] = dy_i[row + j] * dx;
1270d1c53f1SHong Zhang           }
128f8c2866eSHong Zhang           dy_i += nxloc; /* points to dy+i*nxloc */
1290d1c53f1SHong Zhang         }
1300d1c53f1SHong Zhang       }
1310d1c53f1SHong Zhang     } else { /* htype == 'ds' */
1320d1c53f1SHong Zhang       for (l = 0; l < nrows_k; l++) {
1330d1c53f1SHong Zhang         row     = bs * Jentry[nz].row; /* local row index */
1340d1c53f1SHong Zhang         col     = bs * Jentry[nz].col; /* local column index */
135684f2004SHong Zhang         valaddr = Jentry[nz++].valaddr;
136f8c2866eSHong Zhang         spidx   = 0;
137f8c2866eSHong Zhang         dy_i    = dy;
138f8c2866eSHong Zhang         for (i = 0; i < bs; i++) {   /* column of the block */
139f8c2866eSHong Zhang           for (j = 0; j < bs; j++) { /* row of the block */
140f8c2866eSHong Zhang             valaddr[spidx++] = dy_i[row + j] * vscale_array[col + i];
141f8c2866eSHong Zhang           }
142f8c2866eSHong Zhang           dy_i += nxloc; /* points to dy+i*nxloc */
143f8c2866eSHong Zhang         }
1440d1c53f1SHong Zhang       }
1450d1c53f1SHong Zhang     }
1460d1c53f1SHong Zhang   }
1479566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
1489566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
149*48a46eb9SPierre Jolivet   if (vscale) PetscCall(VecRestoreArray(vscale, &vscale_array));
1500d1c53f1SHong Zhang 
1510d1c53f1SHong Zhang   coloring->currentcolor = -1;
1529566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(x1, PETSC_FALSE));
1530d1c53f1SHong Zhang   PetscFunctionReturn(0);
1540d1c53f1SHong Zhang }
155a64fbb32SBarry Smith 
156531c7667SBarry Smith /* this is declared PETSC_EXTERN because it is used by MatFDColoringUseDM() which is in the DM library */
1579371c9d4SSatish Balay PetscErrorCode MatFDColoringApply_AIJ(Mat J, MatFDColoring coloring, Vec x1, void *sctx) {
158fcd7ac73SHong Zhang   PetscErrorCode (*f)(void *, Vec, Vec, void *) = (PetscErrorCode(*)(void *, Vec, Vec, void *))coloring->f;
159a2f2d239SHong Zhang   PetscInt           k, cstart, cend, l, row, col, nz;
1605edff71fSBarry Smith   PetscScalar        dx = 0.0, *y, *w3_array;
1615edff71fSBarry Smith   const PetscScalar *xx;
162fcd7ac73SHong Zhang   PetscScalar       *vscale_array;
163fcd7ac73SHong Zhang   PetscReal          epsilon = coloring->error_rel, umin = coloring->umin, unorm;
1648bc97078SHong Zhang   Vec                w1 = coloring->w1, w2 = coloring->w2, w3, vscale = coloring->vscale;
165fcd7ac73SHong Zhang   void              *fctx  = coloring->fctx;
16691e7fa0fSBarry Smith   ISColoringType     ctype = coloring->ctype;
16791e7fa0fSBarry Smith   PetscInt           nxloc, nrows_k;
168573f477fSHong Zhang   MatEntry          *Jentry  = coloring->matentry;
1690df34763SHong Zhang   MatEntry2         *Jentry2 = coloring->matentry2;
1708bc97078SHong Zhang   const PetscInt     ncolors = coloring->ncolors, *ncolumns = coloring->ncolumns, *nrows = coloring->nrows;
171b294de21SRichard Tran Mills   PetscBool          alreadyboundtocpu;
172fcd7ac73SHong Zhang 
173fcd7ac73SHong Zhang   PetscFunctionBegin;
1749566063dSJacob Faibussowitsch   PetscCall(VecBoundToCPU(x1, &alreadyboundtocpu));
1759566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(x1, PETSC_TRUE));
17608401ef6SPierre Jolivet   PetscCheck(!(ctype == IS_COLORING_LOCAL) || !(J->ops->fdcoloringapply == MatFDColoringApply_AIJ), PetscObjectComm((PetscObject)J), PETSC_ERR_SUP, "Must call MatColoringUseDM() with IS_COLORING_LOCAL");
1778bc97078SHong Zhang   /* (1) Set w1 = F(x1) */
178fcd7ac73SHong Zhang   if (!coloring->fset) {
1799566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(MAT_FDColoringFunction, 0, 0, 0, 0));
1809566063dSJacob Faibussowitsch     PetscCall((*f)(sctx, x1, w1, fctx));
1819566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(MAT_FDColoringFunction, 0, 0, 0, 0));
182fcd7ac73SHong Zhang   } else {
183fcd7ac73SHong Zhang     coloring->fset = PETSC_FALSE;
184fcd7ac73SHong Zhang   }
185fcd7ac73SHong Zhang 
1868bc97078SHong Zhang   /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */
187f6af9589SHong Zhang   if (coloring->htype[0] == 'w') {
188684f2004SHong Zhang     /* vscale = 1./dx is a constant scalar */
1899566063dSJacob Faibussowitsch     PetscCall(VecNorm(x1, NORM_2, &unorm));
190c53567a0SHong Zhang     dx = 1.0 / (PetscSqrtReal(1.0 + unorm) * epsilon);
19170e7395fSHong Zhang   } else {
1929566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(x1, &nxloc));
1939566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(x1, &xx));
1949566063dSJacob Faibussowitsch     PetscCall(VecGetArray(vscale, &vscale_array));
19574d3cef9SHong Zhang     for (col = 0; col < nxloc; col++) {
196fcd7ac73SHong Zhang       dx = xx[col];
19774d3cef9SHong Zhang       if (PetscAbsScalar(dx) < umin) {
19874d3cef9SHong Zhang         if (PetscRealPart(dx) >= 0.0) dx = umin;
19974d3cef9SHong Zhang         else if (PetscRealPart(dx) < 0.0) dx = -umin;
20074d3cef9SHong Zhang       }
201fcd7ac73SHong Zhang       dx *= epsilon;
20274d3cef9SHong Zhang       vscale_array[col] = 1.0 / dx;
203f6af9589SHong Zhang     }
2049566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(x1, &xx));
2059566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(vscale, &vscale_array));
20670e7395fSHong Zhang   }
207684f2004SHong Zhang   if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') {
2089566063dSJacob Faibussowitsch     PetscCall(VecGhostUpdateBegin(vscale, INSERT_VALUES, SCATTER_FORWARD));
2099566063dSJacob Faibussowitsch     PetscCall(VecGhostUpdateEnd(vscale, INSERT_VALUES, SCATTER_FORWARD));
210fcd7ac73SHong Zhang   }
211fcd7ac73SHong Zhang 
2128bc97078SHong Zhang   /* (3) Loop over each color */
2138bc97078SHong Zhang   if (!coloring->w3) {
2149566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x1, &coloring->w3));
2159566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectParent((PetscObject)coloring, (PetscObject)coloring->w3));
2168bc97078SHong Zhang   }
2178bc97078SHong Zhang   w3 = coloring->w3;
218fcd7ac73SHong Zhang 
2199566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(x1, &cstart, &cend)); /* used by ghosted vscale */
220*48a46eb9SPierre Jolivet   if (vscale) PetscCall(VecGetArray(vscale, &vscale_array));
2218bc97078SHong Zhang   nz = 0;
222e2c857f8SHong Zhang 
223b3e1f37bSHong Zhang   if (coloring->bcols > 1) { /* use blocked insertion of Jentry */
224b3e1f37bSHong Zhang     PetscInt     i, m = J->rmap->n, nbcols, bcols = coloring->bcols;
225e2c857f8SHong Zhang     PetscScalar *dy = coloring->dy, *dy_k;
226e2c857f8SHong Zhang 
227e2c857f8SHong Zhang     nbcols = 0;
228e2c857f8SHong Zhang     for (k = 0; k < ncolors; k += bcols) {
229e2c857f8SHong Zhang       /*
230e2c857f8SHong Zhang        (3-1) Loop over each column associated with color
231e2c857f8SHong Zhang        adding the perturbation to the vector w3 = x1 + dx.
232e2c857f8SHong Zhang        */
233e2c857f8SHong Zhang 
234e2c857f8SHong Zhang       dy_k = dy;
235e2c857f8SHong Zhang       if (k + bcols > ncolors) bcols = ncolors - k;
236e2c857f8SHong Zhang       for (i = 0; i < bcols; i++) {
237ceef8a49SBarry Smith         coloring->currentcolor = k + i;
238e2c857f8SHong Zhang 
2399566063dSJacob Faibussowitsch         PetscCall(VecCopy(x1, w3));
2409566063dSJacob Faibussowitsch         PetscCall(VecGetArray(w3, &w3_array));
241e2c857f8SHong Zhang         if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */
242e2c857f8SHong Zhang         if (coloring->htype[0] == 'w') {
243e2c857f8SHong Zhang           for (l = 0; l < ncolumns[k + i]; l++) {
244e2c857f8SHong Zhang             col = coloring->columns[k + i][l]; /* local column (in global index!) of the matrix we are probing for */
245e2c857f8SHong Zhang             w3_array[col] += 1.0 / dx;
246e2c857f8SHong Zhang           }
247e2c857f8SHong Zhang         } else {                  /* htype == 'ds' */
248e2c857f8SHong Zhang           vscale_array -= cstart; /* shift pointer so global index can be used */
249e2c857f8SHong Zhang           for (l = 0; l < ncolumns[k + i]; l++) {
250e2c857f8SHong Zhang             col = coloring->columns[k + i][l]; /* local column (in global index!) of the matrix we are probing for */
251e2c857f8SHong Zhang             w3_array[col] += 1.0 / vscale_array[col];
252e2c857f8SHong Zhang           }
253e2c857f8SHong Zhang           vscale_array += cstart;
254e2c857f8SHong Zhang         }
255e2c857f8SHong Zhang         if (ctype == IS_COLORING_GLOBAL) w3_array += cstart;
2569566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(w3, &w3_array));
257e2c857f8SHong Zhang 
258e2c857f8SHong Zhang         /*
259e2c857f8SHong Zhang          (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
260e2c857f8SHong Zhang                            w2 = F(x1 + dx) - F(x1)
261e2c857f8SHong Zhang          */
2629566063dSJacob Faibussowitsch         PetscCall(PetscLogEventBegin(MAT_FDColoringFunction, 0, 0, 0, 0));
2639566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(w2, dy_k)); /* place w2 to the array dy_i */
2649566063dSJacob Faibussowitsch         PetscCall((*f)(sctx, w3, w2, fctx));
2659566063dSJacob Faibussowitsch         PetscCall(PetscLogEventEnd(MAT_FDColoringFunction, 0, 0, 0, 0));
2669566063dSJacob Faibussowitsch         PetscCall(VecAXPY(w2, -1.0, w1));
2679566063dSJacob Faibussowitsch         PetscCall(VecResetArray(w2));
268e2c857f8SHong Zhang         dy_k += m; /* points to dy+i*nxloc */
269e2c857f8SHong Zhang       }
270e2c857f8SHong Zhang 
271e2c857f8SHong Zhang       /*
272e2c857f8SHong Zhang        (3-3) Loop over block rows of vector, putting results into Jacobian matrix
273e2c857f8SHong Zhang        */
274d880da65SHong Zhang       nrows_k = nrows[nbcols++];
275d880da65SHong Zhang 
276e2c857f8SHong Zhang       if (coloring->htype[0] == 'w') {
277e2c857f8SHong Zhang         for (l = 0; l < nrows_k; l++) {
2780df34763SHong Zhang           row = Jentry2[nz].row; /* local row index */
2795059b303SJunchao 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
2805059b303SJunchao Zhang              another way, and it seems work. See https://lists.mcs.anl.gov/pipermail/petsc-users/2021-December/045158.html
2815059b303SJunchao Zhang            */
2825059b303SJunchao Zhang #if defined(PETSC_USE_COMPLEX)
2835059b303SJunchao Zhang           PetscScalar *tmp = Jentry2[nz].valaddr;
2845059b303SJunchao Zhang           *tmp             = dy[row] * dx;
2855059b303SJunchao Zhang #else
2865059b303SJunchao Zhang           *(Jentry2[nz].valaddr) = dy[row] * dx;
2875059b303SJunchao Zhang #endif
2885059b303SJunchao Zhang           nz++;
289e2c857f8SHong Zhang         }
290e2c857f8SHong Zhang       } else { /* htype == 'ds' */
291e2c857f8SHong Zhang         for (l = 0; l < nrows_k; l++) {
292e2c857f8SHong Zhang           row = Jentry[nz].row; /* local row index */
2935059b303SJunchao Zhang #if defined(PETSC_USE_COMPLEX)  /* See https://lists.mcs.anl.gov/pipermail/petsc-users/2021-December/045158.html */
2945059b303SJunchao Zhang           PetscScalar *tmp = Jentry[nz].valaddr;
2955059b303SJunchao Zhang           *tmp             = dy[row] * vscale_array[Jentry[nz].col];
2965059b303SJunchao Zhang #else
297d880da65SHong Zhang           *(Jentry[nz].valaddr)  = dy[row] * vscale_array[Jentry[nz].col];
2985059b303SJunchao Zhang #endif
299e2c857f8SHong Zhang           nz++;
300e2c857f8SHong Zhang         }
301e2c857f8SHong Zhang       }
302e2c857f8SHong Zhang     }
303b3e1f37bSHong Zhang   } else { /* bcols == 1 */
3048bc97078SHong Zhang     for (k = 0; k < ncolors; k++) {
3058bc97078SHong Zhang       coloring->currentcolor = k;
3069e917edbSHong Zhang 
307fcd7ac73SHong Zhang       /*
3088bc97078SHong Zhang        (3-1) Loop over each column associated with color
309f6af9589SHong Zhang        adding the perturbation to the vector w3 = x1 + dx.
310fcd7ac73SHong Zhang        */
3119566063dSJacob Faibussowitsch       PetscCall(VecCopy(x1, w3));
3129566063dSJacob Faibussowitsch       PetscCall(VecGetArray(w3, &w3_array));
313a2f2d239SHong Zhang       if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */
314b039d6c4SHong Zhang       if (coloring->htype[0] == 'w') {
3158bc97078SHong Zhang         for (l = 0; l < ncolumns[k]; l++) {
316f6af9589SHong Zhang           col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
3179e917edbSHong Zhang           w3_array[col] += 1.0 / dx;
3189e917edbSHong Zhang         }
319b039d6c4SHong Zhang       } else {                  /* htype == 'ds' */
320a2f2d239SHong Zhang         vscale_array -= cstart; /* shift pointer so global index can be used */
321b039d6c4SHong Zhang         for (l = 0; l < ncolumns[k]; l++) {
322b039d6c4SHong Zhang           col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
323a2f2d239SHong Zhang           w3_array[col] += 1.0 / vscale_array[col];
32470e7395fSHong Zhang         }
325a2f2d239SHong Zhang         vscale_array += cstart;
326fcd7ac73SHong Zhang       }
327a2f2d239SHong Zhang       if (ctype == IS_COLORING_GLOBAL) w3_array += cstart;
3289566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(w3, &w3_array));
3299e917edbSHong Zhang 
330fcd7ac73SHong Zhang       /*
3318bc97078SHong Zhang        (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
332fcd7ac73SHong Zhang                            w2 = F(x1 + dx) - F(x1)
333fcd7ac73SHong Zhang        */
3349566063dSJacob Faibussowitsch       PetscCall(PetscLogEventBegin(MAT_FDColoringFunction, 0, 0, 0, 0));
3359566063dSJacob Faibussowitsch       PetscCall((*f)(sctx, w3, w2, fctx));
3369566063dSJacob Faibussowitsch       PetscCall(PetscLogEventEnd(MAT_FDColoringFunction, 0, 0, 0, 0));
3379566063dSJacob Faibussowitsch       PetscCall(VecAXPY(w2, -1.0, w1));
3389e917edbSHong Zhang 
3398bc97078SHong Zhang       /*
3408bc97078SHong Zhang        (3-3) Loop over rows of vector, putting results into Jacobian matrix
3418bc97078SHong Zhang        */
342c53567a0SHong Zhang       nrows_k = nrows[k];
3439566063dSJacob Faibussowitsch       PetscCall(VecGetArray(w2, &y));
344b039d6c4SHong Zhang       if (coloring->htype[0] == 'w') {
345c53567a0SHong Zhang         for (l = 0; l < nrows_k; l++) {
3460df34763SHong Zhang           row = Jentry2[nz].row; /* local row index */
3475059b303SJunchao Zhang #if defined(PETSC_USE_COMPLEX)   /* See https://lists.mcs.anl.gov/pipermail/petsc-users/2021-December/045158.html */
3485059b303SJunchao Zhang           PetscScalar *tmp = Jentry2[nz].valaddr;
3495059b303SJunchao Zhang           *tmp             = y[row] * dx;
3505059b303SJunchao Zhang #else
3515059b303SJunchao Zhang           *(Jentry2[nz].valaddr) = y[row] * dx;
3525059b303SJunchao Zhang #endif
3535059b303SJunchao Zhang           nz++;
3549e917edbSHong Zhang         }
355b039d6c4SHong Zhang       } else { /* htype == 'ds' */
356c53567a0SHong Zhang         for (l = 0; l < nrows_k; l++) {
357b039d6c4SHong Zhang           row = Jentry[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 = Jentry[nz].valaddr;
3605059b303SJunchao Zhang           *tmp             = y[row] * vscale_array[Jentry[nz].col];
3615059b303SJunchao Zhang #else
362b039d6c4SHong Zhang           *(Jentry[nz].valaddr)  = y[row] * vscale_array[Jentry[nz].col];
3635059b303SJunchao Zhang #endif
364573f477fSHong Zhang           nz++;
365fcd7ac73SHong Zhang         }
366b039d6c4SHong Zhang       }
3679566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(w2, &y));
3688bc97078SHong Zhang     }
369b3e1f37bSHong Zhang   }
370b3e1f37bSHong Zhang 
371e2cf4d64SStefano Zampini #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
372c70f7ee4SJunchao Zhang   if (J->offloadmask != PETSC_OFFLOAD_UNALLOCATED) J->offloadmask = PETSC_OFFLOAD_CPU;
373e2cf4d64SStefano Zampini #endif
3749566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
3759566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
376*48a46eb9SPierre Jolivet   if (vscale) PetscCall(VecRestoreArray(vscale, &vscale_array));
377fcd7ac73SHong Zhang   coloring->currentcolor = -1;
3789566063dSJacob Faibussowitsch   if (!alreadyboundtocpu) PetscCall(VecBindToCPU(x1, PETSC_FALSE));
379fcd7ac73SHong Zhang   PetscFunctionReturn(0);
380fcd7ac73SHong Zhang }
381fcd7ac73SHong Zhang 
3829371c9d4SSatish Balay PetscErrorCode MatFDColoringSetUp_MPIXAIJ(Mat mat, ISColoring iscoloring, MatFDColoring c) {
383fcd7ac73SHong Zhang   PetscMPIInt            size, *ncolsonproc, *disp, nn;
384e3225e9dSHong Zhang   PetscInt               i, n, nrows, nrows_i, j, k, m, ncols, col, *rowhit, cstart, cend, colb;
38572c15787SHong Zhang   const PetscInt        *is, *A_ci, *A_cj, *B_ci, *B_cj, *row = NULL, *ltog = NULL;
386071fcb05SBarry Smith   PetscInt               nis = iscoloring->n, nctot, *cols, tmp = 0;
387fcd7ac73SHong Zhang   ISLocalToGlobalMapping map   = mat->cmap->mapping;
388a8971b87SHong Zhang   PetscInt               ctype = c->ctype, *spidxA, *spidxB, nz, bs, bs2, spidx;
3890d1c53f1SHong Zhang   Mat                    A, B;
390e3225e9dSHong Zhang   PetscScalar           *A_val, *B_val, **valaddrhit;
391573f477fSHong Zhang   MatEntry              *Jentry;
3920df34763SHong Zhang   MatEntry2             *Jentry2;
393d4002b98SHong Zhang   PetscBool              isBAIJ, isSELL;
394531e53bdSHong Zhang   PetscInt               bcols = c->bcols;
3950d1c53f1SHong Zhang #if defined(PETSC_USE_CTABLE)
3960d1c53f1SHong Zhang   PetscTable colmap = NULL;
3970d1c53f1SHong Zhang #else
3980d1c53f1SHong Zhang   PetscInt *colmap = NULL;      /* local col number of off-diag col */
3990d1c53f1SHong Zhang #endif
400fcd7ac73SHong Zhang 
401fcd7ac73SHong Zhang   PetscFunctionBegin;
4025bdb020cSBarry Smith   if (ctype == IS_COLORING_LOCAL) {
40328b400f6SJacob Faibussowitsch     PetscCheck(map, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_INCOMP, "When using ghosted differencing matrix must have local to global mapping provided with MatSetLocalToGlobalMapping");
4049566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetIndices(map, &ltog));
40572c15787SHong Zhang   }
406fcd7ac73SHong Zhang 
4079566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(mat, &bs));
4089566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompare((PetscObject)mat, MATMPIBAIJ, &isBAIJ));
4099566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATMPISELL, &isSELL));
410e3225e9dSHong Zhang   if (isBAIJ) {
411195687c5SHong Zhang     Mat_MPIBAIJ *baij = (Mat_MPIBAIJ *)mat->data;
4126c145f34SHong Zhang     Mat_SeqBAIJ *spA, *spB;
4139371c9d4SSatish Balay     A     = baij->A;
4149371c9d4SSatish Balay     spA   = (Mat_SeqBAIJ *)A->data;
4159371c9d4SSatish Balay     A_val = spA->a;
4169371c9d4SSatish Balay     B     = baij->B;
4179371c9d4SSatish Balay     spB   = (Mat_SeqBAIJ *)B->data;
4189371c9d4SSatish Balay     B_val = spB->a;
4190d1c53f1SHong Zhang     nz    = spA->nz + spB->nz; /* total nonzero entries of mat */
420*48a46eb9SPierre Jolivet     if (!baij->colmap) PetscCall(MatCreateColmap_MPIBAIJ_Private(mat));
421097e26fcSJed Brown     colmap = baij->colmap;
4229566063dSJacob Faibussowitsch     PetscCall(MatGetColumnIJ_SeqBAIJ_Color(A, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &A_ci, &A_cj, &spidxA, NULL));
4239566063dSJacob Faibussowitsch     PetscCall(MatGetColumnIJ_SeqBAIJ_Color(B, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &B_ci, &B_cj, &spidxB, NULL));
424195687c5SHong Zhang 
425195687c5SHong Zhang     if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') { /* create vscale for storing dx */
426195687c5SHong Zhang       PetscInt *garray;
4279566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(B->cmap->n, &garray));
428195687c5SHong Zhang       for (i = 0; i < baij->B->cmap->n / bs; i++) {
4299371c9d4SSatish Balay         for (j = 0; j < bs; j++) { garray[i * bs + j] = bs * baij->garray[i] + j; }
430195687c5SHong Zhang       }
4319566063dSJacob Faibussowitsch       PetscCall(VecCreateGhost(PetscObjectComm((PetscObject)mat), mat->cmap->n, PETSC_DETERMINE, B->cmap->n, garray, &c->vscale));
4329566063dSJacob Faibussowitsch       PetscCall(VecBindToCPU(c->vscale, PETSC_TRUE));
4339566063dSJacob Faibussowitsch       PetscCall(PetscFree(garray));
434195687c5SHong Zhang     }
435d4002b98SHong Zhang   } else if (isSELL) {
436d4002b98SHong Zhang     Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
437d4002b98SHong Zhang     Mat_SeqSELL *spA, *spB;
4389371c9d4SSatish Balay     A     = sell->A;
4399371c9d4SSatish Balay     spA   = (Mat_SeqSELL *)A->data;
4409371c9d4SSatish Balay     A_val = spA->val;
4419371c9d4SSatish Balay     B     = sell->B;
4429371c9d4SSatish Balay     spB   = (Mat_SeqSELL *)B->data;
4439371c9d4SSatish Balay     B_val = spB->val;
444f4b713efSHong Zhang     nz    = spA->nz + spB->nz; /* total nonzero entries of mat */
445d4002b98SHong Zhang     if (!sell->colmap) {
446f4b713efSHong Zhang       /* Allow access to data structures of local part of matrix
447f4b713efSHong Zhang        - creates aij->colmap which maps global column number to local number in part B */
4489566063dSJacob Faibussowitsch       PetscCall(MatCreateColmap_MPISELL_Private(mat));
449f4b713efSHong Zhang     }
450d4002b98SHong Zhang     colmap = sell->colmap;
4519566063dSJacob Faibussowitsch     PetscCall(MatGetColumnIJ_SeqSELL_Color(A, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &A_ci, &A_cj, &spidxA, NULL));
4529566063dSJacob Faibussowitsch     PetscCall(MatGetColumnIJ_SeqSELL_Color(B, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &B_ci, &B_cj, &spidxB, NULL));
453f4b713efSHong Zhang 
454f4b713efSHong Zhang     bs = 1; /* only bs=1 is supported for non MPIBAIJ matrix */
455f4b713efSHong Zhang 
456f4b713efSHong Zhang     if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') { /* create vscale for storing dx */
4579566063dSJacob Faibussowitsch       PetscCall(VecCreateGhost(PetscObjectComm((PetscObject)mat), mat->cmap->n, PETSC_DETERMINE, B->cmap->n, sell->garray, &c->vscale));
4589566063dSJacob Faibussowitsch       PetscCall(VecBindToCPU(c->vscale, PETSC_TRUE));
459f4b713efSHong Zhang     }
460e3225e9dSHong Zhang   } else {
461e3225e9dSHong Zhang     Mat_MPIAIJ *aij = (Mat_MPIAIJ *)mat->data;
462e3225e9dSHong Zhang     Mat_SeqAIJ *spA, *spB;
4639371c9d4SSatish Balay     A     = aij->A;
4649371c9d4SSatish Balay     spA   = (Mat_SeqAIJ *)A->data;
4659371c9d4SSatish Balay     A_val = spA->a;
4669371c9d4SSatish Balay     B     = aij->B;
4679371c9d4SSatish Balay     spB   = (Mat_SeqAIJ *)B->data;
4689371c9d4SSatish Balay     B_val = spB->a;
469e3225e9dSHong Zhang     nz    = spA->nz + spB->nz; /* total nonzero entries of mat */
470e3225e9dSHong Zhang     if (!aij->colmap) {
471e3225e9dSHong Zhang       /* Allow access to data structures of local part of matrix
472e3225e9dSHong Zhang        - creates aij->colmap which maps global column number to local number in part B */
4739566063dSJacob Faibussowitsch       PetscCall(MatCreateColmap_MPIAIJ_Private(mat));
474e3225e9dSHong Zhang     }
475097e26fcSJed Brown     colmap = aij->colmap;
4769566063dSJacob Faibussowitsch     PetscCall(MatGetColumnIJ_SeqAIJ_Color(A, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &A_ci, &A_cj, &spidxA, NULL));
4779566063dSJacob Faibussowitsch     PetscCall(MatGetColumnIJ_SeqAIJ_Color(B, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &B_ci, &B_cj, &spidxB, NULL));
478e3225e9dSHong Zhang 
479e3225e9dSHong Zhang     bs = 1; /* only bs=1 is supported for non MPIBAIJ matrix */
480195687c5SHong Zhang 
481195687c5SHong Zhang     if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') { /* create vscale for storing dx */
4829566063dSJacob Faibussowitsch       PetscCall(VecCreateGhost(PetscObjectComm((PetscObject)mat), mat->cmap->n, PETSC_DETERMINE, B->cmap->n, aij->garray, &c->vscale));
4839566063dSJacob Faibussowitsch       PetscCall(VecBindToCPU(c->vscale, PETSC_TRUE));
484195687c5SHong Zhang     }
4850d1c53f1SHong Zhang   }
4860d1c53f1SHong Zhang 
4870d1c53f1SHong Zhang   m      = mat->rmap->n / bs;
4880d1c53f1SHong Zhang   cstart = mat->cmap->rstart / bs;
4890d1c53f1SHong Zhang   cend   = mat->cmap->rend / bs;
490fcd7ac73SHong Zhang 
4919566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nis, &c->ncolumns, nis, &c->columns));
4929566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nis, &c->nrows));
4939566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectMemory((PetscObject)c, 3 * nis * sizeof(PetscInt)));
494fcd7ac73SHong Zhang 
4950df34763SHong Zhang   if (c->htype[0] == 'd') {
4969566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz, &Jentry));
4979566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectMemory((PetscObject)c, nz * sizeof(MatEntry)));
4988bc97078SHong Zhang     c->matentry = Jentry;
4990df34763SHong Zhang   } else if (c->htype[0] == 'w') {
5009566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz, &Jentry2));
5019566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectMemory((PetscObject)c, nz * sizeof(MatEntry2)));
5020df34763SHong Zhang     c->matentry2 = Jentry2;
5030df34763SHong Zhang   } else SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "htype is not supported");
504a774a6f1SHong Zhang 
5059566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(m + 1, &rowhit, m + 1, &valaddrhit));
506fcd7ac73SHong Zhang   nz = 0;
5079566063dSJacob Faibussowitsch   PetscCall(ISColoringGetIS(iscoloring, PETSC_OWN_POINTER, PETSC_IGNORE, &c->isa));
508071fcb05SBarry Smith 
509071fcb05SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
5109566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &size));
5119566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(size, &ncolsonproc, size, &disp));
512071fcb05SBarry Smith   }
513071fcb05SBarry Smith 
514d3825b63SHong Zhang   for (i = 0; i < nis; i++) { /* for each local color */
5159566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(c->isa[i], &n));
5169566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(c->isa[i], &is));
517fcd7ac73SHong Zhang 
518fcd7ac73SHong Zhang     c->ncolumns[i] = n; /* local number of columns of this color on this process */
519071fcb05SBarry Smith     c->columns[i]  = (PetscInt *)is;
520fcd7ac73SHong Zhang 
521fcd7ac73SHong Zhang     if (ctype == IS_COLORING_GLOBAL) {
522d3825b63SHong Zhang       /* Determine nctot, the total (parallel) number of columns of this color */
523d3825b63SHong Zhang       /* ncolsonproc[j]: local ncolumns on proc[j] of this color */
5249566063dSJacob Faibussowitsch       PetscCall(PetscMPIIntCast(n, &nn));
5259566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Allgather(&nn, 1, MPI_INT, ncolsonproc, 1, MPI_INT, PetscObjectComm((PetscObject)mat)));
5269371c9d4SSatish Balay       nctot = 0;
5279371c9d4SSatish Balay       for (j = 0; j < size; j++) nctot += ncolsonproc[j];
528*48a46eb9SPierre Jolivet       if (!nctot) PetscCall(PetscInfo(mat, "Coloring of matrix has some unneeded colors with no corresponding rows\n"));
529fcd7ac73SHong Zhang 
530fcd7ac73SHong Zhang       disp[0] = 0;
5319371c9d4SSatish Balay       for (j = 1; j < size; j++) { disp[j] = disp[j - 1] + ncolsonproc[j - 1]; }
532fcd7ac73SHong Zhang 
533d3825b63SHong Zhang       /* Get cols, the complete list of columns for this color on each process */
5349566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nctot + 1, &cols));
5359566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Allgatherv((void *)is, n, MPIU_INT, cols, ncolsonproc, disp, MPIU_INT, PetscObjectComm((PetscObject)mat)));
5365bdb020cSBarry Smith     } else if (ctype == IS_COLORING_LOCAL) {
537fcd7ac73SHong Zhang       /* Determine local number of columns of this color on this process, including ghost points */
538fcd7ac73SHong Zhang       nctot = n;
539071fcb05SBarry Smith       cols  = (PetscInt *)is;
540fcd7ac73SHong Zhang     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not provided for this MatFDColoring type");
541fcd7ac73SHong Zhang 
5421b97d346SHong Zhang     /* Mark all rows affect by these columns */
5439566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(rowhit, m));
5440d1c53f1SHong Zhang     bs2     = bs * bs;
5454b2e90caSHong Zhang     nrows_i = 0;
5461b97d346SHong Zhang     for (j = 0; j < nctot; j++) { /* loop over columns*/
5475bdb020cSBarry Smith       if (ctype == IS_COLORING_LOCAL) {
548fcd7ac73SHong Zhang         col = ltog[cols[j]];
549fcd7ac73SHong Zhang       } else {
550fcd7ac73SHong Zhang         col = cols[j];
551fcd7ac73SHong Zhang       }
552e3225e9dSHong Zhang       if (col >= cstart && col < cend) { /* column is in A, diagonal block of mat */
553071fcb05SBarry Smith         tmp   = A_ci[col - cstart];
554071fcb05SBarry Smith         row   = A_cj + tmp;
555071fcb05SBarry Smith         nrows = A_ci[col - cstart + 1] - tmp;
5564b2e90caSHong Zhang         nrows_i += nrows;
557071fcb05SBarry Smith 
558fcd7ac73SHong Zhang         /* loop over columns of A marking them in rowhit */
559d3825b63SHong Zhang         for (k = 0; k < nrows; k++) {
560d3825b63SHong Zhang           /* set valaddrhit for part A */
561071fcb05SBarry Smith           spidx            = bs2 * spidxA[tmp + k];
562e3225e9dSHong Zhang           valaddrhit[*row] = &A_val[spidx];
563a774a6f1SHong Zhang           rowhit[*row++]   = col - cstart + 1; /* local column index */
564fcd7ac73SHong Zhang         }
565e3225e9dSHong Zhang       } else { /* column is in B, off-diagonal block of mat */
566fcd7ac73SHong Zhang #if defined(PETSC_USE_CTABLE)
5679566063dSJacob Faibussowitsch         PetscCall(PetscTableFind(colmap, col + 1, &colb));
568fcd7ac73SHong Zhang         colb--;
569fcd7ac73SHong Zhang #else
5700d1c53f1SHong Zhang         colb = colmap[col] - 1; /* local column index */
571fcd7ac73SHong Zhang #endif
572fcd7ac73SHong Zhang         if (colb == -1) {
573d3825b63SHong Zhang           nrows = 0;
574fcd7ac73SHong Zhang         } else {
5750d1c53f1SHong Zhang           colb  = colb / bs;
576071fcb05SBarry Smith           tmp   = B_ci[colb];
577071fcb05SBarry Smith           row   = B_cj + tmp;
578071fcb05SBarry Smith           nrows = B_ci[colb + 1] - tmp;
579fcd7ac73SHong Zhang         }
5804b2e90caSHong Zhang         nrows_i += nrows;
581fcd7ac73SHong Zhang         /* loop over columns of B marking them in rowhit */
582d3825b63SHong Zhang         for (k = 0; k < nrows; k++) {
583d3825b63SHong Zhang           /* set valaddrhit for part B */
584071fcb05SBarry Smith           spidx            = bs2 * spidxB[tmp + k];
585e3225e9dSHong Zhang           valaddrhit[*row] = &B_val[spidx];
58670e7395fSHong Zhang           rowhit[*row++]   = colb + 1 + cend - cstart; /* local column index */
587fcd7ac73SHong Zhang         }
588fcd7ac73SHong Zhang       }
589fcd7ac73SHong Zhang     }
5904b2e90caSHong Zhang     c->nrows[i] = nrows_i;
5918bc97078SHong Zhang 
5920df34763SHong Zhang     if (c->htype[0] == 'd') {
593d3825b63SHong Zhang       for (j = 0; j < m; j++) {
594fcd7ac73SHong Zhang         if (rowhit[j]) {
595573f477fSHong Zhang           Jentry[nz].row     = j;             /* local row index */
596573f477fSHong Zhang           Jentry[nz].col     = rowhit[j] - 1; /* local column index */
597573f477fSHong Zhang           Jentry[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */
598573f477fSHong Zhang           nz++;
599fcd7ac73SHong Zhang         }
600fcd7ac73SHong Zhang       }
6010df34763SHong Zhang     } else { /* c->htype == 'wp' */
6020df34763SHong Zhang       for (j = 0; j < m; j++) {
6030df34763SHong Zhang         if (rowhit[j]) {
6040df34763SHong Zhang           Jentry2[nz].row     = j;             /* local row index */
6050df34763SHong Zhang           Jentry2[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */
6060df34763SHong Zhang           nz++;
6070df34763SHong Zhang         }
6080df34763SHong Zhang       }
6090df34763SHong Zhang     }
610*48a46eb9SPierre Jolivet     if (ctype == IS_COLORING_GLOBAL) PetscCall(PetscFree(cols));
611fcd7ac73SHong Zhang   }
612*48a46eb9SPierre Jolivet   if (ctype == IS_COLORING_GLOBAL) PetscCall(PetscFree2(ncolsonproc, disp));
613fcd7ac73SHong Zhang 
614a8971b87SHong Zhang   if (bcols > 1) { /* reorder Jentry for faster MatFDColoringApply() */
6159566063dSJacob Faibussowitsch     PetscCall(MatFDColoringSetUpBlocked_AIJ_Private(mat, c, nz));
616531e53bdSHong Zhang   }
617531e53bdSHong Zhang 
6180d1c53f1SHong Zhang   if (isBAIJ) {
6199566063dSJacob Faibussowitsch     PetscCall(MatRestoreColumnIJ_SeqBAIJ_Color(A, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &A_ci, &A_cj, &spidxA, NULL));
6209566063dSJacob Faibussowitsch     PetscCall(MatRestoreColumnIJ_SeqBAIJ_Color(B, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &B_ci, &B_cj, &spidxB, NULL));
6219566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(bs * mat->rmap->n, &c->dy));
622d4002b98SHong Zhang   } else if (isSELL) {
6239566063dSJacob Faibussowitsch     PetscCall(MatRestoreColumnIJ_SeqSELL_Color(A, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &A_ci, &A_cj, &spidxA, NULL));
6249566063dSJacob Faibussowitsch     PetscCall(MatRestoreColumnIJ_SeqSELL_Color(B, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &B_ci, &B_cj, &spidxB, NULL));
6250d1c53f1SHong Zhang   } else {
6269566063dSJacob Faibussowitsch     PetscCall(MatRestoreColumnIJ_SeqAIJ_Color(A, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &A_ci, &A_cj, &spidxA, NULL));
6279566063dSJacob Faibussowitsch     PetscCall(MatRestoreColumnIJ_SeqAIJ_Color(B, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &B_ci, &B_cj, &spidxB, NULL));
6280d1c53f1SHong Zhang   }
6290d1c53f1SHong Zhang 
6309566063dSJacob Faibussowitsch   PetscCall(ISColoringRestoreIS(iscoloring, PETSC_OWN_POINTER, &c->isa));
6319566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rowhit, valaddrhit));
632a8971b87SHong Zhang 
633*48a46eb9SPierre Jolivet   if (ctype == IS_COLORING_LOCAL) PetscCall(ISLocalToGlobalMappingRestoreIndices(map, &ltog));
6349566063dSJacob Faibussowitsch   PetscCall(PetscInfo(c, "ncolors %" PetscInt_FMT ", brows %" PetscInt_FMT " and bcols %" PetscInt_FMT " are used.\n", c->ncolors, c->brows, c->bcols));
635fcd7ac73SHong Zhang   PetscFunctionReturn(0);
636fcd7ac73SHong Zhang }
63793dfae19SHong Zhang 
6389371c9d4SSatish Balay PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat mat, ISColoring iscoloring, MatFDColoring c) {
639a8971b87SHong Zhang   PetscInt  bs, nis = iscoloring->n, m = mat->rmap->n;
640d4002b98SHong Zhang   PetscBool isBAIJ, isSELL;
641531e53bdSHong Zhang 
64293dfae19SHong Zhang   PetscFunctionBegin;
643531e53bdSHong Zhang   /* set default brows and bcols for speedup inserting the dense matrix into sparse Jacobian;
644531e53bdSHong Zhang    bcols is chosen s.t. dy-array takes 50% of memory space as mat */
6459566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(mat, &bs));
6469566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompare((PetscObject)mat, MATMPIBAIJ, &isBAIJ));
6479566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATMPISELL, &isSELL));
648b8dfa574SBarry Smith   if (isBAIJ || m == 0) {
649a8971b87SHong Zhang     c->brows = m;
650a8971b87SHong Zhang     c->bcols = 1;
651d4002b98SHong Zhang   } else if (isSELL) {
652f4b713efSHong Zhang     /* bcols is chosen s.t. dy-array takes 50% of local memory space as mat */
653d4002b98SHong Zhang     Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
654d4002b98SHong Zhang     Mat_SeqSELL *spA, *spB;
655f4b713efSHong Zhang     Mat          A, B;
656f4b713efSHong Zhang     PetscInt     nz, brows, bcols;
657f4b713efSHong Zhang     PetscReal    mem;
658f4b713efSHong Zhang 
659d4002b98SHong Zhang     bs = 1; /* only bs=1 is supported for MPISELL matrix */
660f4b713efSHong Zhang 
6619371c9d4SSatish Balay     A     = sell->A;
6629371c9d4SSatish Balay     spA   = (Mat_SeqSELL *)A->data;
6639371c9d4SSatish Balay     B     = sell->B;
6649371c9d4SSatish Balay     spB   = (Mat_SeqSELL *)B->data;
665f4b713efSHong Zhang     nz    = spA->nz + spB->nz; /* total local nonzero entries of mat */
666f4b713efSHong Zhang     mem   = nz * (sizeof(PetscScalar) + sizeof(PetscInt)) + 3 * m * sizeof(PetscInt);
667f4b713efSHong Zhang     bcols = (PetscInt)(0.5 * mem / (m * sizeof(PetscScalar)));
668f4b713efSHong Zhang     brows = 1000 / bcols;
669f4b713efSHong Zhang     if (bcols > nis) bcols = nis;
670f4b713efSHong Zhang     if (brows == 0 || brows > m) brows = m;
671f4b713efSHong Zhang     c->brows = brows;
672f4b713efSHong Zhang     c->bcols = bcols;
673531e53bdSHong Zhang   } else { /* mpiaij matrix */
674531e53bdSHong Zhang     /* bcols is chosen s.t. dy-array takes 50% of local memory space as mat */
675531e53bdSHong Zhang     Mat_MPIAIJ *aij = (Mat_MPIAIJ *)mat->data;
676531e53bdSHong Zhang     Mat_SeqAIJ *spA, *spB;
677531e53bdSHong Zhang     Mat         A, B;
678a8971b87SHong Zhang     PetscInt    nz, brows, bcols;
679531e53bdSHong Zhang     PetscReal   mem;
680531e53bdSHong Zhang 
681531e53bdSHong Zhang     bs = 1; /* only bs=1 is supported for MPIAIJ matrix */
682531e53bdSHong Zhang 
6839371c9d4SSatish Balay     A     = aij->A;
6849371c9d4SSatish Balay     spA   = (Mat_SeqAIJ *)A->data;
6859371c9d4SSatish Balay     B     = aij->B;
6869371c9d4SSatish Balay     spB   = (Mat_SeqAIJ *)B->data;
687531e53bdSHong Zhang     nz    = spA->nz + spB->nz; /* total local nonzero entries of mat */
688531e53bdSHong Zhang     mem   = nz * (sizeof(PetscScalar) + sizeof(PetscInt)) + 3 * m * sizeof(PetscInt);
689531e53bdSHong Zhang     bcols = (PetscInt)(0.5 * mem / (m * sizeof(PetscScalar)));
690531e53bdSHong Zhang     brows = 1000 / bcols;
691531e53bdSHong Zhang     if (bcols > nis) bcols = nis;
692531e53bdSHong Zhang     if (brows == 0 || brows > m) brows = m;
693531e53bdSHong Zhang     c->brows = brows;
694531e53bdSHong Zhang     c->bcols = bcols;
695531e53bdSHong Zhang   }
696a8971b87SHong Zhang 
697a8971b87SHong Zhang   c->M       = mat->rmap->N / bs; /* set the global rows and columns and local rows */
698a8971b87SHong Zhang   c->N       = mat->cmap->N / bs;
699a8971b87SHong Zhang   c->m       = mat->rmap->n / bs;
700a8971b87SHong Zhang   c->rstart  = mat->rmap->rstart / bs;
701a8971b87SHong Zhang   c->ncolors = nis;
70293dfae19SHong Zhang   PetscFunctionReturn(0);
70393dfae19SHong Zhang }
704bdaf1daeSBarry Smith 
705bdaf1daeSBarry Smith /*@C
706bdaf1daeSBarry Smith 
707bdaf1daeSBarry Smith     MatFDColoringSetValues - takes a matrix in compressed color format and enters the matrix into a PETSc Mat
708bdaf1daeSBarry Smith 
709bdaf1daeSBarry Smith    Collective on J
710bdaf1daeSBarry Smith 
711bdaf1daeSBarry Smith    Input Parameters:
712bdaf1daeSBarry Smith +    J - the sparse matrix
713bdaf1daeSBarry Smith .    coloring - created with MatFDColoringCreate() and a local coloring
714bdaf1daeSBarry Smith -    y - column major storage of matrix values with one color of values per column, the number of rows of y should match
715bdaf1daeSBarry Smith          the number of local rows of J and the number of columns is the number of colors.
716bdaf1daeSBarry Smith 
717bdaf1daeSBarry Smith    Level: intermediate
718bdaf1daeSBarry Smith 
7196a9b8d82SBarry Smith    Notes: the matrix in compressed color format may come from an Automatic Differentiation code
720bdaf1daeSBarry Smith 
721bdaf1daeSBarry Smith    The code will be slightly faster if MatFDColoringSetBlockSize(coloring,PETSC_DEFAULT,nc); is called immediately after creating the coloring
722bdaf1daeSBarry Smith 
723db781477SPatrick Sanan .seealso: `MatFDColoringCreate()`, `ISColoring`, `ISColoringCreate()`, `ISColoringSetType()`, `IS_COLORING_LOCAL`, `MatFDColoringSetBlockSize()`
724bdaf1daeSBarry Smith 
725bdaf1daeSBarry Smith @*/
7269371c9d4SSatish Balay PetscErrorCode MatFDColoringSetValues(Mat J, MatFDColoring coloring, const PetscScalar *y) {
727bdaf1daeSBarry Smith   MatEntry2      *Jentry2;
728bdaf1daeSBarry Smith   PetscInt        row, i, nrows_k, l, ncolors, nz = 0, bcols, nbcols = 0;
729bdaf1daeSBarry Smith   const PetscInt *nrows;
730bdaf1daeSBarry Smith   PetscBool       eq;
731bdaf1daeSBarry Smith 
732bdaf1daeSBarry Smith   PetscFunctionBegin;
733bdaf1daeSBarry Smith   PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
734bdaf1daeSBarry Smith   PetscValidHeaderSpecific(coloring, MAT_FDCOLORING_CLASSID, 2);
7359566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompareId((PetscObject)J, coloring->matid, &eq));
73628b400f6SJacob Faibussowitsch   PetscCheck(eq, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_WRONG, "Matrix used with MatFDColoringSetValues() must be that used with MatFDColoringCreate()");
737bdaf1daeSBarry Smith   Jentry2 = coloring->matentry2;
738bdaf1daeSBarry Smith   nrows   = coloring->nrows;
739bdaf1daeSBarry Smith   ncolors = coloring->ncolors;
740bdaf1daeSBarry Smith   bcols   = coloring->bcols;
741bdaf1daeSBarry Smith 
742bdaf1daeSBarry Smith   for (i = 0; i < ncolors; i += bcols) {
743bdaf1daeSBarry Smith     nrows_k = nrows[nbcols++];
744bdaf1daeSBarry Smith     for (l = 0; l < nrows_k; l++) {
745bdaf1daeSBarry Smith       row                      = Jentry2[nz].row; /* local row index */
746bdaf1daeSBarry Smith       *(Jentry2[nz++].valaddr) = y[row];
747bdaf1daeSBarry Smith     }
748bdaf1daeSBarry Smith     y += bcols * coloring->m;
749bdaf1daeSBarry Smith   }
7509566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
7519566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
752bdaf1daeSBarry Smith   PetscFunctionReturn(0);
753bdaf1daeSBarry Smith }
754