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 { 11a077d33dSBarry Smith PetscInt nsmooths; // number of smoothing steps to construct prolongation 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; 16a9ccf005SMark Adams PetscBool use_low_mem_filter; 17*6c34c54dSStefano Zampini PetscBool graph_symmetrize; 182d776b49SBarry Smith MatCoarsen crs; 192e68589bSMark F. Adams } PC_GAMG_AGG; 202e68589bSMark F. Adams 212e68589bSMark F. Adams /*@ 22a077d33dSBarry Smith PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical) used to construct the prolongation operator 232e68589bSMark F. Adams 24c3339decSBarry Smith Logically Collective 252e68589bSMark F. Adams 262e68589bSMark F. Adams Input Parameters: 2720f4b53cSBarry Smith + pc - the preconditioner context 2820f4b53cSBarry Smith - n - the number of smooths 292e68589bSMark F. Adams 302e68589bSMark F. Adams Options Database Key: 31a077d33dSBarry Smith . -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use 322e68589bSMark F. Adams 332e68589bSMark F. Adams Level: intermediate 342e68589bSMark F. Adams 35a077d33dSBarry Smith Note: 36a077d33dSBarry Smith This is a different concept from the number smoothing steps used during the linear solution process which 37a077d33dSBarry Smith can be set with `-mg_levels_ksp_max_it` 38a077d33dSBarry Smith 39a077d33dSBarry Smith Developer Note: 40a077d33dSBarry Smith This should be named `PCGAMGAGGSetNSmooths()`. 41a077d33dSBarry Smith 42a077d33dSBarry Smith .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCMG`, `PCGAMG` 432e68589bSMark F. Adams @*/ 44d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n) 45d71ae5a4SJacob Faibussowitsch { 462e68589bSMark F. Adams PetscFunctionBegin; 472e68589bSMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 48d5d25401SBarry Smith PetscValidLogicalCollectiveInt(pc, n, 2); 49cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetNSmooths_C", (PC, PetscInt), (pc, n)); 503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 512e68589bSMark F. Adams } 522e68589bSMark F. Adams 53d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetNSmooths_AGG(PC pc, PetscInt n) 54d71ae5a4SJacob Faibussowitsch { 552e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 562e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 57c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 582e68589bSMark F. Adams 592e68589bSMark F. Adams PetscFunctionBegin; 60c8b0795cSMark F. Adams pc_gamg_agg->nsmooths = n; 613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 62c8b0795cSMark F. Adams } 63c8b0795cSMark F. Adams 64c8b0795cSMark F. Adams /*@ 65f1580f4eSBarry Smith PCGAMGSetAggressiveLevels - Use aggressive coarsening on first n levels 66ef4ad70eSMark F. Adams 67c3339decSBarry Smith Logically Collective 68ef4ad70eSMark F. Adams 69ef4ad70eSMark F. Adams Input Parameters: 70cab9ed1eSBarry Smith + pc - the preconditioner context 71d5d25401SBarry Smith - n - 0, 1 or more 72ef4ad70eSMark F. Adams 73ef4ad70eSMark F. Adams Options Database Key: 74a077d33dSBarry Smith . -pc_gamg_aggressive_coarsening <n,default = 1> - Number of levels on which to square the graph on before aggregating it 75af3c827dSMark Adams 76ef4ad70eSMark F. Adams Level: intermediate 77ef4ad70eSMark F. Adams 78a077d33dSBarry Smith .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGMISkSetAggressive()`, `PCGAMGSetAggressiveSquareGraph()`, `PCGAMGMISkSetMinDegreeOrdering()`, `PCGAMGSetLowMemoryFilter()` 79ef4ad70eSMark F. Adams @*/ 80d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetAggressiveLevels(PC pc, PetscInt n) 81d71ae5a4SJacob Faibussowitsch { 82ef4ad70eSMark F. Adams PetscFunctionBegin; 83ef4ad70eSMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 84d5d25401SBarry Smith PetscValidLogicalCollectiveInt(pc, n, 2); 85bae903cbSmarkadams4 PetscTryMethod(pc, "PCGAMGSetAggressiveLevels_C", (PC, PetscInt), (pc, n)); 863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 87ef4ad70eSMark F. Adams } 88ef4ad70eSMark F. Adams 89d529f056Smarkadams4 /*@ 90d529f056Smarkadams4 PCGAMGMISkSetAggressive - Number (k) distance in MIS coarsening (>2 is 'aggressive') 91d529f056Smarkadams4 92d529f056Smarkadams4 Logically Collective 93d529f056Smarkadams4 94d529f056Smarkadams4 Input Parameters: 95d529f056Smarkadams4 + pc - the preconditioner context 96d529f056Smarkadams4 - n - 1 or more (default = 2) 97d529f056Smarkadams4 98d529f056Smarkadams4 Options Database Key: 99d529f056Smarkadams4 . -pc_gamg_aggressive_mis_k <n,default=2> - Number (k) distance in MIS coarsening (>2 is 'aggressive') 100d529f056Smarkadams4 101d529f056Smarkadams4 Level: intermediate 102d529f056Smarkadams4 103a077d33dSBarry Smith .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`, `PCGAMGSetAggressiveSquareGraph()`, `PCGAMGMISkSetMinDegreeOrdering()`, `PCGAMGSetLowMemoryFilter()` 104d529f056Smarkadams4 @*/ 105d529f056Smarkadams4 PetscErrorCode PCGAMGMISkSetAggressive(PC pc, PetscInt n) 106d529f056Smarkadams4 { 107d529f056Smarkadams4 PetscFunctionBegin; 108d529f056Smarkadams4 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 109d529f056Smarkadams4 PetscValidLogicalCollectiveInt(pc, n, 2); 110d529f056Smarkadams4 PetscTryMethod(pc, "PCGAMGMISkSetAggressive_C", (PC, PetscInt), (pc, n)); 111d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 112d529f056Smarkadams4 } 113d529f056Smarkadams4 114d529f056Smarkadams4 /*@ 115d529f056Smarkadams4 PCGAMGSetAggressiveSquareGraph - Use graph square A'A for aggressive coarsening, old method 116d529f056Smarkadams4 117d529f056Smarkadams4 Logically Collective 118d529f056Smarkadams4 119d529f056Smarkadams4 Input Parameters: 120d529f056Smarkadams4 + pc - the preconditioner context 121d529f056Smarkadams4 - b - default false - MIS-k is faster 122d529f056Smarkadams4 123d529f056Smarkadams4 Options Database Key: 124d529f056Smarkadams4 . -pc_gamg_aggressive_square_graph <bool,default=false> - Use square graph (A'A) or MIS-k (k=2) for aggressive coarsening 125d529f056Smarkadams4 126d529f056Smarkadams4 Level: intermediate 127d529f056Smarkadams4 128a077d33dSBarry Smith .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`, `PCGAMGMISkSetAggressive()`, `PCGAMGMISkSetMinDegreeOrdering()`, `PCGAMGSetLowMemoryFilter()` 129d529f056Smarkadams4 @*/ 130d529f056Smarkadams4 PetscErrorCode PCGAMGSetAggressiveSquareGraph(PC pc, PetscBool b) 131d529f056Smarkadams4 { 132d529f056Smarkadams4 PetscFunctionBegin; 133d529f056Smarkadams4 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 134d529f056Smarkadams4 PetscValidLogicalCollectiveBool(pc, b, 2); 135d529f056Smarkadams4 PetscTryMethod(pc, "PCGAMGSetAggressiveSquareGraph_C", (PC, PetscBool), (pc, b)); 136d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 137d529f056Smarkadams4 } 138d529f056Smarkadams4 139d529f056Smarkadams4 /*@ 140d529f056Smarkadams4 PCGAMGMISkSetMinDegreeOrdering - Use minimum degree ordering in greedy MIS algorithm 141d529f056Smarkadams4 142d529f056Smarkadams4 Logically Collective 143d529f056Smarkadams4 144d529f056Smarkadams4 Input Parameters: 145d529f056Smarkadams4 + pc - the preconditioner context 146d529f056Smarkadams4 - b - default true 147d529f056Smarkadams4 148d529f056Smarkadams4 Options Database Key: 149d529f056Smarkadams4 . -pc_gamg_mis_k_minimum_degree_ordering <bool,default=true> - Use minimum degree ordering in greedy MIS algorithm 150d529f056Smarkadams4 151d529f056Smarkadams4 Level: intermediate 152d529f056Smarkadams4 153a077d33dSBarry Smith .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`, `PCGAMGMISkSetAggressive()`, `PCGAMGSetAggressiveSquareGraph()`, `PCGAMGSetLowMemoryFilter()` 154d529f056Smarkadams4 @*/ 155d529f056Smarkadams4 PetscErrorCode PCGAMGMISkSetMinDegreeOrdering(PC pc, PetscBool b) 156d529f056Smarkadams4 { 157d529f056Smarkadams4 PetscFunctionBegin; 158d529f056Smarkadams4 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 159d529f056Smarkadams4 PetscValidLogicalCollectiveBool(pc, b, 2); 160d529f056Smarkadams4 PetscTryMethod(pc, "PCGAMGMISkSetMinDegreeOrdering_C", (PC, PetscBool), (pc, b)); 161d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 162d529f056Smarkadams4 } 163d529f056Smarkadams4 164a9ccf005SMark Adams /*@ 165a9ccf005SMark Adams PCGAMGSetLowMemoryFilter - Use low memory graph/matrix filter 166a9ccf005SMark Adams 167a9ccf005SMark Adams Logically Collective 168a9ccf005SMark Adams 169a9ccf005SMark Adams Input Parameters: 170a9ccf005SMark Adams + pc - the preconditioner context 171a9ccf005SMark Adams - b - default false 172a9ccf005SMark Adams 173a9ccf005SMark Adams Options Database Key: 174a9ccf005SMark Adams . -pc_gamg_low_memory_threshold_filter <bool,default=false> - Use low memory graph/matrix filter 175a9ccf005SMark Adams 176a9ccf005SMark Adams Level: intermediate 177a9ccf005SMark Adams 178a077d33dSBarry Smith .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`, 179a077d33dSBarry Smith `PCGAMGMISkSetAggressive()`, `PCGAMGSetAggressiveSquareGraph()`, `PCGAMGMISkSetMinDegreeOrdering()` 180a9ccf005SMark Adams @*/ 181a9ccf005SMark Adams PetscErrorCode PCGAMGSetLowMemoryFilter(PC pc, PetscBool b) 182a9ccf005SMark Adams { 183a9ccf005SMark Adams PetscFunctionBegin; 184a9ccf005SMark Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 185a9ccf005SMark Adams PetscValidLogicalCollectiveBool(pc, b, 2); 186a9ccf005SMark Adams PetscTryMethod(pc, "PCGAMGSetLowMemoryFilter_C", (PC, PetscBool), (pc, b)); 187a9ccf005SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 188a9ccf005SMark Adams } 189a9ccf005SMark Adams 190*6c34c54dSStefano Zampini /*@ 191*6c34c54dSStefano Zampini PCGAMGSetGraphSymmetrize - Set the flag to symmetrize the graph used in coarsening 192*6c34c54dSStefano Zampini 193*6c34c54dSStefano Zampini Logically Collective 194*6c34c54dSStefano Zampini 195*6c34c54dSStefano Zampini Input Parameters: 196*6c34c54dSStefano Zampini + pc - the preconditioner context 197*6c34c54dSStefano Zampini - b - default false 198*6c34c54dSStefano Zampini 199*6c34c54dSStefano Zampini Options Database Key: 200*6c34c54dSStefano Zampini . -pc_gamg_graph_symmetrize <bool,default=false> - Symmetrize the graph 201*6c34c54dSStefano Zampini 202*6c34c54dSStefano Zampini Level: intermediate 203*6c34c54dSStefano Zampini 204*6c34c54dSStefano Zampini .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`, 205*6c34c54dSStefano Zampini `PCGAMGMISkSetAggressive()`, `PCGAMGSetAggressiveSquareGraph()`, `PCGAMGMISkSetMinDegreeOrdering()` 206*6c34c54dSStefano Zampini @*/ 207*6c34c54dSStefano Zampini PetscErrorCode PCGAMGSetGraphSymmetrize(PC pc, PetscBool b) 208*6c34c54dSStefano Zampini { 209*6c34c54dSStefano Zampini PetscFunctionBegin; 210*6c34c54dSStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 211*6c34c54dSStefano Zampini PetscValidLogicalCollectiveBool(pc, b, 2); 212*6c34c54dSStefano Zampini PetscTryMethod(pc, "PCGAMGSetGraphSymmetrize_C", (PC, PetscBool), (pc, b)); 213*6c34c54dSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 214*6c34c54dSStefano Zampini } 215*6c34c54dSStefano Zampini 216d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetAggressiveLevels_AGG(PC pc, PetscInt n) 217d71ae5a4SJacob Faibussowitsch { 218ef4ad70eSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 219ef4ad70eSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 220ef4ad70eSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 221ef4ad70eSMark F. Adams 222ef4ad70eSMark F. Adams PetscFunctionBegin; 223bae903cbSmarkadams4 pc_gamg_agg->aggressive_coarsening_levels = n; 2243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 225ef4ad70eSMark F. Adams } 226ef4ad70eSMark F. Adams 227d529f056Smarkadams4 static PetscErrorCode PCGAMGMISkSetAggressive_AGG(PC pc, PetscInt n) 228d529f056Smarkadams4 { 229d529f056Smarkadams4 PC_MG *mg = (PC_MG *)pc->data; 230d529f056Smarkadams4 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 231d529f056Smarkadams4 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 232d529f056Smarkadams4 233d529f056Smarkadams4 PetscFunctionBegin; 234d529f056Smarkadams4 pc_gamg_agg->aggressive_mis_k = n; 235d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 236d529f056Smarkadams4 } 237d529f056Smarkadams4 238d529f056Smarkadams4 static PetscErrorCode PCGAMGSetAggressiveSquareGraph_AGG(PC pc, PetscBool b) 239d529f056Smarkadams4 { 240d529f056Smarkadams4 PC_MG *mg = (PC_MG *)pc->data; 241d529f056Smarkadams4 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 242d529f056Smarkadams4 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 243d529f056Smarkadams4 244d529f056Smarkadams4 PetscFunctionBegin; 245d529f056Smarkadams4 pc_gamg_agg->use_aggressive_square_graph = b; 246d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 247d529f056Smarkadams4 } 248d529f056Smarkadams4 249a9ccf005SMark Adams static PetscErrorCode PCGAMGSetLowMemoryFilter_AGG(PC pc, PetscBool b) 250a9ccf005SMark Adams { 251a9ccf005SMark Adams PC_MG *mg = (PC_MG *)pc->data; 252a9ccf005SMark Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 253a9ccf005SMark Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 254a9ccf005SMark Adams 255a9ccf005SMark Adams PetscFunctionBegin; 256a9ccf005SMark Adams pc_gamg_agg->use_low_mem_filter = b; 257a9ccf005SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 258a9ccf005SMark Adams } 259a9ccf005SMark Adams 260*6c34c54dSStefano Zampini static PetscErrorCode PCGAMGSetGraphSymmetrize_AGG(PC pc, PetscBool b) 261*6c34c54dSStefano Zampini { 262*6c34c54dSStefano Zampini PC_MG *mg = (PC_MG *)pc->data; 263*6c34c54dSStefano Zampini PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 264*6c34c54dSStefano Zampini PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 265*6c34c54dSStefano Zampini 266*6c34c54dSStefano Zampini PetscFunctionBegin; 267*6c34c54dSStefano Zampini pc_gamg_agg->graph_symmetrize = b; 268*6c34c54dSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 269*6c34c54dSStefano Zampini } 270*6c34c54dSStefano Zampini 271d529f056Smarkadams4 static PetscErrorCode PCGAMGMISkSetMinDegreeOrdering_AGG(PC pc, PetscBool b) 272d529f056Smarkadams4 { 273d529f056Smarkadams4 PC_MG *mg = (PC_MG *)pc->data; 274d529f056Smarkadams4 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 275d529f056Smarkadams4 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 276d529f056Smarkadams4 277d529f056Smarkadams4 PetscFunctionBegin; 278d529f056Smarkadams4 pc_gamg_agg->use_minimum_degree_ordering = b; 279d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 280d529f056Smarkadams4 } 281d529f056Smarkadams4 282d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_GAMG_AGG(PC pc, PetscOptionItems *PetscOptionsObject) 283d71ae5a4SJacob Faibussowitsch { 2842e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 2852e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 286c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 287d4adc10fSMark Adams PetscBool n_aggressive_flg, old_sq_provided = PETSC_FALSE, new_sq_provided = PETSC_FALSE, new_sqr_graph = pc_gamg_agg->use_aggressive_square_graph; 288d4adc10fSMark Adams PetscInt nsq_graph_old = 0; 2892e68589bSMark F. Adams 2902e68589bSMark F. Adams PetscFunctionBegin; 291d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "GAMG-AGG options"); 292a077d33dSBarry Smith PetscCall(PetscOptionsInt("-pc_gamg_agg_nsmooths", "number of smoothing steps to construct prolongation, usually 1", "PCGAMGSetNSmooths", pc_gamg_agg->nsmooths, &pc_gamg_agg->nsmooths, NULL)); 293d4adc10fSMark Adams // aggressive coarsening logic with deprecated -pc_gamg_square_graph 294d4adc10fSMark Adams 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, &n_aggressive_flg)); 295d4adc10fSMark Adams if (!n_aggressive_flg) 296d4adc10fSMark Adams PetscCall(PetscOptionsInt("-pc_gamg_square_graph", "Number of aggressive coarsening (MIS-2) levels from finest (deprecated alias for -pc_gamg_aggressive_coarsening)", "PCGAMGSetAggressiveLevels", nsq_graph_old, &nsq_graph_old, &old_sq_provided)); 297d4adc10fSMark Adams PetscCall(PetscOptionsBool("-pc_gamg_aggressive_square_graph", "Use square graph (A'A) or MIS-k (k=2) for aggressive coarsening", "PCGAMGSetAggressiveSquareGraph", new_sqr_graph, &pc_gamg_agg->use_aggressive_square_graph, &new_sq_provided)); 298d4adc10fSMark Adams if (!new_sq_provided && old_sq_provided) { 299d4adc10fSMark Adams pc_gamg_agg->aggressive_coarsening_levels = nsq_graph_old; // could be zero 300d4adc10fSMark Adams pc_gamg_agg->use_aggressive_square_graph = PETSC_TRUE; 301bae903cbSmarkadams4 } 302d4adc10fSMark Adams if (new_sq_provided && old_sq_provided) 303d4adc10fSMark Adams 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)); 304d529f056Smarkadams4 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)); 305a9ccf005SMark Adams PetscCall(PetscOptionsBool("-pc_gamg_low_memory_threshold_filter", "Use the (built-in) low memory graph/matrix filter", "PCGAMGSetLowMemoryFilter", pc_gamg_agg->use_low_mem_filter, &pc_gamg_agg->use_low_mem_filter, NULL)); 306d529f056Smarkadams4 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)); 307*6c34c54dSStefano Zampini PetscCall(PetscOptionsBool("-pc_gamg_graph_symmetrize", "Symmetrize graph for coarsening", "PCGAMGSetGraphSymmetrize", pc_gamg_agg->graph_symmetrize, &pc_gamg_agg->graph_symmetrize, NULL)); 308d0609cedSBarry Smith PetscOptionsHeadEnd(); 3093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3102e68589bSMark F. Adams } 3112e68589bSMark F. Adams 312d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_GAMG_AGG(PC pc) 313d71ae5a4SJacob Faibussowitsch { 3142e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 3152e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 316e0b7e82fSBarry Smith PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 3172e68589bSMark F. Adams 3182e68589bSMark F. Adams PetscFunctionBegin; 319e0b7e82fSBarry Smith PetscCall(MatCoarsenDestroy(&pc_gamg_agg->crs)); 3209566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->subctx)); 3212e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", NULL)); 322bae903cbSmarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", NULL)); 323d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetAggressive_C", NULL)); 324d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetMinDegreeOrdering_C", NULL)); 325a9ccf005SMark Adams PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetLowMemoryFilter_C", NULL)); 326d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveSquareGraph_C", NULL)); 327*6c34c54dSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetGraphSymmetrize_C", NULL)); 328bae903cbSmarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL)); 3293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3302e68589bSMark F. Adams } 3312e68589bSMark F. Adams 3322e68589bSMark F. Adams /* 3332e68589bSMark F. Adams PCSetCoordinates_AGG 33420f4b53cSBarry Smith 33520f4b53cSBarry Smith Collective 3362e68589bSMark F. Adams 3372e68589bSMark F. Adams Input Parameter: 3382e68589bSMark F. Adams . pc - the preconditioner context 339145b44c9SPierre Jolivet . ndm - dimension of data (used for dof/vertex for Stokes) 340302f38e8SMark F. Adams . a_nloc - number of vertices local 341302f38e8SMark F. Adams . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...} 3422e68589bSMark F. Adams */ 3431e6b0712SBarry Smith 344d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords) 345d71ae5a4SJacob Faibussowitsch { 3462e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 3472e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 34869344418SMark F. Adams PetscInt arrsz, kk, ii, jj, nloc, ndatarows, ndf; 349a2f3521dSMark F. Adams Mat mat = pc->pmat; 3502e68589bSMark F. Adams 3512e68589bSMark F. Adams PetscFunctionBegin; 352a2f3521dSMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 353a2f3521dSMark F. Adams PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 354302f38e8SMark F. Adams nloc = a_nloc; 3552e68589bSMark F. Adams 3562e68589bSMark F. Adams /* SA: null space vectors */ 3579566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(mat, &ndf)); /* this does not work for Stokes */ 35869344418SMark F. Adams if (coords && ndf == 1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */ 359a2f3521dSMark F. Adams else if (coords) { 36063a3b9bcSJacob Faibussowitsch PetscCheck(ndm <= ndf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "degrees of motion %" PetscInt_FMT " > block size %" PetscInt_FMT, ndm, ndf); 36169344418SMark F. Adams pc_gamg->data_cell_cols = (ndm == 2 ? 3 : 6); /* displacement elasticity */ 362ad540459SPierre 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); 363806fa848SBarry Smith } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */ 36471959b99SBarry Smith pc_gamg->data_cell_rows = ndatarows = ndf; 36563a3b9bcSJacob 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); 366c8b0795cSMark F. Adams arrsz = nloc * pc_gamg->data_cell_rows * pc_gamg->data_cell_cols; 3672e68589bSMark F. Adams 3687f66b68fSMark Adams if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) { 3699566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 3709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(arrsz + 1, &pc_gamg->data)); 3712e68589bSMark F. Adams } 3725e116b59SBarry Smith /* copy data in - column-oriented */ 3732e68589bSMark F. Adams for (kk = 0; kk < nloc; kk++) { 37469344418SMark F. Adams const PetscInt M = nloc * pc_gamg->data_cell_rows; /* stride into data */ 37569344418SMark F. Adams PetscReal *data = &pc_gamg->data[kk * ndatarows]; /* start of cell */ 376e0b7e82fSBarry Smith 377c8b0795cSMark F. Adams if (pc_gamg->data_cell_cols == 1) *data = 1.0; 3782e68589bSMark F. Adams else { 37969344418SMark F. Adams /* translational modes */ 3802fa5cd67SKarl Rupp for (ii = 0; ii < ndatarows; ii++) { 3812fa5cd67SKarl Rupp for (jj = 0; jj < ndatarows; jj++) { 38269344418SMark F. Adams if (ii == jj) data[ii * M + jj] = 1.0; 3832e68589bSMark F. Adams else data[ii * M + jj] = 0.0; 3842fa5cd67SKarl Rupp } 3852fa5cd67SKarl Rupp } 38669344418SMark F. Adams 38769344418SMark F. Adams /* rotational modes */ 3882e68589bSMark F. Adams if (coords) { 38969344418SMark F. Adams if (ndm == 2) { 3902e68589bSMark F. Adams data += 2 * M; 3912e68589bSMark F. Adams data[0] = -coords[2 * kk + 1]; 3922e68589bSMark F. Adams data[1] = coords[2 * kk]; 393806fa848SBarry Smith } else { 3942e68589bSMark F. Adams data += 3 * M; 3959371c9d4SSatish Balay data[0] = 0.0; 3969371c9d4SSatish Balay data[M + 0] = coords[3 * kk + 2]; 3979371c9d4SSatish Balay data[2 * M + 0] = -coords[3 * kk + 1]; 3989371c9d4SSatish Balay data[1] = -coords[3 * kk + 2]; 3999371c9d4SSatish Balay data[M + 1] = 0.0; 4009371c9d4SSatish Balay data[2 * M + 1] = coords[3 * kk]; 4019371c9d4SSatish Balay data[2] = coords[3 * kk + 1]; 4029371c9d4SSatish Balay data[M + 2] = -coords[3 * kk]; 4039371c9d4SSatish Balay data[2 * M + 2] = 0.0; 4042e68589bSMark F. Adams } 4052e68589bSMark F. Adams } 4062e68589bSMark F. Adams } 4072e68589bSMark F. Adams } 4082e68589bSMark F. Adams pc_gamg->data_sz = arrsz; 4093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4102e68589bSMark F. Adams } 4112e68589bSMark F. Adams 4122e68589bSMark F. Adams /* 413a2f3521dSMark F. Adams PCSetData_AGG - called if data is not set with PCSetCoordinates. 414a2f3521dSMark F. Adams Looks in Mat for near null space. 415a2f3521dSMark F. Adams Does not work for Stokes 4162e68589bSMark F. Adams 4172e68589bSMark F. Adams Input Parameter: 4182e68589bSMark F. Adams . pc - 419a2f3521dSMark F. Adams . a_A - matrix to get (near) null space out of. 4202e68589bSMark F. Adams */ 421d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A) 422d71ae5a4SJacob Faibussowitsch { 423b8cd405aSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 424b8cd405aSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 425b8cd405aSMark F. Adams MatNullSpace mnull; 42666f2ef4bSMark Adams 427b817416eSBarry Smith PetscFunctionBegin; 4289566063dSJacob Faibussowitsch PetscCall(MatGetNearNullSpace(a_A, &mnull)); 429b8cd405aSMark F. Adams if (!mnull) { 43066f2ef4bSMark Adams DM dm; 431e0b7e82fSBarry Smith 4329566063dSJacob Faibussowitsch PetscCall(PCGetDM(pc, &dm)); 43348a46eb9SPierre Jolivet if (!dm) PetscCall(MatGetDM(a_A, &dm)); 43466f2ef4bSMark Adams if (dm) { 43566f2ef4bSMark Adams PetscObject deformation; 436b0d30dd6SMatthew G. Knepley PetscInt Nf; 437b0d30dd6SMatthew G. Knepley 4389566063dSJacob Faibussowitsch PetscCall(DMGetNumFields(dm, &Nf)); 439b0d30dd6SMatthew G. Knepley if (Nf) { 4409566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, 0, NULL, &deformation)); 4419566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)deformation, "nearnullspace", (PetscObject *)&mnull)); 44248a46eb9SPierre Jolivet if (!mnull) PetscCall(PetscObjectQuery((PetscObject)deformation, "nullspace", (PetscObject *)&mnull)); 44366f2ef4bSMark Adams } 44466f2ef4bSMark Adams } 445b0d30dd6SMatthew G. Knepley } 44666f2ef4bSMark Adams 44766f2ef4bSMark Adams if (!mnull) { 448a2f3521dSMark F. Adams PetscInt bs, NN, MM; 449e0b7e82fSBarry Smith 4509566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(a_A, &bs)); 4519566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(a_A, &MM, &NN)); 45263a3b9bcSJacob Faibussowitsch PetscCheck(MM % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MM %" PetscInt_FMT " must be divisible by bs %" PetscInt_FMT, MM, bs); 4539566063dSJacob Faibussowitsch PetscCall(PCSetCoordinates_AGG(pc, bs, MM / bs, NULL)); 454806fa848SBarry Smith } else { 455b8cd405aSMark F. Adams PetscReal *nullvec; 456b8cd405aSMark F. Adams PetscBool has_const; 457b8cd405aSMark F. Adams PetscInt i, j, mlocal, nvec, bs; 458d5d25401SBarry Smith const Vec *vecs; 459d5d25401SBarry Smith const PetscScalar *v; 460b817416eSBarry Smith 4619566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(a_A, &mlocal, NULL)); 4629566063dSJacob Faibussowitsch PetscCall(MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs)); 463d19c4ebbSmarkadams4 for (i = 0; i < nvec; i++) { 464d19c4ebbSmarkadams4 PetscCall(VecGetLocalSize(vecs[i], &j)); 465d19c4ebbSmarkadams4 PetscCheck(j == mlocal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Attached null space vector size %" PetscInt_FMT " != matrix size %" PetscInt_FMT, j, mlocal); 466d19c4ebbSmarkadams4 } 467a0dea87dSMark Adams pc_gamg->data_sz = (nvec + !!has_const) * mlocal; 4689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((nvec + !!has_const) * mlocal, &nullvec)); 4699371c9d4SSatish Balay if (has_const) 4709371c9d4SSatish Balay for (i = 0; i < mlocal; i++) nullvec[i] = 1.0; 471b8cd405aSMark F. Adams for (i = 0; i < nvec; i++) { 4729566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(vecs[i], &v)); 473b8cd405aSMark F. Adams for (j = 0; j < mlocal; j++) nullvec[(i + !!has_const) * mlocal + j] = PetscRealPart(v[j]); 4749566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(vecs[i], &v)); 475b8cd405aSMark F. Adams } 476b8cd405aSMark F. Adams pc_gamg->data = nullvec; 477b8cd405aSMark F. Adams pc_gamg->data_cell_cols = (nvec + !!has_const); 4789566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(a_A, &bs)); 479b8cd405aSMark F. Adams pc_gamg->data_cell_rows = bs; 480b8cd405aSMark F. Adams } 4813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4822e68589bSMark F. Adams } 4832e68589bSMark F. Adams 4842e68589bSMark F. Adams /* 485bae903cbSmarkadams4 formProl0 - collect null space data for each aggregate, do QR, put R in coarse grid data and Q in P_0 4862e68589bSMark F. Adams 4872e68589bSMark F. Adams Input Parameter: 4882adfe9d3SBarry Smith . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices 4892adfe9d3SBarry Smith . bs - row block size 4900cbbd2e1SMark F. Adams . nSAvec - column bs of new P 4910cbbd2e1SMark F. Adams . my0crs - global index of start of locals 4922adfe9d3SBarry Smith . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec] 4932e68589bSMark F. Adams . data_in[data_stride*nSAvec] - local data on fine grid 4942e68589bSMark F. Adams . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist' 495bae903cbSmarkadams4 4962e68589bSMark F. Adams Output Parameter: 4972e68589bSMark F. Adams . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data 4982e68589bSMark F. Adams . a_Prol - prolongation operator 4992e68589bSMark F. Adams */ 500d71ae5a4SJacob 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) 501d71ae5a4SJacob Faibussowitsch { 502bd026e97SJed Brown PetscInt Istart, my0, Iend, nloc, clid, flid = 0, aggID, kk, jj, ii, mm, nSelected, minsz, nghosts, out_data_stride; 5033b4367a7SBarry Smith MPI_Comm comm; 5042e68589bSMark F. Adams PetscReal *out_data; 505539c167fSBarry Smith PetscCDIntNd *pos; 5061943db53SBarry Smith PCGAMGHashTable fgid_flid; 5070cbbd2e1SMark F. Adams 5082e68589bSMark F. Adams PetscFunctionBegin; 5099566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)a_Prol, &comm)); 5109566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(a_Prol, &Istart, &Iend)); 5119371c9d4SSatish Balay nloc = (Iend - Istart) / bs; 5129371c9d4SSatish Balay my0 = Istart / bs; 51363a3b9bcSJacob 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); 5140cbbd2e1SMark F. Adams Iend /= bs; 5150cbbd2e1SMark F. Adams nghosts = data_stride / bs - nloc; 5162e68589bSMark F. Adams 5179566063dSJacob Faibussowitsch PetscCall(PCGAMGHashTableCreate(2 * nghosts + 1, &fgid_flid)); 51848a46eb9SPierre Jolivet for (kk = 0; kk < nghosts; kk++) PetscCall(PCGAMGHashTableAdd(&fgid_flid, flid_fgid[nloc + kk], nloc + kk)); 5192e68589bSMark F. Adams 5200cbbd2e1SMark F. Adams /* count selected -- same as number of cols of P */ 5210cbbd2e1SMark F. Adams for (nSelected = mm = 0; mm < nloc; mm++) { 522e78576d6SMark F. Adams PetscBool ise; 523e0b7e82fSBarry Smith 5245e62d202SMark Adams PetscCall(PetscCDIsEmptyAt(agg_llists, mm, &ise)); 525e78576d6SMark F. Adams if (!ise) nSelected++; 5260cbbd2e1SMark F. Adams } 5279566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(a_Prol, &ii, &jj)); 52863a3b9bcSJacob Faibussowitsch PetscCheck((ii / nSAvec) == my0crs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ii %" PetscInt_FMT " /nSAvec %" PetscInt_FMT " != my0crs %" PetscInt_FMT, ii, nSAvec, my0crs); 52963a3b9bcSJacob 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); 5300cbbd2e1SMark F. Adams 5312e68589bSMark F. Adams /* aloc space for coarse point data (output) */ 5320cbbd2e1SMark F. Adams out_data_stride = nSelected * nSAvec; 5332fa5cd67SKarl Rupp 5349566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(out_data_stride * nSAvec, &out_data)); 53533812677SSatish Balay for (ii = 0; ii < out_data_stride * nSAvec; ii++) out_data[ii] = PETSC_MAX_REAL; 5360cbbd2e1SMark F. Adams *a_data_out = out_data; /* output - stride nSelected*nSAvec */ 5372e68589bSMark F. Adams 5382e68589bSMark F. Adams /* find points and set prolongation */ 539c8b0795cSMark F. Adams minsz = 100; 5400cbbd2e1SMark F. Adams for (mm = clid = 0; mm < nloc; mm++) { 5415e62d202SMark Adams PetscCall(PetscCDCountAt(agg_llists, mm, &jj)); 542e78576d6SMark F. Adams if (jj > 0) { 5430cbbd2e1SMark F. Adams const PetscInt lid = mm, cgid = my0crs + clid; 5440cbbd2e1SMark F. Adams PetscInt cids[100]; /* max bs */ 5450cbbd2e1SMark F. Adams PetscBLASInt asz = jj, M = asz * bs, N = nSAvec, INFO; 5462e68589bSMark F. Adams PetscBLASInt Mdata = M + ((N - M > 0) ? N - M : 0), LDA = Mdata, LWORK = N * bs; 5472e68589bSMark F. Adams PetscScalar *qqc, *qqr, *TAU, *WORK; 5482e68589bSMark F. Adams PetscInt *fids; 54965d7b583SSatish Balay PetscReal *data; 550b817416eSBarry Smith 5510cbbd2e1SMark F. Adams /* count agg */ 5520cbbd2e1SMark F. Adams if (asz < minsz) minsz = asz; 5530cbbd2e1SMark F. Adams 5540cbbd2e1SMark F. Adams /* get block */ 5559566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(Mdata * N, &qqc, M * N, &qqr, N, &TAU, LWORK, &WORK, M, &fids)); 5562e68589bSMark F. Adams 5572e68589bSMark F. Adams aggID = 0; 5589566063dSJacob Faibussowitsch PetscCall(PetscCDGetHeadPos(agg_llists, lid, &pos)); 559e78576d6SMark F. Adams while (pos) { 560e78576d6SMark F. Adams PetscInt gid1; 561e0b7e82fSBarry Smith 5629566063dSJacob Faibussowitsch PetscCall(PetscCDIntNdGetID(pos, &gid1)); 5639566063dSJacob Faibussowitsch PetscCall(PetscCDGetNextPos(agg_llists, lid, &pos)); 564e78576d6SMark F. Adams 5650cbbd2e1SMark F. Adams if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0; 5660cbbd2e1SMark F. Adams else { 5679566063dSJacob Faibussowitsch PetscCall(PCGAMGHashTableFind(&fgid_flid, gid1, &flid)); 56808401ef6SPierre Jolivet PetscCheck(flid >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot find gid1 in table"); 5690cbbd2e1SMark F. Adams } 5705e116b59SBarry Smith /* copy in B_i matrix - column-oriented */ 57165d7b583SSatish Balay data = &data_in[flid * bs]; 5723d3eaba7SBarry Smith for (ii = 0; ii < bs; ii++) { 5732e68589bSMark F. Adams for (jj = 0; jj < N; jj++) { 5740cbbd2e1SMark F. Adams PetscReal d = data[jj * data_stride + ii]; 575e0b7e82fSBarry Smith 5760cbbd2e1SMark F. Adams qqc[jj * Mdata + aggID * bs + ii] = d; 5772e68589bSMark F. Adams } 5782e68589bSMark F. Adams } 5792e68589bSMark F. Adams /* set fine IDs */ 5802e68589bSMark F. Adams for (kk = 0; kk < bs; kk++) fids[aggID * bs + kk] = flid_fgid[flid] * bs + kk; 5812e68589bSMark F. Adams aggID++; 5820cbbd2e1SMark F. Adams } 5832e68589bSMark F. Adams 5842e68589bSMark F. Adams /* pad with zeros */ 5852e68589bSMark F. Adams for (ii = asz * bs; ii < Mdata; ii++) { 586ad540459SPierre Jolivet for (jj = 0; jj < N; jj++, kk++) qqc[jj * Mdata + ii] = .0; 5872e68589bSMark F. Adams } 5882e68589bSMark F. Adams 5892e68589bSMark F. Adams /* QR */ 5909566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 591792fecdfSBarry Smith PetscCallBLAS("LAPACKgeqrf", LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 5929566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPop()); 59308401ef6SPierre Jolivet PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xGEQRF error"); 5945e116b59SBarry Smith /* get R - column-oriented - output B_{i+1} */ 5952e68589bSMark F. Adams { 5962e68589bSMark F. Adams PetscReal *data = &out_data[clid * nSAvec]; 597e0b7e82fSBarry Smith 5982e68589bSMark F. Adams for (jj = 0; jj < nSAvec; jj++) { 5992e68589bSMark F. Adams for (ii = 0; ii < nSAvec; ii++) { 60008401ef6SPierre 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); 6010cbbd2e1SMark F. Adams if (ii <= jj) data[jj * out_data_stride + ii] = PetscRealPart(qqc[jj * Mdata + ii]); 6020cbbd2e1SMark F. Adams else data[jj * out_data_stride + ii] = 0.; 6032e68589bSMark F. Adams } 6042e68589bSMark F. Adams } 6052e68589bSMark F. Adams } 6062e68589bSMark F. Adams 6075e116b59SBarry Smith /* get Q - row-oriented */ 608792fecdfSBarry Smith PetscCallBLAS("LAPACKorgqr", LAPACKorgqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 60963a3b9bcSJacob Faibussowitsch PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xORGQR error arg %" PetscBLASInt_FMT, -INFO); 6102e68589bSMark F. Adams 6112e68589bSMark F. Adams for (ii = 0; ii < M; ii++) { 612ad540459SPierre Jolivet for (jj = 0; jj < N; jj++) qqr[N * ii + jj] = qqc[jj * Mdata + ii]; 6132e68589bSMark F. Adams } 6142e68589bSMark F. Adams 6152e68589bSMark F. Adams /* add diagonal block of P0 */ 6169371c9d4SSatish Balay for (kk = 0; kk < N; kk++) { cids[kk] = N * cgid + kk; /* global col IDs in P0 */ } 6179566063dSJacob Faibussowitsch PetscCall(MatSetValues(a_Prol, M, fids, N, cids, qqr, INSERT_VALUES)); 6189566063dSJacob Faibussowitsch PetscCall(PetscFree5(qqc, qqr, TAU, WORK, fids)); 619b43b03e9SMark F. Adams clid++; 6200cbbd2e1SMark F. Adams } /* coarse agg */ 6210cbbd2e1SMark F. Adams } /* for all fine nodes */ 6229566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(a_Prol, MAT_FINAL_ASSEMBLY)); 6239566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(a_Prol, MAT_FINAL_ASSEMBLY)); 6249566063dSJacob Faibussowitsch PetscCall(PCGAMGHashTableDestroy(&fgid_flid)); 6253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6262e68589bSMark F. Adams } 6272e68589bSMark F. Adams 628d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_GAMG_AGG(PC pc, PetscViewer viewer) 629d71ae5a4SJacob Faibussowitsch { 6305adeb434SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 6315adeb434SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 6325adeb434SBarry Smith PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 6335adeb434SBarry Smith 6345adeb434SBarry Smith PetscFunctionBegin; 6359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " AGG specific options\n")); 636d529f056Smarkadams4 PetscCall(PetscViewerASCIIPrintf(viewer, " Number of levels of aggressive coarsening %d\n", (int)pc_gamg_agg->aggressive_coarsening_levels)); 637d529f056Smarkadams4 if (pc_gamg_agg->aggressive_coarsening_levels > 0) { 638d529f056Smarkadams4 PetscCall(PetscViewerASCIIPrintf(viewer, " %s aggressive coarsening\n", !pc_gamg_agg->use_aggressive_square_graph ? "MIS-k" : "Square graph")); 639d529f056Smarkadams4 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)); 640d529f056Smarkadams4 } 641e0b7e82fSBarry Smith PetscCall(PetscViewerASCIIPushTab(viewer)); 642e0b7e82fSBarry Smith PetscCall(PetscViewerASCIIPushTab(viewer)); 643e0b7e82fSBarry Smith PetscCall(PetscViewerASCIIPushTab(viewer)); 644e0b7e82fSBarry Smith PetscCall(PetscViewerASCIIPushTab(viewer)); 6450acf423eSBarry Smith if (pc_gamg_agg->crs) PetscCall(MatCoarsenView(pc_gamg_agg->crs, viewer)); 6460acf423eSBarry Smith else PetscCall(PetscViewerASCIIPrintf(viewer, "Coarsening algorithm not yet selected\n")); 647e0b7e82fSBarry Smith PetscCall(PetscViewerASCIIPopTab(viewer)); 648e0b7e82fSBarry Smith PetscCall(PetscViewerASCIIPopTab(viewer)); 649e0b7e82fSBarry Smith PetscCall(PetscViewerASCIIPopTab(viewer)); 650e0b7e82fSBarry Smith PetscCall(PetscViewerASCIIPopTab(viewer)); 651e0b7e82fSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, " Number smoothing steps to construct prolongation %d\n", (int)pc_gamg_agg->nsmooths)); 6523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6535adeb434SBarry Smith } 6545adeb434SBarry Smith 6552d776b49SBarry Smith static PetscErrorCode PCGAMGCreateGraph_AGG(PC pc, Mat Amat, Mat *a_Gmat) 656d71ae5a4SJacob Faibussowitsch { 657c8b0795cSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 658c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 659c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 6602d776b49SBarry Smith const PetscReal vfilter = pc_gamg->threshold[pc_gamg->current_level]; 6619c15c1aeSMark Adams PetscBool ishem, ismis; 6622d776b49SBarry Smith const char *prefix; 663d529f056Smarkadams4 MatInfo info0, info1; 664d529f056Smarkadams4 PetscInt bs; 665c8b0795cSMark F. Adams 666c8b0795cSMark F. Adams PetscFunctionBegin; 667a9ccf005SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0)); 6682d776b49SBarry 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 */ 6692d776b49SBarry Smith /* MATCOARSENHEM requires numerical weights for edges so ensure they are computed */ 670e0b7e82fSBarry Smith PetscCall(MatCoarsenDestroy(&pc_gamg_agg->crs)); 6712d776b49SBarry Smith PetscCall(MatCoarsenCreate(PetscObjectComm((PetscObject)pc), &pc_gamg_agg->crs)); 6722d776b49SBarry Smith PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, &prefix)); 6732d776b49SBarry Smith PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc_gamg_agg->crs, prefix)); 674e0b7e82fSBarry Smith PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)pc_gamg_agg->crs, "pc_gamg_")); 6752d776b49SBarry Smith PetscCall(MatCoarsenSetFromOptions(pc_gamg_agg->crs)); 676e02fb3cdSMark Adams PetscCall(MatGetBlockSize(Amat, &bs)); 677e02fb3cdSMark Adams // check for valid indices wrt bs 678e02fb3cdSMark Adams for (int ii = 0; ii < pc_gamg_agg->crs->strength_index_size; ii++) { 679b65aec2dSMark Adams PetscCheck(pc_gamg_agg->crs->strength_index[ii] >= 0 && pc_gamg_agg->crs->strength_index[ii] < bs, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Indices (%d) must be non-negative and < block size (%d), NB, can not use -mat_coarsen_strength_index with -mat_coarsen_strength_index", 680b65aec2dSMark Adams (int)pc_gamg_agg->crs->strength_index[ii], (int)bs); 681e02fb3cdSMark Adams } 6822d776b49SBarry Smith PetscCall(PetscObjectTypeCompare((PetscObject)pc_gamg_agg->crs, MATCOARSENHEM, &ishem)); 6835e62d202SMark Adams if (ishem) { 6849c15c1aeSMark Adams if (pc_gamg_agg->aggressive_coarsening_levels) PetscCall(PetscInfo(pc, "HEM and aggressive coarsening ignored: HEM using %d iterations\n", (int)pc_gamg_agg->crs->max_it)); 6855e62d202SMark Adams pc_gamg_agg->aggressive_coarsening_levels = 0; // aggressive and HEM does not make sense 6865e62d202SMark Adams PetscCall(MatCoarsenSetMaximumIterations(pc_gamg_agg->crs, pc_gamg_agg->crs->max_it)); // for code coverage 6875e62d202SMark Adams PetscCall(MatCoarsenSetThreshold(pc_gamg_agg->crs, vfilter)); // for code coverage 6889c15c1aeSMark Adams } else { 6899c15c1aeSMark Adams PetscCall(PetscObjectTypeCompare((PetscObject)pc_gamg_agg->crs, MATCOARSENMIS, &ismis)); 6909c15c1aeSMark Adams if (ismis && pc_gamg_agg->aggressive_coarsening_levels && !pc_gamg_agg->use_aggressive_square_graph) { 6919c15c1aeSMark Adams PetscCall(PetscInfo(pc, "MIS and aggressive coarsening and no square graph: force square graph\n")); 6929c15c1aeSMark Adams pc_gamg_agg->use_aggressive_square_graph = PETSC_TRUE; 6939c15c1aeSMark Adams } 6945e62d202SMark Adams } 695a9ccf005SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0)); 696a9ccf005SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0)); 697d529f056Smarkadams4 PetscCall(MatGetInfo(Amat, MAT_LOCAL, &info0)); /* global reduction */ 698a9ccf005SMark Adams 699a9ccf005SMark Adams if (ishem || pc_gamg_agg->use_low_mem_filter) { 700*6c34c54dSStefano Zampini PetscCall(MatCreateGraph(Amat, pc_gamg_agg->graph_symmetrize, (vfilter >= 0 || ishem) ? PETSC_TRUE : PETSC_FALSE, vfilter, pc_gamg_agg->crs->strength_index_size, pc_gamg_agg->crs->strength_index, a_Gmat)); 701a9ccf005SMark Adams } else { 702*6c34c54dSStefano Zampini // make scalar graph, symmetrize if not known to be symmetric, scale, but do not filter (expensive) 703*6c34c54dSStefano Zampini PetscCall(MatCreateGraph(Amat, pc_gamg_agg->graph_symmetrize, PETSC_TRUE, -1, pc_gamg_agg->crs->strength_index_size, pc_gamg_agg->crs->strength_index, a_Gmat)); 704a9ccf005SMark Adams if (vfilter >= 0) { 705a9ccf005SMark Adams PetscInt Istart, Iend, ncols, nnz0, nnz1, NN, MM, nloc; 706a9ccf005SMark Adams Mat tGmat, Gmat = *a_Gmat; 707a9ccf005SMark Adams MPI_Comm comm; 708a9ccf005SMark Adams const PetscScalar *vals; 709a9ccf005SMark Adams const PetscInt *idx; 710a9ccf005SMark Adams PetscInt *d_nnz, *o_nnz, kk, *garray = NULL, *AJ, maxcols = 0; 711a9ccf005SMark Adams MatScalar *AA; // this is checked in graph 712a9ccf005SMark Adams PetscBool isseqaij; 713a9ccf005SMark Adams Mat a, b, c; 714a9ccf005SMark Adams MatType jtype; 715a9ccf005SMark Adams 716a9ccf005SMark Adams PetscCall(PetscObjectGetComm((PetscObject)Gmat, &comm)); 717a9ccf005SMark Adams PetscCall(PetscObjectBaseTypeCompare((PetscObject)Gmat, MATSEQAIJ, &isseqaij)); 718a9ccf005SMark Adams PetscCall(MatGetType(Gmat, &jtype)); 719a9ccf005SMark Adams PetscCall(MatCreate(comm, &tGmat)); 720a9ccf005SMark Adams PetscCall(MatSetType(tGmat, jtype)); 721a9ccf005SMark Adams 722a9ccf005SMark Adams /* TODO GPU: this can be called when filter = 0 -> Probably provide MatAIJThresholdCompress that compresses the entries below a threshold? 723a9ccf005SMark Adams Also, if the matrix is symmetric, can we skip this 724a9ccf005SMark Adams operation? It can be very expensive on large matrices. */ 725a9ccf005SMark Adams 726a9ccf005SMark Adams // global sizes 727a9ccf005SMark Adams PetscCall(MatGetSize(Gmat, &MM, &NN)); 728a9ccf005SMark Adams PetscCall(MatGetOwnershipRange(Gmat, &Istart, &Iend)); 729a9ccf005SMark Adams nloc = Iend - Istart; 730a9ccf005SMark Adams PetscCall(PetscMalloc2(nloc, &d_nnz, nloc, &o_nnz)); 731a9ccf005SMark Adams if (isseqaij) { 732a9ccf005SMark Adams a = Gmat; 733a9ccf005SMark Adams b = NULL; 734a9ccf005SMark Adams } else { 735a9ccf005SMark Adams Mat_MPIAIJ *d = (Mat_MPIAIJ *)Gmat->data; 736e0b7e82fSBarry Smith 737a9ccf005SMark Adams a = d->A; 738a9ccf005SMark Adams b = d->B; 739a9ccf005SMark Adams garray = d->garray; 740a9ccf005SMark Adams } 741a9ccf005SMark Adams /* Determine upper bound on non-zeros needed in new filtered matrix */ 742a9ccf005SMark Adams for (PetscInt row = 0; row < nloc; row++) { 743a9ccf005SMark Adams PetscCall(MatGetRow(a, row, &ncols, NULL, NULL)); 744a9ccf005SMark Adams d_nnz[row] = ncols; 745a9ccf005SMark Adams if (ncols > maxcols) maxcols = ncols; 746a9ccf005SMark Adams PetscCall(MatRestoreRow(a, row, &ncols, NULL, NULL)); 747a9ccf005SMark Adams } 748a9ccf005SMark Adams if (b) { 749a9ccf005SMark Adams for (PetscInt row = 0; row < nloc; row++) { 750a9ccf005SMark Adams PetscCall(MatGetRow(b, row, &ncols, NULL, NULL)); 751a9ccf005SMark Adams o_nnz[row] = ncols; 752a9ccf005SMark Adams if (ncols > maxcols) maxcols = ncols; 753a9ccf005SMark Adams PetscCall(MatRestoreRow(b, row, &ncols, NULL, NULL)); 754a9ccf005SMark Adams } 755a9ccf005SMark Adams } 756a9ccf005SMark Adams PetscCall(MatSetSizes(tGmat, nloc, nloc, MM, MM)); 757a9ccf005SMark Adams PetscCall(MatSetBlockSizes(tGmat, 1, 1)); 758a9ccf005SMark Adams PetscCall(MatSeqAIJSetPreallocation(tGmat, 0, d_nnz)); 759a9ccf005SMark Adams PetscCall(MatMPIAIJSetPreallocation(tGmat, 0, d_nnz, 0, o_nnz)); 760a9ccf005SMark Adams PetscCall(MatSetOption(tGmat, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); 761a9ccf005SMark Adams PetscCall(PetscFree2(d_nnz, o_nnz)); 762a9ccf005SMark Adams PetscCall(PetscMalloc2(maxcols, &AA, maxcols, &AJ)); 763a9ccf005SMark Adams nnz0 = nnz1 = 0; 764a9ccf005SMark Adams for (c = a, kk = 0; c && kk < 2; c = b, kk++) { 765a9ccf005SMark Adams for (PetscInt row = 0, grow = Istart, ncol_row, jj; row < nloc; row++, grow++) { 766a9ccf005SMark Adams PetscCall(MatGetRow(c, row, &ncols, &idx, &vals)); 767a9ccf005SMark Adams for (ncol_row = jj = 0; jj < ncols; jj++, nnz0++) { 768a9ccf005SMark Adams PetscScalar sv = PetscAbs(PetscRealPart(vals[jj])); 769a9ccf005SMark Adams if (PetscRealPart(sv) > vfilter) { 770a9ccf005SMark Adams PetscInt cid = idx[jj] + Istart; //diag 771e0b7e82fSBarry Smith 772a9ccf005SMark Adams nnz1++; 773a9ccf005SMark Adams if (c != a) cid = garray[idx[jj]]; 774a9ccf005SMark Adams AA[ncol_row] = vals[jj]; 775a9ccf005SMark Adams AJ[ncol_row] = cid; 776a9ccf005SMark Adams ncol_row++; 777a9ccf005SMark Adams } 778a9ccf005SMark Adams } 779a9ccf005SMark Adams PetscCall(MatRestoreRow(c, row, &ncols, &idx, &vals)); 780a9ccf005SMark Adams PetscCall(MatSetValues(tGmat, 1, &grow, ncol_row, AJ, AA, INSERT_VALUES)); 781a9ccf005SMark Adams } 782a9ccf005SMark Adams } 783a9ccf005SMark Adams PetscCall(PetscFree2(AA, AJ)); 784a9ccf005SMark Adams PetscCall(MatAssemblyBegin(tGmat, MAT_FINAL_ASSEMBLY)); 785a9ccf005SMark Adams PetscCall(MatAssemblyEnd(tGmat, MAT_FINAL_ASSEMBLY)); 786a9ccf005SMark Adams PetscCall(MatPropagateSymmetryOptions(Gmat, tGmat)); /* Normal Mat options are not relevant ? */ 787a9ccf005SMark Adams PetscCall(PetscInfo(pc, "\t %g%% nnz after filtering, with threshold %g, %g nnz ave. (N=%" PetscInt_FMT ", max row size %" PetscInt_FMT "\n", (!nnz0) ? 1. : 100. * (double)nnz1 / (double)nnz0, (double)vfilter, (!nloc) ? 1. : (double)nnz0 / (double)nloc, MM, maxcols)); 788a9ccf005SMark Adams PetscCall(MatViewFromOptions(tGmat, NULL, "-mat_filter_graph_view")); 789a9ccf005SMark Adams PetscCall(MatDestroy(&Gmat)); 790a9ccf005SMark Adams *a_Gmat = tGmat; 791a9ccf005SMark Adams } 792a9ccf005SMark Adams } 793a9ccf005SMark Adams 794d529f056Smarkadams4 PetscCall(MatGetInfo(*a_Gmat, MAT_LOCAL, &info1)); /* global reduction */ 795d529f056Smarkadams4 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)); 796849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0)); 7973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 798c8b0795cSMark F. Adams } 799c8b0795cSMark F. Adams 800d529f056Smarkadams4 typedef PetscInt NState; 801d529f056Smarkadams4 static const NState NOT_DONE = -2; 802d529f056Smarkadams4 static const NState DELETED = -1; 803d529f056Smarkadams4 static const NState REMOVED = -3; 804d529f056Smarkadams4 #define IS_SELECTED(s) (s != DELETED && s != NOT_DONE && s != REMOVED) 805d529f056Smarkadams4 806d529f056Smarkadams4 /* 807d529f056Smarkadams4 fixAggregatesWithSquare - greedy grab of with G1 (unsquared graph) -- AIJ specific -- change to fixAggregatesWithSquare -- TODD 808d529f056Smarkadams4 - AGG-MG specific: clears singletons out of 'selected_2' 809d529f056Smarkadams4 810d529f056Smarkadams4 Input Parameter: 811d529f056Smarkadams4 . Gmat_2 - global matrix of squared graph (data not defined) 812d529f056Smarkadams4 . Gmat_1 - base graph to grab with base graph 813d529f056Smarkadams4 Input/Output Parameter: 814d529f056Smarkadams4 . aggs_2 - linked list of aggs with gids) 815d529f056Smarkadams4 */ 816d529f056Smarkadams4 static PetscErrorCode fixAggregatesWithSquare(PC pc, Mat Gmat_2, Mat Gmat_1, PetscCoarsenData *aggs_2) 817d529f056Smarkadams4 { 818d529f056Smarkadams4 PetscBool isMPI; 819d529f056Smarkadams4 Mat_SeqAIJ *matA_1, *matB_1 = NULL; 820d529f056Smarkadams4 MPI_Comm comm; 821d529f056Smarkadams4 PetscInt lid, *ii, *idx, ix, Iend, my0, kk, n, j; 822d529f056Smarkadams4 Mat_MPIAIJ *mpimat_2 = NULL, *mpimat_1 = NULL; 823d529f056Smarkadams4 const PetscInt nloc = Gmat_2->rmap->n; 824d529f056Smarkadams4 PetscScalar *cpcol_1_state, *cpcol_2_state, *cpcol_2_par_orig, *lid_parent_gid; 825d529f056Smarkadams4 PetscInt *lid_cprowID_1 = NULL; 826d529f056Smarkadams4 NState *lid_state; 827d529f056Smarkadams4 Vec ghost_par_orig2; 828d529f056Smarkadams4 PetscMPIInt rank; 829d529f056Smarkadams4 830d529f056Smarkadams4 PetscFunctionBegin; 831d529f056Smarkadams4 PetscCall(PetscObjectGetComm((PetscObject)Gmat_2, &comm)); 832d529f056Smarkadams4 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 833d529f056Smarkadams4 PetscCall(MatGetOwnershipRange(Gmat_1, &my0, &Iend)); 834d529f056Smarkadams4 835d529f056Smarkadams4 /* get submatrices */ 836d529f056Smarkadams4 PetscCall(PetscStrbeginswith(((PetscObject)Gmat_1)->type_name, MATMPIAIJ, &isMPI)); 837d529f056Smarkadams4 PetscCall(PetscInfo(pc, "isMPI = %s\n", isMPI ? "yes" : "no")); 838d529f056Smarkadams4 PetscCall(PetscMalloc3(nloc, &lid_state, nloc, &lid_parent_gid, nloc, &lid_cprowID_1)); 839d529f056Smarkadams4 for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1; 840d529f056Smarkadams4 if (isMPI) { 841d529f056Smarkadams4 /* grab matrix objects */ 842d529f056Smarkadams4 mpimat_2 = (Mat_MPIAIJ *)Gmat_2->data; 843d529f056Smarkadams4 mpimat_1 = (Mat_MPIAIJ *)Gmat_1->data; 844d529f056Smarkadams4 matA_1 = (Mat_SeqAIJ *)mpimat_1->A->data; 845d529f056Smarkadams4 matB_1 = (Mat_SeqAIJ *)mpimat_1->B->data; 846d529f056Smarkadams4 847d529f056Smarkadams4 /* force compressed row storage for B matrix in AuxMat */ 848d529f056Smarkadams4 PetscCall(MatCheckCompressedRow(mpimat_1->B, matB_1->nonzerorowcnt, &matB_1->compressedrow, matB_1->i, Gmat_1->rmap->n, -1.0)); 849d529f056Smarkadams4 for (ix = 0; ix < matB_1->compressedrow.nrows; ix++) { 850d529f056Smarkadams4 PetscInt lid = matB_1->compressedrow.rindex[ix]; 851e0b7e82fSBarry Smith 852d529f056Smarkadams4 PetscCheck(lid <= nloc && lid >= -1, PETSC_COMM_SELF, PETSC_ERR_USER, "lid %d out of range. nloc = %d", (int)lid, (int)nloc); 853d529f056Smarkadams4 if (lid != -1) lid_cprowID_1[lid] = ix; 854d529f056Smarkadams4 } 855d529f056Smarkadams4 } else { 856d529f056Smarkadams4 PetscBool isAIJ; 857e0b7e82fSBarry Smith 858d529f056Smarkadams4 PetscCall(PetscStrbeginswith(((PetscObject)Gmat_1)->type_name, MATSEQAIJ, &isAIJ)); 859d529f056Smarkadams4 PetscCheck(isAIJ, PETSC_COMM_SELF, PETSC_ERR_USER, "Require AIJ matrix."); 860d529f056Smarkadams4 matA_1 = (Mat_SeqAIJ *)Gmat_1->data; 861d529f056Smarkadams4 } 862d529f056Smarkadams4 if (nloc > 0) { PetscCheck(!matB_1 || matB_1->compressedrow.use, PETSC_COMM_SELF, PETSC_ERR_PLIB, "matB_1 && !matB_1->compressedrow.use: PETSc bug???"); } 863d529f056Smarkadams4 /* get state of locals and selected gid for deleted */ 864d529f056Smarkadams4 for (lid = 0; lid < nloc; lid++) { 865d529f056Smarkadams4 lid_parent_gid[lid] = -1.0; 866d529f056Smarkadams4 lid_state[lid] = DELETED; 867d529f056Smarkadams4 } 868d529f056Smarkadams4 869d529f056Smarkadams4 /* set lid_state */ 870d529f056Smarkadams4 for (lid = 0; lid < nloc; lid++) { 871d529f056Smarkadams4 PetscCDIntNd *pos; 872e0b7e82fSBarry Smith 873d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos)); 874d529f056Smarkadams4 if (pos) { 875d529f056Smarkadams4 PetscInt gid1; 876d529f056Smarkadams4 877d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gid1)); 878d529f056Smarkadams4 PetscCheck(gid1 == lid + my0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "gid1 %d != lid %d + my0 %d", (int)gid1, (int)lid, (int)my0); 879d529f056Smarkadams4 lid_state[lid] = gid1; 880d529f056Smarkadams4 } 881d529f056Smarkadams4 } 882d529f056Smarkadams4 883d529f056Smarkadams4 /* map local to selected local, DELETED means a ghost owns it */ 884d529f056Smarkadams4 for (lid = kk = 0; lid < nloc; lid++) { 885d529f056Smarkadams4 NState state = lid_state[lid]; 886d529f056Smarkadams4 if (IS_SELECTED(state)) { 887d529f056Smarkadams4 PetscCDIntNd *pos; 888e0b7e82fSBarry Smith 889d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos)); 890d529f056Smarkadams4 while (pos) { 891d529f056Smarkadams4 PetscInt gid1; 892e0b7e82fSBarry Smith 893d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gid1)); 894d529f056Smarkadams4 PetscCall(PetscCDGetNextPos(aggs_2, lid, &pos)); 895d529f056Smarkadams4 if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1 - my0] = (PetscScalar)(lid + my0); 896d529f056Smarkadams4 } 897d529f056Smarkadams4 } 898d529f056Smarkadams4 } 899d529f056Smarkadams4 /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */ 900d529f056Smarkadams4 if (isMPI) { 901d529f056Smarkadams4 Vec tempVec; 902d529f056Smarkadams4 /* get 'cpcol_1_state' */ 903d529f056Smarkadams4 PetscCall(MatCreateVecs(Gmat_1, &tempVec, NULL)); 904d529f056Smarkadams4 for (kk = 0, j = my0; kk < nloc; kk++, j++) { 905d529f056Smarkadams4 PetscScalar v = (PetscScalar)lid_state[kk]; 906e0b7e82fSBarry Smith 907d529f056Smarkadams4 PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES)); 908d529f056Smarkadams4 } 909d529f056Smarkadams4 PetscCall(VecAssemblyBegin(tempVec)); 910d529f056Smarkadams4 PetscCall(VecAssemblyEnd(tempVec)); 911d529f056Smarkadams4 PetscCall(VecScatterBegin(mpimat_1->Mvctx, tempVec, mpimat_1->lvec, INSERT_VALUES, SCATTER_FORWARD)); 912d529f056Smarkadams4 PetscCall(VecScatterEnd(mpimat_1->Mvctx, tempVec, mpimat_1->lvec, INSERT_VALUES, SCATTER_FORWARD)); 913d529f056Smarkadams4 PetscCall(VecGetArray(mpimat_1->lvec, &cpcol_1_state)); 914d529f056Smarkadams4 /* get 'cpcol_2_state' */ 915d529f056Smarkadams4 PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, mpimat_2->lvec, INSERT_VALUES, SCATTER_FORWARD)); 916d529f056Smarkadams4 PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, mpimat_2->lvec, INSERT_VALUES, SCATTER_FORWARD)); 917d529f056Smarkadams4 PetscCall(VecGetArray(mpimat_2->lvec, &cpcol_2_state)); 918d529f056Smarkadams4 /* get 'cpcol_2_par_orig' */ 919d529f056Smarkadams4 for (kk = 0, j = my0; kk < nloc; kk++, j++) { 920d529f056Smarkadams4 PetscScalar v = (PetscScalar)lid_parent_gid[kk]; 921e0b7e82fSBarry Smith 922d529f056Smarkadams4 PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES)); 923d529f056Smarkadams4 } 924d529f056Smarkadams4 PetscCall(VecAssemblyBegin(tempVec)); 925d529f056Smarkadams4 PetscCall(VecAssemblyEnd(tempVec)); 926d529f056Smarkadams4 PetscCall(VecDuplicate(mpimat_2->lvec, &ghost_par_orig2)); 927d529f056Smarkadams4 PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghost_par_orig2, INSERT_VALUES, SCATTER_FORWARD)); 928d529f056Smarkadams4 PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghost_par_orig2, INSERT_VALUES, SCATTER_FORWARD)); 929d529f056Smarkadams4 PetscCall(VecGetArray(ghost_par_orig2, &cpcol_2_par_orig)); 930d529f056Smarkadams4 931d529f056Smarkadams4 PetscCall(VecDestroy(&tempVec)); 932d529f056Smarkadams4 } /* ismpi */ 933d529f056Smarkadams4 for (lid = 0; lid < nloc; lid++) { 934d529f056Smarkadams4 NState state = lid_state[lid]; 935e0b7e82fSBarry Smith 936d529f056Smarkadams4 if (IS_SELECTED(state)) { 937d529f056Smarkadams4 /* steal locals */ 938d529f056Smarkadams4 ii = matA_1->i; 939d529f056Smarkadams4 n = ii[lid + 1] - ii[lid]; 940d529f056Smarkadams4 idx = matA_1->j + ii[lid]; 941d529f056Smarkadams4 for (j = 0; j < n; j++) { 942d529f056Smarkadams4 PetscInt lidj = idx[j], sgid; 943d529f056Smarkadams4 NState statej = lid_state[lidj]; 944e0b7e82fSBarry Smith 945d529f056Smarkadams4 if (statej == DELETED && (sgid = (PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid + my0) { /* steal local */ 946d529f056Smarkadams4 lid_parent_gid[lidj] = (PetscScalar)(lid + my0); /* send this if sgid is not local */ 947d529f056Smarkadams4 if (sgid >= my0 && sgid < Iend) { /* I'm stealing this local from a local sgid */ 948d529f056Smarkadams4 PetscInt hav = 0, slid = sgid - my0, gidj = lidj + my0; 949d529f056Smarkadams4 PetscCDIntNd *pos, *last = NULL; 950e0b7e82fSBarry Smith 951d529f056Smarkadams4 /* looking for local from local so id_llist_2 works */ 952d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, slid, &pos)); 953d529f056Smarkadams4 while (pos) { 954d529f056Smarkadams4 PetscInt gid; 955d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gid)); 956d529f056Smarkadams4 if (gid == gidj) { 957d529f056Smarkadams4 PetscCheck(last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "last cannot be null"); 958d529f056Smarkadams4 PetscCall(PetscCDRemoveNextNode(aggs_2, slid, last)); 959d529f056Smarkadams4 PetscCall(PetscCDAppendNode(aggs_2, lid, pos)); 960d529f056Smarkadams4 hav = 1; 961d529f056Smarkadams4 break; 962d529f056Smarkadams4 } else last = pos; 963d529f056Smarkadams4 PetscCall(PetscCDGetNextPos(aggs_2, slid, &pos)); 964d529f056Smarkadams4 } 965d529f056Smarkadams4 if (hav != 1) { 966d529f056Smarkadams4 PetscCheck(hav, PETSC_COMM_SELF, PETSC_ERR_PLIB, "failed to find adj in 'selected' lists - structurally unsymmetric matrix"); 967d529f056Smarkadams4 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "found node %d times???", (int)hav); 968d529f056Smarkadams4 } 969d529f056Smarkadams4 } else { /* I'm stealing this local, owned by a ghost */ 970d529f056Smarkadams4 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.", 971d529f056Smarkadams4 ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : ""); 972d529f056Smarkadams4 PetscCall(PetscCDAppendID(aggs_2, lid, lidj + my0)); 973d529f056Smarkadams4 } 974d529f056Smarkadams4 } 975d529f056Smarkadams4 } /* local neighbors */ 976d529f056Smarkadams4 } else if (state == DELETED /* && lid_cprowID_1 */) { 977d529f056Smarkadams4 PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]); 978e0b7e82fSBarry Smith 979d529f056Smarkadams4 /* see if I have a selected ghost neighbor that will steal me */ 980d529f056Smarkadams4 if ((ix = lid_cprowID_1[lid]) != -1) { 981d529f056Smarkadams4 ii = matB_1->compressedrow.i; 982d529f056Smarkadams4 n = ii[ix + 1] - ii[ix]; 983d529f056Smarkadams4 idx = matB_1->j + ii[ix]; 984d529f056Smarkadams4 for (j = 0; j < n; j++) { 985d529f056Smarkadams4 PetscInt cpid = idx[j]; 986d529f056Smarkadams4 NState statej = (NState)PetscRealPart(cpcol_1_state[cpid]); 987e0b7e82fSBarry Smith 988d529f056Smarkadams4 if (IS_SELECTED(statej) && sgidold != (PetscInt)statej) { /* ghost will steal this, remove from my list */ 989d529f056Smarkadams4 lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */ 990d529f056Smarkadams4 if (sgidold >= my0 && sgidold < Iend) { /* this was mine */ 991d529f056Smarkadams4 PetscInt hav = 0, oldslidj = sgidold - my0; 992d529f056Smarkadams4 PetscCDIntNd *pos, *last = NULL; 993e0b7e82fSBarry Smith 994d529f056Smarkadams4 /* remove from 'oldslidj' list */ 995d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, oldslidj, &pos)); 996d529f056Smarkadams4 while (pos) { 997d529f056Smarkadams4 PetscInt gid; 998e0b7e82fSBarry Smith 999d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gid)); 1000d529f056Smarkadams4 if (lid + my0 == gid) { 1001d529f056Smarkadams4 /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove lid from oldslidj list *\/ */ 1002d529f056Smarkadams4 PetscCheck(last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "last cannot be null"); 1003d529f056Smarkadams4 PetscCall(PetscCDRemoveNextNode(aggs_2, oldslidj, last)); 1004d529f056Smarkadams4 /* ghost (PetscScalar)statej will add this later */ 1005d529f056Smarkadams4 hav = 1; 1006d529f056Smarkadams4 break; 1007d529f056Smarkadams4 } else last = pos; 1008d529f056Smarkadams4 PetscCall(PetscCDGetNextPos(aggs_2, oldslidj, &pos)); 1009d529f056Smarkadams4 } 1010d529f056Smarkadams4 if (hav != 1) { 1011d529f056Smarkadams4 PetscCheck(hav, PETSC_COMM_SELF, PETSC_ERR_PLIB, "failed to find (hav=%d) adj in 'selected' lists - structurally unsymmetric matrix", (int)hav); 1012d529f056Smarkadams4 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "found node %d times???", (int)hav); 1013d529f056Smarkadams4 } 1014d529f056Smarkadams4 } else { 1015d529f056Smarkadams4 /* TODO: ghosts remove this later */ 1016d529f056Smarkadams4 } 1017d529f056Smarkadams4 } 1018d529f056Smarkadams4 } 1019d529f056Smarkadams4 } 1020d529f056Smarkadams4 } /* selected/deleted */ 1021d529f056Smarkadams4 } /* node loop */ 1022d529f056Smarkadams4 1023d529f056Smarkadams4 if (isMPI) { 1024d529f056Smarkadams4 PetscScalar *cpcol_2_parent, *cpcol_2_gid; 1025d529f056Smarkadams4 Vec tempVec, ghostgids2, ghostparents2; 1026d529f056Smarkadams4 PetscInt cpid, nghost_2; 1027d529f056Smarkadams4 PCGAMGHashTable gid_cpid; 1028d529f056Smarkadams4 1029d529f056Smarkadams4 PetscCall(VecGetSize(mpimat_2->lvec, &nghost_2)); 1030d529f056Smarkadams4 PetscCall(MatCreateVecs(Gmat_2, &tempVec, NULL)); 1031d529f056Smarkadams4 1032d529f056Smarkadams4 /* get 'cpcol_2_parent' */ 1033d529f056Smarkadams4 for (kk = 0, j = my0; kk < nloc; kk++, j++) { PetscCall(VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES)); } 1034d529f056Smarkadams4 PetscCall(VecAssemblyBegin(tempVec)); 1035d529f056Smarkadams4 PetscCall(VecAssemblyEnd(tempVec)); 1036d529f056Smarkadams4 PetscCall(VecDuplicate(mpimat_2->lvec, &ghostparents2)); 1037d529f056Smarkadams4 PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghostparents2, INSERT_VALUES, SCATTER_FORWARD)); 1038d529f056Smarkadams4 PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghostparents2, INSERT_VALUES, SCATTER_FORWARD)); 1039d529f056Smarkadams4 PetscCall(VecGetArray(ghostparents2, &cpcol_2_parent)); 1040d529f056Smarkadams4 1041d529f056Smarkadams4 /* get 'cpcol_2_gid' */ 1042d529f056Smarkadams4 for (kk = 0, j = my0; kk < nloc; kk++, j++) { 1043d529f056Smarkadams4 PetscScalar v = (PetscScalar)j; 1044e0b7e82fSBarry Smith 1045d529f056Smarkadams4 PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES)); 1046d529f056Smarkadams4 } 1047d529f056Smarkadams4 PetscCall(VecAssemblyBegin(tempVec)); 1048d529f056Smarkadams4 PetscCall(VecAssemblyEnd(tempVec)); 1049d529f056Smarkadams4 PetscCall(VecDuplicate(mpimat_2->lvec, &ghostgids2)); 1050d529f056Smarkadams4 PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghostgids2, INSERT_VALUES, SCATTER_FORWARD)); 1051d529f056Smarkadams4 PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghostgids2, INSERT_VALUES, SCATTER_FORWARD)); 1052d529f056Smarkadams4 PetscCall(VecGetArray(ghostgids2, &cpcol_2_gid)); 1053d529f056Smarkadams4 PetscCall(VecDestroy(&tempVec)); 1054d529f056Smarkadams4 1055d529f056Smarkadams4 /* look for deleted ghosts and add to table */ 1056d529f056Smarkadams4 PetscCall(PCGAMGHashTableCreate(2 * nghost_2 + 1, &gid_cpid)); 1057d529f056Smarkadams4 for (cpid = 0; cpid < nghost_2; cpid++) { 1058d529f056Smarkadams4 NState state = (NState)PetscRealPart(cpcol_2_state[cpid]); 1059e0b7e82fSBarry Smith 1060d529f056Smarkadams4 if (state == DELETED) { 1061d529f056Smarkadams4 PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 1062d529f056Smarkadams4 PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]); 1063e0b7e82fSBarry Smith 1064d529f056Smarkadams4 if (sgid_old == -1 && sgid_new != -1) { 1065d529f056Smarkadams4 PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 1066e0b7e82fSBarry Smith 1067d529f056Smarkadams4 PetscCall(PCGAMGHashTableAdd(&gid_cpid, gid, cpid)); 1068d529f056Smarkadams4 } 1069d529f056Smarkadams4 } 1070d529f056Smarkadams4 } 1071d529f056Smarkadams4 1072d529f056Smarkadams4 /* look for deleted ghosts and see if they moved - remove it */ 1073d529f056Smarkadams4 for (lid = 0; lid < nloc; lid++) { 1074d529f056Smarkadams4 NState state = lid_state[lid]; 1075d529f056Smarkadams4 if (IS_SELECTED(state)) { 1076d529f056Smarkadams4 PetscCDIntNd *pos, *last = NULL; 1077d529f056Smarkadams4 /* look for deleted ghosts and see if they moved */ 1078d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos)); 1079d529f056Smarkadams4 while (pos) { 1080d529f056Smarkadams4 PetscInt gid; 1081d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gid)); 1082d529f056Smarkadams4 1083d529f056Smarkadams4 if (gid < my0 || gid >= Iend) { 1084d529f056Smarkadams4 PetscCall(PCGAMGHashTableFind(&gid_cpid, gid, &cpid)); 1085d529f056Smarkadams4 if (cpid != -1) { 1086d529f056Smarkadams4 /* a moved ghost - */ 1087d529f056Smarkadams4 /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove 'flid' from list *\/ */ 1088d529f056Smarkadams4 PetscCall(PetscCDRemoveNextNode(aggs_2, lid, last)); 1089d529f056Smarkadams4 } else last = pos; 1090d529f056Smarkadams4 } else last = pos; 1091d529f056Smarkadams4 1092d529f056Smarkadams4 PetscCall(PetscCDGetNextPos(aggs_2, lid, &pos)); 1093d529f056Smarkadams4 } /* loop over list of deleted */ 1094d529f056Smarkadams4 } /* selected */ 1095d529f056Smarkadams4 } 1096d529f056Smarkadams4 PetscCall(PCGAMGHashTableDestroy(&gid_cpid)); 1097d529f056Smarkadams4 1098d529f056Smarkadams4 /* look at ghosts, see if they changed - and it */ 1099d529f056Smarkadams4 for (cpid = 0; cpid < nghost_2; cpid++) { 1100d529f056Smarkadams4 PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 1101e0b7e82fSBarry Smith 1102d529f056Smarkadams4 if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */ 1103d529f056Smarkadams4 PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 1104d529f056Smarkadams4 PetscInt slid_new = sgid_new - my0, hav = 0; 1105d529f056Smarkadams4 PetscCDIntNd *pos; 1106d529f056Smarkadams4 1107d529f056Smarkadams4 /* search for this gid to see if I have it */ 1108d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, slid_new, &pos)); 1109d529f056Smarkadams4 while (pos) { 1110d529f056Smarkadams4 PetscInt gidj; 1111e0b7e82fSBarry Smith 1112d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gidj)); 1113d529f056Smarkadams4 PetscCall(PetscCDGetNextPos(aggs_2, slid_new, &pos)); 1114d529f056Smarkadams4 1115d529f056Smarkadams4 if (gidj == gid) { 1116d529f056Smarkadams4 hav = 1; 1117d529f056Smarkadams4 break; 1118d529f056Smarkadams4 } 1119d529f056Smarkadams4 } 1120d529f056Smarkadams4 if (hav != 1) { 1121d529f056Smarkadams4 /* insert 'flidj' into head of llist */ 1122d529f056Smarkadams4 PetscCall(PetscCDAppendID(aggs_2, slid_new, gid)); 1123d529f056Smarkadams4 } 1124d529f056Smarkadams4 } 1125d529f056Smarkadams4 } 1126d529f056Smarkadams4 PetscCall(VecRestoreArray(mpimat_1->lvec, &cpcol_1_state)); 1127d529f056Smarkadams4 PetscCall(VecRestoreArray(mpimat_2->lvec, &cpcol_2_state)); 1128d529f056Smarkadams4 PetscCall(VecRestoreArray(ghostparents2, &cpcol_2_parent)); 1129d529f056Smarkadams4 PetscCall(VecRestoreArray(ghostgids2, &cpcol_2_gid)); 1130d529f056Smarkadams4 PetscCall(VecDestroy(&ghostgids2)); 1131d529f056Smarkadams4 PetscCall(VecDestroy(&ghostparents2)); 1132d529f056Smarkadams4 PetscCall(VecDestroy(&ghost_par_orig2)); 1133d529f056Smarkadams4 } 1134d529f056Smarkadams4 PetscCall(PetscFree3(lid_state, lid_parent_gid, lid_cprowID_1)); 1135d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 1136d529f056Smarkadams4 } 1137d529f056Smarkadams4 1138c8b0795cSMark F. Adams /* 1139bae903cbSmarkadams4 PCGAMGCoarsen_AGG - supports squaring the graph (deprecated) and new graph for 1140bae903cbSmarkadams4 communication of QR data used with HEM and MISk coarsening 1141c8b0795cSMark F. Adams 1142c8b0795cSMark F. Adams Input Parameter: 1143e0940f08SMark F. Adams . a_pc - this 1144bae903cbSmarkadams4 1145e0940f08SMark F. Adams Input/Output Parameter: 1146bae903cbSmarkadams4 . a_Gmat1 - graph to coarsen (in), graph off processor edges for QR gather scatter (out) 1147bae903cbSmarkadams4 1148c8b0795cSMark F. Adams Output Parameter: 11490cbbd2e1SMark F. Adams . agg_lists - list of aggregates 1150bae903cbSmarkadams4 1151c8b0795cSMark F. Adams */ 1152d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc, Mat *a_Gmat1, PetscCoarsenData **agg_lists) 1153d71ae5a4SJacob Faibussowitsch { 1154e0940f08SMark F. Adams PC_MG *mg = (PC_MG *)a_pc->data; 1155c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1156c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 11578926f930SMark Adams Mat Gmat2, Gmat1 = *a_Gmat1; /* aggressive graph */ 11580cbbd2e1SMark F. Adams IS perm; 1159bae903cbSmarkadams4 PetscInt Istart, Iend, Ii, nloc, bs, nn; 1160bae903cbSmarkadams4 PetscInt *permute, *degree; 1161c8b0795cSMark F. Adams PetscBool *bIndexSet; 1162aea53286SMark Adams PetscReal hashfact; 1163e2c00dcbSBarry Smith PetscInt iSwapIndex; 11643b50add6SMark Adams PetscRandom random; 1165d529f056Smarkadams4 MPI_Comm comm; 1166c8b0795cSMark F. Adams 1167c8b0795cSMark F. Adams PetscFunctionBegin; 1168d529f056Smarkadams4 PetscCall(PetscObjectGetComm((PetscObject)Gmat1, &comm)); 1169849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0)); 1170bae903cbSmarkadams4 PetscCall(MatGetLocalSize(Gmat1, &nn, NULL)); 11719566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Gmat1, &bs)); 117263a3b9bcSJacob Faibussowitsch PetscCheck(bs == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bs %" PetscInt_FMT " must be 1", bs); 1173bae903cbSmarkadams4 nloc = nn / bs; 11745cfd4bd9SMark Adams /* get MIS aggs - randomize */ 1175bae903cbSmarkadams4 PetscCall(PetscMalloc2(nloc, &permute, nloc, °ree)); 11769566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nloc, &bIndexSet)); 11776e09b0e3SStefano Zampini for (Ii = 0; Ii < nloc; Ii++) permute[Ii] = Ii; 11789566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &random)); 11799566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Gmat1, &Istart, &Iend)); 1180c8b0795cSMark F. Adams for (Ii = 0; Ii < nloc; Ii++) { 1181bae903cbSmarkadams4 PetscInt nc; 1182e0b7e82fSBarry Smith 1183bae903cbSmarkadams4 PetscCall(MatGetRow(Gmat1, Istart + Ii, &nc, NULL, NULL)); 1184bae903cbSmarkadams4 degree[Ii] = nc; 1185bae903cbSmarkadams4 PetscCall(MatRestoreRow(Gmat1, Istart + Ii, &nc, NULL, NULL)); 1186bae903cbSmarkadams4 } 1187bae903cbSmarkadams4 for (Ii = 0; Ii < nloc; Ii++) { 11889566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(random, &hashfact)); 1189aea53286SMark Adams iSwapIndex = (PetscInt)(hashfact * nloc) % nloc; 1190c8b0795cSMark F. Adams if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) { 1191c8b0795cSMark F. Adams PetscInt iTemp = permute[iSwapIndex]; 1192e0b7e82fSBarry Smith 1193c8b0795cSMark F. Adams permute[iSwapIndex] = permute[Ii]; 1194c8b0795cSMark F. Adams permute[Ii] = iTemp; 1195bae903cbSmarkadams4 iTemp = degree[iSwapIndex]; 1196bae903cbSmarkadams4 degree[iSwapIndex] = degree[Ii]; 1197bae903cbSmarkadams4 degree[Ii] = iTemp; 1198c8b0795cSMark F. Adams bIndexSet[iSwapIndex] = PETSC_TRUE; 1199c8b0795cSMark F. Adams } 1200c8b0795cSMark F. Adams } 1201d529f056Smarkadams4 // apply minimum degree ordering -- NEW 1202d529f056Smarkadams4 if (pc_gamg_agg->use_minimum_degree_ordering) { PetscCall(PetscSortIntWithArray(nloc, degree, permute)); } 12039566063dSJacob Faibussowitsch PetscCall(PetscFree(bIndexSet)); 12049566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&random)); 12059566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm)); 1206849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0)); 1207d529f056Smarkadams4 // square graph 1208d529f056Smarkadams4 if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels && pc_gamg_agg->use_aggressive_square_graph) { 1209d529f056Smarkadams4 PetscCall(PCGAMGSquareGraph_GAMG(a_pc, Gmat1, &Gmat2)); 1210d529f056Smarkadams4 } else Gmat2 = Gmat1; 1211d529f056Smarkadams4 // switch to old MIS-1 for square graph 1212d529f056Smarkadams4 if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels) { 1213d529f056Smarkadams4 if (!pc_gamg_agg->use_aggressive_square_graph) PetscCall(MatCoarsenMISKSetDistance(pc_gamg_agg->crs, pc_gamg_agg->aggressive_mis_k)); // hardwire to MIS-2 1214d529f056Smarkadams4 else PetscCall(MatCoarsenSetType(pc_gamg_agg->crs, MATCOARSENMIS)); // old MIS -- side effect 1215d529f056Smarkadams4 } else if (pc_gamg_agg->use_aggressive_square_graph && pc_gamg_agg->aggressive_coarsening_levels > 0) { // we reset the MIS 1216d529f056Smarkadams4 const char *prefix; 1217e0b7e82fSBarry Smith 1218d529f056Smarkadams4 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)a_pc, &prefix)); 1219d529f056Smarkadams4 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc_gamg_agg->crs, prefix)); 1220d529f056Smarkadams4 PetscCall(MatCoarsenSetFromOptions(pc_gamg_agg->crs)); // get the default back on non-aggressive levels when square graph switched to old MIS 1221d529f056Smarkadams4 } 1222d529f056Smarkadams4 PetscCall(MatCoarsenSetAdjacency(pc_gamg_agg->crs, Gmat2)); 12232d776b49SBarry Smith PetscCall(MatCoarsenSetStrictAggs(pc_gamg_agg->crs, PETSC_TRUE)); 1224d529f056Smarkadams4 PetscCall(MatCoarsenSetGreedyOrdering(pc_gamg_agg->crs, perm)); 12252d776b49SBarry Smith PetscCall(MatCoarsenApply(pc_gamg_agg->crs)); 12262d776b49SBarry Smith PetscCall(MatCoarsenGetData(pc_gamg_agg->crs, agg_lists)); /* output */ 1227b43b03e9SMark F. Adams 12289566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 1229bae903cbSmarkadams4 PetscCall(PetscFree2(permute, degree)); 1230849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0)); 12319f3f12c8SMark F. Adams 12329c15c1aeSMark Adams if (Gmat2 != Gmat1) { // square graph, we need ghosts for selected 1233d529f056Smarkadams4 PetscCoarsenData *llist = *agg_lists; 1234e0b7e82fSBarry Smith 1235d529f056Smarkadams4 PetscCall(fixAggregatesWithSquare(a_pc, Gmat2, Gmat1, *agg_lists)); 1236d529f056Smarkadams4 PetscCall(MatDestroy(&Gmat1)); 1237d529f056Smarkadams4 *a_Gmat1 = Gmat2; /* output */ 12388926f930SMark Adams PetscCall(PetscCDSetMat(llist, *a_Gmat1)); /* Need a graph with ghosts here */ 12390cbbd2e1SMark F. Adams } 1240849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0)); 12413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1242c8b0795cSMark F. Adams } 1243c8b0795cSMark F. Adams 1244c8b0795cSMark F. Adams /* 1245e0b7e82fSBarry Smith PCGAMGConstructProlongator_AGG 1246c8b0795cSMark F. Adams 1247c8b0795cSMark F. Adams Input Parameter: 1248c8b0795cSMark F. Adams . pc - this 1249c8b0795cSMark F. Adams . Amat - matrix on this fine level 1250c8b0795cSMark F. Adams . Graph - used to get ghost data for nodes in 12510cbbd2e1SMark F. Adams . agg_lists - list of aggregates 1252c8b0795cSMark F. Adams Output Parameter: 1253c8b0795cSMark F. Adams . a_P_out - prolongation operator to the next level 1254c8b0795cSMark F. Adams */ 1255e0b7e82fSBarry Smith static PetscErrorCode PCGAMGConstructProlongator_AGG(PC pc, Mat Amat, PetscCoarsenData *agg_lists, Mat *a_P_out) 1256d71ae5a4SJacob Faibussowitsch { 12572e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 12582e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 12594a2f8832SBarry Smith const PetscInt col_bs = pc_gamg->data_cell_cols; 1260c8b0795cSMark F. Adams PetscInt Istart, Iend, nloc, ii, jj, kk, my0, nLocalSelected, bs; 12618926f930SMark Adams Mat Gmat, Prol; 1262d5d25401SBarry Smith PetscMPIInt size; 12633b4367a7SBarry Smith MPI_Comm comm; 1264c8b0795cSMark F. Adams PetscReal *data_w_ghost; 1265c8b0795cSMark F. Adams PetscInt myCrs0, nbnodes = 0, *flid_fgid; 1266d9558ea9SBarry Smith MatType mtype; 12672e68589bSMark F. Adams 12682e68589bSMark F. Adams PetscFunctionBegin; 12699566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm)); 127008401ef6SPierre Jolivet PetscCheck(col_bs >= 1, comm, PETSC_ERR_PLIB, "Column bs cannot be less than 1"); 1271849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0)); 12729566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 12739566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Amat, &Istart, &Iend)); 12749566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Amat, &bs)); 12759371c9d4SSatish Balay nloc = (Iend - Istart) / bs; 12769371c9d4SSatish Balay my0 = Istart / bs; 127763a3b9bcSJacob 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); 1278d8b4a066SPierre Jolivet PetscCall(PetscCDGetMat(agg_lists, &Gmat)); // get auxiliary matrix for ghost edges for size > 1 12792e68589bSMark F. Adams 12802e68589bSMark F. Adams /* get 'nLocalSelected' */ 12810cbbd2e1SMark F. Adams for (ii = 0, nLocalSelected = 0; ii < nloc; ii++) { 1282e78576d6SMark F. Adams PetscBool ise; 1283e0b7e82fSBarry Smith 12840cbbd2e1SMark F. Adams /* filter out singletons 0 or 1? */ 12855e62d202SMark Adams PetscCall(PetscCDIsEmptyAt(agg_lists, ii, &ise)); 1286e78576d6SMark F. Adams if (!ise) nLocalSelected++; 12872e68589bSMark F. Adams } 12882e68589bSMark F. Adams 12892e68589bSMark F. Adams /* create prolongator, create P matrix */ 12909566063dSJacob Faibussowitsch PetscCall(MatGetType(Amat, &mtype)); 12919566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &Prol)); 12929566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Prol, nloc * bs, nLocalSelected * col_bs, PETSC_DETERMINE, PETSC_DETERMINE)); 129344eff04eSMark Adams PetscCall(MatSetBlockSizes(Prol, bs, col_bs)); // should this be before MatSetSizes? 12949566063dSJacob Faibussowitsch PetscCall(MatSetType(Prol, mtype)); 12953742c8caSstefanozampini #if PetscDefined(HAVE_DEVICE) 12963742c8caSstefanozampini PetscBool flg; 12973742c8caSstefanozampini PetscCall(MatBoundToCPU(Amat, &flg)); 12983742c8caSstefanozampini PetscCall(MatBindToCPU(Prol, flg)); 12993742c8caSstefanozampini if (flg) PetscCall(MatSetBindingPropagates(Prol, PETSC_TRUE)); 13003742c8caSstefanozampini #endif 13019566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(Prol, col_bs, NULL)); 13029566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(Prol, col_bs, NULL, col_bs, NULL)); 13032e68589bSMark F. Adams 13042e68589bSMark F. Adams /* can get all points "removed" */ 13059566063dSJacob Faibussowitsch PetscCall(MatGetSize(Prol, &kk, &ii)); 13067f66b68fSMark Adams if (!ii) { 130763a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: No selected points on coarse grid\n", ((PetscObject)pc)->prefix)); 13089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Prol)); 13090298fd71SBarry Smith *a_P_out = NULL; /* out */ 1310849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0)); 13113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13122e68589bSMark F. Adams } 131363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: New grid %" PetscInt_FMT " nodes\n", ((PetscObject)pc)->prefix, ii / col_bs)); 13149566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk)); 13150cbbd2e1SMark F. Adams 131663a3b9bcSJacob 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); 1317c8b0795cSMark F. Adams myCrs0 = myCrs0 / col_bs; 131863a3b9bcSJacob 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); 13192e68589bSMark F. Adams 13202e68589bSMark F. Adams /* create global vector of data in 'data_w_ghost' */ 1321849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0)); 1322bae903cbSmarkadams4 if (size > 1) { /* get ghost null space data */ 13232e68589bSMark F. Adams PetscReal *tmp_gdata, *tmp_ldata, *tp2; 1324e0b7e82fSBarry Smith 13259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc, &tmp_ldata)); 13264a2f8832SBarry Smith for (jj = 0; jj < col_bs; jj++) { 1327c8b0795cSMark F. Adams for (kk = 0; kk < bs; kk++) { 1328a2f3521dSMark F. Adams PetscInt ii, stride; 13298e3a54c0SPierre Jolivet const PetscReal *tp = PetscSafePointerPlusOffset(pc_gamg->data, jj * bs * nloc + kk); 1330e0b7e82fSBarry Smith 13312fa5cd67SKarl Rupp for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp; 13322fa5cd67SKarl Rupp 13339566063dSJacob Faibussowitsch PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata)); 1334a2f3521dSMark F. Adams 1335bae903cbSmarkadams4 if (!jj && !kk) { /* now I know how many total nodes - allocate TODO: move below and do in one 'col_bs' call */ 13369566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(stride * bs * col_bs, &data_w_ghost)); 1337a2f3521dSMark F. Adams nbnodes = bs * stride; 13382e68589bSMark F. Adams } 13398e3a54c0SPierre Jolivet tp2 = PetscSafePointerPlusOffset(data_w_ghost, jj * bs * stride + kk); 13402fa5cd67SKarl Rupp for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii]; 13419566063dSJacob Faibussowitsch PetscCall(PetscFree(tmp_gdata)); 13422e68589bSMark F. Adams } 13432e68589bSMark F. Adams } 13449566063dSJacob Faibussowitsch PetscCall(PetscFree(tmp_ldata)); 1345806fa848SBarry Smith } else { 1346c8b0795cSMark F. Adams nbnodes = bs * nloc; 1347c8b0795cSMark F. Adams data_w_ghost = (PetscReal *)pc_gamg->data; 13482e68589bSMark F. Adams } 13492e68589bSMark F. Adams 1350bae903cbSmarkadams4 /* get 'flid_fgid' TODO - move up to get 'stride' and do get null space data above in one step (jj loop) */ 1351c5df96a5SBarry Smith if (size > 1) { 13522e68589bSMark F. Adams PetscReal *fid_glid_loc, *fiddata; 1353a2f3521dSMark F. Adams PetscInt stride; 13542e68589bSMark F. Adams 13559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc, &fid_glid_loc)); 13562e68589bSMark F. Adams for (kk = 0; kk < nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0 + kk); 13579566063dSJacob Faibussowitsch PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata)); 1358bae903cbSmarkadams4 PetscCall(PetscMalloc1(stride, &flid_fgid)); /* copy real data to in */ 1359a2f3521dSMark F. Adams for (kk = 0; kk < stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk]; 13609566063dSJacob Faibussowitsch PetscCall(PetscFree(fiddata)); 1361a2f3521dSMark F. Adams 136263a3b9bcSJacob Faibussowitsch PetscCheck(stride == nbnodes / bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "stride %" PetscInt_FMT " != nbnodes %" PetscInt_FMT "/bs %" PetscInt_FMT, stride, nbnodes, bs); 13639566063dSJacob Faibussowitsch PetscCall(PetscFree(fid_glid_loc)); 1364806fa848SBarry Smith } else { 13659566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc, &flid_fgid)); 13662e68589bSMark F. Adams for (kk = 0; kk < nloc; kk++) flid_fgid[kk] = my0 + kk; 13672e68589bSMark F. Adams } 1368849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0)); 1369bae903cbSmarkadams4 /* get P0 */ 1370849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0)); 1371c8b0795cSMark F. Adams { 13720298fd71SBarry Smith PetscReal *data_out = NULL; 1373e0b7e82fSBarry Smith 13749566063dSJacob Faibussowitsch PetscCall(formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes, data_w_ghost, flid_fgid, &data_out, Prol)); 13759566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 13762fa5cd67SKarl Rupp 1377c8b0795cSMark F. Adams pc_gamg->data = data_out; 13784a2f8832SBarry Smith pc_gamg->data_cell_rows = col_bs; 13794a2f8832SBarry Smith pc_gamg->data_sz = col_bs * col_bs * nLocalSelected; 1380c8b0795cSMark F. Adams } 1381849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0)); 138248a46eb9SPierre Jolivet if (size > 1) PetscCall(PetscFree(data_w_ghost)); 13839566063dSJacob Faibussowitsch PetscCall(PetscFree(flid_fgid)); 13842e68589bSMark F. Adams 1385c8b0795cSMark F. Adams *a_P_out = Prol; /* out */ 1386e0b7e82fSBarry Smith PetscCall(MatViewFromOptions(Prol, NULL, "-pc_gamg_agg_view_initial_prolongation")); 1387ed4e82eaSPeter Brune 1388849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0)); 13893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1390c8b0795cSMark F. Adams } 1391c8b0795cSMark F. Adams 1392c8b0795cSMark F. Adams /* 1393e0b7e82fSBarry Smith PCGAMGOptimizeProlongator_AGG - given the initial prolongator optimizes it by smoothed aggregation pc_gamg_agg->nsmooths times 1394c8b0795cSMark F. Adams 1395c8b0795cSMark F. Adams Input Parameter: 1396c8b0795cSMark F. Adams . pc - this 1397c8b0795cSMark F. Adams . Amat - matrix on this fine level 1398c8b0795cSMark F. Adams In/Output Parameter: 13992a808120SBarry Smith . a_P - prolongation operator to the next level 1400c8b0795cSMark F. Adams */ 1401e0b7e82fSBarry Smith static PetscErrorCode PCGAMGOptimizeProlongator_AGG(PC pc, Mat Amat, Mat *a_P) 1402d71ae5a4SJacob Faibussowitsch { 1403c8b0795cSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 1404c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1405c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 1406c8b0795cSMark F. Adams PetscInt jj; 1407c8b0795cSMark F. Adams Mat Prol = *a_P; 14083b4367a7SBarry Smith MPI_Comm comm; 14092a808120SBarry Smith KSP eksp; 14102a808120SBarry Smith Vec bb, xx; 14112a808120SBarry Smith PC epc; 14122a808120SBarry Smith PetscReal alpha, emax, emin; 1413c8b0795cSMark F. Adams 1414c8b0795cSMark F. Adams PetscFunctionBegin; 14159566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm)); 1416849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0)); 1417c8b0795cSMark F. Adams 1418d5d25401SBarry Smith /* compute maximum singular value of operator to be used in smoother */ 14192a808120SBarry Smith if (0 < pc_gamg_agg->nsmooths) { 142018c3aa7eSMark /* get eigen estimates */ 142118c3aa7eSMark if (pc_gamg->emax > 0) { 142218c3aa7eSMark emin = pc_gamg->emin; 142318c3aa7eSMark emax = pc_gamg->emax; 142418c3aa7eSMark } else { 14250ed2132dSStefano Zampini const char *prefix; 14260ed2132dSStefano Zampini 14279566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Amat, &bb, NULL)); 14289566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Amat, &xx, NULL)); 14299566063dSJacob Faibussowitsch PetscCall(KSPSetNoisy_Private(bb)); 1430e696ed0bSMark F. Adams 14319566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm, &eksp)); 14323821be0aSBarry Smith PetscCall(KSPSetNestLevel(eksp, pc->kspnestlevel)); 14339566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 14349566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(eksp, prefix)); 14359566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(eksp, "pc_gamg_esteig_")); 143673f7197eSJed Brown { 1437b94d7dedSBarry Smith PetscBool isset, sflg; 1438e0b7e82fSBarry Smith 1439b94d7dedSBarry Smith PetscCall(MatIsSPDKnown(Amat, &isset, &sflg)); 1440b94d7dedSBarry Smith if (isset && sflg) PetscCall(KSPSetType(eksp, KSPCG)); 1441d24ecf33SMark } 14429566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(eksp, pc->erroriffailure)); 14439566063dSJacob Faibussowitsch PetscCall(KSPSetNormType(eksp, KSP_NORM_NONE)); 1444c2436cacSMark F. Adams 14459566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(eksp, PETSC_FALSE)); 14469566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(eksp, Amat, Amat)); 14472e68589bSMark F. Adams 14489566063dSJacob Faibussowitsch PetscCall(KSPGetPC(eksp, &epc)); 14499566063dSJacob Faibussowitsch PetscCall(PCSetType(epc, PCJACOBI)); /* smoother in smoothed agg. */ 1450c2436cacSMark F. Adams 1451fb842aefSJose E. Roman PetscCall(KSPSetTolerances(eksp, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT, 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 1452074aec55SMark Adams 14539566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(eksp)); 14549566063dSJacob Faibussowitsch PetscCall(KSPSetComputeSingularValues(eksp, PETSC_TRUE)); 14559566063dSJacob Faibussowitsch PetscCall(KSPSolve(eksp, bb, xx)); 14569566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(eksp, pc, xx)); 14572e68589bSMark F. Adams 14589566063dSJacob Faibussowitsch PetscCall(KSPComputeExtremeSingularValues(eksp, &emax, &emin)); 14599566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: Smooth P0: max eigen=%e min=%e PC=%s\n", ((PetscObject)pc)->prefix, (double)emax, (double)emin, PCJACOBI)); 14609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xx)); 14619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bb)); 14629566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&eksp)); 14632e68589bSMark F. Adams } 146418c3aa7eSMark if (pc_gamg->use_sa_esteig) { 146518c3aa7eSMark mg->min_eigen_DinvA[pc_gamg->current_level] = emin; 146618c3aa7eSMark mg->max_eigen_DinvA[pc_gamg->current_level] = emax; 146763a3b9bcSJacob 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)); 146818c3aa7eSMark } else { 146918c3aa7eSMark mg->min_eigen_DinvA[pc_gamg->current_level] = 0; 147018c3aa7eSMark mg->max_eigen_DinvA[pc_gamg->current_level] = 0; 147118c3aa7eSMark } 147218c3aa7eSMark } else { 147318c3aa7eSMark mg->min_eigen_DinvA[pc_gamg->current_level] = 0; 147418c3aa7eSMark mg->max_eigen_DinvA[pc_gamg->current_level] = 0; 147518c3aa7eSMark } 14762e68589bSMark F. Adams 14772a808120SBarry Smith /* smooth P0 */ 1478e0b7e82fSBarry Smith if (pc_gamg_agg->nsmooths > 0) { 14792a808120SBarry Smith Vec diag; 14802a808120SBarry Smith 1481e0b7e82fSBarry Smith /* TODO: Set a PCFailedReason and exit the building of the AMG preconditioner */ 1482e0b7e82fSBarry Smith PetscCheck(emax != 0.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Computed maximum singular value as zero"); 14832a808120SBarry Smith 1484e0b7e82fSBarry Smith PetscCall(MatCreateVecs(Amat, &diag, NULL)); 1485e0b7e82fSBarry Smith PetscCall(MatGetDiagonal(Amat, diag)); /* effectively PCJACOBI */ 1486e0b7e82fSBarry Smith PetscCall(VecReciprocal(diag)); 1487e0b7e82fSBarry Smith 1488e0b7e82fSBarry Smith for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) { 1489e0b7e82fSBarry Smith Mat tMat; 1490e0b7e82fSBarry Smith 1491e0b7e82fSBarry Smith PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0)); 1492e0b7e82fSBarry Smith /* 1493e0b7e82fSBarry Smith Smooth aggregation on the prolongator 1494e0b7e82fSBarry Smith 1495e0b7e82fSBarry Smith P_{i} := (I - 1.4/emax D^{-1}A) P_i\{i-1} 1496e0b7e82fSBarry Smith */ 14979566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0)); 1498fb842aefSJose E. Roman PetscCall(MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_CURRENT, &tMat)); 14999566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0)); 15009566063dSJacob Faibussowitsch PetscCall(MatProductClear(tMat)); 15019566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(tMat, diag, NULL)); 1502b4da3a1bSStefano Zampini 1503d5d25401SBarry Smith /* TODO: Document the 1.4 and don't hardwire it in this routine */ 1504b4ec6429SMark F. Adams alpha = -1.4 / emax; 15059566063dSJacob Faibussowitsch PetscCall(MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN)); 15069566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Prol)); 15072e68589bSMark F. Adams Prol = tMat; 1508849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0)); 15092e68589bSMark F. Adams } 1510e0b7e82fSBarry Smith PetscCall(VecDestroy(&diag)); 1511e0b7e82fSBarry Smith } 1512849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0)); 1513e0b7e82fSBarry Smith PetscCall(MatViewFromOptions(Prol, NULL, "-pc_gamg_agg_view_prolongation")); 1514c8b0795cSMark F. Adams *a_P = Prol; 15153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15162e68589bSMark F. Adams } 15172e68589bSMark F. Adams 1518a077d33dSBarry Smith /*MC 1519a077d33dSBarry Smith PCGAMGAGG - Smooth aggregation, {cite}`vanek1996algebraic`, {cite}`vanek2001convergence`, variant of PETSc's algebraic multigrid (`PCGAMG`) preconditioner 15202e68589bSMark F. Adams 1521a077d33dSBarry Smith Options Database Keys: 1522a077d33dSBarry Smith + -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation to construct prolongation 1523a077d33dSBarry Smith . -pc_gamg_aggressive_coarsening <n,default=1> - number of aggressive coarsening (MIS-2) levels from finest. 1524a077d33dSBarry Smith . -pc_gamg_aggressive_square_graph <bool,default=false> - Use square graph (A'A) or MIS-k (k=2) for aggressive coarsening 1525a077d33dSBarry Smith . -pc_gamg_mis_k_minimum_degree_ordering <bool,default=true> - Use minimum degree ordering in greedy MIS algorithm 1526a077d33dSBarry Smith . -pc_gamg_pc_gamg_asm_hem_aggs <n,default=0> - Number of HEM aggregation steps for ASM smoother 1527a077d33dSBarry Smith - -pc_gamg_aggressive_mis_k <n,default=2> - Number (k) distance in MIS coarsening (>2 is 'aggressive') 1528a077d33dSBarry Smith 1529a077d33dSBarry Smith Level: intermediate 1530a077d33dSBarry Smith 1531a077d33dSBarry Smith Notes: 1532a077d33dSBarry Smith To obtain good performance for `PCGAMG` for vector valued problems you must 1533a077d33dSBarry Smith call `MatSetBlockSize()` to indicate the number of degrees of freedom per grid point. 1534a077d33dSBarry Smith Call `MatSetNearNullSpace()` (or `PCSetCoordinates()` if solving the equations of elasticity) to indicate the near null space of the operator 1535a077d33dSBarry Smith 1536a077d33dSBarry Smith The many options for `PCMG` and `PCGAMG` such as controlling the smoothers on each level etc. also work for `PCGAMGAGG` 1537a077d33dSBarry Smith 1538a077d33dSBarry Smith .seealso: `PCGAMG`, [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCCreate()`, `PCSetType()`, 1539a077d33dSBarry Smith `MatSetBlockSize()`, `PCMGType`, `PCSetCoordinates()`, `MatSetNearNullSpace()`, `PCGAMGSetType()`, 1540a077d33dSBarry Smith `PCGAMGAGG`, `PCGAMGGEO`, `PCGAMGCLASSICAL`, `PCGAMGSetProcEqLim()`, `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRepartition()`, `PCGAMGRegister()`, 1541a077d33dSBarry Smith `PCGAMGSetReuseInterpolation()`, `PCGAMGASMSetUseAggs()`, `PCGAMGSetParallelCoarseGridSolve()`, `PCGAMGSetNlevels()`, `PCGAMGSetThreshold()`, 1542a077d33dSBarry Smith `PCGAMGGetType()`, `PCGAMGSetUseSAEstEig()` 1543a077d33dSBarry Smith M*/ 1544d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCreateGAMG_AGG(PC pc) 1545d71ae5a4SJacob Faibussowitsch { 1546c8b0795cSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 1547c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1548c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg; 15492e68589bSMark F. Adams 1550c8b0795cSMark F. Adams PetscFunctionBegin; 1551c8b0795cSMark F. Adams /* create sub context for SA */ 15524dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&pc_gamg_agg)); 1553c8b0795cSMark F. Adams pc_gamg->subctx = pc_gamg_agg; 1554c8b0795cSMark F. Adams 15551ab5ffc9SJed Brown pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG; 15569b8ffb57SJed Brown pc_gamg->ops->destroy = PCDestroy_GAMG_AGG; 1557c8b0795cSMark F. Adams /* reset does not do anything; setup not virtual */ 1558c8b0795cSMark F. Adams 1559c8b0795cSMark F. Adams /* set internal function pointers */ 15602d776b49SBarry Smith pc_gamg->ops->creategraph = PCGAMGCreateGraph_AGG; 15611ab5ffc9SJed Brown pc_gamg->ops->coarsen = PCGAMGCoarsen_AGG; 1562e0b7e82fSBarry Smith pc_gamg->ops->prolongator = PCGAMGConstructProlongator_AGG; 1563e0b7e82fSBarry Smith pc_gamg->ops->optprolongator = PCGAMGOptimizeProlongator_AGG; 15641ab5ffc9SJed Brown pc_gamg->ops->createdefaultdata = PCSetData_AGG; 15655adeb434SBarry Smith pc_gamg->ops->view = PCView_GAMG_AGG; 1566c8b0795cSMark F. Adams 1567cab9ed1eSBarry Smith pc_gamg_agg->nsmooths = 1; 1568d529f056Smarkadams4 pc_gamg_agg->aggressive_coarsening_levels = 1; 1569d4adc10fSMark Adams pc_gamg_agg->use_aggressive_square_graph = PETSC_TRUE; 1570d529f056Smarkadams4 pc_gamg_agg->use_minimum_degree_ordering = PETSC_FALSE; 1571a9ccf005SMark Adams pc_gamg_agg->use_low_mem_filter = PETSC_FALSE; 1572d529f056Smarkadams4 pc_gamg_agg->aggressive_mis_k = 2; 1573*6c34c54dSStefano Zampini pc_gamg_agg->graph_symmetrize = PETSC_TRUE; 1574cab9ed1eSBarry Smith 15759566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", PCGAMGSetNSmooths_AGG)); 1576bae903cbSmarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", PCGAMGSetAggressiveLevels_AGG)); 1577d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveSquareGraph_C", PCGAMGSetAggressiveSquareGraph_AGG)); 1578d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetMinDegreeOrdering_C", PCGAMGMISkSetMinDegreeOrdering_AGG)); 1579a9ccf005SMark Adams PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetLowMemoryFilter_C", PCGAMGSetLowMemoryFilter_AGG)); 1580d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetAggressive_C", PCGAMGMISkSetAggressive_AGG)); 1581*6c34c54dSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetGraphSymmetrize_C", PCGAMGSetGraphSymmetrize_AGG)); 15829566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_AGG)); 15833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15842e68589bSMark F. Adams } 1585