xref: /petsc/src/mat/impls/aij/mpi/fdmpiaij.c (revision c3339decea92175325d9368fa13196bcd0e0e58b)
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 
6d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringApply_BAIJ(Mat J, MatFDColoring coloring, Vec x1, void *sctx)
7d71ae5a4SJacob Faibussowitsch {
80d1c53f1SHong Zhang   PetscErrorCode (*f)(void *, Vec, Vec, void *) = (PetscErrorCode(*)(void *, Vec, Vec, void *))coloring->f;
90d1c53f1SHong Zhang   PetscInt           k, cstart, cend, l, row, col, nz, spidx, i, j;
105edff71fSBarry Smith   PetscScalar        dx = 0.0, *w3_array, *dy_i, *dy = coloring->dy;
110d1c53f1SHong Zhang   PetscScalar       *vscale_array;
125edff71fSBarry Smith   const PetscScalar *xx;
130d1c53f1SHong Zhang   PetscReal          epsilon = coloring->error_rel, umin = coloring->umin, unorm;
140d1c53f1SHong Zhang   Vec                w1 = coloring->w1, w2 = coloring->w2, w3, vscale = coloring->vscale;
150d1c53f1SHong Zhang   void              *fctx  = coloring->fctx;
160d1c53f1SHong Zhang   PetscInt           ctype = coloring->ctype, nxloc, nrows_k;
170d1c53f1SHong Zhang   PetscScalar       *valaddr;
180d1c53f1SHong Zhang   MatEntry          *Jentry  = coloring->matentry;
190df34763SHong Zhang   MatEntry2         *Jentry2 = coloring->matentry2;
200d1c53f1SHong Zhang   const PetscInt     ncolors = coloring->ncolors, *ncolumns = coloring->ncolumns, *nrows = coloring->nrows;
210d1c53f1SHong Zhang   PetscInt           bs = J->rmap->bs;
220d1c53f1SHong Zhang 
230d1c53f1SHong Zhang   PetscFunctionBegin;
249566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(x1, PETSC_TRUE));
250d1c53f1SHong Zhang   /* (1) Set w1 = F(x1) */
260d1c53f1SHong Zhang   if (!coloring->fset) {
279566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(MAT_FDColoringFunction, coloring, 0, 0, 0));
289566063dSJacob Faibussowitsch     PetscCall((*f)(sctx, x1, w1, fctx));
299566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(MAT_FDColoringFunction, coloring, 0, 0, 0));
300d1c53f1SHong Zhang   } else {
310d1c53f1SHong Zhang     coloring->fset = PETSC_FALSE;
320d1c53f1SHong Zhang   }
330d1c53f1SHong Zhang 
340d1c53f1SHong Zhang   /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */
359566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(x1, &nxloc));
360d1c53f1SHong Zhang   if (coloring->htype[0] == 'w') {
370d1c53f1SHong Zhang     /* vscale = dx is a constant scalar */
389566063dSJacob Faibussowitsch     PetscCall(VecNorm(x1, NORM_2, &unorm));
390d1c53f1SHong Zhang     dx = 1.0 / (PetscSqrtReal(1.0 + unorm) * epsilon);
400d1c53f1SHong Zhang   } else {
419566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(x1, &xx));
429566063dSJacob Faibussowitsch     PetscCall(VecGetArray(vscale, &vscale_array));
430d1c53f1SHong Zhang     for (col = 0; col < nxloc; col++) {
440d1c53f1SHong Zhang       dx = xx[col];
450d1c53f1SHong Zhang       if (PetscAbsScalar(dx) < umin) {
460d1c53f1SHong Zhang         if (PetscRealPart(dx) >= 0.0) dx = umin;
470d1c53f1SHong Zhang         else if (PetscRealPart(dx) < 0.0) dx = -umin;
480d1c53f1SHong Zhang       }
490d1c53f1SHong Zhang       dx *= epsilon;
500d1c53f1SHong Zhang       vscale_array[col] = 1.0 / dx;
510d1c53f1SHong Zhang     }
529566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(x1, &xx));
539566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(vscale, &vscale_array));
540d1c53f1SHong Zhang   }
55684f2004SHong Zhang   if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') {
569566063dSJacob Faibussowitsch     PetscCall(VecGhostUpdateBegin(vscale, INSERT_VALUES, SCATTER_FORWARD));
579566063dSJacob Faibussowitsch     PetscCall(VecGhostUpdateEnd(vscale, INSERT_VALUES, SCATTER_FORWARD));
580d1c53f1SHong Zhang   }
590d1c53f1SHong Zhang 
600d1c53f1SHong Zhang   /* (3) Loop over each color */
610d1c53f1SHong Zhang   if (!coloring->w3) {
629566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x1, &coloring->w3));
63ce911d59SRichard Tran Mills     /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */
649566063dSJacob Faibussowitsch     PetscCall(VecBindToCPU(coloring->w3, PETSC_TRUE));
650d1c53f1SHong Zhang   }
660d1c53f1SHong Zhang   w3 = coloring->w3;
670d1c53f1SHong Zhang 
689566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(x1, &cstart, &cend)); /* used by ghosted vscale */
6948a46eb9SPierre 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));
14948a46eb9SPierre 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 */
157d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringApply_AIJ(Mat J, MatFDColoring coloring, Vec x1, void *sctx)
158d71ae5a4SJacob Faibussowitsch {
159fcd7ac73SHong Zhang   PetscErrorCode (*f)(void *, Vec, Vec, void *) = (PetscErrorCode(*)(void *, Vec, Vec, void *))coloring->f;
160a2f2d239SHong Zhang   PetscInt           k, cstart, cend, l, row, col, nz;
1615edff71fSBarry Smith   PetscScalar        dx = 0.0, *y, *w3_array;
1625edff71fSBarry Smith   const PetscScalar *xx;
163fcd7ac73SHong Zhang   PetscScalar       *vscale_array;
164fcd7ac73SHong Zhang   PetscReal          epsilon = coloring->error_rel, umin = coloring->umin, unorm;
1658bc97078SHong Zhang   Vec                w1 = coloring->w1, w2 = coloring->w2, w3, vscale = coloring->vscale;
166fcd7ac73SHong Zhang   void              *fctx  = coloring->fctx;
16791e7fa0fSBarry Smith   ISColoringType     ctype = coloring->ctype;
16891e7fa0fSBarry Smith   PetscInt           nxloc, nrows_k;
169573f477fSHong Zhang   MatEntry          *Jentry  = coloring->matentry;
1700df34763SHong Zhang   MatEntry2         *Jentry2 = coloring->matentry2;
1718bc97078SHong Zhang   const PetscInt     ncolors = coloring->ncolors, *ncolumns = coloring->ncolumns, *nrows = coloring->nrows;
172b294de21SRichard Tran Mills   PetscBool          alreadyboundtocpu;
173fcd7ac73SHong Zhang 
174fcd7ac73SHong Zhang   PetscFunctionBegin;
1759566063dSJacob Faibussowitsch   PetscCall(VecBoundToCPU(x1, &alreadyboundtocpu));
1769566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(x1, PETSC_TRUE));
17708401ef6SPierre Jolivet   PetscCheck(!(ctype == IS_COLORING_LOCAL) || !(J->ops->fdcoloringapply == MatFDColoringApply_AIJ), PetscObjectComm((PetscObject)J), PETSC_ERR_SUP, "Must call MatColoringUseDM() with IS_COLORING_LOCAL");
1788bc97078SHong Zhang   /* (1) Set w1 = F(x1) */
179fcd7ac73SHong Zhang   if (!coloring->fset) {
1809566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(MAT_FDColoringFunction, 0, 0, 0, 0));
1819566063dSJacob Faibussowitsch     PetscCall((*f)(sctx, x1, w1, fctx));
1829566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(MAT_FDColoringFunction, 0, 0, 0, 0));
183fcd7ac73SHong Zhang   } else {
184fcd7ac73SHong Zhang     coloring->fset = PETSC_FALSE;
185fcd7ac73SHong Zhang   }
186fcd7ac73SHong Zhang 
1878bc97078SHong Zhang   /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */
188f6af9589SHong Zhang   if (coloring->htype[0] == 'w') {
189684f2004SHong Zhang     /* vscale = 1./dx is a constant scalar */
1909566063dSJacob Faibussowitsch     PetscCall(VecNorm(x1, NORM_2, &unorm));
191c53567a0SHong Zhang     dx = 1.0 / (PetscSqrtReal(1.0 + unorm) * epsilon);
19270e7395fSHong Zhang   } else {
1939566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(x1, &nxloc));
1949566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(x1, &xx));
1959566063dSJacob Faibussowitsch     PetscCall(VecGetArray(vscale, &vscale_array));
19674d3cef9SHong Zhang     for (col = 0; col < nxloc; col++) {
197fcd7ac73SHong Zhang       dx = xx[col];
19874d3cef9SHong Zhang       if (PetscAbsScalar(dx) < umin) {
19974d3cef9SHong Zhang         if (PetscRealPart(dx) >= 0.0) dx = umin;
20074d3cef9SHong Zhang         else if (PetscRealPart(dx) < 0.0) dx = -umin;
20174d3cef9SHong Zhang       }
202fcd7ac73SHong Zhang       dx *= epsilon;
20374d3cef9SHong Zhang       vscale_array[col] = 1.0 / dx;
204f6af9589SHong Zhang     }
2059566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(x1, &xx));
2069566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(vscale, &vscale_array));
20770e7395fSHong Zhang   }
208684f2004SHong Zhang   if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') {
2099566063dSJacob Faibussowitsch     PetscCall(VecGhostUpdateBegin(vscale, INSERT_VALUES, SCATTER_FORWARD));
2109566063dSJacob Faibussowitsch     PetscCall(VecGhostUpdateEnd(vscale, INSERT_VALUES, SCATTER_FORWARD));
211fcd7ac73SHong Zhang   }
212fcd7ac73SHong Zhang 
2138bc97078SHong Zhang   /* (3) Loop over each color */
2144dfa11a4SJacob Faibussowitsch   if (!coloring->w3) { PetscCall(VecDuplicate(x1, &coloring->w3)); }
2158bc97078SHong Zhang   w3 = coloring->w3;
216fcd7ac73SHong Zhang 
2179566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(x1, &cstart, &cend)); /* used by ghosted vscale */
21848a46eb9SPierre Jolivet   if (vscale) PetscCall(VecGetArray(vscale, &vscale_array));
2198bc97078SHong Zhang   nz = 0;
220e2c857f8SHong Zhang 
221b3e1f37bSHong Zhang   if (coloring->bcols > 1) { /* use blocked insertion of Jentry */
222b3e1f37bSHong Zhang     PetscInt     i, m = J->rmap->n, nbcols, bcols = coloring->bcols;
223e2c857f8SHong Zhang     PetscScalar *dy = coloring->dy, *dy_k;
224e2c857f8SHong Zhang 
225e2c857f8SHong Zhang     nbcols = 0;
226e2c857f8SHong Zhang     for (k = 0; k < ncolors; k += bcols) {
227e2c857f8SHong Zhang       /*
228e2c857f8SHong Zhang        (3-1) Loop over each column associated with color
229e2c857f8SHong Zhang        adding the perturbation to the vector w3 = x1 + dx.
230e2c857f8SHong Zhang        */
231e2c857f8SHong Zhang 
232e2c857f8SHong Zhang       dy_k = dy;
233e2c857f8SHong Zhang       if (k + bcols > ncolors) bcols = ncolors - k;
234e2c857f8SHong Zhang       for (i = 0; i < bcols; i++) {
235ceef8a49SBarry Smith         coloring->currentcolor = k + i;
236e2c857f8SHong Zhang 
2379566063dSJacob Faibussowitsch         PetscCall(VecCopy(x1, w3));
2389566063dSJacob Faibussowitsch         PetscCall(VecGetArray(w3, &w3_array));
239e2c857f8SHong Zhang         if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */
240e2c857f8SHong Zhang         if (coloring->htype[0] == 'w') {
241e2c857f8SHong Zhang           for (l = 0; l < ncolumns[k + i]; l++) {
242e2c857f8SHong Zhang             col = coloring->columns[k + i][l]; /* local column (in global index!) of the matrix we are probing for */
243e2c857f8SHong Zhang             w3_array[col] += 1.0 / dx;
244e2c857f8SHong Zhang           }
245e2c857f8SHong Zhang         } else {                  /* htype == 'ds' */
246e2c857f8SHong Zhang           vscale_array -= cstart; /* shift pointer so global index can be used */
247e2c857f8SHong Zhang           for (l = 0; l < ncolumns[k + i]; l++) {
248e2c857f8SHong Zhang             col = coloring->columns[k + i][l]; /* local column (in global index!) of the matrix we are probing for */
249e2c857f8SHong Zhang             w3_array[col] += 1.0 / vscale_array[col];
250e2c857f8SHong Zhang           }
251e2c857f8SHong Zhang           vscale_array += cstart;
252e2c857f8SHong Zhang         }
253e2c857f8SHong Zhang         if (ctype == IS_COLORING_GLOBAL) w3_array += cstart;
2549566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(w3, &w3_array));
255e2c857f8SHong Zhang 
256e2c857f8SHong Zhang         /*
257e2c857f8SHong Zhang          (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
258e2c857f8SHong Zhang                            w2 = F(x1 + dx) - F(x1)
259e2c857f8SHong Zhang          */
2609566063dSJacob Faibussowitsch         PetscCall(PetscLogEventBegin(MAT_FDColoringFunction, 0, 0, 0, 0));
2619566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(w2, dy_k)); /* place w2 to the array dy_i */
2629566063dSJacob Faibussowitsch         PetscCall((*f)(sctx, w3, w2, fctx));
2639566063dSJacob Faibussowitsch         PetscCall(PetscLogEventEnd(MAT_FDColoringFunction, 0, 0, 0, 0));
2649566063dSJacob Faibussowitsch         PetscCall(VecAXPY(w2, -1.0, w1));
2659566063dSJacob Faibussowitsch         PetscCall(VecResetArray(w2));
266e2c857f8SHong Zhang         dy_k += m; /* points to dy+i*nxloc */
267e2c857f8SHong Zhang       }
268e2c857f8SHong Zhang 
269e2c857f8SHong Zhang       /*
270e2c857f8SHong Zhang        (3-3) Loop over block rows of vector, putting results into Jacobian matrix
271e2c857f8SHong Zhang        */
272d880da65SHong Zhang       nrows_k = nrows[nbcols++];
273d880da65SHong Zhang 
274e2c857f8SHong Zhang       if (coloring->htype[0] == 'w') {
275e2c857f8SHong Zhang         for (l = 0; l < nrows_k; l++) {
2760df34763SHong Zhang           row = Jentry2[nz].row; /* local row index */
2775059b303SJunchao 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
2785059b303SJunchao Zhang              another way, and it seems work. See https://lists.mcs.anl.gov/pipermail/petsc-users/2021-December/045158.html
2795059b303SJunchao Zhang            */
2805059b303SJunchao Zhang #if defined(PETSC_USE_COMPLEX)
2815059b303SJunchao Zhang           PetscScalar *tmp = Jentry2[nz].valaddr;
2825059b303SJunchao Zhang           *tmp             = dy[row] * dx;
2835059b303SJunchao Zhang #else
2845059b303SJunchao Zhang           *(Jentry2[nz].valaddr) = dy[row] * dx;
2855059b303SJunchao Zhang #endif
2865059b303SJunchao Zhang           nz++;
287e2c857f8SHong Zhang         }
288e2c857f8SHong Zhang       } else { /* htype == 'ds' */
289e2c857f8SHong Zhang         for (l = 0; l < nrows_k; l++) {
290e2c857f8SHong Zhang           row = Jentry[nz].row; /* local row index */
2915059b303SJunchao Zhang #if defined(PETSC_USE_COMPLEX)  /* See https://lists.mcs.anl.gov/pipermail/petsc-users/2021-December/045158.html */
2925059b303SJunchao Zhang           PetscScalar *tmp = Jentry[nz].valaddr;
2935059b303SJunchao Zhang           *tmp             = dy[row] * vscale_array[Jentry[nz].col];
2945059b303SJunchao Zhang #else
295d880da65SHong Zhang           *(Jentry[nz].valaddr)  = dy[row] * vscale_array[Jentry[nz].col];
2965059b303SJunchao Zhang #endif
297e2c857f8SHong Zhang           nz++;
298e2c857f8SHong Zhang         }
299e2c857f8SHong Zhang       }
300e2c857f8SHong Zhang     }
301b3e1f37bSHong Zhang   } else { /* bcols == 1 */
3028bc97078SHong Zhang     for (k = 0; k < ncolors; k++) {
3038bc97078SHong Zhang       coloring->currentcolor = k;
3049e917edbSHong Zhang 
305fcd7ac73SHong Zhang       /*
3068bc97078SHong Zhang        (3-1) Loop over each column associated with color
307f6af9589SHong Zhang        adding the perturbation to the vector w3 = x1 + dx.
308fcd7ac73SHong Zhang        */
3099566063dSJacob Faibussowitsch       PetscCall(VecCopy(x1, w3));
3109566063dSJacob Faibussowitsch       PetscCall(VecGetArray(w3, &w3_array));
311a2f2d239SHong Zhang       if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */
312b039d6c4SHong Zhang       if (coloring->htype[0] == 'w') {
3138bc97078SHong Zhang         for (l = 0; l < ncolumns[k]; l++) {
314f6af9589SHong Zhang           col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
3159e917edbSHong Zhang           w3_array[col] += 1.0 / dx;
3169e917edbSHong Zhang         }
317b039d6c4SHong Zhang       } else {                  /* htype == 'ds' */
318a2f2d239SHong Zhang         vscale_array -= cstart; /* shift pointer so global index can be used */
319b039d6c4SHong Zhang         for (l = 0; l < ncolumns[k]; l++) {
320b039d6c4SHong Zhang           col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
321a2f2d239SHong Zhang           w3_array[col] += 1.0 / vscale_array[col];
32270e7395fSHong Zhang         }
323a2f2d239SHong Zhang         vscale_array += cstart;
324fcd7ac73SHong Zhang       }
325a2f2d239SHong Zhang       if (ctype == IS_COLORING_GLOBAL) w3_array += cstart;
3269566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(w3, &w3_array));
3279e917edbSHong Zhang 
328fcd7ac73SHong Zhang       /*
3298bc97078SHong Zhang        (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
330fcd7ac73SHong Zhang                            w2 = F(x1 + dx) - F(x1)
331fcd7ac73SHong Zhang        */
3329566063dSJacob Faibussowitsch       PetscCall(PetscLogEventBegin(MAT_FDColoringFunction, 0, 0, 0, 0));
3339566063dSJacob Faibussowitsch       PetscCall((*f)(sctx, w3, w2, fctx));
3349566063dSJacob Faibussowitsch       PetscCall(PetscLogEventEnd(MAT_FDColoringFunction, 0, 0, 0, 0));
3359566063dSJacob Faibussowitsch       PetscCall(VecAXPY(w2, -1.0, w1));
3369e917edbSHong Zhang 
3378bc97078SHong Zhang       /*
3388bc97078SHong Zhang        (3-3) Loop over rows of vector, putting results into Jacobian matrix
3398bc97078SHong Zhang        */
340c53567a0SHong Zhang       nrows_k = nrows[k];
3419566063dSJacob Faibussowitsch       PetscCall(VecGetArray(w2, &y));
342b039d6c4SHong Zhang       if (coloring->htype[0] == 'w') {
343c53567a0SHong Zhang         for (l = 0; l < nrows_k; l++) {
3440df34763SHong Zhang           row = Jentry2[nz].row; /* local row index */
3455059b303SJunchao Zhang #if defined(PETSC_USE_COMPLEX)   /* See https://lists.mcs.anl.gov/pipermail/petsc-users/2021-December/045158.html */
3465059b303SJunchao Zhang           PetscScalar *tmp = Jentry2[nz].valaddr;
3475059b303SJunchao Zhang           *tmp             = y[row] * dx;
3485059b303SJunchao Zhang #else
3495059b303SJunchao Zhang           *(Jentry2[nz].valaddr) = y[row] * dx;
3505059b303SJunchao Zhang #endif
3515059b303SJunchao Zhang           nz++;
3529e917edbSHong Zhang         }
353b039d6c4SHong Zhang       } else { /* htype == 'ds' */
354c53567a0SHong Zhang         for (l = 0; l < nrows_k; l++) {
355b039d6c4SHong Zhang           row = Jentry[nz].row; /* local row index */
3565059b303SJunchao Zhang #if defined(PETSC_USE_COMPLEX)  /* See https://lists.mcs.anl.gov/pipermail/petsc-users/2021-December/045158.html */
3575059b303SJunchao Zhang           PetscScalar *tmp = Jentry[nz].valaddr;
3585059b303SJunchao Zhang           *tmp             = y[row] * vscale_array[Jentry[nz].col];
3595059b303SJunchao Zhang #else
360b039d6c4SHong Zhang           *(Jentry[nz].valaddr)  = y[row] * vscale_array[Jentry[nz].col];
3615059b303SJunchao Zhang #endif
362573f477fSHong Zhang           nz++;
363fcd7ac73SHong Zhang         }
364b039d6c4SHong Zhang       }
3659566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(w2, &y));
3668bc97078SHong Zhang     }
367b3e1f37bSHong Zhang   }
368b3e1f37bSHong Zhang 
369e2cf4d64SStefano Zampini #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
370c70f7ee4SJunchao Zhang   if (J->offloadmask != PETSC_OFFLOAD_UNALLOCATED) J->offloadmask = PETSC_OFFLOAD_CPU;
371e2cf4d64SStefano Zampini #endif
3729566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
3739566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
37448a46eb9SPierre Jolivet   if (vscale) PetscCall(VecRestoreArray(vscale, &vscale_array));
375fcd7ac73SHong Zhang   coloring->currentcolor = -1;
3769566063dSJacob Faibussowitsch   if (!alreadyboundtocpu) PetscCall(VecBindToCPU(x1, PETSC_FALSE));
377fcd7ac73SHong Zhang   PetscFunctionReturn(0);
378fcd7ac73SHong Zhang }
379fcd7ac73SHong Zhang 
380d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetUp_MPIXAIJ(Mat mat, ISColoring iscoloring, MatFDColoring c)
381d71ae5a4SJacob Faibussowitsch {
382fcd7ac73SHong Zhang   PetscMPIInt            size, *ncolsonproc, *disp, nn;
383e3225e9dSHong Zhang   PetscInt               i, n, nrows, nrows_i, j, k, m, ncols, col, *rowhit, cstart, cend, colb;
38472c15787SHong Zhang   const PetscInt        *is, *A_ci, *A_cj, *B_ci, *B_cj, *row = NULL, *ltog = NULL;
385071fcb05SBarry Smith   PetscInt               nis = iscoloring->n, nctot, *cols, tmp = 0;
386fcd7ac73SHong Zhang   ISLocalToGlobalMapping map   = mat->cmap->mapping;
387a8971b87SHong Zhang   PetscInt               ctype = c->ctype, *spidxA, *spidxB, nz, bs, bs2, spidx;
3880d1c53f1SHong Zhang   Mat                    A, B;
389e3225e9dSHong Zhang   PetscScalar           *A_val, *B_val, **valaddrhit;
390573f477fSHong Zhang   MatEntry              *Jentry;
3910df34763SHong Zhang   MatEntry2             *Jentry2;
392d4002b98SHong Zhang   PetscBool              isBAIJ, isSELL;
393531e53bdSHong Zhang   PetscInt               bcols = c->bcols;
3940d1c53f1SHong Zhang #if defined(PETSC_USE_CTABLE)
3950d1c53f1SHong Zhang   PetscTable colmap = NULL;
3960d1c53f1SHong Zhang #else
3970d1c53f1SHong Zhang   PetscInt *colmap = NULL;      /* local col number of off-diag col */
3980d1c53f1SHong Zhang #endif
399fcd7ac73SHong Zhang 
400fcd7ac73SHong Zhang   PetscFunctionBegin;
4015bdb020cSBarry Smith   if (ctype == IS_COLORING_LOCAL) {
40228b400f6SJacob Faibussowitsch     PetscCheck(map, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_INCOMP, "When using ghosted differencing matrix must have local to global mapping provided with MatSetLocalToGlobalMapping");
4039566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetIndices(map, &ltog));
40472c15787SHong Zhang   }
405fcd7ac73SHong Zhang 
4069566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(mat, &bs));
4079566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompare((PetscObject)mat, MATMPIBAIJ, &isBAIJ));
4089566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATMPISELL, &isSELL));
409e3225e9dSHong Zhang   if (isBAIJ) {
410195687c5SHong Zhang     Mat_MPIBAIJ *baij = (Mat_MPIBAIJ *)mat->data;
4116c145f34SHong Zhang     Mat_SeqBAIJ *spA, *spB;
4129371c9d4SSatish Balay     A     = baij->A;
4139371c9d4SSatish Balay     spA   = (Mat_SeqBAIJ *)A->data;
4149371c9d4SSatish Balay     A_val = spA->a;
4159371c9d4SSatish Balay     B     = baij->B;
4169371c9d4SSatish Balay     spB   = (Mat_SeqBAIJ *)B->data;
4179371c9d4SSatish Balay     B_val = spB->a;
4180d1c53f1SHong Zhang     nz    = spA->nz + spB->nz; /* total nonzero entries of mat */
41948a46eb9SPierre Jolivet     if (!baij->colmap) PetscCall(MatCreateColmap_MPIBAIJ_Private(mat));
420097e26fcSJed Brown     colmap = baij->colmap;
4219566063dSJacob Faibussowitsch     PetscCall(MatGetColumnIJ_SeqBAIJ_Color(A, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &A_ci, &A_cj, &spidxA, NULL));
4229566063dSJacob Faibussowitsch     PetscCall(MatGetColumnIJ_SeqBAIJ_Color(B, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &B_ci, &B_cj, &spidxB, NULL));
423195687c5SHong Zhang 
424195687c5SHong Zhang     if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') { /* create vscale for storing dx */
425195687c5SHong Zhang       PetscInt *garray;
4269566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(B->cmap->n, &garray));
427195687c5SHong Zhang       for (i = 0; i < baij->B->cmap->n / bs; i++) {
428ad540459SPierre Jolivet         for (j = 0; j < bs; j++) garray[i * bs + j] = bs * baij->garray[i] + j;
429195687c5SHong Zhang       }
4309566063dSJacob Faibussowitsch       PetscCall(VecCreateGhost(PetscObjectComm((PetscObject)mat), mat->cmap->n, PETSC_DETERMINE, B->cmap->n, garray, &c->vscale));
4319566063dSJacob Faibussowitsch       PetscCall(VecBindToCPU(c->vscale, PETSC_TRUE));
4329566063dSJacob Faibussowitsch       PetscCall(PetscFree(garray));
433195687c5SHong Zhang     }
434d4002b98SHong Zhang   } else if (isSELL) {
435d4002b98SHong Zhang     Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
436d4002b98SHong Zhang     Mat_SeqSELL *spA, *spB;
4379371c9d4SSatish Balay     A     = sell->A;
4389371c9d4SSatish Balay     spA   = (Mat_SeqSELL *)A->data;
4399371c9d4SSatish Balay     A_val = spA->val;
4409371c9d4SSatish Balay     B     = sell->B;
4419371c9d4SSatish Balay     spB   = (Mat_SeqSELL *)B->data;
4429371c9d4SSatish Balay     B_val = spB->val;
443f4b713efSHong Zhang     nz    = spA->nz + spB->nz; /* total nonzero entries of mat */
444d4002b98SHong Zhang     if (!sell->colmap) {
445f4b713efSHong Zhang       /* Allow access to data structures of local part of matrix
446f4b713efSHong Zhang        - creates aij->colmap which maps global column number to local number in part B */
4479566063dSJacob Faibussowitsch       PetscCall(MatCreateColmap_MPISELL_Private(mat));
448f4b713efSHong Zhang     }
449d4002b98SHong Zhang     colmap = sell->colmap;
4509566063dSJacob Faibussowitsch     PetscCall(MatGetColumnIJ_SeqSELL_Color(A, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &A_ci, &A_cj, &spidxA, NULL));
4519566063dSJacob Faibussowitsch     PetscCall(MatGetColumnIJ_SeqSELL_Color(B, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &B_ci, &B_cj, &spidxB, NULL));
452f4b713efSHong Zhang 
453f4b713efSHong Zhang     bs = 1; /* only bs=1 is supported for non MPIBAIJ matrix */
454f4b713efSHong Zhang 
455f4b713efSHong Zhang     if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') { /* create vscale for storing dx */
4569566063dSJacob Faibussowitsch       PetscCall(VecCreateGhost(PetscObjectComm((PetscObject)mat), mat->cmap->n, PETSC_DETERMINE, B->cmap->n, sell->garray, &c->vscale));
4579566063dSJacob Faibussowitsch       PetscCall(VecBindToCPU(c->vscale, PETSC_TRUE));
458f4b713efSHong Zhang     }
459e3225e9dSHong Zhang   } else {
460e3225e9dSHong Zhang     Mat_MPIAIJ *aij = (Mat_MPIAIJ *)mat->data;
461e3225e9dSHong Zhang     Mat_SeqAIJ *spA, *spB;
4629371c9d4SSatish Balay     A     = aij->A;
4639371c9d4SSatish Balay     spA   = (Mat_SeqAIJ *)A->data;
4649371c9d4SSatish Balay     A_val = spA->a;
4659371c9d4SSatish Balay     B     = aij->B;
4669371c9d4SSatish Balay     spB   = (Mat_SeqAIJ *)B->data;
4679371c9d4SSatish Balay     B_val = spB->a;
468e3225e9dSHong Zhang     nz    = spA->nz + spB->nz; /* total nonzero entries of mat */
469e3225e9dSHong Zhang     if (!aij->colmap) {
470e3225e9dSHong Zhang       /* Allow access to data structures of local part of matrix
471e3225e9dSHong Zhang        - creates aij->colmap which maps global column number to local number in part B */
4729566063dSJacob Faibussowitsch       PetscCall(MatCreateColmap_MPIAIJ_Private(mat));
473e3225e9dSHong Zhang     }
474097e26fcSJed Brown     colmap = aij->colmap;
4759566063dSJacob Faibussowitsch     PetscCall(MatGetColumnIJ_SeqAIJ_Color(A, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &A_ci, &A_cj, &spidxA, NULL));
4769566063dSJacob Faibussowitsch     PetscCall(MatGetColumnIJ_SeqAIJ_Color(B, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &B_ci, &B_cj, &spidxB, NULL));
477e3225e9dSHong Zhang 
478e3225e9dSHong Zhang     bs = 1; /* only bs=1 is supported for non MPIBAIJ matrix */
479195687c5SHong Zhang 
480195687c5SHong Zhang     if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') { /* create vscale for storing dx */
4819566063dSJacob Faibussowitsch       PetscCall(VecCreateGhost(PetscObjectComm((PetscObject)mat), mat->cmap->n, PETSC_DETERMINE, B->cmap->n, aij->garray, &c->vscale));
4829566063dSJacob Faibussowitsch       PetscCall(VecBindToCPU(c->vscale, PETSC_TRUE));
483195687c5SHong Zhang     }
4840d1c53f1SHong Zhang   }
4850d1c53f1SHong Zhang 
4860d1c53f1SHong Zhang   m      = mat->rmap->n / bs;
4870d1c53f1SHong Zhang   cstart = mat->cmap->rstart / bs;
4880d1c53f1SHong Zhang   cend   = mat->cmap->rend / bs;
489fcd7ac73SHong Zhang 
4909566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nis, &c->ncolumns, nis, &c->columns));
4919566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nis, &c->nrows));
492fcd7ac73SHong Zhang 
4930df34763SHong Zhang   if (c->htype[0] == 'd') {
4949566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz, &Jentry));
4958bc97078SHong Zhang     c->matentry = Jentry;
4960df34763SHong Zhang   } else if (c->htype[0] == 'w') {
4979566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz, &Jentry2));
4980df34763SHong Zhang     c->matentry2 = Jentry2;
4990df34763SHong Zhang   } else SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "htype is not supported");
500a774a6f1SHong Zhang 
5019566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(m + 1, &rowhit, m + 1, &valaddrhit));
502fcd7ac73SHong Zhang   nz = 0;
5039566063dSJacob Faibussowitsch   PetscCall(ISColoringGetIS(iscoloring, PETSC_OWN_POINTER, PETSC_IGNORE, &c->isa));
504071fcb05SBarry Smith 
505071fcb05SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
5069566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &size));
5079566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(size, &ncolsonproc, size, &disp));
508071fcb05SBarry Smith   }
509071fcb05SBarry Smith 
510d3825b63SHong Zhang   for (i = 0; i < nis; i++) { /* for each local color */
5119566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(c->isa[i], &n));
5129566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(c->isa[i], &is));
513fcd7ac73SHong Zhang 
514fcd7ac73SHong Zhang     c->ncolumns[i] = n; /* local number of columns of this color on this process */
515071fcb05SBarry Smith     c->columns[i]  = (PetscInt *)is;
516fcd7ac73SHong Zhang 
517fcd7ac73SHong Zhang     if (ctype == IS_COLORING_GLOBAL) {
518d3825b63SHong Zhang       /* Determine nctot, the total (parallel) number of columns of this color */
519d3825b63SHong Zhang       /* ncolsonproc[j]: local ncolumns on proc[j] of this color */
5209566063dSJacob Faibussowitsch       PetscCall(PetscMPIIntCast(n, &nn));
5219566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Allgather(&nn, 1, MPI_INT, ncolsonproc, 1, MPI_INT, PetscObjectComm((PetscObject)mat)));
5229371c9d4SSatish Balay       nctot = 0;
5239371c9d4SSatish Balay       for (j = 0; j < size; j++) nctot += ncolsonproc[j];
52448a46eb9SPierre Jolivet       if (!nctot) PetscCall(PetscInfo(mat, "Coloring of matrix has some unneeded colors with no corresponding rows\n"));
525fcd7ac73SHong Zhang 
526fcd7ac73SHong Zhang       disp[0] = 0;
527ad540459SPierre Jolivet       for (j = 1; j < size; j++) disp[j] = disp[j - 1] + ncolsonproc[j - 1];
528fcd7ac73SHong Zhang 
529d3825b63SHong Zhang       /* Get cols, the complete list of columns for this color on each process */
5309566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nctot + 1, &cols));
5319566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Allgatherv((void *)is, n, MPIU_INT, cols, ncolsonproc, disp, MPIU_INT, PetscObjectComm((PetscObject)mat)));
5325bdb020cSBarry Smith     } else if (ctype == IS_COLORING_LOCAL) {
533fcd7ac73SHong Zhang       /* Determine local number of columns of this color on this process, including ghost points */
534fcd7ac73SHong Zhang       nctot = n;
535071fcb05SBarry Smith       cols  = (PetscInt *)is;
536fcd7ac73SHong Zhang     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not provided for this MatFDColoring type");
537fcd7ac73SHong Zhang 
5381b97d346SHong Zhang     /* Mark all rows affect by these columns */
5399566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(rowhit, m));
5400d1c53f1SHong Zhang     bs2     = bs * bs;
5414b2e90caSHong Zhang     nrows_i = 0;
5421b97d346SHong Zhang     for (j = 0; j < nctot; j++) { /* loop over columns*/
5435bdb020cSBarry Smith       if (ctype == IS_COLORING_LOCAL) {
544fcd7ac73SHong Zhang         col = ltog[cols[j]];
545fcd7ac73SHong Zhang       } else {
546fcd7ac73SHong Zhang         col = cols[j];
547fcd7ac73SHong Zhang       }
548e3225e9dSHong Zhang       if (col >= cstart && col < cend) { /* column is in A, diagonal block of mat */
549071fcb05SBarry Smith         tmp   = A_ci[col - cstart];
550071fcb05SBarry Smith         row   = A_cj + tmp;
551071fcb05SBarry Smith         nrows = A_ci[col - cstart + 1] - tmp;
5524b2e90caSHong Zhang         nrows_i += nrows;
553071fcb05SBarry Smith 
554fcd7ac73SHong Zhang         /* loop over columns of A marking them in rowhit */
555d3825b63SHong Zhang         for (k = 0; k < nrows; k++) {
556d3825b63SHong Zhang           /* set valaddrhit for part A */
557071fcb05SBarry Smith           spidx            = bs2 * spidxA[tmp + k];
558e3225e9dSHong Zhang           valaddrhit[*row] = &A_val[spidx];
559a774a6f1SHong Zhang           rowhit[*row++]   = col - cstart + 1; /* local column index */
560fcd7ac73SHong Zhang         }
561e3225e9dSHong Zhang       } else { /* column is in B, off-diagonal block of mat */
562fcd7ac73SHong Zhang #if defined(PETSC_USE_CTABLE)
5639566063dSJacob Faibussowitsch         PetscCall(PetscTableFind(colmap, col + 1, &colb));
564fcd7ac73SHong Zhang         colb--;
565fcd7ac73SHong Zhang #else
5660d1c53f1SHong Zhang         colb = colmap[col] - 1; /* local column index */
567fcd7ac73SHong Zhang #endif
568fcd7ac73SHong Zhang         if (colb == -1) {
569d3825b63SHong Zhang           nrows = 0;
570fcd7ac73SHong Zhang         } else {
5710d1c53f1SHong Zhang           colb  = colb / bs;
572071fcb05SBarry Smith           tmp   = B_ci[colb];
573071fcb05SBarry Smith           row   = B_cj + tmp;
574071fcb05SBarry Smith           nrows = B_ci[colb + 1] - tmp;
575fcd7ac73SHong Zhang         }
5764b2e90caSHong Zhang         nrows_i += nrows;
577fcd7ac73SHong Zhang         /* loop over columns of B marking them in rowhit */
578d3825b63SHong Zhang         for (k = 0; k < nrows; k++) {
579d3825b63SHong Zhang           /* set valaddrhit for part B */
580071fcb05SBarry Smith           spidx            = bs2 * spidxB[tmp + k];
581e3225e9dSHong Zhang           valaddrhit[*row] = &B_val[spidx];
58270e7395fSHong Zhang           rowhit[*row++]   = colb + 1 + cend - cstart; /* local column index */
583fcd7ac73SHong Zhang         }
584fcd7ac73SHong Zhang       }
585fcd7ac73SHong Zhang     }
5864b2e90caSHong Zhang     c->nrows[i] = nrows_i;
5878bc97078SHong Zhang 
5880df34763SHong Zhang     if (c->htype[0] == 'd') {
589d3825b63SHong Zhang       for (j = 0; j < m; j++) {
590fcd7ac73SHong Zhang         if (rowhit[j]) {
591573f477fSHong Zhang           Jentry[nz].row     = j;             /* local row index */
592573f477fSHong Zhang           Jentry[nz].col     = rowhit[j] - 1; /* local column index */
593573f477fSHong Zhang           Jentry[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */
594573f477fSHong Zhang           nz++;
595fcd7ac73SHong Zhang         }
596fcd7ac73SHong Zhang       }
5970df34763SHong Zhang     } else { /* c->htype == 'wp' */
5980df34763SHong Zhang       for (j = 0; j < m; j++) {
5990df34763SHong Zhang         if (rowhit[j]) {
6000df34763SHong Zhang           Jentry2[nz].row     = j;             /* local row index */
6010df34763SHong Zhang           Jentry2[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */
6020df34763SHong Zhang           nz++;
6030df34763SHong Zhang         }
6040df34763SHong Zhang       }
6050df34763SHong Zhang     }
60648a46eb9SPierre Jolivet     if (ctype == IS_COLORING_GLOBAL) PetscCall(PetscFree(cols));
607fcd7ac73SHong Zhang   }
60848a46eb9SPierre Jolivet   if (ctype == IS_COLORING_GLOBAL) PetscCall(PetscFree2(ncolsonproc, disp));
609fcd7ac73SHong Zhang 
610a8971b87SHong Zhang   if (bcols > 1) { /* reorder Jentry for faster MatFDColoringApply() */
6119566063dSJacob Faibussowitsch     PetscCall(MatFDColoringSetUpBlocked_AIJ_Private(mat, c, nz));
612531e53bdSHong Zhang   }
613531e53bdSHong Zhang 
6140d1c53f1SHong Zhang   if (isBAIJ) {
6159566063dSJacob Faibussowitsch     PetscCall(MatRestoreColumnIJ_SeqBAIJ_Color(A, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &A_ci, &A_cj, &spidxA, NULL));
6169566063dSJacob Faibussowitsch     PetscCall(MatRestoreColumnIJ_SeqBAIJ_Color(B, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &B_ci, &B_cj, &spidxB, NULL));
6179566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(bs * mat->rmap->n, &c->dy));
618d4002b98SHong Zhang   } else if (isSELL) {
6199566063dSJacob Faibussowitsch     PetscCall(MatRestoreColumnIJ_SeqSELL_Color(A, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &A_ci, &A_cj, &spidxA, NULL));
6209566063dSJacob Faibussowitsch     PetscCall(MatRestoreColumnIJ_SeqSELL_Color(B, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &B_ci, &B_cj, &spidxB, NULL));
6210d1c53f1SHong Zhang   } else {
6229566063dSJacob Faibussowitsch     PetscCall(MatRestoreColumnIJ_SeqAIJ_Color(A, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &A_ci, &A_cj, &spidxA, NULL));
6239566063dSJacob Faibussowitsch     PetscCall(MatRestoreColumnIJ_SeqAIJ_Color(B, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &B_ci, &B_cj, &spidxB, NULL));
6240d1c53f1SHong Zhang   }
6250d1c53f1SHong Zhang 
6269566063dSJacob Faibussowitsch   PetscCall(ISColoringRestoreIS(iscoloring, PETSC_OWN_POINTER, &c->isa));
6279566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rowhit, valaddrhit));
628a8971b87SHong Zhang 
62948a46eb9SPierre Jolivet   if (ctype == IS_COLORING_LOCAL) PetscCall(ISLocalToGlobalMappingRestoreIndices(map, &ltog));
6309566063dSJacob Faibussowitsch   PetscCall(PetscInfo(c, "ncolors %" PetscInt_FMT ", brows %" PetscInt_FMT " and bcols %" PetscInt_FMT " are used.\n", c->ncolors, c->brows, c->bcols));
631fcd7ac73SHong Zhang   PetscFunctionReturn(0);
632fcd7ac73SHong Zhang }
63393dfae19SHong Zhang 
634d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat mat, ISColoring iscoloring, MatFDColoring c)
635d71ae5a4SJacob Faibussowitsch {
636a8971b87SHong Zhang   PetscInt  bs, nis = iscoloring->n, m = mat->rmap->n;
637d4002b98SHong Zhang   PetscBool isBAIJ, isSELL;
638531e53bdSHong Zhang 
63993dfae19SHong Zhang   PetscFunctionBegin;
640531e53bdSHong Zhang   /* set default brows and bcols for speedup inserting the dense matrix into sparse Jacobian;
641531e53bdSHong Zhang    bcols is chosen s.t. dy-array takes 50% of memory space as mat */
6429566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(mat, &bs));
6439566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompare((PetscObject)mat, MATMPIBAIJ, &isBAIJ));
6449566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATMPISELL, &isSELL));
645b8dfa574SBarry Smith   if (isBAIJ || m == 0) {
646a8971b87SHong Zhang     c->brows = m;
647a8971b87SHong Zhang     c->bcols = 1;
648d4002b98SHong Zhang   } else if (isSELL) {
649f4b713efSHong Zhang     /* bcols is chosen s.t. dy-array takes 50% of local memory space as mat */
650d4002b98SHong Zhang     Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
651d4002b98SHong Zhang     Mat_SeqSELL *spA, *spB;
652f4b713efSHong Zhang     Mat          A, B;
653f4b713efSHong Zhang     PetscInt     nz, brows, bcols;
654f4b713efSHong Zhang     PetscReal    mem;
655f4b713efSHong Zhang 
656d4002b98SHong Zhang     bs = 1; /* only bs=1 is supported for MPISELL matrix */
657f4b713efSHong Zhang 
6589371c9d4SSatish Balay     A     = sell->A;
6599371c9d4SSatish Balay     spA   = (Mat_SeqSELL *)A->data;
6609371c9d4SSatish Balay     B     = sell->B;
6619371c9d4SSatish Balay     spB   = (Mat_SeqSELL *)B->data;
662f4b713efSHong Zhang     nz    = spA->nz + spB->nz; /* total local nonzero entries of mat */
663f4b713efSHong Zhang     mem   = nz * (sizeof(PetscScalar) + sizeof(PetscInt)) + 3 * m * sizeof(PetscInt);
664f4b713efSHong Zhang     bcols = (PetscInt)(0.5 * mem / (m * sizeof(PetscScalar)));
665f4b713efSHong Zhang     brows = 1000 / bcols;
666f4b713efSHong Zhang     if (bcols > nis) bcols = nis;
667f4b713efSHong Zhang     if (brows == 0 || brows > m) brows = m;
668f4b713efSHong Zhang     c->brows = brows;
669f4b713efSHong Zhang     c->bcols = bcols;
670531e53bdSHong Zhang   } else { /* mpiaij matrix */
671531e53bdSHong Zhang     /* bcols is chosen s.t. dy-array takes 50% of local memory space as mat */
672531e53bdSHong Zhang     Mat_MPIAIJ *aij = (Mat_MPIAIJ *)mat->data;
673531e53bdSHong Zhang     Mat_SeqAIJ *spA, *spB;
674531e53bdSHong Zhang     Mat         A, B;
675a8971b87SHong Zhang     PetscInt    nz, brows, bcols;
676531e53bdSHong Zhang     PetscReal   mem;
677531e53bdSHong Zhang 
678531e53bdSHong Zhang     bs = 1; /* only bs=1 is supported for MPIAIJ matrix */
679531e53bdSHong Zhang 
6809371c9d4SSatish Balay     A     = aij->A;
6819371c9d4SSatish Balay     spA   = (Mat_SeqAIJ *)A->data;
6829371c9d4SSatish Balay     B     = aij->B;
6839371c9d4SSatish Balay     spB   = (Mat_SeqAIJ *)B->data;
684531e53bdSHong Zhang     nz    = spA->nz + spB->nz; /* total local nonzero entries of mat */
685531e53bdSHong Zhang     mem   = nz * (sizeof(PetscScalar) + sizeof(PetscInt)) + 3 * m * sizeof(PetscInt);
686531e53bdSHong Zhang     bcols = (PetscInt)(0.5 * mem / (m * sizeof(PetscScalar)));
687531e53bdSHong Zhang     brows = 1000 / bcols;
688531e53bdSHong Zhang     if (bcols > nis) bcols = nis;
689531e53bdSHong Zhang     if (brows == 0 || brows > m) brows = m;
690531e53bdSHong Zhang     c->brows = brows;
691531e53bdSHong Zhang     c->bcols = bcols;
692531e53bdSHong Zhang   }
693a8971b87SHong Zhang 
694a8971b87SHong Zhang   c->M       = mat->rmap->N / bs; /* set the global rows and columns and local rows */
695a8971b87SHong Zhang   c->N       = mat->cmap->N / bs;
696a8971b87SHong Zhang   c->m       = mat->rmap->n / bs;
697a8971b87SHong Zhang   c->rstart  = mat->rmap->rstart / bs;
698a8971b87SHong Zhang   c->ncolors = nis;
69993dfae19SHong Zhang   PetscFunctionReturn(0);
70093dfae19SHong Zhang }
701bdaf1daeSBarry Smith 
702bdaf1daeSBarry Smith /*@C
703bdaf1daeSBarry Smith 
70411a5261eSBarry Smith     MatFDColoringSetValues - takes a matrix in compressed color format and enters the matrix into a PETSc `Mat`
705bdaf1daeSBarry Smith 
706*c3339decSBarry Smith    Collective
707bdaf1daeSBarry Smith 
708bdaf1daeSBarry Smith    Input Parameters:
709bdaf1daeSBarry Smith +    J - the sparse matrix
71011a5261eSBarry Smith .    coloring - created with `MatFDColoringCreate()` and a local coloring
711bdaf1daeSBarry Smith -    y - column major storage of matrix values with one color of values per column, the number of rows of y should match
712bdaf1daeSBarry Smith          the number of local rows of J and the number of columns is the number of colors.
713bdaf1daeSBarry Smith 
714bdaf1daeSBarry Smith    Level: intermediate
715bdaf1daeSBarry Smith 
71611a5261eSBarry Smith    Notes: the matrix in compressed color format may come from an automatic differentiation code
717bdaf1daeSBarry Smith 
71811a5261eSBarry Smith    The code will be slightly faster if `MatFDColoringSetBlockSize`(coloring,`PETSC_DEFAULT`,nc); is called immediately after creating the coloring
719bdaf1daeSBarry Smith 
720db781477SPatrick Sanan .seealso: `MatFDColoringCreate()`, `ISColoring`, `ISColoringCreate()`, `ISColoringSetType()`, `IS_COLORING_LOCAL`, `MatFDColoringSetBlockSize()`
721bdaf1daeSBarry Smith @*/
722d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetValues(Mat J, MatFDColoring coloring, const PetscScalar *y)
723d71ae5a4SJacob Faibussowitsch {
724bdaf1daeSBarry Smith   MatEntry2      *Jentry2;
725bdaf1daeSBarry Smith   PetscInt        row, i, nrows_k, l, ncolors, nz = 0, bcols, nbcols = 0;
726bdaf1daeSBarry Smith   const PetscInt *nrows;
727bdaf1daeSBarry Smith   PetscBool       eq;
728bdaf1daeSBarry Smith 
729bdaf1daeSBarry Smith   PetscFunctionBegin;
730bdaf1daeSBarry Smith   PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
731bdaf1daeSBarry Smith   PetscValidHeaderSpecific(coloring, MAT_FDCOLORING_CLASSID, 2);
7329566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompareId((PetscObject)J, coloring->matid, &eq));
73328b400f6SJacob Faibussowitsch   PetscCheck(eq, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_WRONG, "Matrix used with MatFDColoringSetValues() must be that used with MatFDColoringCreate()");
734bdaf1daeSBarry Smith   Jentry2 = coloring->matentry2;
735bdaf1daeSBarry Smith   nrows   = coloring->nrows;
736bdaf1daeSBarry Smith   ncolors = coloring->ncolors;
737bdaf1daeSBarry Smith   bcols   = coloring->bcols;
738bdaf1daeSBarry Smith 
739bdaf1daeSBarry Smith   for (i = 0; i < ncolors; i += bcols) {
740bdaf1daeSBarry Smith     nrows_k = nrows[nbcols++];
741bdaf1daeSBarry Smith     for (l = 0; l < nrows_k; l++) {
742bdaf1daeSBarry Smith       row                      = Jentry2[nz].row; /* local row index */
743bdaf1daeSBarry Smith       *(Jentry2[nz++].valaddr) = y[row];
744bdaf1daeSBarry Smith     }
745bdaf1daeSBarry Smith     y += bcols * coloring->m;
746bdaf1daeSBarry Smith   }
7479566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
7489566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
749bdaf1daeSBarry Smith   PetscFunctionReturn(0);
750bdaf1daeSBarry Smith }
751