12e68589bSMark F. Adams /* 2b817416eSBarry Smith GAMG geometric-algebric multigrid PC - Mark Adams 2011 32e68589bSMark F. Adams */ 42e68589bSMark F. Adams 52e68589bSMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 62e68589bSMark F. Adams #include <petscblaslapack.h> 7539c167fSBarry Smith #include <petscdm.h> 873f7197eSJed Brown #include <petsc/private/kspimpl.h> 92e68589bSMark F. Adams 102e68589bSMark F. Adams typedef struct { 11c8b0795cSMark F. Adams PetscInt nsmooths; 12bae903cbSmarkadams4 PetscInt aggressive_coarsening_levels; // number of aggressive coarsening levels (square or MISk) 13d529f056Smarkadams4 PetscInt aggressive_mis_k; // the k in MIS-k 14d529f056Smarkadams4 PetscBool use_aggressive_square_graph; 15d529f056Smarkadams4 PetscBool use_minimum_degree_ordering; 16a9ccf005SMark Adams PetscBool use_low_mem_filter; 172d776b49SBarry Smith MatCoarsen crs; 182e68589bSMark F. Adams } PC_GAMG_AGG; 192e68589bSMark F. Adams 202e68589bSMark F. Adams /*@ 21f1580f4eSBarry Smith PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical) used for multigrid on all the levels 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: 30cab9ed1eSBarry Smith . -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation 312e68589bSMark F. Adams 322e68589bSMark F. Adams Level: intermediate 332e68589bSMark F. Adams 34562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCGAMG` 352e68589bSMark F. Adams @*/ 36d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n) 37d71ae5a4SJacob Faibussowitsch { 382e68589bSMark F. Adams PetscFunctionBegin; 392e68589bSMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 40d5d25401SBarry Smith PetscValidLogicalCollectiveInt(pc, n, 2); 41cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGSetNSmooths_C", (PC, PetscInt), (pc, n)); 423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 432e68589bSMark F. Adams } 442e68589bSMark F. Adams 45d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetNSmooths_AGG(PC pc, PetscInt n) 46d71ae5a4SJacob Faibussowitsch { 472e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 482e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 49c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 502e68589bSMark F. Adams 512e68589bSMark F. Adams PetscFunctionBegin; 52c8b0795cSMark F. Adams pc_gamg_agg->nsmooths = n; 533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 54c8b0795cSMark F. Adams } 55c8b0795cSMark F. Adams 56c8b0795cSMark F. Adams /*@ 57f1580f4eSBarry Smith PCGAMGSetAggressiveLevels - Use aggressive coarsening on first n levels 58ef4ad70eSMark F. Adams 59c3339decSBarry Smith Logically Collective 60ef4ad70eSMark F. Adams 61ef4ad70eSMark F. Adams Input Parameters: 62cab9ed1eSBarry Smith + pc - the preconditioner context 63d5d25401SBarry Smith - n - 0, 1 or more 64ef4ad70eSMark F. Adams 65ef4ad70eSMark F. Adams Options Database Key: 66bae903cbSmarkadams4 . -pc_gamg_aggressive_coarsening <n,default = 1> - Number of levels to square the graph on before aggregating it 67af3c827dSMark Adams 68ef4ad70eSMark F. Adams Level: intermediate 69ef4ad70eSMark F. Adams 70cbb74892SSatish Balay .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGMISkSetAggressive()`, `PCGAMGSetAggressiveSquareGraph()`, `PCGAMGMISkSetMinDegreeOrdering()`, `PCGAMGSetLowMemoryFilter()` 71ef4ad70eSMark F. Adams @*/ 72d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetAggressiveLevels(PC pc, PetscInt n) 73d71ae5a4SJacob Faibussowitsch { 74ef4ad70eSMark F. Adams PetscFunctionBegin; 75ef4ad70eSMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 76d5d25401SBarry Smith PetscValidLogicalCollectiveInt(pc, n, 2); 77bae903cbSmarkadams4 PetscTryMethod(pc, "PCGAMGSetAggressiveLevels_C", (PC, PetscInt), (pc, n)); 783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 79ef4ad70eSMark F. Adams } 80ef4ad70eSMark F. Adams 81d529f056Smarkadams4 /*@ 82d529f056Smarkadams4 PCGAMGMISkSetAggressive - Number (k) distance in MIS coarsening (>2 is 'aggressive') 83d529f056Smarkadams4 84d529f056Smarkadams4 Logically Collective 85d529f056Smarkadams4 86d529f056Smarkadams4 Input Parameters: 87d529f056Smarkadams4 + pc - the preconditioner context 88d529f056Smarkadams4 - n - 1 or more (default = 2) 89d529f056Smarkadams4 90d529f056Smarkadams4 Options Database Key: 91d529f056Smarkadams4 . -pc_gamg_aggressive_mis_k <n,default=2> - Number (k) distance in MIS coarsening (>2 is 'aggressive') 92d529f056Smarkadams4 93d529f056Smarkadams4 Level: intermediate 94d529f056Smarkadams4 95cbb74892SSatish Balay .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`, `PCGAMGSetAggressiveSquareGraph()`, `PCGAMGMISkSetMinDegreeOrdering()`, `PCGAMGSetLowMemoryFilter()` 96d529f056Smarkadams4 @*/ 97d529f056Smarkadams4 PetscErrorCode PCGAMGMISkSetAggressive(PC pc, PetscInt n) 98d529f056Smarkadams4 { 99d529f056Smarkadams4 PetscFunctionBegin; 100d529f056Smarkadams4 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 101d529f056Smarkadams4 PetscValidLogicalCollectiveInt(pc, n, 2); 102d529f056Smarkadams4 PetscTryMethod(pc, "PCGAMGMISkSetAggressive_C", (PC, PetscInt), (pc, n)); 103d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 104d529f056Smarkadams4 } 105d529f056Smarkadams4 106d529f056Smarkadams4 /*@ 107d529f056Smarkadams4 PCGAMGSetAggressiveSquareGraph - Use graph square A'A for aggressive coarsening, old method 108d529f056Smarkadams4 109d529f056Smarkadams4 Logically Collective 110d529f056Smarkadams4 111d529f056Smarkadams4 Input Parameters: 112d529f056Smarkadams4 + pc - the preconditioner context 113d529f056Smarkadams4 - b - default false - MIS-k is faster 114d529f056Smarkadams4 115d529f056Smarkadams4 Options Database Key: 116d529f056Smarkadams4 . -pc_gamg_aggressive_square_graph <bool,default=false> - Use square graph (A'A) or MIS-k (k=2) for aggressive coarsening 117d529f056Smarkadams4 118d529f056Smarkadams4 Level: intermediate 119d529f056Smarkadams4 120cbb74892SSatish Balay .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`, `PCGAMGMISkSetAggressive()`, `PCGAMGMISkSetMinDegreeOrdering()`, `PCGAMGSetLowMemoryFilter()` 121d529f056Smarkadams4 @*/ 122d529f056Smarkadams4 PetscErrorCode PCGAMGSetAggressiveSquareGraph(PC pc, PetscBool b) 123d529f056Smarkadams4 { 124d529f056Smarkadams4 PetscFunctionBegin; 125d529f056Smarkadams4 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 126d529f056Smarkadams4 PetscValidLogicalCollectiveBool(pc, b, 2); 127d529f056Smarkadams4 PetscTryMethod(pc, "PCGAMGSetAggressiveSquareGraph_C", (PC, PetscBool), (pc, b)); 128d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 129d529f056Smarkadams4 } 130d529f056Smarkadams4 131d529f056Smarkadams4 /*@ 132d529f056Smarkadams4 PCGAMGMISkSetMinDegreeOrdering - Use minimum degree ordering in greedy MIS algorithm 133d529f056Smarkadams4 134d529f056Smarkadams4 Logically Collective 135d529f056Smarkadams4 136d529f056Smarkadams4 Input Parameters: 137d529f056Smarkadams4 + pc - the preconditioner context 138d529f056Smarkadams4 - b - default true 139d529f056Smarkadams4 140d529f056Smarkadams4 Options Database Key: 141d529f056Smarkadams4 . -pc_gamg_mis_k_minimum_degree_ordering <bool,default=true> - Use minimum degree ordering in greedy MIS algorithm 142d529f056Smarkadams4 143d529f056Smarkadams4 Level: intermediate 144d529f056Smarkadams4 145cbb74892SSatish Balay .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`, `PCGAMGMISkSetAggressive()`, `PCGAMGSetAggressiveSquareGraph()`, `PCGAMGSetLowMemoryFilter()` 146d529f056Smarkadams4 @*/ 147d529f056Smarkadams4 PetscErrorCode PCGAMGMISkSetMinDegreeOrdering(PC pc, PetscBool b) 148d529f056Smarkadams4 { 149d529f056Smarkadams4 PetscFunctionBegin; 150d529f056Smarkadams4 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 151d529f056Smarkadams4 PetscValidLogicalCollectiveBool(pc, b, 2); 152d529f056Smarkadams4 PetscTryMethod(pc, "PCGAMGMISkSetMinDegreeOrdering_C", (PC, PetscBool), (pc, b)); 153d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 154d529f056Smarkadams4 } 155d529f056Smarkadams4 156a9ccf005SMark Adams /*@ 157a9ccf005SMark Adams PCGAMGSetLowMemoryFilter - Use low memory graph/matrix filter 158a9ccf005SMark Adams 159a9ccf005SMark Adams Logically Collective 160a9ccf005SMark Adams 161a9ccf005SMark Adams Input Parameters: 162a9ccf005SMark Adams + pc - the preconditioner context 163a9ccf005SMark Adams - b - default false 164a9ccf005SMark Adams 165a9ccf005SMark Adams Options Database Key: 166a9ccf005SMark Adams . -pc_gamg_low_memory_threshold_filter <bool,default=false> - Use low memory graph/matrix filter 167a9ccf005SMark Adams 168a9ccf005SMark Adams Level: intermediate 169a9ccf005SMark Adams 170a9ccf005SMark Adams .seealso: `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`, `PCGAMGMISkSetAggressive()`, `PCGAMGSetAggressiveSquareGraph()`, `PCGAMGMISkSetMinDegreeOrdering()` 171a9ccf005SMark Adams @*/ 172a9ccf005SMark Adams PetscErrorCode PCGAMGSetLowMemoryFilter(PC pc, PetscBool b) 173a9ccf005SMark Adams { 174a9ccf005SMark Adams PetscFunctionBegin; 175a9ccf005SMark Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 176a9ccf005SMark Adams PetscValidLogicalCollectiveBool(pc, b, 2); 177a9ccf005SMark Adams PetscTryMethod(pc, "PCGAMGSetLowMemoryFilter_C", (PC, PetscBool), (pc, b)); 178a9ccf005SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 179a9ccf005SMark Adams } 180a9ccf005SMark Adams 181d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetAggressiveLevels_AGG(PC pc, PetscInt n) 182d71ae5a4SJacob Faibussowitsch { 183ef4ad70eSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 184ef4ad70eSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 185ef4ad70eSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 186ef4ad70eSMark F. Adams 187ef4ad70eSMark F. Adams PetscFunctionBegin; 188bae903cbSmarkadams4 pc_gamg_agg->aggressive_coarsening_levels = n; 1893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 190ef4ad70eSMark F. Adams } 191ef4ad70eSMark F. Adams 192d529f056Smarkadams4 static PetscErrorCode PCGAMGMISkSetAggressive_AGG(PC pc, PetscInt n) 193d529f056Smarkadams4 { 194d529f056Smarkadams4 PC_MG *mg = (PC_MG *)pc->data; 195d529f056Smarkadams4 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 196d529f056Smarkadams4 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 197d529f056Smarkadams4 198d529f056Smarkadams4 PetscFunctionBegin; 199d529f056Smarkadams4 pc_gamg_agg->aggressive_mis_k = n; 200d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 201d529f056Smarkadams4 } 202d529f056Smarkadams4 203d529f056Smarkadams4 static PetscErrorCode PCGAMGSetAggressiveSquareGraph_AGG(PC pc, PetscBool b) 204d529f056Smarkadams4 { 205d529f056Smarkadams4 PC_MG *mg = (PC_MG *)pc->data; 206d529f056Smarkadams4 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 207d529f056Smarkadams4 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 208d529f056Smarkadams4 209d529f056Smarkadams4 PetscFunctionBegin; 210d529f056Smarkadams4 pc_gamg_agg->use_aggressive_square_graph = b; 211d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 212d529f056Smarkadams4 } 213d529f056Smarkadams4 214a9ccf005SMark Adams static PetscErrorCode PCGAMGSetLowMemoryFilter_AGG(PC pc, PetscBool b) 215a9ccf005SMark Adams { 216a9ccf005SMark Adams PC_MG *mg = (PC_MG *)pc->data; 217a9ccf005SMark Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 218a9ccf005SMark Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 219a9ccf005SMark Adams 220a9ccf005SMark Adams PetscFunctionBegin; 221a9ccf005SMark Adams pc_gamg_agg->use_low_mem_filter = b; 222a9ccf005SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 223a9ccf005SMark Adams } 224a9ccf005SMark Adams 225d529f056Smarkadams4 static PetscErrorCode PCGAMGMISkSetMinDegreeOrdering_AGG(PC pc, PetscBool b) 226d529f056Smarkadams4 { 227d529f056Smarkadams4 PC_MG *mg = (PC_MG *)pc->data; 228d529f056Smarkadams4 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 229d529f056Smarkadams4 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 230d529f056Smarkadams4 231d529f056Smarkadams4 PetscFunctionBegin; 232d529f056Smarkadams4 pc_gamg_agg->use_minimum_degree_ordering = b; 233d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 234d529f056Smarkadams4 } 235d529f056Smarkadams4 236d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_GAMG_AGG(PC pc, PetscOptionItems *PetscOptionsObject) 237d71ae5a4SJacob Faibussowitsch { 2382e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 2392e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 240c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 241d4adc10fSMark 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; 242d4adc10fSMark Adams PetscInt nsq_graph_old = 0; 2432e68589bSMark F. Adams 2442e68589bSMark F. Adams PetscFunctionBegin; 245d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "GAMG-AGG options"); 2469566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_gamg_agg_nsmooths", "smoothing steps for smoothed aggregation, usually 1", "PCGAMGSetNSmooths", pc_gamg_agg->nsmooths, &pc_gamg_agg->nsmooths, NULL)); 247d4adc10fSMark Adams // aggressive coarsening logic with deprecated -pc_gamg_square_graph 248d4adc10fSMark 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)); 249d4adc10fSMark Adams if (!n_aggressive_flg) 250d4adc10fSMark 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)); 251d4adc10fSMark 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)); 252d4adc10fSMark Adams if (!new_sq_provided && old_sq_provided) { 253d4adc10fSMark Adams pc_gamg_agg->aggressive_coarsening_levels = nsq_graph_old; // could be zero 254d4adc10fSMark Adams pc_gamg_agg->use_aggressive_square_graph = PETSC_TRUE; 255bae903cbSmarkadams4 } 256d4adc10fSMark Adams if (new_sq_provided && old_sq_provided) 257d4adc10fSMark 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)); 258d529f056Smarkadams4 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)); 259a9ccf005SMark 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)); 260d529f056Smarkadams4 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)); 261d0609cedSBarry Smith PetscOptionsHeadEnd(); 2623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2632e68589bSMark F. Adams } 2642e68589bSMark F. Adams 265d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_GAMG_AGG(PC pc) 266d71ae5a4SJacob Faibussowitsch { 2672e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 2682e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 2692e68589bSMark F. Adams 2702e68589bSMark F. Adams PetscFunctionBegin; 2719566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->subctx)); 2722e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", NULL)); 273bae903cbSmarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", NULL)); 274d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetAggressive_C", NULL)); 275d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetMinDegreeOrdering_C", NULL)); 276a9ccf005SMark Adams PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetLowMemoryFilter_C", NULL)); 277d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveSquareGraph_C", NULL)); 278bae903cbSmarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL)); 2793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2802e68589bSMark F. Adams } 2812e68589bSMark F. Adams 2822e68589bSMark F. Adams /* 2832e68589bSMark F. Adams PCSetCoordinates_AGG 28420f4b53cSBarry Smith 28520f4b53cSBarry Smith Collective 2862e68589bSMark F. Adams 2872e68589bSMark F. Adams Input Parameter: 2882e68589bSMark F. Adams . pc - the preconditioner context 289145b44c9SPierre Jolivet . ndm - dimension of data (used for dof/vertex for Stokes) 290302f38e8SMark F. Adams . a_nloc - number of vertices local 291302f38e8SMark F. Adams . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...} 2922e68589bSMark F. Adams */ 2931e6b0712SBarry Smith 294d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords) 295d71ae5a4SJacob Faibussowitsch { 2962e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 2972e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 29869344418SMark F. Adams PetscInt arrsz, kk, ii, jj, nloc, ndatarows, ndf; 299a2f3521dSMark F. Adams Mat mat = pc->pmat; 3002e68589bSMark F. Adams 3012e68589bSMark F. Adams PetscFunctionBegin; 302a2f3521dSMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 303a2f3521dSMark F. Adams PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 304302f38e8SMark F. Adams nloc = a_nloc; 3052e68589bSMark F. Adams 3062e68589bSMark F. Adams /* SA: null space vectors */ 3079566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(mat, &ndf)); /* this does not work for Stokes */ 30869344418SMark F. Adams if (coords && ndf == 1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */ 309a2f3521dSMark F. Adams else if (coords) { 31063a3b9bcSJacob Faibussowitsch PetscCheck(ndm <= ndf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "degrees of motion %" PetscInt_FMT " > block size %" PetscInt_FMT, ndm, ndf); 31169344418SMark F. Adams pc_gamg->data_cell_cols = (ndm == 2 ? 3 : 6); /* displacement elasticity */ 312ad540459SPierre 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); 313806fa848SBarry Smith } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */ 31471959b99SBarry Smith pc_gamg->data_cell_rows = ndatarows = ndf; 31563a3b9bcSJacob 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); 316c8b0795cSMark F. Adams arrsz = nloc * pc_gamg->data_cell_rows * pc_gamg->data_cell_cols; 3172e68589bSMark F. Adams 3187f66b68fSMark Adams if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) { 3199566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 3209566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(arrsz + 1, &pc_gamg->data)); 3212e68589bSMark F. Adams } 3222e68589bSMark F. Adams /* copy data in - column oriented */ 3232e68589bSMark F. Adams for (kk = 0; kk < nloc; kk++) { 32469344418SMark F. Adams const PetscInt M = nloc * pc_gamg->data_cell_rows; /* stride into data */ 32569344418SMark F. Adams PetscReal *data = &pc_gamg->data[kk * ndatarows]; /* start of cell */ 326c8b0795cSMark F. Adams if (pc_gamg->data_cell_cols == 1) *data = 1.0; 3272e68589bSMark F. Adams else { 32869344418SMark F. Adams /* translational modes */ 3292fa5cd67SKarl Rupp for (ii = 0; ii < ndatarows; ii++) { 3302fa5cd67SKarl Rupp for (jj = 0; jj < ndatarows; jj++) { 33169344418SMark F. Adams if (ii == jj) data[ii * M + jj] = 1.0; 3322e68589bSMark F. Adams else data[ii * M + jj] = 0.0; 3332fa5cd67SKarl Rupp } 3342fa5cd67SKarl Rupp } 33569344418SMark F. Adams 33669344418SMark F. Adams /* rotational modes */ 3372e68589bSMark F. Adams if (coords) { 33869344418SMark F. Adams if (ndm == 2) { 3392e68589bSMark F. Adams data += 2 * M; 3402e68589bSMark F. Adams data[0] = -coords[2 * kk + 1]; 3412e68589bSMark F. Adams data[1] = coords[2 * kk]; 342806fa848SBarry Smith } else { 3432e68589bSMark F. Adams data += 3 * M; 3449371c9d4SSatish Balay data[0] = 0.0; 3459371c9d4SSatish Balay data[M + 0] = coords[3 * kk + 2]; 3469371c9d4SSatish Balay data[2 * M + 0] = -coords[3 * kk + 1]; 3479371c9d4SSatish Balay data[1] = -coords[3 * kk + 2]; 3489371c9d4SSatish Balay data[M + 1] = 0.0; 3499371c9d4SSatish Balay data[2 * M + 1] = coords[3 * kk]; 3509371c9d4SSatish Balay data[2] = coords[3 * kk + 1]; 3519371c9d4SSatish Balay data[M + 2] = -coords[3 * kk]; 3529371c9d4SSatish Balay data[2 * M + 2] = 0.0; 3532e68589bSMark F. Adams } 3542e68589bSMark F. Adams } 3552e68589bSMark F. Adams } 3562e68589bSMark F. Adams } 3572e68589bSMark F. Adams pc_gamg->data_sz = arrsz; 3583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3592e68589bSMark F. Adams } 3602e68589bSMark F. Adams 3612e68589bSMark F. Adams /* 362a2f3521dSMark F. Adams PCSetData_AGG - called if data is not set with PCSetCoordinates. 363a2f3521dSMark F. Adams Looks in Mat for near null space. 364a2f3521dSMark F. Adams Does not work for Stokes 3652e68589bSMark F. Adams 3662e68589bSMark F. Adams Input Parameter: 3672e68589bSMark F. Adams . pc - 368a2f3521dSMark F. Adams . a_A - matrix to get (near) null space out of. 3692e68589bSMark F. Adams */ 370d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A) 371d71ae5a4SJacob Faibussowitsch { 372b8cd405aSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 373b8cd405aSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 374b8cd405aSMark F. Adams MatNullSpace mnull; 37566f2ef4bSMark Adams 376b817416eSBarry Smith PetscFunctionBegin; 3779566063dSJacob Faibussowitsch PetscCall(MatGetNearNullSpace(a_A, &mnull)); 378b8cd405aSMark F. Adams if (!mnull) { 37966f2ef4bSMark Adams DM dm; 3809566063dSJacob Faibussowitsch PetscCall(PCGetDM(pc, &dm)); 38148a46eb9SPierre Jolivet if (!dm) PetscCall(MatGetDM(a_A, &dm)); 38266f2ef4bSMark Adams if (dm) { 38366f2ef4bSMark Adams PetscObject deformation; 384b0d30dd6SMatthew G. Knepley PetscInt Nf; 385b0d30dd6SMatthew G. Knepley 3869566063dSJacob Faibussowitsch PetscCall(DMGetNumFields(dm, &Nf)); 387b0d30dd6SMatthew G. Knepley if (Nf) { 3889566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, 0, NULL, &deformation)); 3899566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)deformation, "nearnullspace", (PetscObject *)&mnull)); 39048a46eb9SPierre Jolivet if (!mnull) PetscCall(PetscObjectQuery((PetscObject)deformation, "nullspace", (PetscObject *)&mnull)); 39166f2ef4bSMark Adams } 39266f2ef4bSMark Adams } 393b0d30dd6SMatthew G. Knepley } 39466f2ef4bSMark Adams 39566f2ef4bSMark Adams if (!mnull) { 396a2f3521dSMark F. Adams PetscInt bs, NN, MM; 3979566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(a_A, &bs)); 3989566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(a_A, &MM, &NN)); 39963a3b9bcSJacob Faibussowitsch PetscCheck(MM % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MM %" PetscInt_FMT " must be divisible by bs %" PetscInt_FMT, MM, bs); 4009566063dSJacob Faibussowitsch PetscCall(PCSetCoordinates_AGG(pc, bs, MM / bs, NULL)); 401806fa848SBarry Smith } else { 402b8cd405aSMark F. Adams PetscReal *nullvec; 403b8cd405aSMark F. Adams PetscBool has_const; 404b8cd405aSMark F. Adams PetscInt i, j, mlocal, nvec, bs; 405d5d25401SBarry Smith const Vec *vecs; 406d5d25401SBarry Smith const PetscScalar *v; 407b817416eSBarry Smith 4089566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(a_A, &mlocal, NULL)); 4099566063dSJacob Faibussowitsch PetscCall(MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs)); 410d19c4ebbSmarkadams4 for (i = 0; i < nvec; i++) { 411d19c4ebbSmarkadams4 PetscCall(VecGetLocalSize(vecs[i], &j)); 412d19c4ebbSmarkadams4 PetscCheck(j == mlocal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Attached null space vector size %" PetscInt_FMT " != matrix size %" PetscInt_FMT, j, mlocal); 413d19c4ebbSmarkadams4 } 414a0dea87dSMark Adams pc_gamg->data_sz = (nvec + !!has_const) * mlocal; 4159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((nvec + !!has_const) * mlocal, &nullvec)); 4169371c9d4SSatish Balay if (has_const) 4179371c9d4SSatish Balay for (i = 0; i < mlocal; i++) nullvec[i] = 1.0; 418b8cd405aSMark F. Adams for (i = 0; i < nvec; i++) { 4199566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(vecs[i], &v)); 420b8cd405aSMark F. Adams for (j = 0; j < mlocal; j++) nullvec[(i + !!has_const) * mlocal + j] = PetscRealPart(v[j]); 4219566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(vecs[i], &v)); 422b8cd405aSMark F. Adams } 423b8cd405aSMark F. Adams pc_gamg->data = nullvec; 424b8cd405aSMark F. Adams pc_gamg->data_cell_cols = (nvec + !!has_const); 4259566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(a_A, &bs)); 426b8cd405aSMark F. Adams pc_gamg->data_cell_rows = bs; 427b8cd405aSMark F. Adams } 4283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4292e68589bSMark F. Adams } 4302e68589bSMark F. Adams 4312e68589bSMark F. Adams /* 432bae903cbSmarkadams4 formProl0 - collect null space data for each aggregate, do QR, put R in coarse grid data and Q in P_0 4332e68589bSMark F. Adams 4342e68589bSMark F. Adams Input Parameter: 4352adfe9d3SBarry Smith . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices 4362adfe9d3SBarry Smith . bs - row block size 4370cbbd2e1SMark F. Adams . nSAvec - column bs of new P 4380cbbd2e1SMark F. Adams . my0crs - global index of start of locals 4392adfe9d3SBarry Smith . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec] 4402e68589bSMark F. Adams . data_in[data_stride*nSAvec] - local data on fine grid 4412e68589bSMark F. Adams . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist' 442bae903cbSmarkadams4 4432e68589bSMark F. Adams Output Parameter: 4442e68589bSMark F. Adams . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data 4452e68589bSMark F. Adams . a_Prol - prolongation operator 4462e68589bSMark F. Adams */ 447d71ae5a4SJacob 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) 448d71ae5a4SJacob Faibussowitsch { 449bd026e97SJed Brown PetscInt Istart, my0, Iend, nloc, clid, flid = 0, aggID, kk, jj, ii, mm, nSelected, minsz, nghosts, out_data_stride; 4503b4367a7SBarry Smith MPI_Comm comm; 4512e68589bSMark F. Adams PetscReal *out_data; 452539c167fSBarry Smith PetscCDIntNd *pos; 4531943db53SBarry Smith PCGAMGHashTable fgid_flid; 4540cbbd2e1SMark F. Adams 4552e68589bSMark F. Adams PetscFunctionBegin; 4569566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)a_Prol, &comm)); 4579566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(a_Prol, &Istart, &Iend)); 4589371c9d4SSatish Balay nloc = (Iend - Istart) / bs; 4599371c9d4SSatish Balay my0 = Istart / bs; 46063a3b9bcSJacob 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); 4610cbbd2e1SMark F. Adams Iend /= bs; 4620cbbd2e1SMark F. Adams nghosts = data_stride / bs - nloc; 4632e68589bSMark F. Adams 4649566063dSJacob Faibussowitsch PetscCall(PCGAMGHashTableCreate(2 * nghosts + 1, &fgid_flid)); 46548a46eb9SPierre Jolivet for (kk = 0; kk < nghosts; kk++) PetscCall(PCGAMGHashTableAdd(&fgid_flid, flid_fgid[nloc + kk], nloc + kk)); 4662e68589bSMark F. Adams 4670cbbd2e1SMark F. Adams /* count selected -- same as number of cols of P */ 4680cbbd2e1SMark F. Adams for (nSelected = mm = 0; mm < nloc; mm++) { 469e78576d6SMark F. Adams PetscBool ise; 4705e62d202SMark Adams PetscCall(PetscCDIsEmptyAt(agg_llists, mm, &ise)); 471e78576d6SMark F. Adams if (!ise) nSelected++; 4720cbbd2e1SMark F. Adams } 4739566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(a_Prol, &ii, &jj)); 47463a3b9bcSJacob Faibussowitsch PetscCheck((ii / nSAvec) == my0crs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ii %" PetscInt_FMT " /nSAvec %" PetscInt_FMT " != my0crs %" PetscInt_FMT, ii, nSAvec, my0crs); 47563a3b9bcSJacob 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); 4760cbbd2e1SMark F. Adams 4772e68589bSMark F. Adams /* aloc space for coarse point data (output) */ 4780cbbd2e1SMark F. Adams out_data_stride = nSelected * nSAvec; 4792fa5cd67SKarl Rupp 4809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(out_data_stride * nSAvec, &out_data)); 48133812677SSatish Balay for (ii = 0; ii < out_data_stride * nSAvec; ii++) out_data[ii] = PETSC_MAX_REAL; 4820cbbd2e1SMark F. Adams *a_data_out = out_data; /* output - stride nSelected*nSAvec */ 4832e68589bSMark F. Adams 4842e68589bSMark F. Adams /* find points and set prolongation */ 485c8b0795cSMark F. Adams minsz = 100; 4860cbbd2e1SMark F. Adams for (mm = clid = 0; mm < nloc; mm++) { 4875e62d202SMark Adams PetscCall(PetscCDCountAt(agg_llists, mm, &jj)); 488e78576d6SMark F. Adams if (jj > 0) { 4890cbbd2e1SMark F. Adams const PetscInt lid = mm, cgid = my0crs + clid; 4900cbbd2e1SMark F. Adams PetscInt cids[100]; /* max bs */ 4910cbbd2e1SMark F. Adams PetscBLASInt asz = jj, M = asz * bs, N = nSAvec, INFO; 4922e68589bSMark F. Adams PetscBLASInt Mdata = M + ((N - M > 0) ? N - M : 0), LDA = Mdata, LWORK = N * bs; 4932e68589bSMark F. Adams PetscScalar *qqc, *qqr, *TAU, *WORK; 4942e68589bSMark F. Adams PetscInt *fids; 49565d7b583SSatish Balay PetscReal *data; 496b817416eSBarry Smith 4970cbbd2e1SMark F. Adams /* count agg */ 4980cbbd2e1SMark F. Adams if (asz < minsz) minsz = asz; 4990cbbd2e1SMark F. Adams 5000cbbd2e1SMark F. Adams /* get block */ 5019566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(Mdata * N, &qqc, M * N, &qqr, N, &TAU, LWORK, &WORK, M, &fids)); 5022e68589bSMark F. Adams 5032e68589bSMark F. Adams aggID = 0; 5049566063dSJacob Faibussowitsch PetscCall(PetscCDGetHeadPos(agg_llists, lid, &pos)); 505e78576d6SMark F. Adams while (pos) { 506e78576d6SMark F. Adams PetscInt gid1; 5079566063dSJacob Faibussowitsch PetscCall(PetscCDIntNdGetID(pos, &gid1)); 5089566063dSJacob Faibussowitsch PetscCall(PetscCDGetNextPos(agg_llists, lid, &pos)); 509e78576d6SMark F. Adams 5100cbbd2e1SMark F. Adams if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0; 5110cbbd2e1SMark F. Adams else { 5129566063dSJacob Faibussowitsch PetscCall(PCGAMGHashTableFind(&fgid_flid, gid1, &flid)); 51308401ef6SPierre Jolivet PetscCheck(flid >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot find gid1 in table"); 5140cbbd2e1SMark F. Adams } 5152e68589bSMark F. Adams /* copy in B_i matrix - column oriented */ 51665d7b583SSatish Balay data = &data_in[flid * bs]; 5173d3eaba7SBarry Smith for (ii = 0; ii < bs; ii++) { 5182e68589bSMark F. Adams for (jj = 0; jj < N; jj++) { 5190cbbd2e1SMark F. Adams PetscReal d = data[jj * data_stride + ii]; 5200cbbd2e1SMark F. Adams qqc[jj * Mdata + aggID * bs + ii] = d; 5212e68589bSMark F. Adams } 5222e68589bSMark F. Adams } 5232e68589bSMark F. Adams /* set fine IDs */ 5242e68589bSMark F. Adams for (kk = 0; kk < bs; kk++) fids[aggID * bs + kk] = flid_fgid[flid] * bs + kk; 5252e68589bSMark F. Adams aggID++; 5260cbbd2e1SMark F. Adams } 5272e68589bSMark F. Adams 5282e68589bSMark F. Adams /* pad with zeros */ 5292e68589bSMark F. Adams for (ii = asz * bs; ii < Mdata; ii++) { 530ad540459SPierre Jolivet for (jj = 0; jj < N; jj++, kk++) qqc[jj * Mdata + ii] = .0; 5312e68589bSMark F. Adams } 5322e68589bSMark F. Adams 5332e68589bSMark F. Adams /* QR */ 5349566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 535792fecdfSBarry Smith PetscCallBLAS("LAPACKgeqrf", LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 5369566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPop()); 53708401ef6SPierre Jolivet PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xGEQRF error"); 5382e68589bSMark F. Adams /* get R - column oriented - output B_{i+1} */ 5392e68589bSMark F. Adams { 5402e68589bSMark F. Adams PetscReal *data = &out_data[clid * nSAvec]; 5412e68589bSMark F. Adams for (jj = 0; jj < nSAvec; jj++) { 5422e68589bSMark F. Adams for (ii = 0; ii < nSAvec; ii++) { 54308401ef6SPierre 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); 5440cbbd2e1SMark F. Adams if (ii <= jj) data[jj * out_data_stride + ii] = PetscRealPart(qqc[jj * Mdata + ii]); 5450cbbd2e1SMark F. Adams else data[jj * out_data_stride + ii] = 0.; 5462e68589bSMark F. Adams } 5472e68589bSMark F. Adams } 5482e68589bSMark F. Adams } 5492e68589bSMark F. Adams 5502e68589bSMark F. Adams /* get Q - row oriented */ 551792fecdfSBarry Smith PetscCallBLAS("LAPACKorgqr", LAPACKorgqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 55263a3b9bcSJacob Faibussowitsch PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xORGQR error arg %" PetscBLASInt_FMT, -INFO); 5532e68589bSMark F. Adams 5542e68589bSMark F. Adams for (ii = 0; ii < M; ii++) { 555ad540459SPierre Jolivet for (jj = 0; jj < N; jj++) qqr[N * ii + jj] = qqc[jj * Mdata + ii]; 5562e68589bSMark F. Adams } 5572e68589bSMark F. Adams 5582e68589bSMark F. Adams /* add diagonal block of P0 */ 5599371c9d4SSatish Balay for (kk = 0; kk < N; kk++) { cids[kk] = N * cgid + kk; /* global col IDs in P0 */ } 5609566063dSJacob Faibussowitsch PetscCall(MatSetValues(a_Prol, M, fids, N, cids, qqr, INSERT_VALUES)); 5619566063dSJacob Faibussowitsch PetscCall(PetscFree5(qqc, qqr, TAU, WORK, fids)); 562b43b03e9SMark F. Adams clid++; 5630cbbd2e1SMark F. Adams } /* coarse agg */ 5640cbbd2e1SMark F. Adams } /* for all fine nodes */ 5659566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(a_Prol, MAT_FINAL_ASSEMBLY)); 5669566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(a_Prol, MAT_FINAL_ASSEMBLY)); 5679566063dSJacob Faibussowitsch PetscCall(PCGAMGHashTableDestroy(&fgid_flid)); 5683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5692e68589bSMark F. Adams } 5702e68589bSMark F. Adams 571d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_GAMG_AGG(PC pc, PetscViewer viewer) 572d71ae5a4SJacob Faibussowitsch { 5735adeb434SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 5745adeb434SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 5755adeb434SBarry Smith PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 5765adeb434SBarry Smith 5775adeb434SBarry Smith PetscFunctionBegin; 5789566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " AGG specific options\n")); 579d529f056Smarkadams4 PetscCall(PetscViewerASCIIPrintf(viewer, " Number of levels of aggressive coarsening %d\n", (int)pc_gamg_agg->aggressive_coarsening_levels)); 580d529f056Smarkadams4 if (pc_gamg_agg->aggressive_coarsening_levels > 0) { 581d529f056Smarkadams4 PetscCall(PetscViewerASCIIPrintf(viewer, " %s aggressive coarsening\n", !pc_gamg_agg->use_aggressive_square_graph ? "MIS-k" : "Square graph")); 582d529f056Smarkadams4 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)); 583d529f056Smarkadams4 } 584bae903cbSmarkadams4 PetscCall(PetscViewerASCIIPrintf(viewer, " Number smoothing steps %d\n", (int)pc_gamg_agg->nsmooths)); 5853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5865adeb434SBarry Smith } 5875adeb434SBarry Smith 5882d776b49SBarry Smith static PetscErrorCode PCGAMGCreateGraph_AGG(PC pc, Mat Amat, Mat *a_Gmat) 589d71ae5a4SJacob Faibussowitsch { 590c8b0795cSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 591c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 592c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 5932d776b49SBarry Smith const PetscReal vfilter = pc_gamg->threshold[pc_gamg->current_level]; 5949c15c1aeSMark Adams PetscBool ishem, ismis; 5952d776b49SBarry Smith const char *prefix; 596d529f056Smarkadams4 MatInfo info0, info1; 597d529f056Smarkadams4 PetscInt bs; 598c8b0795cSMark F. Adams 599c8b0795cSMark F. Adams PetscFunctionBegin; 600a9ccf005SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0)); 6012d776b49SBarry 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 */ 6022d776b49SBarry Smith /* MATCOARSENHEM requires numerical weights for edges so ensure they are computed */ 6032d776b49SBarry Smith PetscCall(MatCoarsenCreate(PetscObjectComm((PetscObject)pc), &pc_gamg_agg->crs)); 6042d776b49SBarry Smith PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, &prefix)); 6052d776b49SBarry Smith PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc_gamg_agg->crs, prefix)); 6062d776b49SBarry Smith PetscCall(MatCoarsenSetFromOptions(pc_gamg_agg->crs)); 607e02fb3cdSMark Adams PetscCall(MatGetBlockSize(Amat, &bs)); 608e02fb3cdSMark Adams // check for valid indices wrt bs 609e02fb3cdSMark Adams for (int ii = 0; ii < pc_gamg_agg->crs->strength_index_size; ii++) { 610e02fb3cdSMark 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)", (int)pc_gamg_agg->crs->strength_index[ii], (int)bs); 611e02fb3cdSMark Adams } 6122d776b49SBarry Smith PetscCall(PetscObjectTypeCompare((PetscObject)pc_gamg_agg->crs, MATCOARSENHEM, &ishem)); 6135e62d202SMark Adams if (ishem) { 6149c15c1aeSMark 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)); 6155e62d202SMark Adams pc_gamg_agg->aggressive_coarsening_levels = 0; // aggressive and HEM does not make sense 6165e62d202SMark Adams PetscCall(MatCoarsenSetMaximumIterations(pc_gamg_agg->crs, pc_gamg_agg->crs->max_it)); // for code coverage 6175e62d202SMark Adams PetscCall(MatCoarsenSetThreshold(pc_gamg_agg->crs, vfilter)); // for code coverage 6189c15c1aeSMark Adams } else { 6199c15c1aeSMark Adams PetscCall(PetscObjectTypeCompare((PetscObject)pc_gamg_agg->crs, MATCOARSENMIS, &ismis)); 6209c15c1aeSMark Adams if (ismis && pc_gamg_agg->aggressive_coarsening_levels && !pc_gamg_agg->use_aggressive_square_graph) { 6219c15c1aeSMark Adams PetscCall(PetscInfo(pc, "MIS and aggressive coarsening and no square graph: force square graph\n")); 6229c15c1aeSMark Adams pc_gamg_agg->use_aggressive_square_graph = PETSC_TRUE; 6239c15c1aeSMark Adams } 6245e62d202SMark Adams } 625a9ccf005SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0)); 626a9ccf005SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0)); 627d529f056Smarkadams4 PetscCall(MatGetInfo(Amat, MAT_LOCAL, &info0)); /* global reduction */ 628a9ccf005SMark Adams 629a9ccf005SMark Adams if (ishem || pc_gamg_agg->use_low_mem_filter) { 630e02fb3cdSMark 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)); 631a9ccf005SMark Adams } else { 632a9ccf005SMark Adams // make scalar graph, symetrize if not know to be symetric, scale, but do not filter (expensive) 633e02fb3cdSMark Adams PetscCall(MatCreateGraph(Amat, PETSC_TRUE, PETSC_TRUE, -1, pc_gamg_agg->crs->strength_index_size, pc_gamg_agg->crs->strength_index, a_Gmat)); 634a9ccf005SMark Adams if (vfilter >= 0) { 635a9ccf005SMark Adams PetscInt Istart, Iend, ncols, nnz0, nnz1, NN, MM, nloc; 636a9ccf005SMark Adams Mat tGmat, Gmat = *a_Gmat; 637a9ccf005SMark Adams MPI_Comm comm; 638a9ccf005SMark Adams const PetscScalar *vals; 639a9ccf005SMark Adams const PetscInt *idx; 640a9ccf005SMark Adams PetscInt *d_nnz, *o_nnz, kk, *garray = NULL, *AJ, maxcols = 0; 641a9ccf005SMark Adams MatScalar *AA; // this is checked in graph 642a9ccf005SMark Adams PetscBool isseqaij; 643a9ccf005SMark Adams Mat a, b, c; 644a9ccf005SMark Adams MatType jtype; 645a9ccf005SMark Adams 646a9ccf005SMark Adams PetscCall(PetscObjectGetComm((PetscObject)Gmat, &comm)); 647a9ccf005SMark Adams PetscCall(PetscObjectBaseTypeCompare((PetscObject)Gmat, MATSEQAIJ, &isseqaij)); 648a9ccf005SMark Adams PetscCall(MatGetType(Gmat, &jtype)); 649a9ccf005SMark Adams PetscCall(MatCreate(comm, &tGmat)); 650a9ccf005SMark Adams PetscCall(MatSetType(tGmat, jtype)); 651a9ccf005SMark Adams 652a9ccf005SMark Adams /* TODO GPU: this can be called when filter = 0 -> Probably provide MatAIJThresholdCompress that compresses the entries below a threshold? 653a9ccf005SMark Adams Also, if the matrix is symmetric, can we skip this 654a9ccf005SMark Adams operation? It can be very expensive on large matrices. */ 655a9ccf005SMark Adams 656a9ccf005SMark Adams // global sizes 657a9ccf005SMark Adams PetscCall(MatGetSize(Gmat, &MM, &NN)); 658a9ccf005SMark Adams PetscCall(MatGetOwnershipRange(Gmat, &Istart, &Iend)); 659a9ccf005SMark Adams nloc = Iend - Istart; 660a9ccf005SMark Adams PetscCall(PetscMalloc2(nloc, &d_nnz, nloc, &o_nnz)); 661a9ccf005SMark Adams if (isseqaij) { 662a9ccf005SMark Adams a = Gmat; 663a9ccf005SMark Adams b = NULL; 664a9ccf005SMark Adams } else { 665a9ccf005SMark Adams Mat_MPIAIJ *d = (Mat_MPIAIJ *)Gmat->data; 666a9ccf005SMark Adams a = d->A; 667a9ccf005SMark Adams b = d->B; 668a9ccf005SMark Adams garray = d->garray; 669a9ccf005SMark Adams } 670a9ccf005SMark Adams /* Determine upper bound on non-zeros needed in new filtered matrix */ 671a9ccf005SMark Adams for (PetscInt row = 0; row < nloc; row++) { 672a9ccf005SMark Adams PetscCall(MatGetRow(a, row, &ncols, NULL, NULL)); 673a9ccf005SMark Adams d_nnz[row] = ncols; 674a9ccf005SMark Adams if (ncols > maxcols) maxcols = ncols; 675a9ccf005SMark Adams PetscCall(MatRestoreRow(a, row, &ncols, NULL, NULL)); 676a9ccf005SMark Adams } 677a9ccf005SMark Adams if (b) { 678a9ccf005SMark Adams for (PetscInt row = 0; row < nloc; row++) { 679a9ccf005SMark Adams PetscCall(MatGetRow(b, row, &ncols, NULL, NULL)); 680a9ccf005SMark Adams o_nnz[row] = ncols; 681a9ccf005SMark Adams if (ncols > maxcols) maxcols = ncols; 682a9ccf005SMark Adams PetscCall(MatRestoreRow(b, row, &ncols, NULL, NULL)); 683a9ccf005SMark Adams } 684a9ccf005SMark Adams } 685a9ccf005SMark Adams PetscCall(MatSetSizes(tGmat, nloc, nloc, MM, MM)); 686a9ccf005SMark Adams PetscCall(MatSetBlockSizes(tGmat, 1, 1)); 687a9ccf005SMark Adams PetscCall(MatSeqAIJSetPreallocation(tGmat, 0, d_nnz)); 688a9ccf005SMark Adams PetscCall(MatMPIAIJSetPreallocation(tGmat, 0, d_nnz, 0, o_nnz)); 689a9ccf005SMark Adams PetscCall(MatSetOption(tGmat, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); 690a9ccf005SMark Adams PetscCall(PetscFree2(d_nnz, o_nnz)); 691a9ccf005SMark Adams PetscCall(PetscMalloc2(maxcols, &AA, maxcols, &AJ)); 692a9ccf005SMark Adams nnz0 = nnz1 = 0; 693a9ccf005SMark Adams for (c = a, kk = 0; c && kk < 2; c = b, kk++) { 694a9ccf005SMark Adams for (PetscInt row = 0, grow = Istart, ncol_row, jj; row < nloc; row++, grow++) { 695a9ccf005SMark Adams PetscCall(MatGetRow(c, row, &ncols, &idx, &vals)); 696a9ccf005SMark Adams for (ncol_row = jj = 0; jj < ncols; jj++, nnz0++) { 697a9ccf005SMark Adams PetscScalar sv = PetscAbs(PetscRealPart(vals[jj])); 698a9ccf005SMark Adams if (PetscRealPart(sv) > vfilter) { 699a9ccf005SMark Adams PetscInt cid = idx[jj] + Istart; //diag 700a9ccf005SMark Adams nnz1++; 701a9ccf005SMark Adams if (c != a) cid = garray[idx[jj]]; 702a9ccf005SMark Adams AA[ncol_row] = vals[jj]; 703a9ccf005SMark Adams AJ[ncol_row] = cid; 704a9ccf005SMark Adams ncol_row++; 705a9ccf005SMark Adams } 706a9ccf005SMark Adams } 707a9ccf005SMark Adams PetscCall(MatRestoreRow(c, row, &ncols, &idx, &vals)); 708a9ccf005SMark Adams PetscCall(MatSetValues(tGmat, 1, &grow, ncol_row, AJ, AA, INSERT_VALUES)); 709a9ccf005SMark Adams } 710a9ccf005SMark Adams } 711a9ccf005SMark Adams PetscCall(PetscFree2(AA, AJ)); 712a9ccf005SMark Adams PetscCall(MatAssemblyBegin(tGmat, MAT_FINAL_ASSEMBLY)); 713a9ccf005SMark Adams PetscCall(MatAssemblyEnd(tGmat, MAT_FINAL_ASSEMBLY)); 714a9ccf005SMark Adams PetscCall(MatPropagateSymmetryOptions(Gmat, tGmat)); /* Normal Mat options are not relevant ? */ 715a9ccf005SMark 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)); 716a9ccf005SMark Adams PetscCall(MatViewFromOptions(tGmat, NULL, "-mat_filter_graph_view")); 717a9ccf005SMark Adams PetscCall(MatDestroy(&Gmat)); 718a9ccf005SMark Adams *a_Gmat = tGmat; 719a9ccf005SMark Adams } 720a9ccf005SMark Adams } 721a9ccf005SMark Adams 722d529f056Smarkadams4 PetscCall(MatGetInfo(*a_Gmat, MAT_LOCAL, &info1)); /* global reduction */ 723d529f056Smarkadams4 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)); 724849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0)); 7253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 726c8b0795cSMark F. Adams } 727c8b0795cSMark F. Adams 728d529f056Smarkadams4 typedef PetscInt NState; 729d529f056Smarkadams4 static const NState NOT_DONE = -2; 730d529f056Smarkadams4 static const NState DELETED = -1; 731d529f056Smarkadams4 static const NState REMOVED = -3; 732d529f056Smarkadams4 #define IS_SELECTED(s) (s != DELETED && s != NOT_DONE && s != REMOVED) 733d529f056Smarkadams4 734d529f056Smarkadams4 /* 735d529f056Smarkadams4 fixAggregatesWithSquare - greedy grab of with G1 (unsquared graph) -- AIJ specific -- change to fixAggregatesWithSquare -- TODD 736d529f056Smarkadams4 - AGG-MG specific: clears singletons out of 'selected_2' 737d529f056Smarkadams4 738d529f056Smarkadams4 Input Parameter: 739d529f056Smarkadams4 . Gmat_2 - global matrix of squared graph (data not defined) 740d529f056Smarkadams4 . Gmat_1 - base graph to grab with base graph 741d529f056Smarkadams4 Input/Output Parameter: 742d529f056Smarkadams4 . aggs_2 - linked list of aggs with gids) 743d529f056Smarkadams4 */ 744d529f056Smarkadams4 static PetscErrorCode fixAggregatesWithSquare(PC pc, Mat Gmat_2, Mat Gmat_1, PetscCoarsenData *aggs_2) 745d529f056Smarkadams4 { 746d529f056Smarkadams4 PetscBool isMPI; 747d529f056Smarkadams4 Mat_SeqAIJ *matA_1, *matB_1 = NULL; 748d529f056Smarkadams4 MPI_Comm comm; 749d529f056Smarkadams4 PetscInt lid, *ii, *idx, ix, Iend, my0, kk, n, j; 750d529f056Smarkadams4 Mat_MPIAIJ *mpimat_2 = NULL, *mpimat_1 = NULL; 751d529f056Smarkadams4 const PetscInt nloc = Gmat_2->rmap->n; 752d529f056Smarkadams4 PetscScalar *cpcol_1_state, *cpcol_2_state, *cpcol_2_par_orig, *lid_parent_gid; 753d529f056Smarkadams4 PetscInt *lid_cprowID_1 = NULL; 754d529f056Smarkadams4 NState *lid_state; 755d529f056Smarkadams4 Vec ghost_par_orig2; 756d529f056Smarkadams4 PetscMPIInt rank; 757d529f056Smarkadams4 758d529f056Smarkadams4 PetscFunctionBegin; 759d529f056Smarkadams4 PetscCall(PetscObjectGetComm((PetscObject)Gmat_2, &comm)); 760d529f056Smarkadams4 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 761d529f056Smarkadams4 PetscCall(MatGetOwnershipRange(Gmat_1, &my0, &Iend)); 762d529f056Smarkadams4 763d529f056Smarkadams4 /* get submatrices */ 764d529f056Smarkadams4 PetscCall(PetscStrbeginswith(((PetscObject)Gmat_1)->type_name, MATMPIAIJ, &isMPI)); 765d529f056Smarkadams4 PetscCall(PetscInfo(pc, "isMPI = %s\n", isMPI ? "yes" : "no")); 766d529f056Smarkadams4 PetscCall(PetscMalloc3(nloc, &lid_state, nloc, &lid_parent_gid, nloc, &lid_cprowID_1)); 767d529f056Smarkadams4 for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1; 768d529f056Smarkadams4 if (isMPI) { 769d529f056Smarkadams4 /* grab matrix objects */ 770d529f056Smarkadams4 mpimat_2 = (Mat_MPIAIJ *)Gmat_2->data; 771d529f056Smarkadams4 mpimat_1 = (Mat_MPIAIJ *)Gmat_1->data; 772d529f056Smarkadams4 matA_1 = (Mat_SeqAIJ *)mpimat_1->A->data; 773d529f056Smarkadams4 matB_1 = (Mat_SeqAIJ *)mpimat_1->B->data; 774d529f056Smarkadams4 775d529f056Smarkadams4 /* force compressed row storage for B matrix in AuxMat */ 776d529f056Smarkadams4 PetscCall(MatCheckCompressedRow(mpimat_1->B, matB_1->nonzerorowcnt, &matB_1->compressedrow, matB_1->i, Gmat_1->rmap->n, -1.0)); 777d529f056Smarkadams4 for (ix = 0; ix < matB_1->compressedrow.nrows; ix++) { 778d529f056Smarkadams4 PetscInt lid = matB_1->compressedrow.rindex[ix]; 779d529f056Smarkadams4 PetscCheck(lid <= nloc && lid >= -1, PETSC_COMM_SELF, PETSC_ERR_USER, "lid %d out of range. nloc = %d", (int)lid, (int)nloc); 780d529f056Smarkadams4 if (lid != -1) lid_cprowID_1[lid] = ix; 781d529f056Smarkadams4 } 782d529f056Smarkadams4 } else { 783d529f056Smarkadams4 PetscBool isAIJ; 784d529f056Smarkadams4 PetscCall(PetscStrbeginswith(((PetscObject)Gmat_1)->type_name, MATSEQAIJ, &isAIJ)); 785d529f056Smarkadams4 PetscCheck(isAIJ, PETSC_COMM_SELF, PETSC_ERR_USER, "Require AIJ matrix."); 786d529f056Smarkadams4 matA_1 = (Mat_SeqAIJ *)Gmat_1->data; 787d529f056Smarkadams4 } 788d529f056Smarkadams4 if (nloc > 0) { PetscCheck(!matB_1 || matB_1->compressedrow.use, PETSC_COMM_SELF, PETSC_ERR_PLIB, "matB_1 && !matB_1->compressedrow.use: PETSc bug???"); } 789d529f056Smarkadams4 /* get state of locals and selected gid for deleted */ 790d529f056Smarkadams4 for (lid = 0; lid < nloc; lid++) { 791d529f056Smarkadams4 lid_parent_gid[lid] = -1.0; 792d529f056Smarkadams4 lid_state[lid] = DELETED; 793d529f056Smarkadams4 } 794d529f056Smarkadams4 795d529f056Smarkadams4 /* set lid_state */ 796d529f056Smarkadams4 for (lid = 0; lid < nloc; lid++) { 797d529f056Smarkadams4 PetscCDIntNd *pos; 798d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos)); 799d529f056Smarkadams4 if (pos) { 800d529f056Smarkadams4 PetscInt gid1; 801d529f056Smarkadams4 802d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gid1)); 803d529f056Smarkadams4 PetscCheck(gid1 == lid + my0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "gid1 %d != lid %d + my0 %d", (int)gid1, (int)lid, (int)my0); 804d529f056Smarkadams4 lid_state[lid] = gid1; 805d529f056Smarkadams4 } 806d529f056Smarkadams4 } 807d529f056Smarkadams4 808d529f056Smarkadams4 /* map local to selected local, DELETED means a ghost owns it */ 809d529f056Smarkadams4 for (lid = kk = 0; lid < nloc; lid++) { 810d529f056Smarkadams4 NState state = lid_state[lid]; 811d529f056Smarkadams4 if (IS_SELECTED(state)) { 812d529f056Smarkadams4 PetscCDIntNd *pos; 813d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos)); 814d529f056Smarkadams4 while (pos) { 815d529f056Smarkadams4 PetscInt gid1; 816d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gid1)); 817d529f056Smarkadams4 PetscCall(PetscCDGetNextPos(aggs_2, lid, &pos)); 818d529f056Smarkadams4 if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1 - my0] = (PetscScalar)(lid + my0); 819d529f056Smarkadams4 } 820d529f056Smarkadams4 } 821d529f056Smarkadams4 } 822d529f056Smarkadams4 /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */ 823d529f056Smarkadams4 if (isMPI) { 824d529f056Smarkadams4 Vec tempVec; 825d529f056Smarkadams4 /* get 'cpcol_1_state' */ 826d529f056Smarkadams4 PetscCall(MatCreateVecs(Gmat_1, &tempVec, NULL)); 827d529f056Smarkadams4 for (kk = 0, j = my0; kk < nloc; kk++, j++) { 828d529f056Smarkadams4 PetscScalar v = (PetscScalar)lid_state[kk]; 829d529f056Smarkadams4 PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES)); 830d529f056Smarkadams4 } 831d529f056Smarkadams4 PetscCall(VecAssemblyBegin(tempVec)); 832d529f056Smarkadams4 PetscCall(VecAssemblyEnd(tempVec)); 833d529f056Smarkadams4 PetscCall(VecScatterBegin(mpimat_1->Mvctx, tempVec, mpimat_1->lvec, INSERT_VALUES, SCATTER_FORWARD)); 834d529f056Smarkadams4 PetscCall(VecScatterEnd(mpimat_1->Mvctx, tempVec, mpimat_1->lvec, INSERT_VALUES, SCATTER_FORWARD)); 835d529f056Smarkadams4 PetscCall(VecGetArray(mpimat_1->lvec, &cpcol_1_state)); 836d529f056Smarkadams4 /* get 'cpcol_2_state' */ 837d529f056Smarkadams4 PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, mpimat_2->lvec, INSERT_VALUES, SCATTER_FORWARD)); 838d529f056Smarkadams4 PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, mpimat_2->lvec, INSERT_VALUES, SCATTER_FORWARD)); 839d529f056Smarkadams4 PetscCall(VecGetArray(mpimat_2->lvec, &cpcol_2_state)); 840d529f056Smarkadams4 /* get 'cpcol_2_par_orig' */ 841d529f056Smarkadams4 for (kk = 0, j = my0; kk < nloc; kk++, j++) { 842d529f056Smarkadams4 PetscScalar v = (PetscScalar)lid_parent_gid[kk]; 843d529f056Smarkadams4 PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES)); 844d529f056Smarkadams4 } 845d529f056Smarkadams4 PetscCall(VecAssemblyBegin(tempVec)); 846d529f056Smarkadams4 PetscCall(VecAssemblyEnd(tempVec)); 847d529f056Smarkadams4 PetscCall(VecDuplicate(mpimat_2->lvec, &ghost_par_orig2)); 848d529f056Smarkadams4 PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghost_par_orig2, INSERT_VALUES, SCATTER_FORWARD)); 849d529f056Smarkadams4 PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghost_par_orig2, INSERT_VALUES, SCATTER_FORWARD)); 850d529f056Smarkadams4 PetscCall(VecGetArray(ghost_par_orig2, &cpcol_2_par_orig)); 851d529f056Smarkadams4 852d529f056Smarkadams4 PetscCall(VecDestroy(&tempVec)); 853d529f056Smarkadams4 } /* ismpi */ 854d529f056Smarkadams4 for (lid = 0; lid < nloc; lid++) { 855d529f056Smarkadams4 NState state = lid_state[lid]; 856d529f056Smarkadams4 if (IS_SELECTED(state)) { 857d529f056Smarkadams4 /* steal locals */ 858d529f056Smarkadams4 ii = matA_1->i; 859d529f056Smarkadams4 n = ii[lid + 1] - ii[lid]; 860d529f056Smarkadams4 idx = matA_1->j + ii[lid]; 861d529f056Smarkadams4 for (j = 0; j < n; j++) { 862d529f056Smarkadams4 PetscInt lidj = idx[j], sgid; 863d529f056Smarkadams4 NState statej = lid_state[lidj]; 864d529f056Smarkadams4 if (statej == DELETED && (sgid = (PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid + my0) { /* steal local */ 865d529f056Smarkadams4 lid_parent_gid[lidj] = (PetscScalar)(lid + my0); /* send this if sgid is not local */ 866d529f056Smarkadams4 if (sgid >= my0 && sgid < Iend) { /* I'm stealing this local from a local sgid */ 867d529f056Smarkadams4 PetscInt hav = 0, slid = sgid - my0, gidj = lidj + my0; 868d529f056Smarkadams4 PetscCDIntNd *pos, *last = NULL; 869d529f056Smarkadams4 /* looking for local from local so id_llist_2 works */ 870d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, slid, &pos)); 871d529f056Smarkadams4 while (pos) { 872d529f056Smarkadams4 PetscInt gid; 873d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gid)); 874d529f056Smarkadams4 if (gid == gidj) { 875d529f056Smarkadams4 PetscCheck(last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "last cannot be null"); 876d529f056Smarkadams4 PetscCall(PetscCDRemoveNextNode(aggs_2, slid, last)); 877d529f056Smarkadams4 PetscCall(PetscCDAppendNode(aggs_2, lid, pos)); 878d529f056Smarkadams4 hav = 1; 879d529f056Smarkadams4 break; 880d529f056Smarkadams4 } else last = pos; 881d529f056Smarkadams4 PetscCall(PetscCDGetNextPos(aggs_2, slid, &pos)); 882d529f056Smarkadams4 } 883d529f056Smarkadams4 if (hav != 1) { 884d529f056Smarkadams4 PetscCheck(hav, PETSC_COMM_SELF, PETSC_ERR_PLIB, "failed to find adj in 'selected' lists - structurally unsymmetric matrix"); 885d529f056Smarkadams4 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "found node %d times???", (int)hav); 886d529f056Smarkadams4 } 887d529f056Smarkadams4 } else { /* I'm stealing this local, owned by a ghost */ 888d529f056Smarkadams4 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.", 889d529f056Smarkadams4 ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : ""); 890d529f056Smarkadams4 PetscCall(PetscCDAppendID(aggs_2, lid, lidj + my0)); 891d529f056Smarkadams4 } 892d529f056Smarkadams4 } 893d529f056Smarkadams4 } /* local neighbors */ 894d529f056Smarkadams4 } else if (state == DELETED /* && lid_cprowID_1 */) { 895d529f056Smarkadams4 PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]); 896d529f056Smarkadams4 /* see if I have a selected ghost neighbor that will steal me */ 897d529f056Smarkadams4 if ((ix = lid_cprowID_1[lid]) != -1) { 898d529f056Smarkadams4 ii = matB_1->compressedrow.i; 899d529f056Smarkadams4 n = ii[ix + 1] - ii[ix]; 900d529f056Smarkadams4 idx = matB_1->j + ii[ix]; 901d529f056Smarkadams4 for (j = 0; j < n; j++) { 902d529f056Smarkadams4 PetscInt cpid = idx[j]; 903d529f056Smarkadams4 NState statej = (NState)PetscRealPart(cpcol_1_state[cpid]); 904d529f056Smarkadams4 if (IS_SELECTED(statej) && sgidold != (PetscInt)statej) { /* ghost will steal this, remove from my list */ 905d529f056Smarkadams4 lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */ 906d529f056Smarkadams4 if (sgidold >= my0 && sgidold < Iend) { /* this was mine */ 907d529f056Smarkadams4 PetscInt hav = 0, oldslidj = sgidold - my0; 908d529f056Smarkadams4 PetscCDIntNd *pos, *last = NULL; 909d529f056Smarkadams4 /* remove from 'oldslidj' list */ 910d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, oldslidj, &pos)); 911d529f056Smarkadams4 while (pos) { 912d529f056Smarkadams4 PetscInt gid; 913d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gid)); 914d529f056Smarkadams4 if (lid + my0 == gid) { 915d529f056Smarkadams4 /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove lid from oldslidj list *\/ */ 916d529f056Smarkadams4 PetscCheck(last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "last cannot be null"); 917d529f056Smarkadams4 PetscCall(PetscCDRemoveNextNode(aggs_2, oldslidj, last)); 918d529f056Smarkadams4 /* ghost (PetscScalar)statej will add this later */ 919d529f056Smarkadams4 hav = 1; 920d529f056Smarkadams4 break; 921d529f056Smarkadams4 } else last = pos; 922d529f056Smarkadams4 PetscCall(PetscCDGetNextPos(aggs_2, oldslidj, &pos)); 923d529f056Smarkadams4 } 924d529f056Smarkadams4 if (hav != 1) { 925d529f056Smarkadams4 PetscCheck(hav, PETSC_COMM_SELF, PETSC_ERR_PLIB, "failed to find (hav=%d) adj in 'selected' lists - structurally unsymmetric matrix", (int)hav); 926d529f056Smarkadams4 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "found node %d times???", (int)hav); 927d529f056Smarkadams4 } 928d529f056Smarkadams4 } else { 929d529f056Smarkadams4 /* TODO: ghosts remove this later */ 930d529f056Smarkadams4 } 931d529f056Smarkadams4 } 932d529f056Smarkadams4 } 933d529f056Smarkadams4 } 934d529f056Smarkadams4 } /* selected/deleted */ 935d529f056Smarkadams4 } /* node loop */ 936d529f056Smarkadams4 937d529f056Smarkadams4 if (isMPI) { 938d529f056Smarkadams4 PetscScalar *cpcol_2_parent, *cpcol_2_gid; 939d529f056Smarkadams4 Vec tempVec, ghostgids2, ghostparents2; 940d529f056Smarkadams4 PetscInt cpid, nghost_2; 941d529f056Smarkadams4 PCGAMGHashTable gid_cpid; 942d529f056Smarkadams4 943d529f056Smarkadams4 PetscCall(VecGetSize(mpimat_2->lvec, &nghost_2)); 944d529f056Smarkadams4 PetscCall(MatCreateVecs(Gmat_2, &tempVec, NULL)); 945d529f056Smarkadams4 946d529f056Smarkadams4 /* get 'cpcol_2_parent' */ 947d529f056Smarkadams4 for (kk = 0, j = my0; kk < nloc; kk++, j++) { PetscCall(VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES)); } 948d529f056Smarkadams4 PetscCall(VecAssemblyBegin(tempVec)); 949d529f056Smarkadams4 PetscCall(VecAssemblyEnd(tempVec)); 950d529f056Smarkadams4 PetscCall(VecDuplicate(mpimat_2->lvec, &ghostparents2)); 951d529f056Smarkadams4 PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghostparents2, INSERT_VALUES, SCATTER_FORWARD)); 952d529f056Smarkadams4 PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghostparents2, INSERT_VALUES, SCATTER_FORWARD)); 953d529f056Smarkadams4 PetscCall(VecGetArray(ghostparents2, &cpcol_2_parent)); 954d529f056Smarkadams4 955d529f056Smarkadams4 /* get 'cpcol_2_gid' */ 956d529f056Smarkadams4 for (kk = 0, j = my0; kk < nloc; kk++, j++) { 957d529f056Smarkadams4 PetscScalar v = (PetscScalar)j; 958d529f056Smarkadams4 PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES)); 959d529f056Smarkadams4 } 960d529f056Smarkadams4 PetscCall(VecAssemblyBegin(tempVec)); 961d529f056Smarkadams4 PetscCall(VecAssemblyEnd(tempVec)); 962d529f056Smarkadams4 PetscCall(VecDuplicate(mpimat_2->lvec, &ghostgids2)); 963d529f056Smarkadams4 PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghostgids2, INSERT_VALUES, SCATTER_FORWARD)); 964d529f056Smarkadams4 PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghostgids2, INSERT_VALUES, SCATTER_FORWARD)); 965d529f056Smarkadams4 PetscCall(VecGetArray(ghostgids2, &cpcol_2_gid)); 966d529f056Smarkadams4 PetscCall(VecDestroy(&tempVec)); 967d529f056Smarkadams4 968d529f056Smarkadams4 /* look for deleted ghosts and add to table */ 969d529f056Smarkadams4 PetscCall(PCGAMGHashTableCreate(2 * nghost_2 + 1, &gid_cpid)); 970d529f056Smarkadams4 for (cpid = 0; cpid < nghost_2; cpid++) { 971d529f056Smarkadams4 NState state = (NState)PetscRealPart(cpcol_2_state[cpid]); 972d529f056Smarkadams4 if (state == DELETED) { 973d529f056Smarkadams4 PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 974d529f056Smarkadams4 PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]); 975d529f056Smarkadams4 if (sgid_old == -1 && sgid_new != -1) { 976d529f056Smarkadams4 PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 977d529f056Smarkadams4 PetscCall(PCGAMGHashTableAdd(&gid_cpid, gid, cpid)); 978d529f056Smarkadams4 } 979d529f056Smarkadams4 } 980d529f056Smarkadams4 } 981d529f056Smarkadams4 982d529f056Smarkadams4 /* look for deleted ghosts and see if they moved - remove it */ 983d529f056Smarkadams4 for (lid = 0; lid < nloc; lid++) { 984d529f056Smarkadams4 NState state = lid_state[lid]; 985d529f056Smarkadams4 if (IS_SELECTED(state)) { 986d529f056Smarkadams4 PetscCDIntNd *pos, *last = NULL; 987d529f056Smarkadams4 /* look for deleted ghosts and see if they moved */ 988d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos)); 989d529f056Smarkadams4 while (pos) { 990d529f056Smarkadams4 PetscInt gid; 991d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gid)); 992d529f056Smarkadams4 993d529f056Smarkadams4 if (gid < my0 || gid >= Iend) { 994d529f056Smarkadams4 PetscCall(PCGAMGHashTableFind(&gid_cpid, gid, &cpid)); 995d529f056Smarkadams4 if (cpid != -1) { 996d529f056Smarkadams4 /* a moved ghost - */ 997d529f056Smarkadams4 /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove 'flid' from list *\/ */ 998d529f056Smarkadams4 PetscCall(PetscCDRemoveNextNode(aggs_2, lid, last)); 999d529f056Smarkadams4 } else last = pos; 1000d529f056Smarkadams4 } else last = pos; 1001d529f056Smarkadams4 1002d529f056Smarkadams4 PetscCall(PetscCDGetNextPos(aggs_2, lid, &pos)); 1003d529f056Smarkadams4 } /* loop over list of deleted */ 1004d529f056Smarkadams4 } /* selected */ 1005d529f056Smarkadams4 } 1006d529f056Smarkadams4 PetscCall(PCGAMGHashTableDestroy(&gid_cpid)); 1007d529f056Smarkadams4 1008d529f056Smarkadams4 /* look at ghosts, see if they changed - and it */ 1009d529f056Smarkadams4 for (cpid = 0; cpid < nghost_2; cpid++) { 1010d529f056Smarkadams4 PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 1011d529f056Smarkadams4 if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */ 1012d529f056Smarkadams4 PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 1013d529f056Smarkadams4 PetscInt slid_new = sgid_new - my0, hav = 0; 1014d529f056Smarkadams4 PetscCDIntNd *pos; 1015d529f056Smarkadams4 1016d529f056Smarkadams4 /* search for this gid to see if I have it */ 1017d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, slid_new, &pos)); 1018d529f056Smarkadams4 while (pos) { 1019d529f056Smarkadams4 PetscInt gidj; 1020d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gidj)); 1021d529f056Smarkadams4 PetscCall(PetscCDGetNextPos(aggs_2, slid_new, &pos)); 1022d529f056Smarkadams4 1023d529f056Smarkadams4 if (gidj == gid) { 1024d529f056Smarkadams4 hav = 1; 1025d529f056Smarkadams4 break; 1026d529f056Smarkadams4 } 1027d529f056Smarkadams4 } 1028d529f056Smarkadams4 if (hav != 1) { 1029d529f056Smarkadams4 /* insert 'flidj' into head of llist */ 1030d529f056Smarkadams4 PetscCall(PetscCDAppendID(aggs_2, slid_new, gid)); 1031d529f056Smarkadams4 } 1032d529f056Smarkadams4 } 1033d529f056Smarkadams4 } 1034d529f056Smarkadams4 PetscCall(VecRestoreArray(mpimat_1->lvec, &cpcol_1_state)); 1035d529f056Smarkadams4 PetscCall(VecRestoreArray(mpimat_2->lvec, &cpcol_2_state)); 1036d529f056Smarkadams4 PetscCall(VecRestoreArray(ghostparents2, &cpcol_2_parent)); 1037d529f056Smarkadams4 PetscCall(VecRestoreArray(ghostgids2, &cpcol_2_gid)); 1038d529f056Smarkadams4 PetscCall(VecDestroy(&ghostgids2)); 1039d529f056Smarkadams4 PetscCall(VecDestroy(&ghostparents2)); 1040d529f056Smarkadams4 PetscCall(VecDestroy(&ghost_par_orig2)); 1041d529f056Smarkadams4 } 1042d529f056Smarkadams4 PetscCall(PetscFree3(lid_state, lid_parent_gid, lid_cprowID_1)); 1043d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 1044d529f056Smarkadams4 } 1045d529f056Smarkadams4 1046c8b0795cSMark F. Adams /* 1047bae903cbSmarkadams4 PCGAMGCoarsen_AGG - supports squaring the graph (deprecated) and new graph for 1048bae903cbSmarkadams4 communication of QR data used with HEM and MISk coarsening 1049c8b0795cSMark F. Adams 1050c8b0795cSMark F. Adams Input Parameter: 1051e0940f08SMark F. Adams . a_pc - this 1052bae903cbSmarkadams4 1053e0940f08SMark F. Adams Input/Output Parameter: 1054bae903cbSmarkadams4 . a_Gmat1 - graph to coarsen (in), graph off processor edges for QR gather scatter (out) 1055bae903cbSmarkadams4 1056c8b0795cSMark F. Adams Output Parameter: 10570cbbd2e1SMark F. Adams . agg_lists - list of aggregates 1058bae903cbSmarkadams4 1059c8b0795cSMark F. Adams */ 1060d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc, Mat *a_Gmat1, PetscCoarsenData **agg_lists) 1061d71ae5a4SJacob Faibussowitsch { 1062e0940f08SMark F. Adams PC_MG *mg = (PC_MG *)a_pc->data; 1063c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1064c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 10658926f930SMark Adams Mat Gmat2, Gmat1 = *a_Gmat1; /* aggressive graph */ 10660cbbd2e1SMark F. Adams IS perm; 1067bae903cbSmarkadams4 PetscInt Istart, Iend, Ii, nloc, bs, nn; 1068bae903cbSmarkadams4 PetscInt *permute, *degree; 1069c8b0795cSMark F. Adams PetscBool *bIndexSet; 1070aea53286SMark Adams PetscReal hashfact; 1071e2c00dcbSBarry Smith PetscInt iSwapIndex; 10723b50add6SMark Adams PetscRandom random; 1073d529f056Smarkadams4 MPI_Comm comm; 1074c8b0795cSMark F. Adams 1075c8b0795cSMark F. Adams PetscFunctionBegin; 1076d529f056Smarkadams4 PetscCall(PetscObjectGetComm((PetscObject)Gmat1, &comm)); 1077849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0)); 1078bae903cbSmarkadams4 PetscCall(MatGetLocalSize(Gmat1, &nn, NULL)); 10799566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Gmat1, &bs)); 108063a3b9bcSJacob Faibussowitsch PetscCheck(bs == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bs %" PetscInt_FMT " must be 1", bs); 1081bae903cbSmarkadams4 nloc = nn / bs; 10825cfd4bd9SMark Adams /* get MIS aggs - randomize */ 1083bae903cbSmarkadams4 PetscCall(PetscMalloc2(nloc, &permute, nloc, °ree)); 10849566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nloc, &bIndexSet)); 10856e09b0e3SStefano Zampini for (Ii = 0; Ii < nloc; Ii++) permute[Ii] = Ii; 10869566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &random)); 10879566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Gmat1, &Istart, &Iend)); 1088c8b0795cSMark F. Adams for (Ii = 0; Ii < nloc; Ii++) { 1089bae903cbSmarkadams4 PetscInt nc; 1090bae903cbSmarkadams4 PetscCall(MatGetRow(Gmat1, Istart + Ii, &nc, NULL, NULL)); 1091bae903cbSmarkadams4 degree[Ii] = nc; 1092bae903cbSmarkadams4 PetscCall(MatRestoreRow(Gmat1, Istart + Ii, &nc, NULL, NULL)); 1093bae903cbSmarkadams4 } 1094bae903cbSmarkadams4 for (Ii = 0; Ii < nloc; Ii++) { 10959566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(random, &hashfact)); 1096aea53286SMark Adams iSwapIndex = (PetscInt)(hashfact * nloc) % nloc; 1097c8b0795cSMark F. Adams if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) { 1098c8b0795cSMark F. Adams PetscInt iTemp = permute[iSwapIndex]; 1099c8b0795cSMark F. Adams permute[iSwapIndex] = permute[Ii]; 1100c8b0795cSMark F. Adams permute[Ii] = iTemp; 1101bae903cbSmarkadams4 iTemp = degree[iSwapIndex]; 1102bae903cbSmarkadams4 degree[iSwapIndex] = degree[Ii]; 1103bae903cbSmarkadams4 degree[Ii] = iTemp; 1104c8b0795cSMark F. Adams bIndexSet[iSwapIndex] = PETSC_TRUE; 1105c8b0795cSMark F. Adams } 1106c8b0795cSMark F. Adams } 1107d529f056Smarkadams4 // apply minimum degree ordering -- NEW 1108d529f056Smarkadams4 if (pc_gamg_agg->use_minimum_degree_ordering) { PetscCall(PetscSortIntWithArray(nloc, degree, permute)); } 11099566063dSJacob Faibussowitsch PetscCall(PetscFree(bIndexSet)); 11109566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&random)); 11119566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm)); 1112849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0)); 1113d529f056Smarkadams4 // square graph 1114d529f056Smarkadams4 if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels && pc_gamg_agg->use_aggressive_square_graph) { 1115d529f056Smarkadams4 PetscCall(PCGAMGSquareGraph_GAMG(a_pc, Gmat1, &Gmat2)); 1116d529f056Smarkadams4 } else Gmat2 = Gmat1; 1117d529f056Smarkadams4 // switch to old MIS-1 for square graph 1118d529f056Smarkadams4 if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels) { 1119d529f056Smarkadams4 if (!pc_gamg_agg->use_aggressive_square_graph) PetscCall(MatCoarsenMISKSetDistance(pc_gamg_agg->crs, pc_gamg_agg->aggressive_mis_k)); // hardwire to MIS-2 1120d529f056Smarkadams4 else PetscCall(MatCoarsenSetType(pc_gamg_agg->crs, MATCOARSENMIS)); // old MIS -- side effect 1121d529f056Smarkadams4 } else if (pc_gamg_agg->use_aggressive_square_graph && pc_gamg_agg->aggressive_coarsening_levels > 0) { // we reset the MIS 1122d529f056Smarkadams4 const char *prefix; 1123d529f056Smarkadams4 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)a_pc, &prefix)); 1124d529f056Smarkadams4 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc_gamg_agg->crs, prefix)); 1125d529f056Smarkadams4 PetscCall(MatCoarsenSetFromOptions(pc_gamg_agg->crs)); // get the default back on non-aggressive levels when square graph switched to old MIS 1126d529f056Smarkadams4 } 1127d529f056Smarkadams4 PetscCall(MatCoarsenSetAdjacency(pc_gamg_agg->crs, Gmat2)); 11282d776b49SBarry Smith PetscCall(MatCoarsenSetStrictAggs(pc_gamg_agg->crs, PETSC_TRUE)); 1129d529f056Smarkadams4 PetscCall(MatCoarsenSetGreedyOrdering(pc_gamg_agg->crs, perm)); 11302d776b49SBarry Smith PetscCall(MatCoarsenApply(pc_gamg_agg->crs)); 11312d776b49SBarry Smith PetscCall(MatCoarsenViewFromOptions(pc_gamg_agg->crs, NULL, "-mat_coarsen_view")); 11322d776b49SBarry Smith PetscCall(MatCoarsenGetData(pc_gamg_agg->crs, agg_lists)); /* output */ 11332d776b49SBarry Smith PetscCall(MatCoarsenDestroy(&pc_gamg_agg->crs)); 1134b43b03e9SMark F. Adams 11359566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 1136bae903cbSmarkadams4 PetscCall(PetscFree2(permute, degree)); 1137849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0)); 11389f3f12c8SMark F. Adams 11399c15c1aeSMark Adams if (Gmat2 != Gmat1) { // square graph, we need ghosts for selected 1140d529f056Smarkadams4 PetscCoarsenData *llist = *agg_lists; 1141d529f056Smarkadams4 PetscCall(fixAggregatesWithSquare(a_pc, Gmat2, Gmat1, *agg_lists)); 1142d529f056Smarkadams4 PetscCall(MatDestroy(&Gmat1)); 1143d529f056Smarkadams4 *a_Gmat1 = Gmat2; /* output */ 11448926f930SMark Adams PetscCall(PetscCDSetMat(llist, *a_Gmat1)); /* Need a graph with ghosts here */ 11450cbbd2e1SMark F. Adams } 1146849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0)); 11473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1148c8b0795cSMark F. Adams } 1149c8b0795cSMark F. Adams 1150c8b0795cSMark F. Adams /* 11510cbbd2e1SMark F. Adams PCGAMGProlongator_AGG 1152c8b0795cSMark F. Adams 1153c8b0795cSMark F. Adams Input Parameter: 1154c8b0795cSMark F. Adams . pc - this 1155c8b0795cSMark F. Adams . Amat - matrix on this fine level 1156c8b0795cSMark F. Adams . Graph - used to get ghost data for nodes in 11570cbbd2e1SMark F. Adams . agg_lists - list of aggregates 1158c8b0795cSMark F. Adams Output Parameter: 1159c8b0795cSMark F. Adams . a_P_out - prolongation operator to the next level 1160c8b0795cSMark F. Adams */ 11618926f930SMark Adams static PetscErrorCode PCGAMGProlongator_AGG(PC pc, Mat Amat, PetscCoarsenData *agg_lists, Mat *a_P_out) 1162d71ae5a4SJacob Faibussowitsch { 11632e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 11642e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 11654a2f8832SBarry Smith const PetscInt col_bs = pc_gamg->data_cell_cols; 1166c8b0795cSMark F. Adams PetscInt Istart, Iend, nloc, ii, jj, kk, my0, nLocalSelected, bs; 11678926f930SMark Adams Mat Gmat, Prol; 1168d5d25401SBarry Smith PetscMPIInt size; 11693b4367a7SBarry Smith MPI_Comm comm; 1170c8b0795cSMark F. Adams PetscReal *data_w_ghost; 1171c8b0795cSMark F. Adams PetscInt myCrs0, nbnodes = 0, *flid_fgid; 1172d9558ea9SBarry Smith MatType mtype; 11732e68589bSMark F. Adams 11742e68589bSMark F. Adams PetscFunctionBegin; 11759566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm)); 117608401ef6SPierre Jolivet PetscCheck(col_bs >= 1, comm, PETSC_ERR_PLIB, "Column bs cannot be less than 1"); 1177849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0)); 11789566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 11799566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Amat, &Istart, &Iend)); 11809566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Amat, &bs)); 11819371c9d4SSatish Balay nloc = (Iend - Istart) / bs; 11829371c9d4SSatish Balay my0 = Istart / bs; 118363a3b9bcSJacob 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); 11848926f930SMark Adams PetscCall(PetscCDGetMat(agg_lists, &Gmat)); // get auxilary matrix for ghost edges for size > 1 11852e68589bSMark F. Adams 11862e68589bSMark F. Adams /* get 'nLocalSelected' */ 11870cbbd2e1SMark F. Adams for (ii = 0, nLocalSelected = 0; ii < nloc; ii++) { 1188e78576d6SMark F. Adams PetscBool ise; 11890cbbd2e1SMark F. Adams /* filter out singletons 0 or 1? */ 11905e62d202SMark Adams PetscCall(PetscCDIsEmptyAt(agg_lists, ii, &ise)); 1191e78576d6SMark F. Adams if (!ise) nLocalSelected++; 11922e68589bSMark F. Adams } 11932e68589bSMark F. Adams 11942e68589bSMark F. Adams /* create prolongator, create P matrix */ 11959566063dSJacob Faibussowitsch PetscCall(MatGetType(Amat, &mtype)); 11969566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &Prol)); 11979566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Prol, nloc * bs, nLocalSelected * col_bs, PETSC_DETERMINE, PETSC_DETERMINE)); 11989566063dSJacob Faibussowitsch PetscCall(MatSetBlockSizes(Prol, bs, col_bs)); 11999566063dSJacob Faibussowitsch PetscCall(MatSetType(Prol, mtype)); 12003742c8caSstefanozampini #if PetscDefined(HAVE_DEVICE) 12013742c8caSstefanozampini PetscBool flg; 12023742c8caSstefanozampini PetscCall(MatBoundToCPU(Amat, &flg)); 12033742c8caSstefanozampini PetscCall(MatBindToCPU(Prol, flg)); 12043742c8caSstefanozampini if (flg) PetscCall(MatSetBindingPropagates(Prol, PETSC_TRUE)); 12053742c8caSstefanozampini #endif 12069566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(Prol, col_bs, NULL)); 12079566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(Prol, col_bs, NULL, col_bs, NULL)); 12082e68589bSMark F. Adams 12092e68589bSMark F. Adams /* can get all points "removed" */ 12109566063dSJacob Faibussowitsch PetscCall(MatGetSize(Prol, &kk, &ii)); 12117f66b68fSMark Adams if (!ii) { 121263a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: No selected points on coarse grid\n", ((PetscObject)pc)->prefix)); 12139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Prol)); 12140298fd71SBarry Smith *a_P_out = NULL; /* out */ 1215849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0)); 12163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12172e68589bSMark F. Adams } 121863a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: New grid %" PetscInt_FMT " nodes\n", ((PetscObject)pc)->prefix, ii / col_bs)); 12199566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk)); 12200cbbd2e1SMark F. Adams 122163a3b9bcSJacob 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); 1222c8b0795cSMark F. Adams myCrs0 = myCrs0 / col_bs; 122363a3b9bcSJacob 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); 12242e68589bSMark F. Adams 12252e68589bSMark F. Adams /* create global vector of data in 'data_w_ghost' */ 1226849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0)); 1227bae903cbSmarkadams4 if (size > 1) { /* get ghost null space data */ 12282e68589bSMark F. Adams PetscReal *tmp_gdata, *tmp_ldata, *tp2; 12299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc, &tmp_ldata)); 12304a2f8832SBarry Smith for (jj = 0; jj < col_bs; jj++) { 1231c8b0795cSMark F. Adams for (kk = 0; kk < bs; kk++) { 1232a2f3521dSMark F. Adams PetscInt ii, stride; 1233*8e3a54c0SPierre Jolivet const PetscReal *tp = PetscSafePointerPlusOffset(pc_gamg->data, jj * bs * nloc + kk); 12342fa5cd67SKarl Rupp for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp; 12352fa5cd67SKarl Rupp 12369566063dSJacob Faibussowitsch PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata)); 1237a2f3521dSMark F. Adams 1238bae903cbSmarkadams4 if (!jj && !kk) { /* now I know how many total nodes - allocate TODO: move below and do in one 'col_bs' call */ 12399566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(stride * bs * col_bs, &data_w_ghost)); 1240a2f3521dSMark F. Adams nbnodes = bs * stride; 12412e68589bSMark F. Adams } 1242*8e3a54c0SPierre Jolivet tp2 = PetscSafePointerPlusOffset(data_w_ghost, jj * bs * stride + kk); 12432fa5cd67SKarl Rupp for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii]; 12449566063dSJacob Faibussowitsch PetscCall(PetscFree(tmp_gdata)); 12452e68589bSMark F. Adams } 12462e68589bSMark F. Adams } 12479566063dSJacob Faibussowitsch PetscCall(PetscFree(tmp_ldata)); 1248806fa848SBarry Smith } else { 1249c8b0795cSMark F. Adams nbnodes = bs * nloc; 1250c8b0795cSMark F. Adams data_w_ghost = (PetscReal *)pc_gamg->data; 12512e68589bSMark F. Adams } 12522e68589bSMark F. Adams 1253bae903cbSmarkadams4 /* get 'flid_fgid' TODO - move up to get 'stride' and do get null space data above in one step (jj loop) */ 1254c5df96a5SBarry Smith if (size > 1) { 12552e68589bSMark F. Adams PetscReal *fid_glid_loc, *fiddata; 1256a2f3521dSMark F. Adams PetscInt stride; 12572e68589bSMark F. Adams 12589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc, &fid_glid_loc)); 12592e68589bSMark F. Adams for (kk = 0; kk < nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0 + kk); 12609566063dSJacob Faibussowitsch PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata)); 1261bae903cbSmarkadams4 PetscCall(PetscMalloc1(stride, &flid_fgid)); /* copy real data to in */ 1262a2f3521dSMark F. Adams for (kk = 0; kk < stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk]; 12639566063dSJacob Faibussowitsch PetscCall(PetscFree(fiddata)); 1264a2f3521dSMark F. Adams 126563a3b9bcSJacob Faibussowitsch PetscCheck(stride == nbnodes / bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "stride %" PetscInt_FMT " != nbnodes %" PetscInt_FMT "/bs %" PetscInt_FMT, stride, nbnodes, bs); 12669566063dSJacob Faibussowitsch PetscCall(PetscFree(fid_glid_loc)); 1267806fa848SBarry Smith } else { 12689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc, &flid_fgid)); 12692e68589bSMark F. Adams for (kk = 0; kk < nloc; kk++) flid_fgid[kk] = my0 + kk; 12702e68589bSMark F. Adams } 1271849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0)); 1272bae903cbSmarkadams4 /* get P0 */ 1273849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0)); 1274c8b0795cSMark F. Adams { 12750298fd71SBarry Smith PetscReal *data_out = NULL; 12769566063dSJacob Faibussowitsch PetscCall(formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes, data_w_ghost, flid_fgid, &data_out, Prol)); 12779566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 12782fa5cd67SKarl Rupp 1279c8b0795cSMark F. Adams pc_gamg->data = data_out; 12804a2f8832SBarry Smith pc_gamg->data_cell_rows = col_bs; 12814a2f8832SBarry Smith pc_gamg->data_sz = col_bs * col_bs * nLocalSelected; 1282c8b0795cSMark F. Adams } 1283849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0)); 128448a46eb9SPierre Jolivet if (size > 1) PetscCall(PetscFree(data_w_ghost)); 12859566063dSJacob Faibussowitsch PetscCall(PetscFree(flid_fgid)); 12862e68589bSMark F. Adams 1287c8b0795cSMark F. Adams *a_P_out = Prol; /* out */ 12885e62d202SMark Adams PetscCall(MatViewFromOptions(Prol, NULL, "-view_P")); 1289ed4e82eaSPeter Brune 1290849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0)); 12913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1292c8b0795cSMark F. Adams } 1293c8b0795cSMark F. Adams 1294c8b0795cSMark F. Adams /* 1295fd1112cbSBarry Smith PCGAMGOptProlongator_AGG 1296c8b0795cSMark F. Adams 1297c8b0795cSMark F. Adams Input Parameter: 1298c8b0795cSMark F. Adams . pc - this 1299c8b0795cSMark F. Adams . Amat - matrix on this fine level 1300c8b0795cSMark F. Adams In/Output Parameter: 13012a808120SBarry Smith . a_P - prolongation operator to the next level 1302c8b0795cSMark F. Adams */ 1303d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGOptProlongator_AGG(PC pc, Mat Amat, Mat *a_P) 1304d71ae5a4SJacob Faibussowitsch { 1305c8b0795cSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 1306c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1307c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 1308c8b0795cSMark F. Adams PetscInt jj; 1309c8b0795cSMark F. Adams Mat Prol = *a_P; 13103b4367a7SBarry Smith MPI_Comm comm; 13112a808120SBarry Smith KSP eksp; 13122a808120SBarry Smith Vec bb, xx; 13132a808120SBarry Smith PC epc; 13142a808120SBarry Smith PetscReal alpha, emax, emin; 1315c8b0795cSMark F. Adams 1316c8b0795cSMark F. Adams PetscFunctionBegin; 13179566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm)); 1318849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0)); 1319c8b0795cSMark F. Adams 1320d5d25401SBarry Smith /* compute maximum singular value of operator to be used in smoother */ 13212a808120SBarry Smith if (0 < pc_gamg_agg->nsmooths) { 132218c3aa7eSMark /* get eigen estimates */ 132318c3aa7eSMark if (pc_gamg->emax > 0) { 132418c3aa7eSMark emin = pc_gamg->emin; 132518c3aa7eSMark emax = pc_gamg->emax; 132618c3aa7eSMark } else { 13270ed2132dSStefano Zampini const char *prefix; 13280ed2132dSStefano Zampini 13299566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Amat, &bb, NULL)); 13309566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Amat, &xx, NULL)); 13319566063dSJacob Faibussowitsch PetscCall(KSPSetNoisy_Private(bb)); 1332e696ed0bSMark F. Adams 13339566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm, &eksp)); 13343821be0aSBarry Smith PetscCall(KSPSetNestLevel(eksp, pc->kspnestlevel)); 13359566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 13369566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(eksp, prefix)); 13379566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(eksp, "pc_gamg_esteig_")); 133873f7197eSJed Brown { 1339b94d7dedSBarry Smith PetscBool isset, sflg; 1340b94d7dedSBarry Smith PetscCall(MatIsSPDKnown(Amat, &isset, &sflg)); 1341b94d7dedSBarry Smith if (isset && sflg) PetscCall(KSPSetType(eksp, KSPCG)); 1342d24ecf33SMark } 13439566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(eksp, pc->erroriffailure)); 13449566063dSJacob Faibussowitsch PetscCall(KSPSetNormType(eksp, KSP_NORM_NONE)); 1345c2436cacSMark F. Adams 13469566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(eksp, PETSC_FALSE)); 13479566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(eksp, Amat, Amat)); 13482e68589bSMark F. Adams 13499566063dSJacob Faibussowitsch PetscCall(KSPGetPC(eksp, &epc)); 13509566063dSJacob Faibussowitsch PetscCall(PCSetType(epc, PCJACOBI)); /* smoother in smoothed agg. */ 1351c2436cacSMark F. Adams 13529566063dSJacob 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 1353074aec55SMark Adams 13549566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(eksp)); 13559566063dSJacob Faibussowitsch PetscCall(KSPSetComputeSingularValues(eksp, PETSC_TRUE)); 13569566063dSJacob Faibussowitsch PetscCall(KSPSolve(eksp, bb, xx)); 13579566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(eksp, pc, xx)); 13582e68589bSMark F. Adams 13599566063dSJacob Faibussowitsch PetscCall(KSPComputeExtremeSingularValues(eksp, &emax, &emin)); 13609566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: Smooth P0: max eigen=%e min=%e PC=%s\n", ((PetscObject)pc)->prefix, (double)emax, (double)emin, PCJACOBI)); 13619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xx)); 13629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bb)); 13639566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&eksp)); 13642e68589bSMark F. Adams } 136518c3aa7eSMark if (pc_gamg->use_sa_esteig) { 136618c3aa7eSMark mg->min_eigen_DinvA[pc_gamg->current_level] = emin; 136718c3aa7eSMark mg->max_eigen_DinvA[pc_gamg->current_level] = emax; 136863a3b9bcSJacob 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)); 136918c3aa7eSMark } else { 137018c3aa7eSMark mg->min_eigen_DinvA[pc_gamg->current_level] = 0; 137118c3aa7eSMark mg->max_eigen_DinvA[pc_gamg->current_level] = 0; 137218c3aa7eSMark } 137318c3aa7eSMark } else { 137418c3aa7eSMark mg->min_eigen_DinvA[pc_gamg->current_level] = 0; 137518c3aa7eSMark mg->max_eigen_DinvA[pc_gamg->current_level] = 0; 137618c3aa7eSMark } 13772e68589bSMark F. Adams 13782a808120SBarry Smith /* smooth P0 */ 13792a808120SBarry Smith for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) { 13802a808120SBarry Smith Mat tMat; 13812a808120SBarry Smith Vec diag; 13822a808120SBarry Smith 1383849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0)); 13842a808120SBarry Smith 13852e68589bSMark F. Adams /* smooth P1 := (I - omega/lam D^{-1}A)P0 */ 13869566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0)); 13879566063dSJacob Faibussowitsch PetscCall(MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat)); 13889566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0)); 13899566063dSJacob Faibussowitsch PetscCall(MatProductClear(tMat)); 13909566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Amat, &diag, NULL)); 13919566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(Amat, diag)); /* effectively PCJACOBI */ 13929566063dSJacob Faibussowitsch PetscCall(VecReciprocal(diag)); 13939566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(tMat, diag, NULL)); 13949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&diag)); 1395b4da3a1bSStefano Zampini 1396d5d25401SBarry Smith /* TODO: Set a PCFailedReason and exit the building of the AMG preconditioner */ 139708401ef6SPierre Jolivet PetscCheck(emax != 0.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Computed maximum singular value as zero"); 1398d5d25401SBarry Smith /* TODO: Document the 1.4 and don't hardwire it in this routine */ 1399b4ec6429SMark F. Adams alpha = -1.4 / emax; 1400b4da3a1bSStefano Zampini 14019566063dSJacob Faibussowitsch PetscCall(MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN)); 14029566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Prol)); 14032e68589bSMark F. Adams Prol = tMat; 1404849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0)); 14052e68589bSMark F. Adams } 1406849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0)); 1407c8b0795cSMark F. Adams *a_P = Prol; 14083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14092e68589bSMark F. Adams } 14102e68589bSMark F. Adams 1411c8b0795cSMark F. Adams /* 1412c8b0795cSMark F. Adams PCCreateGAMG_AGG 14132e68589bSMark F. Adams 1414c8b0795cSMark F. Adams Input Parameter: 1415c8b0795cSMark F. Adams . pc - 1416c8b0795cSMark F. Adams */ 1417d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCreateGAMG_AGG(PC pc) 1418d71ae5a4SJacob Faibussowitsch { 1419c8b0795cSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 1420c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1421c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg; 14222e68589bSMark F. Adams 1423c8b0795cSMark F. Adams PetscFunctionBegin; 1424c8b0795cSMark F. Adams /* create sub context for SA */ 14254dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&pc_gamg_agg)); 1426c8b0795cSMark F. Adams pc_gamg->subctx = pc_gamg_agg; 1427c8b0795cSMark F. Adams 14281ab5ffc9SJed Brown pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG; 14299b8ffb57SJed Brown pc_gamg->ops->destroy = PCDestroy_GAMG_AGG; 1430c8b0795cSMark F. Adams /* reset does not do anything; setup not virtual */ 1431c8b0795cSMark F. Adams 1432c8b0795cSMark F. Adams /* set internal function pointers */ 14332d776b49SBarry Smith pc_gamg->ops->creategraph = PCGAMGCreateGraph_AGG; 14341ab5ffc9SJed Brown pc_gamg->ops->coarsen = PCGAMGCoarsen_AGG; 14351ab5ffc9SJed Brown pc_gamg->ops->prolongator = PCGAMGProlongator_AGG; 1436fd1112cbSBarry Smith pc_gamg->ops->optprolongator = PCGAMGOptProlongator_AGG; 14371ab5ffc9SJed Brown pc_gamg->ops->createdefaultdata = PCSetData_AGG; 14385adeb434SBarry Smith pc_gamg->ops->view = PCView_GAMG_AGG; 1439c8b0795cSMark F. Adams 1440cab9ed1eSBarry Smith pc_gamg_agg->nsmooths = 1; 1441d529f056Smarkadams4 pc_gamg_agg->aggressive_coarsening_levels = 1; 1442d4adc10fSMark Adams pc_gamg_agg->use_aggressive_square_graph = PETSC_TRUE; 1443d529f056Smarkadams4 pc_gamg_agg->use_minimum_degree_ordering = PETSC_FALSE; 1444a9ccf005SMark Adams pc_gamg_agg->use_low_mem_filter = PETSC_FALSE; 1445d529f056Smarkadams4 pc_gamg_agg->aggressive_mis_k = 2; 1446cab9ed1eSBarry Smith 14479566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", PCGAMGSetNSmooths_AGG)); 1448bae903cbSmarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", PCGAMGSetAggressiveLevels_AGG)); 1449d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveSquareGraph_C", PCGAMGSetAggressiveSquareGraph_AGG)); 1450d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetMinDegreeOrdering_C", PCGAMGMISkSetMinDegreeOrdering_AGG)); 1451a9ccf005SMark Adams PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetLowMemoryFilter_C", PCGAMGSetLowMemoryFilter_AGG)); 1452d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetAggressive_C", PCGAMGMISkSetAggressive_AGG)); 14539566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_AGG)); 14543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14552e68589bSMark F. Adams } 1456