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 6589ce454SStefano Zampini static PetscErrorCode MatFDColoringMarkHost_AIJ(Mat J) 7589ce454SStefano Zampini { 862cf5dd6SStefano Zampini PetscBool isseqAIJ, ismpiAIJ, issell; 9589ce454SStefano Zampini PetscScalar *v; 10589ce454SStefano Zampini 11589ce454SStefano Zampini PetscFunctionBegin; 12589ce454SStefano Zampini PetscCall(PetscObjectBaseTypeCompare((PetscObject)J, MATMPIAIJ, &ismpiAIJ)); 13589ce454SStefano Zampini PetscCall(PetscObjectBaseTypeCompare((PetscObject)J, MATSEQAIJ, &isseqAIJ)); 1462cf5dd6SStefano Zampini PetscCall(PetscObjectTypeCompareAny((PetscObject)J, &issell, MATSEQSELLCUDA, MATMPISELLCUDA, "")); 1562cf5dd6SStefano Zampini PetscCheck(!issell, PETSC_COMM_SELF, PETSC_ERR_SUP, "Not coded for %s. Send an email to petsc-dev@mcs.anl.gov to request this feature", ((PetscObject)J)->type_name); 16589ce454SStefano Zampini if (isseqAIJ) { 17589ce454SStefano Zampini PetscCall(MatSeqAIJGetArrayWrite(J, &v)); 18589ce454SStefano Zampini PetscCall(MatSeqAIJRestoreArrayWrite(J, &v)); 19589ce454SStefano Zampini } else if (ismpiAIJ) { 20589ce454SStefano Zampini Mat dJ, oJ; 21589ce454SStefano Zampini 22589ce454SStefano Zampini PetscCall(MatMPIAIJGetSeqAIJ(J, &dJ, &oJ, NULL)); 23589ce454SStefano Zampini PetscCall(MatSeqAIJGetArrayWrite(dJ, &v)); 24589ce454SStefano Zampini PetscCall(MatSeqAIJRestoreArrayWrite(dJ, &v)); 25589ce454SStefano Zampini PetscCall(MatSeqAIJGetArrayWrite(oJ, &v)); 26589ce454SStefano Zampini PetscCall(MatSeqAIJRestoreArrayWrite(oJ, &v)); 27589ce454SStefano Zampini } 28589ce454SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 29589ce454SStefano Zampini } 30589ce454SStefano Zampini 31d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringApply_BAIJ(Mat J, MatFDColoring coloring, Vec x1, void *sctx) 32d71ae5a4SJacob Faibussowitsch { 330d1c53f1SHong Zhang PetscErrorCode (*f)(void *, Vec, Vec, void *) = (PetscErrorCode(*)(void *, Vec, Vec, void *))coloring->f; 340d1c53f1SHong Zhang PetscInt k, cstart, cend, l, row, col, nz, spidx, i, j; 355edff71fSBarry Smith PetscScalar dx = 0.0, *w3_array, *dy_i, *dy = coloring->dy; 360d1c53f1SHong Zhang PetscScalar *vscale_array; 375edff71fSBarry Smith const PetscScalar *xx; 380d1c53f1SHong Zhang PetscReal epsilon = coloring->error_rel, umin = coloring->umin, unorm; 390d1c53f1SHong Zhang Vec w1 = coloring->w1, w2 = coloring->w2, w3, vscale = coloring->vscale; 400d1c53f1SHong Zhang void *fctx = coloring->fctx; 410d1c53f1SHong Zhang PetscInt ctype = coloring->ctype, nxloc, nrows_k; 420d1c53f1SHong Zhang PetscScalar *valaddr; 430d1c53f1SHong Zhang MatEntry *Jentry = coloring->matentry; 440df34763SHong Zhang MatEntry2 *Jentry2 = coloring->matentry2; 450d1c53f1SHong Zhang const PetscInt ncolors = coloring->ncolors, *ncolumns = coloring->ncolumns, *nrows = coloring->nrows; 460d1c53f1SHong Zhang PetscInt bs = J->rmap->bs; 470d1c53f1SHong Zhang 480d1c53f1SHong Zhang PetscFunctionBegin; 499566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(x1, PETSC_TRUE)); 500d1c53f1SHong Zhang /* (1) Set w1 = F(x1) */ 510d1c53f1SHong Zhang if (!coloring->fset) { 529566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MAT_FDColoringFunction, coloring, 0, 0, 0)); 539566063dSJacob Faibussowitsch PetscCall((*f)(sctx, x1, w1, fctx)); 549566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MAT_FDColoringFunction, coloring, 0, 0, 0)); 550d1c53f1SHong Zhang } else { 560d1c53f1SHong Zhang coloring->fset = PETSC_FALSE; 570d1c53f1SHong Zhang } 580d1c53f1SHong Zhang 590d1c53f1SHong Zhang /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */ 609566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(x1, &nxloc)); 610d1c53f1SHong Zhang if (coloring->htype[0] == 'w') { 620d1c53f1SHong Zhang /* vscale = dx is a constant scalar */ 639566063dSJacob Faibussowitsch PetscCall(VecNorm(x1, NORM_2, &unorm)); 640d1c53f1SHong Zhang dx = 1.0 / (PetscSqrtReal(1.0 + unorm) * epsilon); 650d1c53f1SHong Zhang } else { 669566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x1, &xx)); 679566063dSJacob Faibussowitsch PetscCall(VecGetArray(vscale, &vscale_array)); 680d1c53f1SHong Zhang for (col = 0; col < nxloc; col++) { 690d1c53f1SHong Zhang dx = xx[col]; 700d1c53f1SHong Zhang if (PetscAbsScalar(dx) < umin) { 710d1c53f1SHong Zhang if (PetscRealPart(dx) >= 0.0) dx = umin; 720d1c53f1SHong Zhang else if (PetscRealPart(dx) < 0.0) dx = -umin; 730d1c53f1SHong Zhang } 740d1c53f1SHong Zhang dx *= epsilon; 750d1c53f1SHong Zhang vscale_array[col] = 1.0 / dx; 760d1c53f1SHong Zhang } 779566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x1, &xx)); 789566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vscale, &vscale_array)); 790d1c53f1SHong Zhang } 80684f2004SHong Zhang if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') { 819566063dSJacob Faibussowitsch PetscCall(VecGhostUpdateBegin(vscale, INSERT_VALUES, SCATTER_FORWARD)); 829566063dSJacob Faibussowitsch PetscCall(VecGhostUpdateEnd(vscale, INSERT_VALUES, SCATTER_FORWARD)); 830d1c53f1SHong Zhang } 840d1c53f1SHong Zhang 850d1c53f1SHong Zhang /* (3) Loop over each color */ 860d1c53f1SHong Zhang if (!coloring->w3) { 879566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x1, &coloring->w3)); 88ce911d59SRichard Tran Mills /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */ 899566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(coloring->w3, PETSC_TRUE)); 900d1c53f1SHong Zhang } 910d1c53f1SHong Zhang w3 = coloring->w3; 920d1c53f1SHong Zhang 939566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x1, &cstart, &cend)); /* used by ghosted vscale */ 9448a46eb9SPierre Jolivet if (vscale) PetscCall(VecGetArray(vscale, &vscale_array)); 950d1c53f1SHong Zhang nz = 0; 960d1c53f1SHong Zhang for (k = 0; k < ncolors; k++) { 970d1c53f1SHong Zhang coloring->currentcolor = k; 980d1c53f1SHong Zhang 990d1c53f1SHong Zhang /* 1000d1c53f1SHong Zhang (3-1) Loop over each column associated with color 1010d1c53f1SHong Zhang adding the perturbation to the vector w3 = x1 + dx. 1020d1c53f1SHong Zhang */ 1039566063dSJacob Faibussowitsch PetscCall(VecCopy(x1, w3)); 1040d1c53f1SHong Zhang dy_i = dy; 1050d1c53f1SHong Zhang for (i = 0; i < bs; i++) { /* Loop over a block of columns */ 1069566063dSJacob Faibussowitsch PetscCall(VecGetArray(w3, &w3_array)); 1070d1c53f1SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */ 1080d1c53f1SHong Zhang if (coloring->htype[0] == 'w') { 1090d1c53f1SHong Zhang for (l = 0; l < ncolumns[k]; l++) { 1100d1c53f1SHong Zhang col = i + bs * coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 1110d1c53f1SHong Zhang w3_array[col] += 1.0 / dx; 112684f2004SHong Zhang if (i) w3_array[col - 1] -= 1.0 / dx; /* resume original w3[col-1] */ 1130d1c53f1SHong Zhang } 1140d1c53f1SHong Zhang } else { /* htype == 'ds' */ 1150d1c53f1SHong Zhang vscale_array -= cstart; /* shift pointer so global index can be used */ 1160d1c53f1SHong Zhang for (l = 0; l < ncolumns[k]; l++) { 117f8c2866eSHong Zhang col = i + bs * coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 1180d1c53f1SHong Zhang w3_array[col] += 1.0 / vscale_array[col]; 119684f2004SHong Zhang if (i) w3_array[col - 1] -= 1.0 / vscale_array[col - 1]; /* resume original w3[col-1] */ 1200d1c53f1SHong Zhang } 1210d1c53f1SHong Zhang vscale_array += cstart; 1220d1c53f1SHong Zhang } 1230d1c53f1SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array += cstart; 1249566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(w3, &w3_array)); 1250d1c53f1SHong Zhang 1260d1c53f1SHong Zhang /* 1270d1c53f1SHong Zhang (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 1280d1c53f1SHong Zhang w2 = F(x1 + dx) - F(x1) 1290d1c53f1SHong Zhang */ 1309566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MAT_FDColoringFunction, 0, 0, 0, 0)); 1319566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(w2, dy_i)); /* place w2 to the array dy_i */ 1329566063dSJacob Faibussowitsch PetscCall((*f)(sctx, w3, w2, fctx)); 1339566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MAT_FDColoringFunction, 0, 0, 0, 0)); 1349566063dSJacob Faibussowitsch PetscCall(VecAXPY(w2, -1.0, w1)); 1359566063dSJacob Faibussowitsch PetscCall(VecResetArray(w2)); 1360d1c53f1SHong Zhang dy_i += nxloc; /* points to dy+i*nxloc */ 1370d1c53f1SHong Zhang } 1380d1c53f1SHong Zhang 1390d1c53f1SHong Zhang /* 1400d1c53f1SHong Zhang (3-3) Loop over rows of vector, putting results into Jacobian matrix 1410d1c53f1SHong Zhang */ 1420d1c53f1SHong Zhang nrows_k = nrows[k]; 1430d1c53f1SHong Zhang if (coloring->htype[0] == 'w') { 1440d1c53f1SHong Zhang for (l = 0; l < nrows_k; l++) { 1450df34763SHong Zhang row = bs * Jentry2[nz].row; /* local row index */ 1460df34763SHong Zhang valaddr = Jentry2[nz++].valaddr; 1470d1c53f1SHong Zhang spidx = 0; 1480d1c53f1SHong Zhang dy_i = dy; 1490d1c53f1SHong Zhang for (i = 0; i < bs; i++) { /* column of the block */ 1500d1c53f1SHong Zhang for (j = 0; j < bs; j++) { /* row of the block */ 1510d1c53f1SHong Zhang valaddr[spidx++] = dy_i[row + j] * dx; 1520d1c53f1SHong Zhang } 153f8c2866eSHong Zhang dy_i += nxloc; /* points to dy+i*nxloc */ 1540d1c53f1SHong Zhang } 1550d1c53f1SHong Zhang } 1560d1c53f1SHong Zhang } else { /* htype == 'ds' */ 1570d1c53f1SHong Zhang for (l = 0; l < nrows_k; l++) { 1580d1c53f1SHong Zhang row = bs * Jentry[nz].row; /* local row index */ 1590d1c53f1SHong Zhang col = bs * Jentry[nz].col; /* local column index */ 160684f2004SHong Zhang valaddr = Jentry[nz++].valaddr; 161f8c2866eSHong Zhang spidx = 0; 162f8c2866eSHong Zhang dy_i = dy; 163f8c2866eSHong Zhang for (i = 0; i < bs; i++) { /* column of the block */ 164f8c2866eSHong Zhang for (j = 0; j < bs; j++) { /* row of the block */ 165f8c2866eSHong Zhang valaddr[spidx++] = dy_i[row + j] * vscale_array[col + i]; 166f8c2866eSHong Zhang } 167f8c2866eSHong Zhang dy_i += nxloc; /* points to dy+i*nxloc */ 168f8c2866eSHong Zhang } 1690d1c53f1SHong Zhang } 1700d1c53f1SHong Zhang } 1710d1c53f1SHong Zhang } 1729566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 1739566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 17448a46eb9SPierre Jolivet if (vscale) PetscCall(VecRestoreArray(vscale, &vscale_array)); 1750d1c53f1SHong Zhang 1760d1c53f1SHong Zhang coloring->currentcolor = -1; 1779566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(x1, PETSC_FALSE)); 1783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1790d1c53f1SHong Zhang } 180a64fbb32SBarry Smith 181531c7667SBarry Smith /* this is declared PETSC_EXTERN because it is used by MatFDColoringUseDM() which is in the DM library */ 182d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringApply_AIJ(Mat J, MatFDColoring coloring, Vec x1, void *sctx) 183d71ae5a4SJacob Faibussowitsch { 184fcd7ac73SHong Zhang PetscErrorCode (*f)(void *, Vec, Vec, void *) = (PetscErrorCode(*)(void *, Vec, Vec, void *))coloring->f; 185a2f2d239SHong Zhang PetscInt k, cstart, cend, l, row, col, nz; 1865edff71fSBarry Smith PetscScalar dx = 0.0, *y, *w3_array; 1875edff71fSBarry Smith const PetscScalar *xx; 188fcd7ac73SHong Zhang PetscScalar *vscale_array; 189fcd7ac73SHong Zhang PetscReal epsilon = coloring->error_rel, umin = coloring->umin, unorm; 1908bc97078SHong Zhang Vec w1 = coloring->w1, w2 = coloring->w2, w3, vscale = coloring->vscale; 191fcd7ac73SHong Zhang void *fctx = coloring->fctx; 19291e7fa0fSBarry Smith ISColoringType ctype = coloring->ctype; 19391e7fa0fSBarry Smith PetscInt nxloc, nrows_k; 194573f477fSHong Zhang MatEntry *Jentry = coloring->matentry; 1950df34763SHong Zhang MatEntry2 *Jentry2 = coloring->matentry2; 1968bc97078SHong Zhang const PetscInt ncolors = coloring->ncolors, *ncolumns = coloring->ncolumns, *nrows = coloring->nrows; 197b294de21SRichard Tran Mills PetscBool alreadyboundtocpu; 198fcd7ac73SHong Zhang 199fcd7ac73SHong Zhang PetscFunctionBegin; 200589ce454SStefano Zampini PetscCall(MatFDColoringMarkHost_AIJ(J)); 2019566063dSJacob Faibussowitsch PetscCall(VecBoundToCPU(x1, &alreadyboundtocpu)); 2029566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(x1, PETSC_TRUE)); 20308401ef6SPierre Jolivet PetscCheck(!(ctype == IS_COLORING_LOCAL) || !(J->ops->fdcoloringapply == MatFDColoringApply_AIJ), PetscObjectComm((PetscObject)J), PETSC_ERR_SUP, "Must call MatColoringUseDM() with IS_COLORING_LOCAL"); 2048bc97078SHong Zhang /* (1) Set w1 = F(x1) */ 205fcd7ac73SHong Zhang if (!coloring->fset) { 2069566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MAT_FDColoringFunction, 0, 0, 0, 0)); 2079566063dSJacob Faibussowitsch PetscCall((*f)(sctx, x1, w1, fctx)); 2089566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MAT_FDColoringFunction, 0, 0, 0, 0)); 209fcd7ac73SHong Zhang } else { 210fcd7ac73SHong Zhang coloring->fset = PETSC_FALSE; 211fcd7ac73SHong Zhang } 212fcd7ac73SHong Zhang 2138bc97078SHong Zhang /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */ 214f6af9589SHong Zhang if (coloring->htype[0] == 'w') { 215684f2004SHong Zhang /* vscale = 1./dx is a constant scalar */ 2169566063dSJacob Faibussowitsch PetscCall(VecNorm(x1, NORM_2, &unorm)); 217c53567a0SHong Zhang dx = 1.0 / (PetscSqrtReal(1.0 + unorm) * epsilon); 21870e7395fSHong Zhang } else { 2199566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(x1, &nxloc)); 2209566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x1, &xx)); 2219566063dSJacob Faibussowitsch PetscCall(VecGetArray(vscale, &vscale_array)); 22274d3cef9SHong Zhang for (col = 0; col < nxloc; col++) { 223fcd7ac73SHong Zhang dx = xx[col]; 22474d3cef9SHong Zhang if (PetscAbsScalar(dx) < umin) { 22574d3cef9SHong Zhang if (PetscRealPart(dx) >= 0.0) dx = umin; 22674d3cef9SHong Zhang else if (PetscRealPart(dx) < 0.0) dx = -umin; 22774d3cef9SHong Zhang } 228fcd7ac73SHong Zhang dx *= epsilon; 22974d3cef9SHong Zhang vscale_array[col] = 1.0 / dx; 230f6af9589SHong Zhang } 2319566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x1, &xx)); 2329566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vscale, &vscale_array)); 23370e7395fSHong Zhang } 234684f2004SHong Zhang if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') { 2359566063dSJacob Faibussowitsch PetscCall(VecGhostUpdateBegin(vscale, INSERT_VALUES, SCATTER_FORWARD)); 2369566063dSJacob Faibussowitsch PetscCall(VecGhostUpdateEnd(vscale, INSERT_VALUES, SCATTER_FORWARD)); 237fcd7ac73SHong Zhang } 238fcd7ac73SHong Zhang 2398bc97078SHong Zhang /* (3) Loop over each color */ 240aa624791SPierre Jolivet if (!coloring->w3) PetscCall(VecDuplicate(x1, &coloring->w3)); 2418bc97078SHong Zhang w3 = coloring->w3; 242fcd7ac73SHong Zhang 2439566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x1, &cstart, &cend)); /* used by ghosted vscale */ 24448a46eb9SPierre Jolivet if (vscale) PetscCall(VecGetArray(vscale, &vscale_array)); 2458bc97078SHong Zhang nz = 0; 246e2c857f8SHong Zhang 247b3e1f37bSHong Zhang if (coloring->bcols > 1) { /* use blocked insertion of Jentry */ 248b3e1f37bSHong Zhang PetscInt i, m = J->rmap->n, nbcols, bcols = coloring->bcols; 249e2c857f8SHong Zhang PetscScalar *dy = coloring->dy, *dy_k; 250e2c857f8SHong Zhang 251e2c857f8SHong Zhang nbcols = 0; 252e2c857f8SHong Zhang for (k = 0; k < ncolors; k += bcols) { 253e2c857f8SHong Zhang /* 254e2c857f8SHong Zhang (3-1) Loop over each column associated with color 255e2c857f8SHong Zhang adding the perturbation to the vector w3 = x1 + dx. 256e2c857f8SHong Zhang */ 257e2c857f8SHong Zhang 258e2c857f8SHong Zhang dy_k = dy; 259e2c857f8SHong Zhang if (k + bcols > ncolors) bcols = ncolors - k; 260e2c857f8SHong Zhang for (i = 0; i < bcols; i++) { 261ceef8a49SBarry Smith coloring->currentcolor = k + i; 262e2c857f8SHong Zhang 2639566063dSJacob Faibussowitsch PetscCall(VecCopy(x1, w3)); 2649566063dSJacob Faibussowitsch PetscCall(VecGetArray(w3, &w3_array)); 265e2c857f8SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */ 266e2c857f8SHong Zhang if (coloring->htype[0] == 'w') { 267e2c857f8SHong Zhang for (l = 0; l < ncolumns[k + i]; l++) { 268e2c857f8SHong Zhang col = coloring->columns[k + i][l]; /* local column (in global index!) of the matrix we are probing for */ 269e2c857f8SHong Zhang w3_array[col] += 1.0 / dx; 270e2c857f8SHong Zhang } 271e2c857f8SHong Zhang } else { /* htype == 'ds' */ 272e2c857f8SHong Zhang vscale_array -= cstart; /* shift pointer so global index can be used */ 273e2c857f8SHong Zhang for (l = 0; l < ncolumns[k + i]; l++) { 274e2c857f8SHong Zhang col = coloring->columns[k + i][l]; /* local column (in global index!) of the matrix we are probing for */ 275e2c857f8SHong Zhang w3_array[col] += 1.0 / vscale_array[col]; 276e2c857f8SHong Zhang } 277e2c857f8SHong Zhang vscale_array += cstart; 278e2c857f8SHong Zhang } 279e2c857f8SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array += cstart; 2809566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(w3, &w3_array)); 281e2c857f8SHong Zhang 282e2c857f8SHong Zhang /* 283e2c857f8SHong Zhang (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 284e2c857f8SHong Zhang w2 = F(x1 + dx) - F(x1) 285e2c857f8SHong Zhang */ 2869566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MAT_FDColoringFunction, 0, 0, 0, 0)); 2879566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(w2, dy_k)); /* place w2 to the array dy_i */ 2889566063dSJacob Faibussowitsch PetscCall((*f)(sctx, w3, w2, fctx)); 2899566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MAT_FDColoringFunction, 0, 0, 0, 0)); 2909566063dSJacob Faibussowitsch PetscCall(VecAXPY(w2, -1.0, w1)); 2919566063dSJacob Faibussowitsch PetscCall(VecResetArray(w2)); 292e2c857f8SHong Zhang dy_k += m; /* points to dy+i*nxloc */ 293e2c857f8SHong Zhang } 294e2c857f8SHong Zhang 295e2c857f8SHong Zhang /* 296e2c857f8SHong Zhang (3-3) Loop over block rows of vector, putting results into Jacobian matrix 297e2c857f8SHong Zhang */ 298d880da65SHong Zhang nrows_k = nrows[nbcols++]; 299d880da65SHong Zhang 300e2c857f8SHong Zhang if (coloring->htype[0] == 'w') { 301e2c857f8SHong Zhang for (l = 0; l < nrows_k; l++) { 3020df34763SHong Zhang row = Jentry2[nz].row; /* local row index */ 3035059b303SJunchao 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 3045059b303SJunchao Zhang another way, and it seems work. See https://lists.mcs.anl.gov/pipermail/petsc-users/2021-December/045158.html 3055059b303SJunchao Zhang */ 3065059b303SJunchao Zhang #if defined(PETSC_USE_COMPLEX) 3075059b303SJunchao Zhang PetscScalar *tmp = Jentry2[nz].valaddr; 3085059b303SJunchao Zhang *tmp = dy[row] * dx; 3095059b303SJunchao Zhang #else 310*f4f49eeaSPierre Jolivet *Jentry2[nz].valaddr = dy[row] * dx; 3115059b303SJunchao Zhang #endif 3125059b303SJunchao Zhang nz++; 313e2c857f8SHong Zhang } 314e2c857f8SHong Zhang } else { /* htype == 'ds' */ 315e2c857f8SHong Zhang for (l = 0; l < nrows_k; l++) { 316e2c857f8SHong Zhang row = Jentry[nz].row; /* local row index */ 3175059b303SJunchao Zhang #if defined(PETSC_USE_COMPLEX) /* See https://lists.mcs.anl.gov/pipermail/petsc-users/2021-December/045158.html */ 3185059b303SJunchao Zhang PetscScalar *tmp = Jentry[nz].valaddr; 3195059b303SJunchao Zhang *tmp = dy[row] * vscale_array[Jentry[nz].col]; 3205059b303SJunchao Zhang #else 321*f4f49eeaSPierre Jolivet *Jentry[nz].valaddr = dy[row] * vscale_array[Jentry[nz].col]; 3225059b303SJunchao Zhang #endif 323e2c857f8SHong Zhang nz++; 324e2c857f8SHong Zhang } 325e2c857f8SHong Zhang } 326e2c857f8SHong Zhang } 327b3e1f37bSHong Zhang } else { /* bcols == 1 */ 3288bc97078SHong Zhang for (k = 0; k < ncolors; k++) { 3298bc97078SHong Zhang coloring->currentcolor = k; 3309e917edbSHong Zhang 331fcd7ac73SHong Zhang /* 3328bc97078SHong Zhang (3-1) Loop over each column associated with color 333f6af9589SHong Zhang adding the perturbation to the vector w3 = x1 + dx. 334fcd7ac73SHong Zhang */ 3359566063dSJacob Faibussowitsch PetscCall(VecCopy(x1, w3)); 3369566063dSJacob Faibussowitsch PetscCall(VecGetArray(w3, &w3_array)); 337a2f2d239SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */ 338b039d6c4SHong Zhang if (coloring->htype[0] == 'w') { 3398bc97078SHong Zhang for (l = 0; l < ncolumns[k]; l++) { 340f6af9589SHong Zhang col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 3419e917edbSHong Zhang w3_array[col] += 1.0 / dx; 3429e917edbSHong Zhang } 343b039d6c4SHong Zhang } else { /* htype == 'ds' */ 344a2f2d239SHong Zhang vscale_array -= cstart; /* shift pointer so global index can be used */ 345b039d6c4SHong Zhang for (l = 0; l < ncolumns[k]; l++) { 346b039d6c4SHong Zhang col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 347a2f2d239SHong Zhang w3_array[col] += 1.0 / vscale_array[col]; 34870e7395fSHong Zhang } 349a2f2d239SHong Zhang vscale_array += cstart; 350fcd7ac73SHong Zhang } 351a2f2d239SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array += cstart; 3529566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(w3, &w3_array)); 3539e917edbSHong Zhang 354fcd7ac73SHong Zhang /* 3558bc97078SHong Zhang (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 356fcd7ac73SHong Zhang w2 = F(x1 + dx) - F(x1) 357fcd7ac73SHong Zhang */ 3589566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MAT_FDColoringFunction, 0, 0, 0, 0)); 3599566063dSJacob Faibussowitsch PetscCall((*f)(sctx, w3, w2, fctx)); 3609566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MAT_FDColoringFunction, 0, 0, 0, 0)); 3619566063dSJacob Faibussowitsch PetscCall(VecAXPY(w2, -1.0, w1)); 3629e917edbSHong Zhang 3638bc97078SHong Zhang /* 3648bc97078SHong Zhang (3-3) Loop over rows of vector, putting results into Jacobian matrix 3658bc97078SHong Zhang */ 366c53567a0SHong Zhang nrows_k = nrows[k]; 3679566063dSJacob Faibussowitsch PetscCall(VecGetArray(w2, &y)); 368b039d6c4SHong Zhang if (coloring->htype[0] == 'w') { 369c53567a0SHong Zhang for (l = 0; l < nrows_k; l++) { 3700df34763SHong Zhang row = Jentry2[nz].row; /* local row index */ 3715059b303SJunchao Zhang #if defined(PETSC_USE_COMPLEX) /* See https://lists.mcs.anl.gov/pipermail/petsc-users/2021-December/045158.html */ 3725059b303SJunchao Zhang PetscScalar *tmp = Jentry2[nz].valaddr; 3735059b303SJunchao Zhang *tmp = y[row] * dx; 3745059b303SJunchao Zhang #else 375*f4f49eeaSPierre Jolivet *Jentry2[nz].valaddr = y[row] * dx; 3765059b303SJunchao Zhang #endif 3775059b303SJunchao Zhang nz++; 3789e917edbSHong Zhang } 379b039d6c4SHong Zhang } else { /* htype == 'ds' */ 380c53567a0SHong Zhang for (l = 0; l < nrows_k; l++) { 381b039d6c4SHong Zhang row = Jentry[nz].row; /* local row index */ 3825059b303SJunchao Zhang #if defined(PETSC_USE_COMPLEX) /* See https://lists.mcs.anl.gov/pipermail/petsc-users/2021-December/045158.html */ 3835059b303SJunchao Zhang PetscScalar *tmp = Jentry[nz].valaddr; 3845059b303SJunchao Zhang *tmp = y[row] * vscale_array[Jentry[nz].col]; 3855059b303SJunchao Zhang #else 386*f4f49eeaSPierre Jolivet *Jentry[nz].valaddr = y[row] * vscale_array[Jentry[nz].col]; 3875059b303SJunchao Zhang #endif 388573f477fSHong Zhang nz++; 389fcd7ac73SHong Zhang } 390b039d6c4SHong Zhang } 3919566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(w2, &y)); 3928bc97078SHong Zhang } 393b3e1f37bSHong Zhang } 394b3e1f37bSHong Zhang 3959566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 3969566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 39748a46eb9SPierre Jolivet if (vscale) PetscCall(VecRestoreArray(vscale, &vscale_array)); 398fcd7ac73SHong Zhang coloring->currentcolor = -1; 3999566063dSJacob Faibussowitsch if (!alreadyboundtocpu) PetscCall(VecBindToCPU(x1, PETSC_FALSE)); 4003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 401fcd7ac73SHong Zhang } 402fcd7ac73SHong Zhang 403d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetUp_MPIXAIJ(Mat mat, ISColoring iscoloring, MatFDColoring c) 404d71ae5a4SJacob Faibussowitsch { 405fcd7ac73SHong Zhang PetscMPIInt size, *ncolsonproc, *disp, nn; 406e3225e9dSHong Zhang PetscInt i, n, nrows, nrows_i, j, k, m, ncols, col, *rowhit, cstart, cend, colb; 40772c15787SHong Zhang const PetscInt *is, *A_ci, *A_cj, *B_ci, *B_cj, *row = NULL, *ltog = NULL; 408071fcb05SBarry Smith PetscInt nis = iscoloring->n, nctot, *cols, tmp = 0; 409fcd7ac73SHong Zhang ISLocalToGlobalMapping map = mat->cmap->mapping; 410a8971b87SHong Zhang PetscInt ctype = c->ctype, *spidxA, *spidxB, nz, bs, bs2, spidx; 4110d1c53f1SHong Zhang Mat A, B; 412e3225e9dSHong Zhang PetscScalar *A_val, *B_val, **valaddrhit; 413573f477fSHong Zhang MatEntry *Jentry; 4140df34763SHong Zhang MatEntry2 *Jentry2; 415d4002b98SHong Zhang PetscBool isBAIJ, isSELL; 416531e53bdSHong Zhang PetscInt bcols = c->bcols; 4170d1c53f1SHong Zhang #if defined(PETSC_USE_CTABLE) 418eec179cfSJacob Faibussowitsch PetscHMapI colmap = NULL; 4190d1c53f1SHong Zhang #else 4200d1c53f1SHong Zhang PetscInt *colmap = NULL; /* local col number of off-diag col */ 4210d1c53f1SHong Zhang #endif 422fcd7ac73SHong Zhang 423fcd7ac73SHong Zhang PetscFunctionBegin; 4245bdb020cSBarry Smith if (ctype == IS_COLORING_LOCAL) { 42528b400f6SJacob Faibussowitsch PetscCheck(map, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_INCOMP, "When using ghosted differencing matrix must have local to global mapping provided with MatSetLocalToGlobalMapping"); 4269566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(map, <og)); 42772c15787SHong Zhang } 428fcd7ac73SHong Zhang 4299566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(mat, &bs)); 4309566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)mat, MATMPIBAIJ, &isBAIJ)); 4319566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATMPISELL, &isSELL)); 432e3225e9dSHong Zhang if (isBAIJ) { 433195687c5SHong Zhang Mat_MPIBAIJ *baij = (Mat_MPIBAIJ *)mat->data; 4346c145f34SHong Zhang Mat_SeqBAIJ *spA, *spB; 4359371c9d4SSatish Balay A = baij->A; 4369371c9d4SSatish Balay spA = (Mat_SeqBAIJ *)A->data; 4379371c9d4SSatish Balay A_val = spA->a; 4389371c9d4SSatish Balay B = baij->B; 4399371c9d4SSatish Balay spB = (Mat_SeqBAIJ *)B->data; 4409371c9d4SSatish Balay B_val = spB->a; 4410d1c53f1SHong Zhang nz = spA->nz + spB->nz; /* total nonzero entries of mat */ 44248a46eb9SPierre Jolivet if (!baij->colmap) PetscCall(MatCreateColmap_MPIBAIJ_Private(mat)); 443097e26fcSJed Brown colmap = baij->colmap; 4449566063dSJacob Faibussowitsch PetscCall(MatGetColumnIJ_SeqBAIJ_Color(A, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &A_ci, &A_cj, &spidxA, NULL)); 4459566063dSJacob Faibussowitsch PetscCall(MatGetColumnIJ_SeqBAIJ_Color(B, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &B_ci, &B_cj, &spidxB, NULL)); 446195687c5SHong Zhang 447195687c5SHong Zhang if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') { /* create vscale for storing dx */ 448195687c5SHong Zhang PetscInt *garray; 4499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(B->cmap->n, &garray)); 450195687c5SHong Zhang for (i = 0; i < baij->B->cmap->n / bs; i++) { 451ad540459SPierre Jolivet for (j = 0; j < bs; j++) garray[i * bs + j] = bs * baij->garray[i] + j; 452195687c5SHong Zhang } 4539566063dSJacob Faibussowitsch PetscCall(VecCreateGhost(PetscObjectComm((PetscObject)mat), mat->cmap->n, PETSC_DETERMINE, B->cmap->n, garray, &c->vscale)); 4549566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(c->vscale, PETSC_TRUE)); 4559566063dSJacob Faibussowitsch PetscCall(PetscFree(garray)); 456195687c5SHong Zhang } 457d4002b98SHong Zhang } else if (isSELL) { 458d4002b98SHong Zhang Mat_MPISELL *sell = (Mat_MPISELL *)mat->data; 459d4002b98SHong Zhang Mat_SeqSELL *spA, *spB; 4609371c9d4SSatish Balay A = sell->A; 4619371c9d4SSatish Balay spA = (Mat_SeqSELL *)A->data; 4629371c9d4SSatish Balay A_val = spA->val; 4639371c9d4SSatish Balay B = sell->B; 4649371c9d4SSatish Balay spB = (Mat_SeqSELL *)B->data; 4659371c9d4SSatish Balay B_val = spB->val; 466f4b713efSHong Zhang nz = spA->nz + spB->nz; /* total nonzero entries of mat */ 467d4002b98SHong Zhang if (!sell->colmap) { 468f4b713efSHong Zhang /* Allow access to data structures of local part of matrix 469f4b713efSHong Zhang - creates aij->colmap which maps global column number to local number in part B */ 4709566063dSJacob Faibussowitsch PetscCall(MatCreateColmap_MPISELL_Private(mat)); 471f4b713efSHong Zhang } 472d4002b98SHong Zhang colmap = sell->colmap; 4739566063dSJacob Faibussowitsch PetscCall(MatGetColumnIJ_SeqSELL_Color(A, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &A_ci, &A_cj, &spidxA, NULL)); 4749566063dSJacob Faibussowitsch PetscCall(MatGetColumnIJ_SeqSELL_Color(B, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &B_ci, &B_cj, &spidxB, NULL)); 475f4b713efSHong Zhang 476f4b713efSHong Zhang bs = 1; /* only bs=1 is supported for non MPIBAIJ matrix */ 477f4b713efSHong Zhang 478f4b713efSHong Zhang if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') { /* create vscale for storing dx */ 4799566063dSJacob Faibussowitsch PetscCall(VecCreateGhost(PetscObjectComm((PetscObject)mat), mat->cmap->n, PETSC_DETERMINE, B->cmap->n, sell->garray, &c->vscale)); 4809566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(c->vscale, PETSC_TRUE)); 481f4b713efSHong Zhang } 482e3225e9dSHong Zhang } else { 483e3225e9dSHong Zhang Mat_MPIAIJ *aij = (Mat_MPIAIJ *)mat->data; 484e3225e9dSHong Zhang Mat_SeqAIJ *spA, *spB; 4859371c9d4SSatish Balay A = aij->A; 4869371c9d4SSatish Balay spA = (Mat_SeqAIJ *)A->data; 4879371c9d4SSatish Balay A_val = spA->a; 4889371c9d4SSatish Balay B = aij->B; 4899371c9d4SSatish Balay spB = (Mat_SeqAIJ *)B->data; 4909371c9d4SSatish Balay B_val = spB->a; 491e3225e9dSHong Zhang nz = spA->nz + spB->nz; /* total nonzero entries of mat */ 492e3225e9dSHong Zhang if (!aij->colmap) { 493e3225e9dSHong Zhang /* Allow access to data structures of local part of matrix 494e3225e9dSHong Zhang - creates aij->colmap which maps global column number to local number in part B */ 4959566063dSJacob Faibussowitsch PetscCall(MatCreateColmap_MPIAIJ_Private(mat)); 496e3225e9dSHong Zhang } 497097e26fcSJed Brown colmap = aij->colmap; 4989566063dSJacob Faibussowitsch PetscCall(MatGetColumnIJ_SeqAIJ_Color(A, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &A_ci, &A_cj, &spidxA, NULL)); 4999566063dSJacob Faibussowitsch PetscCall(MatGetColumnIJ_SeqAIJ_Color(B, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &B_ci, &B_cj, &spidxB, NULL)); 500e3225e9dSHong Zhang 501e3225e9dSHong Zhang bs = 1; /* only bs=1 is supported for non MPIBAIJ matrix */ 502195687c5SHong Zhang 503195687c5SHong Zhang if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') { /* create vscale for storing dx */ 5049566063dSJacob Faibussowitsch PetscCall(VecCreateGhost(PetscObjectComm((PetscObject)mat), mat->cmap->n, PETSC_DETERMINE, B->cmap->n, aij->garray, &c->vscale)); 5059566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(c->vscale, PETSC_TRUE)); 506195687c5SHong Zhang } 5070d1c53f1SHong Zhang } 5080d1c53f1SHong Zhang 5090d1c53f1SHong Zhang m = mat->rmap->n / bs; 5100d1c53f1SHong Zhang cstart = mat->cmap->rstart / bs; 5110d1c53f1SHong Zhang cend = mat->cmap->rend / bs; 512fcd7ac73SHong Zhang 5139566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nis, &c->ncolumns, nis, &c->columns)); 5149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nis, &c->nrows)); 515fcd7ac73SHong Zhang 5160df34763SHong Zhang if (c->htype[0] == 'd') { 5179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nz, &Jentry)); 5188bc97078SHong Zhang c->matentry = Jentry; 5190df34763SHong Zhang } else if (c->htype[0] == 'w') { 5209566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nz, &Jentry2)); 5210df34763SHong Zhang c->matentry2 = Jentry2; 5220df34763SHong Zhang } else SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "htype is not supported"); 523a774a6f1SHong Zhang 5249566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(m + 1, &rowhit, m + 1, &valaddrhit)); 525fcd7ac73SHong Zhang nz = 0; 5269566063dSJacob Faibussowitsch PetscCall(ISColoringGetIS(iscoloring, PETSC_OWN_POINTER, PETSC_IGNORE, &c->isa)); 527071fcb05SBarry Smith 528071fcb05SBarry Smith if (ctype == IS_COLORING_GLOBAL) { 5299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &size)); 5309566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size, &ncolsonproc, size, &disp)); 531071fcb05SBarry Smith } 532071fcb05SBarry Smith 533d3825b63SHong Zhang for (i = 0; i < nis; i++) { /* for each local color */ 5349566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(c->isa[i], &n)); 5359566063dSJacob Faibussowitsch PetscCall(ISGetIndices(c->isa[i], &is)); 536fcd7ac73SHong Zhang 537fcd7ac73SHong Zhang c->ncolumns[i] = n; /* local number of columns of this color on this process */ 538071fcb05SBarry Smith c->columns[i] = (PetscInt *)is; 539fcd7ac73SHong Zhang 540fcd7ac73SHong Zhang if (ctype == IS_COLORING_GLOBAL) { 541d3825b63SHong Zhang /* Determine nctot, the total (parallel) number of columns of this color */ 542d3825b63SHong Zhang /* ncolsonproc[j]: local ncolumns on proc[j] of this color */ 5439566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(n, &nn)); 5449566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&nn, 1, MPI_INT, ncolsonproc, 1, MPI_INT, PetscObjectComm((PetscObject)mat))); 5459371c9d4SSatish Balay nctot = 0; 5469371c9d4SSatish Balay for (j = 0; j < size; j++) nctot += ncolsonproc[j]; 54748a46eb9SPierre Jolivet if (!nctot) PetscCall(PetscInfo(mat, "Coloring of matrix has some unneeded colors with no corresponding rows\n")); 548fcd7ac73SHong Zhang 549fcd7ac73SHong Zhang disp[0] = 0; 550ad540459SPierre Jolivet for (j = 1; j < size; j++) disp[j] = disp[j - 1] + ncolsonproc[j - 1]; 551fcd7ac73SHong Zhang 552d3825b63SHong Zhang /* Get cols, the complete list of columns for this color on each process */ 5539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nctot + 1, &cols)); 5549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv((void *)is, n, MPIU_INT, cols, ncolsonproc, disp, MPIU_INT, PetscObjectComm((PetscObject)mat))); 5555bdb020cSBarry Smith } else if (ctype == IS_COLORING_LOCAL) { 556fcd7ac73SHong Zhang /* Determine local number of columns of this color on this process, including ghost points */ 557fcd7ac73SHong Zhang nctot = n; 558071fcb05SBarry Smith cols = (PetscInt *)is; 559fcd7ac73SHong Zhang } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not provided for this MatFDColoring type"); 560fcd7ac73SHong Zhang 5611b97d346SHong Zhang /* Mark all rows affect by these columns */ 5629566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(rowhit, m)); 5630d1c53f1SHong Zhang bs2 = bs * bs; 5644b2e90caSHong Zhang nrows_i = 0; 5651b97d346SHong Zhang for (j = 0; j < nctot; j++) { /* loop over columns*/ 5665bdb020cSBarry Smith if (ctype == IS_COLORING_LOCAL) { 567fcd7ac73SHong Zhang col = ltog[cols[j]]; 568fcd7ac73SHong Zhang } else { 569fcd7ac73SHong Zhang col = cols[j]; 570fcd7ac73SHong Zhang } 571e3225e9dSHong Zhang if (col >= cstart && col < cend) { /* column is in A, diagonal block of mat */ 572071fcb05SBarry Smith tmp = A_ci[col - cstart]; 573071fcb05SBarry Smith row = A_cj + tmp; 574071fcb05SBarry Smith nrows = A_ci[col - cstart + 1] - tmp; 5754b2e90caSHong Zhang nrows_i += nrows; 576071fcb05SBarry Smith 577fcd7ac73SHong Zhang /* loop over columns of A marking them in rowhit */ 578d3825b63SHong Zhang for (k = 0; k < nrows; k++) { 579d3825b63SHong Zhang /* set valaddrhit for part A */ 580071fcb05SBarry Smith spidx = bs2 * spidxA[tmp + k]; 581e3225e9dSHong Zhang valaddrhit[*row] = &A_val[spidx]; 582a774a6f1SHong Zhang rowhit[*row++] = col - cstart + 1; /* local column index */ 583fcd7ac73SHong Zhang } 584e3225e9dSHong Zhang } else { /* column is in B, off-diagonal block of mat */ 585fcd7ac73SHong Zhang #if defined(PETSC_USE_CTABLE) 586eec179cfSJacob Faibussowitsch PetscCall(PetscHMapIGetWithDefault(colmap, col + 1, 0, &colb)); 587fcd7ac73SHong Zhang colb--; 588fcd7ac73SHong Zhang #else 5890d1c53f1SHong Zhang colb = colmap[col] - 1; /* local column index */ 590fcd7ac73SHong Zhang #endif 591fcd7ac73SHong Zhang if (colb == -1) { 592d3825b63SHong Zhang nrows = 0; 593fcd7ac73SHong Zhang } else { 5940d1c53f1SHong Zhang colb = colb / bs; 595071fcb05SBarry Smith tmp = B_ci[colb]; 596071fcb05SBarry Smith row = B_cj + tmp; 597071fcb05SBarry Smith nrows = B_ci[colb + 1] - tmp; 598fcd7ac73SHong Zhang } 5994b2e90caSHong Zhang nrows_i += nrows; 600fcd7ac73SHong Zhang /* loop over columns of B marking them in rowhit */ 601d3825b63SHong Zhang for (k = 0; k < nrows; k++) { 602d3825b63SHong Zhang /* set valaddrhit for part B */ 603071fcb05SBarry Smith spidx = bs2 * spidxB[tmp + k]; 604e3225e9dSHong Zhang valaddrhit[*row] = &B_val[spidx]; 60570e7395fSHong Zhang rowhit[*row++] = colb + 1 + cend - cstart; /* local column index */ 606fcd7ac73SHong Zhang } 607fcd7ac73SHong Zhang } 608fcd7ac73SHong Zhang } 6094b2e90caSHong Zhang c->nrows[i] = nrows_i; 6108bc97078SHong Zhang 6110df34763SHong Zhang if (c->htype[0] == 'd') { 612d3825b63SHong Zhang for (j = 0; j < m; j++) { 613fcd7ac73SHong Zhang if (rowhit[j]) { 614573f477fSHong Zhang Jentry[nz].row = j; /* local row index */ 615573f477fSHong Zhang Jentry[nz].col = rowhit[j] - 1; /* local column index */ 616573f477fSHong Zhang Jentry[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */ 617573f477fSHong Zhang nz++; 618fcd7ac73SHong Zhang } 619fcd7ac73SHong Zhang } 6200df34763SHong Zhang } else { /* c->htype == 'wp' */ 6210df34763SHong Zhang for (j = 0; j < m; j++) { 6220df34763SHong Zhang if (rowhit[j]) { 6230df34763SHong Zhang Jentry2[nz].row = j; /* local row index */ 6240df34763SHong Zhang Jentry2[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */ 6250df34763SHong Zhang nz++; 6260df34763SHong Zhang } 6270df34763SHong Zhang } 6280df34763SHong Zhang } 62948a46eb9SPierre Jolivet if (ctype == IS_COLORING_GLOBAL) PetscCall(PetscFree(cols)); 630fcd7ac73SHong Zhang } 63148a46eb9SPierre Jolivet if (ctype == IS_COLORING_GLOBAL) PetscCall(PetscFree2(ncolsonproc, disp)); 632fcd7ac73SHong Zhang 633a8971b87SHong Zhang if (bcols > 1) { /* reorder Jentry for faster MatFDColoringApply() */ 6349566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetUpBlocked_AIJ_Private(mat, c, nz)); 635531e53bdSHong Zhang } 636531e53bdSHong Zhang 6370d1c53f1SHong Zhang if (isBAIJ) { 6389566063dSJacob Faibussowitsch PetscCall(MatRestoreColumnIJ_SeqBAIJ_Color(A, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &A_ci, &A_cj, &spidxA, NULL)); 6399566063dSJacob Faibussowitsch PetscCall(MatRestoreColumnIJ_SeqBAIJ_Color(B, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &B_ci, &B_cj, &spidxB, NULL)); 6409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs * mat->rmap->n, &c->dy)); 641d4002b98SHong Zhang } else if (isSELL) { 6429566063dSJacob Faibussowitsch PetscCall(MatRestoreColumnIJ_SeqSELL_Color(A, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &A_ci, &A_cj, &spidxA, NULL)); 6439566063dSJacob Faibussowitsch PetscCall(MatRestoreColumnIJ_SeqSELL_Color(B, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &B_ci, &B_cj, &spidxB, NULL)); 6440d1c53f1SHong Zhang } else { 6459566063dSJacob Faibussowitsch PetscCall(MatRestoreColumnIJ_SeqAIJ_Color(A, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &A_ci, &A_cj, &spidxA, NULL)); 6469566063dSJacob Faibussowitsch PetscCall(MatRestoreColumnIJ_SeqAIJ_Color(B, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &B_ci, &B_cj, &spidxB, NULL)); 6470d1c53f1SHong Zhang } 6480d1c53f1SHong Zhang 6499566063dSJacob Faibussowitsch PetscCall(ISColoringRestoreIS(iscoloring, PETSC_OWN_POINTER, &c->isa)); 6509566063dSJacob Faibussowitsch PetscCall(PetscFree2(rowhit, valaddrhit)); 651a8971b87SHong Zhang 65248a46eb9SPierre Jolivet if (ctype == IS_COLORING_LOCAL) PetscCall(ISLocalToGlobalMappingRestoreIndices(map, <og)); 6539566063dSJacob Faibussowitsch PetscCall(PetscInfo(c, "ncolors %" PetscInt_FMT ", brows %" PetscInt_FMT " and bcols %" PetscInt_FMT " are used.\n", c->ncolors, c->brows, c->bcols)); 6543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 655fcd7ac73SHong Zhang } 65693dfae19SHong Zhang 657d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat mat, ISColoring iscoloring, MatFDColoring c) 658d71ae5a4SJacob Faibussowitsch { 659a8971b87SHong Zhang PetscInt bs, nis = iscoloring->n, m = mat->rmap->n; 660d4002b98SHong Zhang PetscBool isBAIJ, isSELL; 661531e53bdSHong Zhang 66293dfae19SHong Zhang PetscFunctionBegin; 663531e53bdSHong Zhang /* set default brows and bcols for speedup inserting the dense matrix into sparse Jacobian; 664531e53bdSHong Zhang bcols is chosen s.t. dy-array takes 50% of memory space as mat */ 6659566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(mat, &bs)); 6669566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)mat, MATMPIBAIJ, &isBAIJ)); 6679566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATMPISELL, &isSELL)); 668b8dfa574SBarry Smith if (isBAIJ || m == 0) { 669a8971b87SHong Zhang c->brows = m; 670a8971b87SHong Zhang c->bcols = 1; 671d4002b98SHong Zhang } else if (isSELL) { 672f4b713efSHong Zhang /* bcols is chosen s.t. dy-array takes 50% of local memory space as mat */ 673d4002b98SHong Zhang Mat_MPISELL *sell = (Mat_MPISELL *)mat->data; 674d4002b98SHong Zhang Mat_SeqSELL *spA, *spB; 675f4b713efSHong Zhang Mat A, B; 676f4b713efSHong Zhang PetscInt nz, brows, bcols; 677f4b713efSHong Zhang PetscReal mem; 678f4b713efSHong Zhang 679d4002b98SHong Zhang bs = 1; /* only bs=1 is supported for MPISELL matrix */ 680f4b713efSHong Zhang 6819371c9d4SSatish Balay A = sell->A; 6829371c9d4SSatish Balay spA = (Mat_SeqSELL *)A->data; 6839371c9d4SSatish Balay B = sell->B; 6849371c9d4SSatish Balay spB = (Mat_SeqSELL *)B->data; 685f4b713efSHong Zhang nz = spA->nz + spB->nz; /* total local nonzero entries of mat */ 686f4b713efSHong Zhang mem = nz * (sizeof(PetscScalar) + sizeof(PetscInt)) + 3 * m * sizeof(PetscInt); 687f4b713efSHong Zhang bcols = (PetscInt)(0.5 * mem / (m * sizeof(PetscScalar))); 688f4b713efSHong Zhang brows = 1000 / bcols; 689f4b713efSHong Zhang if (bcols > nis) bcols = nis; 690f4b713efSHong Zhang if (brows == 0 || brows > m) brows = m; 691f4b713efSHong Zhang c->brows = brows; 692f4b713efSHong Zhang c->bcols = bcols; 693531e53bdSHong Zhang } else { /* mpiaij matrix */ 694531e53bdSHong Zhang /* bcols is chosen s.t. dy-array takes 50% of local memory space as mat */ 695531e53bdSHong Zhang Mat_MPIAIJ *aij = (Mat_MPIAIJ *)mat->data; 696531e53bdSHong Zhang Mat_SeqAIJ *spA, *spB; 697531e53bdSHong Zhang Mat A, B; 698a8971b87SHong Zhang PetscInt nz, brows, bcols; 699531e53bdSHong Zhang PetscReal mem; 700531e53bdSHong Zhang 701531e53bdSHong Zhang bs = 1; /* only bs=1 is supported for MPIAIJ matrix */ 702531e53bdSHong Zhang 7039371c9d4SSatish Balay A = aij->A; 7049371c9d4SSatish Balay spA = (Mat_SeqAIJ *)A->data; 7059371c9d4SSatish Balay B = aij->B; 7069371c9d4SSatish Balay spB = (Mat_SeqAIJ *)B->data; 707531e53bdSHong Zhang nz = spA->nz + spB->nz; /* total local nonzero entries of mat */ 708531e53bdSHong Zhang mem = nz * (sizeof(PetscScalar) + sizeof(PetscInt)) + 3 * m * sizeof(PetscInt); 709531e53bdSHong Zhang bcols = (PetscInt)(0.5 * mem / (m * sizeof(PetscScalar))); 710531e53bdSHong Zhang brows = 1000 / bcols; 711531e53bdSHong Zhang if (bcols > nis) bcols = nis; 712531e53bdSHong Zhang if (brows == 0 || brows > m) brows = m; 713531e53bdSHong Zhang c->brows = brows; 714531e53bdSHong Zhang c->bcols = bcols; 715531e53bdSHong Zhang } 716a8971b87SHong Zhang 717a8971b87SHong Zhang c->M = mat->rmap->N / bs; /* set the global rows and columns and local rows */ 718a8971b87SHong Zhang c->N = mat->cmap->N / bs; 719a8971b87SHong Zhang c->m = mat->rmap->n / bs; 720a8971b87SHong Zhang c->rstart = mat->rmap->rstart / bs; 721a8971b87SHong Zhang c->ncolors = nis; 7223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 72393dfae19SHong Zhang } 724bdaf1daeSBarry Smith 725bdaf1daeSBarry Smith /*@C 726bdaf1daeSBarry Smith 72711a5261eSBarry Smith MatFDColoringSetValues - takes a matrix in compressed color format and enters the matrix into a PETSc `Mat` 728bdaf1daeSBarry Smith 729c3339decSBarry Smith Collective 730bdaf1daeSBarry Smith 731bdaf1daeSBarry Smith Input Parameters: 732bdaf1daeSBarry Smith + J - the sparse matrix 73311a5261eSBarry Smith . coloring - created with `MatFDColoringCreate()` and a local coloring 734bdaf1daeSBarry Smith - y - column major storage of matrix values with one color of values per column, the number of rows of y should match 7352ef1f0ffSBarry Smith the number of local rows of `J` and the number of columns is the number of colors. 736bdaf1daeSBarry Smith 737bdaf1daeSBarry Smith Level: intermediate 738bdaf1daeSBarry Smith 7392ef1f0ffSBarry Smith Notes: 7402ef1f0ffSBarry Smith The matrix in compressed color format may come from an automatic differentiation code 741bdaf1daeSBarry Smith 74211a5261eSBarry Smith The code will be slightly faster if `MatFDColoringSetBlockSize`(coloring,`PETSC_DEFAULT`,nc); is called immediately after creating the coloring 743bdaf1daeSBarry Smith 7441cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatFDColoringCreate()`, `ISColoring`, `ISColoringCreate()`, `ISColoringSetType()`, `IS_COLORING_LOCAL`, `MatFDColoringSetBlockSize()` 745bdaf1daeSBarry Smith @*/ 746d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetValues(Mat J, MatFDColoring coloring, const PetscScalar *y) 747d71ae5a4SJacob Faibussowitsch { 748bdaf1daeSBarry Smith MatEntry2 *Jentry2; 749bdaf1daeSBarry Smith PetscInt row, i, nrows_k, l, ncolors, nz = 0, bcols, nbcols = 0; 750bdaf1daeSBarry Smith const PetscInt *nrows; 751bdaf1daeSBarry Smith PetscBool eq; 752bdaf1daeSBarry Smith 753bdaf1daeSBarry Smith PetscFunctionBegin; 754bdaf1daeSBarry Smith PetscValidHeaderSpecific(J, MAT_CLASSID, 1); 755bdaf1daeSBarry Smith PetscValidHeaderSpecific(coloring, MAT_FDCOLORING_CLASSID, 2); 7569566063dSJacob Faibussowitsch PetscCall(PetscObjectCompareId((PetscObject)J, coloring->matid, &eq)); 75728b400f6SJacob Faibussowitsch PetscCheck(eq, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_WRONG, "Matrix used with MatFDColoringSetValues() must be that used with MatFDColoringCreate()"); 758bdaf1daeSBarry Smith Jentry2 = coloring->matentry2; 759bdaf1daeSBarry Smith nrows = coloring->nrows; 760bdaf1daeSBarry Smith ncolors = coloring->ncolors; 761bdaf1daeSBarry Smith bcols = coloring->bcols; 762bdaf1daeSBarry Smith 763bdaf1daeSBarry Smith for (i = 0; i < ncolors; i += bcols) { 764bdaf1daeSBarry Smith nrows_k = nrows[nbcols++]; 765bdaf1daeSBarry Smith for (l = 0; l < nrows_k; l++) { 766bdaf1daeSBarry Smith row = Jentry2[nz].row; /* local row index */ 767*f4f49eeaSPierre Jolivet *Jentry2[nz++].valaddr = y[row]; 768bdaf1daeSBarry Smith } 769bdaf1daeSBarry Smith y += bcols * coloring->m; 770bdaf1daeSBarry Smith } 7719566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 7729566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 7733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 774bdaf1daeSBarry Smith } 775