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 { 11*a077d33dSBarry 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; 172d776b49SBarry Smith MatCoarsen crs; 182e68589bSMark F. Adams } PC_GAMG_AGG; 192e68589bSMark F. Adams 202e68589bSMark F. Adams /*@ 21*a077d33dSBarry Smith PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical) used to construct the prolongation operator 222e68589bSMark F. Adams 23c3339decSBarry Smith Logically Collective 242e68589bSMark F. Adams 252e68589bSMark F. Adams Input Parameters: 2620f4b53cSBarry Smith + pc - the preconditioner context 2720f4b53cSBarry Smith - n - the number of smooths 282e68589bSMark F. Adams 292e68589bSMark F. Adams Options Database Key: 30*a077d33dSBarry Smith . -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use 312e68589bSMark F. Adams 322e68589bSMark F. Adams Level: intermediate 332e68589bSMark F. Adams 34*a077d33dSBarry Smith Note: 35*a077d33dSBarry Smith This is a different concept from the number smoothing steps used during the linear solution process which 36*a077d33dSBarry Smith can be set with `-mg_levels_ksp_max_it` 37*a077d33dSBarry Smith 38*a077d33dSBarry Smith Developer Note: 39*a077d33dSBarry Smith This should be named `PCGAMGAGGSetNSmooths()`. 40*a077d33dSBarry Smith 41*a077d33dSBarry Smith .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCMG`, `PCGAMG` 422e68589bSMark F. Adams @*/ 43d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n) 44d71ae5a4SJacob Faibussowitsch { 452e68589bSMark F. Adams PetscFunctionBegin; 462e68589bSMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 47d5d25401SBarry Smith PetscValidLogicalCollectiveInt(pc, n, 2); 48cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetNSmooths_C", (PC, PetscInt), (pc, n)); 493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 502e68589bSMark F. Adams } 512e68589bSMark F. Adams 52d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetNSmooths_AGG(PC pc, PetscInt n) 53d71ae5a4SJacob Faibussowitsch { 542e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 552e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 56c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 572e68589bSMark F. Adams 582e68589bSMark F. Adams PetscFunctionBegin; 59c8b0795cSMark F. Adams pc_gamg_agg->nsmooths = n; 603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 61c8b0795cSMark F. Adams } 62c8b0795cSMark F. Adams 63c8b0795cSMark F. Adams /*@ 64f1580f4eSBarry Smith PCGAMGSetAggressiveLevels - Use aggressive coarsening on first n levels 65ef4ad70eSMark F. Adams 66c3339decSBarry Smith Logically Collective 67ef4ad70eSMark F. Adams 68ef4ad70eSMark F. Adams Input Parameters: 69cab9ed1eSBarry Smith + pc - the preconditioner context 70d5d25401SBarry Smith - n - 0, 1 or more 71ef4ad70eSMark F. Adams 72ef4ad70eSMark F. Adams Options Database Key: 73*a077d33dSBarry Smith . -pc_gamg_aggressive_coarsening <n,default = 1> - Number of levels on which to square the graph on before aggregating it 74af3c827dSMark Adams 75ef4ad70eSMark F. Adams Level: intermediate 76ef4ad70eSMark F. Adams 77*a077d33dSBarry 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()` 78ef4ad70eSMark F. Adams @*/ 79d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetAggressiveLevels(PC pc, PetscInt n) 80d71ae5a4SJacob Faibussowitsch { 81ef4ad70eSMark F. Adams PetscFunctionBegin; 82ef4ad70eSMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 83d5d25401SBarry Smith PetscValidLogicalCollectiveInt(pc, n, 2); 84bae903cbSmarkadams4 PetscTryMethod(pc, "PCGAMGSetAggressiveLevels_C", (PC, PetscInt), (pc, n)); 853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 86ef4ad70eSMark F. Adams } 87ef4ad70eSMark F. Adams 88d529f056Smarkadams4 /*@ 89d529f056Smarkadams4 PCGAMGMISkSetAggressive - Number (k) distance in MIS coarsening (>2 is 'aggressive') 90d529f056Smarkadams4 91d529f056Smarkadams4 Logically Collective 92d529f056Smarkadams4 93d529f056Smarkadams4 Input Parameters: 94d529f056Smarkadams4 + pc - the preconditioner context 95d529f056Smarkadams4 - n - 1 or more (default = 2) 96d529f056Smarkadams4 97d529f056Smarkadams4 Options Database Key: 98d529f056Smarkadams4 . -pc_gamg_aggressive_mis_k <n,default=2> - Number (k) distance in MIS coarsening (>2 is 'aggressive') 99d529f056Smarkadams4 100d529f056Smarkadams4 Level: intermediate 101d529f056Smarkadams4 102*a077d33dSBarry 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()` 103d529f056Smarkadams4 @*/ 104d529f056Smarkadams4 PetscErrorCode PCGAMGMISkSetAggressive(PC pc, PetscInt n) 105d529f056Smarkadams4 { 106d529f056Smarkadams4 PetscFunctionBegin; 107d529f056Smarkadams4 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 108d529f056Smarkadams4 PetscValidLogicalCollectiveInt(pc, n, 2); 109d529f056Smarkadams4 PetscTryMethod(pc, "PCGAMGMISkSetAggressive_C", (PC, PetscInt), (pc, n)); 110d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 111d529f056Smarkadams4 } 112d529f056Smarkadams4 113d529f056Smarkadams4 /*@ 114d529f056Smarkadams4 PCGAMGSetAggressiveSquareGraph - Use graph square A'A for aggressive coarsening, old method 115d529f056Smarkadams4 116d529f056Smarkadams4 Logically Collective 117d529f056Smarkadams4 118d529f056Smarkadams4 Input Parameters: 119d529f056Smarkadams4 + pc - the preconditioner context 120d529f056Smarkadams4 - b - default false - MIS-k is faster 121d529f056Smarkadams4 122d529f056Smarkadams4 Options Database Key: 123d529f056Smarkadams4 . -pc_gamg_aggressive_square_graph <bool,default=false> - Use square graph (A'A) or MIS-k (k=2) for aggressive coarsening 124d529f056Smarkadams4 125d529f056Smarkadams4 Level: intermediate 126d529f056Smarkadams4 127*a077d33dSBarry 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()` 128d529f056Smarkadams4 @*/ 129d529f056Smarkadams4 PetscErrorCode PCGAMGSetAggressiveSquareGraph(PC pc, PetscBool b) 130d529f056Smarkadams4 { 131d529f056Smarkadams4 PetscFunctionBegin; 132d529f056Smarkadams4 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 133d529f056Smarkadams4 PetscValidLogicalCollectiveBool(pc, b, 2); 134d529f056Smarkadams4 PetscTryMethod(pc, "PCGAMGSetAggressiveSquareGraph_C", (PC, PetscBool), (pc, b)); 135d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 136d529f056Smarkadams4 } 137d529f056Smarkadams4 138d529f056Smarkadams4 /*@ 139d529f056Smarkadams4 PCGAMGMISkSetMinDegreeOrdering - Use minimum degree ordering in greedy MIS algorithm 140d529f056Smarkadams4 141d529f056Smarkadams4 Logically Collective 142d529f056Smarkadams4 143d529f056Smarkadams4 Input Parameters: 144d529f056Smarkadams4 + pc - the preconditioner context 145d529f056Smarkadams4 - b - default true 146d529f056Smarkadams4 147d529f056Smarkadams4 Options Database Key: 148d529f056Smarkadams4 . -pc_gamg_mis_k_minimum_degree_ordering <bool,default=true> - Use minimum degree ordering in greedy MIS algorithm 149d529f056Smarkadams4 150d529f056Smarkadams4 Level: intermediate 151d529f056Smarkadams4 152*a077d33dSBarry 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()` 153d529f056Smarkadams4 @*/ 154d529f056Smarkadams4 PetscErrorCode PCGAMGMISkSetMinDegreeOrdering(PC pc, PetscBool b) 155d529f056Smarkadams4 { 156d529f056Smarkadams4 PetscFunctionBegin; 157d529f056Smarkadams4 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 158d529f056Smarkadams4 PetscValidLogicalCollectiveBool(pc, b, 2); 159d529f056Smarkadams4 PetscTryMethod(pc, "PCGAMGMISkSetMinDegreeOrdering_C", (PC, PetscBool), (pc, b)); 160d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 161d529f056Smarkadams4 } 162d529f056Smarkadams4 163a9ccf005SMark Adams /*@ 164a9ccf005SMark Adams PCGAMGSetLowMemoryFilter - Use low memory graph/matrix filter 165a9ccf005SMark Adams 166a9ccf005SMark Adams Logically Collective 167a9ccf005SMark Adams 168a9ccf005SMark Adams Input Parameters: 169a9ccf005SMark Adams + pc - the preconditioner context 170a9ccf005SMark Adams - b - default false 171a9ccf005SMark Adams 172a9ccf005SMark Adams Options Database Key: 173a9ccf005SMark Adams . -pc_gamg_low_memory_threshold_filter <bool,default=false> - Use low memory graph/matrix filter 174a9ccf005SMark Adams 175a9ccf005SMark Adams Level: intermediate 176a9ccf005SMark Adams 177*a077d33dSBarry Smith .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`, 178*a077d33dSBarry Smith `PCGAMGMISkSetAggressive()`, `PCGAMGSetAggressiveSquareGraph()`, `PCGAMGMISkSetMinDegreeOrdering()` 179a9ccf005SMark Adams @*/ 180a9ccf005SMark Adams PetscErrorCode PCGAMGSetLowMemoryFilter(PC pc, PetscBool b) 181a9ccf005SMark Adams { 182a9ccf005SMark Adams PetscFunctionBegin; 183a9ccf005SMark Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 184a9ccf005SMark Adams PetscValidLogicalCollectiveBool(pc, b, 2); 185a9ccf005SMark Adams PetscTryMethod(pc, "PCGAMGSetLowMemoryFilter_C", (PC, PetscBool), (pc, b)); 186a9ccf005SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 187a9ccf005SMark Adams } 188a9ccf005SMark Adams 189d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetAggressiveLevels_AGG(PC pc, PetscInt n) 190d71ae5a4SJacob Faibussowitsch { 191ef4ad70eSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 192ef4ad70eSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 193ef4ad70eSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 194ef4ad70eSMark F. Adams 195ef4ad70eSMark F. Adams PetscFunctionBegin; 196bae903cbSmarkadams4 pc_gamg_agg->aggressive_coarsening_levels = n; 1973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 198ef4ad70eSMark F. Adams } 199ef4ad70eSMark F. Adams 200d529f056Smarkadams4 static PetscErrorCode PCGAMGMISkSetAggressive_AGG(PC pc, PetscInt n) 201d529f056Smarkadams4 { 202d529f056Smarkadams4 PC_MG *mg = (PC_MG *)pc->data; 203d529f056Smarkadams4 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 204d529f056Smarkadams4 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 205d529f056Smarkadams4 206d529f056Smarkadams4 PetscFunctionBegin; 207d529f056Smarkadams4 pc_gamg_agg->aggressive_mis_k = n; 208d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 209d529f056Smarkadams4 } 210d529f056Smarkadams4 211d529f056Smarkadams4 static PetscErrorCode PCGAMGSetAggressiveSquareGraph_AGG(PC pc, PetscBool b) 212d529f056Smarkadams4 { 213d529f056Smarkadams4 PC_MG *mg = (PC_MG *)pc->data; 214d529f056Smarkadams4 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 215d529f056Smarkadams4 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 216d529f056Smarkadams4 217d529f056Smarkadams4 PetscFunctionBegin; 218d529f056Smarkadams4 pc_gamg_agg->use_aggressive_square_graph = b; 219d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 220d529f056Smarkadams4 } 221d529f056Smarkadams4 222a9ccf005SMark Adams static PetscErrorCode PCGAMGSetLowMemoryFilter_AGG(PC pc, PetscBool b) 223a9ccf005SMark Adams { 224a9ccf005SMark Adams PC_MG *mg = (PC_MG *)pc->data; 225a9ccf005SMark Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 226a9ccf005SMark Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 227a9ccf005SMark Adams 228a9ccf005SMark Adams PetscFunctionBegin; 229a9ccf005SMark Adams pc_gamg_agg->use_low_mem_filter = b; 230a9ccf005SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 231a9ccf005SMark Adams } 232a9ccf005SMark Adams 233d529f056Smarkadams4 static PetscErrorCode PCGAMGMISkSetMinDegreeOrdering_AGG(PC pc, PetscBool b) 234d529f056Smarkadams4 { 235d529f056Smarkadams4 PC_MG *mg = (PC_MG *)pc->data; 236d529f056Smarkadams4 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 237d529f056Smarkadams4 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 238d529f056Smarkadams4 239d529f056Smarkadams4 PetscFunctionBegin; 240d529f056Smarkadams4 pc_gamg_agg->use_minimum_degree_ordering = b; 241d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 242d529f056Smarkadams4 } 243d529f056Smarkadams4 244d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_GAMG_AGG(PC pc, PetscOptionItems *PetscOptionsObject) 245d71ae5a4SJacob Faibussowitsch { 2462e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 2472e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 248c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 249d4adc10fSMark 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; 250d4adc10fSMark Adams PetscInt nsq_graph_old = 0; 2512e68589bSMark F. Adams 2522e68589bSMark F. Adams PetscFunctionBegin; 253d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "GAMG-AGG options"); 254*a077d33dSBarry 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)); 255d4adc10fSMark Adams // aggressive coarsening logic with deprecated -pc_gamg_square_graph 256d4adc10fSMark 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)); 257d4adc10fSMark Adams if (!n_aggressive_flg) 258d4adc10fSMark 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)); 259d4adc10fSMark 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)); 260d4adc10fSMark Adams if (!new_sq_provided && old_sq_provided) { 261d4adc10fSMark Adams pc_gamg_agg->aggressive_coarsening_levels = nsq_graph_old; // could be zero 262d4adc10fSMark Adams pc_gamg_agg->use_aggressive_square_graph = PETSC_TRUE; 263bae903cbSmarkadams4 } 264d4adc10fSMark Adams if (new_sq_provided && old_sq_provided) 265d4adc10fSMark 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)); 266d529f056Smarkadams4 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)); 267a9ccf005SMark 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)); 268d529f056Smarkadams4 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)); 269d0609cedSBarry Smith PetscOptionsHeadEnd(); 2703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2712e68589bSMark F. Adams } 2722e68589bSMark F. Adams 273d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_GAMG_AGG(PC pc) 274d71ae5a4SJacob Faibussowitsch { 2752e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 2762e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 2772e68589bSMark F. Adams 2782e68589bSMark F. Adams PetscFunctionBegin; 2799566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->subctx)); 2802e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", NULL)); 281bae903cbSmarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", NULL)); 282d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetAggressive_C", NULL)); 283d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetMinDegreeOrdering_C", NULL)); 284a9ccf005SMark Adams PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetLowMemoryFilter_C", NULL)); 285d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveSquareGraph_C", NULL)); 286bae903cbSmarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL)); 2873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2882e68589bSMark F. Adams } 2892e68589bSMark F. Adams 2902e68589bSMark F. Adams /* 2912e68589bSMark F. Adams PCSetCoordinates_AGG 29220f4b53cSBarry Smith 29320f4b53cSBarry Smith Collective 2942e68589bSMark F. Adams 2952e68589bSMark F. Adams Input Parameter: 2962e68589bSMark F. Adams . pc - the preconditioner context 297145b44c9SPierre Jolivet . ndm - dimension of data (used for dof/vertex for Stokes) 298302f38e8SMark F. Adams . a_nloc - number of vertices local 299302f38e8SMark F. Adams . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...} 3002e68589bSMark F. Adams */ 3011e6b0712SBarry Smith 302d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords) 303d71ae5a4SJacob Faibussowitsch { 3042e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 3052e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 30669344418SMark F. Adams PetscInt arrsz, kk, ii, jj, nloc, ndatarows, ndf; 307a2f3521dSMark F. Adams Mat mat = pc->pmat; 3082e68589bSMark F. Adams 3092e68589bSMark F. Adams PetscFunctionBegin; 310a2f3521dSMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 311a2f3521dSMark F. Adams PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 312302f38e8SMark F. Adams nloc = a_nloc; 3132e68589bSMark F. Adams 3142e68589bSMark F. Adams /* SA: null space vectors */ 3159566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(mat, &ndf)); /* this does not work for Stokes */ 31669344418SMark F. Adams if (coords && ndf == 1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */ 317a2f3521dSMark F. Adams else if (coords) { 31863a3b9bcSJacob Faibussowitsch PetscCheck(ndm <= ndf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "degrees of motion %" PetscInt_FMT " > block size %" PetscInt_FMT, ndm, ndf); 31969344418SMark F. Adams pc_gamg->data_cell_cols = (ndm == 2 ? 3 : 6); /* displacement elasticity */ 320ad540459SPierre 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); 321806fa848SBarry Smith } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */ 32271959b99SBarry Smith pc_gamg->data_cell_rows = ndatarows = ndf; 32363a3b9bcSJacob 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); 324c8b0795cSMark F. Adams arrsz = nloc * pc_gamg->data_cell_rows * pc_gamg->data_cell_cols; 3252e68589bSMark F. Adams 3267f66b68fSMark Adams if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) { 3279566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 3289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(arrsz + 1, &pc_gamg->data)); 3292e68589bSMark F. Adams } 3305e116b59SBarry Smith /* copy data in - column-oriented */ 3312e68589bSMark F. Adams for (kk = 0; kk < nloc; kk++) { 33269344418SMark F. Adams const PetscInt M = nloc * pc_gamg->data_cell_rows; /* stride into data */ 33369344418SMark F. Adams PetscReal *data = &pc_gamg->data[kk * ndatarows]; /* start of cell */ 334c8b0795cSMark F. Adams if (pc_gamg->data_cell_cols == 1) *data = 1.0; 3352e68589bSMark F. Adams else { 33669344418SMark F. Adams /* translational modes */ 3372fa5cd67SKarl Rupp for (ii = 0; ii < ndatarows; ii++) { 3382fa5cd67SKarl Rupp for (jj = 0; jj < ndatarows; jj++) { 33969344418SMark F. Adams if (ii == jj) data[ii * M + jj] = 1.0; 3402e68589bSMark F. Adams else data[ii * M + jj] = 0.0; 3412fa5cd67SKarl Rupp } 3422fa5cd67SKarl Rupp } 34369344418SMark F. Adams 34469344418SMark F. Adams /* rotational modes */ 3452e68589bSMark F. Adams if (coords) { 34669344418SMark F. Adams if (ndm == 2) { 3472e68589bSMark F. Adams data += 2 * M; 3482e68589bSMark F. Adams data[0] = -coords[2 * kk + 1]; 3492e68589bSMark F. Adams data[1] = coords[2 * kk]; 350806fa848SBarry Smith } else { 3512e68589bSMark F. Adams data += 3 * M; 3529371c9d4SSatish Balay data[0] = 0.0; 3539371c9d4SSatish Balay data[M + 0] = coords[3 * kk + 2]; 3549371c9d4SSatish Balay data[2 * M + 0] = -coords[3 * kk + 1]; 3559371c9d4SSatish Balay data[1] = -coords[3 * kk + 2]; 3569371c9d4SSatish Balay data[M + 1] = 0.0; 3579371c9d4SSatish Balay data[2 * M + 1] = coords[3 * kk]; 3589371c9d4SSatish Balay data[2] = coords[3 * kk + 1]; 3599371c9d4SSatish Balay data[M + 2] = -coords[3 * kk]; 3609371c9d4SSatish Balay data[2 * M + 2] = 0.0; 3612e68589bSMark F. Adams } 3622e68589bSMark F. Adams } 3632e68589bSMark F. Adams } 3642e68589bSMark F. Adams } 3652e68589bSMark F. Adams pc_gamg->data_sz = arrsz; 3663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3672e68589bSMark F. Adams } 3682e68589bSMark F. Adams 3692e68589bSMark F. Adams /* 370a2f3521dSMark F. Adams PCSetData_AGG - called if data is not set with PCSetCoordinates. 371a2f3521dSMark F. Adams Looks in Mat for near null space. 372a2f3521dSMark F. Adams Does not work for Stokes 3732e68589bSMark F. Adams 3742e68589bSMark F. Adams Input Parameter: 3752e68589bSMark F. Adams . pc - 376a2f3521dSMark F. Adams . a_A - matrix to get (near) null space out of. 3772e68589bSMark F. Adams */ 378d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A) 379d71ae5a4SJacob Faibussowitsch { 380b8cd405aSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 381b8cd405aSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 382b8cd405aSMark F. Adams MatNullSpace mnull; 38366f2ef4bSMark Adams 384b817416eSBarry Smith PetscFunctionBegin; 3859566063dSJacob Faibussowitsch PetscCall(MatGetNearNullSpace(a_A, &mnull)); 386b8cd405aSMark F. Adams if (!mnull) { 38766f2ef4bSMark Adams DM dm; 3889566063dSJacob Faibussowitsch PetscCall(PCGetDM(pc, &dm)); 38948a46eb9SPierre Jolivet if (!dm) PetscCall(MatGetDM(a_A, &dm)); 39066f2ef4bSMark Adams if (dm) { 39166f2ef4bSMark Adams PetscObject deformation; 392b0d30dd6SMatthew G. Knepley PetscInt Nf; 393b0d30dd6SMatthew G. Knepley 3949566063dSJacob Faibussowitsch PetscCall(DMGetNumFields(dm, &Nf)); 395b0d30dd6SMatthew G. Knepley if (Nf) { 3969566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, 0, NULL, &deformation)); 3979566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)deformation, "nearnullspace", (PetscObject *)&mnull)); 39848a46eb9SPierre Jolivet if (!mnull) PetscCall(PetscObjectQuery((PetscObject)deformation, "nullspace", (PetscObject *)&mnull)); 39966f2ef4bSMark Adams } 40066f2ef4bSMark Adams } 401b0d30dd6SMatthew G. Knepley } 40266f2ef4bSMark Adams 40366f2ef4bSMark Adams if (!mnull) { 404a2f3521dSMark F. Adams PetscInt bs, NN, MM; 4059566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(a_A, &bs)); 4069566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(a_A, &MM, &NN)); 40763a3b9bcSJacob Faibussowitsch PetscCheck(MM % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MM %" PetscInt_FMT " must be divisible by bs %" PetscInt_FMT, MM, bs); 4089566063dSJacob Faibussowitsch PetscCall(PCSetCoordinates_AGG(pc, bs, MM / bs, NULL)); 409806fa848SBarry Smith } else { 410b8cd405aSMark F. Adams PetscReal *nullvec; 411b8cd405aSMark F. Adams PetscBool has_const; 412b8cd405aSMark F. Adams PetscInt i, j, mlocal, nvec, bs; 413d5d25401SBarry Smith const Vec *vecs; 414d5d25401SBarry Smith const PetscScalar *v; 415b817416eSBarry Smith 4169566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(a_A, &mlocal, NULL)); 4179566063dSJacob Faibussowitsch PetscCall(MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs)); 418d19c4ebbSmarkadams4 for (i = 0; i < nvec; i++) { 419d19c4ebbSmarkadams4 PetscCall(VecGetLocalSize(vecs[i], &j)); 420d19c4ebbSmarkadams4 PetscCheck(j == mlocal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Attached null space vector size %" PetscInt_FMT " != matrix size %" PetscInt_FMT, j, mlocal); 421d19c4ebbSmarkadams4 } 422a0dea87dSMark Adams pc_gamg->data_sz = (nvec + !!has_const) * mlocal; 4239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((nvec + !!has_const) * mlocal, &nullvec)); 4249371c9d4SSatish Balay if (has_const) 4259371c9d4SSatish Balay for (i = 0; i < mlocal; i++) nullvec[i] = 1.0; 426b8cd405aSMark F. Adams for (i = 0; i < nvec; i++) { 4279566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(vecs[i], &v)); 428b8cd405aSMark F. Adams for (j = 0; j < mlocal; j++) nullvec[(i + !!has_const) * mlocal + j] = PetscRealPart(v[j]); 4299566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(vecs[i], &v)); 430b8cd405aSMark F. Adams } 431b8cd405aSMark F. Adams pc_gamg->data = nullvec; 432b8cd405aSMark F. Adams pc_gamg->data_cell_cols = (nvec + !!has_const); 4339566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(a_A, &bs)); 434b8cd405aSMark F. Adams pc_gamg->data_cell_rows = bs; 435b8cd405aSMark F. Adams } 4363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4372e68589bSMark F. Adams } 4382e68589bSMark F. Adams 4392e68589bSMark F. Adams /* 440bae903cbSmarkadams4 formProl0 - collect null space data for each aggregate, do QR, put R in coarse grid data and Q in P_0 4412e68589bSMark F. Adams 4422e68589bSMark F. Adams Input Parameter: 4432adfe9d3SBarry Smith . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices 4442adfe9d3SBarry Smith . bs - row block size 4450cbbd2e1SMark F. Adams . nSAvec - column bs of new P 4460cbbd2e1SMark F. Adams . my0crs - global index of start of locals 4472adfe9d3SBarry Smith . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec] 4482e68589bSMark F. Adams . data_in[data_stride*nSAvec] - local data on fine grid 4492e68589bSMark F. Adams . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist' 450bae903cbSmarkadams4 4512e68589bSMark F. Adams Output Parameter: 4522e68589bSMark F. Adams . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data 4532e68589bSMark F. Adams . a_Prol - prolongation operator 4542e68589bSMark F. Adams */ 455d71ae5a4SJacob 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) 456d71ae5a4SJacob Faibussowitsch { 457bd026e97SJed Brown PetscInt Istart, my0, Iend, nloc, clid, flid = 0, aggID, kk, jj, ii, mm, nSelected, minsz, nghosts, out_data_stride; 4583b4367a7SBarry Smith MPI_Comm comm; 4592e68589bSMark F. Adams PetscReal *out_data; 460539c167fSBarry Smith PetscCDIntNd *pos; 4611943db53SBarry Smith PCGAMGHashTable fgid_flid; 4620cbbd2e1SMark F. Adams 4632e68589bSMark F. Adams PetscFunctionBegin; 4649566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)a_Prol, &comm)); 4659566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(a_Prol, &Istart, &Iend)); 4669371c9d4SSatish Balay nloc = (Iend - Istart) / bs; 4679371c9d4SSatish Balay my0 = Istart / bs; 46863a3b9bcSJacob 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); 4690cbbd2e1SMark F. Adams Iend /= bs; 4700cbbd2e1SMark F. Adams nghosts = data_stride / bs - nloc; 4712e68589bSMark F. Adams 4729566063dSJacob Faibussowitsch PetscCall(PCGAMGHashTableCreate(2 * nghosts + 1, &fgid_flid)); 47348a46eb9SPierre Jolivet for (kk = 0; kk < nghosts; kk++) PetscCall(PCGAMGHashTableAdd(&fgid_flid, flid_fgid[nloc + kk], nloc + kk)); 4742e68589bSMark F. Adams 4750cbbd2e1SMark F. Adams /* count selected -- same as number of cols of P */ 4760cbbd2e1SMark F. Adams for (nSelected = mm = 0; mm < nloc; mm++) { 477e78576d6SMark F. Adams PetscBool ise; 4785e62d202SMark Adams PetscCall(PetscCDIsEmptyAt(agg_llists, mm, &ise)); 479e78576d6SMark F. Adams if (!ise) nSelected++; 4800cbbd2e1SMark F. Adams } 4819566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(a_Prol, &ii, &jj)); 48263a3b9bcSJacob Faibussowitsch PetscCheck((ii / nSAvec) == my0crs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ii %" PetscInt_FMT " /nSAvec %" PetscInt_FMT " != my0crs %" PetscInt_FMT, ii, nSAvec, my0crs); 48363a3b9bcSJacob 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); 4840cbbd2e1SMark F. Adams 4852e68589bSMark F. Adams /* aloc space for coarse point data (output) */ 4860cbbd2e1SMark F. Adams out_data_stride = nSelected * nSAvec; 4872fa5cd67SKarl Rupp 4889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(out_data_stride * nSAvec, &out_data)); 48933812677SSatish Balay for (ii = 0; ii < out_data_stride * nSAvec; ii++) out_data[ii] = PETSC_MAX_REAL; 4900cbbd2e1SMark F. Adams *a_data_out = out_data; /* output - stride nSelected*nSAvec */ 4912e68589bSMark F. Adams 4922e68589bSMark F. Adams /* find points and set prolongation */ 493c8b0795cSMark F. Adams minsz = 100; 4940cbbd2e1SMark F. Adams for (mm = clid = 0; mm < nloc; mm++) { 4955e62d202SMark Adams PetscCall(PetscCDCountAt(agg_llists, mm, &jj)); 496e78576d6SMark F. Adams if (jj > 0) { 4970cbbd2e1SMark F. Adams const PetscInt lid = mm, cgid = my0crs + clid; 4980cbbd2e1SMark F. Adams PetscInt cids[100]; /* max bs */ 4990cbbd2e1SMark F. Adams PetscBLASInt asz = jj, M = asz * bs, N = nSAvec, INFO; 5002e68589bSMark F. Adams PetscBLASInt Mdata = M + ((N - M > 0) ? N - M : 0), LDA = Mdata, LWORK = N * bs; 5012e68589bSMark F. Adams PetscScalar *qqc, *qqr, *TAU, *WORK; 5022e68589bSMark F. Adams PetscInt *fids; 50365d7b583SSatish Balay PetscReal *data; 504b817416eSBarry Smith 5050cbbd2e1SMark F. Adams /* count agg */ 5060cbbd2e1SMark F. Adams if (asz < minsz) minsz = asz; 5070cbbd2e1SMark F. Adams 5080cbbd2e1SMark F. Adams /* get block */ 5099566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(Mdata * N, &qqc, M * N, &qqr, N, &TAU, LWORK, &WORK, M, &fids)); 5102e68589bSMark F. Adams 5112e68589bSMark F. Adams aggID = 0; 5129566063dSJacob Faibussowitsch PetscCall(PetscCDGetHeadPos(agg_llists, lid, &pos)); 513e78576d6SMark F. Adams while (pos) { 514e78576d6SMark F. Adams PetscInt gid1; 5159566063dSJacob Faibussowitsch PetscCall(PetscCDIntNdGetID(pos, &gid1)); 5169566063dSJacob Faibussowitsch PetscCall(PetscCDGetNextPos(agg_llists, lid, &pos)); 517e78576d6SMark F. Adams 5180cbbd2e1SMark F. Adams if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0; 5190cbbd2e1SMark F. Adams else { 5209566063dSJacob Faibussowitsch PetscCall(PCGAMGHashTableFind(&fgid_flid, gid1, &flid)); 52108401ef6SPierre Jolivet PetscCheck(flid >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot find gid1 in table"); 5220cbbd2e1SMark F. Adams } 5235e116b59SBarry Smith /* copy in B_i matrix - column-oriented */ 52465d7b583SSatish Balay data = &data_in[flid * bs]; 5253d3eaba7SBarry Smith for (ii = 0; ii < bs; ii++) { 5262e68589bSMark F. Adams for (jj = 0; jj < N; jj++) { 5270cbbd2e1SMark F. Adams PetscReal d = data[jj * data_stride + ii]; 5280cbbd2e1SMark F. Adams qqc[jj * Mdata + aggID * bs + ii] = d; 5292e68589bSMark F. Adams } 5302e68589bSMark F. Adams } 5312e68589bSMark F. Adams /* set fine IDs */ 5322e68589bSMark F. Adams for (kk = 0; kk < bs; kk++) fids[aggID * bs + kk] = flid_fgid[flid] * bs + kk; 5332e68589bSMark F. Adams aggID++; 5340cbbd2e1SMark F. Adams } 5352e68589bSMark F. Adams 5362e68589bSMark F. Adams /* pad with zeros */ 5372e68589bSMark F. Adams for (ii = asz * bs; ii < Mdata; ii++) { 538ad540459SPierre Jolivet for (jj = 0; jj < N; jj++, kk++) qqc[jj * Mdata + ii] = .0; 5392e68589bSMark F. Adams } 5402e68589bSMark F. Adams 5412e68589bSMark F. Adams /* QR */ 5429566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 543792fecdfSBarry Smith PetscCallBLAS("LAPACKgeqrf", LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 5449566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPop()); 54508401ef6SPierre Jolivet PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xGEQRF error"); 5465e116b59SBarry Smith /* get R - column-oriented - output B_{i+1} */ 5472e68589bSMark F. Adams { 5482e68589bSMark F. Adams PetscReal *data = &out_data[clid * nSAvec]; 5492e68589bSMark F. Adams for (jj = 0; jj < nSAvec; jj++) { 5502e68589bSMark F. Adams for (ii = 0; ii < nSAvec; ii++) { 55108401ef6SPierre 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); 5520cbbd2e1SMark F. Adams if (ii <= jj) data[jj * out_data_stride + ii] = PetscRealPart(qqc[jj * Mdata + ii]); 5530cbbd2e1SMark F. Adams else data[jj * out_data_stride + ii] = 0.; 5542e68589bSMark F. Adams } 5552e68589bSMark F. Adams } 5562e68589bSMark F. Adams } 5572e68589bSMark F. Adams 5585e116b59SBarry Smith /* get Q - row-oriented */ 559792fecdfSBarry Smith PetscCallBLAS("LAPACKorgqr", LAPACKorgqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 56063a3b9bcSJacob Faibussowitsch PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xORGQR error arg %" PetscBLASInt_FMT, -INFO); 5612e68589bSMark F. Adams 5622e68589bSMark F. Adams for (ii = 0; ii < M; ii++) { 563ad540459SPierre Jolivet for (jj = 0; jj < N; jj++) qqr[N * ii + jj] = qqc[jj * Mdata + ii]; 5642e68589bSMark F. Adams } 5652e68589bSMark F. Adams 5662e68589bSMark F. Adams /* add diagonal block of P0 */ 5679371c9d4SSatish Balay for (kk = 0; kk < N; kk++) { cids[kk] = N * cgid + kk; /* global col IDs in P0 */ } 5689566063dSJacob Faibussowitsch PetscCall(MatSetValues(a_Prol, M, fids, N, cids, qqr, INSERT_VALUES)); 5699566063dSJacob Faibussowitsch PetscCall(PetscFree5(qqc, qqr, TAU, WORK, fids)); 570b43b03e9SMark F. Adams clid++; 5710cbbd2e1SMark F. Adams } /* coarse agg */ 5720cbbd2e1SMark F. Adams } /* for all fine nodes */ 5739566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(a_Prol, MAT_FINAL_ASSEMBLY)); 5749566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(a_Prol, MAT_FINAL_ASSEMBLY)); 5759566063dSJacob Faibussowitsch PetscCall(PCGAMGHashTableDestroy(&fgid_flid)); 5763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5772e68589bSMark F. Adams } 5782e68589bSMark F. Adams 579d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_GAMG_AGG(PC pc, PetscViewer viewer) 580d71ae5a4SJacob Faibussowitsch { 5815adeb434SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 5825adeb434SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 5835adeb434SBarry Smith PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 5845adeb434SBarry Smith 5855adeb434SBarry Smith PetscFunctionBegin; 5869566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " AGG specific options\n")); 587d529f056Smarkadams4 PetscCall(PetscViewerASCIIPrintf(viewer, " Number of levels of aggressive coarsening %d\n", (int)pc_gamg_agg->aggressive_coarsening_levels)); 588d529f056Smarkadams4 if (pc_gamg_agg->aggressive_coarsening_levels > 0) { 589d529f056Smarkadams4 PetscCall(PetscViewerASCIIPrintf(viewer, " %s aggressive coarsening\n", !pc_gamg_agg->use_aggressive_square_graph ? "MIS-k" : "Square graph")); 590d529f056Smarkadams4 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)); 591d529f056Smarkadams4 } 592bae903cbSmarkadams4 PetscCall(PetscViewerASCIIPrintf(viewer, " Number smoothing steps %d\n", (int)pc_gamg_agg->nsmooths)); 5933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5945adeb434SBarry Smith } 5955adeb434SBarry Smith 5962d776b49SBarry Smith static PetscErrorCode PCGAMGCreateGraph_AGG(PC pc, Mat Amat, Mat *a_Gmat) 597d71ae5a4SJacob Faibussowitsch { 598c8b0795cSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 599c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 600c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 6012d776b49SBarry Smith const PetscReal vfilter = pc_gamg->threshold[pc_gamg->current_level]; 6029c15c1aeSMark Adams PetscBool ishem, ismis; 6032d776b49SBarry Smith const char *prefix; 604d529f056Smarkadams4 MatInfo info0, info1; 605d529f056Smarkadams4 PetscInt bs; 606c8b0795cSMark F. Adams 607c8b0795cSMark F. Adams PetscFunctionBegin; 608a9ccf005SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0)); 6092d776b49SBarry 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 */ 6102d776b49SBarry Smith /* MATCOARSENHEM requires numerical weights for edges so ensure they are computed */ 6112d776b49SBarry Smith PetscCall(MatCoarsenCreate(PetscObjectComm((PetscObject)pc), &pc_gamg_agg->crs)); 6122d776b49SBarry Smith PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, &prefix)); 6132d776b49SBarry Smith PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc_gamg_agg->crs, prefix)); 6142d776b49SBarry Smith PetscCall(MatCoarsenSetFromOptions(pc_gamg_agg->crs)); 615e02fb3cdSMark Adams PetscCall(MatGetBlockSize(Amat, &bs)); 616e02fb3cdSMark Adams // check for valid indices wrt bs 617e02fb3cdSMark Adams for (int ii = 0; ii < pc_gamg_agg->crs->strength_index_size; ii++) { 618b65aec2dSMark 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", 619b65aec2dSMark Adams (int)pc_gamg_agg->crs->strength_index[ii], (int)bs); 620e02fb3cdSMark Adams } 6212d776b49SBarry Smith PetscCall(PetscObjectTypeCompare((PetscObject)pc_gamg_agg->crs, MATCOARSENHEM, &ishem)); 6225e62d202SMark Adams if (ishem) { 6239c15c1aeSMark 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)); 6245e62d202SMark Adams pc_gamg_agg->aggressive_coarsening_levels = 0; // aggressive and HEM does not make sense 6255e62d202SMark Adams PetscCall(MatCoarsenSetMaximumIterations(pc_gamg_agg->crs, pc_gamg_agg->crs->max_it)); // for code coverage 6265e62d202SMark Adams PetscCall(MatCoarsenSetThreshold(pc_gamg_agg->crs, vfilter)); // for code coverage 6279c15c1aeSMark Adams } else { 6289c15c1aeSMark Adams PetscCall(PetscObjectTypeCompare((PetscObject)pc_gamg_agg->crs, MATCOARSENMIS, &ismis)); 6299c15c1aeSMark Adams if (ismis && pc_gamg_agg->aggressive_coarsening_levels && !pc_gamg_agg->use_aggressive_square_graph) { 6309c15c1aeSMark Adams PetscCall(PetscInfo(pc, "MIS and aggressive coarsening and no square graph: force square graph\n")); 6319c15c1aeSMark Adams pc_gamg_agg->use_aggressive_square_graph = PETSC_TRUE; 6329c15c1aeSMark Adams } 6335e62d202SMark Adams } 634a9ccf005SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0)); 635a9ccf005SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0)); 636d529f056Smarkadams4 PetscCall(MatGetInfo(Amat, MAT_LOCAL, &info0)); /* global reduction */ 637a9ccf005SMark Adams 638a9ccf005SMark Adams if (ishem || pc_gamg_agg->use_low_mem_filter) { 639e02fb3cdSMark Adams PetscCall(MatCreateGraph(Amat, PETSC_TRUE, (vfilter >= 0 || ishem) ? PETSC_TRUE : PETSC_FALSE, vfilter, pc_gamg_agg->crs->strength_index_size, pc_gamg_agg->crs->strength_index, a_Gmat)); 640a9ccf005SMark Adams } else { 641d8b4a066SPierre Jolivet // make scalar graph, symetrize if not know to be symmetric, scale, but do not filter (expensive) 642e02fb3cdSMark Adams PetscCall(MatCreateGraph(Amat, PETSC_TRUE, PETSC_TRUE, -1, pc_gamg_agg->crs->strength_index_size, pc_gamg_agg->crs->strength_index, a_Gmat)); 643a9ccf005SMark Adams if (vfilter >= 0) { 644a9ccf005SMark Adams PetscInt Istart, Iend, ncols, nnz0, nnz1, NN, MM, nloc; 645a9ccf005SMark Adams Mat tGmat, Gmat = *a_Gmat; 646a9ccf005SMark Adams MPI_Comm comm; 647a9ccf005SMark Adams const PetscScalar *vals; 648a9ccf005SMark Adams const PetscInt *idx; 649a9ccf005SMark Adams PetscInt *d_nnz, *o_nnz, kk, *garray = NULL, *AJ, maxcols = 0; 650a9ccf005SMark Adams MatScalar *AA; // this is checked in graph 651a9ccf005SMark Adams PetscBool isseqaij; 652a9ccf005SMark Adams Mat a, b, c; 653a9ccf005SMark Adams MatType jtype; 654a9ccf005SMark Adams 655a9ccf005SMark Adams PetscCall(PetscObjectGetComm((PetscObject)Gmat, &comm)); 656a9ccf005SMark Adams PetscCall(PetscObjectBaseTypeCompare((PetscObject)Gmat, MATSEQAIJ, &isseqaij)); 657a9ccf005SMark Adams PetscCall(MatGetType(Gmat, &jtype)); 658a9ccf005SMark Adams PetscCall(MatCreate(comm, &tGmat)); 659a9ccf005SMark Adams PetscCall(MatSetType(tGmat, jtype)); 660a9ccf005SMark Adams 661a9ccf005SMark Adams /* TODO GPU: this can be called when filter = 0 -> Probably provide MatAIJThresholdCompress that compresses the entries below a threshold? 662a9ccf005SMark Adams Also, if the matrix is symmetric, can we skip this 663a9ccf005SMark Adams operation? It can be very expensive on large matrices. */ 664a9ccf005SMark Adams 665a9ccf005SMark Adams // global sizes 666a9ccf005SMark Adams PetscCall(MatGetSize(Gmat, &MM, &NN)); 667a9ccf005SMark Adams PetscCall(MatGetOwnershipRange(Gmat, &Istart, &Iend)); 668a9ccf005SMark Adams nloc = Iend - Istart; 669a9ccf005SMark Adams PetscCall(PetscMalloc2(nloc, &d_nnz, nloc, &o_nnz)); 670a9ccf005SMark Adams if (isseqaij) { 671a9ccf005SMark Adams a = Gmat; 672a9ccf005SMark Adams b = NULL; 673a9ccf005SMark Adams } else { 674a9ccf005SMark Adams Mat_MPIAIJ *d = (Mat_MPIAIJ *)Gmat->data; 675a9ccf005SMark Adams a = d->A; 676a9ccf005SMark Adams b = d->B; 677a9ccf005SMark Adams garray = d->garray; 678a9ccf005SMark Adams } 679a9ccf005SMark Adams /* Determine upper bound on non-zeros needed in new filtered matrix */ 680a9ccf005SMark Adams for (PetscInt row = 0; row < nloc; row++) { 681a9ccf005SMark Adams PetscCall(MatGetRow(a, row, &ncols, NULL, NULL)); 682a9ccf005SMark Adams d_nnz[row] = ncols; 683a9ccf005SMark Adams if (ncols > maxcols) maxcols = ncols; 684a9ccf005SMark Adams PetscCall(MatRestoreRow(a, row, &ncols, NULL, NULL)); 685a9ccf005SMark Adams } 686a9ccf005SMark Adams if (b) { 687a9ccf005SMark Adams for (PetscInt row = 0; row < nloc; row++) { 688a9ccf005SMark Adams PetscCall(MatGetRow(b, row, &ncols, NULL, NULL)); 689a9ccf005SMark Adams o_nnz[row] = ncols; 690a9ccf005SMark Adams if (ncols > maxcols) maxcols = ncols; 691a9ccf005SMark Adams PetscCall(MatRestoreRow(b, row, &ncols, NULL, NULL)); 692a9ccf005SMark Adams } 693a9ccf005SMark Adams } 694a9ccf005SMark Adams PetscCall(MatSetSizes(tGmat, nloc, nloc, MM, MM)); 695a9ccf005SMark Adams PetscCall(MatSetBlockSizes(tGmat, 1, 1)); 696a9ccf005SMark Adams PetscCall(MatSeqAIJSetPreallocation(tGmat, 0, d_nnz)); 697a9ccf005SMark Adams PetscCall(MatMPIAIJSetPreallocation(tGmat, 0, d_nnz, 0, o_nnz)); 698a9ccf005SMark Adams PetscCall(MatSetOption(tGmat, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); 699a9ccf005SMark Adams PetscCall(PetscFree2(d_nnz, o_nnz)); 700a9ccf005SMark Adams PetscCall(PetscMalloc2(maxcols, &AA, maxcols, &AJ)); 701a9ccf005SMark Adams nnz0 = nnz1 = 0; 702a9ccf005SMark Adams for (c = a, kk = 0; c && kk < 2; c = b, kk++) { 703a9ccf005SMark Adams for (PetscInt row = 0, grow = Istart, ncol_row, jj; row < nloc; row++, grow++) { 704a9ccf005SMark Adams PetscCall(MatGetRow(c, row, &ncols, &idx, &vals)); 705a9ccf005SMark Adams for (ncol_row = jj = 0; jj < ncols; jj++, nnz0++) { 706a9ccf005SMark Adams PetscScalar sv = PetscAbs(PetscRealPart(vals[jj])); 707a9ccf005SMark Adams if (PetscRealPart(sv) > vfilter) { 708a9ccf005SMark Adams PetscInt cid = idx[jj] + Istart; //diag 709a9ccf005SMark Adams nnz1++; 710a9ccf005SMark Adams if (c != a) cid = garray[idx[jj]]; 711a9ccf005SMark Adams AA[ncol_row] = vals[jj]; 712a9ccf005SMark Adams AJ[ncol_row] = cid; 713a9ccf005SMark Adams ncol_row++; 714a9ccf005SMark Adams } 715a9ccf005SMark Adams } 716a9ccf005SMark Adams PetscCall(MatRestoreRow(c, row, &ncols, &idx, &vals)); 717a9ccf005SMark Adams PetscCall(MatSetValues(tGmat, 1, &grow, ncol_row, AJ, AA, INSERT_VALUES)); 718a9ccf005SMark Adams } 719a9ccf005SMark Adams } 720a9ccf005SMark Adams PetscCall(PetscFree2(AA, AJ)); 721a9ccf005SMark Adams PetscCall(MatAssemblyBegin(tGmat, MAT_FINAL_ASSEMBLY)); 722a9ccf005SMark Adams PetscCall(MatAssemblyEnd(tGmat, MAT_FINAL_ASSEMBLY)); 723a9ccf005SMark Adams PetscCall(MatPropagateSymmetryOptions(Gmat, tGmat)); /* Normal Mat options are not relevant ? */ 724a9ccf005SMark 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)); 725a9ccf005SMark Adams PetscCall(MatViewFromOptions(tGmat, NULL, "-mat_filter_graph_view")); 726a9ccf005SMark Adams PetscCall(MatDestroy(&Gmat)); 727a9ccf005SMark Adams *a_Gmat = tGmat; 728a9ccf005SMark Adams } 729a9ccf005SMark Adams } 730a9ccf005SMark Adams 731d529f056Smarkadams4 PetscCall(MatGetInfo(*a_Gmat, MAT_LOCAL, &info1)); /* global reduction */ 732d529f056Smarkadams4 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)); 733849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0)); 7343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 735c8b0795cSMark F. Adams } 736c8b0795cSMark F. Adams 737d529f056Smarkadams4 typedef PetscInt NState; 738d529f056Smarkadams4 static const NState NOT_DONE = -2; 739d529f056Smarkadams4 static const NState DELETED = -1; 740d529f056Smarkadams4 static const NState REMOVED = -3; 741d529f056Smarkadams4 #define IS_SELECTED(s) (s != DELETED && s != NOT_DONE && s != REMOVED) 742d529f056Smarkadams4 743d529f056Smarkadams4 /* 744d529f056Smarkadams4 fixAggregatesWithSquare - greedy grab of with G1 (unsquared graph) -- AIJ specific -- change to fixAggregatesWithSquare -- TODD 745d529f056Smarkadams4 - AGG-MG specific: clears singletons out of 'selected_2' 746d529f056Smarkadams4 747d529f056Smarkadams4 Input Parameter: 748d529f056Smarkadams4 . Gmat_2 - global matrix of squared graph (data not defined) 749d529f056Smarkadams4 . Gmat_1 - base graph to grab with base graph 750d529f056Smarkadams4 Input/Output Parameter: 751d529f056Smarkadams4 . aggs_2 - linked list of aggs with gids) 752d529f056Smarkadams4 */ 753d529f056Smarkadams4 static PetscErrorCode fixAggregatesWithSquare(PC pc, Mat Gmat_2, Mat Gmat_1, PetscCoarsenData *aggs_2) 754d529f056Smarkadams4 { 755d529f056Smarkadams4 PetscBool isMPI; 756d529f056Smarkadams4 Mat_SeqAIJ *matA_1, *matB_1 = NULL; 757d529f056Smarkadams4 MPI_Comm comm; 758d529f056Smarkadams4 PetscInt lid, *ii, *idx, ix, Iend, my0, kk, n, j; 759d529f056Smarkadams4 Mat_MPIAIJ *mpimat_2 = NULL, *mpimat_1 = NULL; 760d529f056Smarkadams4 const PetscInt nloc = Gmat_2->rmap->n; 761d529f056Smarkadams4 PetscScalar *cpcol_1_state, *cpcol_2_state, *cpcol_2_par_orig, *lid_parent_gid; 762d529f056Smarkadams4 PetscInt *lid_cprowID_1 = NULL; 763d529f056Smarkadams4 NState *lid_state; 764d529f056Smarkadams4 Vec ghost_par_orig2; 765d529f056Smarkadams4 PetscMPIInt rank; 766d529f056Smarkadams4 767d529f056Smarkadams4 PetscFunctionBegin; 768d529f056Smarkadams4 PetscCall(PetscObjectGetComm((PetscObject)Gmat_2, &comm)); 769d529f056Smarkadams4 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 770d529f056Smarkadams4 PetscCall(MatGetOwnershipRange(Gmat_1, &my0, &Iend)); 771d529f056Smarkadams4 772d529f056Smarkadams4 /* get submatrices */ 773d529f056Smarkadams4 PetscCall(PetscStrbeginswith(((PetscObject)Gmat_1)->type_name, MATMPIAIJ, &isMPI)); 774d529f056Smarkadams4 PetscCall(PetscInfo(pc, "isMPI = %s\n", isMPI ? "yes" : "no")); 775d529f056Smarkadams4 PetscCall(PetscMalloc3(nloc, &lid_state, nloc, &lid_parent_gid, nloc, &lid_cprowID_1)); 776d529f056Smarkadams4 for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1; 777d529f056Smarkadams4 if (isMPI) { 778d529f056Smarkadams4 /* grab matrix objects */ 779d529f056Smarkadams4 mpimat_2 = (Mat_MPIAIJ *)Gmat_2->data; 780d529f056Smarkadams4 mpimat_1 = (Mat_MPIAIJ *)Gmat_1->data; 781d529f056Smarkadams4 matA_1 = (Mat_SeqAIJ *)mpimat_1->A->data; 782d529f056Smarkadams4 matB_1 = (Mat_SeqAIJ *)mpimat_1->B->data; 783d529f056Smarkadams4 784d529f056Smarkadams4 /* force compressed row storage for B matrix in AuxMat */ 785d529f056Smarkadams4 PetscCall(MatCheckCompressedRow(mpimat_1->B, matB_1->nonzerorowcnt, &matB_1->compressedrow, matB_1->i, Gmat_1->rmap->n, -1.0)); 786d529f056Smarkadams4 for (ix = 0; ix < matB_1->compressedrow.nrows; ix++) { 787d529f056Smarkadams4 PetscInt lid = matB_1->compressedrow.rindex[ix]; 788d529f056Smarkadams4 PetscCheck(lid <= nloc && lid >= -1, PETSC_COMM_SELF, PETSC_ERR_USER, "lid %d out of range. nloc = %d", (int)lid, (int)nloc); 789d529f056Smarkadams4 if (lid != -1) lid_cprowID_1[lid] = ix; 790d529f056Smarkadams4 } 791d529f056Smarkadams4 } else { 792d529f056Smarkadams4 PetscBool isAIJ; 793d529f056Smarkadams4 PetscCall(PetscStrbeginswith(((PetscObject)Gmat_1)->type_name, MATSEQAIJ, &isAIJ)); 794d529f056Smarkadams4 PetscCheck(isAIJ, PETSC_COMM_SELF, PETSC_ERR_USER, "Require AIJ matrix."); 795d529f056Smarkadams4 matA_1 = (Mat_SeqAIJ *)Gmat_1->data; 796d529f056Smarkadams4 } 797d529f056Smarkadams4 if (nloc > 0) { PetscCheck(!matB_1 || matB_1->compressedrow.use, PETSC_COMM_SELF, PETSC_ERR_PLIB, "matB_1 && !matB_1->compressedrow.use: PETSc bug???"); } 798d529f056Smarkadams4 /* get state of locals and selected gid for deleted */ 799d529f056Smarkadams4 for (lid = 0; lid < nloc; lid++) { 800d529f056Smarkadams4 lid_parent_gid[lid] = -1.0; 801d529f056Smarkadams4 lid_state[lid] = DELETED; 802d529f056Smarkadams4 } 803d529f056Smarkadams4 804d529f056Smarkadams4 /* set lid_state */ 805d529f056Smarkadams4 for (lid = 0; lid < nloc; lid++) { 806d529f056Smarkadams4 PetscCDIntNd *pos; 807d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos)); 808d529f056Smarkadams4 if (pos) { 809d529f056Smarkadams4 PetscInt gid1; 810d529f056Smarkadams4 811d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gid1)); 812d529f056Smarkadams4 PetscCheck(gid1 == lid + my0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "gid1 %d != lid %d + my0 %d", (int)gid1, (int)lid, (int)my0); 813d529f056Smarkadams4 lid_state[lid] = gid1; 814d529f056Smarkadams4 } 815d529f056Smarkadams4 } 816d529f056Smarkadams4 817d529f056Smarkadams4 /* map local to selected local, DELETED means a ghost owns it */ 818d529f056Smarkadams4 for (lid = kk = 0; lid < nloc; lid++) { 819d529f056Smarkadams4 NState state = lid_state[lid]; 820d529f056Smarkadams4 if (IS_SELECTED(state)) { 821d529f056Smarkadams4 PetscCDIntNd *pos; 822d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos)); 823d529f056Smarkadams4 while (pos) { 824d529f056Smarkadams4 PetscInt gid1; 825d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gid1)); 826d529f056Smarkadams4 PetscCall(PetscCDGetNextPos(aggs_2, lid, &pos)); 827d529f056Smarkadams4 if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1 - my0] = (PetscScalar)(lid + my0); 828d529f056Smarkadams4 } 829d529f056Smarkadams4 } 830d529f056Smarkadams4 } 831d529f056Smarkadams4 /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */ 832d529f056Smarkadams4 if (isMPI) { 833d529f056Smarkadams4 Vec tempVec; 834d529f056Smarkadams4 /* get 'cpcol_1_state' */ 835d529f056Smarkadams4 PetscCall(MatCreateVecs(Gmat_1, &tempVec, NULL)); 836d529f056Smarkadams4 for (kk = 0, j = my0; kk < nloc; kk++, j++) { 837d529f056Smarkadams4 PetscScalar v = (PetscScalar)lid_state[kk]; 838d529f056Smarkadams4 PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES)); 839d529f056Smarkadams4 } 840d529f056Smarkadams4 PetscCall(VecAssemblyBegin(tempVec)); 841d529f056Smarkadams4 PetscCall(VecAssemblyEnd(tempVec)); 842d529f056Smarkadams4 PetscCall(VecScatterBegin(mpimat_1->Mvctx, tempVec, mpimat_1->lvec, INSERT_VALUES, SCATTER_FORWARD)); 843d529f056Smarkadams4 PetscCall(VecScatterEnd(mpimat_1->Mvctx, tempVec, mpimat_1->lvec, INSERT_VALUES, SCATTER_FORWARD)); 844d529f056Smarkadams4 PetscCall(VecGetArray(mpimat_1->lvec, &cpcol_1_state)); 845d529f056Smarkadams4 /* get 'cpcol_2_state' */ 846d529f056Smarkadams4 PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, mpimat_2->lvec, INSERT_VALUES, SCATTER_FORWARD)); 847d529f056Smarkadams4 PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, mpimat_2->lvec, INSERT_VALUES, SCATTER_FORWARD)); 848d529f056Smarkadams4 PetscCall(VecGetArray(mpimat_2->lvec, &cpcol_2_state)); 849d529f056Smarkadams4 /* get 'cpcol_2_par_orig' */ 850d529f056Smarkadams4 for (kk = 0, j = my0; kk < nloc; kk++, j++) { 851d529f056Smarkadams4 PetscScalar v = (PetscScalar)lid_parent_gid[kk]; 852d529f056Smarkadams4 PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES)); 853d529f056Smarkadams4 } 854d529f056Smarkadams4 PetscCall(VecAssemblyBegin(tempVec)); 855d529f056Smarkadams4 PetscCall(VecAssemblyEnd(tempVec)); 856d529f056Smarkadams4 PetscCall(VecDuplicate(mpimat_2->lvec, &ghost_par_orig2)); 857d529f056Smarkadams4 PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghost_par_orig2, INSERT_VALUES, SCATTER_FORWARD)); 858d529f056Smarkadams4 PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghost_par_orig2, INSERT_VALUES, SCATTER_FORWARD)); 859d529f056Smarkadams4 PetscCall(VecGetArray(ghost_par_orig2, &cpcol_2_par_orig)); 860d529f056Smarkadams4 861d529f056Smarkadams4 PetscCall(VecDestroy(&tempVec)); 862d529f056Smarkadams4 } /* ismpi */ 863d529f056Smarkadams4 for (lid = 0; lid < nloc; lid++) { 864d529f056Smarkadams4 NState state = lid_state[lid]; 865d529f056Smarkadams4 if (IS_SELECTED(state)) { 866d529f056Smarkadams4 /* steal locals */ 867d529f056Smarkadams4 ii = matA_1->i; 868d529f056Smarkadams4 n = ii[lid + 1] - ii[lid]; 869d529f056Smarkadams4 idx = matA_1->j + ii[lid]; 870d529f056Smarkadams4 for (j = 0; j < n; j++) { 871d529f056Smarkadams4 PetscInt lidj = idx[j], sgid; 872d529f056Smarkadams4 NState statej = lid_state[lidj]; 873d529f056Smarkadams4 if (statej == DELETED && (sgid = (PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid + my0) { /* steal local */ 874d529f056Smarkadams4 lid_parent_gid[lidj] = (PetscScalar)(lid + my0); /* send this if sgid is not local */ 875d529f056Smarkadams4 if (sgid >= my0 && sgid < Iend) { /* I'm stealing this local from a local sgid */ 876d529f056Smarkadams4 PetscInt hav = 0, slid = sgid - my0, gidj = lidj + my0; 877d529f056Smarkadams4 PetscCDIntNd *pos, *last = NULL; 878d529f056Smarkadams4 /* looking for local from local so id_llist_2 works */ 879d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, slid, &pos)); 880d529f056Smarkadams4 while (pos) { 881d529f056Smarkadams4 PetscInt gid; 882d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gid)); 883d529f056Smarkadams4 if (gid == gidj) { 884d529f056Smarkadams4 PetscCheck(last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "last cannot be null"); 885d529f056Smarkadams4 PetscCall(PetscCDRemoveNextNode(aggs_2, slid, last)); 886d529f056Smarkadams4 PetscCall(PetscCDAppendNode(aggs_2, lid, pos)); 887d529f056Smarkadams4 hav = 1; 888d529f056Smarkadams4 break; 889d529f056Smarkadams4 } else last = pos; 890d529f056Smarkadams4 PetscCall(PetscCDGetNextPos(aggs_2, slid, &pos)); 891d529f056Smarkadams4 } 892d529f056Smarkadams4 if (hav != 1) { 893d529f056Smarkadams4 PetscCheck(hav, PETSC_COMM_SELF, PETSC_ERR_PLIB, "failed to find adj in 'selected' lists - structurally unsymmetric matrix"); 894d529f056Smarkadams4 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "found node %d times???", (int)hav); 895d529f056Smarkadams4 } 896d529f056Smarkadams4 } else { /* I'm stealing this local, owned by a ghost */ 897d529f056Smarkadams4 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.", 898d529f056Smarkadams4 ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : ""); 899d529f056Smarkadams4 PetscCall(PetscCDAppendID(aggs_2, lid, lidj + my0)); 900d529f056Smarkadams4 } 901d529f056Smarkadams4 } 902d529f056Smarkadams4 } /* local neighbors */ 903d529f056Smarkadams4 } else if (state == DELETED /* && lid_cprowID_1 */) { 904d529f056Smarkadams4 PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]); 905d529f056Smarkadams4 /* see if I have a selected ghost neighbor that will steal me */ 906d529f056Smarkadams4 if ((ix = lid_cprowID_1[lid]) != -1) { 907d529f056Smarkadams4 ii = matB_1->compressedrow.i; 908d529f056Smarkadams4 n = ii[ix + 1] - ii[ix]; 909d529f056Smarkadams4 idx = matB_1->j + ii[ix]; 910d529f056Smarkadams4 for (j = 0; j < n; j++) { 911d529f056Smarkadams4 PetscInt cpid = idx[j]; 912d529f056Smarkadams4 NState statej = (NState)PetscRealPart(cpcol_1_state[cpid]); 913d529f056Smarkadams4 if (IS_SELECTED(statej) && sgidold != (PetscInt)statej) { /* ghost will steal this, remove from my list */ 914d529f056Smarkadams4 lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */ 915d529f056Smarkadams4 if (sgidold >= my0 && sgidold < Iend) { /* this was mine */ 916d529f056Smarkadams4 PetscInt hav = 0, oldslidj = sgidold - my0; 917d529f056Smarkadams4 PetscCDIntNd *pos, *last = NULL; 918d529f056Smarkadams4 /* remove from 'oldslidj' list */ 919d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, oldslidj, &pos)); 920d529f056Smarkadams4 while (pos) { 921d529f056Smarkadams4 PetscInt gid; 922d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gid)); 923d529f056Smarkadams4 if (lid + my0 == gid) { 924d529f056Smarkadams4 /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove lid from oldslidj list *\/ */ 925d529f056Smarkadams4 PetscCheck(last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "last cannot be null"); 926d529f056Smarkadams4 PetscCall(PetscCDRemoveNextNode(aggs_2, oldslidj, last)); 927d529f056Smarkadams4 /* ghost (PetscScalar)statej will add this later */ 928d529f056Smarkadams4 hav = 1; 929d529f056Smarkadams4 break; 930d529f056Smarkadams4 } else last = pos; 931d529f056Smarkadams4 PetscCall(PetscCDGetNextPos(aggs_2, oldslidj, &pos)); 932d529f056Smarkadams4 } 933d529f056Smarkadams4 if (hav != 1) { 934d529f056Smarkadams4 PetscCheck(hav, PETSC_COMM_SELF, PETSC_ERR_PLIB, "failed to find (hav=%d) adj in 'selected' lists - structurally unsymmetric matrix", (int)hav); 935d529f056Smarkadams4 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "found node %d times???", (int)hav); 936d529f056Smarkadams4 } 937d529f056Smarkadams4 } else { 938d529f056Smarkadams4 /* TODO: ghosts remove this later */ 939d529f056Smarkadams4 } 940d529f056Smarkadams4 } 941d529f056Smarkadams4 } 942d529f056Smarkadams4 } 943d529f056Smarkadams4 } /* selected/deleted */ 944d529f056Smarkadams4 } /* node loop */ 945d529f056Smarkadams4 946d529f056Smarkadams4 if (isMPI) { 947d529f056Smarkadams4 PetscScalar *cpcol_2_parent, *cpcol_2_gid; 948d529f056Smarkadams4 Vec tempVec, ghostgids2, ghostparents2; 949d529f056Smarkadams4 PetscInt cpid, nghost_2; 950d529f056Smarkadams4 PCGAMGHashTable gid_cpid; 951d529f056Smarkadams4 952d529f056Smarkadams4 PetscCall(VecGetSize(mpimat_2->lvec, &nghost_2)); 953d529f056Smarkadams4 PetscCall(MatCreateVecs(Gmat_2, &tempVec, NULL)); 954d529f056Smarkadams4 955d529f056Smarkadams4 /* get 'cpcol_2_parent' */ 956d529f056Smarkadams4 for (kk = 0, j = my0; kk < nloc; kk++, j++) { PetscCall(VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES)); } 957d529f056Smarkadams4 PetscCall(VecAssemblyBegin(tempVec)); 958d529f056Smarkadams4 PetscCall(VecAssemblyEnd(tempVec)); 959d529f056Smarkadams4 PetscCall(VecDuplicate(mpimat_2->lvec, &ghostparents2)); 960d529f056Smarkadams4 PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghostparents2, INSERT_VALUES, SCATTER_FORWARD)); 961d529f056Smarkadams4 PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghostparents2, INSERT_VALUES, SCATTER_FORWARD)); 962d529f056Smarkadams4 PetscCall(VecGetArray(ghostparents2, &cpcol_2_parent)); 963d529f056Smarkadams4 964d529f056Smarkadams4 /* get 'cpcol_2_gid' */ 965d529f056Smarkadams4 for (kk = 0, j = my0; kk < nloc; kk++, j++) { 966d529f056Smarkadams4 PetscScalar v = (PetscScalar)j; 967d529f056Smarkadams4 PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES)); 968d529f056Smarkadams4 } 969d529f056Smarkadams4 PetscCall(VecAssemblyBegin(tempVec)); 970d529f056Smarkadams4 PetscCall(VecAssemblyEnd(tempVec)); 971d529f056Smarkadams4 PetscCall(VecDuplicate(mpimat_2->lvec, &ghostgids2)); 972d529f056Smarkadams4 PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghostgids2, INSERT_VALUES, SCATTER_FORWARD)); 973d529f056Smarkadams4 PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghostgids2, INSERT_VALUES, SCATTER_FORWARD)); 974d529f056Smarkadams4 PetscCall(VecGetArray(ghostgids2, &cpcol_2_gid)); 975d529f056Smarkadams4 PetscCall(VecDestroy(&tempVec)); 976d529f056Smarkadams4 977d529f056Smarkadams4 /* look for deleted ghosts and add to table */ 978d529f056Smarkadams4 PetscCall(PCGAMGHashTableCreate(2 * nghost_2 + 1, &gid_cpid)); 979d529f056Smarkadams4 for (cpid = 0; cpid < nghost_2; cpid++) { 980d529f056Smarkadams4 NState state = (NState)PetscRealPart(cpcol_2_state[cpid]); 981d529f056Smarkadams4 if (state == DELETED) { 982d529f056Smarkadams4 PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 983d529f056Smarkadams4 PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]); 984d529f056Smarkadams4 if (sgid_old == -1 && sgid_new != -1) { 985d529f056Smarkadams4 PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 986d529f056Smarkadams4 PetscCall(PCGAMGHashTableAdd(&gid_cpid, gid, cpid)); 987d529f056Smarkadams4 } 988d529f056Smarkadams4 } 989d529f056Smarkadams4 } 990d529f056Smarkadams4 991d529f056Smarkadams4 /* look for deleted ghosts and see if they moved - remove it */ 992d529f056Smarkadams4 for (lid = 0; lid < nloc; lid++) { 993d529f056Smarkadams4 NState state = lid_state[lid]; 994d529f056Smarkadams4 if (IS_SELECTED(state)) { 995d529f056Smarkadams4 PetscCDIntNd *pos, *last = NULL; 996d529f056Smarkadams4 /* look for deleted ghosts and see if they moved */ 997d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos)); 998d529f056Smarkadams4 while (pos) { 999d529f056Smarkadams4 PetscInt gid; 1000d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gid)); 1001d529f056Smarkadams4 1002d529f056Smarkadams4 if (gid < my0 || gid >= Iend) { 1003d529f056Smarkadams4 PetscCall(PCGAMGHashTableFind(&gid_cpid, gid, &cpid)); 1004d529f056Smarkadams4 if (cpid != -1) { 1005d529f056Smarkadams4 /* a moved ghost - */ 1006d529f056Smarkadams4 /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove 'flid' from list *\/ */ 1007d529f056Smarkadams4 PetscCall(PetscCDRemoveNextNode(aggs_2, lid, last)); 1008d529f056Smarkadams4 } else last = pos; 1009d529f056Smarkadams4 } else last = pos; 1010d529f056Smarkadams4 1011d529f056Smarkadams4 PetscCall(PetscCDGetNextPos(aggs_2, lid, &pos)); 1012d529f056Smarkadams4 } /* loop over list of deleted */ 1013d529f056Smarkadams4 } /* selected */ 1014d529f056Smarkadams4 } 1015d529f056Smarkadams4 PetscCall(PCGAMGHashTableDestroy(&gid_cpid)); 1016d529f056Smarkadams4 1017d529f056Smarkadams4 /* look at ghosts, see if they changed - and it */ 1018d529f056Smarkadams4 for (cpid = 0; cpid < nghost_2; cpid++) { 1019d529f056Smarkadams4 PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 1020d529f056Smarkadams4 if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */ 1021d529f056Smarkadams4 PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 1022d529f056Smarkadams4 PetscInt slid_new = sgid_new - my0, hav = 0; 1023d529f056Smarkadams4 PetscCDIntNd *pos; 1024d529f056Smarkadams4 1025d529f056Smarkadams4 /* search for this gid to see if I have it */ 1026d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, slid_new, &pos)); 1027d529f056Smarkadams4 while (pos) { 1028d529f056Smarkadams4 PetscInt gidj; 1029d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gidj)); 1030d529f056Smarkadams4 PetscCall(PetscCDGetNextPos(aggs_2, slid_new, &pos)); 1031d529f056Smarkadams4 1032d529f056Smarkadams4 if (gidj == gid) { 1033d529f056Smarkadams4 hav = 1; 1034d529f056Smarkadams4 break; 1035d529f056Smarkadams4 } 1036d529f056Smarkadams4 } 1037d529f056Smarkadams4 if (hav != 1) { 1038d529f056Smarkadams4 /* insert 'flidj' into head of llist */ 1039d529f056Smarkadams4 PetscCall(PetscCDAppendID(aggs_2, slid_new, gid)); 1040d529f056Smarkadams4 } 1041d529f056Smarkadams4 } 1042d529f056Smarkadams4 } 1043d529f056Smarkadams4 PetscCall(VecRestoreArray(mpimat_1->lvec, &cpcol_1_state)); 1044d529f056Smarkadams4 PetscCall(VecRestoreArray(mpimat_2->lvec, &cpcol_2_state)); 1045d529f056Smarkadams4 PetscCall(VecRestoreArray(ghostparents2, &cpcol_2_parent)); 1046d529f056Smarkadams4 PetscCall(VecRestoreArray(ghostgids2, &cpcol_2_gid)); 1047d529f056Smarkadams4 PetscCall(VecDestroy(&ghostgids2)); 1048d529f056Smarkadams4 PetscCall(VecDestroy(&ghostparents2)); 1049d529f056Smarkadams4 PetscCall(VecDestroy(&ghost_par_orig2)); 1050d529f056Smarkadams4 } 1051d529f056Smarkadams4 PetscCall(PetscFree3(lid_state, lid_parent_gid, lid_cprowID_1)); 1052d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 1053d529f056Smarkadams4 } 1054d529f056Smarkadams4 1055c8b0795cSMark F. Adams /* 1056bae903cbSmarkadams4 PCGAMGCoarsen_AGG - supports squaring the graph (deprecated) and new graph for 1057bae903cbSmarkadams4 communication of QR data used with HEM and MISk coarsening 1058c8b0795cSMark F. Adams 1059c8b0795cSMark F. Adams Input Parameter: 1060e0940f08SMark F. Adams . a_pc - this 1061bae903cbSmarkadams4 1062e0940f08SMark F. Adams Input/Output Parameter: 1063bae903cbSmarkadams4 . a_Gmat1 - graph to coarsen (in), graph off processor edges for QR gather scatter (out) 1064bae903cbSmarkadams4 1065c8b0795cSMark F. Adams Output Parameter: 10660cbbd2e1SMark F. Adams . agg_lists - list of aggregates 1067bae903cbSmarkadams4 1068c8b0795cSMark F. Adams */ 1069d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc, Mat *a_Gmat1, PetscCoarsenData **agg_lists) 1070d71ae5a4SJacob Faibussowitsch { 1071e0940f08SMark F. Adams PC_MG *mg = (PC_MG *)a_pc->data; 1072c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1073c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 10748926f930SMark Adams Mat Gmat2, Gmat1 = *a_Gmat1; /* aggressive graph */ 10750cbbd2e1SMark F. Adams IS perm; 1076bae903cbSmarkadams4 PetscInt Istart, Iend, Ii, nloc, bs, nn; 1077bae903cbSmarkadams4 PetscInt *permute, *degree; 1078c8b0795cSMark F. Adams PetscBool *bIndexSet; 1079aea53286SMark Adams PetscReal hashfact; 1080e2c00dcbSBarry Smith PetscInt iSwapIndex; 10813b50add6SMark Adams PetscRandom random; 1082d529f056Smarkadams4 MPI_Comm comm; 1083c8b0795cSMark F. Adams 1084c8b0795cSMark F. Adams PetscFunctionBegin; 1085d529f056Smarkadams4 PetscCall(PetscObjectGetComm((PetscObject)Gmat1, &comm)); 1086849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0)); 1087bae903cbSmarkadams4 PetscCall(MatGetLocalSize(Gmat1, &nn, NULL)); 10889566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Gmat1, &bs)); 108963a3b9bcSJacob Faibussowitsch PetscCheck(bs == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bs %" PetscInt_FMT " must be 1", bs); 1090bae903cbSmarkadams4 nloc = nn / bs; 10915cfd4bd9SMark Adams /* get MIS aggs - randomize */ 1092bae903cbSmarkadams4 PetscCall(PetscMalloc2(nloc, &permute, nloc, °ree)); 10939566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nloc, &bIndexSet)); 10946e09b0e3SStefano Zampini for (Ii = 0; Ii < nloc; Ii++) permute[Ii] = Ii; 10959566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &random)); 10969566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Gmat1, &Istart, &Iend)); 1097c8b0795cSMark F. Adams for (Ii = 0; Ii < nloc; Ii++) { 1098bae903cbSmarkadams4 PetscInt nc; 1099bae903cbSmarkadams4 PetscCall(MatGetRow(Gmat1, Istart + Ii, &nc, NULL, NULL)); 1100bae903cbSmarkadams4 degree[Ii] = nc; 1101bae903cbSmarkadams4 PetscCall(MatRestoreRow(Gmat1, Istart + Ii, &nc, NULL, NULL)); 1102bae903cbSmarkadams4 } 1103bae903cbSmarkadams4 for (Ii = 0; Ii < nloc; Ii++) { 11049566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(random, &hashfact)); 1105aea53286SMark Adams iSwapIndex = (PetscInt)(hashfact * nloc) % nloc; 1106c8b0795cSMark F. Adams if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) { 1107c8b0795cSMark F. Adams PetscInt iTemp = permute[iSwapIndex]; 1108c8b0795cSMark F. Adams permute[iSwapIndex] = permute[Ii]; 1109c8b0795cSMark F. Adams permute[Ii] = iTemp; 1110bae903cbSmarkadams4 iTemp = degree[iSwapIndex]; 1111bae903cbSmarkadams4 degree[iSwapIndex] = degree[Ii]; 1112bae903cbSmarkadams4 degree[Ii] = iTemp; 1113c8b0795cSMark F. Adams bIndexSet[iSwapIndex] = PETSC_TRUE; 1114c8b0795cSMark F. Adams } 1115c8b0795cSMark F. Adams } 1116d529f056Smarkadams4 // apply minimum degree ordering -- NEW 1117d529f056Smarkadams4 if (pc_gamg_agg->use_minimum_degree_ordering) { PetscCall(PetscSortIntWithArray(nloc, degree, permute)); } 11189566063dSJacob Faibussowitsch PetscCall(PetscFree(bIndexSet)); 11199566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&random)); 11209566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm)); 1121849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0)); 1122d529f056Smarkadams4 // square graph 1123d529f056Smarkadams4 if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels && pc_gamg_agg->use_aggressive_square_graph) { 1124d529f056Smarkadams4 PetscCall(PCGAMGSquareGraph_GAMG(a_pc, Gmat1, &Gmat2)); 1125d529f056Smarkadams4 } else Gmat2 = Gmat1; 1126d529f056Smarkadams4 // switch to old MIS-1 for square graph 1127d529f056Smarkadams4 if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels) { 1128d529f056Smarkadams4 if (!pc_gamg_agg->use_aggressive_square_graph) PetscCall(MatCoarsenMISKSetDistance(pc_gamg_agg->crs, pc_gamg_agg->aggressive_mis_k)); // hardwire to MIS-2 1129d529f056Smarkadams4 else PetscCall(MatCoarsenSetType(pc_gamg_agg->crs, MATCOARSENMIS)); // old MIS -- side effect 1130d529f056Smarkadams4 } else if (pc_gamg_agg->use_aggressive_square_graph && pc_gamg_agg->aggressive_coarsening_levels > 0) { // we reset the MIS 1131d529f056Smarkadams4 const char *prefix; 1132d529f056Smarkadams4 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)a_pc, &prefix)); 1133d529f056Smarkadams4 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc_gamg_agg->crs, prefix)); 1134d529f056Smarkadams4 PetscCall(MatCoarsenSetFromOptions(pc_gamg_agg->crs)); // get the default back on non-aggressive levels when square graph switched to old MIS 1135d529f056Smarkadams4 } 1136d529f056Smarkadams4 PetscCall(MatCoarsenSetAdjacency(pc_gamg_agg->crs, Gmat2)); 11372d776b49SBarry Smith PetscCall(MatCoarsenSetStrictAggs(pc_gamg_agg->crs, PETSC_TRUE)); 1138d529f056Smarkadams4 PetscCall(MatCoarsenSetGreedyOrdering(pc_gamg_agg->crs, perm)); 11392d776b49SBarry Smith PetscCall(MatCoarsenApply(pc_gamg_agg->crs)); 11402d776b49SBarry Smith PetscCall(MatCoarsenViewFromOptions(pc_gamg_agg->crs, NULL, "-mat_coarsen_view")); 11412d776b49SBarry Smith PetscCall(MatCoarsenGetData(pc_gamg_agg->crs, agg_lists)); /* output */ 11422d776b49SBarry Smith PetscCall(MatCoarsenDestroy(&pc_gamg_agg->crs)); 1143b43b03e9SMark F. Adams 11449566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 1145bae903cbSmarkadams4 PetscCall(PetscFree2(permute, degree)); 1146849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0)); 11479f3f12c8SMark F. Adams 11489c15c1aeSMark Adams if (Gmat2 != Gmat1) { // square graph, we need ghosts for selected 1149d529f056Smarkadams4 PetscCoarsenData *llist = *agg_lists; 1150d529f056Smarkadams4 PetscCall(fixAggregatesWithSquare(a_pc, Gmat2, Gmat1, *agg_lists)); 1151d529f056Smarkadams4 PetscCall(MatDestroy(&Gmat1)); 1152d529f056Smarkadams4 *a_Gmat1 = Gmat2; /* output */ 11538926f930SMark Adams PetscCall(PetscCDSetMat(llist, *a_Gmat1)); /* Need a graph with ghosts here */ 11540cbbd2e1SMark F. Adams } 1155849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0)); 11563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1157c8b0795cSMark F. Adams } 1158c8b0795cSMark F. Adams 1159c8b0795cSMark F. Adams /* 11600cbbd2e1SMark F. Adams PCGAMGProlongator_AGG 1161c8b0795cSMark F. Adams 1162c8b0795cSMark F. Adams Input Parameter: 1163c8b0795cSMark F. Adams . pc - this 1164c8b0795cSMark F. Adams . Amat - matrix on this fine level 1165c8b0795cSMark F. Adams . Graph - used to get ghost data for nodes in 11660cbbd2e1SMark F. Adams . agg_lists - list of aggregates 1167c8b0795cSMark F. Adams Output Parameter: 1168c8b0795cSMark F. Adams . a_P_out - prolongation operator to the next level 1169c8b0795cSMark F. Adams */ 11708926f930SMark Adams static PetscErrorCode PCGAMGProlongator_AGG(PC pc, Mat Amat, PetscCoarsenData *agg_lists, Mat *a_P_out) 1171d71ae5a4SJacob Faibussowitsch { 11722e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 11732e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 11744a2f8832SBarry Smith const PetscInt col_bs = pc_gamg->data_cell_cols; 1175c8b0795cSMark F. Adams PetscInt Istart, Iend, nloc, ii, jj, kk, my0, nLocalSelected, bs; 11768926f930SMark Adams Mat Gmat, Prol; 1177d5d25401SBarry Smith PetscMPIInt size; 11783b4367a7SBarry Smith MPI_Comm comm; 1179c8b0795cSMark F. Adams PetscReal *data_w_ghost; 1180c8b0795cSMark F. Adams PetscInt myCrs0, nbnodes = 0, *flid_fgid; 1181d9558ea9SBarry Smith MatType mtype; 11822e68589bSMark F. Adams 11832e68589bSMark F. Adams PetscFunctionBegin; 11849566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm)); 118508401ef6SPierre Jolivet PetscCheck(col_bs >= 1, comm, PETSC_ERR_PLIB, "Column bs cannot be less than 1"); 1186849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0)); 11879566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 11889566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Amat, &Istart, &Iend)); 11899566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Amat, &bs)); 11909371c9d4SSatish Balay nloc = (Iend - Istart) / bs; 11919371c9d4SSatish Balay my0 = Istart / bs; 119263a3b9bcSJacob 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); 1193d8b4a066SPierre Jolivet PetscCall(PetscCDGetMat(agg_lists, &Gmat)); // get auxiliary matrix for ghost edges for size > 1 11942e68589bSMark F. Adams 11952e68589bSMark F. Adams /* get 'nLocalSelected' */ 11960cbbd2e1SMark F. Adams for (ii = 0, nLocalSelected = 0; ii < nloc; ii++) { 1197e78576d6SMark F. Adams PetscBool ise; 11980cbbd2e1SMark F. Adams /* filter out singletons 0 or 1? */ 11995e62d202SMark Adams PetscCall(PetscCDIsEmptyAt(agg_lists, ii, &ise)); 1200e78576d6SMark F. Adams if (!ise) nLocalSelected++; 12012e68589bSMark F. Adams } 12022e68589bSMark F. Adams 12032e68589bSMark F. Adams /* create prolongator, create P matrix */ 12049566063dSJacob Faibussowitsch PetscCall(MatGetType(Amat, &mtype)); 12059566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &Prol)); 12069566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Prol, nloc * bs, nLocalSelected * col_bs, PETSC_DETERMINE, PETSC_DETERMINE)); 120744eff04eSMark Adams PetscCall(MatSetBlockSizes(Prol, bs, col_bs)); // should this be before MatSetSizes? 12089566063dSJacob Faibussowitsch PetscCall(MatSetType(Prol, mtype)); 12093742c8caSstefanozampini #if PetscDefined(HAVE_DEVICE) 12103742c8caSstefanozampini PetscBool flg; 12113742c8caSstefanozampini PetscCall(MatBoundToCPU(Amat, &flg)); 12123742c8caSstefanozampini PetscCall(MatBindToCPU(Prol, flg)); 12133742c8caSstefanozampini if (flg) PetscCall(MatSetBindingPropagates(Prol, PETSC_TRUE)); 12143742c8caSstefanozampini #endif 12159566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(Prol, col_bs, NULL)); 12169566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(Prol, col_bs, NULL, col_bs, NULL)); 12172e68589bSMark F. Adams 12182e68589bSMark F. Adams /* can get all points "removed" */ 12199566063dSJacob Faibussowitsch PetscCall(MatGetSize(Prol, &kk, &ii)); 12207f66b68fSMark Adams if (!ii) { 122163a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: No selected points on coarse grid\n", ((PetscObject)pc)->prefix)); 12229566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Prol)); 12230298fd71SBarry Smith *a_P_out = NULL; /* out */ 1224849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0)); 12253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12262e68589bSMark F. Adams } 122763a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: New grid %" PetscInt_FMT " nodes\n", ((PetscObject)pc)->prefix, ii / col_bs)); 12289566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk)); 12290cbbd2e1SMark F. Adams 123063a3b9bcSJacob 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); 1231c8b0795cSMark F. Adams myCrs0 = myCrs0 / col_bs; 123263a3b9bcSJacob 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); 12332e68589bSMark F. Adams 12342e68589bSMark F. Adams /* create global vector of data in 'data_w_ghost' */ 1235849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0)); 1236bae903cbSmarkadams4 if (size > 1) { /* get ghost null space data */ 12372e68589bSMark F. Adams PetscReal *tmp_gdata, *tmp_ldata, *tp2; 12389566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc, &tmp_ldata)); 12394a2f8832SBarry Smith for (jj = 0; jj < col_bs; jj++) { 1240c8b0795cSMark F. Adams for (kk = 0; kk < bs; kk++) { 1241a2f3521dSMark F. Adams PetscInt ii, stride; 12428e3a54c0SPierre Jolivet const PetscReal *tp = PetscSafePointerPlusOffset(pc_gamg->data, jj * bs * nloc + kk); 12432fa5cd67SKarl Rupp for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp; 12442fa5cd67SKarl Rupp 12459566063dSJacob Faibussowitsch PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata)); 1246a2f3521dSMark F. Adams 1247bae903cbSmarkadams4 if (!jj && !kk) { /* now I know how many total nodes - allocate TODO: move below and do in one 'col_bs' call */ 12489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(stride * bs * col_bs, &data_w_ghost)); 1249a2f3521dSMark F. Adams nbnodes = bs * stride; 12502e68589bSMark F. Adams } 12518e3a54c0SPierre Jolivet tp2 = PetscSafePointerPlusOffset(data_w_ghost, jj * bs * stride + kk); 12522fa5cd67SKarl Rupp for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii]; 12539566063dSJacob Faibussowitsch PetscCall(PetscFree(tmp_gdata)); 12542e68589bSMark F. Adams } 12552e68589bSMark F. Adams } 12569566063dSJacob Faibussowitsch PetscCall(PetscFree(tmp_ldata)); 1257806fa848SBarry Smith } else { 1258c8b0795cSMark F. Adams nbnodes = bs * nloc; 1259c8b0795cSMark F. Adams data_w_ghost = (PetscReal *)pc_gamg->data; 12602e68589bSMark F. Adams } 12612e68589bSMark F. Adams 1262bae903cbSmarkadams4 /* get 'flid_fgid' TODO - move up to get 'stride' and do get null space data above in one step (jj loop) */ 1263c5df96a5SBarry Smith if (size > 1) { 12642e68589bSMark F. Adams PetscReal *fid_glid_loc, *fiddata; 1265a2f3521dSMark F. Adams PetscInt stride; 12662e68589bSMark F. Adams 12679566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc, &fid_glid_loc)); 12682e68589bSMark F. Adams for (kk = 0; kk < nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0 + kk); 12699566063dSJacob Faibussowitsch PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata)); 1270bae903cbSmarkadams4 PetscCall(PetscMalloc1(stride, &flid_fgid)); /* copy real data to in */ 1271a2f3521dSMark F. Adams for (kk = 0; kk < stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk]; 12729566063dSJacob Faibussowitsch PetscCall(PetscFree(fiddata)); 1273a2f3521dSMark F. Adams 127463a3b9bcSJacob Faibussowitsch PetscCheck(stride == nbnodes / bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "stride %" PetscInt_FMT " != nbnodes %" PetscInt_FMT "/bs %" PetscInt_FMT, stride, nbnodes, bs); 12759566063dSJacob Faibussowitsch PetscCall(PetscFree(fid_glid_loc)); 1276806fa848SBarry Smith } else { 12779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc, &flid_fgid)); 12782e68589bSMark F. Adams for (kk = 0; kk < nloc; kk++) flid_fgid[kk] = my0 + kk; 12792e68589bSMark F. Adams } 1280849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0)); 1281bae903cbSmarkadams4 /* get P0 */ 1282849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0)); 1283c8b0795cSMark F. Adams { 12840298fd71SBarry Smith PetscReal *data_out = NULL; 12859566063dSJacob Faibussowitsch PetscCall(formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes, data_w_ghost, flid_fgid, &data_out, Prol)); 12869566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 12872fa5cd67SKarl Rupp 1288c8b0795cSMark F. Adams pc_gamg->data = data_out; 12894a2f8832SBarry Smith pc_gamg->data_cell_rows = col_bs; 12904a2f8832SBarry Smith pc_gamg->data_sz = col_bs * col_bs * nLocalSelected; 1291c8b0795cSMark F. Adams } 1292849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0)); 129348a46eb9SPierre Jolivet if (size > 1) PetscCall(PetscFree(data_w_ghost)); 12949566063dSJacob Faibussowitsch PetscCall(PetscFree(flid_fgid)); 12952e68589bSMark F. Adams 1296c8b0795cSMark F. Adams *a_P_out = Prol; /* out */ 12975e62d202SMark Adams PetscCall(MatViewFromOptions(Prol, NULL, "-view_P")); 1298ed4e82eaSPeter Brune 1299849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0)); 13003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1301c8b0795cSMark F. Adams } 1302c8b0795cSMark F. Adams 1303c8b0795cSMark F. Adams /* 1304fd1112cbSBarry Smith PCGAMGOptProlongator_AGG 1305c8b0795cSMark F. Adams 1306c8b0795cSMark F. Adams Input Parameter: 1307c8b0795cSMark F. Adams . pc - this 1308c8b0795cSMark F. Adams . Amat - matrix on this fine level 1309c8b0795cSMark F. Adams In/Output Parameter: 13102a808120SBarry Smith . a_P - prolongation operator to the next level 1311c8b0795cSMark F. Adams */ 1312d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGOptProlongator_AGG(PC pc, Mat Amat, Mat *a_P) 1313d71ae5a4SJacob Faibussowitsch { 1314c8b0795cSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 1315c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1316c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 1317c8b0795cSMark F. Adams PetscInt jj; 1318c8b0795cSMark F. Adams Mat Prol = *a_P; 13193b4367a7SBarry Smith MPI_Comm comm; 13202a808120SBarry Smith KSP eksp; 13212a808120SBarry Smith Vec bb, xx; 13222a808120SBarry Smith PC epc; 13232a808120SBarry Smith PetscReal alpha, emax, emin; 1324c8b0795cSMark F. Adams 1325c8b0795cSMark F. Adams PetscFunctionBegin; 13269566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm)); 1327849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0)); 1328c8b0795cSMark F. Adams 1329d5d25401SBarry Smith /* compute maximum singular value of operator to be used in smoother */ 13302a808120SBarry Smith if (0 < pc_gamg_agg->nsmooths) { 133118c3aa7eSMark /* get eigen estimates */ 133218c3aa7eSMark if (pc_gamg->emax > 0) { 133318c3aa7eSMark emin = pc_gamg->emin; 133418c3aa7eSMark emax = pc_gamg->emax; 133518c3aa7eSMark } else { 13360ed2132dSStefano Zampini const char *prefix; 13370ed2132dSStefano Zampini 13389566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Amat, &bb, NULL)); 13399566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Amat, &xx, NULL)); 13409566063dSJacob Faibussowitsch PetscCall(KSPSetNoisy_Private(bb)); 1341e696ed0bSMark F. Adams 13429566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm, &eksp)); 13433821be0aSBarry Smith PetscCall(KSPSetNestLevel(eksp, pc->kspnestlevel)); 13449566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 13459566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(eksp, prefix)); 13469566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(eksp, "pc_gamg_esteig_")); 134773f7197eSJed Brown { 1348b94d7dedSBarry Smith PetscBool isset, sflg; 1349b94d7dedSBarry Smith PetscCall(MatIsSPDKnown(Amat, &isset, &sflg)); 1350b94d7dedSBarry Smith if (isset && sflg) PetscCall(KSPSetType(eksp, KSPCG)); 1351d24ecf33SMark } 13529566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(eksp, pc->erroriffailure)); 13539566063dSJacob Faibussowitsch PetscCall(KSPSetNormType(eksp, KSP_NORM_NONE)); 1354c2436cacSMark F. Adams 13559566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(eksp, PETSC_FALSE)); 13569566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(eksp, Amat, Amat)); 13572e68589bSMark F. Adams 13589566063dSJacob Faibussowitsch PetscCall(KSPGetPC(eksp, &epc)); 13599566063dSJacob Faibussowitsch PetscCall(PCSetType(epc, PCJACOBI)); /* smoother in smoothed agg. */ 1360c2436cacSMark F. Adams 13619566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10)); // 10 is safer, but 5 is often fine, can override with -pc_gamg_esteig_ksp_max_it -mg_levels_ksp_chebyshev_esteig 0,0.25,0,1.2 1362074aec55SMark Adams 13639566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(eksp)); 13649566063dSJacob Faibussowitsch PetscCall(KSPSetComputeSingularValues(eksp, PETSC_TRUE)); 13659566063dSJacob Faibussowitsch PetscCall(KSPSolve(eksp, bb, xx)); 13669566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(eksp, pc, xx)); 13672e68589bSMark F. Adams 13689566063dSJacob Faibussowitsch PetscCall(KSPComputeExtremeSingularValues(eksp, &emax, &emin)); 13699566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: Smooth P0: max eigen=%e min=%e PC=%s\n", ((PetscObject)pc)->prefix, (double)emax, (double)emin, PCJACOBI)); 13709566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xx)); 13719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bb)); 13729566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&eksp)); 13732e68589bSMark F. Adams } 137418c3aa7eSMark if (pc_gamg->use_sa_esteig) { 137518c3aa7eSMark mg->min_eigen_DinvA[pc_gamg->current_level] = emin; 137618c3aa7eSMark mg->max_eigen_DinvA[pc_gamg->current_level] = emax; 137763a3b9bcSJacob 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)); 137818c3aa7eSMark } else { 137918c3aa7eSMark mg->min_eigen_DinvA[pc_gamg->current_level] = 0; 138018c3aa7eSMark mg->max_eigen_DinvA[pc_gamg->current_level] = 0; 138118c3aa7eSMark } 138218c3aa7eSMark } else { 138318c3aa7eSMark mg->min_eigen_DinvA[pc_gamg->current_level] = 0; 138418c3aa7eSMark mg->max_eigen_DinvA[pc_gamg->current_level] = 0; 138518c3aa7eSMark } 13862e68589bSMark F. Adams 13872a808120SBarry Smith /* smooth P0 */ 13882a808120SBarry Smith for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) { 13892a808120SBarry Smith Mat tMat; 13902a808120SBarry Smith Vec diag; 13912a808120SBarry Smith 1392849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0)); 13932a808120SBarry Smith 13942e68589bSMark F. Adams /* smooth P1 := (I - omega/lam D^{-1}A)P0 */ 13959566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0)); 13969566063dSJacob Faibussowitsch PetscCall(MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat)); 13979566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0)); 13989566063dSJacob Faibussowitsch PetscCall(MatProductClear(tMat)); 13999566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Amat, &diag, NULL)); 14009566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(Amat, diag)); /* effectively PCJACOBI */ 14019566063dSJacob Faibussowitsch PetscCall(VecReciprocal(diag)); 14029566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(tMat, diag, NULL)); 14039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&diag)); 1404b4da3a1bSStefano Zampini 1405d5d25401SBarry Smith /* TODO: Set a PCFailedReason and exit the building of the AMG preconditioner */ 140608401ef6SPierre Jolivet PetscCheck(emax != 0.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Computed maximum singular value as zero"); 1407d5d25401SBarry Smith /* TODO: Document the 1.4 and don't hardwire it in this routine */ 1408b4ec6429SMark F. Adams alpha = -1.4 / emax; 1409b4da3a1bSStefano Zampini 14109566063dSJacob Faibussowitsch PetscCall(MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN)); 14119566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Prol)); 14122e68589bSMark F. Adams Prol = tMat; 1413849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0)); 14142e68589bSMark F. Adams } 1415849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0)); 1416c8b0795cSMark F. Adams *a_P = Prol; 14173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14182e68589bSMark F. Adams } 14192e68589bSMark F. Adams 1420*a077d33dSBarry Smith /*MC 1421*a077d33dSBarry Smith PCGAMGAGG - Smooth aggregation, {cite}`vanek1996algebraic`, {cite}`vanek2001convergence`, variant of PETSc's algebraic multigrid (`PCGAMG`) preconditioner 14222e68589bSMark F. Adams 1423*a077d33dSBarry Smith Options Database Keys: 1424*a077d33dSBarry Smith + -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation to construct prolongation 1425*a077d33dSBarry Smith . -pc_gamg_aggressive_coarsening <n,default=1> - number of aggressive coarsening (MIS-2) levels from finest. 1426*a077d33dSBarry Smith . -pc_gamg_aggressive_square_graph <bool,default=false> - Use square graph (A'A) or MIS-k (k=2) for aggressive coarsening 1427*a077d33dSBarry Smith . -pc_gamg_mis_k_minimum_degree_ordering <bool,default=true> - Use minimum degree ordering in greedy MIS algorithm 1428*a077d33dSBarry Smith . -pc_gamg_pc_gamg_asm_hem_aggs <n,default=0> - Number of HEM aggregation steps for ASM smoother 1429*a077d33dSBarry Smith - -pc_gamg_aggressive_mis_k <n,default=2> - Number (k) distance in MIS coarsening (>2 is 'aggressive') 1430*a077d33dSBarry Smith 1431*a077d33dSBarry Smith Level: intermediate 1432*a077d33dSBarry Smith 1433*a077d33dSBarry Smith Notes: 1434*a077d33dSBarry Smith To obtain good performance for `PCGAMG` for vector valued problems you must 1435*a077d33dSBarry Smith call `MatSetBlockSize()` to indicate the number of degrees of freedom per grid point. 1436*a077d33dSBarry Smith Call `MatSetNearNullSpace()` (or `PCSetCoordinates()` if solving the equations of elasticity) to indicate the near null space of the operator 1437*a077d33dSBarry Smith 1438*a077d33dSBarry Smith The many options for `PCMG` and `PCGAMG` such as controlling the smoothers on each level etc. also work for `PCGAMGAGG` 1439*a077d33dSBarry Smith 1440*a077d33dSBarry Smith .seealso: `PCGAMG`, [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCCreate()`, `PCSetType()`, 1441*a077d33dSBarry Smith `MatSetBlockSize()`, `PCMGType`, `PCSetCoordinates()`, `MatSetNearNullSpace()`, `PCGAMGSetType()`, 1442*a077d33dSBarry Smith `PCGAMGAGG`, `PCGAMGGEO`, `PCGAMGCLASSICAL`, `PCGAMGSetProcEqLim()`, `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRepartition()`, `PCGAMGRegister()`, 1443*a077d33dSBarry Smith `PCGAMGSetReuseInterpolation()`, `PCGAMGASMSetUseAggs()`, `PCGAMGSetParallelCoarseGridSolve()`, `PCGAMGSetNlevels()`, `PCGAMGSetThreshold()`, 1444*a077d33dSBarry Smith `PCGAMGGetType()`, `PCGAMGSetUseSAEstEig()` 1445*a077d33dSBarry Smith M*/ 1446d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCreateGAMG_AGG(PC pc) 1447d71ae5a4SJacob Faibussowitsch { 1448c8b0795cSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 1449c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1450c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg; 14512e68589bSMark F. Adams 1452c8b0795cSMark F. Adams PetscFunctionBegin; 1453c8b0795cSMark F. Adams /* create sub context for SA */ 14544dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&pc_gamg_agg)); 1455c8b0795cSMark F. Adams pc_gamg->subctx = pc_gamg_agg; 1456c8b0795cSMark F. Adams 14571ab5ffc9SJed Brown pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG; 14589b8ffb57SJed Brown pc_gamg->ops->destroy = PCDestroy_GAMG_AGG; 1459c8b0795cSMark F. Adams /* reset does not do anything; setup not virtual */ 1460c8b0795cSMark F. Adams 1461c8b0795cSMark F. Adams /* set internal function pointers */ 14622d776b49SBarry Smith pc_gamg->ops->creategraph = PCGAMGCreateGraph_AGG; 14631ab5ffc9SJed Brown pc_gamg->ops->coarsen = PCGAMGCoarsen_AGG; 14641ab5ffc9SJed Brown pc_gamg->ops->prolongator = PCGAMGProlongator_AGG; 1465fd1112cbSBarry Smith pc_gamg->ops->optprolongator = PCGAMGOptProlongator_AGG; 14661ab5ffc9SJed Brown pc_gamg->ops->createdefaultdata = PCSetData_AGG; 14675adeb434SBarry Smith pc_gamg->ops->view = PCView_GAMG_AGG; 1468c8b0795cSMark F. Adams 1469cab9ed1eSBarry Smith pc_gamg_agg->nsmooths = 1; 1470d529f056Smarkadams4 pc_gamg_agg->aggressive_coarsening_levels = 1; 1471d4adc10fSMark Adams pc_gamg_agg->use_aggressive_square_graph = PETSC_TRUE; 1472d529f056Smarkadams4 pc_gamg_agg->use_minimum_degree_ordering = PETSC_FALSE; 1473a9ccf005SMark Adams pc_gamg_agg->use_low_mem_filter = PETSC_FALSE; 1474d529f056Smarkadams4 pc_gamg_agg->aggressive_mis_k = 2; 1475cab9ed1eSBarry Smith 14769566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", PCGAMGSetNSmooths_AGG)); 1477bae903cbSmarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", PCGAMGSetAggressiveLevels_AGG)); 1478d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveSquareGraph_C", PCGAMGSetAggressiveSquareGraph_AGG)); 1479d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetMinDegreeOrdering_C", PCGAMGMISkSetMinDegreeOrdering_AGG)); 1480a9ccf005SMark Adams PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetLowMemoryFilter_C", PCGAMGSetLowMemoryFilter_AGG)); 1481d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetAggressive_C", PCGAMGMISkSetAggressive_AGG)); 14829566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_AGG)); 14833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14842e68589bSMark F. Adams } 1485