12e68589bSMark F. Adams /* 2b817416eSBarry Smith GAMG geometric-algebric multigrid PC - Mark Adams 2011 32e68589bSMark F. Adams */ 42e68589bSMark F. Adams 52e68589bSMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 62e68589bSMark F. Adams #include <petscblaslapack.h> 7539c167fSBarry Smith #include <petscdm.h> 873f7197eSJed Brown #include <petsc/private/kspimpl.h> 92e68589bSMark F. Adams 102e68589bSMark F. Adams typedef struct { 11c8b0795cSMark F. Adams PetscInt nsmooths; 12bae903cbSmarkadams4 PetscInt aggressive_coarsening_levels; // number of aggressive coarsening levels (square or MISk) 132d776b49SBarry Smith MatCoarsen crs; 142e68589bSMark F. Adams } PC_GAMG_AGG; 152e68589bSMark F. Adams 162e68589bSMark F. Adams /*@ 17f1580f4eSBarry Smith PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical) used for multigrid on all the levels 182e68589bSMark F. Adams 19c3339decSBarry Smith Logically Collective 202e68589bSMark F. Adams 212e68589bSMark F. Adams Input Parameters: 22*20f4b53cSBarry Smith + pc - the preconditioner context 23*20f4b53cSBarry Smith - n - the number of smooths 242e68589bSMark F. Adams 252e68589bSMark F. Adams Options Database Key: 26cab9ed1eSBarry Smith . -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation 272e68589bSMark F. Adams 282e68589bSMark F. Adams Level: intermediate 292e68589bSMark F. Adams 30f1580f4eSBarry Smith .seealso: `PCMG`, `PCGAMG` 312e68589bSMark F. Adams @*/ 32d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n) 33d71ae5a4SJacob Faibussowitsch { 342e68589bSMark F. Adams PetscFunctionBegin; 352e68589bSMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 36d5d25401SBarry Smith PetscValidLogicalCollectiveInt(pc, n, 2); 37cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetNSmooths_C", (PC, PetscInt), (pc, n)); 383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 392e68589bSMark F. Adams } 402e68589bSMark F. Adams 41d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetNSmooths_AGG(PC pc, PetscInt n) 42d71ae5a4SJacob Faibussowitsch { 432e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 442e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 45c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 462e68589bSMark F. Adams 472e68589bSMark F. Adams PetscFunctionBegin; 48c8b0795cSMark F. Adams pc_gamg_agg->nsmooths = n; 493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 50c8b0795cSMark F. Adams } 51c8b0795cSMark F. Adams 52c8b0795cSMark F. Adams /*@ 53f1580f4eSBarry Smith PCGAMGSetAggressiveLevels - Use aggressive coarsening on first n levels 54ef4ad70eSMark F. Adams 55c3339decSBarry Smith Logically Collective 56ef4ad70eSMark F. Adams 57ef4ad70eSMark F. Adams Input Parameters: 58cab9ed1eSBarry Smith + pc - the preconditioner context 59d5d25401SBarry Smith - n - 0, 1 or more 60ef4ad70eSMark F. Adams 61ef4ad70eSMark F. Adams Options Database Key: 62bae903cbSmarkadams4 . -pc_gamg_aggressive_coarsening <n,default = 1> - Number of levels to square the graph on before aggregating it 63af3c827dSMark Adams 64ef4ad70eSMark F. Adams Level: intermediate 65ef4ad70eSMark F. Adams 662d776b49SBarry Smith .seealso: `PCGAMG`, `PCGAMGSetThreshold()` 67ef4ad70eSMark F. Adams @*/ 68d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetAggressiveLevels(PC pc, PetscInt n) 69d71ae5a4SJacob Faibussowitsch { 70ef4ad70eSMark F. Adams PetscFunctionBegin; 71ef4ad70eSMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 72d5d25401SBarry Smith PetscValidLogicalCollectiveInt(pc, n, 2); 73bae903cbSmarkadams4 PetscTryMethod(pc, "PCGAMGSetAggressiveLevels_C", (PC, PetscInt), (pc, n)); 743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 75ef4ad70eSMark F. Adams } 76ef4ad70eSMark F. Adams 77d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetAggressiveLevels_AGG(PC pc, PetscInt n) 78d71ae5a4SJacob Faibussowitsch { 79ef4ad70eSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 80ef4ad70eSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 81ef4ad70eSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 82ef4ad70eSMark F. Adams 83ef4ad70eSMark F. Adams PetscFunctionBegin; 84bae903cbSmarkadams4 pc_gamg_agg->aggressive_coarsening_levels = n; 853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 86ef4ad70eSMark F. Adams } 87ef4ad70eSMark F. Adams 88d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_GAMG_AGG(PC pc, PetscOptionItems *PetscOptionsObject) 89d71ae5a4SJacob Faibussowitsch { 902e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 912e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 92c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 932e68589bSMark F. Adams 942e68589bSMark F. Adams PetscFunctionBegin; 95d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "GAMG-AGG options"); 962e68589bSMark F. Adams { 97bae903cbSmarkadams4 PetscBool flg; 989566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_gamg_agg_nsmooths", "smoothing steps for smoothed aggregation, usually 1", "PCGAMGSetNSmooths", pc_gamg_agg->nsmooths, &pc_gamg_agg->nsmooths, NULL)); 999371c9d4SSatish Balay PetscCall( 1009371c9d4SSatish Balay PetscOptionsInt("-pc_gamg_square_graph", "Number of aggressive coarsening (MIS-2) levels from finest (alias for -pc_gamg_aggressive_coarsening, deprecated)", "PCGAMGSetAggressiveLevels", pc_gamg_agg->aggressive_coarsening_levels, &pc_gamg_agg->aggressive_coarsening_levels, &flg)); 101bae903cbSmarkadams4 if (!flg) { 102bae903cbSmarkadams4 PetscCall(PetscOptionsInt("-pc_gamg_aggressive_coarsening", "Number of aggressive coarsening (MIS-2) levels from finest", "PCGAMGSetAggressiveLevels", pc_gamg_agg->aggressive_coarsening_levels, &pc_gamg_agg->aggressive_coarsening_levels, NULL)); 103bae903cbSmarkadams4 } else { 104bae903cbSmarkadams4 PetscCall(PetscOptionsInt("-pc_gamg_aggressive_coarsening", "Number of aggressive coarsening (MIS-2) levels from finest", "PCGAMGSetAggressiveLevels", pc_gamg_agg->aggressive_coarsening_levels, &pc_gamg_agg->aggressive_coarsening_levels, &flg)); 105bae903cbSmarkadams4 if (flg) PetscCall(PetscInfo(pc, "Warning: both -pc_gamg_square_graph and -pc_gamg_aggressive_coarsening are used. -pc_gamg_square_graph is deprecated, Number of aggressive levels is %d\n", (int)pc_gamg_agg->aggressive_coarsening_levels)); 106bae903cbSmarkadams4 } 1072e68589bSMark F. Adams } 108d0609cedSBarry Smith PetscOptionsHeadEnd(); 1093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1102e68589bSMark F. Adams } 1112e68589bSMark F. Adams 112d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_GAMG_AGG(PC pc) 113d71ae5a4SJacob Faibussowitsch { 1142e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 1152e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1162e68589bSMark F. Adams 1172e68589bSMark F. Adams PetscFunctionBegin; 1189566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->subctx)); 1192e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", NULL)); 120bae903cbSmarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", NULL)); 121bae903cbSmarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL)); 1223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1232e68589bSMark F. Adams } 1242e68589bSMark F. Adams 1252e68589bSMark F. Adams /* 1262e68589bSMark F. Adams PCSetCoordinates_AGG 127*20f4b53cSBarry Smith 128*20f4b53cSBarry Smith Collective 1292e68589bSMark F. Adams 1302e68589bSMark F. Adams Input Parameter: 1312e68589bSMark F. Adams . pc - the preconditioner context 132a2f3521dSMark F. Adams . ndm - dimesion of data (used for dof/vertex for Stokes) 133302f38e8SMark F. Adams . a_nloc - number of vertices local 134302f38e8SMark F. Adams . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...} 1352e68589bSMark F. Adams */ 1361e6b0712SBarry Smith 137d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords) 138d71ae5a4SJacob Faibussowitsch { 1392e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 1402e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 14169344418SMark F. Adams PetscInt arrsz, kk, ii, jj, nloc, ndatarows, ndf; 142a2f3521dSMark F. Adams Mat mat = pc->pmat; 1432e68589bSMark F. Adams 1442e68589bSMark F. Adams PetscFunctionBegin; 145a2f3521dSMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 146a2f3521dSMark F. Adams PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 147302f38e8SMark F. Adams nloc = a_nloc; 1482e68589bSMark F. Adams 1492e68589bSMark F. Adams /* SA: null space vectors */ 1509566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(mat, &ndf)); /* this does not work for Stokes */ 15169344418SMark F. Adams if (coords && ndf == 1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */ 152a2f3521dSMark F. Adams else if (coords) { 15363a3b9bcSJacob Faibussowitsch PetscCheck(ndm <= ndf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "degrees of motion %" PetscInt_FMT " > block size %" PetscInt_FMT, ndm, ndf); 15469344418SMark F. Adams pc_gamg->data_cell_cols = (ndm == 2 ? 3 : 6); /* displacement elasticity */ 155ad540459SPierre Jolivet if (ndm != ndf) PetscCheck(pc_gamg->data_cell_cols == ndf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Don't know how to create null space for ndm=%" PetscInt_FMT ", ndf=%" PetscInt_FMT ". Use MatSetNearNullSpace().", ndm, ndf); 156806fa848SBarry Smith } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */ 15771959b99SBarry Smith pc_gamg->data_cell_rows = ndatarows = ndf; 15863a3b9bcSJacob Faibussowitsch PetscCheck(pc_gamg->data_cell_cols > 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "pc_gamg->data_cell_cols %" PetscInt_FMT " <= 0", pc_gamg->data_cell_cols); 159c8b0795cSMark F. Adams arrsz = nloc * pc_gamg->data_cell_rows * pc_gamg->data_cell_cols; 1602e68589bSMark F. Adams 1617f66b68fSMark Adams if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) { 1629566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 1639566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(arrsz + 1, &pc_gamg->data)); 1642e68589bSMark F. Adams } 1652e68589bSMark F. Adams /* copy data in - column oriented */ 1662e68589bSMark F. Adams for (kk = 0; kk < nloc; kk++) { 16769344418SMark F. Adams const PetscInt M = nloc * pc_gamg->data_cell_rows; /* stride into data */ 16869344418SMark F. Adams PetscReal *data = &pc_gamg->data[kk * ndatarows]; /* start of cell */ 169c8b0795cSMark F. Adams if (pc_gamg->data_cell_cols == 1) *data = 1.0; 1702e68589bSMark F. Adams else { 17169344418SMark F. Adams /* translational modes */ 1722fa5cd67SKarl Rupp for (ii = 0; ii < ndatarows; ii++) { 1732fa5cd67SKarl Rupp for (jj = 0; jj < ndatarows; jj++) { 17469344418SMark F. Adams if (ii == jj) data[ii * M + jj] = 1.0; 1752e68589bSMark F. Adams else data[ii * M + jj] = 0.0; 1762fa5cd67SKarl Rupp } 1772fa5cd67SKarl Rupp } 17869344418SMark F. Adams 17969344418SMark F. Adams /* rotational modes */ 1802e68589bSMark F. Adams if (coords) { 18169344418SMark F. Adams if (ndm == 2) { 1822e68589bSMark F. Adams data += 2 * M; 1832e68589bSMark F. Adams data[0] = -coords[2 * kk + 1]; 1842e68589bSMark F. Adams data[1] = coords[2 * kk]; 185806fa848SBarry Smith } else { 1862e68589bSMark F. Adams data += 3 * M; 1879371c9d4SSatish Balay data[0] = 0.0; 1889371c9d4SSatish Balay data[M + 0] = coords[3 * kk + 2]; 1899371c9d4SSatish Balay data[2 * M + 0] = -coords[3 * kk + 1]; 1909371c9d4SSatish Balay data[1] = -coords[3 * kk + 2]; 1919371c9d4SSatish Balay data[M + 1] = 0.0; 1929371c9d4SSatish Balay data[2 * M + 1] = coords[3 * kk]; 1939371c9d4SSatish Balay data[2] = coords[3 * kk + 1]; 1949371c9d4SSatish Balay data[M + 2] = -coords[3 * kk]; 1959371c9d4SSatish Balay data[2 * M + 2] = 0.0; 1962e68589bSMark F. Adams } 1972e68589bSMark F. Adams } 1982e68589bSMark F. Adams } 1992e68589bSMark F. Adams } 2002e68589bSMark F. Adams pc_gamg->data_sz = arrsz; 2013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2022e68589bSMark F. Adams } 2032e68589bSMark F. Adams 2042e68589bSMark F. Adams /* 205a2f3521dSMark F. Adams PCSetData_AGG - called if data is not set with PCSetCoordinates. 206a2f3521dSMark F. Adams Looks in Mat for near null space. 207a2f3521dSMark F. Adams Does not work for Stokes 2082e68589bSMark F. Adams 2092e68589bSMark F. Adams Input Parameter: 2102e68589bSMark F. Adams . pc - 211a2f3521dSMark F. Adams . a_A - matrix to get (near) null space out of. 2122e68589bSMark F. Adams */ 213d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A) 214d71ae5a4SJacob Faibussowitsch { 215b8cd405aSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 216b8cd405aSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 217b8cd405aSMark F. Adams MatNullSpace mnull; 21866f2ef4bSMark Adams 219b817416eSBarry Smith PetscFunctionBegin; 2209566063dSJacob Faibussowitsch PetscCall(MatGetNearNullSpace(a_A, &mnull)); 221b8cd405aSMark F. Adams if (!mnull) { 22266f2ef4bSMark Adams DM dm; 2239566063dSJacob Faibussowitsch PetscCall(PCGetDM(pc, &dm)); 22448a46eb9SPierre Jolivet if (!dm) PetscCall(MatGetDM(a_A, &dm)); 22566f2ef4bSMark Adams if (dm) { 22666f2ef4bSMark Adams PetscObject deformation; 227b0d30dd6SMatthew G. Knepley PetscInt Nf; 228b0d30dd6SMatthew G. Knepley 2299566063dSJacob Faibussowitsch PetscCall(DMGetNumFields(dm, &Nf)); 230b0d30dd6SMatthew G. Knepley if (Nf) { 2319566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, 0, NULL, &deformation)); 2329566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)deformation, "nearnullspace", (PetscObject *)&mnull)); 23348a46eb9SPierre Jolivet if (!mnull) PetscCall(PetscObjectQuery((PetscObject)deformation, "nullspace", (PetscObject *)&mnull)); 23466f2ef4bSMark Adams } 23566f2ef4bSMark Adams } 236b0d30dd6SMatthew G. Knepley } 23766f2ef4bSMark Adams 23866f2ef4bSMark Adams if (!mnull) { 239a2f3521dSMark F. Adams PetscInt bs, NN, MM; 2409566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(a_A, &bs)); 2419566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(a_A, &MM, &NN)); 24263a3b9bcSJacob Faibussowitsch PetscCheck(MM % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MM %" PetscInt_FMT " must be divisible by bs %" PetscInt_FMT, MM, bs); 2439566063dSJacob Faibussowitsch PetscCall(PCSetCoordinates_AGG(pc, bs, MM / bs, NULL)); 244806fa848SBarry Smith } else { 245b8cd405aSMark F. Adams PetscReal *nullvec; 246b8cd405aSMark F. Adams PetscBool has_const; 247b8cd405aSMark F. Adams PetscInt i, j, mlocal, nvec, bs; 248d5d25401SBarry Smith const Vec *vecs; 249d5d25401SBarry Smith const PetscScalar *v; 250b817416eSBarry Smith 2519566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(a_A, &mlocal, NULL)); 2529566063dSJacob Faibussowitsch PetscCall(MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs)); 253d19c4ebbSmarkadams4 for (i = 0; i < nvec; i++) { 254d19c4ebbSmarkadams4 PetscCall(VecGetLocalSize(vecs[i], &j)); 255d19c4ebbSmarkadams4 PetscCheck(j == mlocal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Attached null space vector size %" PetscInt_FMT " != matrix size %" PetscInt_FMT, j, mlocal); 256d19c4ebbSmarkadams4 } 257a0dea87dSMark Adams pc_gamg->data_sz = (nvec + !!has_const) * mlocal; 2589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((nvec + !!has_const) * mlocal, &nullvec)); 2599371c9d4SSatish Balay if (has_const) 2609371c9d4SSatish Balay for (i = 0; i < mlocal; i++) nullvec[i] = 1.0; 261b8cd405aSMark F. Adams for (i = 0; i < nvec; i++) { 2629566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(vecs[i], &v)); 263b8cd405aSMark F. Adams for (j = 0; j < mlocal; j++) nullvec[(i + !!has_const) * mlocal + j] = PetscRealPart(v[j]); 2649566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(vecs[i], &v)); 265b8cd405aSMark F. Adams } 266b8cd405aSMark F. Adams pc_gamg->data = nullvec; 267b8cd405aSMark F. Adams pc_gamg->data_cell_cols = (nvec + !!has_const); 2689566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(a_A, &bs)); 269b8cd405aSMark F. Adams pc_gamg->data_cell_rows = bs; 270b8cd405aSMark F. Adams } 2713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2722e68589bSMark F. Adams } 2732e68589bSMark F. Adams 2742e68589bSMark F. Adams /* 275bae903cbSmarkadams4 formProl0 - collect null space data for each aggregate, do QR, put R in coarse grid data and Q in P_0 2762e68589bSMark F. Adams 2772e68589bSMark F. Adams Input Parameter: 2782adfe9d3SBarry Smith . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices 2792adfe9d3SBarry Smith . bs - row block size 2800cbbd2e1SMark F. Adams . nSAvec - column bs of new P 2810cbbd2e1SMark F. Adams . my0crs - global index of start of locals 2822adfe9d3SBarry Smith . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec] 2832e68589bSMark F. Adams . data_in[data_stride*nSAvec] - local data on fine grid 2842e68589bSMark F. Adams . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist' 285bae903cbSmarkadams4 2862e68589bSMark F. Adams Output Parameter: 2872e68589bSMark F. Adams . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data 2882e68589bSMark F. Adams . a_Prol - prolongation operator 2892e68589bSMark F. Adams */ 290d71ae5a4SJacob Faibussowitsch static PetscErrorCode formProl0(PetscCoarsenData *agg_llists, PetscInt bs, PetscInt nSAvec, PetscInt my0crs, PetscInt data_stride, PetscReal data_in[], const PetscInt flid_fgid[], PetscReal **a_data_out, Mat a_Prol) 291d71ae5a4SJacob Faibussowitsch { 292bd026e97SJed Brown PetscInt Istart, my0, Iend, nloc, clid, flid = 0, aggID, kk, jj, ii, mm, nSelected, minsz, nghosts, out_data_stride; 2933b4367a7SBarry Smith MPI_Comm comm; 2942e68589bSMark F. Adams PetscReal *out_data; 295539c167fSBarry Smith PetscCDIntNd *pos; 2961943db53SBarry Smith PCGAMGHashTable fgid_flid; 2970cbbd2e1SMark F. Adams 2982e68589bSMark F. Adams PetscFunctionBegin; 2999566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)a_Prol, &comm)); 3009566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(a_Prol, &Istart, &Iend)); 3019371c9d4SSatish Balay nloc = (Iend - Istart) / bs; 3029371c9d4SSatish Balay my0 = Istart / bs; 30363a3b9bcSJacob Faibussowitsch PetscCheck((Iend - Istart) % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Iend %" PetscInt_FMT " - Istart %" PetscInt_FMT " must be divisible by bs %" PetscInt_FMT, Iend, Istart, bs); 3040cbbd2e1SMark F. Adams Iend /= bs; 3050cbbd2e1SMark F. Adams nghosts = data_stride / bs - nloc; 3062e68589bSMark F. Adams 3079566063dSJacob Faibussowitsch PetscCall(PCGAMGHashTableCreate(2 * nghosts + 1, &fgid_flid)); 30848a46eb9SPierre Jolivet for (kk = 0; kk < nghosts; kk++) PetscCall(PCGAMGHashTableAdd(&fgid_flid, flid_fgid[nloc + kk], nloc + kk)); 3092e68589bSMark F. Adams 3100cbbd2e1SMark F. Adams /* count selected -- same as number of cols of P */ 3110cbbd2e1SMark F. Adams for (nSelected = mm = 0; mm < nloc; mm++) { 312e78576d6SMark F. Adams PetscBool ise; 3139566063dSJacob Faibussowitsch PetscCall(PetscCDEmptyAt(agg_llists, mm, &ise)); 314e78576d6SMark F. Adams if (!ise) nSelected++; 3150cbbd2e1SMark F. Adams } 3169566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(a_Prol, &ii, &jj)); 31763a3b9bcSJacob Faibussowitsch PetscCheck((ii / nSAvec) == my0crs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ii %" PetscInt_FMT " /nSAvec %" PetscInt_FMT " != my0crs %" PetscInt_FMT, ii, nSAvec, my0crs); 31863a3b9bcSJacob Faibussowitsch PetscCheck(nSelected == (jj - ii) / nSAvec, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nSelected %" PetscInt_FMT " != (jj %" PetscInt_FMT " - ii %" PetscInt_FMT ")/nSAvec %" PetscInt_FMT, nSelected, jj, ii, nSAvec); 3190cbbd2e1SMark F. Adams 3202e68589bSMark F. Adams /* aloc space for coarse point data (output) */ 3210cbbd2e1SMark F. Adams out_data_stride = nSelected * nSAvec; 3222fa5cd67SKarl Rupp 3239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(out_data_stride * nSAvec, &out_data)); 32433812677SSatish Balay for (ii = 0; ii < out_data_stride * nSAvec; ii++) out_data[ii] = PETSC_MAX_REAL; 3250cbbd2e1SMark F. Adams *a_data_out = out_data; /* output - stride nSelected*nSAvec */ 3262e68589bSMark F. Adams 3272e68589bSMark F. Adams /* find points and set prolongation */ 328c8b0795cSMark F. Adams minsz = 100; 3290cbbd2e1SMark F. Adams for (mm = clid = 0; mm < nloc; mm++) { 3309566063dSJacob Faibussowitsch PetscCall(PetscCDSizeAt(agg_llists, mm, &jj)); 331e78576d6SMark F. Adams if (jj > 0) { 3320cbbd2e1SMark F. Adams const PetscInt lid = mm, cgid = my0crs + clid; 3330cbbd2e1SMark F. Adams PetscInt cids[100]; /* max bs */ 3340cbbd2e1SMark F. Adams PetscBLASInt asz = jj, M = asz * bs, N = nSAvec, INFO; 3352e68589bSMark F. Adams PetscBLASInt Mdata = M + ((N - M > 0) ? N - M : 0), LDA = Mdata, LWORK = N * bs; 3362e68589bSMark F. Adams PetscScalar *qqc, *qqr, *TAU, *WORK; 3372e68589bSMark F. Adams PetscInt *fids; 33865d7b583SSatish Balay PetscReal *data; 339b817416eSBarry Smith 3400cbbd2e1SMark F. Adams /* count agg */ 3410cbbd2e1SMark F. Adams if (asz < minsz) minsz = asz; 3420cbbd2e1SMark F. Adams 3430cbbd2e1SMark F. Adams /* get block */ 3449566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(Mdata * N, &qqc, M * N, &qqr, N, &TAU, LWORK, &WORK, M, &fids)); 3452e68589bSMark F. Adams 3462e68589bSMark F. Adams aggID = 0; 3479566063dSJacob Faibussowitsch PetscCall(PetscCDGetHeadPos(agg_llists, lid, &pos)); 348e78576d6SMark F. Adams while (pos) { 349e78576d6SMark F. Adams PetscInt gid1; 3509566063dSJacob Faibussowitsch PetscCall(PetscCDIntNdGetID(pos, &gid1)); 3519566063dSJacob Faibussowitsch PetscCall(PetscCDGetNextPos(agg_llists, lid, &pos)); 352e78576d6SMark F. Adams 3530cbbd2e1SMark F. Adams if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0; 3540cbbd2e1SMark F. Adams else { 3559566063dSJacob Faibussowitsch PetscCall(PCGAMGHashTableFind(&fgid_flid, gid1, &flid)); 35608401ef6SPierre Jolivet PetscCheck(flid >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot find gid1 in table"); 3570cbbd2e1SMark F. Adams } 3582e68589bSMark F. Adams /* copy in B_i matrix - column oriented */ 35965d7b583SSatish Balay data = &data_in[flid * bs]; 3603d3eaba7SBarry Smith for (ii = 0; ii < bs; ii++) { 3612e68589bSMark F. Adams for (jj = 0; jj < N; jj++) { 3620cbbd2e1SMark F. Adams PetscReal d = data[jj * data_stride + ii]; 3630cbbd2e1SMark F. Adams qqc[jj * Mdata + aggID * bs + ii] = d; 3642e68589bSMark F. Adams } 3652e68589bSMark F. Adams } 3662e68589bSMark F. Adams /* set fine IDs */ 3672e68589bSMark F. Adams for (kk = 0; kk < bs; kk++) fids[aggID * bs + kk] = flid_fgid[flid] * bs + kk; 3682e68589bSMark F. Adams aggID++; 3690cbbd2e1SMark F. Adams } 3702e68589bSMark F. Adams 3712e68589bSMark F. Adams /* pad with zeros */ 3722e68589bSMark F. Adams for (ii = asz * bs; ii < Mdata; ii++) { 373ad540459SPierre Jolivet for (jj = 0; jj < N; jj++, kk++) qqc[jj * Mdata + ii] = .0; 3742e68589bSMark F. Adams } 3752e68589bSMark F. Adams 3762e68589bSMark F. Adams /* QR */ 3779566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 378792fecdfSBarry Smith PetscCallBLAS("LAPACKgeqrf", LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 3799566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPop()); 38008401ef6SPierre Jolivet PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xGEQRF error"); 3812e68589bSMark F. Adams /* get R - column oriented - output B_{i+1} */ 3822e68589bSMark F. Adams { 3832e68589bSMark F. Adams PetscReal *data = &out_data[clid * nSAvec]; 3842e68589bSMark F. Adams for (jj = 0; jj < nSAvec; jj++) { 3852e68589bSMark F. Adams for (ii = 0; ii < nSAvec; ii++) { 38608401ef6SPierre Jolivet PetscCheck(data[jj * out_data_stride + ii] == PETSC_MAX_REAL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "data[jj*out_data_stride + ii] != %e", (double)PETSC_MAX_REAL); 3870cbbd2e1SMark F. Adams if (ii <= jj) data[jj * out_data_stride + ii] = PetscRealPart(qqc[jj * Mdata + ii]); 3880cbbd2e1SMark F. Adams else data[jj * out_data_stride + ii] = 0.; 3892e68589bSMark F. Adams } 3902e68589bSMark F. Adams } 3912e68589bSMark F. Adams } 3922e68589bSMark F. Adams 3932e68589bSMark F. Adams /* get Q - row oriented */ 394792fecdfSBarry Smith PetscCallBLAS("LAPACKorgqr", LAPACKorgqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 39563a3b9bcSJacob Faibussowitsch PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xORGQR error arg %" PetscBLASInt_FMT, -INFO); 3962e68589bSMark F. Adams 3972e68589bSMark F. Adams for (ii = 0; ii < M; ii++) { 398ad540459SPierre Jolivet for (jj = 0; jj < N; jj++) qqr[N * ii + jj] = qqc[jj * Mdata + ii]; 3992e68589bSMark F. Adams } 4002e68589bSMark F. Adams 4012e68589bSMark F. Adams /* add diagonal block of P0 */ 4029371c9d4SSatish Balay for (kk = 0; kk < N; kk++) { cids[kk] = N * cgid + kk; /* global col IDs in P0 */ } 4039566063dSJacob Faibussowitsch PetscCall(MatSetValues(a_Prol, M, fids, N, cids, qqr, INSERT_VALUES)); 4049566063dSJacob Faibussowitsch PetscCall(PetscFree5(qqc, qqr, TAU, WORK, fids)); 405b43b03e9SMark F. Adams clid++; 4060cbbd2e1SMark F. Adams } /* coarse agg */ 4070cbbd2e1SMark F. Adams } /* for all fine nodes */ 4089566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(a_Prol, MAT_FINAL_ASSEMBLY)); 4099566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(a_Prol, MAT_FINAL_ASSEMBLY)); 4109566063dSJacob Faibussowitsch PetscCall(PCGAMGHashTableDestroy(&fgid_flid)); 4113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4122e68589bSMark F. Adams } 4132e68589bSMark F. Adams 414d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_GAMG_AGG(PC pc, PetscViewer viewer) 415d71ae5a4SJacob Faibussowitsch { 4165adeb434SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 4175adeb434SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 4185adeb434SBarry Smith PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 4195adeb434SBarry Smith 4205adeb434SBarry Smith PetscFunctionBegin; 4219566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " AGG specific options\n")); 422bae903cbSmarkadams4 PetscCall(PetscViewerASCIIPrintf(viewer, " Number of levels to square graph %d\n", (int)pc_gamg_agg->aggressive_coarsening_levels)); 423bae903cbSmarkadams4 PetscCall(PetscViewerASCIIPrintf(viewer, " Number smoothing steps %d\n", (int)pc_gamg_agg->nsmooths)); 4243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4255adeb434SBarry Smith } 4265adeb434SBarry Smith 4272d776b49SBarry Smith static PetscErrorCode PCGAMGCreateGraph_AGG(PC pc, Mat Amat, Mat *a_Gmat) 428d71ae5a4SJacob Faibussowitsch { 429c8b0795cSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 430c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 431c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 4322d776b49SBarry Smith const PetscReal vfilter = pc_gamg->threshold[pc_gamg->current_level]; 4332d776b49SBarry Smith PetscBool ishem; 4342d776b49SBarry Smith const char *prefix; 435c8b0795cSMark F. Adams 436c8b0795cSMark F. Adams PetscFunctionBegin; 437849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0)); 4382d776b49SBarry Smith /* Note: depending on the algorithm that will be used for computing the coarse grid points this should pass PETSC_TRUE or PETSC_FALSE as the first argument */ 4392d776b49SBarry Smith 4402d776b49SBarry Smith /* MATCOARSENHEM requires numerical weights for edges so ensure they are computed */ 4412d776b49SBarry Smith PetscCall(MatCoarsenCreate(PetscObjectComm((PetscObject)pc), &pc_gamg_agg->crs)); 4422d776b49SBarry Smith PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, &prefix)); 4432d776b49SBarry Smith PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc_gamg_agg->crs, prefix)); 4442d776b49SBarry Smith PetscCall(MatCoarsenSetFromOptions(pc_gamg_agg->crs)); 4452d776b49SBarry Smith PetscCall(PetscObjectTypeCompare((PetscObject)pc_gamg_agg->crs, MATCOARSENHEM, &ishem)); 4462d776b49SBarry Smith 4472d776b49SBarry Smith PetscCall(MatCreateGraph(Amat, PETSC_TRUE, (vfilter >= 0 || ishem) ? PETSC_TRUE : PETSC_FALSE, vfilter, a_Gmat)); 448849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0)); 4493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 450c8b0795cSMark F. Adams } 451c8b0795cSMark F. Adams 452c8b0795cSMark F. Adams /* 453bae903cbSmarkadams4 PCGAMGCoarsen_AGG - supports squaring the graph (deprecated) and new graph for 454bae903cbSmarkadams4 communication of QR data used with HEM and MISk coarsening 455c8b0795cSMark F. Adams 456c8b0795cSMark F. Adams Input Parameter: 457e0940f08SMark F. Adams . a_pc - this 458bae903cbSmarkadams4 459e0940f08SMark F. Adams Input/Output Parameter: 460bae903cbSmarkadams4 . a_Gmat1 - graph to coarsen (in), graph off processor edges for QR gather scatter (out) 461bae903cbSmarkadams4 462c8b0795cSMark F. Adams Output Parameter: 4630cbbd2e1SMark F. Adams . agg_lists - list of aggregates 464bae903cbSmarkadams4 465c8b0795cSMark F. Adams */ 466d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc, Mat *a_Gmat1, PetscCoarsenData **agg_lists) 467d71ae5a4SJacob Faibussowitsch { 468e0940f08SMark F. Adams PC_MG *mg = (PC_MG *)a_pc->data; 469c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 470c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 471bae903cbSmarkadams4 Mat mat, Gmat1 = *a_Gmat1; /* aggressive graph */ 4720cbbd2e1SMark F. Adams IS perm; 473bae903cbSmarkadams4 PetscInt Istart, Iend, Ii, nloc, bs, nn; 474bae903cbSmarkadams4 PetscInt *permute, *degree; 475c8b0795cSMark F. Adams PetscBool *bIndexSet; 476aea53286SMark Adams PetscReal hashfact; 477e2c00dcbSBarry Smith PetscInt iSwapIndex; 4783b50add6SMark Adams PetscRandom random; 479c8b0795cSMark F. Adams 480c8b0795cSMark F. Adams PetscFunctionBegin; 481849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0)); 482bae903cbSmarkadams4 PetscCall(MatGetLocalSize(Gmat1, &nn, NULL)); 4839566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Gmat1, &bs)); 48463a3b9bcSJacob Faibussowitsch PetscCheck(bs == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bs %" PetscInt_FMT " must be 1", bs); 485bae903cbSmarkadams4 nloc = nn / bs; 486c8b0795cSMark F. Adams 4875cfd4bd9SMark Adams /* get MIS aggs - randomize */ 488bae903cbSmarkadams4 PetscCall(PetscMalloc2(nloc, &permute, nloc, °ree)); 4899566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nloc, &bIndexSet)); 4906e09b0e3SStefano Zampini for (Ii = 0; Ii < nloc; Ii++) permute[Ii] = Ii; 4919566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &random)); 4929566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Gmat1, &Istart, &Iend)); 493c8b0795cSMark F. Adams for (Ii = 0; Ii < nloc; Ii++) { 494bae903cbSmarkadams4 PetscInt nc; 495bae903cbSmarkadams4 PetscCall(MatGetRow(Gmat1, Istart + Ii, &nc, NULL, NULL)); 496bae903cbSmarkadams4 degree[Ii] = nc; 497bae903cbSmarkadams4 PetscCall(MatRestoreRow(Gmat1, Istart + Ii, &nc, NULL, NULL)); 498bae903cbSmarkadams4 } 499bae903cbSmarkadams4 for (Ii = 0; Ii < nloc; Ii++) { 5009566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(random, &hashfact)); 501aea53286SMark Adams iSwapIndex = (PetscInt)(hashfact * nloc) % nloc; 502c8b0795cSMark F. Adams if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) { 503c8b0795cSMark F. Adams PetscInt iTemp = permute[iSwapIndex]; 504c8b0795cSMark F. Adams permute[iSwapIndex] = permute[Ii]; 505c8b0795cSMark F. Adams permute[Ii] = iTemp; 506bae903cbSmarkadams4 iTemp = degree[iSwapIndex]; 507bae903cbSmarkadams4 degree[iSwapIndex] = degree[Ii]; 508bae903cbSmarkadams4 degree[Ii] = iTemp; 509c8b0795cSMark F. Adams bIndexSet[iSwapIndex] = PETSC_TRUE; 510c8b0795cSMark F. Adams } 511c8b0795cSMark F. Adams } 512bae903cbSmarkadams4 // create minimum degree ordering 513bae903cbSmarkadams4 PetscCall(PetscSortIntWithArray(nloc, degree, permute)); 514bae903cbSmarkadams4 5159566063dSJacob Faibussowitsch PetscCall(PetscFree(bIndexSet)); 5169566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&random)); 5179566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm)); 518849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0)); 5192d776b49SBarry Smith PetscCall(MatCoarsenSetGreedyOrdering(pc_gamg_agg->crs, perm)); 5202d776b49SBarry Smith PetscCall(MatCoarsenSetAdjacency(pc_gamg_agg->crs, Gmat1)); 5212d776b49SBarry Smith PetscCall(MatCoarsenSetStrictAggs(pc_gamg_agg->crs, PETSC_TRUE)); 5222d776b49SBarry Smith if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels) PetscCall(MatCoarsenMISKSetDistance(pc_gamg_agg->crs, 2)); // hardwire to MIS-2 5232d776b49SBarry Smith else PetscCall(MatCoarsenMISKSetDistance(pc_gamg_agg->crs, 1)); // MIS 5242d776b49SBarry Smith PetscCall(MatCoarsenApply(pc_gamg_agg->crs)); 5252d776b49SBarry Smith PetscCall(MatCoarsenViewFromOptions(pc_gamg_agg->crs, NULL, "-mat_coarsen_view")); 5262d776b49SBarry Smith PetscCall(MatCoarsenGetData(pc_gamg_agg->crs, agg_lists)); /* output */ 5272d776b49SBarry Smith PetscCall(MatCoarsenDestroy(&pc_gamg_agg->crs)); 528b43b03e9SMark F. Adams 5299566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 530bae903cbSmarkadams4 PetscCall(PetscFree2(permute, degree)); 531849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0)); 5329f3f12c8SMark F. Adams 533bae903cbSmarkadams4 { 534bae903cbSmarkadams4 PetscCoarsenData *llist = *agg_lists; 5359ab59c8bSMark Adams /* see if we have a matrix that takes precedence (returned from MatCoarsenApply) */ 5369566063dSJacob Faibussowitsch PetscCall(PetscCDGetMat(llist, &mat)); 5370cbbd2e1SMark F. Adams if (mat) { 5389566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Gmat1)); 5390cbbd2e1SMark F. Adams *a_Gmat1 = mat; /* output */ 5400cbbd2e1SMark F. Adams } 5410cbbd2e1SMark F. Adams } 542849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0)); 5433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 544c8b0795cSMark F. Adams } 545c8b0795cSMark F. Adams 546c8b0795cSMark F. Adams /* 5470cbbd2e1SMark F. Adams PCGAMGProlongator_AGG 548c8b0795cSMark F. Adams 549c8b0795cSMark F. Adams Input Parameter: 550c8b0795cSMark F. Adams . pc - this 551c8b0795cSMark F. Adams . Amat - matrix on this fine level 552c8b0795cSMark F. Adams . Graph - used to get ghost data for nodes in 5530cbbd2e1SMark F. Adams . agg_lists - list of aggregates 554c8b0795cSMark F. Adams Output Parameter: 555c8b0795cSMark F. Adams . a_P_out - prolongation operator to the next level 556c8b0795cSMark F. Adams */ 557d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGProlongator_AGG(PC pc, Mat Amat, Mat Gmat, PetscCoarsenData *agg_lists, Mat *a_P_out) 558d71ae5a4SJacob Faibussowitsch { 5592e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 5602e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 5614a2f8832SBarry Smith const PetscInt col_bs = pc_gamg->data_cell_cols; 562c8b0795cSMark F. Adams PetscInt Istart, Iend, nloc, ii, jj, kk, my0, nLocalSelected, bs; 563c8b0795cSMark F. Adams Mat Prol; 564d5d25401SBarry Smith PetscMPIInt size; 5653b4367a7SBarry Smith MPI_Comm comm; 566c8b0795cSMark F. Adams PetscReal *data_w_ghost; 567c8b0795cSMark F. Adams PetscInt myCrs0, nbnodes = 0, *flid_fgid; 568d9558ea9SBarry Smith MatType mtype; 5692e68589bSMark F. Adams 5702e68589bSMark F. Adams PetscFunctionBegin; 5719566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm)); 57208401ef6SPierre Jolivet PetscCheck(col_bs >= 1, comm, PETSC_ERR_PLIB, "Column bs cannot be less than 1"); 573849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0)); 5749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 5759566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Amat, &Istart, &Iend)); 5769566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Amat, &bs)); 5779371c9d4SSatish Balay nloc = (Iend - Istart) / bs; 5789371c9d4SSatish Balay my0 = Istart / bs; 57963a3b9bcSJacob Faibussowitsch PetscCheck((Iend - Istart) % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "(Iend %" PetscInt_FMT " - Istart %" PetscInt_FMT ") not divisible by bs %" PetscInt_FMT, Iend, Istart, bs); 5802e68589bSMark F. Adams 5812e68589bSMark F. Adams /* get 'nLocalSelected' */ 5820cbbd2e1SMark F. Adams for (ii = 0, nLocalSelected = 0; ii < nloc; ii++) { 583e78576d6SMark F. Adams PetscBool ise; 5840cbbd2e1SMark F. Adams /* filter out singletons 0 or 1? */ 5859566063dSJacob Faibussowitsch PetscCall(PetscCDEmptyAt(agg_lists, ii, &ise)); 586e78576d6SMark F. Adams if (!ise) nLocalSelected++; 5872e68589bSMark F. Adams } 5882e68589bSMark F. Adams 5892e68589bSMark F. Adams /* create prolongator, create P matrix */ 5909566063dSJacob Faibussowitsch PetscCall(MatGetType(Amat, &mtype)); 5919566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &Prol)); 5929566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Prol, nloc * bs, nLocalSelected * col_bs, PETSC_DETERMINE, PETSC_DETERMINE)); 5939566063dSJacob Faibussowitsch PetscCall(MatSetBlockSizes(Prol, bs, col_bs)); 5949566063dSJacob Faibussowitsch PetscCall(MatSetType(Prol, mtype)); 5959566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(Prol, col_bs, NULL)); 5969566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(Prol, col_bs, NULL, col_bs, NULL)); 5972e68589bSMark F. Adams 5982e68589bSMark F. Adams /* can get all points "removed" */ 5999566063dSJacob Faibussowitsch PetscCall(MatGetSize(Prol, &kk, &ii)); 6007f66b68fSMark Adams if (!ii) { 60163a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: No selected points on coarse grid\n", ((PetscObject)pc)->prefix)); 6029566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Prol)); 6030298fd71SBarry Smith *a_P_out = NULL; /* out */ 604849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0)); 6053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6062e68589bSMark F. Adams } 60763a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: New grid %" PetscInt_FMT " nodes\n", ((PetscObject)pc)->prefix, ii / col_bs)); 6089566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk)); 6090cbbd2e1SMark F. Adams 61063a3b9bcSJacob Faibussowitsch PetscCheck((kk - myCrs0) % col_bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "(kk %" PetscInt_FMT " -myCrs0 %" PetscInt_FMT ") not divisible by col_bs %" PetscInt_FMT, kk, myCrs0, col_bs); 611c8b0795cSMark F. Adams myCrs0 = myCrs0 / col_bs; 61263a3b9bcSJacob Faibussowitsch PetscCheck((kk / col_bs - myCrs0) == nLocalSelected, PETSC_COMM_SELF, PETSC_ERR_PLIB, "(kk %" PetscInt_FMT "/col_bs %" PetscInt_FMT " - myCrs0 %" PetscInt_FMT ") != nLocalSelected %" PetscInt_FMT ")", kk, col_bs, myCrs0, nLocalSelected); 6132e68589bSMark F. Adams 6142e68589bSMark F. Adams /* create global vector of data in 'data_w_ghost' */ 615849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0)); 616bae903cbSmarkadams4 if (size > 1) { /* get ghost null space data */ 6172e68589bSMark F. Adams PetscReal *tmp_gdata, *tmp_ldata, *tp2; 6189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc, &tmp_ldata)); 6194a2f8832SBarry Smith for (jj = 0; jj < col_bs; jj++) { 620c8b0795cSMark F. Adams for (kk = 0; kk < bs; kk++) { 621a2f3521dSMark F. Adams PetscInt ii, stride; 622c8b0795cSMark F. Adams const PetscReal *tp = pc_gamg->data + jj * bs * nloc + kk; 6232fa5cd67SKarl Rupp for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp; 6242fa5cd67SKarl Rupp 6259566063dSJacob Faibussowitsch PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata)); 626a2f3521dSMark F. Adams 627bae903cbSmarkadams4 if (!jj && !kk) { /* now I know how many total nodes - allocate TODO: move below and do in one 'col_bs' call */ 6289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(stride * bs * col_bs, &data_w_ghost)); 629a2f3521dSMark F. Adams nbnodes = bs * stride; 6302e68589bSMark F. Adams } 631a2f3521dSMark F. Adams tp2 = data_w_ghost + jj * bs * stride + kk; 6322fa5cd67SKarl Rupp for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii]; 6339566063dSJacob Faibussowitsch PetscCall(PetscFree(tmp_gdata)); 6342e68589bSMark F. Adams } 6352e68589bSMark F. Adams } 6369566063dSJacob Faibussowitsch PetscCall(PetscFree(tmp_ldata)); 637806fa848SBarry Smith } else { 638c8b0795cSMark F. Adams nbnodes = bs * nloc; 639c8b0795cSMark F. Adams data_w_ghost = (PetscReal *)pc_gamg->data; 6402e68589bSMark F. Adams } 6412e68589bSMark F. Adams 642bae903cbSmarkadams4 /* get 'flid_fgid' TODO - move up to get 'stride' and do get null space data above in one step (jj loop) */ 643c5df96a5SBarry Smith if (size > 1) { 6442e68589bSMark F. Adams PetscReal *fid_glid_loc, *fiddata; 645a2f3521dSMark F. Adams PetscInt stride; 6462e68589bSMark F. Adams 6479566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc, &fid_glid_loc)); 6482e68589bSMark F. Adams for (kk = 0; kk < nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0 + kk); 6499566063dSJacob Faibussowitsch PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata)); 650bae903cbSmarkadams4 PetscCall(PetscMalloc1(stride, &flid_fgid)); /* copy real data to in */ 651a2f3521dSMark F. Adams for (kk = 0; kk < stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk]; 6529566063dSJacob Faibussowitsch PetscCall(PetscFree(fiddata)); 653a2f3521dSMark F. Adams 65463a3b9bcSJacob Faibussowitsch PetscCheck(stride == nbnodes / bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "stride %" PetscInt_FMT " != nbnodes %" PetscInt_FMT "/bs %" PetscInt_FMT, stride, nbnodes, bs); 6559566063dSJacob Faibussowitsch PetscCall(PetscFree(fid_glid_loc)); 656806fa848SBarry Smith } else { 6579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc, &flid_fgid)); 6582e68589bSMark F. Adams for (kk = 0; kk < nloc; kk++) flid_fgid[kk] = my0 + kk; 6592e68589bSMark F. Adams } 660849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0)); 661bae903cbSmarkadams4 /* get P0 */ 662849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0)); 663c8b0795cSMark F. Adams { 6640298fd71SBarry Smith PetscReal *data_out = NULL; 6659566063dSJacob Faibussowitsch PetscCall(formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes, data_w_ghost, flid_fgid, &data_out, Prol)); 6669566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 6672fa5cd67SKarl Rupp 668c8b0795cSMark F. Adams pc_gamg->data = data_out; 6694a2f8832SBarry Smith pc_gamg->data_cell_rows = col_bs; 6704a2f8832SBarry Smith pc_gamg->data_sz = col_bs * col_bs * nLocalSelected; 671c8b0795cSMark F. Adams } 672849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0)); 67348a46eb9SPierre Jolivet if (size > 1) PetscCall(PetscFree(data_w_ghost)); 6749566063dSJacob Faibussowitsch PetscCall(PetscFree(flid_fgid)); 6752e68589bSMark F. Adams 676c8b0795cSMark F. Adams *a_P_out = Prol; /* out */ 677ed4e82eaSPeter Brune 678849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0)); 6793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 680c8b0795cSMark F. Adams } 681c8b0795cSMark F. Adams 682c8b0795cSMark F. Adams /* 683fd1112cbSBarry Smith PCGAMGOptProlongator_AGG 684c8b0795cSMark F. Adams 685c8b0795cSMark F. Adams Input Parameter: 686c8b0795cSMark F. Adams . pc - this 687c8b0795cSMark F. Adams . Amat - matrix on this fine level 688c8b0795cSMark F. Adams In/Output Parameter: 6892a808120SBarry Smith . a_P - prolongation operator to the next level 690c8b0795cSMark F. Adams */ 691d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGOptProlongator_AGG(PC pc, Mat Amat, Mat *a_P) 692d71ae5a4SJacob Faibussowitsch { 693c8b0795cSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 694c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 695c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 696c8b0795cSMark F. Adams PetscInt jj; 697c8b0795cSMark F. Adams Mat Prol = *a_P; 6983b4367a7SBarry Smith MPI_Comm comm; 6992a808120SBarry Smith KSP eksp; 7002a808120SBarry Smith Vec bb, xx; 7012a808120SBarry Smith PC epc; 7022a808120SBarry Smith PetscReal alpha, emax, emin; 703c8b0795cSMark F. Adams 704c8b0795cSMark F. Adams PetscFunctionBegin; 7059566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm)); 706849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0)); 707c8b0795cSMark F. Adams 708d5d25401SBarry Smith /* compute maximum singular value of operator to be used in smoother */ 7092a808120SBarry Smith if (0 < pc_gamg_agg->nsmooths) { 71018c3aa7eSMark /* get eigen estimates */ 71118c3aa7eSMark if (pc_gamg->emax > 0) { 71218c3aa7eSMark emin = pc_gamg->emin; 71318c3aa7eSMark emax = pc_gamg->emax; 71418c3aa7eSMark } else { 7150ed2132dSStefano Zampini const char *prefix; 7160ed2132dSStefano Zampini 7179566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Amat, &bb, NULL)); 7189566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Amat, &xx, NULL)); 7199566063dSJacob Faibussowitsch PetscCall(KSPSetNoisy_Private(bb)); 720e696ed0bSMark F. Adams 7219566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm, &eksp)); 7229566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 7239566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(eksp, prefix)); 7249566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(eksp, "pc_gamg_esteig_")); 72573f7197eSJed Brown { 726b94d7dedSBarry Smith PetscBool isset, sflg; 727b94d7dedSBarry Smith PetscCall(MatIsSPDKnown(Amat, &isset, &sflg)); 728b94d7dedSBarry Smith if (isset && sflg) PetscCall(KSPSetType(eksp, KSPCG)); 729d24ecf33SMark } 7309566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(eksp, pc->erroriffailure)); 7319566063dSJacob Faibussowitsch PetscCall(KSPSetNormType(eksp, KSP_NORM_NONE)); 732c2436cacSMark F. Adams 7339566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(eksp, PETSC_FALSE)); 7349566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(eksp, Amat, Amat)); 7352e68589bSMark F. Adams 7369566063dSJacob Faibussowitsch PetscCall(KSPGetPC(eksp, &epc)); 7379566063dSJacob Faibussowitsch PetscCall(PCSetType(epc, PCJACOBI)); /* smoother in smoothed agg. */ 738c2436cacSMark F. Adams 7399566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10)); // 10 is safer, but 5 is often fine, can override with -pc_gamg_esteig_ksp_max_it -mg_levels_ksp_chebyshev_esteig 0,0.25,0,1.2 740074aec55SMark Adams 7419566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(eksp)); 7429566063dSJacob Faibussowitsch PetscCall(KSPSetComputeSingularValues(eksp, PETSC_TRUE)); 7439566063dSJacob Faibussowitsch PetscCall(KSPSolve(eksp, bb, xx)); 7449566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(eksp, pc, xx)); 7452e68589bSMark F. Adams 7469566063dSJacob Faibussowitsch PetscCall(KSPComputeExtremeSingularValues(eksp, &emax, &emin)); 7479566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: Smooth P0: max eigen=%e min=%e PC=%s\n", ((PetscObject)pc)->prefix, (double)emax, (double)emin, PCJACOBI)); 7489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xx)); 7499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bb)); 7509566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&eksp)); 7512e68589bSMark F. Adams } 75218c3aa7eSMark if (pc_gamg->use_sa_esteig) { 75318c3aa7eSMark mg->min_eigen_DinvA[pc_gamg->current_level] = emin; 75418c3aa7eSMark mg->max_eigen_DinvA[pc_gamg->current_level] = emax; 75563a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: Smooth P0: level %" PetscInt_FMT ", cache spectra %g %g\n", ((PetscObject)pc)->prefix, pc_gamg->current_level, (double)emin, (double)emax)); 75618c3aa7eSMark } else { 75718c3aa7eSMark mg->min_eigen_DinvA[pc_gamg->current_level] = 0; 75818c3aa7eSMark mg->max_eigen_DinvA[pc_gamg->current_level] = 0; 75918c3aa7eSMark } 76018c3aa7eSMark } else { 76118c3aa7eSMark mg->min_eigen_DinvA[pc_gamg->current_level] = 0; 76218c3aa7eSMark mg->max_eigen_DinvA[pc_gamg->current_level] = 0; 76318c3aa7eSMark } 7642e68589bSMark F. Adams 7652a808120SBarry Smith /* smooth P0 */ 7662a808120SBarry Smith for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) { 7672a808120SBarry Smith Mat tMat; 7682a808120SBarry Smith Vec diag; 7692a808120SBarry Smith 770849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0)); 7712a808120SBarry Smith 7722e68589bSMark F. Adams /* smooth P1 := (I - omega/lam D^{-1}A)P0 */ 7739566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0)); 7749566063dSJacob Faibussowitsch PetscCall(MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat)); 7759566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0)); 7769566063dSJacob Faibussowitsch PetscCall(MatProductClear(tMat)); 7779566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Amat, &diag, NULL)); 7789566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(Amat, diag)); /* effectively PCJACOBI */ 7799566063dSJacob Faibussowitsch PetscCall(VecReciprocal(diag)); 7809566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(tMat, diag, NULL)); 7819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&diag)); 782b4da3a1bSStefano Zampini 783d5d25401SBarry Smith /* TODO: Set a PCFailedReason and exit the building of the AMG preconditioner */ 78408401ef6SPierre Jolivet PetscCheck(emax != 0.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Computed maximum singular value as zero"); 785d5d25401SBarry Smith /* TODO: Document the 1.4 and don't hardwire it in this routine */ 786b4ec6429SMark F. Adams alpha = -1.4 / emax; 787b4da3a1bSStefano Zampini 7889566063dSJacob Faibussowitsch PetscCall(MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN)); 7899566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Prol)); 7902e68589bSMark F. Adams Prol = tMat; 791849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0)); 7922e68589bSMark F. Adams } 793849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0)); 794c8b0795cSMark F. Adams *a_P = Prol; 7953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7962e68589bSMark F. Adams } 7972e68589bSMark F. Adams 798c8b0795cSMark F. Adams /* 799c8b0795cSMark F. Adams PCCreateGAMG_AGG 8002e68589bSMark F. Adams 801c8b0795cSMark F. Adams Input Parameter: 802c8b0795cSMark F. Adams . pc - 803c8b0795cSMark F. Adams */ 804d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCreateGAMG_AGG(PC pc) 805d71ae5a4SJacob Faibussowitsch { 806c8b0795cSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 807c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 808c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg; 8092e68589bSMark F. Adams 810c8b0795cSMark F. Adams PetscFunctionBegin; 811c8b0795cSMark F. Adams /* create sub context for SA */ 8124dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&pc_gamg_agg)); 813c8b0795cSMark F. Adams pc_gamg->subctx = pc_gamg_agg; 814c8b0795cSMark F. Adams 8151ab5ffc9SJed Brown pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG; 8169b8ffb57SJed Brown pc_gamg->ops->destroy = PCDestroy_GAMG_AGG; 817c8b0795cSMark F. Adams /* reset does not do anything; setup not virtual */ 818c8b0795cSMark F. Adams 819c8b0795cSMark F. Adams /* set internal function pointers */ 8202d776b49SBarry Smith pc_gamg->ops->creategraph = PCGAMGCreateGraph_AGG; 8211ab5ffc9SJed Brown pc_gamg->ops->coarsen = PCGAMGCoarsen_AGG; 8221ab5ffc9SJed Brown pc_gamg->ops->prolongator = PCGAMGProlongator_AGG; 823fd1112cbSBarry Smith pc_gamg->ops->optprolongator = PCGAMGOptProlongator_AGG; 8241ab5ffc9SJed Brown pc_gamg->ops->createdefaultdata = PCSetData_AGG; 8255adeb434SBarry Smith pc_gamg->ops->view = PCView_GAMG_AGG; 826c8b0795cSMark F. Adams 827599482c2Smarkadams4 pc_gamg_agg->aggressive_coarsening_levels = 1; 828cab9ed1eSBarry Smith pc_gamg_agg->nsmooths = 1; 829cab9ed1eSBarry Smith 8309566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", PCGAMGSetNSmooths_AGG)); 831bae903cbSmarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", PCGAMGSetAggressiveLevels_AGG)); 8329566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_AGG)); 8333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8342e68589bSMark F. Adams } 835