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)); 640d1c53f1SHong Zhang } 650d1c53f1SHong Zhang w3 = coloring->w3; 660d1c53f1SHong Zhang 679566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x1, &cstart, &cend)); /* used by ghosted vscale */ 6848a46eb9SPierre Jolivet if (vscale) PetscCall(VecGetArray(vscale, &vscale_array)); 690d1c53f1SHong Zhang nz = 0; 700d1c53f1SHong Zhang for (k = 0; k < ncolors; k++) { 710d1c53f1SHong Zhang coloring->currentcolor = k; 720d1c53f1SHong Zhang 730d1c53f1SHong Zhang /* 740d1c53f1SHong Zhang (3-1) Loop over each column associated with color 750d1c53f1SHong Zhang adding the perturbation to the vector w3 = x1 + dx. 760d1c53f1SHong Zhang */ 779566063dSJacob Faibussowitsch PetscCall(VecCopy(x1, w3)); 780d1c53f1SHong Zhang dy_i = dy; 790d1c53f1SHong Zhang for (i = 0; i < bs; i++) { /* Loop over a block of columns */ 809566063dSJacob Faibussowitsch PetscCall(VecGetArray(w3, &w3_array)); 810d1c53f1SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */ 820d1c53f1SHong Zhang if (coloring->htype[0] == 'w') { 830d1c53f1SHong Zhang for (l = 0; l < ncolumns[k]; l++) { 840d1c53f1SHong Zhang col = i + bs * coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 850d1c53f1SHong Zhang w3_array[col] += 1.0 / dx; 86684f2004SHong Zhang if (i) w3_array[col - 1] -= 1.0 / dx; /* resume original w3[col-1] */ 870d1c53f1SHong Zhang } 880d1c53f1SHong Zhang } else { /* htype == 'ds' */ 890d1c53f1SHong Zhang vscale_array -= cstart; /* shift pointer so global index can be used */ 900d1c53f1SHong Zhang for (l = 0; l < ncolumns[k]; l++) { 91f8c2866eSHong Zhang col = i + bs * coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 920d1c53f1SHong Zhang w3_array[col] += 1.0 / vscale_array[col]; 93684f2004SHong Zhang if (i) w3_array[col - 1] -= 1.0 / vscale_array[col - 1]; /* resume original w3[col-1] */ 940d1c53f1SHong Zhang } 950d1c53f1SHong Zhang vscale_array += cstart; 960d1c53f1SHong Zhang } 970d1c53f1SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array += cstart; 989566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(w3, &w3_array)); 990d1c53f1SHong Zhang 1000d1c53f1SHong Zhang /* 1010d1c53f1SHong Zhang (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 1020d1c53f1SHong Zhang w2 = F(x1 + dx) - F(x1) 1030d1c53f1SHong Zhang */ 1049566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MAT_FDColoringFunction, 0, 0, 0, 0)); 1059566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(w2, dy_i)); /* place w2 to the array dy_i */ 1069566063dSJacob Faibussowitsch PetscCall((*f)(sctx, w3, w2, fctx)); 1079566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MAT_FDColoringFunction, 0, 0, 0, 0)); 1089566063dSJacob Faibussowitsch PetscCall(VecAXPY(w2, -1.0, w1)); 1099566063dSJacob Faibussowitsch PetscCall(VecResetArray(w2)); 1100d1c53f1SHong Zhang dy_i += nxloc; /* points to dy+i*nxloc */ 1110d1c53f1SHong Zhang } 1120d1c53f1SHong Zhang 1130d1c53f1SHong Zhang /* 1140d1c53f1SHong Zhang (3-3) Loop over rows of vector, putting results into Jacobian matrix 1150d1c53f1SHong Zhang */ 1160d1c53f1SHong Zhang nrows_k = nrows[k]; 1170d1c53f1SHong Zhang if (coloring->htype[0] == 'w') { 1180d1c53f1SHong Zhang for (l = 0; l < nrows_k; l++) { 1190df34763SHong Zhang row = bs * Jentry2[nz].row; /* local row index */ 1200df34763SHong Zhang valaddr = Jentry2[nz++].valaddr; 1210d1c53f1SHong Zhang spidx = 0; 1220d1c53f1SHong Zhang dy_i = dy; 1230d1c53f1SHong Zhang for (i = 0; i < bs; i++) { /* column of the block */ 1240d1c53f1SHong Zhang for (j = 0; j < bs; j++) { /* row of the block */ 1250d1c53f1SHong Zhang valaddr[spidx++] = dy_i[row + j] * dx; 1260d1c53f1SHong Zhang } 127f8c2866eSHong Zhang dy_i += nxloc; /* points to dy+i*nxloc */ 1280d1c53f1SHong Zhang } 1290d1c53f1SHong Zhang } 1300d1c53f1SHong Zhang } else { /* htype == 'ds' */ 1310d1c53f1SHong Zhang for (l = 0; l < nrows_k; l++) { 1320d1c53f1SHong Zhang row = bs * Jentry[nz].row; /* local row index */ 1330d1c53f1SHong Zhang col = bs * Jentry[nz].col; /* local column index */ 134684f2004SHong Zhang valaddr = Jentry[nz++].valaddr; 135f8c2866eSHong Zhang spidx = 0; 136f8c2866eSHong Zhang dy_i = dy; 137f8c2866eSHong Zhang for (i = 0; i < bs; i++) { /* column of the block */ 138f8c2866eSHong Zhang for (j = 0; j < bs; j++) { /* row of the block */ 139f8c2866eSHong Zhang valaddr[spidx++] = dy_i[row + j] * vscale_array[col + i]; 140f8c2866eSHong Zhang } 141f8c2866eSHong Zhang dy_i += nxloc; /* points to dy+i*nxloc */ 142f8c2866eSHong Zhang } 1430d1c53f1SHong Zhang } 1440d1c53f1SHong Zhang } 1450d1c53f1SHong Zhang } 1469566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 1479566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 14848a46eb9SPierre Jolivet if (vscale) PetscCall(VecRestoreArray(vscale, &vscale_array)); 1490d1c53f1SHong Zhang 1500d1c53f1SHong Zhang coloring->currentcolor = -1; 1519566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(x1, PETSC_FALSE)); 1520d1c53f1SHong Zhang PetscFunctionReturn(0); 1530d1c53f1SHong Zhang } 154a64fbb32SBarry Smith 155531c7667SBarry Smith /* this is declared PETSC_EXTERN because it is used by MatFDColoringUseDM() which is in the DM library */ 1569371c9d4SSatish Balay PetscErrorCode MatFDColoringApply_AIJ(Mat J, MatFDColoring coloring, Vec x1, void *sctx) { 157fcd7ac73SHong Zhang PetscErrorCode (*f)(void *, Vec, Vec, void *) = (PetscErrorCode(*)(void *, Vec, Vec, void *))coloring->f; 158a2f2d239SHong Zhang PetscInt k, cstart, cend, l, row, col, nz; 1595edff71fSBarry Smith PetscScalar dx = 0.0, *y, *w3_array; 1605edff71fSBarry Smith const PetscScalar *xx; 161fcd7ac73SHong Zhang PetscScalar *vscale_array; 162fcd7ac73SHong Zhang PetscReal epsilon = coloring->error_rel, umin = coloring->umin, unorm; 1638bc97078SHong Zhang Vec w1 = coloring->w1, w2 = coloring->w2, w3, vscale = coloring->vscale; 164fcd7ac73SHong Zhang void *fctx = coloring->fctx; 16591e7fa0fSBarry Smith ISColoringType ctype = coloring->ctype; 16691e7fa0fSBarry Smith PetscInt nxloc, nrows_k; 167573f477fSHong Zhang MatEntry *Jentry = coloring->matentry; 1680df34763SHong Zhang MatEntry2 *Jentry2 = coloring->matentry2; 1698bc97078SHong Zhang const PetscInt ncolors = coloring->ncolors, *ncolumns = coloring->ncolumns, *nrows = coloring->nrows; 170b294de21SRichard Tran Mills PetscBool alreadyboundtocpu; 171fcd7ac73SHong Zhang 172fcd7ac73SHong Zhang PetscFunctionBegin; 1739566063dSJacob Faibussowitsch PetscCall(VecBoundToCPU(x1, &alreadyboundtocpu)); 1749566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(x1, PETSC_TRUE)); 17508401ef6SPierre Jolivet PetscCheck(!(ctype == IS_COLORING_LOCAL) || !(J->ops->fdcoloringapply == MatFDColoringApply_AIJ), PetscObjectComm((PetscObject)J), PETSC_ERR_SUP, "Must call MatColoringUseDM() with IS_COLORING_LOCAL"); 1768bc97078SHong Zhang /* (1) Set w1 = F(x1) */ 177fcd7ac73SHong Zhang if (!coloring->fset) { 1789566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MAT_FDColoringFunction, 0, 0, 0, 0)); 1799566063dSJacob Faibussowitsch PetscCall((*f)(sctx, x1, w1, fctx)); 1809566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MAT_FDColoringFunction, 0, 0, 0, 0)); 181fcd7ac73SHong Zhang } else { 182fcd7ac73SHong Zhang coloring->fset = PETSC_FALSE; 183fcd7ac73SHong Zhang } 184fcd7ac73SHong Zhang 1858bc97078SHong Zhang /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */ 186f6af9589SHong Zhang if (coloring->htype[0] == 'w') { 187684f2004SHong Zhang /* vscale = 1./dx is a constant scalar */ 1889566063dSJacob Faibussowitsch PetscCall(VecNorm(x1, NORM_2, &unorm)); 189c53567a0SHong Zhang dx = 1.0 / (PetscSqrtReal(1.0 + unorm) * epsilon); 19070e7395fSHong Zhang } else { 1919566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(x1, &nxloc)); 1929566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x1, &xx)); 1939566063dSJacob Faibussowitsch PetscCall(VecGetArray(vscale, &vscale_array)); 19474d3cef9SHong Zhang for (col = 0; col < nxloc; col++) { 195fcd7ac73SHong Zhang dx = xx[col]; 19674d3cef9SHong Zhang if (PetscAbsScalar(dx) < umin) { 19774d3cef9SHong Zhang if (PetscRealPart(dx) >= 0.0) dx = umin; 19874d3cef9SHong Zhang else if (PetscRealPart(dx) < 0.0) dx = -umin; 19974d3cef9SHong Zhang } 200fcd7ac73SHong Zhang dx *= epsilon; 20174d3cef9SHong Zhang vscale_array[col] = 1.0 / dx; 202f6af9589SHong Zhang } 2039566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x1, &xx)); 2049566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vscale, &vscale_array)); 20570e7395fSHong Zhang } 206684f2004SHong Zhang if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') { 2079566063dSJacob Faibussowitsch PetscCall(VecGhostUpdateBegin(vscale, INSERT_VALUES, SCATTER_FORWARD)); 2089566063dSJacob Faibussowitsch PetscCall(VecGhostUpdateEnd(vscale, INSERT_VALUES, SCATTER_FORWARD)); 209fcd7ac73SHong Zhang } 210fcd7ac73SHong Zhang 2118bc97078SHong Zhang /* (3) Loop over each color */ 212*4dfa11a4SJacob Faibussowitsch if (!coloring->w3) { PetscCall(VecDuplicate(x1, &coloring->w3)); } 2138bc97078SHong Zhang w3 = coloring->w3; 214fcd7ac73SHong Zhang 2159566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x1, &cstart, &cend)); /* used by ghosted vscale */ 21648a46eb9SPierre Jolivet if (vscale) PetscCall(VecGetArray(vscale, &vscale_array)); 2178bc97078SHong Zhang nz = 0; 218e2c857f8SHong Zhang 219b3e1f37bSHong Zhang if (coloring->bcols > 1) { /* use blocked insertion of Jentry */ 220b3e1f37bSHong Zhang PetscInt i, m = J->rmap->n, nbcols, bcols = coloring->bcols; 221e2c857f8SHong Zhang PetscScalar *dy = coloring->dy, *dy_k; 222e2c857f8SHong Zhang 223e2c857f8SHong Zhang nbcols = 0; 224e2c857f8SHong Zhang for (k = 0; k < ncolors; k += bcols) { 225e2c857f8SHong Zhang /* 226e2c857f8SHong Zhang (3-1) Loop over each column associated with color 227e2c857f8SHong Zhang adding the perturbation to the vector w3 = x1 + dx. 228e2c857f8SHong Zhang */ 229e2c857f8SHong Zhang 230e2c857f8SHong Zhang dy_k = dy; 231e2c857f8SHong Zhang if (k + bcols > ncolors) bcols = ncolors - k; 232e2c857f8SHong Zhang for (i = 0; i < bcols; i++) { 233ceef8a49SBarry Smith coloring->currentcolor = k + i; 234e2c857f8SHong Zhang 2359566063dSJacob Faibussowitsch PetscCall(VecCopy(x1, w3)); 2369566063dSJacob Faibussowitsch PetscCall(VecGetArray(w3, &w3_array)); 237e2c857f8SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */ 238e2c857f8SHong Zhang if (coloring->htype[0] == 'w') { 239e2c857f8SHong Zhang for (l = 0; l < ncolumns[k + i]; l++) { 240e2c857f8SHong Zhang col = coloring->columns[k + i][l]; /* local column (in global index!) of the matrix we are probing for */ 241e2c857f8SHong Zhang w3_array[col] += 1.0 / dx; 242e2c857f8SHong Zhang } 243e2c857f8SHong Zhang } else { /* htype == 'ds' */ 244e2c857f8SHong Zhang vscale_array -= cstart; /* shift pointer so global index can be used */ 245e2c857f8SHong Zhang for (l = 0; l < ncolumns[k + i]; l++) { 246e2c857f8SHong Zhang col = coloring->columns[k + i][l]; /* local column (in global index!) of the matrix we are probing for */ 247e2c857f8SHong Zhang w3_array[col] += 1.0 / vscale_array[col]; 248e2c857f8SHong Zhang } 249e2c857f8SHong Zhang vscale_array += cstart; 250e2c857f8SHong Zhang } 251e2c857f8SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array += cstart; 2529566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(w3, &w3_array)); 253e2c857f8SHong Zhang 254e2c857f8SHong Zhang /* 255e2c857f8SHong Zhang (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 256e2c857f8SHong Zhang w2 = F(x1 + dx) - F(x1) 257e2c857f8SHong Zhang */ 2589566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MAT_FDColoringFunction, 0, 0, 0, 0)); 2599566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(w2, dy_k)); /* place w2 to the array dy_i */ 2609566063dSJacob Faibussowitsch PetscCall((*f)(sctx, w3, w2, fctx)); 2619566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MAT_FDColoringFunction, 0, 0, 0, 0)); 2629566063dSJacob Faibussowitsch PetscCall(VecAXPY(w2, -1.0, w1)); 2639566063dSJacob Faibussowitsch PetscCall(VecResetArray(w2)); 264e2c857f8SHong Zhang dy_k += m; /* points to dy+i*nxloc */ 265e2c857f8SHong Zhang } 266e2c857f8SHong Zhang 267e2c857f8SHong Zhang /* 268e2c857f8SHong Zhang (3-3) Loop over block rows of vector, putting results into Jacobian matrix 269e2c857f8SHong Zhang */ 270d880da65SHong Zhang nrows_k = nrows[nbcols++]; 271d880da65SHong Zhang 272e2c857f8SHong Zhang if (coloring->htype[0] == 'w') { 273e2c857f8SHong Zhang for (l = 0; l < nrows_k; l++) { 2740df34763SHong Zhang row = Jentry2[nz].row; /* local row index */ 2755059b303SJunchao 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 2765059b303SJunchao Zhang another way, and it seems work. See https://lists.mcs.anl.gov/pipermail/petsc-users/2021-December/045158.html 2775059b303SJunchao Zhang */ 2785059b303SJunchao Zhang #if defined(PETSC_USE_COMPLEX) 2795059b303SJunchao Zhang PetscScalar *tmp = Jentry2[nz].valaddr; 2805059b303SJunchao Zhang *tmp = dy[row] * dx; 2815059b303SJunchao Zhang #else 2825059b303SJunchao Zhang *(Jentry2[nz].valaddr) = dy[row] * dx; 2835059b303SJunchao Zhang #endif 2845059b303SJunchao Zhang nz++; 285e2c857f8SHong Zhang } 286e2c857f8SHong Zhang } else { /* htype == 'ds' */ 287e2c857f8SHong Zhang for (l = 0; l < nrows_k; l++) { 288e2c857f8SHong Zhang row = Jentry[nz].row; /* local row index */ 2895059b303SJunchao Zhang #if defined(PETSC_USE_COMPLEX) /* See https://lists.mcs.anl.gov/pipermail/petsc-users/2021-December/045158.html */ 2905059b303SJunchao Zhang PetscScalar *tmp = Jentry[nz].valaddr; 2915059b303SJunchao Zhang *tmp = dy[row] * vscale_array[Jentry[nz].col]; 2925059b303SJunchao Zhang #else 293d880da65SHong Zhang *(Jentry[nz].valaddr) = dy[row] * vscale_array[Jentry[nz].col]; 2945059b303SJunchao Zhang #endif 295e2c857f8SHong Zhang nz++; 296e2c857f8SHong Zhang } 297e2c857f8SHong Zhang } 298e2c857f8SHong Zhang } 299b3e1f37bSHong Zhang } else { /* bcols == 1 */ 3008bc97078SHong Zhang for (k = 0; k < ncolors; k++) { 3018bc97078SHong Zhang coloring->currentcolor = k; 3029e917edbSHong Zhang 303fcd7ac73SHong Zhang /* 3048bc97078SHong Zhang (3-1) Loop over each column associated with color 305f6af9589SHong Zhang adding the perturbation to the vector w3 = x1 + dx. 306fcd7ac73SHong Zhang */ 3079566063dSJacob Faibussowitsch PetscCall(VecCopy(x1, w3)); 3089566063dSJacob Faibussowitsch PetscCall(VecGetArray(w3, &w3_array)); 309a2f2d239SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */ 310b039d6c4SHong Zhang if (coloring->htype[0] == 'w') { 3118bc97078SHong Zhang for (l = 0; l < ncolumns[k]; l++) { 312f6af9589SHong Zhang col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 3139e917edbSHong Zhang w3_array[col] += 1.0 / dx; 3149e917edbSHong Zhang } 315b039d6c4SHong Zhang } else { /* htype == 'ds' */ 316a2f2d239SHong Zhang vscale_array -= cstart; /* shift pointer so global index can be used */ 317b039d6c4SHong Zhang for (l = 0; l < ncolumns[k]; l++) { 318b039d6c4SHong Zhang col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 319a2f2d239SHong Zhang w3_array[col] += 1.0 / vscale_array[col]; 32070e7395fSHong Zhang } 321a2f2d239SHong Zhang vscale_array += cstart; 322fcd7ac73SHong Zhang } 323a2f2d239SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array += cstart; 3249566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(w3, &w3_array)); 3259e917edbSHong Zhang 326fcd7ac73SHong Zhang /* 3278bc97078SHong Zhang (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 328fcd7ac73SHong Zhang w2 = F(x1 + dx) - F(x1) 329fcd7ac73SHong Zhang */ 3309566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MAT_FDColoringFunction, 0, 0, 0, 0)); 3319566063dSJacob Faibussowitsch PetscCall((*f)(sctx, w3, w2, fctx)); 3329566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MAT_FDColoringFunction, 0, 0, 0, 0)); 3339566063dSJacob Faibussowitsch PetscCall(VecAXPY(w2, -1.0, w1)); 3349e917edbSHong Zhang 3358bc97078SHong Zhang /* 3368bc97078SHong Zhang (3-3) Loop over rows of vector, putting results into Jacobian matrix 3378bc97078SHong Zhang */ 338c53567a0SHong Zhang nrows_k = nrows[k]; 3399566063dSJacob Faibussowitsch PetscCall(VecGetArray(w2, &y)); 340b039d6c4SHong Zhang if (coloring->htype[0] == 'w') { 341c53567a0SHong Zhang for (l = 0; l < nrows_k; l++) { 3420df34763SHong Zhang row = Jentry2[nz].row; /* local row index */ 3435059b303SJunchao Zhang #if defined(PETSC_USE_COMPLEX) /* See https://lists.mcs.anl.gov/pipermail/petsc-users/2021-December/045158.html */ 3445059b303SJunchao Zhang PetscScalar *tmp = Jentry2[nz].valaddr; 3455059b303SJunchao Zhang *tmp = y[row] * dx; 3465059b303SJunchao Zhang #else 3475059b303SJunchao Zhang *(Jentry2[nz].valaddr) = y[row] * dx; 3485059b303SJunchao Zhang #endif 3495059b303SJunchao Zhang nz++; 3509e917edbSHong Zhang } 351b039d6c4SHong Zhang } else { /* htype == 'ds' */ 352c53567a0SHong Zhang for (l = 0; l < nrows_k; l++) { 353b039d6c4SHong Zhang row = Jentry[nz].row; /* local row index */ 3545059b303SJunchao Zhang #if defined(PETSC_USE_COMPLEX) /* See https://lists.mcs.anl.gov/pipermail/petsc-users/2021-December/045158.html */ 3555059b303SJunchao Zhang PetscScalar *tmp = Jentry[nz].valaddr; 3565059b303SJunchao Zhang *tmp = y[row] * vscale_array[Jentry[nz].col]; 3575059b303SJunchao Zhang #else 358b039d6c4SHong Zhang *(Jentry[nz].valaddr) = y[row] * vscale_array[Jentry[nz].col]; 3595059b303SJunchao Zhang #endif 360573f477fSHong Zhang nz++; 361fcd7ac73SHong Zhang } 362b039d6c4SHong Zhang } 3639566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(w2, &y)); 3648bc97078SHong Zhang } 365b3e1f37bSHong Zhang } 366b3e1f37bSHong Zhang 367e2cf4d64SStefano Zampini #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 368c70f7ee4SJunchao Zhang if (J->offloadmask != PETSC_OFFLOAD_UNALLOCATED) J->offloadmask = PETSC_OFFLOAD_CPU; 369e2cf4d64SStefano Zampini #endif 3709566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 3719566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 37248a46eb9SPierre Jolivet if (vscale) PetscCall(VecRestoreArray(vscale, &vscale_array)); 373fcd7ac73SHong Zhang coloring->currentcolor = -1; 3749566063dSJacob Faibussowitsch if (!alreadyboundtocpu) PetscCall(VecBindToCPU(x1, PETSC_FALSE)); 375fcd7ac73SHong Zhang PetscFunctionReturn(0); 376fcd7ac73SHong Zhang } 377fcd7ac73SHong Zhang 3789371c9d4SSatish Balay PetscErrorCode MatFDColoringSetUp_MPIXAIJ(Mat mat, ISColoring iscoloring, MatFDColoring c) { 379fcd7ac73SHong Zhang PetscMPIInt size, *ncolsonproc, *disp, nn; 380e3225e9dSHong Zhang PetscInt i, n, nrows, nrows_i, j, k, m, ncols, col, *rowhit, cstart, cend, colb; 38172c15787SHong Zhang const PetscInt *is, *A_ci, *A_cj, *B_ci, *B_cj, *row = NULL, *ltog = NULL; 382071fcb05SBarry Smith PetscInt nis = iscoloring->n, nctot, *cols, tmp = 0; 383fcd7ac73SHong Zhang ISLocalToGlobalMapping map = mat->cmap->mapping; 384a8971b87SHong Zhang PetscInt ctype = c->ctype, *spidxA, *spidxB, nz, bs, bs2, spidx; 3850d1c53f1SHong Zhang Mat A, B; 386e3225e9dSHong Zhang PetscScalar *A_val, *B_val, **valaddrhit; 387573f477fSHong Zhang MatEntry *Jentry; 3880df34763SHong Zhang MatEntry2 *Jentry2; 389d4002b98SHong Zhang PetscBool isBAIJ, isSELL; 390531e53bdSHong Zhang PetscInt bcols = c->bcols; 3910d1c53f1SHong Zhang #if defined(PETSC_USE_CTABLE) 3920d1c53f1SHong Zhang PetscTable colmap = NULL; 3930d1c53f1SHong Zhang #else 3940d1c53f1SHong Zhang PetscInt *colmap = NULL; /* local col number of off-diag col */ 3950d1c53f1SHong Zhang #endif 396fcd7ac73SHong Zhang 397fcd7ac73SHong Zhang PetscFunctionBegin; 3985bdb020cSBarry Smith if (ctype == IS_COLORING_LOCAL) { 39928b400f6SJacob Faibussowitsch PetscCheck(map, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_INCOMP, "When using ghosted differencing matrix must have local to global mapping provided with MatSetLocalToGlobalMapping"); 4009566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(map, <og)); 40172c15787SHong Zhang } 402fcd7ac73SHong Zhang 4039566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(mat, &bs)); 4049566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)mat, MATMPIBAIJ, &isBAIJ)); 4059566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATMPISELL, &isSELL)); 406e3225e9dSHong Zhang if (isBAIJ) { 407195687c5SHong Zhang Mat_MPIBAIJ *baij = (Mat_MPIBAIJ *)mat->data; 4086c145f34SHong Zhang Mat_SeqBAIJ *spA, *spB; 4099371c9d4SSatish Balay A = baij->A; 4109371c9d4SSatish Balay spA = (Mat_SeqBAIJ *)A->data; 4119371c9d4SSatish Balay A_val = spA->a; 4129371c9d4SSatish Balay B = baij->B; 4139371c9d4SSatish Balay spB = (Mat_SeqBAIJ *)B->data; 4149371c9d4SSatish Balay B_val = spB->a; 4150d1c53f1SHong Zhang nz = spA->nz + spB->nz; /* total nonzero entries of mat */ 41648a46eb9SPierre Jolivet if (!baij->colmap) PetscCall(MatCreateColmap_MPIBAIJ_Private(mat)); 417097e26fcSJed Brown colmap = baij->colmap; 4189566063dSJacob Faibussowitsch PetscCall(MatGetColumnIJ_SeqBAIJ_Color(A, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &A_ci, &A_cj, &spidxA, NULL)); 4199566063dSJacob Faibussowitsch PetscCall(MatGetColumnIJ_SeqBAIJ_Color(B, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &B_ci, &B_cj, &spidxB, NULL)); 420195687c5SHong Zhang 421195687c5SHong Zhang if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') { /* create vscale for storing dx */ 422195687c5SHong Zhang PetscInt *garray; 4239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(B->cmap->n, &garray)); 424195687c5SHong Zhang for (i = 0; i < baij->B->cmap->n / bs; i++) { 425ad540459SPierre Jolivet for (j = 0; j < bs; j++) garray[i * bs + j] = bs * baij->garray[i] + j; 426195687c5SHong Zhang } 4279566063dSJacob Faibussowitsch PetscCall(VecCreateGhost(PetscObjectComm((PetscObject)mat), mat->cmap->n, PETSC_DETERMINE, B->cmap->n, garray, &c->vscale)); 4289566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(c->vscale, PETSC_TRUE)); 4299566063dSJacob Faibussowitsch PetscCall(PetscFree(garray)); 430195687c5SHong Zhang } 431d4002b98SHong Zhang } else if (isSELL) { 432d4002b98SHong Zhang Mat_MPISELL *sell = (Mat_MPISELL *)mat->data; 433d4002b98SHong Zhang Mat_SeqSELL *spA, *spB; 4349371c9d4SSatish Balay A = sell->A; 4359371c9d4SSatish Balay spA = (Mat_SeqSELL *)A->data; 4369371c9d4SSatish Balay A_val = spA->val; 4379371c9d4SSatish Balay B = sell->B; 4389371c9d4SSatish Balay spB = (Mat_SeqSELL *)B->data; 4399371c9d4SSatish Balay B_val = spB->val; 440f4b713efSHong Zhang nz = spA->nz + spB->nz; /* total nonzero entries of mat */ 441d4002b98SHong Zhang if (!sell->colmap) { 442f4b713efSHong Zhang /* Allow access to data structures of local part of matrix 443f4b713efSHong Zhang - creates aij->colmap which maps global column number to local number in part B */ 4449566063dSJacob Faibussowitsch PetscCall(MatCreateColmap_MPISELL_Private(mat)); 445f4b713efSHong Zhang } 446d4002b98SHong Zhang colmap = sell->colmap; 4479566063dSJacob Faibussowitsch PetscCall(MatGetColumnIJ_SeqSELL_Color(A, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &A_ci, &A_cj, &spidxA, NULL)); 4489566063dSJacob Faibussowitsch PetscCall(MatGetColumnIJ_SeqSELL_Color(B, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &B_ci, &B_cj, &spidxB, NULL)); 449f4b713efSHong Zhang 450f4b713efSHong Zhang bs = 1; /* only bs=1 is supported for non MPIBAIJ matrix */ 451f4b713efSHong Zhang 452f4b713efSHong Zhang if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') { /* create vscale for storing dx */ 4539566063dSJacob Faibussowitsch PetscCall(VecCreateGhost(PetscObjectComm((PetscObject)mat), mat->cmap->n, PETSC_DETERMINE, B->cmap->n, sell->garray, &c->vscale)); 4549566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(c->vscale, PETSC_TRUE)); 455f4b713efSHong Zhang } 456e3225e9dSHong Zhang } else { 457e3225e9dSHong Zhang Mat_MPIAIJ *aij = (Mat_MPIAIJ *)mat->data; 458e3225e9dSHong Zhang Mat_SeqAIJ *spA, *spB; 4599371c9d4SSatish Balay A = aij->A; 4609371c9d4SSatish Balay spA = (Mat_SeqAIJ *)A->data; 4619371c9d4SSatish Balay A_val = spA->a; 4629371c9d4SSatish Balay B = aij->B; 4639371c9d4SSatish Balay spB = (Mat_SeqAIJ *)B->data; 4649371c9d4SSatish Balay B_val = spB->a; 465e3225e9dSHong Zhang nz = spA->nz + spB->nz; /* total nonzero entries of mat */ 466e3225e9dSHong Zhang if (!aij->colmap) { 467e3225e9dSHong Zhang /* Allow access to data structures of local part of matrix 468e3225e9dSHong Zhang - creates aij->colmap which maps global column number to local number in part B */ 4699566063dSJacob Faibussowitsch PetscCall(MatCreateColmap_MPIAIJ_Private(mat)); 470e3225e9dSHong Zhang } 471097e26fcSJed Brown colmap = aij->colmap; 4729566063dSJacob Faibussowitsch PetscCall(MatGetColumnIJ_SeqAIJ_Color(A, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &A_ci, &A_cj, &spidxA, NULL)); 4739566063dSJacob Faibussowitsch PetscCall(MatGetColumnIJ_SeqAIJ_Color(B, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &B_ci, &B_cj, &spidxB, NULL)); 474e3225e9dSHong Zhang 475e3225e9dSHong Zhang bs = 1; /* only bs=1 is supported for non MPIBAIJ matrix */ 476195687c5SHong Zhang 477195687c5SHong Zhang if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') { /* create vscale for storing dx */ 4789566063dSJacob Faibussowitsch PetscCall(VecCreateGhost(PetscObjectComm((PetscObject)mat), mat->cmap->n, PETSC_DETERMINE, B->cmap->n, aij->garray, &c->vscale)); 4799566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(c->vscale, PETSC_TRUE)); 480195687c5SHong Zhang } 4810d1c53f1SHong Zhang } 4820d1c53f1SHong Zhang 4830d1c53f1SHong Zhang m = mat->rmap->n / bs; 4840d1c53f1SHong Zhang cstart = mat->cmap->rstart / bs; 4850d1c53f1SHong Zhang cend = mat->cmap->rend / bs; 486fcd7ac73SHong Zhang 4879566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nis, &c->ncolumns, nis, &c->columns)); 4889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nis, &c->nrows)); 489fcd7ac73SHong Zhang 4900df34763SHong Zhang if (c->htype[0] == 'd') { 4919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nz, &Jentry)); 4928bc97078SHong Zhang c->matentry = Jentry; 4930df34763SHong Zhang } else if (c->htype[0] == 'w') { 4949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nz, &Jentry2)); 4950df34763SHong Zhang c->matentry2 = Jentry2; 4960df34763SHong Zhang } else SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "htype is not supported"); 497a774a6f1SHong Zhang 4989566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(m + 1, &rowhit, m + 1, &valaddrhit)); 499fcd7ac73SHong Zhang nz = 0; 5009566063dSJacob Faibussowitsch PetscCall(ISColoringGetIS(iscoloring, PETSC_OWN_POINTER, PETSC_IGNORE, &c->isa)); 501071fcb05SBarry Smith 502071fcb05SBarry Smith if (ctype == IS_COLORING_GLOBAL) { 5039566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &size)); 5049566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size, &ncolsonproc, size, &disp)); 505071fcb05SBarry Smith } 506071fcb05SBarry Smith 507d3825b63SHong Zhang for (i = 0; i < nis; i++) { /* for each local color */ 5089566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(c->isa[i], &n)); 5099566063dSJacob Faibussowitsch PetscCall(ISGetIndices(c->isa[i], &is)); 510fcd7ac73SHong Zhang 511fcd7ac73SHong Zhang c->ncolumns[i] = n; /* local number of columns of this color on this process */ 512071fcb05SBarry Smith c->columns[i] = (PetscInt *)is; 513fcd7ac73SHong Zhang 514fcd7ac73SHong Zhang if (ctype == IS_COLORING_GLOBAL) { 515d3825b63SHong Zhang /* Determine nctot, the total (parallel) number of columns of this color */ 516d3825b63SHong Zhang /* ncolsonproc[j]: local ncolumns on proc[j] of this color */ 5179566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(n, &nn)); 5189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&nn, 1, MPI_INT, ncolsonproc, 1, MPI_INT, PetscObjectComm((PetscObject)mat))); 5199371c9d4SSatish Balay nctot = 0; 5209371c9d4SSatish Balay for (j = 0; j < size; j++) nctot += ncolsonproc[j]; 52148a46eb9SPierre Jolivet if (!nctot) PetscCall(PetscInfo(mat, "Coloring of matrix has some unneeded colors with no corresponding rows\n")); 522fcd7ac73SHong Zhang 523fcd7ac73SHong Zhang disp[0] = 0; 524ad540459SPierre Jolivet for (j = 1; j < size; j++) disp[j] = disp[j - 1] + ncolsonproc[j - 1]; 525fcd7ac73SHong Zhang 526d3825b63SHong Zhang /* Get cols, the complete list of columns for this color on each process */ 5279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nctot + 1, &cols)); 5289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv((void *)is, n, MPIU_INT, cols, ncolsonproc, disp, MPIU_INT, PetscObjectComm((PetscObject)mat))); 5295bdb020cSBarry Smith } else if (ctype == IS_COLORING_LOCAL) { 530fcd7ac73SHong Zhang /* Determine local number of columns of this color on this process, including ghost points */ 531fcd7ac73SHong Zhang nctot = n; 532071fcb05SBarry Smith cols = (PetscInt *)is; 533fcd7ac73SHong Zhang } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not provided for this MatFDColoring type"); 534fcd7ac73SHong Zhang 5351b97d346SHong Zhang /* Mark all rows affect by these columns */ 5369566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(rowhit, m)); 5370d1c53f1SHong Zhang bs2 = bs * bs; 5384b2e90caSHong Zhang nrows_i = 0; 5391b97d346SHong Zhang for (j = 0; j < nctot; j++) { /* loop over columns*/ 5405bdb020cSBarry Smith if (ctype == IS_COLORING_LOCAL) { 541fcd7ac73SHong Zhang col = ltog[cols[j]]; 542fcd7ac73SHong Zhang } else { 543fcd7ac73SHong Zhang col = cols[j]; 544fcd7ac73SHong Zhang } 545e3225e9dSHong Zhang if (col >= cstart && col < cend) { /* column is in A, diagonal block of mat */ 546071fcb05SBarry Smith tmp = A_ci[col - cstart]; 547071fcb05SBarry Smith row = A_cj + tmp; 548071fcb05SBarry Smith nrows = A_ci[col - cstart + 1] - tmp; 5494b2e90caSHong Zhang nrows_i += nrows; 550071fcb05SBarry Smith 551fcd7ac73SHong Zhang /* loop over columns of A marking them in rowhit */ 552d3825b63SHong Zhang for (k = 0; k < nrows; k++) { 553d3825b63SHong Zhang /* set valaddrhit for part A */ 554071fcb05SBarry Smith spidx = bs2 * spidxA[tmp + k]; 555e3225e9dSHong Zhang valaddrhit[*row] = &A_val[spidx]; 556a774a6f1SHong Zhang rowhit[*row++] = col - cstart + 1; /* local column index */ 557fcd7ac73SHong Zhang } 558e3225e9dSHong Zhang } else { /* column is in B, off-diagonal block of mat */ 559fcd7ac73SHong Zhang #if defined(PETSC_USE_CTABLE) 5609566063dSJacob Faibussowitsch PetscCall(PetscTableFind(colmap, col + 1, &colb)); 561fcd7ac73SHong Zhang colb--; 562fcd7ac73SHong Zhang #else 5630d1c53f1SHong Zhang colb = colmap[col] - 1; /* local column index */ 564fcd7ac73SHong Zhang #endif 565fcd7ac73SHong Zhang if (colb == -1) { 566d3825b63SHong Zhang nrows = 0; 567fcd7ac73SHong Zhang } else { 5680d1c53f1SHong Zhang colb = colb / bs; 569071fcb05SBarry Smith tmp = B_ci[colb]; 570071fcb05SBarry Smith row = B_cj + tmp; 571071fcb05SBarry Smith nrows = B_ci[colb + 1] - tmp; 572fcd7ac73SHong Zhang } 5734b2e90caSHong Zhang nrows_i += nrows; 574fcd7ac73SHong Zhang /* loop over columns of B marking them in rowhit */ 575d3825b63SHong Zhang for (k = 0; k < nrows; k++) { 576d3825b63SHong Zhang /* set valaddrhit for part B */ 577071fcb05SBarry Smith spidx = bs2 * spidxB[tmp + k]; 578e3225e9dSHong Zhang valaddrhit[*row] = &B_val[spidx]; 57970e7395fSHong Zhang rowhit[*row++] = colb + 1 + cend - cstart; /* local column index */ 580fcd7ac73SHong Zhang } 581fcd7ac73SHong Zhang } 582fcd7ac73SHong Zhang } 5834b2e90caSHong Zhang c->nrows[i] = nrows_i; 5848bc97078SHong Zhang 5850df34763SHong Zhang if (c->htype[0] == 'd') { 586d3825b63SHong Zhang for (j = 0; j < m; j++) { 587fcd7ac73SHong Zhang if (rowhit[j]) { 588573f477fSHong Zhang Jentry[nz].row = j; /* local row index */ 589573f477fSHong Zhang Jentry[nz].col = rowhit[j] - 1; /* local column index */ 590573f477fSHong Zhang Jentry[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */ 591573f477fSHong Zhang nz++; 592fcd7ac73SHong Zhang } 593fcd7ac73SHong Zhang } 5940df34763SHong Zhang } else { /* c->htype == 'wp' */ 5950df34763SHong Zhang for (j = 0; j < m; j++) { 5960df34763SHong Zhang if (rowhit[j]) { 5970df34763SHong Zhang Jentry2[nz].row = j; /* local row index */ 5980df34763SHong Zhang Jentry2[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */ 5990df34763SHong Zhang nz++; 6000df34763SHong Zhang } 6010df34763SHong Zhang } 6020df34763SHong Zhang } 60348a46eb9SPierre Jolivet if (ctype == IS_COLORING_GLOBAL) PetscCall(PetscFree(cols)); 604fcd7ac73SHong Zhang } 60548a46eb9SPierre Jolivet if (ctype == IS_COLORING_GLOBAL) PetscCall(PetscFree2(ncolsonproc, disp)); 606fcd7ac73SHong Zhang 607a8971b87SHong Zhang if (bcols > 1) { /* reorder Jentry for faster MatFDColoringApply() */ 6089566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetUpBlocked_AIJ_Private(mat, c, nz)); 609531e53bdSHong Zhang } 610531e53bdSHong Zhang 6110d1c53f1SHong Zhang if (isBAIJ) { 6129566063dSJacob Faibussowitsch PetscCall(MatRestoreColumnIJ_SeqBAIJ_Color(A, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &A_ci, &A_cj, &spidxA, NULL)); 6139566063dSJacob Faibussowitsch PetscCall(MatRestoreColumnIJ_SeqBAIJ_Color(B, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &B_ci, &B_cj, &spidxB, NULL)); 6149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs * mat->rmap->n, &c->dy)); 615d4002b98SHong Zhang } else if (isSELL) { 6169566063dSJacob Faibussowitsch PetscCall(MatRestoreColumnIJ_SeqSELL_Color(A, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &A_ci, &A_cj, &spidxA, NULL)); 6179566063dSJacob Faibussowitsch PetscCall(MatRestoreColumnIJ_SeqSELL_Color(B, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &B_ci, &B_cj, &spidxB, NULL)); 6180d1c53f1SHong Zhang } else { 6199566063dSJacob Faibussowitsch PetscCall(MatRestoreColumnIJ_SeqAIJ_Color(A, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &A_ci, &A_cj, &spidxA, NULL)); 6209566063dSJacob Faibussowitsch PetscCall(MatRestoreColumnIJ_SeqAIJ_Color(B, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &B_ci, &B_cj, &spidxB, NULL)); 6210d1c53f1SHong Zhang } 6220d1c53f1SHong Zhang 6239566063dSJacob Faibussowitsch PetscCall(ISColoringRestoreIS(iscoloring, PETSC_OWN_POINTER, &c->isa)); 6249566063dSJacob Faibussowitsch PetscCall(PetscFree2(rowhit, valaddrhit)); 625a8971b87SHong Zhang 62648a46eb9SPierre Jolivet if (ctype == IS_COLORING_LOCAL) PetscCall(ISLocalToGlobalMappingRestoreIndices(map, <og)); 6279566063dSJacob Faibussowitsch PetscCall(PetscInfo(c, "ncolors %" PetscInt_FMT ", brows %" PetscInt_FMT " and bcols %" PetscInt_FMT " are used.\n", c->ncolors, c->brows, c->bcols)); 628fcd7ac73SHong Zhang PetscFunctionReturn(0); 629fcd7ac73SHong Zhang } 63093dfae19SHong Zhang 6319371c9d4SSatish Balay PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat mat, ISColoring iscoloring, MatFDColoring c) { 632a8971b87SHong Zhang PetscInt bs, nis = iscoloring->n, m = mat->rmap->n; 633d4002b98SHong Zhang PetscBool isBAIJ, isSELL; 634531e53bdSHong Zhang 63593dfae19SHong Zhang PetscFunctionBegin; 636531e53bdSHong Zhang /* set default brows and bcols for speedup inserting the dense matrix into sparse Jacobian; 637531e53bdSHong Zhang bcols is chosen s.t. dy-array takes 50% of memory space as mat */ 6389566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(mat, &bs)); 6399566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)mat, MATMPIBAIJ, &isBAIJ)); 6409566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATMPISELL, &isSELL)); 641b8dfa574SBarry Smith if (isBAIJ || m == 0) { 642a8971b87SHong Zhang c->brows = m; 643a8971b87SHong Zhang c->bcols = 1; 644d4002b98SHong Zhang } else if (isSELL) { 645f4b713efSHong Zhang /* bcols is chosen s.t. dy-array takes 50% of local memory space as mat */ 646d4002b98SHong Zhang Mat_MPISELL *sell = (Mat_MPISELL *)mat->data; 647d4002b98SHong Zhang Mat_SeqSELL *spA, *spB; 648f4b713efSHong Zhang Mat A, B; 649f4b713efSHong Zhang PetscInt nz, brows, bcols; 650f4b713efSHong Zhang PetscReal mem; 651f4b713efSHong Zhang 652d4002b98SHong Zhang bs = 1; /* only bs=1 is supported for MPISELL matrix */ 653f4b713efSHong Zhang 6549371c9d4SSatish Balay A = sell->A; 6559371c9d4SSatish Balay spA = (Mat_SeqSELL *)A->data; 6569371c9d4SSatish Balay B = sell->B; 6579371c9d4SSatish Balay spB = (Mat_SeqSELL *)B->data; 658f4b713efSHong Zhang nz = spA->nz + spB->nz; /* total local nonzero entries of mat */ 659f4b713efSHong Zhang mem = nz * (sizeof(PetscScalar) + sizeof(PetscInt)) + 3 * m * sizeof(PetscInt); 660f4b713efSHong Zhang bcols = (PetscInt)(0.5 * mem / (m * sizeof(PetscScalar))); 661f4b713efSHong Zhang brows = 1000 / bcols; 662f4b713efSHong Zhang if (bcols > nis) bcols = nis; 663f4b713efSHong Zhang if (brows == 0 || brows > m) brows = m; 664f4b713efSHong Zhang c->brows = brows; 665f4b713efSHong Zhang c->bcols = bcols; 666531e53bdSHong Zhang } else { /* mpiaij matrix */ 667531e53bdSHong Zhang /* bcols is chosen s.t. dy-array takes 50% of local memory space as mat */ 668531e53bdSHong Zhang Mat_MPIAIJ *aij = (Mat_MPIAIJ *)mat->data; 669531e53bdSHong Zhang Mat_SeqAIJ *spA, *spB; 670531e53bdSHong Zhang Mat A, B; 671a8971b87SHong Zhang PetscInt nz, brows, bcols; 672531e53bdSHong Zhang PetscReal mem; 673531e53bdSHong Zhang 674531e53bdSHong Zhang bs = 1; /* only bs=1 is supported for MPIAIJ matrix */ 675531e53bdSHong Zhang 6769371c9d4SSatish Balay A = aij->A; 6779371c9d4SSatish Balay spA = (Mat_SeqAIJ *)A->data; 6789371c9d4SSatish Balay B = aij->B; 6799371c9d4SSatish Balay spB = (Mat_SeqAIJ *)B->data; 680531e53bdSHong Zhang nz = spA->nz + spB->nz; /* total local nonzero entries of mat */ 681531e53bdSHong Zhang mem = nz * (sizeof(PetscScalar) + sizeof(PetscInt)) + 3 * m * sizeof(PetscInt); 682531e53bdSHong Zhang bcols = (PetscInt)(0.5 * mem / (m * sizeof(PetscScalar))); 683531e53bdSHong Zhang brows = 1000 / bcols; 684531e53bdSHong Zhang if (bcols > nis) bcols = nis; 685531e53bdSHong Zhang if (brows == 0 || brows > m) brows = m; 686531e53bdSHong Zhang c->brows = brows; 687531e53bdSHong Zhang c->bcols = bcols; 688531e53bdSHong Zhang } 689a8971b87SHong Zhang 690a8971b87SHong Zhang c->M = mat->rmap->N / bs; /* set the global rows and columns and local rows */ 691a8971b87SHong Zhang c->N = mat->cmap->N / bs; 692a8971b87SHong Zhang c->m = mat->rmap->n / bs; 693a8971b87SHong Zhang c->rstart = mat->rmap->rstart / bs; 694a8971b87SHong Zhang c->ncolors = nis; 69593dfae19SHong Zhang PetscFunctionReturn(0); 69693dfae19SHong Zhang } 697bdaf1daeSBarry Smith 698bdaf1daeSBarry Smith /*@C 699bdaf1daeSBarry Smith 70011a5261eSBarry Smith MatFDColoringSetValues - takes a matrix in compressed color format and enters the matrix into a PETSc `Mat` 701bdaf1daeSBarry Smith 702bdaf1daeSBarry Smith Collective on J 703bdaf1daeSBarry Smith 704bdaf1daeSBarry Smith Input Parameters: 705bdaf1daeSBarry Smith + J - the sparse matrix 70611a5261eSBarry Smith . coloring - created with `MatFDColoringCreate()` and a local coloring 707bdaf1daeSBarry Smith - y - column major storage of matrix values with one color of values per column, the number of rows of y should match 708bdaf1daeSBarry Smith the number of local rows of J and the number of columns is the number of colors. 709bdaf1daeSBarry Smith 710bdaf1daeSBarry Smith Level: intermediate 711bdaf1daeSBarry Smith 71211a5261eSBarry Smith Notes: the matrix in compressed color format may come from an automatic differentiation code 713bdaf1daeSBarry Smith 71411a5261eSBarry Smith The code will be slightly faster if `MatFDColoringSetBlockSize`(coloring,`PETSC_DEFAULT`,nc); is called immediately after creating the coloring 715bdaf1daeSBarry Smith 716db781477SPatrick Sanan .seealso: `MatFDColoringCreate()`, `ISColoring`, `ISColoringCreate()`, `ISColoringSetType()`, `IS_COLORING_LOCAL`, `MatFDColoringSetBlockSize()` 717bdaf1daeSBarry Smith @*/ 7189371c9d4SSatish Balay PetscErrorCode MatFDColoringSetValues(Mat J, MatFDColoring coloring, const PetscScalar *y) { 719bdaf1daeSBarry Smith MatEntry2 *Jentry2; 720bdaf1daeSBarry Smith PetscInt row, i, nrows_k, l, ncolors, nz = 0, bcols, nbcols = 0; 721bdaf1daeSBarry Smith const PetscInt *nrows; 722bdaf1daeSBarry Smith PetscBool eq; 723bdaf1daeSBarry Smith 724bdaf1daeSBarry Smith PetscFunctionBegin; 725bdaf1daeSBarry Smith PetscValidHeaderSpecific(J, MAT_CLASSID, 1); 726bdaf1daeSBarry Smith PetscValidHeaderSpecific(coloring, MAT_FDCOLORING_CLASSID, 2); 7279566063dSJacob Faibussowitsch PetscCall(PetscObjectCompareId((PetscObject)J, coloring->matid, &eq)); 72828b400f6SJacob Faibussowitsch PetscCheck(eq, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_WRONG, "Matrix used with MatFDColoringSetValues() must be that used with MatFDColoringCreate()"); 729bdaf1daeSBarry Smith Jentry2 = coloring->matentry2; 730bdaf1daeSBarry Smith nrows = coloring->nrows; 731bdaf1daeSBarry Smith ncolors = coloring->ncolors; 732bdaf1daeSBarry Smith bcols = coloring->bcols; 733bdaf1daeSBarry Smith 734bdaf1daeSBarry Smith for (i = 0; i < ncolors; i += bcols) { 735bdaf1daeSBarry Smith nrows_k = nrows[nbcols++]; 736bdaf1daeSBarry Smith for (l = 0; l < nrows_k; l++) { 737bdaf1daeSBarry Smith row = Jentry2[nz].row; /* local row index */ 738bdaf1daeSBarry Smith *(Jentry2[nz++].valaddr) = y[row]; 739bdaf1daeSBarry Smith } 740bdaf1daeSBarry Smith y += bcols * coloring->m; 741bdaf1daeSBarry Smith } 7429566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 7439566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 744bdaf1daeSBarry Smith PetscFunctionReturn(0); 745bdaf1daeSBarry Smith } 746