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; 16*a9ccf005SMark 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 34f1580f4eSBarry Smith .seealso: `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 70*a9ccf005SMark Adams .seealso: `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 95*a9ccf005SMark Adams .seealso: `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 120*a9ccf005SMark Adams .seealso: `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 145*a9ccf005SMark Adams .seealso: `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 156*a9ccf005SMark Adams /*@ 157*a9ccf005SMark Adams PCGAMGSetLowMemoryFilter - Use low memory graph/matrix filter 158*a9ccf005SMark Adams 159*a9ccf005SMark Adams Logically Collective 160*a9ccf005SMark Adams 161*a9ccf005SMark Adams Input Parameters: 162*a9ccf005SMark Adams + pc - the preconditioner context 163*a9ccf005SMark Adams - b - default false 164*a9ccf005SMark Adams 165*a9ccf005SMark Adams Options Database Key: 166*a9ccf005SMark Adams . -pc_gamg_low_memory_threshold_filter <bool,default=false> - Use low memory graph/matrix filter 167*a9ccf005SMark Adams 168*a9ccf005SMark Adams Level: intermediate 169*a9ccf005SMark Adams 170*a9ccf005SMark Adams .seealso: `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`, `PCGAMGMISkSetAggressive()`, `PCGAMGSetAggressiveSquareGraph()`, `PCGAMGMISkSetMinDegreeOrdering()` 171*a9ccf005SMark Adams @*/ 172*a9ccf005SMark Adams PetscErrorCode PCGAMGSetLowMemoryFilter(PC pc, PetscBool b) 173*a9ccf005SMark Adams { 174*a9ccf005SMark Adams PetscFunctionBegin; 175*a9ccf005SMark Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 176*a9ccf005SMark Adams PetscValidLogicalCollectiveBool(pc, b, 2); 177*a9ccf005SMark Adams PetscTryMethod(pc, "PCGAMGSetLowMemoryFilter_C", (PC, PetscBool), (pc, b)); 178*a9ccf005SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 179*a9ccf005SMark Adams } 180*a9ccf005SMark 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 214*a9ccf005SMark Adams static PetscErrorCode PCGAMGSetLowMemoryFilter_AGG(PC pc, PetscBool b) 215*a9ccf005SMark Adams { 216*a9ccf005SMark Adams PC_MG *mg = (PC_MG *)pc->data; 217*a9ccf005SMark Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 218*a9ccf005SMark Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 219*a9ccf005SMark Adams 220*a9ccf005SMark Adams PetscFunctionBegin; 221*a9ccf005SMark Adams pc_gamg_agg->use_low_mem_filter = b; 222*a9ccf005SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 223*a9ccf005SMark Adams } 224*a9ccf005SMark 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; 2412e68589bSMark F. Adams 2422e68589bSMark F. Adams PetscFunctionBegin; 243d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "GAMG-AGG options"); 2442e68589bSMark F. Adams { 245bae903cbSmarkadams4 PetscBool flg; 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)); 247d529f056Smarkadams4 PetscCall(PetscOptionsInt("-pc_gamg_aggressive_coarsening", "Number of aggressive coarsening (MIS-2) levels from finest", "PCGAMGSetAggressiveLevels", pc_gamg_agg->aggressive_coarsening_levels, &pc_gamg_agg->aggressive_coarsening_levels, &flg)); 248d529f056Smarkadams4 if (!flg) { 249d529f056Smarkadams4 PetscCall( 250d529f056Smarkadams4 PetscOptionsInt("-pc_gamg_square_graph", "Number of aggressive coarsening (MIS-2) levels from finest (deprecated alias for -pc_gamg_aggressive_coarsening)", "PCGAMGSetAggressiveLevels", pc_gamg_agg->aggressive_coarsening_levels, &pc_gamg_agg->aggressive_coarsening_levels, NULL)); 251d529f056Smarkadams4 } else { 2529371c9d4SSatish Balay PetscCall( 2539371c9d4SSatish Balay PetscOptionsInt("-pc_gamg_square_graph", "Number of aggressive coarsening (MIS-2) levels from finest (alias for -pc_gamg_aggressive_coarsening, deprecated)", "PCGAMGSetAggressiveLevels", pc_gamg_agg->aggressive_coarsening_levels, &pc_gamg_agg->aggressive_coarsening_levels, &flg)); 254bae903cbSmarkadams4 if (flg) PetscCall(PetscInfo(pc, "Warning: both -pc_gamg_square_graph and -pc_gamg_aggressive_coarsening are used. -pc_gamg_square_graph is deprecated, Number of aggressive levels is %d\n", (int)pc_gamg_agg->aggressive_coarsening_levels)); 255bae903cbSmarkadams4 } 256d529f056Smarkadams4 if (pc_gamg_agg->aggressive_coarsening_levels > 0) { 257d529f056Smarkadams4 PetscCall(PetscOptionsBool("-pc_gamg_aggressive_square_graph", "Use square graph (A'A) or MIS-k (k=2) for aggressive coarsening", "PCGAMGSetAggressiveSquareGraph", pc_gamg_agg->use_aggressive_square_graph, &pc_gamg_agg->use_aggressive_square_graph, NULL)); 258d529f056Smarkadams4 } 259d529f056Smarkadams4 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)); 260*a9ccf005SMark 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)); 261d529f056Smarkadams4 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)); 2622e68589bSMark F. Adams } 263d0609cedSBarry Smith PetscOptionsHeadEnd(); 2643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2652e68589bSMark F. Adams } 2662e68589bSMark F. Adams 267d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_GAMG_AGG(PC pc) 268d71ae5a4SJacob Faibussowitsch { 2692e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 2702e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 2712e68589bSMark F. Adams 2722e68589bSMark F. Adams PetscFunctionBegin; 2739566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->subctx)); 2742e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", NULL)); 275bae903cbSmarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", NULL)); 276d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetAggressive_C", NULL)); 277d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetMinDegreeOrdering_C", NULL)); 278*a9ccf005SMark Adams PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetLowMemoryFilter_C", NULL)); 279d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveSquareGraph_C", NULL)); 280bae903cbSmarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL)); 2813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2822e68589bSMark F. Adams } 2832e68589bSMark F. Adams 2842e68589bSMark F. Adams /* 2852e68589bSMark F. Adams PCSetCoordinates_AGG 28620f4b53cSBarry Smith 28720f4b53cSBarry Smith Collective 2882e68589bSMark F. Adams 2892e68589bSMark F. Adams Input Parameter: 2902e68589bSMark F. Adams . pc - the preconditioner context 291145b44c9SPierre Jolivet . ndm - dimension of data (used for dof/vertex for Stokes) 292302f38e8SMark F. Adams . a_nloc - number of vertices local 293302f38e8SMark F. Adams . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...} 2942e68589bSMark F. Adams */ 2951e6b0712SBarry Smith 296d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords) 297d71ae5a4SJacob Faibussowitsch { 2982e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 2992e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 30069344418SMark F. Adams PetscInt arrsz, kk, ii, jj, nloc, ndatarows, ndf; 301a2f3521dSMark F. Adams Mat mat = pc->pmat; 3022e68589bSMark F. Adams 3032e68589bSMark F. Adams PetscFunctionBegin; 304a2f3521dSMark F. Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 305a2f3521dSMark F. Adams PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 306302f38e8SMark F. Adams nloc = a_nloc; 3072e68589bSMark F. Adams 3082e68589bSMark F. Adams /* SA: null space vectors */ 3099566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(mat, &ndf)); /* this does not work for Stokes */ 31069344418SMark F. Adams if (coords && ndf == 1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */ 311a2f3521dSMark F. Adams else if (coords) { 31263a3b9bcSJacob Faibussowitsch PetscCheck(ndm <= ndf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "degrees of motion %" PetscInt_FMT " > block size %" PetscInt_FMT, ndm, ndf); 31369344418SMark F. Adams pc_gamg->data_cell_cols = (ndm == 2 ? 3 : 6); /* displacement elasticity */ 314ad540459SPierre 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); 315806fa848SBarry Smith } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */ 31671959b99SBarry Smith pc_gamg->data_cell_rows = ndatarows = ndf; 31763a3b9bcSJacob 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); 318c8b0795cSMark F. Adams arrsz = nloc * pc_gamg->data_cell_rows * pc_gamg->data_cell_cols; 3192e68589bSMark F. Adams 3207f66b68fSMark Adams if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) { 3219566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 3229566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(arrsz + 1, &pc_gamg->data)); 3232e68589bSMark F. Adams } 3242e68589bSMark F. Adams /* copy data in - column oriented */ 3252e68589bSMark F. Adams for (kk = 0; kk < nloc; kk++) { 32669344418SMark F. Adams const PetscInt M = nloc * pc_gamg->data_cell_rows; /* stride into data */ 32769344418SMark F. Adams PetscReal *data = &pc_gamg->data[kk * ndatarows]; /* start of cell */ 328c8b0795cSMark F. Adams if (pc_gamg->data_cell_cols == 1) *data = 1.0; 3292e68589bSMark F. Adams else { 33069344418SMark F. Adams /* translational modes */ 3312fa5cd67SKarl Rupp for (ii = 0; ii < ndatarows; ii++) { 3322fa5cd67SKarl Rupp for (jj = 0; jj < ndatarows; jj++) { 33369344418SMark F. Adams if (ii == jj) data[ii * M + jj] = 1.0; 3342e68589bSMark F. Adams else data[ii * M + jj] = 0.0; 3352fa5cd67SKarl Rupp } 3362fa5cd67SKarl Rupp } 33769344418SMark F. Adams 33869344418SMark F. Adams /* rotational modes */ 3392e68589bSMark F. Adams if (coords) { 34069344418SMark F. Adams if (ndm == 2) { 3412e68589bSMark F. Adams data += 2 * M; 3422e68589bSMark F. Adams data[0] = -coords[2 * kk + 1]; 3432e68589bSMark F. Adams data[1] = coords[2 * kk]; 344806fa848SBarry Smith } else { 3452e68589bSMark F. Adams data += 3 * M; 3469371c9d4SSatish Balay data[0] = 0.0; 3479371c9d4SSatish Balay data[M + 0] = coords[3 * kk + 2]; 3489371c9d4SSatish Balay data[2 * M + 0] = -coords[3 * kk + 1]; 3499371c9d4SSatish Balay data[1] = -coords[3 * kk + 2]; 3509371c9d4SSatish Balay data[M + 1] = 0.0; 3519371c9d4SSatish Balay data[2 * M + 1] = coords[3 * kk]; 3529371c9d4SSatish Balay data[2] = coords[3 * kk + 1]; 3539371c9d4SSatish Balay data[M + 2] = -coords[3 * kk]; 3549371c9d4SSatish Balay data[2 * M + 2] = 0.0; 3552e68589bSMark F. Adams } 3562e68589bSMark F. Adams } 3572e68589bSMark F. Adams } 3582e68589bSMark F. Adams } 3592e68589bSMark F. Adams pc_gamg->data_sz = arrsz; 3603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3612e68589bSMark F. Adams } 3622e68589bSMark F. Adams 3632e68589bSMark F. Adams /* 364a2f3521dSMark F. Adams PCSetData_AGG - called if data is not set with PCSetCoordinates. 365a2f3521dSMark F. Adams Looks in Mat for near null space. 366a2f3521dSMark F. Adams Does not work for Stokes 3672e68589bSMark F. Adams 3682e68589bSMark F. Adams Input Parameter: 3692e68589bSMark F. Adams . pc - 370a2f3521dSMark F. Adams . a_A - matrix to get (near) null space out of. 3712e68589bSMark F. Adams */ 372d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A) 373d71ae5a4SJacob Faibussowitsch { 374b8cd405aSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 375b8cd405aSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 376b8cd405aSMark F. Adams MatNullSpace mnull; 37766f2ef4bSMark Adams 378b817416eSBarry Smith PetscFunctionBegin; 3799566063dSJacob Faibussowitsch PetscCall(MatGetNearNullSpace(a_A, &mnull)); 380b8cd405aSMark F. Adams if (!mnull) { 38166f2ef4bSMark Adams DM dm; 3829566063dSJacob Faibussowitsch PetscCall(PCGetDM(pc, &dm)); 38348a46eb9SPierre Jolivet if (!dm) PetscCall(MatGetDM(a_A, &dm)); 38466f2ef4bSMark Adams if (dm) { 38566f2ef4bSMark Adams PetscObject deformation; 386b0d30dd6SMatthew G. Knepley PetscInt Nf; 387b0d30dd6SMatthew G. Knepley 3889566063dSJacob Faibussowitsch PetscCall(DMGetNumFields(dm, &Nf)); 389b0d30dd6SMatthew G. Knepley if (Nf) { 3909566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, 0, NULL, &deformation)); 3919566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)deformation, "nearnullspace", (PetscObject *)&mnull)); 39248a46eb9SPierre Jolivet if (!mnull) PetscCall(PetscObjectQuery((PetscObject)deformation, "nullspace", (PetscObject *)&mnull)); 39366f2ef4bSMark Adams } 39466f2ef4bSMark Adams } 395b0d30dd6SMatthew G. Knepley } 39666f2ef4bSMark Adams 39766f2ef4bSMark Adams if (!mnull) { 398a2f3521dSMark F. Adams PetscInt bs, NN, MM; 3999566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(a_A, &bs)); 4009566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(a_A, &MM, &NN)); 40163a3b9bcSJacob Faibussowitsch PetscCheck(MM % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MM %" PetscInt_FMT " must be divisible by bs %" PetscInt_FMT, MM, bs); 4029566063dSJacob Faibussowitsch PetscCall(PCSetCoordinates_AGG(pc, bs, MM / bs, NULL)); 403806fa848SBarry Smith } else { 404b8cd405aSMark F. Adams PetscReal *nullvec; 405b8cd405aSMark F. Adams PetscBool has_const; 406b8cd405aSMark F. Adams PetscInt i, j, mlocal, nvec, bs; 407d5d25401SBarry Smith const Vec *vecs; 408d5d25401SBarry Smith const PetscScalar *v; 409b817416eSBarry Smith 4109566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(a_A, &mlocal, NULL)); 4119566063dSJacob Faibussowitsch PetscCall(MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs)); 412d19c4ebbSmarkadams4 for (i = 0; i < nvec; i++) { 413d19c4ebbSmarkadams4 PetscCall(VecGetLocalSize(vecs[i], &j)); 414d19c4ebbSmarkadams4 PetscCheck(j == mlocal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Attached null space vector size %" PetscInt_FMT " != matrix size %" PetscInt_FMT, j, mlocal); 415d19c4ebbSmarkadams4 } 416a0dea87dSMark Adams pc_gamg->data_sz = (nvec + !!has_const) * mlocal; 4179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((nvec + !!has_const) * mlocal, &nullvec)); 4189371c9d4SSatish Balay if (has_const) 4199371c9d4SSatish Balay for (i = 0; i < mlocal; i++) nullvec[i] = 1.0; 420b8cd405aSMark F. Adams for (i = 0; i < nvec; i++) { 4219566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(vecs[i], &v)); 422b8cd405aSMark F. Adams for (j = 0; j < mlocal; j++) nullvec[(i + !!has_const) * mlocal + j] = PetscRealPart(v[j]); 4239566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(vecs[i], &v)); 424b8cd405aSMark F. Adams } 425b8cd405aSMark F. Adams pc_gamg->data = nullvec; 426b8cd405aSMark F. Adams pc_gamg->data_cell_cols = (nvec + !!has_const); 4279566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(a_A, &bs)); 428b8cd405aSMark F. Adams pc_gamg->data_cell_rows = bs; 429b8cd405aSMark F. Adams } 4303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4312e68589bSMark F. Adams } 4322e68589bSMark F. Adams 4332e68589bSMark F. Adams /* 434bae903cbSmarkadams4 formProl0 - collect null space data for each aggregate, do QR, put R in coarse grid data and Q in P_0 4352e68589bSMark F. Adams 4362e68589bSMark F. Adams Input Parameter: 4372adfe9d3SBarry Smith . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices 4382adfe9d3SBarry Smith . bs - row block size 4390cbbd2e1SMark F. Adams . nSAvec - column bs of new P 4400cbbd2e1SMark F. Adams . my0crs - global index of start of locals 4412adfe9d3SBarry Smith . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec] 4422e68589bSMark F. Adams . data_in[data_stride*nSAvec] - local data on fine grid 4432e68589bSMark F. Adams . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist' 444bae903cbSmarkadams4 4452e68589bSMark F. Adams Output Parameter: 4462e68589bSMark F. Adams . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data 4472e68589bSMark F. Adams . a_Prol - prolongation operator 4482e68589bSMark F. Adams */ 449d71ae5a4SJacob 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) 450d71ae5a4SJacob Faibussowitsch { 451bd026e97SJed Brown PetscInt Istart, my0, Iend, nloc, clid, flid = 0, aggID, kk, jj, ii, mm, nSelected, minsz, nghosts, out_data_stride; 4523b4367a7SBarry Smith MPI_Comm comm; 4532e68589bSMark F. Adams PetscReal *out_data; 454539c167fSBarry Smith PetscCDIntNd *pos; 4551943db53SBarry Smith PCGAMGHashTable fgid_flid; 4560cbbd2e1SMark F. Adams 4572e68589bSMark F. Adams PetscFunctionBegin; 4589566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)a_Prol, &comm)); 4599566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(a_Prol, &Istart, &Iend)); 4609371c9d4SSatish Balay nloc = (Iend - Istart) / bs; 4619371c9d4SSatish Balay my0 = Istart / bs; 46263a3b9bcSJacob 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); 4630cbbd2e1SMark F. Adams Iend /= bs; 4640cbbd2e1SMark F. Adams nghosts = data_stride / bs - nloc; 4652e68589bSMark F. Adams 4669566063dSJacob Faibussowitsch PetscCall(PCGAMGHashTableCreate(2 * nghosts + 1, &fgid_flid)); 46748a46eb9SPierre Jolivet for (kk = 0; kk < nghosts; kk++) PetscCall(PCGAMGHashTableAdd(&fgid_flid, flid_fgid[nloc + kk], nloc + kk)); 4682e68589bSMark F. Adams 4690cbbd2e1SMark F. Adams /* count selected -- same as number of cols of P */ 4700cbbd2e1SMark F. Adams for (nSelected = mm = 0; mm < nloc; mm++) { 471e78576d6SMark F. Adams PetscBool ise; 4729566063dSJacob Faibussowitsch PetscCall(PetscCDEmptyAt(agg_llists, mm, &ise)); 473e78576d6SMark F. Adams if (!ise) nSelected++; 4740cbbd2e1SMark F. Adams } 4759566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(a_Prol, &ii, &jj)); 47663a3b9bcSJacob Faibussowitsch PetscCheck((ii / nSAvec) == my0crs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ii %" PetscInt_FMT " /nSAvec %" PetscInt_FMT " != my0crs %" PetscInt_FMT, ii, nSAvec, my0crs); 47763a3b9bcSJacob 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); 4780cbbd2e1SMark F. Adams 4792e68589bSMark F. Adams /* aloc space for coarse point data (output) */ 4800cbbd2e1SMark F. Adams out_data_stride = nSelected * nSAvec; 4812fa5cd67SKarl Rupp 4829566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(out_data_stride * nSAvec, &out_data)); 48333812677SSatish Balay for (ii = 0; ii < out_data_stride * nSAvec; ii++) out_data[ii] = PETSC_MAX_REAL; 4840cbbd2e1SMark F. Adams *a_data_out = out_data; /* output - stride nSelected*nSAvec */ 4852e68589bSMark F. Adams 4862e68589bSMark F. Adams /* find points and set prolongation */ 487c8b0795cSMark F. Adams minsz = 100; 4880cbbd2e1SMark F. Adams for (mm = clid = 0; mm < nloc; mm++) { 4899566063dSJacob Faibussowitsch PetscCall(PetscCDSizeAt(agg_llists, mm, &jj)); 490e78576d6SMark F. Adams if (jj > 0) { 4910cbbd2e1SMark F. Adams const PetscInt lid = mm, cgid = my0crs + clid; 4920cbbd2e1SMark F. Adams PetscInt cids[100]; /* max bs */ 4930cbbd2e1SMark F. Adams PetscBLASInt asz = jj, M = asz * bs, N = nSAvec, INFO; 4942e68589bSMark F. Adams PetscBLASInt Mdata = M + ((N - M > 0) ? N - M : 0), LDA = Mdata, LWORK = N * bs; 4952e68589bSMark F. Adams PetscScalar *qqc, *qqr, *TAU, *WORK; 4962e68589bSMark F. Adams PetscInt *fids; 49765d7b583SSatish Balay PetscReal *data; 498b817416eSBarry Smith 4990cbbd2e1SMark F. Adams /* count agg */ 5000cbbd2e1SMark F. Adams if (asz < minsz) minsz = asz; 5010cbbd2e1SMark F. Adams 5020cbbd2e1SMark F. Adams /* get block */ 5039566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(Mdata * N, &qqc, M * N, &qqr, N, &TAU, LWORK, &WORK, M, &fids)); 5042e68589bSMark F. Adams 5052e68589bSMark F. Adams aggID = 0; 5069566063dSJacob Faibussowitsch PetscCall(PetscCDGetHeadPos(agg_llists, lid, &pos)); 507e78576d6SMark F. Adams while (pos) { 508e78576d6SMark F. Adams PetscInt gid1; 5099566063dSJacob Faibussowitsch PetscCall(PetscCDIntNdGetID(pos, &gid1)); 5109566063dSJacob Faibussowitsch PetscCall(PetscCDGetNextPos(agg_llists, lid, &pos)); 511e78576d6SMark F. Adams 5120cbbd2e1SMark F. Adams if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0; 5130cbbd2e1SMark F. Adams else { 5149566063dSJacob Faibussowitsch PetscCall(PCGAMGHashTableFind(&fgid_flid, gid1, &flid)); 51508401ef6SPierre Jolivet PetscCheck(flid >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot find gid1 in table"); 5160cbbd2e1SMark F. Adams } 5172e68589bSMark F. Adams /* copy in B_i matrix - column oriented */ 51865d7b583SSatish Balay data = &data_in[flid * bs]; 5193d3eaba7SBarry Smith for (ii = 0; ii < bs; ii++) { 5202e68589bSMark F. Adams for (jj = 0; jj < N; jj++) { 5210cbbd2e1SMark F. Adams PetscReal d = data[jj * data_stride + ii]; 5220cbbd2e1SMark F. Adams qqc[jj * Mdata + aggID * bs + ii] = d; 5232e68589bSMark F. Adams } 5242e68589bSMark F. Adams } 5252e68589bSMark F. Adams /* set fine IDs */ 5262e68589bSMark F. Adams for (kk = 0; kk < bs; kk++) fids[aggID * bs + kk] = flid_fgid[flid] * bs + kk; 5272e68589bSMark F. Adams aggID++; 5280cbbd2e1SMark F. Adams } 5292e68589bSMark F. Adams 5302e68589bSMark F. Adams /* pad with zeros */ 5312e68589bSMark F. Adams for (ii = asz * bs; ii < Mdata; ii++) { 532ad540459SPierre Jolivet for (jj = 0; jj < N; jj++, kk++) qqc[jj * Mdata + ii] = .0; 5332e68589bSMark F. Adams } 5342e68589bSMark F. Adams 5352e68589bSMark F. Adams /* QR */ 5369566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 537792fecdfSBarry Smith PetscCallBLAS("LAPACKgeqrf", LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 5389566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPop()); 53908401ef6SPierre Jolivet PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xGEQRF error"); 5402e68589bSMark F. Adams /* get R - column oriented - output B_{i+1} */ 5412e68589bSMark F. Adams { 5422e68589bSMark F. Adams PetscReal *data = &out_data[clid * nSAvec]; 5432e68589bSMark F. Adams for (jj = 0; jj < nSAvec; jj++) { 5442e68589bSMark F. Adams for (ii = 0; ii < nSAvec; ii++) { 54508401ef6SPierre 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); 5460cbbd2e1SMark F. Adams if (ii <= jj) data[jj * out_data_stride + ii] = PetscRealPart(qqc[jj * Mdata + ii]); 5470cbbd2e1SMark F. Adams else data[jj * out_data_stride + ii] = 0.; 5482e68589bSMark F. Adams } 5492e68589bSMark F. Adams } 5502e68589bSMark F. Adams } 5512e68589bSMark F. Adams 5522e68589bSMark F. Adams /* get Q - row oriented */ 553792fecdfSBarry Smith PetscCallBLAS("LAPACKorgqr", LAPACKorgqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 55463a3b9bcSJacob Faibussowitsch PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xORGQR error arg %" PetscBLASInt_FMT, -INFO); 5552e68589bSMark F. Adams 5562e68589bSMark F. Adams for (ii = 0; ii < M; ii++) { 557ad540459SPierre Jolivet for (jj = 0; jj < N; jj++) qqr[N * ii + jj] = qqc[jj * Mdata + ii]; 5582e68589bSMark F. Adams } 5592e68589bSMark F. Adams 5602e68589bSMark F. Adams /* add diagonal block of P0 */ 5619371c9d4SSatish Balay for (kk = 0; kk < N; kk++) { cids[kk] = N * cgid + kk; /* global col IDs in P0 */ } 5629566063dSJacob Faibussowitsch PetscCall(MatSetValues(a_Prol, M, fids, N, cids, qqr, INSERT_VALUES)); 5639566063dSJacob Faibussowitsch PetscCall(PetscFree5(qqc, qqr, TAU, WORK, fids)); 564b43b03e9SMark F. Adams clid++; 5650cbbd2e1SMark F. Adams } /* coarse agg */ 5660cbbd2e1SMark F. Adams } /* for all fine nodes */ 5679566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(a_Prol, MAT_FINAL_ASSEMBLY)); 5689566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(a_Prol, MAT_FINAL_ASSEMBLY)); 5699566063dSJacob Faibussowitsch PetscCall(PCGAMGHashTableDestroy(&fgid_flid)); 5703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5712e68589bSMark F. Adams } 5722e68589bSMark F. Adams 573d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_GAMG_AGG(PC pc, PetscViewer viewer) 574d71ae5a4SJacob Faibussowitsch { 5755adeb434SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 5765adeb434SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 5775adeb434SBarry Smith PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 5785adeb434SBarry Smith 5795adeb434SBarry Smith PetscFunctionBegin; 5809566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " AGG specific options\n")); 581d529f056Smarkadams4 PetscCall(PetscViewerASCIIPrintf(viewer, " Number of levels of aggressive coarsening %d\n", (int)pc_gamg_agg->aggressive_coarsening_levels)); 582d529f056Smarkadams4 if (pc_gamg_agg->aggressive_coarsening_levels > 0) { 583d529f056Smarkadams4 PetscCall(PetscViewerASCIIPrintf(viewer, " %s aggressive coarsening\n", !pc_gamg_agg->use_aggressive_square_graph ? "MIS-k" : "Square graph")); 584d529f056Smarkadams4 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)); 585d529f056Smarkadams4 } 586bae903cbSmarkadams4 PetscCall(PetscViewerASCIIPrintf(viewer, " Number smoothing steps %d\n", (int)pc_gamg_agg->nsmooths)); 5873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5885adeb434SBarry Smith } 5895adeb434SBarry Smith 5902d776b49SBarry Smith static PetscErrorCode PCGAMGCreateGraph_AGG(PC pc, Mat Amat, Mat *a_Gmat) 591d71ae5a4SJacob Faibussowitsch { 592c8b0795cSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 593c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 594c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 5952d776b49SBarry Smith const PetscReal vfilter = pc_gamg->threshold[pc_gamg->current_level]; 5962d776b49SBarry Smith PetscBool ishem; 5972d776b49SBarry Smith const char *prefix; 598d529f056Smarkadams4 MatInfo info0, info1; 599d529f056Smarkadams4 PetscInt bs; 600c8b0795cSMark F. Adams 601c8b0795cSMark F. Adams PetscFunctionBegin; 602*a9ccf005SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0)); 6032d776b49SBarry 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 */ 6042d776b49SBarry Smith /* MATCOARSENHEM requires numerical weights for edges so ensure they are computed */ 6052d776b49SBarry Smith PetscCall(MatCoarsenCreate(PetscObjectComm((PetscObject)pc), &pc_gamg_agg->crs)); 6062d776b49SBarry Smith PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, &prefix)); 6072d776b49SBarry Smith PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc_gamg_agg->crs, prefix)); 6082d776b49SBarry Smith PetscCall(MatCoarsenSetFromOptions(pc_gamg_agg->crs)); 6092d776b49SBarry Smith PetscCall(PetscObjectTypeCompare((PetscObject)pc_gamg_agg->crs, MATCOARSENHEM, &ishem)); 610d529f056Smarkadams4 if (ishem) pc_gamg_agg->aggressive_coarsening_levels = 0; // aggressive and HEM does not make sense 611*a9ccf005SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0)); 612*a9ccf005SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0)); 613d529f056Smarkadams4 PetscCall(MatGetInfo(Amat, MAT_LOCAL, &info0)); /* global reduction */ 614*a9ccf005SMark Adams 615*a9ccf005SMark Adams if (ishem || pc_gamg_agg->use_low_mem_filter) { 6162d776b49SBarry Smith PetscCall(MatCreateGraph(Amat, PETSC_TRUE, (vfilter >= 0 || ishem) ? PETSC_TRUE : PETSC_FALSE, vfilter, a_Gmat)); 617*a9ccf005SMark Adams } else { 618*a9ccf005SMark Adams // make scalar graph, symetrize if not know to be symetric, scale, but do not filter (expensive) 619*a9ccf005SMark Adams PetscCall(MatCreateGraph(Amat, PETSC_TRUE, PETSC_TRUE, -1, a_Gmat)); 620*a9ccf005SMark Adams if (vfilter >= 0) { 621*a9ccf005SMark Adams PetscInt Istart, Iend, ncols, nnz0, nnz1, NN, MM, nloc; 622*a9ccf005SMark Adams Mat tGmat, Gmat = *a_Gmat; 623*a9ccf005SMark Adams MPI_Comm comm; 624*a9ccf005SMark Adams const PetscScalar *vals; 625*a9ccf005SMark Adams const PetscInt *idx; 626*a9ccf005SMark Adams PetscInt *d_nnz, *o_nnz, kk, *garray = NULL, *AJ, maxcols = 0; 627*a9ccf005SMark Adams MatScalar *AA; // this is checked in graph 628*a9ccf005SMark Adams PetscBool isseqaij; 629*a9ccf005SMark Adams Mat a, b, c; 630*a9ccf005SMark Adams MatType jtype; 631*a9ccf005SMark Adams 632*a9ccf005SMark Adams PetscCall(PetscObjectGetComm((PetscObject)Gmat, &comm)); 633*a9ccf005SMark Adams PetscCall(PetscObjectBaseTypeCompare((PetscObject)Gmat, MATSEQAIJ, &isseqaij)); 634*a9ccf005SMark Adams PetscCall(MatGetType(Gmat, &jtype)); 635*a9ccf005SMark Adams PetscCall(MatCreate(comm, &tGmat)); 636*a9ccf005SMark Adams PetscCall(MatSetType(tGmat, jtype)); 637*a9ccf005SMark Adams 638*a9ccf005SMark Adams /* TODO GPU: this can be called when filter = 0 -> Probably provide MatAIJThresholdCompress that compresses the entries below a threshold? 639*a9ccf005SMark Adams Also, if the matrix is symmetric, can we skip this 640*a9ccf005SMark Adams operation? It can be very expensive on large matrices. */ 641*a9ccf005SMark Adams 642*a9ccf005SMark Adams // global sizes 643*a9ccf005SMark Adams PetscCall(MatGetSize(Gmat, &MM, &NN)); 644*a9ccf005SMark Adams PetscCall(MatGetOwnershipRange(Gmat, &Istart, &Iend)); 645*a9ccf005SMark Adams nloc = Iend - Istart; 646*a9ccf005SMark Adams PetscCall(PetscMalloc2(nloc, &d_nnz, nloc, &o_nnz)); 647*a9ccf005SMark Adams if (isseqaij) { 648*a9ccf005SMark Adams a = Gmat; 649*a9ccf005SMark Adams b = NULL; 650*a9ccf005SMark Adams } else { 651*a9ccf005SMark Adams Mat_MPIAIJ *d = (Mat_MPIAIJ *)Gmat->data; 652*a9ccf005SMark Adams a = d->A; 653*a9ccf005SMark Adams b = d->B; 654*a9ccf005SMark Adams garray = d->garray; 655*a9ccf005SMark Adams } 656*a9ccf005SMark Adams /* Determine upper bound on non-zeros needed in new filtered matrix */ 657*a9ccf005SMark Adams for (PetscInt row = 0; row < nloc; row++) { 658*a9ccf005SMark Adams PetscCall(MatGetRow(a, row, &ncols, NULL, NULL)); 659*a9ccf005SMark Adams d_nnz[row] = ncols; 660*a9ccf005SMark Adams if (ncols > maxcols) maxcols = ncols; 661*a9ccf005SMark Adams PetscCall(MatRestoreRow(a, row, &ncols, NULL, NULL)); 662*a9ccf005SMark Adams } 663*a9ccf005SMark Adams if (b) { 664*a9ccf005SMark Adams for (PetscInt row = 0; row < nloc; row++) { 665*a9ccf005SMark Adams PetscCall(MatGetRow(b, row, &ncols, NULL, NULL)); 666*a9ccf005SMark Adams o_nnz[row] = ncols; 667*a9ccf005SMark Adams if (ncols > maxcols) maxcols = ncols; 668*a9ccf005SMark Adams PetscCall(MatRestoreRow(b, row, &ncols, NULL, NULL)); 669*a9ccf005SMark Adams } 670*a9ccf005SMark Adams } 671*a9ccf005SMark Adams PetscCall(MatSetSizes(tGmat, nloc, nloc, MM, MM)); 672*a9ccf005SMark Adams PetscCall(MatSetBlockSizes(tGmat, 1, 1)); 673*a9ccf005SMark Adams PetscCall(MatSeqAIJSetPreallocation(tGmat, 0, d_nnz)); 674*a9ccf005SMark Adams PetscCall(MatMPIAIJSetPreallocation(tGmat, 0, d_nnz, 0, o_nnz)); 675*a9ccf005SMark Adams PetscCall(MatSetOption(tGmat, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); 676*a9ccf005SMark Adams PetscCall(PetscFree2(d_nnz, o_nnz)); 677*a9ccf005SMark Adams PetscCall(PetscMalloc2(maxcols, &AA, maxcols, &AJ)); 678*a9ccf005SMark Adams nnz0 = nnz1 = 0; 679*a9ccf005SMark Adams for (c = a, kk = 0; c && kk < 2; c = b, kk++) { 680*a9ccf005SMark Adams for (PetscInt row = 0, grow = Istart, ncol_row, jj; row < nloc; row++, grow++) { 681*a9ccf005SMark Adams PetscCall(MatGetRow(c, row, &ncols, &idx, &vals)); 682*a9ccf005SMark Adams for (ncol_row = jj = 0; jj < ncols; jj++, nnz0++) { 683*a9ccf005SMark Adams PetscScalar sv = PetscAbs(PetscRealPart(vals[jj])); 684*a9ccf005SMark Adams if (PetscRealPart(sv) > vfilter) { 685*a9ccf005SMark Adams PetscInt cid = idx[jj] + Istart; //diag 686*a9ccf005SMark Adams nnz1++; 687*a9ccf005SMark Adams if (c != a) cid = garray[idx[jj]]; 688*a9ccf005SMark Adams AA[ncol_row] = vals[jj]; 689*a9ccf005SMark Adams AJ[ncol_row] = cid; 690*a9ccf005SMark Adams ncol_row++; 691*a9ccf005SMark Adams } 692*a9ccf005SMark Adams } 693*a9ccf005SMark Adams PetscCall(MatRestoreRow(c, row, &ncols, &idx, &vals)); 694*a9ccf005SMark Adams PetscCall(MatSetValues(tGmat, 1, &grow, ncol_row, AJ, AA, INSERT_VALUES)); 695*a9ccf005SMark Adams } 696*a9ccf005SMark Adams } 697*a9ccf005SMark Adams PetscCall(PetscFree2(AA, AJ)); 698*a9ccf005SMark Adams PetscCall(MatAssemblyBegin(tGmat, MAT_FINAL_ASSEMBLY)); 699*a9ccf005SMark Adams PetscCall(MatAssemblyEnd(tGmat, MAT_FINAL_ASSEMBLY)); 700*a9ccf005SMark Adams PetscCall(MatPropagateSymmetryOptions(Gmat, tGmat)); /* Normal Mat options are not relevant ? */ 701*a9ccf005SMark 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)); 702*a9ccf005SMark Adams PetscCall(MatViewFromOptions(tGmat, NULL, "-mat_filter_graph_view")); 703*a9ccf005SMark Adams PetscCall(MatDestroy(&Gmat)); 704*a9ccf005SMark Adams *a_Gmat = tGmat; 705*a9ccf005SMark Adams } 706*a9ccf005SMark Adams } 707*a9ccf005SMark Adams 708d529f056Smarkadams4 PetscCall(MatGetInfo(*a_Gmat, MAT_LOCAL, &info1)); /* global reduction */ 709d529f056Smarkadams4 PetscCall(MatGetBlockSize(Amat, &bs)); 710d529f056Smarkadams4 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)); 711849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0)); 7123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 713c8b0795cSMark F. Adams } 714c8b0795cSMark F. Adams 715d529f056Smarkadams4 typedef PetscInt NState; 716d529f056Smarkadams4 static const NState NOT_DONE = -2; 717d529f056Smarkadams4 static const NState DELETED = -1; 718d529f056Smarkadams4 static const NState REMOVED = -3; 719d529f056Smarkadams4 #define IS_SELECTED(s) (s != DELETED && s != NOT_DONE && s != REMOVED) 720d529f056Smarkadams4 721d529f056Smarkadams4 /* 722d529f056Smarkadams4 fixAggregatesWithSquare - greedy grab of with G1 (unsquared graph) -- AIJ specific -- change to fixAggregatesWithSquare -- TODD 723d529f056Smarkadams4 - AGG-MG specific: clears singletons out of 'selected_2' 724d529f056Smarkadams4 725d529f056Smarkadams4 Input Parameter: 726d529f056Smarkadams4 . Gmat_2 - global matrix of squared graph (data not defined) 727d529f056Smarkadams4 . Gmat_1 - base graph to grab with base graph 728d529f056Smarkadams4 Input/Output Parameter: 729d529f056Smarkadams4 . aggs_2 - linked list of aggs with gids) 730d529f056Smarkadams4 */ 731d529f056Smarkadams4 static PetscErrorCode fixAggregatesWithSquare(PC pc, Mat Gmat_2, Mat Gmat_1, PetscCoarsenData *aggs_2) 732d529f056Smarkadams4 { 733d529f056Smarkadams4 PetscBool isMPI; 734d529f056Smarkadams4 Mat_SeqAIJ *matA_1, *matB_1 = NULL; 735d529f056Smarkadams4 MPI_Comm comm; 736d529f056Smarkadams4 PetscInt lid, *ii, *idx, ix, Iend, my0, kk, n, j; 737d529f056Smarkadams4 Mat_MPIAIJ *mpimat_2 = NULL, *mpimat_1 = NULL; 738d529f056Smarkadams4 const PetscInt nloc = Gmat_2->rmap->n; 739d529f056Smarkadams4 PetscScalar *cpcol_1_state, *cpcol_2_state, *cpcol_2_par_orig, *lid_parent_gid; 740d529f056Smarkadams4 PetscInt *lid_cprowID_1 = NULL; 741d529f056Smarkadams4 NState *lid_state; 742d529f056Smarkadams4 Vec ghost_par_orig2; 743d529f056Smarkadams4 PetscMPIInt rank; 744d529f056Smarkadams4 745d529f056Smarkadams4 PetscFunctionBegin; 746d529f056Smarkadams4 PetscCall(PetscObjectGetComm((PetscObject)Gmat_2, &comm)); 747d529f056Smarkadams4 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 748d529f056Smarkadams4 PetscCall(MatGetOwnershipRange(Gmat_1, &my0, &Iend)); 749d529f056Smarkadams4 750d529f056Smarkadams4 /* get submatrices */ 751d529f056Smarkadams4 PetscCall(PetscStrbeginswith(((PetscObject)Gmat_1)->type_name, MATMPIAIJ, &isMPI)); 752d529f056Smarkadams4 PetscCall(PetscInfo(pc, "isMPI = %s\n", isMPI ? "yes" : "no")); 753d529f056Smarkadams4 PetscCall(PetscMalloc3(nloc, &lid_state, nloc, &lid_parent_gid, nloc, &lid_cprowID_1)); 754d529f056Smarkadams4 for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1; 755d529f056Smarkadams4 if (isMPI) { 756d529f056Smarkadams4 /* grab matrix objects */ 757d529f056Smarkadams4 mpimat_2 = (Mat_MPIAIJ *)Gmat_2->data; 758d529f056Smarkadams4 mpimat_1 = (Mat_MPIAIJ *)Gmat_1->data; 759d529f056Smarkadams4 matA_1 = (Mat_SeqAIJ *)mpimat_1->A->data; 760d529f056Smarkadams4 matB_1 = (Mat_SeqAIJ *)mpimat_1->B->data; 761d529f056Smarkadams4 762d529f056Smarkadams4 /* force compressed row storage for B matrix in AuxMat */ 763d529f056Smarkadams4 PetscCall(MatCheckCompressedRow(mpimat_1->B, matB_1->nonzerorowcnt, &matB_1->compressedrow, matB_1->i, Gmat_1->rmap->n, -1.0)); 764d529f056Smarkadams4 for (ix = 0; ix < matB_1->compressedrow.nrows; ix++) { 765d529f056Smarkadams4 PetscInt lid = matB_1->compressedrow.rindex[ix]; 766d529f056Smarkadams4 PetscCheck(lid <= nloc && lid >= -1, PETSC_COMM_SELF, PETSC_ERR_USER, "lid %d out of range. nloc = %d", (int)lid, (int)nloc); 767d529f056Smarkadams4 if (lid != -1) lid_cprowID_1[lid] = ix; 768d529f056Smarkadams4 } 769d529f056Smarkadams4 } else { 770d529f056Smarkadams4 PetscBool isAIJ; 771d529f056Smarkadams4 PetscCall(PetscStrbeginswith(((PetscObject)Gmat_1)->type_name, MATSEQAIJ, &isAIJ)); 772d529f056Smarkadams4 PetscCheck(isAIJ, PETSC_COMM_SELF, PETSC_ERR_USER, "Require AIJ matrix."); 773d529f056Smarkadams4 matA_1 = (Mat_SeqAIJ *)Gmat_1->data; 774d529f056Smarkadams4 } 775d529f056Smarkadams4 if (nloc > 0) { PetscCheck(!matB_1 || matB_1->compressedrow.use, PETSC_COMM_SELF, PETSC_ERR_PLIB, "matB_1 && !matB_1->compressedrow.use: PETSc bug???"); } 776d529f056Smarkadams4 /* get state of locals and selected gid for deleted */ 777d529f056Smarkadams4 for (lid = 0; lid < nloc; lid++) { 778d529f056Smarkadams4 lid_parent_gid[lid] = -1.0; 779d529f056Smarkadams4 lid_state[lid] = DELETED; 780d529f056Smarkadams4 } 781d529f056Smarkadams4 782d529f056Smarkadams4 /* set lid_state */ 783d529f056Smarkadams4 for (lid = 0; lid < nloc; lid++) { 784d529f056Smarkadams4 PetscCDIntNd *pos; 785d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos)); 786d529f056Smarkadams4 if (pos) { 787d529f056Smarkadams4 PetscInt gid1; 788d529f056Smarkadams4 789d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gid1)); 790d529f056Smarkadams4 PetscCheck(gid1 == lid + my0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "gid1 %d != lid %d + my0 %d", (int)gid1, (int)lid, (int)my0); 791d529f056Smarkadams4 lid_state[lid] = gid1; 792d529f056Smarkadams4 } 793d529f056Smarkadams4 } 794d529f056Smarkadams4 795d529f056Smarkadams4 /* map local to selected local, DELETED means a ghost owns it */ 796d529f056Smarkadams4 for (lid = kk = 0; lid < nloc; lid++) { 797d529f056Smarkadams4 NState state = lid_state[lid]; 798d529f056Smarkadams4 if (IS_SELECTED(state)) { 799d529f056Smarkadams4 PetscCDIntNd *pos; 800d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos)); 801d529f056Smarkadams4 while (pos) { 802d529f056Smarkadams4 PetscInt gid1; 803d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gid1)); 804d529f056Smarkadams4 PetscCall(PetscCDGetNextPos(aggs_2, lid, &pos)); 805d529f056Smarkadams4 if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1 - my0] = (PetscScalar)(lid + my0); 806d529f056Smarkadams4 } 807d529f056Smarkadams4 } 808d529f056Smarkadams4 } 809d529f056Smarkadams4 /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */ 810d529f056Smarkadams4 if (isMPI) { 811d529f056Smarkadams4 Vec tempVec; 812d529f056Smarkadams4 /* get 'cpcol_1_state' */ 813d529f056Smarkadams4 PetscCall(MatCreateVecs(Gmat_1, &tempVec, NULL)); 814d529f056Smarkadams4 for (kk = 0, j = my0; kk < nloc; kk++, j++) { 815d529f056Smarkadams4 PetscScalar v = (PetscScalar)lid_state[kk]; 816d529f056Smarkadams4 PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES)); 817d529f056Smarkadams4 } 818d529f056Smarkadams4 PetscCall(VecAssemblyBegin(tempVec)); 819d529f056Smarkadams4 PetscCall(VecAssemblyEnd(tempVec)); 820d529f056Smarkadams4 PetscCall(VecScatterBegin(mpimat_1->Mvctx, tempVec, mpimat_1->lvec, INSERT_VALUES, SCATTER_FORWARD)); 821d529f056Smarkadams4 PetscCall(VecScatterEnd(mpimat_1->Mvctx, tempVec, mpimat_1->lvec, INSERT_VALUES, SCATTER_FORWARD)); 822d529f056Smarkadams4 PetscCall(VecGetArray(mpimat_1->lvec, &cpcol_1_state)); 823d529f056Smarkadams4 /* get 'cpcol_2_state' */ 824d529f056Smarkadams4 PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, mpimat_2->lvec, INSERT_VALUES, SCATTER_FORWARD)); 825d529f056Smarkadams4 PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, mpimat_2->lvec, INSERT_VALUES, SCATTER_FORWARD)); 826d529f056Smarkadams4 PetscCall(VecGetArray(mpimat_2->lvec, &cpcol_2_state)); 827d529f056Smarkadams4 /* get 'cpcol_2_par_orig' */ 828d529f056Smarkadams4 for (kk = 0, j = my0; kk < nloc; kk++, j++) { 829d529f056Smarkadams4 PetscScalar v = (PetscScalar)lid_parent_gid[kk]; 830d529f056Smarkadams4 PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES)); 831d529f056Smarkadams4 } 832d529f056Smarkadams4 PetscCall(VecAssemblyBegin(tempVec)); 833d529f056Smarkadams4 PetscCall(VecAssemblyEnd(tempVec)); 834d529f056Smarkadams4 PetscCall(VecDuplicate(mpimat_2->lvec, &ghost_par_orig2)); 835d529f056Smarkadams4 PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghost_par_orig2, INSERT_VALUES, SCATTER_FORWARD)); 836d529f056Smarkadams4 PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghost_par_orig2, INSERT_VALUES, SCATTER_FORWARD)); 837d529f056Smarkadams4 PetscCall(VecGetArray(ghost_par_orig2, &cpcol_2_par_orig)); 838d529f056Smarkadams4 839d529f056Smarkadams4 PetscCall(VecDestroy(&tempVec)); 840d529f056Smarkadams4 } /* ismpi */ 841d529f056Smarkadams4 for (lid = 0; lid < nloc; lid++) { 842d529f056Smarkadams4 NState state = lid_state[lid]; 843d529f056Smarkadams4 if (IS_SELECTED(state)) { 844d529f056Smarkadams4 /* steal locals */ 845d529f056Smarkadams4 ii = matA_1->i; 846d529f056Smarkadams4 n = ii[lid + 1] - ii[lid]; 847d529f056Smarkadams4 idx = matA_1->j + ii[lid]; 848d529f056Smarkadams4 for (j = 0; j < n; j++) { 849d529f056Smarkadams4 PetscInt lidj = idx[j], sgid; 850d529f056Smarkadams4 NState statej = lid_state[lidj]; 851d529f056Smarkadams4 if (statej == DELETED && (sgid = (PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid + my0) { /* steal local */ 852d529f056Smarkadams4 lid_parent_gid[lidj] = (PetscScalar)(lid + my0); /* send this if sgid is not local */ 853d529f056Smarkadams4 if (sgid >= my0 && sgid < Iend) { /* I'm stealing this local from a local sgid */ 854d529f056Smarkadams4 PetscInt hav = 0, slid = sgid - my0, gidj = lidj + my0; 855d529f056Smarkadams4 PetscCDIntNd *pos, *last = NULL; 856d529f056Smarkadams4 /* looking for local from local so id_llist_2 works */ 857d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, slid, &pos)); 858d529f056Smarkadams4 while (pos) { 859d529f056Smarkadams4 PetscInt gid; 860d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gid)); 861d529f056Smarkadams4 if (gid == gidj) { 862d529f056Smarkadams4 PetscCheck(last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "last cannot be null"); 863d529f056Smarkadams4 PetscCall(PetscCDRemoveNextNode(aggs_2, slid, last)); 864d529f056Smarkadams4 PetscCall(PetscCDAppendNode(aggs_2, lid, pos)); 865d529f056Smarkadams4 hav = 1; 866d529f056Smarkadams4 break; 867d529f056Smarkadams4 } else last = pos; 868d529f056Smarkadams4 PetscCall(PetscCDGetNextPos(aggs_2, slid, &pos)); 869d529f056Smarkadams4 } 870d529f056Smarkadams4 if (hav != 1) { 871d529f056Smarkadams4 PetscCheck(hav, PETSC_COMM_SELF, PETSC_ERR_PLIB, "failed to find adj in 'selected' lists - structurally unsymmetric matrix"); 872d529f056Smarkadams4 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "found node %d times???", (int)hav); 873d529f056Smarkadams4 } 874d529f056Smarkadams4 } else { /* I'm stealing this local, owned by a ghost */ 875d529f056Smarkadams4 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.", 876d529f056Smarkadams4 ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : ""); 877d529f056Smarkadams4 PetscCall(PetscCDAppendID(aggs_2, lid, lidj + my0)); 878d529f056Smarkadams4 } 879d529f056Smarkadams4 } 880d529f056Smarkadams4 } /* local neighbors */ 881d529f056Smarkadams4 } else if (state == DELETED /* && lid_cprowID_1 */) { 882d529f056Smarkadams4 PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]); 883d529f056Smarkadams4 /* see if I have a selected ghost neighbor that will steal me */ 884d529f056Smarkadams4 if ((ix = lid_cprowID_1[lid]) != -1) { 885d529f056Smarkadams4 ii = matB_1->compressedrow.i; 886d529f056Smarkadams4 n = ii[ix + 1] - ii[ix]; 887d529f056Smarkadams4 idx = matB_1->j + ii[ix]; 888d529f056Smarkadams4 for (j = 0; j < n; j++) { 889d529f056Smarkadams4 PetscInt cpid = idx[j]; 890d529f056Smarkadams4 NState statej = (NState)PetscRealPart(cpcol_1_state[cpid]); 891d529f056Smarkadams4 if (IS_SELECTED(statej) && sgidold != (PetscInt)statej) { /* ghost will steal this, remove from my list */ 892d529f056Smarkadams4 lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */ 893d529f056Smarkadams4 if (sgidold >= my0 && sgidold < Iend) { /* this was mine */ 894d529f056Smarkadams4 PetscInt hav = 0, oldslidj = sgidold - my0; 895d529f056Smarkadams4 PetscCDIntNd *pos, *last = NULL; 896d529f056Smarkadams4 /* remove from 'oldslidj' list */ 897d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, oldslidj, &pos)); 898d529f056Smarkadams4 while (pos) { 899d529f056Smarkadams4 PetscInt gid; 900d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gid)); 901d529f056Smarkadams4 if (lid + my0 == gid) { 902d529f056Smarkadams4 /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove lid from oldslidj list *\/ */ 903d529f056Smarkadams4 PetscCheck(last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "last cannot be null"); 904d529f056Smarkadams4 PetscCall(PetscCDRemoveNextNode(aggs_2, oldslidj, last)); 905d529f056Smarkadams4 /* ghost (PetscScalar)statej will add this later */ 906d529f056Smarkadams4 hav = 1; 907d529f056Smarkadams4 break; 908d529f056Smarkadams4 } else last = pos; 909d529f056Smarkadams4 PetscCall(PetscCDGetNextPos(aggs_2, oldslidj, &pos)); 910d529f056Smarkadams4 } 911d529f056Smarkadams4 if (hav != 1) { 912d529f056Smarkadams4 PetscCheck(hav, PETSC_COMM_SELF, PETSC_ERR_PLIB, "failed to find (hav=%d) adj in 'selected' lists - structurally unsymmetric matrix", (int)hav); 913d529f056Smarkadams4 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "found node %d times???", (int)hav); 914d529f056Smarkadams4 } 915d529f056Smarkadams4 } else { 916d529f056Smarkadams4 /* TODO: ghosts remove this later */ 917d529f056Smarkadams4 } 918d529f056Smarkadams4 } 919d529f056Smarkadams4 } 920d529f056Smarkadams4 } 921d529f056Smarkadams4 } /* selected/deleted */ 922d529f056Smarkadams4 } /* node loop */ 923d529f056Smarkadams4 924d529f056Smarkadams4 if (isMPI) { 925d529f056Smarkadams4 PetscScalar *cpcol_2_parent, *cpcol_2_gid; 926d529f056Smarkadams4 Vec tempVec, ghostgids2, ghostparents2; 927d529f056Smarkadams4 PetscInt cpid, nghost_2; 928d529f056Smarkadams4 PCGAMGHashTable gid_cpid; 929d529f056Smarkadams4 930d529f056Smarkadams4 PetscCall(VecGetSize(mpimat_2->lvec, &nghost_2)); 931d529f056Smarkadams4 PetscCall(MatCreateVecs(Gmat_2, &tempVec, NULL)); 932d529f056Smarkadams4 933d529f056Smarkadams4 /* get 'cpcol_2_parent' */ 934d529f056Smarkadams4 for (kk = 0, j = my0; kk < nloc; kk++, j++) { PetscCall(VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES)); } 935d529f056Smarkadams4 PetscCall(VecAssemblyBegin(tempVec)); 936d529f056Smarkadams4 PetscCall(VecAssemblyEnd(tempVec)); 937d529f056Smarkadams4 PetscCall(VecDuplicate(mpimat_2->lvec, &ghostparents2)); 938d529f056Smarkadams4 PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghostparents2, INSERT_VALUES, SCATTER_FORWARD)); 939d529f056Smarkadams4 PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghostparents2, INSERT_VALUES, SCATTER_FORWARD)); 940d529f056Smarkadams4 PetscCall(VecGetArray(ghostparents2, &cpcol_2_parent)); 941d529f056Smarkadams4 942d529f056Smarkadams4 /* get 'cpcol_2_gid' */ 943d529f056Smarkadams4 for (kk = 0, j = my0; kk < nloc; kk++, j++) { 944d529f056Smarkadams4 PetscScalar v = (PetscScalar)j; 945d529f056Smarkadams4 PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES)); 946d529f056Smarkadams4 } 947d529f056Smarkadams4 PetscCall(VecAssemblyBegin(tempVec)); 948d529f056Smarkadams4 PetscCall(VecAssemblyEnd(tempVec)); 949d529f056Smarkadams4 PetscCall(VecDuplicate(mpimat_2->lvec, &ghostgids2)); 950d529f056Smarkadams4 PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghostgids2, INSERT_VALUES, SCATTER_FORWARD)); 951d529f056Smarkadams4 PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghostgids2, INSERT_VALUES, SCATTER_FORWARD)); 952d529f056Smarkadams4 PetscCall(VecGetArray(ghostgids2, &cpcol_2_gid)); 953d529f056Smarkadams4 PetscCall(VecDestroy(&tempVec)); 954d529f056Smarkadams4 955d529f056Smarkadams4 /* look for deleted ghosts and add to table */ 956d529f056Smarkadams4 PetscCall(PCGAMGHashTableCreate(2 * nghost_2 + 1, &gid_cpid)); 957d529f056Smarkadams4 for (cpid = 0; cpid < nghost_2; cpid++) { 958d529f056Smarkadams4 NState state = (NState)PetscRealPart(cpcol_2_state[cpid]); 959d529f056Smarkadams4 if (state == DELETED) { 960d529f056Smarkadams4 PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 961d529f056Smarkadams4 PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]); 962d529f056Smarkadams4 if (sgid_old == -1 && sgid_new != -1) { 963d529f056Smarkadams4 PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 964d529f056Smarkadams4 PetscCall(PCGAMGHashTableAdd(&gid_cpid, gid, cpid)); 965d529f056Smarkadams4 } 966d529f056Smarkadams4 } 967d529f056Smarkadams4 } 968d529f056Smarkadams4 969d529f056Smarkadams4 /* look for deleted ghosts and see if they moved - remove it */ 970d529f056Smarkadams4 for (lid = 0; lid < nloc; lid++) { 971d529f056Smarkadams4 NState state = lid_state[lid]; 972d529f056Smarkadams4 if (IS_SELECTED(state)) { 973d529f056Smarkadams4 PetscCDIntNd *pos, *last = NULL; 974d529f056Smarkadams4 /* look for deleted ghosts and see if they moved */ 975d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos)); 976d529f056Smarkadams4 while (pos) { 977d529f056Smarkadams4 PetscInt gid; 978d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gid)); 979d529f056Smarkadams4 980d529f056Smarkadams4 if (gid < my0 || gid >= Iend) { 981d529f056Smarkadams4 PetscCall(PCGAMGHashTableFind(&gid_cpid, gid, &cpid)); 982d529f056Smarkadams4 if (cpid != -1) { 983d529f056Smarkadams4 /* a moved ghost - */ 984d529f056Smarkadams4 /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove 'flid' from list *\/ */ 985d529f056Smarkadams4 PetscCall(PetscCDRemoveNextNode(aggs_2, lid, last)); 986d529f056Smarkadams4 } else last = pos; 987d529f056Smarkadams4 } else last = pos; 988d529f056Smarkadams4 989d529f056Smarkadams4 PetscCall(PetscCDGetNextPos(aggs_2, lid, &pos)); 990d529f056Smarkadams4 } /* loop over list of deleted */ 991d529f056Smarkadams4 } /* selected */ 992d529f056Smarkadams4 } 993d529f056Smarkadams4 PetscCall(PCGAMGHashTableDestroy(&gid_cpid)); 994d529f056Smarkadams4 995d529f056Smarkadams4 /* look at ghosts, see if they changed - and it */ 996d529f056Smarkadams4 for (cpid = 0; cpid < nghost_2; cpid++) { 997d529f056Smarkadams4 PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 998d529f056Smarkadams4 if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */ 999d529f056Smarkadams4 PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 1000d529f056Smarkadams4 PetscInt slid_new = sgid_new - my0, hav = 0; 1001d529f056Smarkadams4 PetscCDIntNd *pos; 1002d529f056Smarkadams4 1003d529f056Smarkadams4 /* search for this gid to see if I have it */ 1004d529f056Smarkadams4 PetscCall(PetscCDGetHeadPos(aggs_2, slid_new, &pos)); 1005d529f056Smarkadams4 while (pos) { 1006d529f056Smarkadams4 PetscInt gidj; 1007d529f056Smarkadams4 PetscCall(PetscCDIntNdGetID(pos, &gidj)); 1008d529f056Smarkadams4 PetscCall(PetscCDGetNextPos(aggs_2, slid_new, &pos)); 1009d529f056Smarkadams4 1010d529f056Smarkadams4 if (gidj == gid) { 1011d529f056Smarkadams4 hav = 1; 1012d529f056Smarkadams4 break; 1013d529f056Smarkadams4 } 1014d529f056Smarkadams4 } 1015d529f056Smarkadams4 if (hav != 1) { 1016d529f056Smarkadams4 /* insert 'flidj' into head of llist */ 1017d529f056Smarkadams4 PetscCall(PetscCDAppendID(aggs_2, slid_new, gid)); 1018d529f056Smarkadams4 } 1019d529f056Smarkadams4 } 1020d529f056Smarkadams4 } 1021d529f056Smarkadams4 PetscCall(VecRestoreArray(mpimat_1->lvec, &cpcol_1_state)); 1022d529f056Smarkadams4 PetscCall(VecRestoreArray(mpimat_2->lvec, &cpcol_2_state)); 1023d529f056Smarkadams4 PetscCall(VecRestoreArray(ghostparents2, &cpcol_2_parent)); 1024d529f056Smarkadams4 PetscCall(VecRestoreArray(ghostgids2, &cpcol_2_gid)); 1025d529f056Smarkadams4 PetscCall(VecDestroy(&ghostgids2)); 1026d529f056Smarkadams4 PetscCall(VecDestroy(&ghostparents2)); 1027d529f056Smarkadams4 PetscCall(VecDestroy(&ghost_par_orig2)); 1028d529f056Smarkadams4 } 1029d529f056Smarkadams4 PetscCall(PetscFree3(lid_state, lid_parent_gid, lid_cprowID_1)); 1030d529f056Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS); 1031d529f056Smarkadams4 } 1032d529f056Smarkadams4 1033c8b0795cSMark F. Adams /* 1034bae903cbSmarkadams4 PCGAMGCoarsen_AGG - supports squaring the graph (deprecated) and new graph for 1035bae903cbSmarkadams4 communication of QR data used with HEM and MISk coarsening 1036c8b0795cSMark F. Adams 1037c8b0795cSMark F. Adams Input Parameter: 1038e0940f08SMark F. Adams . a_pc - this 1039bae903cbSmarkadams4 1040e0940f08SMark F. Adams Input/Output Parameter: 1041bae903cbSmarkadams4 . a_Gmat1 - graph to coarsen (in), graph off processor edges for QR gather scatter (out) 1042bae903cbSmarkadams4 1043c8b0795cSMark F. Adams Output Parameter: 10440cbbd2e1SMark F. Adams . agg_lists - list of aggregates 1045bae903cbSmarkadams4 1046c8b0795cSMark F. Adams */ 1047d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc, Mat *a_Gmat1, PetscCoarsenData **agg_lists) 1048d71ae5a4SJacob Faibussowitsch { 1049e0940f08SMark F. Adams PC_MG *mg = (PC_MG *)a_pc->data; 1050c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1051c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 1052d529f056Smarkadams4 Mat mat, Gmat2, Gmat1 = *a_Gmat1; /* aggressive graph */ 10530cbbd2e1SMark F. Adams IS perm; 1054bae903cbSmarkadams4 PetscInt Istart, Iend, Ii, nloc, bs, nn; 1055bae903cbSmarkadams4 PetscInt *permute, *degree; 1056c8b0795cSMark F. Adams PetscBool *bIndexSet; 1057aea53286SMark Adams PetscReal hashfact; 1058e2c00dcbSBarry Smith PetscInt iSwapIndex; 10593b50add6SMark Adams PetscRandom random; 1060d529f056Smarkadams4 MPI_Comm comm; 1061c8b0795cSMark F. Adams 1062c8b0795cSMark F. Adams PetscFunctionBegin; 1063d529f056Smarkadams4 PetscCall(PetscObjectGetComm((PetscObject)Gmat1, &comm)); 1064849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0)); 1065bae903cbSmarkadams4 PetscCall(MatGetLocalSize(Gmat1, &nn, NULL)); 10669566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Gmat1, &bs)); 106763a3b9bcSJacob Faibussowitsch PetscCheck(bs == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bs %" PetscInt_FMT " must be 1", bs); 1068bae903cbSmarkadams4 nloc = nn / bs; 10695cfd4bd9SMark Adams /* get MIS aggs - randomize */ 1070bae903cbSmarkadams4 PetscCall(PetscMalloc2(nloc, &permute, nloc, °ree)); 10719566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nloc, &bIndexSet)); 10726e09b0e3SStefano Zampini for (Ii = 0; Ii < nloc; Ii++) permute[Ii] = Ii; 10739566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &random)); 10749566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Gmat1, &Istart, &Iend)); 1075c8b0795cSMark F. Adams for (Ii = 0; Ii < nloc; Ii++) { 1076bae903cbSmarkadams4 PetscInt nc; 1077bae903cbSmarkadams4 PetscCall(MatGetRow(Gmat1, Istart + Ii, &nc, NULL, NULL)); 1078bae903cbSmarkadams4 degree[Ii] = nc; 1079bae903cbSmarkadams4 PetscCall(MatRestoreRow(Gmat1, Istart + Ii, &nc, NULL, NULL)); 1080bae903cbSmarkadams4 } 1081bae903cbSmarkadams4 for (Ii = 0; Ii < nloc; Ii++) { 10829566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(random, &hashfact)); 1083aea53286SMark Adams iSwapIndex = (PetscInt)(hashfact * nloc) % nloc; 1084c8b0795cSMark F. Adams if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) { 1085c8b0795cSMark F. Adams PetscInt iTemp = permute[iSwapIndex]; 1086c8b0795cSMark F. Adams permute[iSwapIndex] = permute[Ii]; 1087c8b0795cSMark F. Adams permute[Ii] = iTemp; 1088bae903cbSmarkadams4 iTemp = degree[iSwapIndex]; 1089bae903cbSmarkadams4 degree[iSwapIndex] = degree[Ii]; 1090bae903cbSmarkadams4 degree[Ii] = iTemp; 1091c8b0795cSMark F. Adams bIndexSet[iSwapIndex] = PETSC_TRUE; 1092c8b0795cSMark F. Adams } 1093c8b0795cSMark F. Adams } 1094d529f056Smarkadams4 // apply minimum degree ordering -- NEW 1095d529f056Smarkadams4 if (pc_gamg_agg->use_minimum_degree_ordering) { PetscCall(PetscSortIntWithArray(nloc, degree, permute)); } 10969566063dSJacob Faibussowitsch PetscCall(PetscFree(bIndexSet)); 10979566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&random)); 10989566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm)); 1099849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0)); 1100d529f056Smarkadams4 // square graph 1101d529f056Smarkadams4 if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels && pc_gamg_agg->use_aggressive_square_graph) { 1102d529f056Smarkadams4 PetscCall(PCGAMGSquareGraph_GAMG(a_pc, Gmat1, &Gmat2)); 1103d529f056Smarkadams4 } else Gmat2 = Gmat1; 1104d529f056Smarkadams4 // switch to old MIS-1 for square graph 1105d529f056Smarkadams4 if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels) { 1106d529f056Smarkadams4 if (!pc_gamg_agg->use_aggressive_square_graph) PetscCall(MatCoarsenMISKSetDistance(pc_gamg_agg->crs, pc_gamg_agg->aggressive_mis_k)); // hardwire to MIS-2 1107d529f056Smarkadams4 else PetscCall(MatCoarsenSetType(pc_gamg_agg->crs, MATCOARSENMIS)); // old MIS -- side effect 1108d529f056Smarkadams4 } else if (pc_gamg_agg->use_aggressive_square_graph && pc_gamg_agg->aggressive_coarsening_levels > 0) { // we reset the MIS 1109d529f056Smarkadams4 const char *prefix; 1110d529f056Smarkadams4 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)a_pc, &prefix)); 1111d529f056Smarkadams4 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc_gamg_agg->crs, prefix)); 1112d529f056Smarkadams4 PetscCall(MatCoarsenSetFromOptions(pc_gamg_agg->crs)); // get the default back on non-aggressive levels when square graph switched to old MIS 1113d529f056Smarkadams4 } 1114d529f056Smarkadams4 PetscCall(MatCoarsenSetAdjacency(pc_gamg_agg->crs, Gmat2)); 11152d776b49SBarry Smith PetscCall(MatCoarsenSetStrictAggs(pc_gamg_agg->crs, PETSC_TRUE)); 1116d529f056Smarkadams4 PetscCall(MatCoarsenSetGreedyOrdering(pc_gamg_agg->crs, perm)); 11172d776b49SBarry Smith PetscCall(MatCoarsenApply(pc_gamg_agg->crs)); 11182d776b49SBarry Smith PetscCall(MatCoarsenViewFromOptions(pc_gamg_agg->crs, NULL, "-mat_coarsen_view")); 11192d776b49SBarry Smith PetscCall(MatCoarsenGetData(pc_gamg_agg->crs, agg_lists)); /* output */ 11202d776b49SBarry Smith PetscCall(MatCoarsenDestroy(&pc_gamg_agg->crs)); 1121b43b03e9SMark F. Adams 11229566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 1123bae903cbSmarkadams4 PetscCall(PetscFree2(permute, degree)); 1124849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0)); 11259f3f12c8SMark F. Adams 1126d529f056Smarkadams4 if (Gmat2 != Gmat1) { 1127d529f056Smarkadams4 PetscCoarsenData *llist = *agg_lists; 1128d529f056Smarkadams4 PetscCall(fixAggregatesWithSquare(a_pc, Gmat2, Gmat1, *agg_lists)); 1129d529f056Smarkadams4 PetscCall(MatDestroy(&Gmat1)); 1130d529f056Smarkadams4 *a_Gmat1 = Gmat2; /* output */ 1131d529f056Smarkadams4 PetscCall(PetscCDGetMat(llist, &mat)); 1132baca6076SPierre Jolivet PetscCheck(!mat, comm, PETSC_ERR_ARG_WRONG, "Unexpected auxiliary matrix with squared graph"); 1133d529f056Smarkadams4 } else { 1134bae903cbSmarkadams4 PetscCoarsenData *llist = *agg_lists; 11359ab59c8bSMark Adams /* see if we have a matrix that takes precedence (returned from MatCoarsenApply) */ 11369566063dSJacob Faibussowitsch PetscCall(PetscCDGetMat(llist, &mat)); 11370cbbd2e1SMark F. Adams if (mat) { 1138d529f056Smarkadams4 PetscCall(MatDestroy(a_Gmat1)); 11390cbbd2e1SMark F. Adams *a_Gmat1 = mat; /* output */ 11400cbbd2e1SMark F. Adams } 11410cbbd2e1SMark F. Adams } 1142849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0)); 11433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1144c8b0795cSMark F. Adams } 1145c8b0795cSMark F. Adams 1146c8b0795cSMark F. Adams /* 11470cbbd2e1SMark F. Adams PCGAMGProlongator_AGG 1148c8b0795cSMark F. Adams 1149c8b0795cSMark F. Adams Input Parameter: 1150c8b0795cSMark F. Adams . pc - this 1151c8b0795cSMark F. Adams . Amat - matrix on this fine level 1152c8b0795cSMark F. Adams . Graph - used to get ghost data for nodes in 11530cbbd2e1SMark F. Adams . agg_lists - list of aggregates 1154c8b0795cSMark F. Adams Output Parameter: 1155c8b0795cSMark F. Adams . a_P_out - prolongation operator to the next level 1156c8b0795cSMark F. Adams */ 1157d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGProlongator_AGG(PC pc, Mat Amat, Mat Gmat, PetscCoarsenData *agg_lists, Mat *a_P_out) 1158d71ae5a4SJacob Faibussowitsch { 11592e68589bSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 11602e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 11614a2f8832SBarry Smith const PetscInt col_bs = pc_gamg->data_cell_cols; 1162c8b0795cSMark F. Adams PetscInt Istart, Iend, nloc, ii, jj, kk, my0, nLocalSelected, bs; 1163c8b0795cSMark F. Adams Mat Prol; 1164d5d25401SBarry Smith PetscMPIInt size; 11653b4367a7SBarry Smith MPI_Comm comm; 1166c8b0795cSMark F. Adams PetscReal *data_w_ghost; 1167c8b0795cSMark F. Adams PetscInt myCrs0, nbnodes = 0, *flid_fgid; 1168d9558ea9SBarry Smith MatType mtype; 11692e68589bSMark F. Adams 11702e68589bSMark F. Adams PetscFunctionBegin; 11719566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm)); 117208401ef6SPierre Jolivet PetscCheck(col_bs >= 1, comm, PETSC_ERR_PLIB, "Column bs cannot be less than 1"); 1173849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0)); 11749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 11759566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Amat, &Istart, &Iend)); 11769566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(Amat, &bs)); 11779371c9d4SSatish Balay nloc = (Iend - Istart) / bs; 11789371c9d4SSatish Balay my0 = Istart / bs; 117963a3b9bcSJacob 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); 11802e68589bSMark F. Adams 11812e68589bSMark F. Adams /* get 'nLocalSelected' */ 11820cbbd2e1SMark F. Adams for (ii = 0, nLocalSelected = 0; ii < nloc; ii++) { 1183e78576d6SMark F. Adams PetscBool ise; 11840cbbd2e1SMark F. Adams /* filter out singletons 0 or 1? */ 11859566063dSJacob Faibussowitsch PetscCall(PetscCDEmptyAt(agg_lists, ii, &ise)); 1186e78576d6SMark F. Adams if (!ise) nLocalSelected++; 11872e68589bSMark F. Adams } 11882e68589bSMark F. Adams 11892e68589bSMark F. Adams /* create prolongator, create P matrix */ 11909566063dSJacob Faibussowitsch PetscCall(MatGetType(Amat, &mtype)); 11919566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &Prol)); 11929566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Prol, nloc * bs, nLocalSelected * col_bs, PETSC_DETERMINE, PETSC_DETERMINE)); 11939566063dSJacob Faibussowitsch PetscCall(MatSetBlockSizes(Prol, bs, col_bs)); 11949566063dSJacob Faibussowitsch PetscCall(MatSetType(Prol, mtype)); 11953742c8caSstefanozampini #if PetscDefined(HAVE_DEVICE) 11963742c8caSstefanozampini PetscBool flg; 11973742c8caSstefanozampini PetscCall(MatBoundToCPU(Amat, &flg)); 11983742c8caSstefanozampini PetscCall(MatBindToCPU(Prol, flg)); 11993742c8caSstefanozampini if (flg) PetscCall(MatSetBindingPropagates(Prol, PETSC_TRUE)); 12003742c8caSstefanozampini #endif 12019566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(Prol, col_bs, NULL)); 12029566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(Prol, col_bs, NULL, col_bs, NULL)); 12032e68589bSMark F. Adams 12042e68589bSMark F. Adams /* can get all points "removed" */ 12059566063dSJacob Faibussowitsch PetscCall(MatGetSize(Prol, &kk, &ii)); 12067f66b68fSMark Adams if (!ii) { 120763a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: No selected points on coarse grid\n", ((PetscObject)pc)->prefix)); 12089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Prol)); 12090298fd71SBarry Smith *a_P_out = NULL; /* out */ 1210849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0)); 12113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12122e68589bSMark F. Adams } 121363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: New grid %" PetscInt_FMT " nodes\n", ((PetscObject)pc)->prefix, ii / col_bs)); 12149566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk)); 12150cbbd2e1SMark F. Adams 121663a3b9bcSJacob 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); 1217c8b0795cSMark F. Adams myCrs0 = myCrs0 / col_bs; 121863a3b9bcSJacob 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); 12192e68589bSMark F. Adams 12202e68589bSMark F. Adams /* create global vector of data in 'data_w_ghost' */ 1221849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0)); 1222bae903cbSmarkadams4 if (size > 1) { /* get ghost null space data */ 12232e68589bSMark F. Adams PetscReal *tmp_gdata, *tmp_ldata, *tp2; 12249566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc, &tmp_ldata)); 12254a2f8832SBarry Smith for (jj = 0; jj < col_bs; jj++) { 1226c8b0795cSMark F. Adams for (kk = 0; kk < bs; kk++) { 1227a2f3521dSMark F. Adams PetscInt ii, stride; 1228c8b0795cSMark F. Adams const PetscReal *tp = pc_gamg->data + jj * bs * nloc + kk; 12292fa5cd67SKarl Rupp for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp; 12302fa5cd67SKarl Rupp 12319566063dSJacob Faibussowitsch PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata)); 1232a2f3521dSMark F. Adams 1233bae903cbSmarkadams4 if (!jj && !kk) { /* now I know how many total nodes - allocate TODO: move below and do in one 'col_bs' call */ 12349566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(stride * bs * col_bs, &data_w_ghost)); 1235a2f3521dSMark F. Adams nbnodes = bs * stride; 12362e68589bSMark F. Adams } 1237a2f3521dSMark F. Adams tp2 = data_w_ghost + jj * bs * stride + kk; 12382fa5cd67SKarl Rupp for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii]; 12399566063dSJacob Faibussowitsch PetscCall(PetscFree(tmp_gdata)); 12402e68589bSMark F. Adams } 12412e68589bSMark F. Adams } 12429566063dSJacob Faibussowitsch PetscCall(PetscFree(tmp_ldata)); 1243806fa848SBarry Smith } else { 1244c8b0795cSMark F. Adams nbnodes = bs * nloc; 1245c8b0795cSMark F. Adams data_w_ghost = (PetscReal *)pc_gamg->data; 12462e68589bSMark F. Adams } 12472e68589bSMark F. Adams 1248bae903cbSmarkadams4 /* get 'flid_fgid' TODO - move up to get 'stride' and do get null space data above in one step (jj loop) */ 1249c5df96a5SBarry Smith if (size > 1) { 12502e68589bSMark F. Adams PetscReal *fid_glid_loc, *fiddata; 1251a2f3521dSMark F. Adams PetscInt stride; 12522e68589bSMark F. Adams 12539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc, &fid_glid_loc)); 12542e68589bSMark F. Adams for (kk = 0; kk < nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0 + kk); 12559566063dSJacob Faibussowitsch PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata)); 1256bae903cbSmarkadams4 PetscCall(PetscMalloc1(stride, &flid_fgid)); /* copy real data to in */ 1257a2f3521dSMark F. Adams for (kk = 0; kk < stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk]; 12589566063dSJacob Faibussowitsch PetscCall(PetscFree(fiddata)); 1259a2f3521dSMark F. Adams 126063a3b9bcSJacob Faibussowitsch PetscCheck(stride == nbnodes / bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "stride %" PetscInt_FMT " != nbnodes %" PetscInt_FMT "/bs %" PetscInt_FMT, stride, nbnodes, bs); 12619566063dSJacob Faibussowitsch PetscCall(PetscFree(fid_glid_loc)); 1262806fa848SBarry Smith } else { 12639566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc, &flid_fgid)); 12642e68589bSMark F. Adams for (kk = 0; kk < nloc; kk++) flid_fgid[kk] = my0 + kk; 12652e68589bSMark F. Adams } 1266849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0)); 1267bae903cbSmarkadams4 /* get P0 */ 1268849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0)); 1269c8b0795cSMark F. Adams { 12700298fd71SBarry Smith PetscReal *data_out = NULL; 12719566063dSJacob Faibussowitsch PetscCall(formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes, data_w_ghost, flid_fgid, &data_out, Prol)); 12729566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->data)); 12732fa5cd67SKarl Rupp 1274c8b0795cSMark F. Adams pc_gamg->data = data_out; 12754a2f8832SBarry Smith pc_gamg->data_cell_rows = col_bs; 12764a2f8832SBarry Smith pc_gamg->data_sz = col_bs * col_bs * nLocalSelected; 1277c8b0795cSMark F. Adams } 1278849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0)); 127948a46eb9SPierre Jolivet if (size > 1) PetscCall(PetscFree(data_w_ghost)); 12809566063dSJacob Faibussowitsch PetscCall(PetscFree(flid_fgid)); 12812e68589bSMark F. Adams 1282c8b0795cSMark F. Adams *a_P_out = Prol; /* out */ 1283ed4e82eaSPeter Brune 1284849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0)); 12853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1286c8b0795cSMark F. Adams } 1287c8b0795cSMark F. Adams 1288c8b0795cSMark F. Adams /* 1289fd1112cbSBarry Smith PCGAMGOptProlongator_AGG 1290c8b0795cSMark F. Adams 1291c8b0795cSMark F. Adams Input Parameter: 1292c8b0795cSMark F. Adams . pc - this 1293c8b0795cSMark F. Adams . Amat - matrix on this fine level 1294c8b0795cSMark F. Adams In/Output Parameter: 12952a808120SBarry Smith . a_P - prolongation operator to the next level 1296c8b0795cSMark F. Adams */ 1297d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGOptProlongator_AGG(PC pc, Mat Amat, Mat *a_P) 1298d71ae5a4SJacob Faibussowitsch { 1299c8b0795cSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 1300c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1301c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 1302c8b0795cSMark F. Adams PetscInt jj; 1303c8b0795cSMark F. Adams Mat Prol = *a_P; 13043b4367a7SBarry Smith MPI_Comm comm; 13052a808120SBarry Smith KSP eksp; 13062a808120SBarry Smith Vec bb, xx; 13072a808120SBarry Smith PC epc; 13082a808120SBarry Smith PetscReal alpha, emax, emin; 1309c8b0795cSMark F. Adams 1310c8b0795cSMark F. Adams PetscFunctionBegin; 13119566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm)); 1312849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0)); 1313c8b0795cSMark F. Adams 1314d5d25401SBarry Smith /* compute maximum singular value of operator to be used in smoother */ 13152a808120SBarry Smith if (0 < pc_gamg_agg->nsmooths) { 131618c3aa7eSMark /* get eigen estimates */ 131718c3aa7eSMark if (pc_gamg->emax > 0) { 131818c3aa7eSMark emin = pc_gamg->emin; 131918c3aa7eSMark emax = pc_gamg->emax; 132018c3aa7eSMark } else { 13210ed2132dSStefano Zampini const char *prefix; 13220ed2132dSStefano Zampini 13239566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Amat, &bb, NULL)); 13249566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Amat, &xx, NULL)); 13259566063dSJacob Faibussowitsch PetscCall(KSPSetNoisy_Private(bb)); 1326e696ed0bSMark F. Adams 13279566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm, &eksp)); 13283821be0aSBarry Smith PetscCall(KSPSetNestLevel(eksp, pc->kspnestlevel)); 13299566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 13309566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(eksp, prefix)); 13319566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(eksp, "pc_gamg_esteig_")); 133273f7197eSJed Brown { 1333b94d7dedSBarry Smith PetscBool isset, sflg; 1334b94d7dedSBarry Smith PetscCall(MatIsSPDKnown(Amat, &isset, &sflg)); 1335b94d7dedSBarry Smith if (isset && sflg) PetscCall(KSPSetType(eksp, KSPCG)); 1336d24ecf33SMark } 13379566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(eksp, pc->erroriffailure)); 13389566063dSJacob Faibussowitsch PetscCall(KSPSetNormType(eksp, KSP_NORM_NONE)); 1339c2436cacSMark F. Adams 13409566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(eksp, PETSC_FALSE)); 13419566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(eksp, Amat, Amat)); 13422e68589bSMark F. Adams 13439566063dSJacob Faibussowitsch PetscCall(KSPGetPC(eksp, &epc)); 13449566063dSJacob Faibussowitsch PetscCall(PCSetType(epc, PCJACOBI)); /* smoother in smoothed agg. */ 1345c2436cacSMark F. Adams 13469566063dSJacob 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 1347074aec55SMark Adams 13489566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(eksp)); 13499566063dSJacob Faibussowitsch PetscCall(KSPSetComputeSingularValues(eksp, PETSC_TRUE)); 13509566063dSJacob Faibussowitsch PetscCall(KSPSolve(eksp, bb, xx)); 13519566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(eksp, pc, xx)); 13522e68589bSMark F. Adams 13539566063dSJacob Faibussowitsch PetscCall(KSPComputeExtremeSingularValues(eksp, &emax, &emin)); 13549566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "%s: Smooth P0: max eigen=%e min=%e PC=%s\n", ((PetscObject)pc)->prefix, (double)emax, (double)emin, PCJACOBI)); 13559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xx)); 13569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bb)); 13579566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&eksp)); 13582e68589bSMark F. Adams } 135918c3aa7eSMark if (pc_gamg->use_sa_esteig) { 136018c3aa7eSMark mg->min_eigen_DinvA[pc_gamg->current_level] = emin; 136118c3aa7eSMark mg->max_eigen_DinvA[pc_gamg->current_level] = emax; 136263a3b9bcSJacob 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)); 136318c3aa7eSMark } else { 136418c3aa7eSMark mg->min_eigen_DinvA[pc_gamg->current_level] = 0; 136518c3aa7eSMark mg->max_eigen_DinvA[pc_gamg->current_level] = 0; 136618c3aa7eSMark } 136718c3aa7eSMark } else { 136818c3aa7eSMark mg->min_eigen_DinvA[pc_gamg->current_level] = 0; 136918c3aa7eSMark mg->max_eigen_DinvA[pc_gamg->current_level] = 0; 137018c3aa7eSMark } 13712e68589bSMark F. Adams 13722a808120SBarry Smith /* smooth P0 */ 13732a808120SBarry Smith for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) { 13742a808120SBarry Smith Mat tMat; 13752a808120SBarry Smith Vec diag; 13762a808120SBarry Smith 1377849bee69SMark Adams PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0)); 13782a808120SBarry Smith 13792e68589bSMark F. Adams /* smooth P1 := (I - omega/lam D^{-1}A)P0 */ 13809566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0)); 13819566063dSJacob Faibussowitsch PetscCall(MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat)); 13829566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0)); 13839566063dSJacob Faibussowitsch PetscCall(MatProductClear(tMat)); 13849566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Amat, &diag, NULL)); 13859566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(Amat, diag)); /* effectively PCJACOBI */ 13869566063dSJacob Faibussowitsch PetscCall(VecReciprocal(diag)); 13879566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(tMat, diag, NULL)); 13889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&diag)); 1389b4da3a1bSStefano Zampini 1390d5d25401SBarry Smith /* TODO: Set a PCFailedReason and exit the building of the AMG preconditioner */ 139108401ef6SPierre Jolivet PetscCheck(emax != 0.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Computed maximum singular value as zero"); 1392d5d25401SBarry Smith /* TODO: Document the 1.4 and don't hardwire it in this routine */ 1393b4ec6429SMark F. Adams alpha = -1.4 / emax; 1394b4da3a1bSStefano Zampini 13959566063dSJacob Faibussowitsch PetscCall(MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN)); 13969566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Prol)); 13972e68589bSMark F. Adams Prol = tMat; 1398849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0)); 13992e68589bSMark F. Adams } 1400849bee69SMark Adams PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0)); 1401c8b0795cSMark F. Adams *a_P = Prol; 14023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14032e68589bSMark F. Adams } 14042e68589bSMark F. Adams 1405c8b0795cSMark F. Adams /* 1406c8b0795cSMark F. Adams PCCreateGAMG_AGG 14072e68589bSMark F. Adams 1408c8b0795cSMark F. Adams Input Parameter: 1409c8b0795cSMark F. Adams . pc - 1410c8b0795cSMark F. Adams */ 1411d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCreateGAMG_AGG(PC pc) 1412d71ae5a4SJacob Faibussowitsch { 1413c8b0795cSMark F. Adams PC_MG *mg = (PC_MG *)pc->data; 1414c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1415c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg; 14162e68589bSMark F. Adams 1417c8b0795cSMark F. Adams PetscFunctionBegin; 1418c8b0795cSMark F. Adams /* create sub context for SA */ 14194dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&pc_gamg_agg)); 1420c8b0795cSMark F. Adams pc_gamg->subctx = pc_gamg_agg; 1421c8b0795cSMark F. Adams 14221ab5ffc9SJed Brown pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG; 14239b8ffb57SJed Brown pc_gamg->ops->destroy = PCDestroy_GAMG_AGG; 1424c8b0795cSMark F. Adams /* reset does not do anything; setup not virtual */ 1425c8b0795cSMark F. Adams 1426c8b0795cSMark F. Adams /* set internal function pointers */ 14272d776b49SBarry Smith pc_gamg->ops->creategraph = PCGAMGCreateGraph_AGG; 14281ab5ffc9SJed Brown pc_gamg->ops->coarsen = PCGAMGCoarsen_AGG; 14291ab5ffc9SJed Brown pc_gamg->ops->prolongator = PCGAMGProlongator_AGG; 1430fd1112cbSBarry Smith pc_gamg->ops->optprolongator = PCGAMGOptProlongator_AGG; 14311ab5ffc9SJed Brown pc_gamg->ops->createdefaultdata = PCSetData_AGG; 14325adeb434SBarry Smith pc_gamg->ops->view = PCView_GAMG_AGG; 1433c8b0795cSMark F. Adams 1434cab9ed1eSBarry Smith pc_gamg_agg->nsmooths = 1; 1435d529f056Smarkadams4 pc_gamg_agg->aggressive_coarsening_levels = 1; 1436d529f056Smarkadams4 pc_gamg_agg->use_aggressive_square_graph = PETSC_FALSE; 1437d529f056Smarkadams4 pc_gamg_agg->use_minimum_degree_ordering = PETSC_FALSE; 1438*a9ccf005SMark Adams pc_gamg_agg->use_low_mem_filter = PETSC_FALSE; 1439d529f056Smarkadams4 pc_gamg_agg->aggressive_mis_k = 2; 1440cab9ed1eSBarry Smith 14419566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", PCGAMGSetNSmooths_AGG)); 1442bae903cbSmarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", PCGAMGSetAggressiveLevels_AGG)); 1443d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveSquareGraph_C", PCGAMGSetAggressiveSquareGraph_AGG)); 1444d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetMinDegreeOrdering_C", PCGAMGMISkSetMinDegreeOrdering_AGG)); 1445*a9ccf005SMark Adams PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetLowMemoryFilter_C", PCGAMGSetLowMemoryFilter_AGG)); 1446d529f056Smarkadams4 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetAggressive_C", PCGAMGMISkSetAggressive_AGG)); 14479566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_AGG)); 14483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14492e68589bSMark F. Adams } 1450