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) 13d529f056Smarkadams4 PetscInt aggressive_mis_k; // the k in MIS-k 14d529f056Smarkadams4 PetscBool use_aggressive_square_graph; 15d529f056Smarkadams4 PetscBool use_minimum_degree_ordering; 162d776b49SBarry Smith MatCoarsen crs; 172e68589bSMark F. Adams } PC_GAMG_AGG; 182e68589bSMark F. Adams 192e68589bSMark F. Adams /*@ 20f1580f4eSBarry Smith PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical) used for multigrid on all the levels 212e68589bSMark F. Adams 22c3339decSBarry Smith Logically Collective 232e68589bSMark F. Adams 242e68589bSMark F. Adams Input Parameters: 2520f4b53cSBarry Smith + pc - the preconditioner context 2620f4b53cSBarry Smith - n - the number of smooths 272e68589bSMark F. Adams 282e68589bSMark F. Adams Options Database Key: 29cab9ed1eSBarry Smith . -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation 302e68589bSMark F. Adams 312e68589bSMark F. Adams Level: intermediate 322e68589bSMark F. Adams 33*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCGAMG` 342e68589bSMark F. Adams @*/ 35d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n) 36d71ae5a4SJacob Faibussowitsch { 372e68589bSMark F. Adams PetscFunctionBegin; 382e68589bSMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 39d5d25401SBarry Smith PetscValidLogicalCollectiveInt(pc, n, 2); 40cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetNSmooths_C", (PC, PetscInt), (pc, n)); 413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 422e68589bSMark F. Adams } 432e68589bSMark F. Adams 44d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetNSmooths_AGG(PC pc, PetscInt n) 45d71ae5a4SJacob Faibussowitsch { 462e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 472e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 48c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 492e68589bSMark F. Adams 502e68589bSMark F. Adams PetscFunctionBegin; 51c8b0795cSMark F. Adams pc_gamg_agg->nsmooths = n; 523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 53c8b0795cSMark F. Adams } 54c8b0795cSMark F. Adams 55c8b0795cSMark F. Adams /*@ 56f1580f4eSBarry Smith PCGAMGSetAggressiveLevels - Use aggressive coarsening on first n levels 57ef4ad70eSMark F. Adams 58c3339decSBarry Smith Logically Collective 59ef4ad70eSMark F. Adams 60ef4ad70eSMark F. Adams Input Parameters: 61cab9ed1eSBarry Smith + pc - the preconditioner context 62d5d25401SBarry Smith - n - 0, 1 or more 63ef4ad70eSMark F. Adams 64ef4ad70eSMark F. Adams Options Database Key: 65bae903cbSmarkadams4 . -pc_gamg_aggressive_coarsening <n,default = 1> - Number of levels to square the graph on before aggregating it 66af3c827dSMark Adams 67ef4ad70eSMark F. Adams Level: intermediate 68ef4ad70eSMark F. Adams 69*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGMISkSetAggressive()`, `PCGAMGSetAggressiveSquareGraph()`, `PCGAMGMISkSetMinDegreeOrdering()` 70ef4ad70eSMark F. Adams @*/ 71d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetAggressiveLevels(PC pc, PetscInt n) 72d71ae5a4SJacob Faibussowitsch { 73ef4ad70eSMark F. Adams PetscFunctionBegin; 74ef4ad70eSMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 75d5d25401SBarry Smith PetscValidLogicalCollectiveInt(pc, n, 2); 76bae903cbSmarkadams4 PetscTryMethod(pc, "PCGAMGSetAggressiveLevels_C", (PC, PetscInt), (pc, n)); 773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 78ef4ad70eSMark F. Adams } 79ef4ad70eSMark F. Adams 80d529f056Smarkadams4 /*@ 81d529f056Smarkadams4 PCGAMGMISkSetAggressive - Number (k) distance in MIS coarsening (>2 is 'aggressive') 82d529f056Smarkadams4 83d529f056Smarkadams4 Logically Collective 84d529f056Smarkadams4 85d529f056Smarkadams4 Input Parameters: 86d529f056Smarkadams4 + pc - the preconditioner context 87d529f056Smarkadams4 - n - 1 or more (default = 2) 88d529f056Smarkadams4 89d529f056Smarkadams4 Options Database Key: 90d529f056Smarkadams4 . -pc_gamg_aggressive_mis_k <n,default=2> - Number (k) distance in MIS coarsening (>2 is 'aggressive') 91d529f056Smarkadams4 92d529f056Smarkadams4 Level: intermediate 93d529f056Smarkadams4 94*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`, `PCGAMGSetAggressiveSquareGraph()`, `PCGAMGMISkSetMinDegreeOrdering()` 95d529f056Smarkadams4 @*/ 96d529f056Smarkadams4 PetscErrorCode PCGAMGMISkSetAggressive(PC pc, PetscInt n) 97d529f056Smarkadams4 { 98d529f056Smarkadams4 PetscFunctionBegin; 99d529f056Smarkadams4 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 100d529f056Smarkadams4 PetscValidLogicalCollectiveInt(pc, n, 2); 101d529f056Smarkadams4 PetscTryMethod(pc, "PCGAMGMISkSetAggressive_C", (PC, PetscInt), (pc, n)); 102d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 103d529f056Smarkadams4 } 104d529f056Smarkadams4 105d529f056Smarkadams4 /*@ 106d529f056Smarkadams4 PCGAMGSetAggressiveSquareGraph - Use graph square A'A for aggressive coarsening, old method 107d529f056Smarkadams4 108d529f056Smarkadams4 Logically Collective 109d529f056Smarkadams4 110d529f056Smarkadams4 Input Parameters: 111d529f056Smarkadams4 + pc - the preconditioner context 112d529f056Smarkadams4 - b - default false - MIS-k is faster 113d529f056Smarkadams4 114d529f056Smarkadams4 Options Database Key: 115d529f056Smarkadams4 . -pc_gamg_aggressive_square_graph <bool,default=false> - Use square graph (A'A) or MIS-k (k=2) for aggressive coarsening 116d529f056Smarkadams4 117d529f056Smarkadams4 Level: intermediate 118d529f056Smarkadams4 119*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`, `PCGAMGMISkSetAggressive()`, `PCGAMGMISkSetMinDegreeOrdering()` 120d529f056Smarkadams4 @*/ 121d529f056Smarkadams4 PetscErrorCode PCGAMGSetAggressiveSquareGraph(PC pc, PetscBool b) 122d529f056Smarkadams4 { 123d529f056Smarkadams4 PetscFunctionBegin; 124d529f056Smarkadams4 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 125d529f056Smarkadams4 PetscValidLogicalCollectiveBool(pc, b, 2); 126d529f056Smarkadams4 PetscTryMethod(pc, "PCGAMGSetAggressiveSquareGraph_C", (PC, PetscBool), (pc, b)); 127d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 128d529f056Smarkadams4 } 129d529f056Smarkadams4 130d529f056Smarkadams4 /*@ 131d529f056Smarkadams4 PCGAMGMISkSetMinDegreeOrdering - Use minimum degree ordering in greedy MIS algorithm 132d529f056Smarkadams4 133d529f056Smarkadams4 Logically Collective 134d529f056Smarkadams4 135d529f056Smarkadams4 Input Parameters: 136d529f056Smarkadams4 + pc - the preconditioner context 137d529f056Smarkadams4 - b - default true 138d529f056Smarkadams4 139d529f056Smarkadams4 Options Database Key: 140d529f056Smarkadams4 . -pc_gamg_mis_k_minimum_degree_ordering <bool,default=true> - Use minimum degree ordering in greedy MIS algorithm 141d529f056Smarkadams4 142d529f056Smarkadams4 Level: intermediate 143d529f056Smarkadams4 144*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`, `PCGAMGMISkSetAggressive()`, `PCGAMGSetAggressiveSquareGraph()` 145d529f056Smarkadams4 @*/ 146d529f056Smarkadams4 PetscErrorCode PCGAMGMISkSetMinDegreeOrdering(PC pc, PetscBool b) 147d529f056Smarkadams4 { 148d529f056Smarkadams4 PetscFunctionBegin; 149d529f056Smarkadams4 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 150d529f056Smarkadams4 PetscValidLogicalCollectiveBool(pc, b, 2); 151d529f056Smarkadams4 PetscTryMethod(pc, "PCGAMGMISkSetMinDegreeOrdering_C", (PC, PetscBool), (pc, b)); 152d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 153d529f056Smarkadams4 } 154d529f056Smarkadams4 155d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetAggressiveLevels_AGG(PC pc, PetscInt n) 156d71ae5a4SJacob Faibussowitsch { 157ef4ad70eSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 158ef4ad70eSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 159ef4ad70eSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 160ef4ad70eSMark F. Adams 161ef4ad70eSMark F. Adams PetscFunctionBegin; 162bae903cbSmarkadams4 pc_gamg_agg->aggressive_coarsening_levels = n; 1633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 164ef4ad70eSMark F. Adams } 165ef4ad70eSMark F. Adams 166d529f056Smarkadams4 static PetscErrorCode PCGAMGMISkSetAggressive_AGG(PC pc, PetscInt n) 167d529f056Smarkadams4 { 168d529f056Smarkadams4 PC_MG *mg = (PC_MG *)pc->data; 169d529f056Smarkadams4 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 170d529f056Smarkadams4 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 171d529f056Smarkadams4 172d529f056Smarkadams4 PetscFunctionBegin; 173d529f056Smarkadams4 pc_gamg_agg->aggressive_mis_k = n; 174d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 175d529f056Smarkadams4 } 176d529f056Smarkadams4 177d529f056Smarkadams4 static PetscErrorCode PCGAMGSetAggressiveSquareGraph_AGG(PC pc, PetscBool b) 178d529f056Smarkadams4 { 179d529f056Smarkadams4 PC_MG *mg = (PC_MG *)pc->data; 180d529f056Smarkadams4 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 181d529f056Smarkadams4 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 182d529f056Smarkadams4 183d529f056Smarkadams4 PetscFunctionBegin; 184d529f056Smarkadams4 pc_gamg_agg->use_aggressive_square_graph = b; 185d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 186d529f056Smarkadams4 } 187d529f056Smarkadams4 188d529f056Smarkadams4 static PetscErrorCode PCGAMGMISkSetMinDegreeOrdering_AGG(PC pc, PetscBool b) 189d529f056Smarkadams4 { 190d529f056Smarkadams4 PC_MG *mg = (PC_MG *)pc->data; 191d529f056Smarkadams4 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 192d529f056Smarkadams4 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 193d529f056Smarkadams4 194d529f056Smarkadams4 PetscFunctionBegin; 195d529f056Smarkadams4 pc_gamg_agg->use_minimum_degree_ordering = b; 196d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 197d529f056Smarkadams4 } 198d529f056Smarkadams4 199d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_GAMG_AGG(PC pc, PetscOptionItems *PetscOptionsObject) 200d71ae5a4SJacob Faibussowitsch { 2012e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 2022e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 203c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 2042e68589bSMark F. Adams 2052e68589bSMark F. Adams PetscFunctionBegin; 206d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "GAMG-AGG options"); 2072e68589bSMark F. Adams { 208bae903cbSmarkadams4 PetscBool flg; 2099566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_gamg_agg_nsmooths", "smoothing steps for smoothed aggregation, usually 1", "PCGAMGSetNSmooths", pc_gamg_agg->nsmooths, &pc_gamg_agg->nsmooths, NULL)); 210d529f056Smarkadams4 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)); 211d529f056Smarkadams4 if (!flg) { 212d529f056Smarkadams4 PetscCall( 213d529f056Smarkadams4 PetscOptionsInt("-pc_gamg_square_graph", "Number of aggressive coarsening (MIS-2) levels from finest (deprecated alias for -pc_gamg_aggressive_coarsening)", "PCGAMGSetAggressiveLevels", pc_gamg_agg->aggressive_coarsening_levels, &pc_gamg_agg->aggressive_coarsening_levels, NULL)); 214d529f056Smarkadams4 } else { 2159371c9d4SSatish Balay PetscCall( 2169371c9d4SSatish 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)); 217bae903cbSmarkadams4 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)); 218bae903cbSmarkadams4 } 219d529f056Smarkadams4 if (pc_gamg_agg->aggressive_coarsening_levels > 0) { 220d529f056Smarkadams4 PetscCall(PetscOptionsBool("-pc_gamg_aggressive_square_graph", "Use square graph (A'A) or MIS-k (k=2) for aggressive coarsening", "PCGAMGSetAggressiveSquareGraph", pc_gamg_agg->use_aggressive_square_graph, &pc_gamg_agg->use_aggressive_square_graph, NULL)); 221d529f056Smarkadams4 } 222d529f056Smarkadams4 PetscCall(PetscOptionsBool("-pc_gamg_mis_k_minimum_degree_ordering", "Use minimum degree ordering for greedy MIS", "PCGAMGMISkSetMinDegreeOrdering", pc_gamg_agg->use_minimum_degree_ordering, &pc_gamg_agg->use_minimum_degree_ordering, NULL)); 223d529f056Smarkadams4 PetscCall(PetscOptionsInt("-pc_gamg_aggressive_mis_k", "Number of levels of multigrid to use.", "PCGAMGMISkSetAggressive", pc_gamg_agg->aggressive_mis_k, &pc_gamg_agg->aggressive_mis_k, NULL)); 2242e68589bSMark F. Adams } 225d0609cedSBarry Smith PetscOptionsHeadEnd(); 2263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2272e68589bSMark F. Adams } 2282e68589bSMark F. Adams 229d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_GAMG_AGG(PC pc) 230d71ae5a4SJacob Faibussowitsch { 2312e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 2322e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 2332e68589bSMark F. Adams 2342e68589bSMark F. Adams PetscFunctionBegin; 2359566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->subctx)); 2362e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", NULL)); 237bae903cbSmarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", NULL)); 238d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetAggressive_C", NULL)); 239d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetMinDegreeOrdering_C", NULL)); 240d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveSquareGraph_C", NULL)); 241bae903cbSmarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL)); 2423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2432e68589bSMark F. Adams } 2442e68589bSMark F. Adams 2452e68589bSMark F. Adams /* 2462e68589bSMark F. Adams PCSetCoordinates_AGG 24720f4b53cSBarry Smith 24820f4b53cSBarry Smith Collective 2492e68589bSMark F. Adams 2502e68589bSMark F. Adams Input Parameter: 2512e68589bSMark F. Adams . pc - the preconditioner context 252145b44c9SPierre Jolivet . ndm - dimension of data (used for dof/vertex for Stokes) 253302f38e8SMark F. Adams . a_nloc - number of vertices local 254302f38e8SMark F. Adams . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...} 2552e68589bSMark F. Adams */ 2561e6b0712SBarry Smith 257d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords) 258d71ae5a4SJacob Faibussowitsch { 2592e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 2602e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 26169344418SMark F. Adams PetscInt arrsz, kk, ii, jj, nloc, ndatarows, ndf; 262a2f3521dSMark F. Adams Mat mat = pc->pmat; 2632e68589bSMark F. Adams 2642e68589bSMark F. Adams PetscFunctionBegin; 265a2f3521dSMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 266a2f3521dSMark F. Adams PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 267302f38e8SMark F. Adams nloc = a_nloc; 2682e68589bSMark F. Adams 2692e68589bSMark F. Adams /* SA: null space vectors */ 2709566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(mat, &ndf)); /* this does not work for Stokes */ 27169344418SMark F. Adams if (coords && ndf == 1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */ 272a2f3521dSMark F. Adams else if (coords) { 27363a3b9bcSJacob Faibussowitsch PetscCheck(ndm <= ndf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "degrees of motion %" PetscInt_FMT " > block size %" PetscInt_FMT, ndm, ndf); 27469344418SMark F. Adams pc_gamg->data_cell_cols = (ndm == 2 ? 3 : 6); /* displacement elasticity */ 275ad540459SPierre 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); 276806fa848SBarry Smith } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */ 27771959b99SBarry Smith pc_gamg->data_cell_rows = ndatarows = ndf; 27863a3b9bcSJacob 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); 279c8b0795cSMark F. Adams arrsz = nloc * pc_gamg->data_cell_rows * pc_gamg->data_cell_cols; 2802e68589bSMark F. Adams 2817f66b68fSMark Adams if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) { 2829566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 2839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(arrsz + 1, &pc_gamg->data)); 2842e68589bSMark F. Adams } 2852e68589bSMark F. Adams /* copy data in - column oriented */ 2862e68589bSMark F. Adams for (kk = 0; kk < nloc; kk++) { 28769344418SMark F. Adams const PetscInt M = nloc * pc_gamg->data_cell_rows; /* stride into data */ 28869344418SMark F. Adams PetscReal *data = &pc_gamg->data[kk * ndatarows]; /* start of cell */ 289c8b0795cSMark F. Adams if (pc_gamg->data_cell_cols == 1) *data = 1.0; 2902e68589bSMark F. Adams else { 29169344418SMark F. Adams /* translational modes */ 2922fa5cd67SKarl Rupp for (ii = 0; ii < ndatarows; ii++) { 2932fa5cd67SKarl Rupp for (jj = 0; jj < ndatarows; jj++) { 29469344418SMark F. Adams if (ii == jj) data[ii * M + jj] = 1.0; 2952e68589bSMark F. Adams else data[ii * M + jj] = 0.0; 2962fa5cd67SKarl Rupp } 2972fa5cd67SKarl Rupp } 29869344418SMark F. Adams 29969344418SMark F. Adams /* rotational modes */ 3002e68589bSMark F. Adams if (coords) { 30169344418SMark F. Adams if (ndm == 2) { 3022e68589bSMark F. Adams data += 2 * M; 3032e68589bSMark F. Adams data[0] = -coords[2 * kk + 1]; 3042e68589bSMark F. Adams data[1] = coords[2 * kk]; 305806fa848SBarry Smith } else { 3062e68589bSMark F. Adams data += 3 * M; 3079371c9d4SSatish Balay data[0] = 0.0; 3089371c9d4SSatish Balay data[M + 0] = coords[3 * kk + 2]; 3099371c9d4SSatish Balay data[2 * M + 0] = -coords[3 * kk + 1]; 3109371c9d4SSatish Balay data[1] = -coords[3 * kk + 2]; 3119371c9d4SSatish Balay data[M + 1] = 0.0; 3129371c9d4SSatish Balay data[2 * M + 1] = coords[3 * kk]; 3139371c9d4SSatish Balay data[2] = coords[3 * kk + 1]; 3149371c9d4SSatish Balay data[M + 2] = -coords[3 * kk]; 3159371c9d4SSatish Balay data[2 * M + 2] = 0.0; 3162e68589bSMark F. Adams } 3172e68589bSMark F. Adams } 3182e68589bSMark F. Adams } 3192e68589bSMark F. Adams } 3202e68589bSMark F. Adams pc_gamg->data_sz = arrsz; 3213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3222e68589bSMark F. Adams } 3232e68589bSMark F. Adams 3242e68589bSMark F. Adams /* 325a2f3521dSMark F. Adams PCSetData_AGG - called if data is not set with PCSetCoordinates. 326a2f3521dSMark F. Adams Looks in Mat for near null space. 327a2f3521dSMark F. Adams Does not work for Stokes 3282e68589bSMark F. Adams 3292e68589bSMark F. Adams Input Parameter: 3302e68589bSMark F. Adams . pc - 331a2f3521dSMark F. Adams . a_A - matrix to get (near) null space out of. 3322e68589bSMark F. Adams */ 333d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A) 334d71ae5a4SJacob Faibussowitsch { 335b8cd405aSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 336b8cd405aSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 337b8cd405aSMark F. Adams MatNullSpace mnull; 33866f2ef4bSMark Adams 339b817416eSBarry Smith PetscFunctionBegin; 3409566063dSJacob Faibussowitsch PetscCall(MatGetNearNullSpace(a_A, &mnull)); 341b8cd405aSMark F. Adams if (!mnull) { 34266f2ef4bSMark Adams DM dm; 3439566063dSJacob Faibussowitsch PetscCall(PCGetDM(pc, &dm)); 34448a46eb9SPierre Jolivet if (!dm) PetscCall(MatGetDM(a_A, &dm)); 34566f2ef4bSMark Adams if (dm) { 34666f2ef4bSMark Adams PetscObject deformation; 347b0d30dd6SMatthew G. Knepley PetscInt Nf; 348b0d30dd6SMatthew G. Knepley 3499566063dSJacob Faibussowitsch PetscCall(DMGetNumFields(dm, &Nf)); 350b0d30dd6SMatthew G. Knepley if (Nf) { 3519566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, 0, NULL, &deformation)); 3529566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)deformation, "nearnullspace", (PetscObject *)&mnull)); 35348a46eb9SPierre Jolivet if (!mnull) PetscCall(PetscObjectQuery((PetscObject)deformation, "nullspace", (PetscObject *)&mnull)); 35466f2ef4bSMark Adams } 35566f2ef4bSMark Adams } 356b0d30dd6SMatthew G. Knepley } 35766f2ef4bSMark Adams 35866f2ef4bSMark Adams if (!mnull) { 359a2f3521dSMark F. Adams PetscInt bs, NN, MM; 3609566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(a_A, &bs)); 3619566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(a_A, &MM, &NN)); 36263a3b9bcSJacob Faibussowitsch PetscCheck(MM % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MM %" PetscInt_FMT " must be divisible by bs %" PetscInt_FMT, MM, bs); 3639566063dSJacob Faibussowitsch PetscCall(PCSetCoordinates_AGG(pc, bs, MM / bs, NULL)); 364806fa848SBarry Smith } else { 365b8cd405aSMark F. Adams PetscReal *nullvec; 366b8cd405aSMark F. Adams PetscBool has_const; 367b8cd405aSMark F. Adams PetscInt i, j, mlocal, nvec, bs; 368d5d25401SBarry Smith const Vec *vecs; 369d5d25401SBarry Smith const PetscScalar *v; 370b817416eSBarry Smith 3719566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(a_A, &mlocal, NULL)); 3729566063dSJacob Faibussowitsch PetscCall(MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs)); 373d19c4ebbSmarkadams4 for (i = 0; i < nvec; i++) { 374d19c4ebbSmarkadams4 PetscCall(VecGetLocalSize(vecs[i], &j)); 375d19c4ebbSmarkadams4 PetscCheck(j == mlocal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Attached null space vector size %" PetscInt_FMT " != matrix size %" PetscInt_FMT, j, mlocal); 376d19c4ebbSmarkadams4 } 377a0dea87dSMark Adams pc_gamg->data_sz = (nvec + !!has_const) * mlocal; 3789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((nvec + !!has_const) * mlocal, &nullvec)); 3799371c9d4SSatish Balay if (has_const) 3809371c9d4SSatish Balay for (i = 0; i < mlocal; i++) nullvec[i] = 1.0; 381b8cd405aSMark F. Adams for (i = 0; i < nvec; i++) { 3829566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(vecs[i], &v)); 383b8cd405aSMark F. Adams for (j = 0; j < mlocal; j++) nullvec[(i + !!has_const) * mlocal + j] = PetscRealPart(v[j]); 3849566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(vecs[i], &v)); 385b8cd405aSMark F. Adams } 386b8cd405aSMark F. Adams pc_gamg->data = nullvec; 387b8cd405aSMark F. Adams pc_gamg->data_cell_cols = (nvec + !!has_const); 3889566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(a_A, &bs)); 389b8cd405aSMark F. Adams pc_gamg->data_cell_rows = bs; 390b8cd405aSMark F. Adams } 3913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3922e68589bSMark F. Adams } 3932e68589bSMark F. Adams 3942e68589bSMark F. Adams /* 395bae903cbSmarkadams4 formProl0 - collect null space data for each aggregate, do QR, put R in coarse grid data and Q in P_0 3962e68589bSMark F. Adams 3972e68589bSMark F. Adams Input Parameter: 3982adfe9d3SBarry Smith . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices 3992adfe9d3SBarry Smith . bs - row block size 4000cbbd2e1SMark F. Adams . nSAvec - column bs of new P 4010cbbd2e1SMark F. Adams . my0crs - global index of start of locals 4022adfe9d3SBarry Smith . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec] 4032e68589bSMark F. Adams . data_in[data_stride*nSAvec] - local data on fine grid 4042e68589bSMark F. Adams . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist' 405bae903cbSmarkadams4 4062e68589bSMark F. Adams Output Parameter: 4072e68589bSMark F. Adams . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data 4082e68589bSMark F. Adams . a_Prol - prolongation operator 4092e68589bSMark F. Adams */ 410d71ae5a4SJacob 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) 411d71ae5a4SJacob Faibussowitsch { 412bd026e97SJed Brown PetscInt Istart, my0, Iend, nloc, clid, flid = 0, aggID, kk, jj, ii, mm, nSelected, minsz, nghosts, out_data_stride; 4133b4367a7SBarry Smith MPI_Comm comm; 4142e68589bSMark F. Adams PetscReal *out_data; 415539c167fSBarry Smith PetscCDIntNd *pos; 4161943db53SBarry Smith PCGAMGHashTable fgid_flid; 4170cbbd2e1SMark F. Adams 4182e68589bSMark F. Adams PetscFunctionBegin; 4199566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)a_Prol, &comm)); 4209566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(a_Prol, &Istart, &Iend)); 4219371c9d4SSatish Balay nloc = (Iend - Istart) / bs; 4229371c9d4SSatish Balay my0 = Istart / bs; 42363a3b9bcSJacob 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); 4240cbbd2e1SMark F. Adams Iend /= bs; 4250cbbd2e1SMark F. Adams nghosts = data_stride / bs - nloc; 4262e68589bSMark F. Adams 4279566063dSJacob Faibussowitsch PetscCall(PCGAMGHashTableCreate(2 * nghosts + 1, &fgid_flid)); 42848a46eb9SPierre Jolivet for (kk = 0; kk < nghosts; kk++) PetscCall(PCGAMGHashTableAdd(&fgid_flid, flid_fgid[nloc + kk], nloc + kk)); 4292e68589bSMark F. Adams 4300cbbd2e1SMark F. Adams /* count selected -- same as number of cols of P */ 4310cbbd2e1SMark F. Adams for (nSelected = mm = 0; mm < nloc; mm++) { 432e78576d6SMark F. Adams PetscBool ise; 4339566063dSJacob Faibussowitsch PetscCall(PetscCDEmptyAt(agg_llists, mm, &ise)); 434e78576d6SMark F. Adams if (!ise) nSelected++; 4350cbbd2e1SMark F. Adams } 4369566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(a_Prol, &ii, &jj)); 43763a3b9bcSJacob Faibussowitsch PetscCheck((ii / nSAvec) == my0crs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ii %" PetscInt_FMT " /nSAvec %" PetscInt_FMT " != my0crs %" PetscInt_FMT, ii, nSAvec, my0crs); 43863a3b9bcSJacob 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); 4390cbbd2e1SMark F. Adams 4402e68589bSMark F. Adams /* aloc space for coarse point data (output) */ 4410cbbd2e1SMark F. Adams out_data_stride = nSelected * nSAvec; 4422fa5cd67SKarl Rupp 4439566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(out_data_stride * nSAvec, &out_data)); 44433812677SSatish Balay for (ii = 0; ii < out_data_stride * nSAvec; ii++) out_data[ii] = PETSC_MAX_REAL; 4450cbbd2e1SMark F. Adams *a_data_out = out_data; /* output - stride nSelected*nSAvec */ 4462e68589bSMark F. Adams 4472e68589bSMark F. Adams /* find points and set prolongation */ 448c8b0795cSMark F. Adams minsz = 100; 4490cbbd2e1SMark F. Adams for (mm = clid = 0; mm < nloc; mm++) { 4509566063dSJacob Faibussowitsch PetscCall(PetscCDSizeAt(agg_llists, mm, &jj)); 451e78576d6SMark F. Adams if (jj > 0) { 4520cbbd2e1SMark F. Adams const PetscInt lid = mm, cgid = my0crs + clid; 4530cbbd2e1SMark F. Adams PetscInt cids[100]; /* max bs */ 4540cbbd2e1SMark F. Adams PetscBLASInt asz = jj, M = asz * bs, N = nSAvec, INFO; 4552e68589bSMark F. Adams PetscBLASInt Mdata = M + ((N - M > 0) ? N - M : 0), LDA = Mdata, LWORK = N * bs; 4562e68589bSMark F. Adams PetscScalar *qqc, *qqr, *TAU, *WORK; 4572e68589bSMark F. Adams PetscInt *fids; 45865d7b583SSatish Balay PetscReal *data; 459b817416eSBarry Smith 4600cbbd2e1SMark F. Adams /* count agg */ 4610cbbd2e1SMark F. Adams if (asz < minsz) minsz = asz; 4620cbbd2e1SMark F. Adams 4630cbbd2e1SMark F. Adams /* get block */ 4649566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(Mdata * N, &qqc, M * N, &qqr, N, &TAU, LWORK, &WORK, M, &fids)); 4652e68589bSMark F. Adams 4662e68589bSMark F. Adams aggID = 0; 4679566063dSJacob Faibussowitsch PetscCall(PetscCDGetHeadPos(agg_llists, lid, &pos)); 468e78576d6SMark F. Adams while (pos) { 469e78576d6SMark F. Adams PetscInt gid1; 4709566063dSJacob Faibussowitsch PetscCall(PetscCDIntNdGetID(pos, &gid1)); 4719566063dSJacob Faibussowitsch PetscCall(PetscCDGetNextPos(agg_llists, lid, &pos)); 472e78576d6SMark F. Adams 4730cbbd2e1SMark F. Adams if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0; 4740cbbd2e1SMark F. Adams else { 4759566063dSJacob Faibussowitsch PetscCall(PCGAMGHashTableFind(&fgid_flid, gid1, &flid)); 47608401ef6SPierre Jolivet PetscCheck(flid >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot find gid1 in table"); 4770cbbd2e1SMark F. Adams } 4782e68589bSMark F. Adams /* copy in B_i matrix - column oriented */ 47965d7b583SSatish Balay data = &data_in[flid * bs]; 4803d3eaba7SBarry Smith for (ii = 0; ii < bs; ii++) { 4812e68589bSMark F. Adams for (jj = 0; jj < N; jj++) { 4820cbbd2e1SMark F. Adams PetscReal d = data[jj * data_stride + ii]; 4830cbbd2e1SMark F. Adams qqc[jj * Mdata + aggID * bs + ii] = d; 4842e68589bSMark F. Adams } 4852e68589bSMark F. Adams } 4862e68589bSMark F. Adams /* set fine IDs */ 4872e68589bSMark F. Adams for (kk = 0; kk < bs; kk++) fids[aggID * bs + kk] = flid_fgid[flid] * bs + kk; 4882e68589bSMark F. Adams aggID++; 4890cbbd2e1SMark F. Adams } 4902e68589bSMark F. Adams 4912e68589bSMark F. Adams /* pad with zeros */ 4922e68589bSMark F. Adams for (ii = asz * bs; ii < Mdata; ii++) { 493ad540459SPierre Jolivet for (jj = 0; jj < N; jj++, kk++) qqc[jj * Mdata + ii] = .0; 4942e68589bSMark F. Adams } 4952e68589bSMark F. Adams 4962e68589bSMark F. Adams /* QR */ 4979566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 498792fecdfSBarry Smith PetscCallBLAS("LAPACKgeqrf", LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 4999566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPop()); 50008401ef6SPierre Jolivet PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xGEQRF error"); 5012e68589bSMark F. Adams /* get R - column oriented - output B_{i+1} */ 5022e68589bSMark F. Adams { 5032e68589bSMark F. Adams PetscReal *data = &out_data[clid * nSAvec]; 5042e68589bSMark F. Adams for (jj = 0; jj < nSAvec; jj++) { 5052e68589bSMark F. Adams for (ii = 0; ii < nSAvec; ii++) { 50608401ef6SPierre 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); 5070cbbd2e1SMark F. Adams if (ii <= jj) data[jj * out_data_stride + ii] = PetscRealPart(qqc[jj * Mdata + ii]); 5080cbbd2e1SMark F. Adams else data[jj * out_data_stride + ii] = 0.; 5092e68589bSMark F. Adams } 5102e68589bSMark F. Adams } 5112e68589bSMark F. Adams } 5122e68589bSMark F. Adams 5132e68589bSMark F. Adams /* get Q - row oriented */ 514792fecdfSBarry Smith PetscCallBLAS("LAPACKorgqr", LAPACKorgqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 51563a3b9bcSJacob Faibussowitsch PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xORGQR error arg %" PetscBLASInt_FMT, -INFO); 5162e68589bSMark F. Adams 5172e68589bSMark F. Adams for (ii = 0; ii < M; ii++) { 518ad540459SPierre Jolivet for (jj = 0; jj < N; jj++) qqr[N * ii + jj] = qqc[jj * Mdata + ii]; 5192e68589bSMark F. Adams } 5202e68589bSMark F. Adams 5212e68589bSMark F. Adams /* add diagonal block of P0 */ 5229371c9d4SSatish Balay for (kk = 0; kk < N; kk++) { cids[kk] = N * cgid + kk; /* global col IDs in P0 */ } 5239566063dSJacob Faibussowitsch PetscCall(MatSetValues(a_Prol, M, fids, N, cids, qqr, INSERT_VALUES)); 5249566063dSJacob Faibussowitsch PetscCall(PetscFree5(qqc, qqr, TAU, WORK, fids)); 525b43b03e9SMark F. Adams clid++; 5260cbbd2e1SMark F. Adams } /* coarse agg */ 5270cbbd2e1SMark F. Adams } /* for all fine nodes */ 5289566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(a_Prol, MAT_FINAL_ASSEMBLY)); 5299566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(a_Prol, MAT_FINAL_ASSEMBLY)); 5309566063dSJacob Faibussowitsch PetscCall(PCGAMGHashTableDestroy(&fgid_flid)); 5313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5322e68589bSMark F. Adams } 5332e68589bSMark F. Adams 534d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_GAMG_AGG(PC pc, PetscViewer viewer) 535d71ae5a4SJacob Faibussowitsch { 5365adeb434SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 5375adeb434SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 5385adeb434SBarry Smith PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 5395adeb434SBarry Smith 5405adeb434SBarry Smith PetscFunctionBegin; 5419566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " AGG specific options\n")); 542d529f056Smarkadams4 PetscCall(PetscViewerASCIIPrintf(viewer, " Number of levels of aggressive coarsening %d\n", (int)pc_gamg_agg->aggressive_coarsening_levels)); 543d529f056Smarkadams4 if (pc_gamg_agg->aggressive_coarsening_levels > 0) { 544d529f056Smarkadams4 PetscCall(PetscViewerASCIIPrintf(viewer, " %s aggressive coarsening\n", !pc_gamg_agg->use_aggressive_square_graph ? "MIS-k" : "Square graph")); 545d529f056Smarkadams4 if (!pc_gamg_agg->use_aggressive_square_graph) PetscCall(PetscViewerASCIIPrintf(viewer, " MIS-%d coarsening on aggressive levels\n", (int)pc_gamg_agg->aggressive_mis_k)); 546d529f056Smarkadams4 } 547bae903cbSmarkadams4 PetscCall(PetscViewerASCIIPrintf(viewer, " Number smoothing steps %d\n", (int)pc_gamg_agg->nsmooths)); 5483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5495adeb434SBarry Smith } 5505adeb434SBarry Smith 5512d776b49SBarry Smith static PetscErrorCode PCGAMGCreateGraph_AGG(PC pc, Mat Amat, Mat *a_Gmat) 552d71ae5a4SJacob Faibussowitsch { 553c8b0795cSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 554c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 555c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 5562d776b49SBarry Smith const PetscReal vfilter = pc_gamg->threshold[pc_gamg->current_level]; 5572d776b49SBarry Smith PetscBool ishem; 5582d776b49SBarry Smith const char *prefix; 559d529f056Smarkadams4 MatInfo info0, info1; 560d529f056Smarkadams4 PetscInt bs; 561c8b0795cSMark F. Adams 562c8b0795cSMark F. Adams PetscFunctionBegin; 563849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0)); 5642d776b49SBarry 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 */ 5652d776b49SBarry Smith 5662d776b49SBarry Smith /* MATCOARSENHEM requires numerical weights for edges so ensure they are computed */ 5672d776b49SBarry Smith PetscCall(MatCoarsenCreate(PetscObjectComm((PetscObject)pc), &pc_gamg_agg->crs)); 5682d776b49SBarry Smith PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, &prefix)); 5692d776b49SBarry Smith PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc_gamg_agg->crs, prefix)); 5702d776b49SBarry Smith PetscCall(MatCoarsenSetFromOptions(pc_gamg_agg->crs)); 5712d776b49SBarry Smith PetscCall(PetscObjectTypeCompare((PetscObject)pc_gamg_agg->crs, MATCOARSENHEM, &ishem)); 572d529f056Smarkadams4 if (ishem) pc_gamg_agg->aggressive_coarsening_levels = 0; // aggressive and HEM does not make sense 573d529f056Smarkadams4 PetscCall(MatGetInfo(Amat, MAT_LOCAL, &info0)); /* global reduction */ 5742d776b49SBarry Smith PetscCall(MatCreateGraph(Amat, PETSC_TRUE, (vfilter >= 0 || ishem) ? PETSC_TRUE : PETSC_FALSE, vfilter, a_Gmat)); 575d529f056Smarkadams4 PetscCall(MatGetInfo(*a_Gmat, MAT_LOCAL, &info1)); /* global reduction */ 576d529f056Smarkadams4 PetscCall(MatGetBlockSize(Amat, &bs)); 577d529f056Smarkadams4 if (info0.nz_used > 0) PetscCall(PetscInfo(pc, "Filtering left %g %% edges in graph (%e %e)\n", 100.0 * info1.nz_used * (double)(bs * bs) / info0.nz_used, info0.nz_used, info1.nz_used)); 578849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0)); 5793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 580c8b0795cSMark F. Adams } 581c8b0795cSMark F. Adams 582d529f056Smarkadams4 typedef PetscInt NState; 583d529f056Smarkadams4 static const NState NOT_DONE = -2; 584d529f056Smarkadams4 static const NState DELETED = -1; 585d529f056Smarkadams4 static const NState REMOVED = -3; 586d529f056Smarkadams4 #define IS_SELECTED(s) (s != DELETED && s != NOT_DONE && s != REMOVED) 587d529f056Smarkadams4 588d529f056Smarkadams4 /* 589d529f056Smarkadams4 fixAggregatesWithSquare - greedy grab of with G1 (unsquared graph) -- AIJ specific -- change to fixAggregatesWithSquare -- TODD 590d529f056Smarkadams4 - AGG-MG specific: clears singletons out of 'selected_2' 591d529f056Smarkadams4 592d529f056Smarkadams4 Input Parameter: 593d529f056Smarkadams4 . Gmat_2 - global matrix of squared graph (data not defined) 594d529f056Smarkadams4 . Gmat_1 - base graph to grab with base graph 595d529f056Smarkadams4 Input/Output Parameter: 596d529f056Smarkadams4 . aggs_2 - linked list of aggs with gids) 597d529f056Smarkadams4 */ 598d529f056Smarkadams4 static PetscErrorCode fixAggregatesWithSquare(PC pc, Mat Gmat_2, Mat Gmat_1, PetscCoarsenData *aggs_2) 599d529f056Smarkadams4 { 600d529f056Smarkadams4 PetscBool isMPI; 601d529f056Smarkadams4 Mat_SeqAIJ *matA_1, *matB_1 = NULL; 602d529f056Smarkadams4 MPI_Comm comm; 603d529f056Smarkadams4 PetscInt lid, *ii, *idx, ix, Iend, my0, kk, n, j; 604d529f056Smarkadams4 Mat_MPIAIJ *mpimat_2 = NULL, *mpimat_1 = NULL; 605d529f056Smarkadams4 const PetscInt nloc = Gmat_2->rmap->n; 606d529f056Smarkadams4 PetscScalar *cpcol_1_state, *cpcol_2_state, *cpcol_2_par_orig, *lid_parent_gid; 607d529f056Smarkadams4 PetscInt *lid_cprowID_1 = NULL; 608d529f056Smarkadams4 NState *lid_state; 609d529f056Smarkadams4 Vec ghost_par_orig2; 610d529f056Smarkadams4 PetscMPIInt rank; 611d529f056Smarkadams4 612d529f056Smarkadams4 PetscFunctionBegin; 613d529f056Smarkadams4 PetscCall(PetscObjectGetComm((PetscObject)Gmat_2, &comm)); 614d529f056Smarkadams4 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 615d529f056Smarkadams4 PetscCall(MatGetOwnershipRange(Gmat_1, &my0, &Iend)); 616d529f056Smarkadams4 617d529f056Smarkadams4 /* get submatrices */ 618d529f056Smarkadams4 PetscCall(PetscStrbeginswith(((PetscObject)Gmat_1)->type_name, MATMPIAIJ, &isMPI)); 619d529f056Smarkadams4 PetscCall(PetscInfo(pc, "isMPI = %s\n", isMPI ? "yes" : "no")); 620d529f056Smarkadams4 PetscCall(PetscMalloc3(nloc, &lid_state, nloc, &lid_parent_gid, nloc, &lid_cprowID_1)); 621d529f056Smarkadams4 for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1; 622d529f056Smarkadams4 if (isMPI) { 623d529f056Smarkadams4 /* grab matrix objects */ 624d529f056Smarkadams4 mpimat_2 = (Mat_MPIAIJ *)Gmat_2->data; 625d529f056Smarkadams4 mpimat_1 = (Mat_MPIAIJ *)Gmat_1->data; 626d529f056Smarkadams4 matA_1 = (Mat_SeqAIJ *)mpimat_1->A->data; 627d529f056Smarkadams4 matB_1 = (Mat_SeqAIJ *)mpimat_1->B->data; 628d529f056Smarkadams4 629d529f056Smarkadams4 /* force compressed row storage for B matrix in AuxMat */ 630d529f056Smarkadams4 PetscCall(MatCheckCompressedRow(mpimat_1->B, matB_1->nonzerorowcnt, &matB_1->compressedrow, matB_1->i, Gmat_1->rmap->n, -1.0)); 631d529f056Smarkadams4 for (ix = 0; ix < matB_1->compressedrow.nrows; ix++) { 632d529f056Smarkadams4 PetscInt lid = matB_1->compressedrow.rindex[ix]; 633d529f056Smarkadams4 PetscCheck(lid <= nloc && lid >= -1, PETSC_COMM_SELF, PETSC_ERR_USER, "lid %d out of range. nloc = %d", (int)lid, (int)nloc); 634d529f056Smarkadams4 if (lid != -1) lid_cprowID_1[lid] = ix; 635d529f056Smarkadams4 } 636d529f056Smarkadams4 } else { 637d529f056Smarkadams4 PetscBool isAIJ; 638d529f056Smarkadams4 PetscCall(PetscStrbeginswith(((PetscObject)Gmat_1)->type_name, MATSEQAIJ, &isAIJ)); 639d529f056Smarkadams4 PetscCheck(isAIJ, PETSC_COMM_SELF, PETSC_ERR_USER, "Require AIJ matrix."); 640d529f056Smarkadams4 matA_1 = (Mat_SeqAIJ *)Gmat_1->data; 641d529f056Smarkadams4 } 642d529f056Smarkadams4 if (nloc > 0) { PetscCheck(!matB_1 || matB_1->compressedrow.use, PETSC_COMM_SELF, PETSC_ERR_PLIB, "matB_1 && !matB_1->compressedrow.use: PETSc bug???"); } 643d529f056Smarkadams4 /* get state of locals and selected gid for deleted */ 644d529f056Smarkadams4 for (lid = 0; lid < nloc; lid++) { 645d529f056Smarkadams4 lid_parent_gid[lid] = -1.0; 646d529f056Smarkadams4 lid_state[lid] = DELETED; 647d529f056Smarkadams4 } 648d529f056Smarkadams4 649d529f056Smarkadams4 /* set lid_state */ 650d529f056Smarkadams4 for (lid = 0; lid < nloc; lid++) { 651d529f056Smarkadams4 PetscCDIntNd *pos; 652d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos)); 653d529f056Smarkadams4 if (pos) { 654d529f056Smarkadams4 PetscInt gid1; 655d529f056Smarkadams4 656d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gid1)); 657d529f056Smarkadams4 PetscCheck(gid1 == lid + my0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "gid1 %d != lid %d + my0 %d", (int)gid1, (int)lid, (int)my0); 658d529f056Smarkadams4 lid_state[lid] = gid1; 659d529f056Smarkadams4 } 660d529f056Smarkadams4 } 661d529f056Smarkadams4 662d529f056Smarkadams4 /* map local to selected local, DELETED means a ghost owns it */ 663d529f056Smarkadams4 for (lid = kk = 0; lid < nloc; lid++) { 664d529f056Smarkadams4 NState state = lid_state[lid]; 665d529f056Smarkadams4 if (IS_SELECTED(state)) { 666d529f056Smarkadams4 PetscCDIntNd *pos; 667d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos)); 668d529f056Smarkadams4 while (pos) { 669d529f056Smarkadams4 PetscInt gid1; 670d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gid1)); 671d529f056Smarkadams4 PetscCall(PetscCDGetNextPos(aggs_2, lid, &pos)); 672d529f056Smarkadams4 if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1 - my0] = (PetscScalar)(lid + my0); 673d529f056Smarkadams4 } 674d529f056Smarkadams4 } 675d529f056Smarkadams4 } 676d529f056Smarkadams4 /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */ 677d529f056Smarkadams4 if (isMPI) { 678d529f056Smarkadams4 Vec tempVec; 679d529f056Smarkadams4 /* get 'cpcol_1_state' */ 680d529f056Smarkadams4 PetscCall(MatCreateVecs(Gmat_1, &tempVec, NULL)); 681d529f056Smarkadams4 for (kk = 0, j = my0; kk < nloc; kk++, j++) { 682d529f056Smarkadams4 PetscScalar v = (PetscScalar)lid_state[kk]; 683d529f056Smarkadams4 PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES)); 684d529f056Smarkadams4 } 685d529f056Smarkadams4 PetscCall(VecAssemblyBegin(tempVec)); 686d529f056Smarkadams4 PetscCall(VecAssemblyEnd(tempVec)); 687d529f056Smarkadams4 PetscCall(VecScatterBegin(mpimat_1->Mvctx, tempVec, mpimat_1->lvec, INSERT_VALUES, SCATTER_FORWARD)); 688d529f056Smarkadams4 PetscCall(VecScatterEnd(mpimat_1->Mvctx, tempVec, mpimat_1->lvec, INSERT_VALUES, SCATTER_FORWARD)); 689d529f056Smarkadams4 PetscCall(VecGetArray(mpimat_1->lvec, &cpcol_1_state)); 690d529f056Smarkadams4 /* get 'cpcol_2_state' */ 691d529f056Smarkadams4 PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, mpimat_2->lvec, INSERT_VALUES, SCATTER_FORWARD)); 692d529f056Smarkadams4 PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, mpimat_2->lvec, INSERT_VALUES, SCATTER_FORWARD)); 693d529f056Smarkadams4 PetscCall(VecGetArray(mpimat_2->lvec, &cpcol_2_state)); 694d529f056Smarkadams4 /* get 'cpcol_2_par_orig' */ 695d529f056Smarkadams4 for (kk = 0, j = my0; kk < nloc; kk++, j++) { 696d529f056Smarkadams4 PetscScalar v = (PetscScalar)lid_parent_gid[kk]; 697d529f056Smarkadams4 PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES)); 698d529f056Smarkadams4 } 699d529f056Smarkadams4 PetscCall(VecAssemblyBegin(tempVec)); 700d529f056Smarkadams4 PetscCall(VecAssemblyEnd(tempVec)); 701d529f056Smarkadams4 PetscCall(VecDuplicate(mpimat_2->lvec, &ghost_par_orig2)); 702d529f056Smarkadams4 PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghost_par_orig2, INSERT_VALUES, SCATTER_FORWARD)); 703d529f056Smarkadams4 PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghost_par_orig2, INSERT_VALUES, SCATTER_FORWARD)); 704d529f056Smarkadams4 PetscCall(VecGetArray(ghost_par_orig2, &cpcol_2_par_orig)); 705d529f056Smarkadams4 706d529f056Smarkadams4 PetscCall(VecDestroy(&tempVec)); 707d529f056Smarkadams4 } /* ismpi */ 708d529f056Smarkadams4 for (lid = 0; lid < nloc; lid++) { 709d529f056Smarkadams4 NState state = lid_state[lid]; 710d529f056Smarkadams4 if (IS_SELECTED(state)) { 711d529f056Smarkadams4 /* steal locals */ 712d529f056Smarkadams4 ii = matA_1->i; 713d529f056Smarkadams4 n = ii[lid + 1] - ii[lid]; 714d529f056Smarkadams4 idx = matA_1->j + ii[lid]; 715d529f056Smarkadams4 for (j = 0; j < n; j++) { 716d529f056Smarkadams4 PetscInt lidj = idx[j], sgid; 717d529f056Smarkadams4 NState statej = lid_state[lidj]; 718d529f056Smarkadams4 if (statej == DELETED && (sgid = (PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid + my0) { /* steal local */ 719d529f056Smarkadams4 lid_parent_gid[lidj] = (PetscScalar)(lid + my0); /* send this if sgid is not local */ 720d529f056Smarkadams4 if (sgid >= my0 && sgid < Iend) { /* I'm stealing this local from a local sgid */ 721d529f056Smarkadams4 PetscInt hav = 0, slid = sgid - my0, gidj = lidj + my0; 722d529f056Smarkadams4 PetscCDIntNd *pos, *last = NULL; 723d529f056Smarkadams4 /* looking for local from local so id_llist_2 works */ 724d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, slid, &pos)); 725d529f056Smarkadams4 while (pos) { 726d529f056Smarkadams4 PetscInt gid; 727d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gid)); 728d529f056Smarkadams4 if (gid == gidj) { 729d529f056Smarkadams4 PetscCheck(last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "last cannot be null"); 730d529f056Smarkadams4 PetscCall(PetscCDRemoveNextNode(aggs_2, slid, last)); 731d529f056Smarkadams4 PetscCall(PetscCDAppendNode(aggs_2, lid, pos)); 732d529f056Smarkadams4 hav = 1; 733d529f056Smarkadams4 break; 734d529f056Smarkadams4 } else last = pos; 735d529f056Smarkadams4 PetscCall(PetscCDGetNextPos(aggs_2, slid, &pos)); 736d529f056Smarkadams4 } 737d529f056Smarkadams4 if (hav != 1) { 738d529f056Smarkadams4 PetscCheck(hav, PETSC_COMM_SELF, PETSC_ERR_PLIB, "failed to find adj in 'selected' lists - structurally unsymmetric matrix"); 739d529f056Smarkadams4 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "found node %d times???", (int)hav); 740d529f056Smarkadams4 } 741d529f056Smarkadams4 } else { /* I'm stealing this local, owned by a ghost */ 742d529f056Smarkadams4 PetscCheck(sgid == -1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Mat has an un-symmetric graph. Use '-%spc_gamg_sym_graph true' to symmetrize the graph or '-%spc_gamg_threshold -1' if the matrix is structurally symmetric.", 743d529f056Smarkadams4 ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : ""); 744d529f056Smarkadams4 PetscCall(PetscCDAppendID(aggs_2, lid, lidj + my0)); 745d529f056Smarkadams4 } 746d529f056Smarkadams4 } 747d529f056Smarkadams4 } /* local neighbors */ 748d529f056Smarkadams4 } else if (state == DELETED /* && lid_cprowID_1 */) { 749d529f056Smarkadams4 PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]); 750d529f056Smarkadams4 /* see if I have a selected ghost neighbor that will steal me */ 751d529f056Smarkadams4 if ((ix = lid_cprowID_1[lid]) != -1) { 752d529f056Smarkadams4 ii = matB_1->compressedrow.i; 753d529f056Smarkadams4 n = ii[ix + 1] - ii[ix]; 754d529f056Smarkadams4 idx = matB_1->j + ii[ix]; 755d529f056Smarkadams4 for (j = 0; j < n; j++) { 756d529f056Smarkadams4 PetscInt cpid = idx[j]; 757d529f056Smarkadams4 NState statej = (NState)PetscRealPart(cpcol_1_state[cpid]); 758d529f056Smarkadams4 if (IS_SELECTED(statej) && sgidold != (PetscInt)statej) { /* ghost will steal this, remove from my list */ 759d529f056Smarkadams4 lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */ 760d529f056Smarkadams4 if (sgidold >= my0 && sgidold < Iend) { /* this was mine */ 761d529f056Smarkadams4 PetscInt hav = 0, oldslidj = sgidold - my0; 762d529f056Smarkadams4 PetscCDIntNd *pos, *last = NULL; 763d529f056Smarkadams4 /* remove from 'oldslidj' list */ 764d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, oldslidj, &pos)); 765d529f056Smarkadams4 while (pos) { 766d529f056Smarkadams4 PetscInt gid; 767d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gid)); 768d529f056Smarkadams4 if (lid + my0 == gid) { 769d529f056Smarkadams4 /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove lid from oldslidj list *\/ */ 770d529f056Smarkadams4 PetscCheck(last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "last cannot be null"); 771d529f056Smarkadams4 PetscCall(PetscCDRemoveNextNode(aggs_2, oldslidj, last)); 772d529f056Smarkadams4 /* ghost (PetscScalar)statej will add this later */ 773d529f056Smarkadams4 hav = 1; 774d529f056Smarkadams4 break; 775d529f056Smarkadams4 } else last = pos; 776d529f056Smarkadams4 PetscCall(PetscCDGetNextPos(aggs_2, oldslidj, &pos)); 777d529f056Smarkadams4 } 778d529f056Smarkadams4 if (hav != 1) { 779d529f056Smarkadams4 PetscCheck(hav, PETSC_COMM_SELF, PETSC_ERR_PLIB, "failed to find (hav=%d) adj in 'selected' lists - structurally unsymmetric matrix", (int)hav); 780d529f056Smarkadams4 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "found node %d times???", (int)hav); 781d529f056Smarkadams4 } 782d529f056Smarkadams4 } else { 783d529f056Smarkadams4 /* TODO: ghosts remove this later */ 784d529f056Smarkadams4 } 785d529f056Smarkadams4 } 786d529f056Smarkadams4 } 787d529f056Smarkadams4 } 788d529f056Smarkadams4 } /* selected/deleted */ 789d529f056Smarkadams4 } /* node loop */ 790d529f056Smarkadams4 791d529f056Smarkadams4 if (isMPI) { 792d529f056Smarkadams4 PetscScalar *cpcol_2_parent, *cpcol_2_gid; 793d529f056Smarkadams4 Vec tempVec, ghostgids2, ghostparents2; 794d529f056Smarkadams4 PetscInt cpid, nghost_2; 795d529f056Smarkadams4 PCGAMGHashTable gid_cpid; 796d529f056Smarkadams4 797d529f056Smarkadams4 PetscCall(VecGetSize(mpimat_2->lvec, &nghost_2)); 798d529f056Smarkadams4 PetscCall(MatCreateVecs(Gmat_2, &tempVec, NULL)); 799d529f056Smarkadams4 800d529f056Smarkadams4 /* get 'cpcol_2_parent' */ 801d529f056Smarkadams4 for (kk = 0, j = my0; kk < nloc; kk++, j++) { PetscCall(VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES)); } 802d529f056Smarkadams4 PetscCall(VecAssemblyBegin(tempVec)); 803d529f056Smarkadams4 PetscCall(VecAssemblyEnd(tempVec)); 804d529f056Smarkadams4 PetscCall(VecDuplicate(mpimat_2->lvec, &ghostparents2)); 805d529f056Smarkadams4 PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghostparents2, INSERT_VALUES, SCATTER_FORWARD)); 806d529f056Smarkadams4 PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghostparents2, INSERT_VALUES, SCATTER_FORWARD)); 807d529f056Smarkadams4 PetscCall(VecGetArray(ghostparents2, &cpcol_2_parent)); 808d529f056Smarkadams4 809d529f056Smarkadams4 /* get 'cpcol_2_gid' */ 810d529f056Smarkadams4 for (kk = 0, j = my0; kk < nloc; kk++, j++) { 811d529f056Smarkadams4 PetscScalar v = (PetscScalar)j; 812d529f056Smarkadams4 PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES)); 813d529f056Smarkadams4 } 814d529f056Smarkadams4 PetscCall(VecAssemblyBegin(tempVec)); 815d529f056Smarkadams4 PetscCall(VecAssemblyEnd(tempVec)); 816d529f056Smarkadams4 PetscCall(VecDuplicate(mpimat_2->lvec, &ghostgids2)); 817d529f056Smarkadams4 PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghostgids2, INSERT_VALUES, SCATTER_FORWARD)); 818d529f056Smarkadams4 PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghostgids2, INSERT_VALUES, SCATTER_FORWARD)); 819d529f056Smarkadams4 PetscCall(VecGetArray(ghostgids2, &cpcol_2_gid)); 820d529f056Smarkadams4 PetscCall(VecDestroy(&tempVec)); 821d529f056Smarkadams4 822d529f056Smarkadams4 /* look for deleted ghosts and add to table */ 823d529f056Smarkadams4 PetscCall(PCGAMGHashTableCreate(2 * nghost_2 + 1, &gid_cpid)); 824d529f056Smarkadams4 for (cpid = 0; cpid < nghost_2; cpid++) { 825d529f056Smarkadams4 NState state = (NState)PetscRealPart(cpcol_2_state[cpid]); 826d529f056Smarkadams4 if (state == DELETED) { 827d529f056Smarkadams4 PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 828d529f056Smarkadams4 PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]); 829d529f056Smarkadams4 if (sgid_old == -1 && sgid_new != -1) { 830d529f056Smarkadams4 PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 831d529f056Smarkadams4 PetscCall(PCGAMGHashTableAdd(&gid_cpid, gid, cpid)); 832d529f056Smarkadams4 } 833d529f056Smarkadams4 } 834d529f056Smarkadams4 } 835d529f056Smarkadams4 836d529f056Smarkadams4 /* look for deleted ghosts and see if they moved - remove it */ 837d529f056Smarkadams4 for (lid = 0; lid < nloc; lid++) { 838d529f056Smarkadams4 NState state = lid_state[lid]; 839d529f056Smarkadams4 if (IS_SELECTED(state)) { 840d529f056Smarkadams4 PetscCDIntNd *pos, *last = NULL; 841d529f056Smarkadams4 /* look for deleted ghosts and see if they moved */ 842d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos)); 843d529f056Smarkadams4 while (pos) { 844d529f056Smarkadams4 PetscInt gid; 845d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gid)); 846d529f056Smarkadams4 847d529f056Smarkadams4 if (gid < my0 || gid >= Iend) { 848d529f056Smarkadams4 PetscCall(PCGAMGHashTableFind(&gid_cpid, gid, &cpid)); 849d529f056Smarkadams4 if (cpid != -1) { 850d529f056Smarkadams4 /* a moved ghost - */ 851d529f056Smarkadams4 /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove 'flid' from list *\/ */ 852d529f056Smarkadams4 PetscCall(PetscCDRemoveNextNode(aggs_2, lid, last)); 853d529f056Smarkadams4 } else last = pos; 854d529f056Smarkadams4 } else last = pos; 855d529f056Smarkadams4 856d529f056Smarkadams4 PetscCall(PetscCDGetNextPos(aggs_2, lid, &pos)); 857d529f056Smarkadams4 } /* loop over list of deleted */ 858d529f056Smarkadams4 } /* selected */ 859d529f056Smarkadams4 } 860d529f056Smarkadams4 PetscCall(PCGAMGHashTableDestroy(&gid_cpid)); 861d529f056Smarkadams4 862d529f056Smarkadams4 /* look at ghosts, see if they changed - and it */ 863d529f056Smarkadams4 for (cpid = 0; cpid < nghost_2; cpid++) { 864d529f056Smarkadams4 PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 865d529f056Smarkadams4 if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */ 866d529f056Smarkadams4 PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 867d529f056Smarkadams4 PetscInt slid_new = sgid_new - my0, hav = 0; 868d529f056Smarkadams4 PetscCDIntNd *pos; 869d529f056Smarkadams4 870d529f056Smarkadams4 /* search for this gid to see if I have it */ 871d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, slid_new, &pos)); 872d529f056Smarkadams4 while (pos) { 873d529f056Smarkadams4 PetscInt gidj; 874d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gidj)); 875d529f056Smarkadams4 PetscCall(PetscCDGetNextPos(aggs_2, slid_new, &pos)); 876d529f056Smarkadams4 877d529f056Smarkadams4 if (gidj == gid) { 878d529f056Smarkadams4 hav = 1; 879d529f056Smarkadams4 break; 880d529f056Smarkadams4 } 881d529f056Smarkadams4 } 882d529f056Smarkadams4 if (hav != 1) { 883d529f056Smarkadams4 /* insert 'flidj' into head of llist */ 884d529f056Smarkadams4 PetscCall(PetscCDAppendID(aggs_2, slid_new, gid)); 885d529f056Smarkadams4 } 886d529f056Smarkadams4 } 887d529f056Smarkadams4 } 888d529f056Smarkadams4 PetscCall(VecRestoreArray(mpimat_1->lvec, &cpcol_1_state)); 889d529f056Smarkadams4 PetscCall(VecRestoreArray(mpimat_2->lvec, &cpcol_2_state)); 890d529f056Smarkadams4 PetscCall(VecRestoreArray(ghostparents2, &cpcol_2_parent)); 891d529f056Smarkadams4 PetscCall(VecRestoreArray(ghostgids2, &cpcol_2_gid)); 892d529f056Smarkadams4 PetscCall(VecDestroy(&ghostgids2)); 893d529f056Smarkadams4 PetscCall(VecDestroy(&ghostparents2)); 894d529f056Smarkadams4 PetscCall(VecDestroy(&ghost_par_orig2)); 895d529f056Smarkadams4 } 896d529f056Smarkadams4 PetscCall(PetscFree3(lid_state, lid_parent_gid, lid_cprowID_1)); 897d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 898d529f056Smarkadams4 } 899d529f056Smarkadams4 900c8b0795cSMark F. Adams /* 901bae903cbSmarkadams4 PCGAMGCoarsen_AGG - supports squaring the graph (deprecated) and new graph for 902bae903cbSmarkadams4 communication of QR data used with HEM and MISk coarsening 903c8b0795cSMark F. Adams 904c8b0795cSMark F. Adams Input Parameter: 905e0940f08SMark F. Adams . a_pc - this 906bae903cbSmarkadams4 907e0940f08SMark F. Adams Input/Output Parameter: 908bae903cbSmarkadams4 . a_Gmat1 - graph to coarsen (in), graph off processor edges for QR gather scatter (out) 909bae903cbSmarkadams4 910c8b0795cSMark F. Adams Output Parameter: 9110cbbd2e1SMark F. Adams . agg_lists - list of aggregates 912bae903cbSmarkadams4 913c8b0795cSMark F. Adams */ 914d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc, Mat *a_Gmat1, PetscCoarsenData **agg_lists) 915d71ae5a4SJacob Faibussowitsch { 916e0940f08SMark F. Adams PC_MG *mg = (PC_MG *)a_pc->data; 917c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 918c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 919d529f056Smarkadams4 Mat mat, Gmat2, Gmat1 = *a_Gmat1; /* aggressive graph */ 9200cbbd2e1SMark F. Adams IS perm; 921bae903cbSmarkadams4 PetscInt Istart, Iend, Ii, nloc, bs, nn; 922bae903cbSmarkadams4 PetscInt *permute, *degree; 923c8b0795cSMark F. Adams PetscBool *bIndexSet; 924aea53286SMark Adams PetscReal hashfact; 925e2c00dcbSBarry Smith PetscInt iSwapIndex; 9263b50add6SMark Adams PetscRandom random; 927d529f056Smarkadams4 MPI_Comm comm; 928c8b0795cSMark F. Adams 929c8b0795cSMark F. Adams PetscFunctionBegin; 930d529f056Smarkadams4 PetscCall(PetscObjectGetComm((PetscObject)Gmat1, &comm)); 931849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0)); 932bae903cbSmarkadams4 PetscCall(MatGetLocalSize(Gmat1, &nn, NULL)); 9339566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Gmat1, &bs)); 93463a3b9bcSJacob Faibussowitsch PetscCheck(bs == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bs %" PetscInt_FMT " must be 1", bs); 935bae903cbSmarkadams4 nloc = nn / bs; 9365cfd4bd9SMark Adams /* get MIS aggs - randomize */ 937bae903cbSmarkadams4 PetscCall(PetscMalloc2(nloc, &permute, nloc, °ree)); 9389566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nloc, &bIndexSet)); 9396e09b0e3SStefano Zampini for (Ii = 0; Ii < nloc; Ii++) permute[Ii] = Ii; 9409566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &random)); 9419566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Gmat1, &Istart, &Iend)); 942c8b0795cSMark F. Adams for (Ii = 0; Ii < nloc; Ii++) { 943bae903cbSmarkadams4 PetscInt nc; 944bae903cbSmarkadams4 PetscCall(MatGetRow(Gmat1, Istart + Ii, &nc, NULL, NULL)); 945bae903cbSmarkadams4 degree[Ii] = nc; 946bae903cbSmarkadams4 PetscCall(MatRestoreRow(Gmat1, Istart + Ii, &nc, NULL, NULL)); 947bae903cbSmarkadams4 } 948bae903cbSmarkadams4 for (Ii = 0; Ii < nloc; Ii++) { 9499566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(random, &hashfact)); 950aea53286SMark Adams iSwapIndex = (PetscInt)(hashfact * nloc) % nloc; 951c8b0795cSMark F. Adams if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) { 952c8b0795cSMark F. Adams PetscInt iTemp = permute[iSwapIndex]; 953c8b0795cSMark F. Adams permute[iSwapIndex] = permute[Ii]; 954c8b0795cSMark F. Adams permute[Ii] = iTemp; 955bae903cbSmarkadams4 iTemp = degree[iSwapIndex]; 956bae903cbSmarkadams4 degree[iSwapIndex] = degree[Ii]; 957bae903cbSmarkadams4 degree[Ii] = iTemp; 958c8b0795cSMark F. Adams bIndexSet[iSwapIndex] = PETSC_TRUE; 959c8b0795cSMark F. Adams } 960c8b0795cSMark F. Adams } 961d529f056Smarkadams4 // apply minimum degree ordering -- NEW 962d529f056Smarkadams4 if (pc_gamg_agg->use_minimum_degree_ordering) { PetscCall(PetscSortIntWithArray(nloc, degree, permute)); } 9639566063dSJacob Faibussowitsch PetscCall(PetscFree(bIndexSet)); 9649566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&random)); 9659566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm)); 966849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0)); 967d529f056Smarkadams4 // square graph 968d529f056Smarkadams4 if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels && pc_gamg_agg->use_aggressive_square_graph) { 969d529f056Smarkadams4 PetscCall(PCGAMGSquareGraph_GAMG(a_pc, Gmat1, &Gmat2)); 970d529f056Smarkadams4 } else Gmat2 = Gmat1; 971d529f056Smarkadams4 // switch to old MIS-1 for square graph 972d529f056Smarkadams4 if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels) { 973d529f056Smarkadams4 if (!pc_gamg_agg->use_aggressive_square_graph) PetscCall(MatCoarsenMISKSetDistance(pc_gamg_agg->crs, pc_gamg_agg->aggressive_mis_k)); // hardwire to MIS-2 974d529f056Smarkadams4 else PetscCall(MatCoarsenSetType(pc_gamg_agg->crs, MATCOARSENMIS)); // old MIS -- side effect 975d529f056Smarkadams4 } else if (pc_gamg_agg->use_aggressive_square_graph && pc_gamg_agg->aggressive_coarsening_levels > 0) { // we reset the MIS 976d529f056Smarkadams4 const char *prefix; 977d529f056Smarkadams4 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)a_pc, &prefix)); 978d529f056Smarkadams4 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc_gamg_agg->crs, prefix)); 979d529f056Smarkadams4 PetscCall(MatCoarsenSetFromOptions(pc_gamg_agg->crs)); // get the default back on non-aggressive levels when square graph switched to old MIS 980d529f056Smarkadams4 } 981d529f056Smarkadams4 PetscCall(MatCoarsenSetAdjacency(pc_gamg_agg->crs, Gmat2)); 9822d776b49SBarry Smith PetscCall(MatCoarsenSetStrictAggs(pc_gamg_agg->crs, PETSC_TRUE)); 983d529f056Smarkadams4 PetscCall(MatCoarsenSetGreedyOrdering(pc_gamg_agg->crs, perm)); 9842d776b49SBarry Smith PetscCall(MatCoarsenApply(pc_gamg_agg->crs)); 9852d776b49SBarry Smith PetscCall(MatCoarsenViewFromOptions(pc_gamg_agg->crs, NULL, "-mat_coarsen_view")); 9862d776b49SBarry Smith PetscCall(MatCoarsenGetData(pc_gamg_agg->crs, agg_lists)); /* output */ 9872d776b49SBarry Smith PetscCall(MatCoarsenDestroy(&pc_gamg_agg->crs)); 988b43b03e9SMark F. Adams 9899566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 990bae903cbSmarkadams4 PetscCall(PetscFree2(permute, degree)); 991849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0)); 9929f3f12c8SMark F. Adams 993d529f056Smarkadams4 if (Gmat2 != Gmat1) { 994d529f056Smarkadams4 PetscCoarsenData *llist = *agg_lists; 995d529f056Smarkadams4 PetscCall(fixAggregatesWithSquare(a_pc, Gmat2, Gmat1, *agg_lists)); 996d529f056Smarkadams4 PetscCall(MatDestroy(&Gmat1)); 997d529f056Smarkadams4 *a_Gmat1 = Gmat2; /* output */ 998d529f056Smarkadams4 PetscCall(PetscCDGetMat(llist, &mat)); 999baca6076SPierre Jolivet PetscCheck(!mat, comm, PETSC_ERR_ARG_WRONG, "Unexpected auxiliary matrix with squared graph"); 1000d529f056Smarkadams4 } else { 1001bae903cbSmarkadams4 PetscCoarsenData *llist = *agg_lists; 10029ab59c8bSMark Adams /* see if we have a matrix that takes precedence (returned from MatCoarsenApply) */ 10039566063dSJacob Faibussowitsch PetscCall(PetscCDGetMat(llist, &mat)); 10040cbbd2e1SMark F. Adams if (mat) { 1005d529f056Smarkadams4 PetscCall(MatDestroy(a_Gmat1)); 10060cbbd2e1SMark F. Adams *a_Gmat1 = mat; /* output */ 10070cbbd2e1SMark F. Adams } 10080cbbd2e1SMark F. Adams } 1009849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0)); 10103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1011c8b0795cSMark F. Adams } 1012c8b0795cSMark F. Adams 1013c8b0795cSMark F. Adams /* 10140cbbd2e1SMark F. Adams PCGAMGProlongator_AGG 1015c8b0795cSMark F. Adams 1016c8b0795cSMark F. Adams Input Parameter: 1017c8b0795cSMark F. Adams . pc - this 1018c8b0795cSMark F. Adams . Amat - matrix on this fine level 1019c8b0795cSMark F. Adams . Graph - used to get ghost data for nodes in 10200cbbd2e1SMark F. Adams . agg_lists - list of aggregates 1021c8b0795cSMark F. Adams Output Parameter: 1022c8b0795cSMark F. Adams . a_P_out - prolongation operator to the next level 1023c8b0795cSMark F. Adams */ 1024d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGProlongator_AGG(PC pc, Mat Amat, Mat Gmat, PetscCoarsenData *agg_lists, Mat *a_P_out) 1025d71ae5a4SJacob Faibussowitsch { 10262e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 10272e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 10284a2f8832SBarry Smith const PetscInt col_bs = pc_gamg->data_cell_cols; 1029c8b0795cSMark F. Adams PetscInt Istart, Iend, nloc, ii, jj, kk, my0, nLocalSelected, bs; 1030c8b0795cSMark F. Adams Mat Prol; 1031d5d25401SBarry Smith PetscMPIInt size; 10323b4367a7SBarry Smith MPI_Comm comm; 1033c8b0795cSMark F. Adams PetscReal *data_w_ghost; 1034c8b0795cSMark F. Adams PetscInt myCrs0, nbnodes = 0, *flid_fgid; 1035d9558ea9SBarry Smith MatType mtype; 10362e68589bSMark F. Adams 10372e68589bSMark F. Adams PetscFunctionBegin; 10389566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm)); 103908401ef6SPierre Jolivet PetscCheck(col_bs >= 1, comm, PETSC_ERR_PLIB, "Column bs cannot be less than 1"); 1040849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0)); 10419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 10429566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Amat, &Istart, &Iend)); 10439566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Amat, &bs)); 10449371c9d4SSatish Balay nloc = (Iend - Istart) / bs; 10459371c9d4SSatish Balay my0 = Istart / bs; 104663a3b9bcSJacob 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); 10472e68589bSMark F. Adams 10482e68589bSMark F. Adams /* get 'nLocalSelected' */ 10490cbbd2e1SMark F. Adams for (ii = 0, nLocalSelected = 0; ii < nloc; ii++) { 1050e78576d6SMark F. Adams PetscBool ise; 10510cbbd2e1SMark F. Adams /* filter out singletons 0 or 1? */ 10529566063dSJacob Faibussowitsch PetscCall(PetscCDEmptyAt(agg_lists, ii, &ise)); 1053e78576d6SMark F. Adams if (!ise) nLocalSelected++; 10542e68589bSMark F. Adams } 10552e68589bSMark F. Adams 10562e68589bSMark F. Adams /* create prolongator, create P matrix */ 10579566063dSJacob Faibussowitsch PetscCall(MatGetType(Amat, &mtype)); 10589566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &Prol)); 10599566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Prol, nloc * bs, nLocalSelected * col_bs, PETSC_DETERMINE, PETSC_DETERMINE)); 10609566063dSJacob Faibussowitsch PetscCall(MatSetBlockSizes(Prol, bs, col_bs)); 10619566063dSJacob Faibussowitsch PetscCall(MatSetType(Prol, mtype)); 10623742c8caSstefanozampini #if PetscDefined(HAVE_DEVICE) 10633742c8caSstefanozampini PetscBool flg; 10643742c8caSstefanozampini PetscCall(MatBoundToCPU(Amat, &flg)); 10653742c8caSstefanozampini PetscCall(MatBindToCPU(Prol, flg)); 10663742c8caSstefanozampini if (flg) PetscCall(MatSetBindingPropagates(Prol, PETSC_TRUE)); 10673742c8caSstefanozampini #endif 10689566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(Prol, col_bs, NULL)); 10699566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(Prol, col_bs, NULL, col_bs, NULL)); 10702e68589bSMark F. Adams 10712e68589bSMark F. Adams /* can get all points "removed" */ 10729566063dSJacob Faibussowitsch PetscCall(MatGetSize(Prol, &kk, &ii)); 10737f66b68fSMark Adams if (!ii) { 107463a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: No selected points on coarse grid\n", ((PetscObject)pc)->prefix)); 10759566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Prol)); 10760298fd71SBarry Smith *a_P_out = NULL; /* out */ 1077849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0)); 10783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10792e68589bSMark F. Adams } 108063a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: New grid %" PetscInt_FMT " nodes\n", ((PetscObject)pc)->prefix, ii / col_bs)); 10819566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk)); 10820cbbd2e1SMark F. Adams 108363a3b9bcSJacob 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); 1084c8b0795cSMark F. Adams myCrs0 = myCrs0 / col_bs; 108563a3b9bcSJacob 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); 10862e68589bSMark F. Adams 10872e68589bSMark F. Adams /* create global vector of data in 'data_w_ghost' */ 1088849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0)); 1089bae903cbSmarkadams4 if (size > 1) { /* get ghost null space data */ 10902e68589bSMark F. Adams PetscReal *tmp_gdata, *tmp_ldata, *tp2; 10919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc, &tmp_ldata)); 10924a2f8832SBarry Smith for (jj = 0; jj < col_bs; jj++) { 1093c8b0795cSMark F. Adams for (kk = 0; kk < bs; kk++) { 1094a2f3521dSMark F. Adams PetscInt ii, stride; 1095c8b0795cSMark F. Adams const PetscReal *tp = pc_gamg->data + jj * bs * nloc + kk; 10962fa5cd67SKarl Rupp for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp; 10972fa5cd67SKarl Rupp 10989566063dSJacob Faibussowitsch PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata)); 1099a2f3521dSMark F. Adams 1100bae903cbSmarkadams4 if (!jj && !kk) { /* now I know how many total nodes - allocate TODO: move below and do in one 'col_bs' call */ 11019566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(stride * bs * col_bs, &data_w_ghost)); 1102a2f3521dSMark F. Adams nbnodes = bs * stride; 11032e68589bSMark F. Adams } 1104a2f3521dSMark F. Adams tp2 = data_w_ghost + jj * bs * stride + kk; 11052fa5cd67SKarl Rupp for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii]; 11069566063dSJacob Faibussowitsch PetscCall(PetscFree(tmp_gdata)); 11072e68589bSMark F. Adams } 11082e68589bSMark F. Adams } 11099566063dSJacob Faibussowitsch PetscCall(PetscFree(tmp_ldata)); 1110806fa848SBarry Smith } else { 1111c8b0795cSMark F. Adams nbnodes = bs * nloc; 1112c8b0795cSMark F. Adams data_w_ghost = (PetscReal *)pc_gamg->data; 11132e68589bSMark F. Adams } 11142e68589bSMark F. Adams 1115bae903cbSmarkadams4 /* get 'flid_fgid' TODO - move up to get 'stride' and do get null space data above in one step (jj loop) */ 1116c5df96a5SBarry Smith if (size > 1) { 11172e68589bSMark F. Adams PetscReal *fid_glid_loc, *fiddata; 1118a2f3521dSMark F. Adams PetscInt stride; 11192e68589bSMark F. Adams 11209566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc, &fid_glid_loc)); 11212e68589bSMark F. Adams for (kk = 0; kk < nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0 + kk); 11229566063dSJacob Faibussowitsch PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata)); 1123bae903cbSmarkadams4 PetscCall(PetscMalloc1(stride, &flid_fgid)); /* copy real data to in */ 1124a2f3521dSMark F. Adams for (kk = 0; kk < stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk]; 11259566063dSJacob Faibussowitsch PetscCall(PetscFree(fiddata)); 1126a2f3521dSMark F. Adams 112763a3b9bcSJacob Faibussowitsch PetscCheck(stride == nbnodes / bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "stride %" PetscInt_FMT " != nbnodes %" PetscInt_FMT "/bs %" PetscInt_FMT, stride, nbnodes, bs); 11289566063dSJacob Faibussowitsch PetscCall(PetscFree(fid_glid_loc)); 1129806fa848SBarry Smith } else { 11309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc, &flid_fgid)); 11312e68589bSMark F. Adams for (kk = 0; kk < nloc; kk++) flid_fgid[kk] = my0 + kk; 11322e68589bSMark F. Adams } 1133849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0)); 1134bae903cbSmarkadams4 /* get P0 */ 1135849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0)); 1136c8b0795cSMark F. Adams { 11370298fd71SBarry Smith PetscReal *data_out = NULL; 11389566063dSJacob Faibussowitsch PetscCall(formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes, data_w_ghost, flid_fgid, &data_out, Prol)); 11399566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 11402fa5cd67SKarl Rupp 1141c8b0795cSMark F. Adams pc_gamg->data = data_out; 11424a2f8832SBarry Smith pc_gamg->data_cell_rows = col_bs; 11434a2f8832SBarry Smith pc_gamg->data_sz = col_bs * col_bs * nLocalSelected; 1144c8b0795cSMark F. Adams } 1145849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0)); 114648a46eb9SPierre Jolivet if (size > 1) PetscCall(PetscFree(data_w_ghost)); 11479566063dSJacob Faibussowitsch PetscCall(PetscFree(flid_fgid)); 11482e68589bSMark F. Adams 1149c8b0795cSMark F. Adams *a_P_out = Prol; /* out */ 1150ed4e82eaSPeter Brune 1151849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0)); 11523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1153c8b0795cSMark F. Adams } 1154c8b0795cSMark F. Adams 1155c8b0795cSMark F. Adams /* 1156fd1112cbSBarry Smith PCGAMGOptProlongator_AGG 1157c8b0795cSMark F. Adams 1158c8b0795cSMark F. Adams Input Parameter: 1159c8b0795cSMark F. Adams . pc - this 1160c8b0795cSMark F. Adams . Amat - matrix on this fine level 1161c8b0795cSMark F. Adams In/Output Parameter: 11622a808120SBarry Smith . a_P - prolongation operator to the next level 1163c8b0795cSMark F. Adams */ 1164d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGOptProlongator_AGG(PC pc, Mat Amat, Mat *a_P) 1165d71ae5a4SJacob Faibussowitsch { 1166c8b0795cSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 1167c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1168c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 1169c8b0795cSMark F. Adams PetscInt jj; 1170c8b0795cSMark F. Adams Mat Prol = *a_P; 11713b4367a7SBarry Smith MPI_Comm comm; 11722a808120SBarry Smith KSP eksp; 11732a808120SBarry Smith Vec bb, xx; 11742a808120SBarry Smith PC epc; 11752a808120SBarry Smith PetscReal alpha, emax, emin; 1176c8b0795cSMark F. Adams 1177c8b0795cSMark F. Adams PetscFunctionBegin; 11789566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm)); 1179849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0)); 1180c8b0795cSMark F. Adams 1181d5d25401SBarry Smith /* compute maximum singular value of operator to be used in smoother */ 11822a808120SBarry Smith if (0 < pc_gamg_agg->nsmooths) { 118318c3aa7eSMark /* get eigen estimates */ 118418c3aa7eSMark if (pc_gamg->emax > 0) { 118518c3aa7eSMark emin = pc_gamg->emin; 118618c3aa7eSMark emax = pc_gamg->emax; 118718c3aa7eSMark } else { 11880ed2132dSStefano Zampini const char *prefix; 11890ed2132dSStefano Zampini 11909566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Amat, &bb, NULL)); 11919566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Amat, &xx, NULL)); 11929566063dSJacob Faibussowitsch PetscCall(KSPSetNoisy_Private(bb)); 1193e696ed0bSMark F. Adams 11949566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm, &eksp)); 11953821be0aSBarry Smith PetscCall(KSPSetNestLevel(eksp, pc->kspnestlevel)); 11969566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 11979566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(eksp, prefix)); 11989566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(eksp, "pc_gamg_esteig_")); 119973f7197eSJed Brown { 1200b94d7dedSBarry Smith PetscBool isset, sflg; 1201b94d7dedSBarry Smith PetscCall(MatIsSPDKnown(Amat, &isset, &sflg)); 1202b94d7dedSBarry Smith if (isset && sflg) PetscCall(KSPSetType(eksp, KSPCG)); 1203d24ecf33SMark } 12049566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(eksp, pc->erroriffailure)); 12059566063dSJacob Faibussowitsch PetscCall(KSPSetNormType(eksp, KSP_NORM_NONE)); 1206c2436cacSMark F. Adams 12079566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(eksp, PETSC_FALSE)); 12089566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(eksp, Amat, Amat)); 12092e68589bSMark F. Adams 12109566063dSJacob Faibussowitsch PetscCall(KSPGetPC(eksp, &epc)); 12119566063dSJacob Faibussowitsch PetscCall(PCSetType(epc, PCJACOBI)); /* smoother in smoothed agg. */ 1212c2436cacSMark F. Adams 12139566063dSJacob 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 1214074aec55SMark Adams 12159566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(eksp)); 12169566063dSJacob Faibussowitsch PetscCall(KSPSetComputeSingularValues(eksp, PETSC_TRUE)); 12179566063dSJacob Faibussowitsch PetscCall(KSPSolve(eksp, bb, xx)); 12189566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(eksp, pc, xx)); 12192e68589bSMark F. Adams 12209566063dSJacob Faibussowitsch PetscCall(KSPComputeExtremeSingularValues(eksp, &emax, &emin)); 12219566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: Smooth P0: max eigen=%e min=%e PC=%s\n", ((PetscObject)pc)->prefix, (double)emax, (double)emin, PCJACOBI)); 12229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xx)); 12239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bb)); 12249566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&eksp)); 12252e68589bSMark F. Adams } 122618c3aa7eSMark if (pc_gamg->use_sa_esteig) { 122718c3aa7eSMark mg->min_eigen_DinvA[pc_gamg->current_level] = emin; 122818c3aa7eSMark mg->max_eigen_DinvA[pc_gamg->current_level] = emax; 122963a3b9bcSJacob 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)); 123018c3aa7eSMark } else { 123118c3aa7eSMark mg->min_eigen_DinvA[pc_gamg->current_level] = 0; 123218c3aa7eSMark mg->max_eigen_DinvA[pc_gamg->current_level] = 0; 123318c3aa7eSMark } 123418c3aa7eSMark } else { 123518c3aa7eSMark mg->min_eigen_DinvA[pc_gamg->current_level] = 0; 123618c3aa7eSMark mg->max_eigen_DinvA[pc_gamg->current_level] = 0; 123718c3aa7eSMark } 12382e68589bSMark F. Adams 12392a808120SBarry Smith /* smooth P0 */ 12402a808120SBarry Smith for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) { 12412a808120SBarry Smith Mat tMat; 12422a808120SBarry Smith Vec diag; 12432a808120SBarry Smith 1244849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0)); 12452a808120SBarry Smith 12462e68589bSMark F. Adams /* smooth P1 := (I - omega/lam D^{-1}A)P0 */ 12479566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0)); 12489566063dSJacob Faibussowitsch PetscCall(MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat)); 12499566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0)); 12509566063dSJacob Faibussowitsch PetscCall(MatProductClear(tMat)); 12519566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Amat, &diag, NULL)); 12529566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(Amat, diag)); /* effectively PCJACOBI */ 12539566063dSJacob Faibussowitsch PetscCall(VecReciprocal(diag)); 12549566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(tMat, diag, NULL)); 12559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&diag)); 1256b4da3a1bSStefano Zampini 1257d5d25401SBarry Smith /* TODO: Set a PCFailedReason and exit the building of the AMG preconditioner */ 125808401ef6SPierre Jolivet PetscCheck(emax != 0.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Computed maximum singular value as zero"); 1259d5d25401SBarry Smith /* TODO: Document the 1.4 and don't hardwire it in this routine */ 1260b4ec6429SMark F. Adams alpha = -1.4 / emax; 1261b4da3a1bSStefano Zampini 12629566063dSJacob Faibussowitsch PetscCall(MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN)); 12639566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Prol)); 12642e68589bSMark F. Adams Prol = tMat; 1265849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0)); 12662e68589bSMark F. Adams } 1267849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0)); 1268c8b0795cSMark F. Adams *a_P = Prol; 12693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12702e68589bSMark F. Adams } 12712e68589bSMark F. Adams 1272c8b0795cSMark F. Adams /* 1273c8b0795cSMark F. Adams PCCreateGAMG_AGG 12742e68589bSMark F. Adams 1275c8b0795cSMark F. Adams Input Parameter: 1276c8b0795cSMark F. Adams . pc - 1277c8b0795cSMark F. Adams */ 1278d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCreateGAMG_AGG(PC pc) 1279d71ae5a4SJacob Faibussowitsch { 1280c8b0795cSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 1281c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1282c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg; 12832e68589bSMark F. Adams 1284c8b0795cSMark F. Adams PetscFunctionBegin; 1285c8b0795cSMark F. Adams /* create sub context for SA */ 12864dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&pc_gamg_agg)); 1287c8b0795cSMark F. Adams pc_gamg->subctx = pc_gamg_agg; 1288c8b0795cSMark F. Adams 12891ab5ffc9SJed Brown pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG; 12909b8ffb57SJed Brown pc_gamg->ops->destroy = PCDestroy_GAMG_AGG; 1291c8b0795cSMark F. Adams /* reset does not do anything; setup not virtual */ 1292c8b0795cSMark F. Adams 1293c8b0795cSMark F. Adams /* set internal function pointers */ 12942d776b49SBarry Smith pc_gamg->ops->creategraph = PCGAMGCreateGraph_AGG; 12951ab5ffc9SJed Brown pc_gamg->ops->coarsen = PCGAMGCoarsen_AGG; 12961ab5ffc9SJed Brown pc_gamg->ops->prolongator = PCGAMGProlongator_AGG; 1297fd1112cbSBarry Smith pc_gamg->ops->optprolongator = PCGAMGOptProlongator_AGG; 12981ab5ffc9SJed Brown pc_gamg->ops->createdefaultdata = PCSetData_AGG; 12995adeb434SBarry Smith pc_gamg->ops->view = PCView_GAMG_AGG; 1300c8b0795cSMark F. Adams 1301cab9ed1eSBarry Smith pc_gamg_agg->nsmooths = 1; 1302d529f056Smarkadams4 pc_gamg_agg->aggressive_coarsening_levels = 1; 1303d529f056Smarkadams4 pc_gamg_agg->use_aggressive_square_graph = PETSC_FALSE; 1304d529f056Smarkadams4 pc_gamg_agg->use_minimum_degree_ordering = PETSC_FALSE; 1305d529f056Smarkadams4 pc_gamg_agg->aggressive_mis_k = 2; 1306cab9ed1eSBarry Smith 13079566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", PCGAMGSetNSmooths_AGG)); 1308bae903cbSmarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", PCGAMGSetAggressiveLevels_AGG)); 1309d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveSquareGraph_C", PCGAMGSetAggressiveSquareGraph_AGG)); 1310d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetMinDegreeOrdering_C", PCGAMGMISkSetMinDegreeOrdering_AGG)); 1311d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetAggressive_C", PCGAMGMISkSetAggressive_AGG)); 13129566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_AGG)); 13133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13142e68589bSMark F. Adams } 1315