18be712e4SBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 28be712e4SBarry Smith #include <../src/mat/impls/aij/seq/aij.h> 38be712e4SBarry Smith #include <../src/mat/impls/aij/mpi/mpiaij.h> 4452acf8bSMark Adams #include <petscdm.h> 58be712e4SBarry Smith 68be712e4SBarry Smith /* linked list methods 78be712e4SBarry Smith * 88be712e4SBarry Smith * PetscCDCreate 98be712e4SBarry Smith */ 108be712e4SBarry Smith PetscErrorCode PetscCDCreate(PetscInt a_size, PetscCoarsenData **a_out) 118be712e4SBarry Smith { 128be712e4SBarry Smith PetscCoarsenData *ail; 138be712e4SBarry Smith 148be712e4SBarry Smith PetscFunctionBegin; 158be712e4SBarry Smith /* allocate pool, partially */ 168be712e4SBarry Smith PetscCall(PetscNew(&ail)); 178be712e4SBarry Smith *a_out = ail; 188be712e4SBarry Smith ail->pool_list.next = NULL; 198be712e4SBarry Smith ail->pool_list.array = NULL; 208be712e4SBarry Smith ail->chk_sz = 0; 218be712e4SBarry Smith /* allocate array */ 228be712e4SBarry Smith ail->size = a_size; 238be712e4SBarry Smith PetscCall(PetscCalloc1(a_size, &ail->array)); 248be712e4SBarry Smith ail->extra_nodes = NULL; 258be712e4SBarry Smith ail->mat = NULL; 268be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 278be712e4SBarry Smith } 288be712e4SBarry Smith 298be712e4SBarry Smith /* PetscCDDestroy 308be712e4SBarry Smith */ 318be712e4SBarry Smith PetscErrorCode PetscCDDestroy(PetscCoarsenData *ail) 328be712e4SBarry Smith { 338be712e4SBarry Smith PetscCDArrNd *n = &ail->pool_list; 348be712e4SBarry Smith 358be712e4SBarry Smith PetscFunctionBegin; 368be712e4SBarry Smith n = n->next; 378be712e4SBarry Smith while (n) { 388be712e4SBarry Smith PetscCDArrNd *lstn = n; 39*c376f201SBarry Smith 408be712e4SBarry Smith n = n->next; 418be712e4SBarry Smith PetscCall(PetscFree(lstn)); 428be712e4SBarry Smith } 438be712e4SBarry Smith if (ail->pool_list.array) PetscCall(PetscFree(ail->pool_list.array)); 448be712e4SBarry Smith PetscCall(PetscFree(ail->array)); 458be712e4SBarry Smith if (ail->mat) PetscCall(MatDestroy(&ail->mat)); 468be712e4SBarry Smith /* delete this (+agg+pool array) */ 478be712e4SBarry Smith PetscCall(PetscFree(ail)); 488be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 498be712e4SBarry Smith } 508be712e4SBarry Smith 518be712e4SBarry Smith /* PetscCDSetChunkSize 528be712e4SBarry Smith */ 538be712e4SBarry Smith PetscErrorCode PetscCDSetChunkSize(PetscCoarsenData *ail, PetscInt a_sz) 548be712e4SBarry Smith { 558be712e4SBarry Smith PetscFunctionBegin; 568be712e4SBarry Smith ail->chk_sz = a_sz; 578be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 588be712e4SBarry Smith } 598be712e4SBarry Smith 608be712e4SBarry Smith /* PetscCDGetNewNode 618be712e4SBarry Smith */ 628be712e4SBarry Smith static PetscErrorCode PetscCDGetNewNode(PetscCoarsenData *ail, PetscCDIntNd **a_out, PetscInt a_id) 638be712e4SBarry Smith { 648be712e4SBarry Smith PetscFunctionBegin; 658be712e4SBarry Smith *a_out = NULL; /* squelch -Wmaybe-uninitialized */ 668be712e4SBarry Smith if (ail->extra_nodes) { 678be712e4SBarry Smith PetscCDIntNd *node = ail->extra_nodes; 68*c376f201SBarry Smith 698be712e4SBarry Smith ail->extra_nodes = node->next; 708be712e4SBarry Smith node->gid = a_id; 718be712e4SBarry Smith node->next = NULL; 728be712e4SBarry Smith *a_out = node; 738be712e4SBarry Smith } else { 748be712e4SBarry Smith if (!ail->pool_list.array) { 758be712e4SBarry Smith if (!ail->chk_sz) ail->chk_sz = 10; /* use a chuck size of ail->size? */ 768be712e4SBarry Smith PetscCall(PetscMalloc1(ail->chk_sz, &ail->pool_list.array)); 778be712e4SBarry Smith ail->new_node = ail->pool_list.array; 788be712e4SBarry Smith ail->new_left = ail->chk_sz; 798be712e4SBarry Smith ail->new_node->next = NULL; 808be712e4SBarry Smith } else if (!ail->new_left) { 818be712e4SBarry Smith PetscCDArrNd *node; 82*c376f201SBarry Smith 838be712e4SBarry Smith PetscCall(PetscMalloc(ail->chk_sz * sizeof(PetscCDIntNd) + sizeof(PetscCDArrNd), &node)); 848be712e4SBarry Smith node->array = (PetscCDIntNd *)(node + 1); 858be712e4SBarry Smith node->next = ail->pool_list.next; 868be712e4SBarry Smith ail->pool_list.next = node; 878be712e4SBarry Smith ail->new_left = ail->chk_sz; 888be712e4SBarry Smith ail->new_node = node->array; 898be712e4SBarry Smith } 908be712e4SBarry Smith ail->new_node->gid = a_id; 918be712e4SBarry Smith ail->new_node->next = NULL; 928be712e4SBarry Smith *a_out = ail->new_node++; 938be712e4SBarry Smith ail->new_left--; 948be712e4SBarry Smith } 958be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 968be712e4SBarry Smith } 978be712e4SBarry Smith 988be712e4SBarry Smith /* PetscCDIntNdSetID 998be712e4SBarry Smith */ 1008be712e4SBarry Smith PetscErrorCode PetscCDIntNdSetID(PetscCDIntNd *a_this, PetscInt a_id) 1018be712e4SBarry Smith { 1028be712e4SBarry Smith PetscFunctionBegin; 1038be712e4SBarry Smith a_this->gid = a_id; 1048be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 1058be712e4SBarry Smith } 1068be712e4SBarry Smith 1078be712e4SBarry Smith /* PetscCDIntNdGetID 1088be712e4SBarry Smith */ 1098be712e4SBarry Smith PetscErrorCode PetscCDIntNdGetID(const PetscCDIntNd *a_this, PetscInt *a_gid) 1108be712e4SBarry Smith { 1118be712e4SBarry Smith PetscFunctionBegin; 1128be712e4SBarry Smith *a_gid = a_this->gid; 1138be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 1148be712e4SBarry Smith } 1158be712e4SBarry Smith 1168be712e4SBarry Smith /* PetscCDGetHeadPos 1178be712e4SBarry Smith */ 1188be712e4SBarry Smith PetscErrorCode PetscCDGetHeadPos(const PetscCoarsenData *ail, PetscInt a_idx, PetscCDIntNd **pos) 1198be712e4SBarry Smith { 1208be712e4SBarry Smith PetscFunctionBegin; 1218be712e4SBarry Smith PetscCheck(a_idx < ail->size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "a_idx >= ail->size: a_idx=%" PetscInt_FMT ".", a_idx); 1228be712e4SBarry Smith *pos = ail->array[a_idx]; 1238be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 1248be712e4SBarry Smith } 1258be712e4SBarry Smith 1268be712e4SBarry Smith /* PetscCDGetNextPos 1278be712e4SBarry Smith */ 1288be712e4SBarry Smith PetscErrorCode PetscCDGetNextPos(const PetscCoarsenData *ail, PetscInt l_idx, PetscCDIntNd **pos) 1298be712e4SBarry Smith { 1308be712e4SBarry Smith PetscFunctionBegin; 131f4f49eeaSPierre Jolivet PetscCheck(*pos, PETSC_COMM_SELF, PETSC_ERR_PLIB, "NULL input position."); 1328be712e4SBarry Smith *pos = (*pos)->next; 1338be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 1348be712e4SBarry Smith } 1358be712e4SBarry Smith 1368be712e4SBarry Smith /* PetscCDAppendID 1378be712e4SBarry Smith */ 1388be712e4SBarry Smith PetscErrorCode PetscCDAppendID(PetscCoarsenData *ail, PetscInt a_idx, PetscInt a_id) 1398be712e4SBarry Smith { 1408be712e4SBarry Smith PetscCDIntNd *n, *n2; 1418be712e4SBarry Smith 1428be712e4SBarry Smith PetscFunctionBegin; 1438be712e4SBarry Smith PetscCall(PetscCDGetNewNode(ail, &n, a_id)); 1448be712e4SBarry Smith PetscCheck(a_idx < ail->size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Index %" PetscInt_FMT " out of range.", a_idx); 1458be712e4SBarry Smith if (!(n2 = ail->array[a_idx])) ail->array[a_idx] = n; 1468be712e4SBarry Smith else { 1478be712e4SBarry Smith do { 1488be712e4SBarry Smith if (!n2->next) { 1498be712e4SBarry Smith n2->next = n; 1508be712e4SBarry Smith PetscCheck(!n->next, PETSC_COMM_SELF, PETSC_ERR_PLIB, "n should not have a next"); 1518be712e4SBarry Smith break; 1528be712e4SBarry Smith } 1538be712e4SBarry Smith n2 = n2->next; 1548be712e4SBarry Smith } while (n2); 1558be712e4SBarry Smith PetscCheck(n2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "n2 should be non-null"); 1568be712e4SBarry Smith } 1578be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 1588be712e4SBarry Smith } 1598be712e4SBarry Smith 1608be712e4SBarry Smith /* PetscCDAppendNode 1618be712e4SBarry Smith */ 1628be712e4SBarry Smith PetscErrorCode PetscCDAppendNode(PetscCoarsenData *ail, PetscInt a_idx, PetscCDIntNd *a_n) 1638be712e4SBarry Smith { 1648be712e4SBarry Smith PetscCDIntNd *n2; 1658be712e4SBarry Smith 1668be712e4SBarry Smith PetscFunctionBegin; 1678be712e4SBarry Smith PetscCheck(a_idx < ail->size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Index %" PetscInt_FMT " out of range.", a_idx); 1688be712e4SBarry Smith if (!(n2 = ail->array[a_idx])) ail->array[a_idx] = a_n; 1698be712e4SBarry Smith else { 1708be712e4SBarry Smith do { 1718be712e4SBarry Smith if (!n2->next) { 1728be712e4SBarry Smith n2->next = a_n; 1738be712e4SBarry Smith a_n->next = NULL; 1748be712e4SBarry Smith break; 1758be712e4SBarry Smith } 1768be712e4SBarry Smith n2 = n2->next; 1778be712e4SBarry Smith } while (n2); 1788be712e4SBarry Smith PetscCheck(n2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "n2 should be non-null"); 1798be712e4SBarry Smith } 1808be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 1818be712e4SBarry Smith } 1828be712e4SBarry Smith 1838be712e4SBarry Smith /* PetscCDRemoveNextNode: a_last->next, this exposes single linked list structure to API (not used) 1848be712e4SBarry Smith */ 1858be712e4SBarry Smith PetscErrorCode PetscCDRemoveNextNode(PetscCoarsenData *ail, PetscInt a_idx, PetscCDIntNd *a_last) 1868be712e4SBarry Smith { 1878be712e4SBarry Smith PetscCDIntNd *del; 1888be712e4SBarry Smith 1898be712e4SBarry Smith PetscFunctionBegin; 1908be712e4SBarry Smith PetscCheck(a_idx < ail->size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Index %" PetscInt_FMT " out of range.", a_idx); 1918be712e4SBarry Smith PetscCheck(a_last->next, PETSC_COMM_SELF, PETSC_ERR_PLIB, "a_last should have a next"); 1928be712e4SBarry Smith del = a_last->next; 1938be712e4SBarry Smith a_last->next = del->next; 1948be712e4SBarry Smith /* del->next = NULL; -- this still used in a iterator so keep it intact -- need to fix this with a double linked list */ 1958be712e4SBarry Smith /* could reuse n2 but PetscCDAppendNode sometimes uses it */ 1968be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 1978be712e4SBarry Smith } 1988be712e4SBarry Smith 1998be712e4SBarry Smith /* PetscCDPrint 2008be712e4SBarry Smith */ 201*c376f201SBarry Smith PetscErrorCode PetscCDPrint(const PetscCoarsenData *ail, PetscInt Istart, MPI_Comm comm) 2028be712e4SBarry Smith { 2038be712e4SBarry Smith PetscCDIntNd *n, *n2; 2048be712e4SBarry Smith PetscInt ii; 2058be712e4SBarry Smith 2068be712e4SBarry Smith PetscFunctionBegin; 2078be712e4SBarry Smith for (ii = 0; ii < ail->size; ii++) { 2088be712e4SBarry Smith n2 = n = ail->array[ii]; 209*c376f201SBarry Smith if (n) PetscCall(PetscSynchronizedPrintf(comm, "list %" PetscInt_FMT ":", ii + Istart)); 2108be712e4SBarry Smith while (n) { 2118be712e4SBarry Smith PetscCall(PetscSynchronizedPrintf(comm, " %" PetscInt_FMT, n->gid)); 2128be712e4SBarry Smith n = n->next; 2138be712e4SBarry Smith } 2148be712e4SBarry Smith if (n2) PetscCall(PetscSynchronizedPrintf(comm, "\n")); 2158be712e4SBarry Smith } 2168be712e4SBarry Smith PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT)); 2178be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 2188be712e4SBarry Smith } 2198be712e4SBarry Smith 2208be712e4SBarry Smith /* PetscCDMoveAppend - take list in a_srcidx and appends to destidx 2218be712e4SBarry Smith */ 2228be712e4SBarry Smith PetscErrorCode PetscCDMoveAppend(PetscCoarsenData *ail, PetscInt a_destidx, PetscInt a_srcidx) 2238be712e4SBarry Smith { 2248be712e4SBarry Smith PetscCDIntNd *n; 2258be712e4SBarry Smith 2268be712e4SBarry Smith PetscFunctionBegin; 2278be712e4SBarry Smith PetscCheck(a_srcidx < ail->size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Index %" PetscInt_FMT " out of range.", a_srcidx); 2288be712e4SBarry Smith PetscCheck(a_destidx < ail->size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Index %" PetscInt_FMT " out of range.", a_destidx); 2298be712e4SBarry Smith PetscCheck(a_destidx != a_srcidx, PETSC_COMM_SELF, PETSC_ERR_PLIB, "a_destidx==a_srcidx %" PetscInt_FMT ".", a_destidx); 2308be712e4SBarry Smith n = ail->array[a_destidx]; 2318be712e4SBarry Smith if (!n) ail->array[a_destidx] = ail->array[a_srcidx]; 2328be712e4SBarry Smith else { 2338be712e4SBarry Smith do { 2348be712e4SBarry Smith if (!n->next) { 2358be712e4SBarry Smith n->next = ail->array[a_srcidx]; // append 2368be712e4SBarry Smith break; 2378be712e4SBarry Smith } 2388be712e4SBarry Smith n = n->next; 2398be712e4SBarry Smith } while (1); 2408be712e4SBarry Smith } 2418be712e4SBarry Smith ail->array[a_srcidx] = NULL; // empty 2428be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 2438be712e4SBarry Smith } 2448be712e4SBarry Smith 2458be712e4SBarry Smith /* PetscCDRemoveAllAt - empty one list and move data to cache 2468be712e4SBarry Smith */ 2478be712e4SBarry Smith PetscErrorCode PetscCDRemoveAllAt(PetscCoarsenData *ail, PetscInt a_idx) 2488be712e4SBarry Smith { 2498be712e4SBarry Smith PetscCDIntNd *rem, *n1; 2508be712e4SBarry Smith 2518be712e4SBarry Smith PetscFunctionBegin; 2528be712e4SBarry Smith PetscCheck(a_idx < ail->size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Index %" PetscInt_FMT " out of range.", a_idx); 2538be712e4SBarry Smith rem = ail->array[a_idx]; 2548be712e4SBarry Smith ail->array[a_idx] = NULL; 2558be712e4SBarry Smith if (!(n1 = ail->extra_nodes)) ail->extra_nodes = rem; 2568be712e4SBarry Smith else { 2578be712e4SBarry Smith while (n1->next) n1 = n1->next; 2588be712e4SBarry Smith n1->next = rem; 2598be712e4SBarry Smith } 2608be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 2618be712e4SBarry Smith } 2628be712e4SBarry Smith 2638be712e4SBarry Smith /* PetscCDCountAt 2648be712e4SBarry Smith */ 2658be712e4SBarry Smith PetscErrorCode PetscCDCountAt(const PetscCoarsenData *ail, PetscInt a_idx, PetscInt *a_sz) 2668be712e4SBarry Smith { 2678be712e4SBarry Smith PetscCDIntNd *n1; 2688be712e4SBarry Smith PetscInt sz = 0; 2698be712e4SBarry Smith 2708be712e4SBarry Smith PetscFunctionBegin; 2718be712e4SBarry Smith PetscCheck(a_idx < ail->size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Index %" PetscInt_FMT " out of range.", a_idx); 2728be712e4SBarry Smith n1 = ail->array[a_idx]; 2738be712e4SBarry Smith while (n1) { 2748be712e4SBarry Smith n1 = n1->next; 2758be712e4SBarry Smith sz++; 2768be712e4SBarry Smith } 2778be712e4SBarry Smith *a_sz = sz; 2788be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 2798be712e4SBarry Smith } 2808be712e4SBarry Smith 2818be712e4SBarry Smith /* PetscCDSize 2828be712e4SBarry Smith */ 2838be712e4SBarry Smith PetscErrorCode PetscCDCount(const PetscCoarsenData *ail, PetscInt *a_sz) 2848be712e4SBarry Smith { 2858be712e4SBarry Smith PetscInt sz = 0; 2868be712e4SBarry Smith 2878be712e4SBarry Smith PetscFunctionBegin; 288*c376f201SBarry Smith for (PetscInt ii = 0; ii < ail->size; ii++) { 2898be712e4SBarry Smith PetscCDIntNd *n1 = ail->array[ii]; 290*c376f201SBarry Smith 2918be712e4SBarry Smith while (n1) { 2928be712e4SBarry Smith n1 = n1->next; 2938be712e4SBarry Smith sz++; 2948be712e4SBarry Smith } 2958be712e4SBarry Smith } 2968be712e4SBarry Smith *a_sz = sz; 2978be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 2988be712e4SBarry Smith } 2998be712e4SBarry Smith 3008be712e4SBarry Smith /* PetscCDIsEmptyAt - Is the list empty? (not used) 3018be712e4SBarry Smith */ 3028be712e4SBarry Smith PetscErrorCode PetscCDIsEmptyAt(const PetscCoarsenData *ail, PetscInt a_idx, PetscBool *a_e) 3038be712e4SBarry Smith { 3048be712e4SBarry Smith PetscFunctionBegin; 3058be712e4SBarry Smith PetscCheck(a_idx < ail->size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Index %" PetscInt_FMT " out of range.", a_idx); 3068be712e4SBarry Smith *a_e = (PetscBool)(ail->array[a_idx] == NULL); 3078be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 3088be712e4SBarry Smith } 3098be712e4SBarry Smith 3108be712e4SBarry Smith /* PetscCDGetNonemptyIS - used for C-F methods 3118be712e4SBarry Smith */ 3128be712e4SBarry Smith PetscErrorCode PetscCDGetNonemptyIS(PetscCoarsenData *ail, IS *a_mis) 3138be712e4SBarry Smith { 3148be712e4SBarry Smith PetscCDIntNd *n; 3158be712e4SBarry Smith PetscInt ii, kk; 3168be712e4SBarry Smith PetscInt *permute; 3178be712e4SBarry Smith 3188be712e4SBarry Smith PetscFunctionBegin; 3198be712e4SBarry Smith for (ii = kk = 0; ii < ail->size; ii++) { 3208be712e4SBarry Smith n = ail->array[ii]; 3218be712e4SBarry Smith if (n) kk++; 3228be712e4SBarry Smith } 3238be712e4SBarry Smith PetscCall(PetscMalloc1(kk, &permute)); 3248be712e4SBarry Smith for (ii = kk = 0; ii < ail->size; ii++) { 3258be712e4SBarry Smith n = ail->array[ii]; 3268be712e4SBarry Smith if (n) permute[kk++] = ii; 3278be712e4SBarry Smith } 3288be712e4SBarry Smith PetscCall(ISCreateGeneral(PETSC_COMM_SELF, kk, permute, PETSC_OWN_POINTER, a_mis)); 3298be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 3308be712e4SBarry Smith } 3318be712e4SBarry Smith 3328be712e4SBarry Smith /* PetscCDGetMat 3338be712e4SBarry Smith */ 3348be712e4SBarry Smith PetscErrorCode PetscCDGetMat(PetscCoarsenData *ail, Mat *a_mat) 3358be712e4SBarry Smith { 3368be712e4SBarry Smith PetscFunctionBegin; 3378be712e4SBarry Smith *a_mat = ail->mat; 3388be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 3398be712e4SBarry Smith } 3408be712e4SBarry Smith 3418be712e4SBarry Smith /* PetscCDSetMat 3428be712e4SBarry Smith */ 3438be712e4SBarry Smith PetscErrorCode PetscCDSetMat(PetscCoarsenData *ail, Mat a_mat) 3448be712e4SBarry Smith { 3458be712e4SBarry Smith PetscFunctionBegin; 3468be712e4SBarry Smith if (ail->mat) { 3478be712e4SBarry Smith PetscCall(MatDestroy(&ail->mat)); //should not happen 3488be712e4SBarry Smith } 3498be712e4SBarry Smith ail->mat = a_mat; 3508be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 3518be712e4SBarry Smith } 3528be712e4SBarry Smith 3538926f930SMark Adams /* PetscCDClearMat 3548926f930SMark Adams */ 3558926f930SMark Adams PetscErrorCode PetscCDClearMat(PetscCoarsenData *ail) 3568926f930SMark Adams { 3578926f930SMark Adams PetscFunctionBegin; 3588926f930SMark Adams ail->mat = NULL; 3598926f930SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 3608926f930SMark Adams } 3618926f930SMark Adams 3628be712e4SBarry Smith /* PetscCDGetASMBlocks - get IS of aggregates for ASM smoothers 3638be712e4SBarry Smith */ 3648be712e4SBarry Smith PetscErrorCode PetscCDGetASMBlocks(const PetscCoarsenData *ail, const PetscInt a_bs, PetscInt *a_sz, IS **a_local_is) 3658be712e4SBarry Smith { 3668be712e4SBarry Smith PetscCDIntNd *n; 3678be712e4SBarry Smith PetscInt lsz, ii, kk, *idxs, jj, gid; 3688be712e4SBarry Smith IS *is_loc = NULL; 3698be712e4SBarry Smith 3708be712e4SBarry Smith PetscFunctionBegin; 3718be712e4SBarry Smith for (ii = kk = 0; ii < ail->size; ii++) { 3728be712e4SBarry Smith if (ail->array[ii]) kk++; 3738be712e4SBarry Smith } 3748be712e4SBarry Smith *a_sz = kk; 3758be712e4SBarry Smith PetscCall(PetscMalloc1(kk, &is_loc)); 3768be712e4SBarry Smith for (ii = kk = 0; ii < ail->size; ii++) { 3778be712e4SBarry Smith for (lsz = 0, n = ail->array[ii]; n; lsz++, n = n->next) /* void */ 3788be712e4SBarry Smith ; 3798be712e4SBarry Smith if (lsz) { 3808be712e4SBarry Smith PetscCall(PetscMalloc1(a_bs * lsz, &idxs)); 3818be712e4SBarry Smith for (lsz = 0, n = ail->array[ii]; n; n = n->next) { 3828be712e4SBarry Smith PetscCall(PetscCDIntNdGetID(n, &gid)); 3838be712e4SBarry Smith for (jj = 0; jj < a_bs; lsz++, jj++) idxs[lsz] = a_bs * gid + jj; 3848be712e4SBarry Smith } 3858be712e4SBarry Smith PetscCall(ISCreateGeneral(PETSC_COMM_SELF, lsz, idxs, PETSC_OWN_POINTER, &is_loc[kk++])); 3868be712e4SBarry Smith } 3878be712e4SBarry Smith } 3888be712e4SBarry Smith PetscCheck(*a_sz == kk, PETSC_COMM_SELF, PETSC_ERR_PLIB, "*a_sz %" PetscInt_FMT " != kk %" PetscInt_FMT, *a_sz, kk); 3898be712e4SBarry Smith *a_local_is = is_loc; /* out */ 3908be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 3918be712e4SBarry Smith } 3928be712e4SBarry Smith 3938be712e4SBarry Smith /* edge for priority queue */ 3948be712e4SBarry Smith typedef struct edge_tag { 3958be712e4SBarry Smith PetscReal weight; 3968be712e4SBarry Smith PetscInt lid0, gid1, ghost1_idx; 3978be712e4SBarry Smith } Edge; 3988be712e4SBarry Smith 3998be712e4SBarry Smith #define MY_MEPS (PETSC_MACHINE_EPSILON * 100) 4008be712e4SBarry Smith static int gamg_hem_compare(const void *a, const void *b) 4018be712e4SBarry Smith { 4028be712e4SBarry Smith PetscReal va = ((Edge *)a)->weight, vb = ((Edge *)b)->weight; 4038be712e4SBarry Smith return (va <= vb - MY_MEPS) ? 1 : (va > vb + MY_MEPS) ? -1 : 0; /* 0 for equal */ 4048be712e4SBarry Smith } 4058be712e4SBarry Smith 4068be712e4SBarry Smith /* 4078be712e4SBarry Smith MatCoarsenApply_HEM_private - parallel heavy edge matching 4088be712e4SBarry Smith 4098be712e4SBarry Smith Input Parameter: 4108be712e4SBarry Smith . a_Gmat - global matrix of the graph 4118be712e4SBarry Smith . n_iter - number of matching iterations 4128be712e4SBarry Smith . threshold - threshold for filtering graphs 4138be712e4SBarry Smith 4148be712e4SBarry Smith Output Parameter: 4158be712e4SBarry Smith . a_locals_llist - array of list of local nodes rooted at local node 4168be712e4SBarry Smith */ 4178be712e4SBarry Smith static PetscErrorCode MatCoarsenApply_HEM_private(Mat a_Gmat, const PetscInt n_iter, const PetscReal threshold, PetscCoarsenData **a_locals_llist) 4188be712e4SBarry Smith { 4198be712e4SBarry Smith #define REQ_BF_SIZE 100 4208be712e4SBarry Smith PetscBool isMPI; 4218be712e4SBarry Smith MPI_Comm comm; 422*c376f201SBarry Smith PetscInt ix, *ii, *aj, Istart, bc_agg = -1, *rbuff = NULL, rbuff_sz = 0; 423*c376f201SBarry Smith PetscMPIInt rank, size, comm_procs[REQ_BF_SIZE], ncomm_procs, *lid_max_pe; 424*c376f201SBarry Smith const PetscInt nloc = a_Gmat->rmap->n, request_size = PetscCeilInt((int)sizeof(MPI_Request), (int)sizeof(PetscInt)); 4258be712e4SBarry Smith PetscInt *lid_cprowID; 4268be712e4SBarry Smith PetscBool *lid_matched; 4278be712e4SBarry Smith Mat_SeqAIJ *matA, *matB = NULL; 4288be712e4SBarry Smith Mat_MPIAIJ *mpimat = NULL; 4298be712e4SBarry Smith PetscScalar one = 1.; 4308be712e4SBarry Smith PetscCoarsenData *agg_llists = NULL, *ghost_deleted_list = NULL, *bc_list = NULL; 4318be712e4SBarry Smith Mat cMat, tMat, P; 4328be712e4SBarry Smith MatScalar *ap; 4338be712e4SBarry Smith IS info_is; 4348be712e4SBarry Smith 4358be712e4SBarry Smith PetscFunctionBegin; 4368be712e4SBarry Smith PetscCall(PetscObjectGetComm((PetscObject)a_Gmat, &comm)); 4378be712e4SBarry Smith PetscCallMPI(MPI_Comm_rank(comm, &rank)); 4388be712e4SBarry Smith PetscCallMPI(MPI_Comm_size(comm, &size)); 439*c376f201SBarry Smith PetscCall(MatGetOwnershipRange(a_Gmat, &Istart, NULL)); 4408be712e4SBarry Smith PetscCall(ISCreate(comm, &info_is)); 4418926f930SMark Adams PetscCall(PetscInfo(info_is, "Start %" PetscInt_FMT " iterations of HEM.\n", n_iter)); 4428be712e4SBarry Smith 443*c376f201SBarry Smith PetscCall(PetscMalloc3(nloc, &lid_matched, nloc, &lid_cprowID, nloc, &lid_max_pe)); 4448be712e4SBarry Smith PetscCall(PetscCDCreate(nloc, &agg_llists)); 4458be712e4SBarry Smith PetscCall(PetscCDSetChunkSize(agg_llists, nloc + 1)); 4468be712e4SBarry Smith *a_locals_llist = agg_llists; 4478be712e4SBarry Smith /* add self to all lists */ 448*c376f201SBarry Smith for (PetscInt kk = 0; kk < nloc; kk++) PetscCall(PetscCDAppendID(agg_llists, kk, Istart + kk)); 4498be712e4SBarry Smith /* make a copy of the graph, this gets destroyed in iterates */ 4508be712e4SBarry Smith PetscCall(MatDuplicate(a_Gmat, MAT_COPY_VALUES, &cMat)); 4518be712e4SBarry Smith PetscCall(MatConvert(cMat, MATAIJ, MAT_INPLACE_MATRIX, &cMat)); 4528be712e4SBarry Smith isMPI = (PetscBool)(size > 1); 4538be712e4SBarry Smith if (isMPI) { 4548be712e4SBarry Smith /* list of deleted ghosts, should compress this */ 4558be712e4SBarry Smith PetscCall(PetscCDCreate(size, &ghost_deleted_list)); 4568be712e4SBarry Smith PetscCall(PetscCDSetChunkSize(ghost_deleted_list, 100)); 4578be712e4SBarry Smith } 458*c376f201SBarry Smith for (PetscInt iter = 0; iter < n_iter; iter++) { 4598be712e4SBarry Smith const PetscScalar *lghost_max_ew, *lid_max_ew; 4608be712e4SBarry Smith PetscBool *lghost_matched; 4618be712e4SBarry Smith PetscMPIInt *lghost_pe, *lghost_max_pe; 4628be712e4SBarry Smith Vec locMaxEdge, ghostMaxEdge, ghostMaxPE, locMaxPE; 4638be712e4SBarry Smith PetscInt *lghost_gid, nEdges, nEdges0, num_ghosts = 0; 4648be712e4SBarry Smith Edge *Edges; 465*c376f201SBarry Smith const PetscInt n_sub_its = 1000; // in case of a bug, stop at some point 466*c376f201SBarry Smith 4678be712e4SBarry Smith /* get submatrices of cMat */ 468*c376f201SBarry Smith for (PetscInt kk = 0; kk < nloc; kk++) lid_cprowID[kk] = -1; 4698be712e4SBarry Smith if (isMPI) { 4708be712e4SBarry Smith mpimat = (Mat_MPIAIJ *)cMat->data; 4718be712e4SBarry Smith matA = (Mat_SeqAIJ *)mpimat->A->data; 4728be712e4SBarry Smith matB = (Mat_SeqAIJ *)mpimat->B->data; 4738be712e4SBarry Smith if (!matB->compressedrow.use) { 4748be712e4SBarry Smith /* force construction of compressed row data structure since code below requires it */ 4758be712e4SBarry Smith PetscCall(MatCheckCompressedRow(mpimat->B, matB->nonzerorowcnt, &matB->compressedrow, matB->i, mpimat->B->rmap->n, -1.0)); 4768be712e4SBarry Smith } 4778be712e4SBarry Smith /* set index into compressed row 'lid_cprowID' */ 4788be712e4SBarry Smith for (ix = 0; ix < matB->compressedrow.nrows; ix++) { 4798be712e4SBarry Smith PetscInt *ridx = matB->compressedrow.rindex, lid = ridx[ix]; 4808be712e4SBarry Smith if (ridx[ix] >= 0) lid_cprowID[lid] = ix; 4818be712e4SBarry Smith } 4828be712e4SBarry Smith } else { 4838be712e4SBarry Smith matA = (Mat_SeqAIJ *)cMat->data; 4848be712e4SBarry Smith } 4858be712e4SBarry Smith /* set matched flags: true for empty list */ 486*c376f201SBarry Smith for (PetscInt kk = 0; kk < nloc; kk++) { 4878be712e4SBarry Smith PetscCall(PetscCDCountAt(agg_llists, kk, &ix)); 4888be712e4SBarry Smith if (ix > 0) lid_matched[kk] = PETSC_FALSE; 4898be712e4SBarry Smith else lid_matched[kk] = PETSC_TRUE; // call deleted gids as matched 4908be712e4SBarry Smith } 4918be712e4SBarry Smith /* max edge and pe vecs */ 4928be712e4SBarry Smith PetscCall(MatCreateVecs(cMat, &locMaxEdge, NULL)); 4938be712e4SBarry Smith PetscCall(MatCreateVecs(cMat, &locMaxPE, NULL)); 4948be712e4SBarry Smith /* get 'lghost_pe' & 'lghost_gid' & init. 'lghost_matched' using 'mpimat->lvec' */ 4958be712e4SBarry Smith if (isMPI) { 4968be712e4SBarry Smith Vec vec; 4978be712e4SBarry Smith PetscScalar vval; 4988be712e4SBarry Smith const PetscScalar *buf; 499*c376f201SBarry Smith 5008be712e4SBarry Smith PetscCall(MatCreateVecs(cMat, &vec, NULL)); 5018be712e4SBarry Smith PetscCall(VecGetLocalSize(mpimat->lvec, &num_ghosts)); 5028be712e4SBarry Smith /* lghost_matched */ 503*c376f201SBarry Smith for (PetscInt kk = 0, gid = Istart; kk < nloc; kk++, gid++) { 5048be712e4SBarry Smith PetscScalar vval = lid_matched[kk] ? 1.0 : 0.0; 505*c376f201SBarry Smith 5068be712e4SBarry Smith PetscCall(VecSetValues(vec, 1, &gid, &vval, INSERT_VALUES)); 5078be712e4SBarry Smith } 5088be712e4SBarry Smith PetscCall(VecAssemblyBegin(vec)); 5098be712e4SBarry Smith PetscCall(VecAssemblyEnd(vec)); 5108be712e4SBarry Smith PetscCall(VecScatterBegin(mpimat->Mvctx, vec, mpimat->lvec, INSERT_VALUES, SCATTER_FORWARD)); 5118be712e4SBarry Smith PetscCall(VecScatterEnd(mpimat->Mvctx, vec, mpimat->lvec, INSERT_VALUES, SCATTER_FORWARD)); 5128be712e4SBarry Smith PetscCall(VecGetArrayRead(mpimat->lvec, &buf)); /* get proc ID in 'buf' */ 513*c376f201SBarry Smith PetscCall(PetscMalloc4(num_ghosts, &lghost_matched, num_ghosts, &lghost_pe, num_ghosts, &lghost_gid, num_ghosts, &lghost_max_pe)); 514*c376f201SBarry Smith 515*c376f201SBarry Smith for (PetscInt kk = 0; kk < num_ghosts; kk++) { 5168be712e4SBarry Smith lghost_matched[kk] = (PetscBool)(PetscRealPart(buf[kk]) != 0); // the proc of the ghost for now 5178be712e4SBarry Smith } 5188be712e4SBarry Smith PetscCall(VecRestoreArrayRead(mpimat->lvec, &buf)); 5198be712e4SBarry Smith /* lghost_pe */ 5208be712e4SBarry Smith vval = (PetscScalar)(rank); 521*c376f201SBarry Smith for (PetscInt kk = 0, gid = Istart; kk < nloc; kk++, gid++) PetscCall(VecSetValues(vec, 1, &gid, &vval, INSERT_VALUES)); /* set with GID */ 5228be712e4SBarry Smith PetscCall(VecAssemblyBegin(vec)); 5238be712e4SBarry Smith PetscCall(VecAssemblyEnd(vec)); 5248be712e4SBarry Smith PetscCall(VecScatterBegin(mpimat->Mvctx, vec, mpimat->lvec, INSERT_VALUES, SCATTER_FORWARD)); 5258be712e4SBarry Smith PetscCall(VecScatterEnd(mpimat->Mvctx, vec, mpimat->lvec, INSERT_VALUES, SCATTER_FORWARD)); 5268be712e4SBarry Smith PetscCall(VecGetArrayRead(mpimat->lvec, &buf)); /* get proc ID in 'buf' */ 527*c376f201SBarry Smith for (PetscInt kk = 0; kk < num_ghosts; kk++) lghost_pe[kk] = (PetscMPIInt)PetscRealPart(buf[kk]); // the proc of the ghost for now 5288be712e4SBarry Smith PetscCall(VecRestoreArrayRead(mpimat->lvec, &buf)); 5298be712e4SBarry Smith /* lghost_gid */ 530*c376f201SBarry Smith for (PetscInt kk = 0, gid = Istart; kk < nloc; kk++, gid++) { 5318be712e4SBarry Smith vval = (PetscScalar)(gid); 532*c376f201SBarry Smith 5338be712e4SBarry Smith PetscCall(VecSetValues(vec, 1, &gid, &vval, INSERT_VALUES)); /* set with GID */ 5348be712e4SBarry Smith } 5358be712e4SBarry Smith PetscCall(VecAssemblyBegin(vec)); 5368be712e4SBarry Smith PetscCall(VecAssemblyEnd(vec)); 5378be712e4SBarry Smith PetscCall(VecScatterBegin(mpimat->Mvctx, vec, mpimat->lvec, INSERT_VALUES, SCATTER_FORWARD)); 5388be712e4SBarry Smith PetscCall(VecScatterEnd(mpimat->Mvctx, vec, mpimat->lvec, INSERT_VALUES, SCATTER_FORWARD)); 5398be712e4SBarry Smith PetscCall(VecDestroy(&vec)); 5408be712e4SBarry Smith PetscCall(VecGetArrayRead(mpimat->lvec, &buf)); /* get proc ID in 'lghost_gid' */ 541*c376f201SBarry Smith for (PetscInt kk = 0; kk < num_ghosts; kk++) lghost_gid[kk] = (PetscInt)PetscRealPart(buf[kk]); 5428be712e4SBarry Smith PetscCall(VecRestoreArrayRead(mpimat->lvec, &buf)); 5438be712e4SBarry Smith } 5448be712e4SBarry Smith // get 'comm_procs' (could hoist) 545*c376f201SBarry Smith for (PetscInt kk = 0; kk < REQ_BF_SIZE; kk++) comm_procs[kk] = -1; 5468be712e4SBarry Smith for (ix = 0, ncomm_procs = 0; ix < num_ghosts; ix++) { 5478be712e4SBarry Smith PetscMPIInt proc = lghost_pe[ix], idx = -1; 548*c376f201SBarry Smith 549*c376f201SBarry Smith for (PetscInt k = 0; k < ncomm_procs && idx == -1; k++) 5508be712e4SBarry Smith if (comm_procs[k] == proc) idx = k; 5518be712e4SBarry Smith if (idx == -1) { comm_procs[ncomm_procs++] = proc; } 552*c376f201SBarry Smith PetscCheck(ncomm_procs != REQ_BF_SIZE, PETSC_COMM_SELF, PETSC_ERR_SUP, "Receive request array too small: %d", ncomm_procs); 5538be712e4SBarry Smith } 5548be712e4SBarry Smith /* count edges, compute initial 'locMaxEdge', 'locMaxPE' */ 5558be712e4SBarry Smith nEdges0 = 0; 556*c376f201SBarry Smith for (PetscInt kk = 0, gid = Istart; kk < nloc; kk++, gid++) { 5578be712e4SBarry Smith PetscReal max_e = 0., tt; 5588be712e4SBarry Smith PetscScalar vval; 5598be712e4SBarry Smith PetscInt lid = kk, max_pe = rank, pe, n; 560*c376f201SBarry Smith 5618be712e4SBarry Smith ii = matA->i; 5628be712e4SBarry Smith n = ii[lid + 1] - ii[lid]; 5638e3a54c0SPierre Jolivet aj = PetscSafePointerPlusOffset(matA->j, ii[lid]); 5648e3a54c0SPierre Jolivet ap = PetscSafePointerPlusOffset(matA->a, ii[lid]); 565*c376f201SBarry Smith for (PetscInt jj = 0; jj < n; jj++) { 5668be712e4SBarry Smith PetscInt lidj = aj[jj]; 567*c376f201SBarry Smith 5688be712e4SBarry Smith if ((tt = PetscRealPart(ap[jj])) > threshold && lidj != lid) { 5698be712e4SBarry Smith if (tt > max_e) max_e = tt; 5708be712e4SBarry Smith if (lidj > lid) nEdges0++; 5718be712e4SBarry Smith } 5728be712e4SBarry Smith } 5738be712e4SBarry Smith if ((ix = lid_cprowID[lid]) != -1) { /* if I have any ghost neighbors */ 5748be712e4SBarry Smith ii = matB->compressedrow.i; 5758be712e4SBarry Smith n = ii[ix + 1] - ii[ix]; 5768be712e4SBarry Smith ap = matB->a + ii[ix]; 5778be712e4SBarry Smith aj = matB->j + ii[ix]; 578*c376f201SBarry Smith for (PetscInt jj = 0; jj < n; jj++) { 5798be712e4SBarry Smith if ((tt = PetscRealPart(ap[jj])) > threshold) { 5808be712e4SBarry Smith if (tt > max_e) max_e = tt; 5818be712e4SBarry Smith nEdges0++; 5828be712e4SBarry Smith if ((pe = lghost_pe[aj[jj]]) > max_pe) max_pe = pe; 5838be712e4SBarry Smith } 5848be712e4SBarry Smith } 5858be712e4SBarry Smith } 5868be712e4SBarry Smith vval = max_e; 5878be712e4SBarry Smith PetscCall(VecSetValues(locMaxEdge, 1, &gid, &vval, INSERT_VALUES)); 5888be712e4SBarry Smith vval = (PetscScalar)max_pe; 5898be712e4SBarry Smith PetscCall(VecSetValues(locMaxPE, 1, &gid, &vval, INSERT_VALUES)); 5908be712e4SBarry Smith if (iter == 0 && max_e <= MY_MEPS) { // add BCs to fake aggregate 5918be712e4SBarry Smith lid_matched[lid] = PETSC_TRUE; 5928be712e4SBarry Smith if (bc_agg == -1) { 5938be712e4SBarry Smith bc_agg = lid; 5948be712e4SBarry Smith PetscCall(PetscCDCreate(1, &bc_list)); 5958be712e4SBarry Smith } 5968be712e4SBarry Smith PetscCall(PetscCDRemoveAllAt(agg_llists, lid)); 597*c376f201SBarry Smith PetscCall(PetscCDAppendID(bc_list, 0, Istart + lid)); 5988be712e4SBarry Smith } 5998be712e4SBarry Smith } 6008be712e4SBarry Smith PetscCall(VecAssemblyBegin(locMaxEdge)); 6018be712e4SBarry Smith PetscCall(VecAssemblyEnd(locMaxEdge)); 6028be712e4SBarry Smith PetscCall(VecAssemblyBegin(locMaxPE)); 6038be712e4SBarry Smith PetscCall(VecAssemblyEnd(locMaxPE)); 6048be712e4SBarry Smith /* make 'ghostMaxEdge_max_ew', 'lghost_max_pe' */ 605*c376f201SBarry Smith if (isMPI) { 6068be712e4SBarry Smith const PetscScalar *buf; 607*c376f201SBarry Smith 6088be712e4SBarry Smith PetscCall(VecDuplicate(mpimat->lvec, &ghostMaxEdge)); 6098be712e4SBarry Smith PetscCall(VecScatterBegin(mpimat->Mvctx, locMaxEdge, ghostMaxEdge, INSERT_VALUES, SCATTER_FORWARD)); 6108be712e4SBarry Smith PetscCall(VecScatterEnd(mpimat->Mvctx, locMaxEdge, ghostMaxEdge, INSERT_VALUES, SCATTER_FORWARD)); 6118be712e4SBarry Smith 6128be712e4SBarry Smith PetscCall(VecDuplicate(mpimat->lvec, &ghostMaxPE)); 6138be712e4SBarry Smith PetscCall(VecScatterBegin(mpimat->Mvctx, locMaxPE, ghostMaxPE, INSERT_VALUES, SCATTER_FORWARD)); 6148be712e4SBarry Smith PetscCall(VecScatterEnd(mpimat->Mvctx, locMaxPE, ghostMaxPE, INSERT_VALUES, SCATTER_FORWARD)); 6158be712e4SBarry Smith PetscCall(VecGetArrayRead(ghostMaxPE, &buf)); 616*c376f201SBarry Smith for (PetscInt kk = 0; kk < num_ghosts; kk++) lghost_max_pe[kk] = (PetscMPIInt)PetscRealPart(buf[kk]); // the MAX proc of the ghost now 6178be712e4SBarry Smith PetscCall(VecRestoreArrayRead(ghostMaxPE, &buf)); 6188be712e4SBarry Smith } 6198be712e4SBarry Smith { // make lid_max_pe 6208be712e4SBarry Smith const PetscScalar *buf; 621*c376f201SBarry Smith 6228be712e4SBarry Smith PetscCall(VecGetArrayRead(locMaxPE, &buf)); 623*c376f201SBarry Smith for (PetscInt kk = 0; kk < nloc; kk++) lid_max_pe[kk] = (PetscMPIInt)PetscRealPart(buf[kk]); // the MAX proc of the ghost now 6248be712e4SBarry Smith PetscCall(VecRestoreArrayRead(locMaxPE, &buf)); 6258be712e4SBarry Smith } 6268be712e4SBarry Smith /* setup sorted list of edges, and make 'Edges' */ 6278be712e4SBarry Smith PetscCall(PetscMalloc1(nEdges0, &Edges)); 6288be712e4SBarry Smith nEdges = 0; 629*c376f201SBarry Smith for (PetscInt kk = 0, n; kk < nloc; kk++) { 6308be712e4SBarry Smith const PetscInt lid = kk; 6318be712e4SBarry Smith PetscReal tt; 632*c376f201SBarry Smith 6338be712e4SBarry Smith ii = matA->i; 6348be712e4SBarry Smith n = ii[lid + 1] - ii[lid]; 6358e3a54c0SPierre Jolivet aj = PetscSafePointerPlusOffset(matA->j, ii[lid]); 6368e3a54c0SPierre Jolivet ap = PetscSafePointerPlusOffset(matA->a, ii[lid]); 637*c376f201SBarry Smith for (PetscInt jj = 0; jj < n; jj++) { 6388be712e4SBarry Smith PetscInt lidj = aj[jj]; 639*c376f201SBarry Smith 6408be712e4SBarry Smith if ((tt = PetscRealPart(ap[jj])) > threshold && lidj != lid) { 6418be712e4SBarry Smith if (lidj > lid) { 6428be712e4SBarry Smith Edges[nEdges].lid0 = lid; 643*c376f201SBarry Smith Edges[nEdges].gid1 = lidj + Istart; 6448be712e4SBarry Smith Edges[nEdges].ghost1_idx = -1; 6458be712e4SBarry Smith Edges[nEdges].weight = tt; 6468be712e4SBarry Smith nEdges++; 6478be712e4SBarry Smith } 6488be712e4SBarry Smith } 6498be712e4SBarry Smith } 6508be712e4SBarry Smith if ((ix = lid_cprowID[lid]) != -1) { /* if I have any ghost neighbor */ 6518be712e4SBarry Smith ii = matB->compressedrow.i; 6528be712e4SBarry Smith n = ii[ix + 1] - ii[ix]; 6538be712e4SBarry Smith ap = matB->a + ii[ix]; 6548be712e4SBarry Smith aj = matB->j + ii[ix]; 655*c376f201SBarry Smith for (PetscInt jj = 0; jj < n; jj++) { 6568be712e4SBarry Smith if ((tt = PetscRealPart(ap[jj])) > threshold) { 6578be712e4SBarry Smith Edges[nEdges].lid0 = lid; 6588be712e4SBarry Smith Edges[nEdges].gid1 = lghost_gid[aj[jj]]; 6598be712e4SBarry Smith Edges[nEdges].ghost1_idx = aj[jj]; 6608be712e4SBarry Smith Edges[nEdges].weight = tt; 6618be712e4SBarry Smith nEdges++; 6628be712e4SBarry Smith } 6638be712e4SBarry Smith } 6648be712e4SBarry Smith } 6658be712e4SBarry Smith } 666*c376f201SBarry Smith PetscCheck(nEdges == nEdges0, PETSC_COMM_SELF, PETSC_ERR_SUP, "nEdges != nEdges0: %" PetscInt_FMT " %" PetscInt_FMT, nEdges0, nEdges); 667810441c8SPierre Jolivet if (Edges) qsort(Edges, nEdges, sizeof(Edge), gamg_hem_compare); 6688be712e4SBarry Smith 669*c376f201SBarry Smith PetscCall(PetscInfo(info_is, "[%d] HEM iteration %" PetscInt_FMT " with %" PetscInt_FMT " edges\n", rank, iter, nEdges)); 6708be712e4SBarry Smith 6718be712e4SBarry Smith /* projection matrix */ 6728be712e4SBarry Smith PetscCall(MatCreate(comm, &P)); 6738be712e4SBarry Smith PetscCall(MatSetType(P, MATAIJ)); 6748be712e4SBarry Smith PetscCall(MatSetSizes(P, nloc, nloc, PETSC_DETERMINE, PETSC_DETERMINE)); 6758be712e4SBarry Smith PetscCall(MatMPIAIJSetPreallocation(P, 1, NULL, 1, NULL)); 6768be712e4SBarry Smith PetscCall(MatSeqAIJSetPreallocation(P, 1, NULL)); 6778be712e4SBarry Smith PetscCall(MatSetUp(P)); 6788be712e4SBarry Smith /* process - communicate - process */ 679*c376f201SBarry Smith for (PetscInt sub_it = 0, old_num_edge = 0; /* sub_it < n_sub_its */; /* sub_it++ */) { 6808be712e4SBarry Smith PetscInt nactive_edges = 0, n_act_n[3], gn_act_n[3]; 6818be712e4SBarry Smith PetscMPIInt tag1, tag2; 682*c376f201SBarry Smith 6838be712e4SBarry Smith PetscCall(VecGetArrayRead(locMaxEdge, &lid_max_ew)); 6848be712e4SBarry Smith if (isMPI) { 6858be712e4SBarry Smith PetscCall(VecGetArrayRead(ghostMaxEdge, &lghost_max_ew)); 6868be712e4SBarry Smith PetscCall(PetscCommGetNewTag(comm, &tag1)); 6878be712e4SBarry Smith PetscCall(PetscCommGetNewTag(comm, &tag2)); 6888be712e4SBarry Smith } 689*c376f201SBarry Smith for (PetscInt kk = 0; kk < nEdges; kk++) { 6908be712e4SBarry Smith const Edge *e = &Edges[kk]; 691*c376f201SBarry Smith const PetscInt lid0 = e->lid0, gid1 = e->gid1, ghost1_idx = e->ghost1_idx, gid0 = lid0 + Istart, lid1 = gid1 - Istart; 6928be712e4SBarry Smith PetscBool isOK = PETSC_TRUE, print = PETSC_FALSE; 693*c376f201SBarry Smith 694*c376f201SBarry Smith if (print) 695*c376f201SBarry Smith PetscCall(PetscSynchronizedPrintf(comm, "\t[%d] edge (%" PetscInt_FMT " %" PetscInt_FMT "), %s %s %s\n", rank, gid0, gid1, lid_matched[lid0] ? "true" : "false", (ghost1_idx != -1 && lghost_matched[ghost1_idx]) ? "true" : "false", (ghost1_idx == -1 && lid_matched[lid1]) ? "true" : "false")); 6968be712e4SBarry Smith /* skip if either vertex is matched already */ 6978be712e4SBarry Smith if (lid_matched[lid0] || (ghost1_idx != -1 && lghost_matched[ghost1_idx]) || (ghost1_idx == -1 && lid_matched[lid1])) continue; 6988be712e4SBarry Smith 6998be712e4SBarry Smith nactive_edges++; 7008be712e4SBarry Smith PetscCheck(PetscRealPart(lid_max_ew[lid0]) >= e->weight - MY_MEPS, PETSC_COMM_SELF, PETSC_ERR_SUP, "edge weight %e > max %e", (double)e->weight, (double)PetscRealPart(lid_max_ew[lid0])); 701*c376f201SBarry Smith if (print) PetscCall(PetscSynchronizedPrintf(comm, "\t[%d] active edge (%" PetscInt_FMT " %" PetscInt_FMT "), diff0 = %10.4e\n", rank, gid0, gid1, (double)(PetscRealPart(lid_max_ew[lid0]) - (double)e->weight))); 7028be712e4SBarry Smith // smaller edge, lid_max_ew get updated - e0 7038be712e4SBarry Smith if (PetscRealPart(lid_max_ew[lid0]) > e->weight + MY_MEPS) { 7048be712e4SBarry Smith if (print) 705*c376f201SBarry Smith PetscCall(PetscSynchronizedPrintf(comm, "\t\t[%d] 1) e0 SKIPPING small edge %20.14e edge (%" PetscInt_FMT " %" PetscInt_FMT "), diff = %10.4e to proc %d. max = %20.14e, w = %20.14e\n", rank, (double)e->weight, gid0, gid1, (double)(PetscRealPart(lid_max_ew[lid0]) - e->weight), ghost1_idx != -1 ? lghost_pe[ghost1_idx] : rank, (double)PetscRealPart(lid_max_ew[lid0]), 7068be712e4SBarry Smith (double)e->weight)); 7078be712e4SBarry Smith continue; // we are basically filter edges here 7088be712e4SBarry Smith } 7098be712e4SBarry Smith // e1 - local 7108be712e4SBarry Smith if (ghost1_idx == -1) { 7118be712e4SBarry Smith if (PetscRealPart(lid_max_ew[lid1]) > e->weight + MY_MEPS) { 7128be712e4SBarry Smith if (print) 713*c376f201SBarry Smith PetscCall(PetscSynchronizedPrintf(comm, "\t\t%c[%d] 2) e1 SKIPPING small local edge %20.14e edge (%" PetscInt_FMT " %" PetscInt_FMT "), diff = %10.4e\n", ghost1_idx != -1 ? '\t' : ' ', rank, (double)e->weight, gid0, gid1, (double)(PetscRealPart(lid_max_ew[lid1]) - e->weight))); 7148be712e4SBarry Smith continue; // we are basically filter edges here 7158be712e4SBarry Smith } 7168be712e4SBarry Smith } else { // e1 - ghost 7178be712e4SBarry Smith /* see if edge might get matched on other proc */ 7188be712e4SBarry Smith PetscReal g_max_e1 = PetscRealPart(lghost_max_ew[ghost1_idx]); 719*c376f201SBarry Smith 7208be712e4SBarry Smith if (print) 721*c376f201SBarry Smith PetscCall(PetscSynchronizedPrintf(comm, "\t\t\t[%d] CHECK GHOST e1, edge (%" PetscInt_FMT " %" PetscInt_FMT "), E0 MAX EDGE WEIGHT = %10.4e, EDGE WEIGHT = %10.4e, diff1 = %10.4e, ghost proc %d with max pe %d on e0 and %d on e1\n", rank, gid0, gid1, (double)PetscRealPart(lid_max_ew[lid0]), 722*c376f201SBarry Smith (double)e->weight, (double)(PetscRealPart(lghost_max_ew[ghost1_idx]) - e->weight), lghost_pe[ghost1_idx], lid_max_pe[lid0], lghost_max_pe[ghost1_idx])); 7238be712e4SBarry Smith if (g_max_e1 > e->weight + MY_MEPS) { 724*c376f201SBarry Smith /* PetscCall(PetscSynchronizedPrintf(comm,"\t\t\t\t[%d] 3) ghost e1 SKIPPING small edge (%d %d), diff = %10.4e from proc %d with max pe %d. max = %20.14e, w = %20.14e\n", rank, gid0, gid1, g_max_e1 - e->weight, lghost_pe[ghost1_idx], lghost_max_pe[ghost1_idx], g_max_e1, e->weight )); */ 7258be712e4SBarry Smith continue; 726452acf8bSMark Adams } else if (g_max_e1 >= e->weight - MY_MEPS && lghost_pe[ghost1_idx] > rank) { // is 'lghost_max_pe[ghost1_idx] > rank' needed? 7278be712e4SBarry Smith /* check for max_ea == to this edge and larger processor that will deal with this */ 7288be712e4SBarry Smith if (print) 729*c376f201SBarry Smith PetscCall(PetscSynchronizedPrintf(comm, "\t\t\t[%d] ghost e1 SKIPPING EQUAL (%" PetscInt_FMT " %" PetscInt_FMT "), diff = %10.4e from larger proc %d with max pe %d. max = %20.14e, w = %20.14e\n", rank, gid0, gid1, 730*c376f201SBarry Smith (double)(PetscRealPart(lid_max_ew[lid0]) - (double)e->weight), lghost_pe[ghost1_idx], lghost_max_pe[ghost1_idx], (double)g_max_e1, (double)e->weight)); 7318be712e4SBarry Smith isOK = PETSC_FALSE; // this guy could delete me 7328be712e4SBarry Smith continue; 7338be712e4SBarry Smith } else { 734*c376f201SBarry Smith /* PetscCall(PetscSynchronizedPrintf(comm,"\t[%d] Edge (%d %d) passes gid0 tests, diff = %10.4e from proc %d with max pe %d. max = %20.14e, w = %20.14e\n", rank, gid0, gid1, g_max_e1 - e->weight, lghost_pe[ghost1_idx], lghost_max_pe[ghost1_idx], g_max_e1, e->weight )); */ 7358be712e4SBarry Smith } 7368be712e4SBarry Smith } 7378be712e4SBarry Smith /* check ghost for v0 */ 7388be712e4SBarry Smith if (isOK) { 7398be712e4SBarry Smith PetscReal max_e, ew; 740*c376f201SBarry Smith 7418be712e4SBarry Smith if ((ix = lid_cprowID[lid0]) != -1) { /* if I have any ghost neighbors */ 742*c376f201SBarry Smith PetscInt n; 743*c376f201SBarry Smith 7448be712e4SBarry Smith ii = matB->compressedrow.i; 7458be712e4SBarry Smith n = ii[ix + 1] - ii[ix]; 7468be712e4SBarry Smith ap = matB->a + ii[ix]; 7478be712e4SBarry Smith aj = matB->j + ii[ix]; 748*c376f201SBarry Smith for (PetscInt jj = 0; jj < n && isOK; jj++) { 7498be712e4SBarry Smith PetscInt lidj = aj[jj]; 750*c376f201SBarry Smith 7518be712e4SBarry Smith if (lghost_matched[lidj]) continue; 7528be712e4SBarry Smith ew = PetscRealPart(ap[jj]); 7538be712e4SBarry Smith if (ew <= threshold) continue; 7548be712e4SBarry Smith max_e = PetscRealPart(lghost_max_ew[lidj]); 755*c376f201SBarry Smith 7568be712e4SBarry Smith /* check for max_e == to this edge and larger processor that will deal with this */ 757452acf8bSMark Adams if (ew >= PetscRealPart(lid_max_ew[lid0]) - MY_MEPS && lghost_max_pe[lidj] > rank) isOK = PETSC_FALSE; 758*c376f201SBarry Smith PetscCheck(ew <= max_e + MY_MEPS, PETSC_COMM_SELF, PETSC_ERR_SUP, "edge weight %e > max %e. ncols = %" PetscInt_FMT ", gid0 = %" PetscInt_FMT ", gid1 = %" PetscInt_FMT, (double)PetscRealPart(ew), (double)PetscRealPart(max_e), n, lid0 + Istart, lghost_gid[lidj]); 7598be712e4SBarry Smith if (print) 760*c376f201SBarry Smith PetscCall(PetscSynchronizedPrintf(comm, "\t\t\t\t[%d] e0: looked at ghost adj (%" PetscInt_FMT " %" PetscInt_FMT "), diff = %10.4e, ghost on proc %d (max %d). isOK = %d, %d %d %d; ew = %e, lid0 max ew = %e, diff = %e, eps = %e\n", rank, gid0, lghost_gid[lidj], (double)(max_e - ew), lghost_pe[lidj], lghost_max_pe[lidj], isOK, (double)(ew) >= (double)(max_e - MY_MEPS), ew >= PetscRealPart(lid_max_ew[lid0]) - MY_MEPS, lghost_pe[lidj] > rank, (double)ew, (double)PetscRealPart(lid_max_ew[lid0]), (double)(ew - PetscRealPart(lid_max_ew[lid0])), (double)MY_MEPS)); 7618be712e4SBarry Smith } 762*c376f201SBarry Smith if (!isOK && print) PetscCall(PetscSynchronizedPrintf(comm, "\t\t[%d] skip edge (%" PetscInt_FMT " %" PetscInt_FMT ") from ghost inspection\n", rank, gid0, gid1)); 7638be712e4SBarry Smith } 7648be712e4SBarry Smith /* check local v1 */ 7658be712e4SBarry Smith if (ghost1_idx == -1) { 7668be712e4SBarry Smith if ((ix = lid_cprowID[lid1]) != -1) { /* if I have any ghost neighbors */ 767*c376f201SBarry Smith PetscInt n; 768*c376f201SBarry Smith 7698be712e4SBarry Smith ii = matB->compressedrow.i; 7708be712e4SBarry Smith n = ii[ix + 1] - ii[ix]; 7718be712e4SBarry Smith ap = matB->a + ii[ix]; 7728be712e4SBarry Smith aj = matB->j + ii[ix]; 773*c376f201SBarry Smith for (PetscInt jj = 0; jj < n && isOK; jj++) { 7748be712e4SBarry Smith PetscInt lidj = aj[jj]; 775*c376f201SBarry Smith 7768be712e4SBarry Smith if (lghost_matched[lidj]) continue; 7778be712e4SBarry Smith ew = PetscRealPart(ap[jj]); 7788be712e4SBarry Smith if (ew <= threshold) continue; 7798be712e4SBarry Smith max_e = PetscRealPart(lghost_max_ew[lidj]); 7808be712e4SBarry Smith /* check for max_e == to this edge and larger processor that will deal with this */ 781452acf8bSMark Adams if (ew >= PetscRealPart(lid_max_ew[lid1]) - MY_MEPS && lghost_max_pe[lidj] > rank) isOK = PETSC_FALSE; 7828be712e4SBarry Smith PetscCheck(ew <= max_e + MY_MEPS, PETSC_COMM_SELF, PETSC_ERR_SUP, "edge weight %e > max %e", (double)PetscRealPart(ew), (double)PetscRealPart(max_e)); 7838be712e4SBarry Smith if (print) 784*c376f201SBarry Smith PetscCall(PetscSynchronizedPrintf(comm, "\t\t\t\t\t[%d] e1: looked at ghost adj (%" PetscInt_FMT " %" PetscInt_FMT "), diff = %10.4e, ghost on proc %d (max %d)\n", rank, gid0, lghost_gid[lidj], (double)(max_e - ew), lghost_pe[lidj], lghost_max_pe[lidj])); 7858be712e4SBarry Smith } 7868be712e4SBarry Smith } 787*c376f201SBarry Smith if (!isOK && print) PetscCall(PetscSynchronizedPrintf(comm, "\t\t[%d] skip edge (%" PetscInt_FMT " %" PetscInt_FMT ") from ghost inspection\n", rank, gid0, gid1)); 7888be712e4SBarry Smith } 7898be712e4SBarry Smith } 7908be712e4SBarry Smith PetscReal e1_max_w = (ghost1_idx == -1 ? PetscRealPart(lid_max_ew[lid0]) : PetscRealPart(lghost_max_ew[ghost1_idx])); 7918be712e4SBarry Smith if (print) 792*c376f201SBarry Smith PetscCall(PetscSynchronizedPrintf(comm, "\t[%d] MATCHING (%" PetscInt_FMT " %" PetscInt_FMT ") e1 max weight = %e, e1 weight diff %e, %s. isOK = %d\n", rank, gid0, gid1, (double)e1_max_w, (double)(e1_max_w - e->weight), ghost1_idx == -1 ? "local" : "ghost", isOK)); 7938be712e4SBarry Smith /* do it */ 7948be712e4SBarry Smith if (isOK) { 7958be712e4SBarry Smith if (ghost1_idx == -1) { 796*c376f201SBarry Smith PetscCheck(!lid_matched[lid1], PETSC_COMM_SELF, PETSC_ERR_SUP, "local %" PetscInt_FMT " is matched", gid1); 7978be712e4SBarry Smith lid_matched[lid1] = PETSC_TRUE; /* keep track of what we've done this round */ 7988be712e4SBarry Smith PetscCall(PetscCDMoveAppend(agg_llists, lid0, lid1)); // takes lid1's list and appends to lid0's 7998be712e4SBarry Smith } else { 8008be712e4SBarry Smith /* add gid1 to list of ghost deleted by me -- I need their children */ 8018be712e4SBarry Smith PetscMPIInt proc = lghost_pe[ghost1_idx]; 802*c376f201SBarry Smith PetscCheck(!lghost_matched[ghost1_idx], PETSC_COMM_SELF, PETSC_ERR_SUP, "ghost %" PetscInt_FMT " is matched", lghost_gid[ghost1_idx]); 8038be712e4SBarry Smith lghost_matched[ghost1_idx] = PETSC_TRUE; 8048be712e4SBarry Smith PetscCall(PetscCDAppendID(ghost_deleted_list, proc, ghost1_idx)); /* cache to send messages */ 8058be712e4SBarry Smith PetscCall(PetscCDAppendID(ghost_deleted_list, proc, lid0)); 8068be712e4SBarry Smith } 8078be712e4SBarry Smith lid_matched[lid0] = PETSC_TRUE; /* keep track of what we've done this round */ 8088be712e4SBarry Smith /* set projection */ 8098be712e4SBarry Smith PetscCall(MatSetValues(P, 1, &gid0, 1, &gid0, &one, INSERT_VALUES)); 8108be712e4SBarry Smith PetscCall(MatSetValues(P, 1, &gid1, 1, &gid0, &one, INSERT_VALUES)); 811*c376f201SBarry Smith //PetscCall(PetscPrintf(comm,"\t %" PetscInt_FMT ".%" PetscInt_FMT ") match active EDGE %" PetscInt_FMT " : (%" PetscInt_FMT " %" PetscInt_FMT ")\n",iter,sub_it, nactive_edges, gid0, gid1)); 8128be712e4SBarry Smith } /* matched */ 8138be712e4SBarry Smith } /* edge loop */ 814452acf8bSMark Adams PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT)); 8158be712e4SBarry Smith if (isMPI) PetscCall(VecRestoreArrayRead(ghostMaxEdge, &lghost_max_ew)); 8168be712e4SBarry Smith PetscCall(VecRestoreArrayRead(locMaxEdge, &lid_max_ew)); 8178be712e4SBarry Smith // count active for test, latter, update deleted ghosts 8188be712e4SBarry Smith n_act_n[0] = nactive_edges; 8198be712e4SBarry Smith if (ghost_deleted_list) PetscCall(PetscCDCount(ghost_deleted_list, &n_act_n[2])); 8208be712e4SBarry Smith else n_act_n[2] = 0; 8218be712e4SBarry Smith PetscCall(PetscCDCount(agg_llists, &n_act_n[1])); 8228be712e4SBarry Smith PetscCall(MPIU_Allreduce(n_act_n, gn_act_n, 3, MPIU_INT, MPI_SUM, comm)); 823*c376f201SBarry Smith PetscCall(PetscInfo(info_is, "[%d] %" PetscInt_FMT ".%" PetscInt_FMT ") nactive edges=%" PetscInt_FMT ", ncomm_procs=%d, nEdges=%" PetscInt_FMT ", %" PetscInt_FMT " deleted ghosts, N=%" PetscInt_FMT "\n", rank, iter, sub_it, gn_act_n[0], ncomm_procs, nEdges, gn_act_n[2], gn_act_n[1])); 8248be712e4SBarry Smith /* deal with deleted ghost */ 8258be712e4SBarry Smith if (isMPI) { 8268be712e4SBarry Smith PetscCDIntNd *pos; 8278be712e4SBarry Smith PetscInt *sbuffs1[REQ_BF_SIZE], ndel; 8288be712e4SBarry Smith PetscInt *sbuffs2[REQ_BF_SIZE]; 8298be712e4SBarry Smith MPI_Status status; 8308be712e4SBarry Smith 8318be712e4SBarry Smith /* send deleted ghosts */ 832*c376f201SBarry Smith for (PetscInt proc_idx = 0; proc_idx < ncomm_procs; proc_idx++) { 8338be712e4SBarry Smith const PetscMPIInt proc = comm_procs[proc_idx]; 8348be712e4SBarry Smith PetscInt *sbuff, *pt, scount; 8358be712e4SBarry Smith MPI_Request *request; 836*c376f201SBarry Smith 8378be712e4SBarry Smith /* count ghosts */ 8388be712e4SBarry Smith PetscCall(PetscCDCountAt(ghost_deleted_list, proc, &ndel)); 8398be712e4SBarry Smith ndel /= 2; // two entries for each proc 8408be712e4SBarry Smith scount = 2 + 2 * ndel; 8418be712e4SBarry Smith PetscCall(PetscMalloc1(scount + request_size, &sbuff)); 8428be712e4SBarry Smith /* save requests */ 8438be712e4SBarry Smith sbuffs1[proc_idx] = sbuff; 8448be712e4SBarry Smith request = (MPI_Request *)sbuff; 8458be712e4SBarry Smith sbuff = pt = (PetscInt *)(sbuff + request_size); 8468be712e4SBarry Smith /* write [ndel, proc, n*[gid1,gid0] */ 8478be712e4SBarry Smith *pt++ = ndel; // number of deleted to send 8488be712e4SBarry Smith *pt++ = rank; // proc (not used) 8498be712e4SBarry Smith PetscCall(PetscCDGetHeadPos(ghost_deleted_list, proc, &pos)); 8508be712e4SBarry Smith while (pos) { 8518be712e4SBarry Smith PetscInt lid0, ghost_idx, gid1; 852*c376f201SBarry Smith 8538be712e4SBarry Smith PetscCall(PetscCDIntNdGetID(pos, &ghost_idx)); 8548be712e4SBarry Smith gid1 = lghost_gid[ghost_idx]; 8558be712e4SBarry Smith PetscCall(PetscCDGetNextPos(ghost_deleted_list, proc, &pos)); 8568be712e4SBarry Smith PetscCall(PetscCDIntNdGetID(pos, &lid0)); 8578be712e4SBarry Smith PetscCall(PetscCDGetNextPos(ghost_deleted_list, proc, &pos)); 8588be712e4SBarry Smith *pt++ = gid1; 859*c376f201SBarry Smith *pt++ = lid0 + Istart; // gid0 8608be712e4SBarry Smith } 861*c376f201SBarry Smith PetscCheck(pt - sbuff == (ptrdiff_t)scount, PETSC_COMM_SELF, PETSC_ERR_SUP, "sbuff-pt != scount: %zu", (pt - sbuff)); 8628be712e4SBarry Smith /* MPI_Isend: tag1 [ndel, proc, n*[gid1,gid0] ] */ 8638be712e4SBarry Smith PetscCallMPI(MPI_Isend(sbuff, scount, MPIU_INT, proc, tag1, comm, request)); 8648be712e4SBarry Smith PetscCall(PetscCDRemoveAllAt(ghost_deleted_list, proc)); // done with this list 8658be712e4SBarry Smith } 8668be712e4SBarry Smith /* receive deleted, send back partial aggregates, clear lists */ 867*c376f201SBarry Smith for (PetscInt proc_idx = 0; proc_idx < ncomm_procs; proc_idx++) { 8688be712e4SBarry Smith PetscCallMPI(MPI_Probe(comm_procs[proc_idx] /* MPI_ANY_SOURCE */, tag1, comm, &status)); 8698be712e4SBarry Smith { 8708be712e4SBarry Smith PetscInt *pt, *pt2, *pt3, *sbuff, tmp; 8718be712e4SBarry Smith MPI_Request *request; 872*c376f201SBarry Smith PetscMPIInt rcount, scount; 8738be712e4SBarry Smith const PetscMPIInt proc = status.MPI_SOURCE; 874*c376f201SBarry Smith 8758be712e4SBarry Smith PetscCallMPI(MPI_Get_count(&status, MPIU_INT, &rcount)); 8768be712e4SBarry Smith if (rcount > rbuff_sz) { 8778be712e4SBarry Smith if (rbuff) PetscCall(PetscFree(rbuff)); 8788be712e4SBarry Smith PetscCall(PetscMalloc1(rcount, &rbuff)); 8798be712e4SBarry Smith rbuff_sz = rcount; 8808be712e4SBarry Smith } 8818be712e4SBarry Smith /* MPI_Recv: tag1 [ndel, proc, ndel*[gid1,gid0] ] */ 8828be712e4SBarry Smith PetscCallMPI(MPI_Recv(rbuff, rcount, MPIU_INT, proc, tag1, comm, &status)); 8838be712e4SBarry Smith /* read and count sends *[lid0, n, n*[gid] ] */ 8848be712e4SBarry Smith pt = rbuff; 8858be712e4SBarry Smith scount = 0; 8868be712e4SBarry Smith ndel = *pt++; // number of deleted to recv 8878be712e4SBarry Smith tmp = *pt++; // proc (not used) 8888be712e4SBarry Smith while (ndel--) { 889*c376f201SBarry Smith PetscInt gid1 = *pt++, lid1 = gid1 - Istart; 8908be712e4SBarry Smith int gh_gid0 = *pt++; // gid on other proc (not used here to count) 891*c376f201SBarry Smith 892*c376f201SBarry Smith PetscCheck(lid1 >= 0 && lid1 < nloc, PETSC_COMM_SELF, PETSC_ERR_SUP, "received ghost deleted %" PetscInt_FMT, gid1); 893*c376f201SBarry Smith PetscCheck(!lid_matched[lid1], PETSC_COMM_SELF, PETSC_ERR_PLIB, "%" PetscInt_FMT ") received matched local gid %" PetscInt_FMT ",%d, with ghost (lid) %" PetscInt_FMT " from proc %d", sub_it, gid1, gh_gid0, tmp, proc); 8948be712e4SBarry Smith lid_matched[lid1] = PETSC_TRUE; /* keep track of what we've done this round */ 8958be712e4SBarry Smith PetscCall(PetscCDCountAt(agg_llists, lid1, &tmp)); // n 896*c376f201SBarry Smith /* PetscCheck(tmp == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "sending %" PetscInt_FMT " (!= 1) size aggregate. gid-0 %" PetscInt_FMT ", from %d (gid-1 %" PetscInt_FMT ")", tmp, gid, proc, gh_gid0); */ 8978be712e4SBarry Smith scount += tmp + 2; // lid0, n, n*[gid] 8988be712e4SBarry Smith } 899*c376f201SBarry Smith PetscCheck((pt - rbuff) == (ptrdiff_t)rcount, PETSC_COMM_SELF, PETSC_ERR_SUP, "receive buffer size != num read: %zu; rcount: %d", pt - rbuff, rcount); 9008be712e4SBarry Smith /* send tag2: *[gid0, n, n*[gid] ] */ 9018be712e4SBarry Smith PetscCall(PetscMalloc1(scount + request_size, &sbuff)); 9028be712e4SBarry Smith sbuffs2[proc_idx] = sbuff; /* cache request */ 9038be712e4SBarry Smith request = (MPI_Request *)sbuff; 9048be712e4SBarry Smith pt2 = sbuff = sbuff + request_size; 9058be712e4SBarry Smith // read again: n, proc, n*[gid1,gid0] 9068be712e4SBarry Smith pt = rbuff; 9078be712e4SBarry Smith ndel = *pt++; 9088be712e4SBarry Smith tmp = *pt++; // proc (not used) 9098be712e4SBarry Smith while (ndel--) { 910*c376f201SBarry Smith PetscInt gid1 = *pt++, lid1 = gid1 - Istart, gh_gid0 = *pt++; 911*c376f201SBarry Smith 9128be712e4SBarry Smith /* write [gid0, aggSz, aggSz[gid] ] */ 9138be712e4SBarry Smith *pt2++ = gh_gid0; 9148be712e4SBarry Smith pt3 = pt2++; /* save pointer for later */ 9158be712e4SBarry Smith PetscCall(PetscCDGetHeadPos(agg_llists, lid1, &pos)); 9168be712e4SBarry Smith while (pos) { 9178be712e4SBarry Smith PetscInt gid; 918*c376f201SBarry Smith 9198be712e4SBarry Smith PetscCall(PetscCDIntNdGetID(pos, &gid)); 9208be712e4SBarry Smith PetscCall(PetscCDGetNextPos(agg_llists, lid1, &pos)); 9218be712e4SBarry Smith *pt2++ = gid; 9228be712e4SBarry Smith } 9238be712e4SBarry Smith *pt3 = (pt2 - pt3) - 1; 9248be712e4SBarry Smith /* clear list */ 9258be712e4SBarry Smith PetscCall(PetscCDRemoveAllAt(agg_llists, lid1)); 9268be712e4SBarry Smith } 927*c376f201SBarry Smith PetscCheck((pt2 - sbuff) == (ptrdiff_t)scount, PETSC_COMM_SELF, PETSC_ERR_SUP, "buffer size != num write: %zu %d", pt2 - sbuff, scount); 9288be712e4SBarry Smith /* MPI_Isend: requested data tag2 *[lid0, n, n*[gid1] ] */ 9298be712e4SBarry Smith PetscCallMPI(MPI_Isend(sbuff, scount, MPIU_INT, proc, tag2, comm, request)); 9308be712e4SBarry Smith } 9318be712e4SBarry Smith } // proc_idx 9328be712e4SBarry Smith /* receive tag2 *[gid0, n, n*[gid] ] */ 933*c376f201SBarry Smith for (PetscMPIInt proc_idx = 0; proc_idx < ncomm_procs; proc_idx++) { 9348be712e4SBarry Smith PetscMPIInt proc; 9358be712e4SBarry Smith PetscInt *pt; 9368be712e4SBarry Smith int rcount; 937*c376f201SBarry Smith 9388be712e4SBarry Smith PetscCallMPI(MPI_Probe(comm_procs[proc_idx] /* MPI_ANY_SOURCE */, tag2, comm, &status)); 9398be712e4SBarry Smith PetscCallMPI(MPI_Get_count(&status, MPIU_INT, &rcount)); 9408be712e4SBarry Smith if (rcount > rbuff_sz) { 9418be712e4SBarry Smith if (rbuff) PetscCall(PetscFree(rbuff)); 9428be712e4SBarry Smith PetscCall(PetscMalloc1(rcount, &rbuff)); 9438be712e4SBarry Smith rbuff_sz = rcount; 9448be712e4SBarry Smith } 9458be712e4SBarry Smith proc = status.MPI_SOURCE; 9468be712e4SBarry Smith /* MPI_Recv: tag1 [n, proc, n*[gid1,lid0] ] */ 9478be712e4SBarry Smith PetscCallMPI(MPI_Recv(rbuff, rcount, MPIU_INT, proc, tag2, comm, &status)); 9488be712e4SBarry Smith pt = rbuff; 9498be712e4SBarry Smith while (pt - rbuff < rcount) { 9508be712e4SBarry Smith PetscInt gid0 = *pt++, n = *pt++; 951*c376f201SBarry Smith 9528be712e4SBarry Smith while (n--) { 9538be712e4SBarry Smith PetscInt gid1 = *pt++; 954*c376f201SBarry Smith 955*c376f201SBarry Smith PetscCall(PetscCDAppendID(agg_llists, gid0 - Istart, gid1)); 9568be712e4SBarry Smith } 9578be712e4SBarry Smith } 958*c376f201SBarry Smith PetscCheck((pt - rbuff) == (ptrdiff_t)rcount, PETSC_COMM_SELF, PETSC_ERR_SUP, "recv buffer size != num read: %zu %d", pt - rbuff, rcount); 9598be712e4SBarry Smith } 9608be712e4SBarry Smith /* wait for tag1 isends */ 961*c376f201SBarry Smith for (PetscMPIInt proc_idx = 0; proc_idx < ncomm_procs; proc_idx++) { 962*c376f201SBarry Smith MPI_Request *request = (MPI_Request *)sbuffs1[proc_idx]; 963*c376f201SBarry Smith 9648be712e4SBarry Smith PetscCallMPI(MPI_Wait(request, &status)); 9658be712e4SBarry Smith PetscCall(PetscFree(sbuffs1[proc_idx])); 9668be712e4SBarry Smith } 9678be712e4SBarry Smith /* wait for tag2 isends */ 968*c376f201SBarry Smith for (PetscMPIInt proc_idx = 0; proc_idx < ncomm_procs; proc_idx++) { 9698be712e4SBarry Smith MPI_Request *request = (MPI_Request *)sbuffs2[proc_idx]; 970*c376f201SBarry Smith 9718be712e4SBarry Smith PetscCallMPI(MPI_Wait(request, &status)); 9728be712e4SBarry Smith PetscCall(PetscFree(sbuffs2[proc_idx])); 9738be712e4SBarry Smith } 9748be712e4SBarry Smith } /* MPI */ 9758be712e4SBarry Smith /* set 'lghost_matched' - use locMaxEdge, ghostMaxEdge (recomputed next) */ 9768be712e4SBarry Smith if (isMPI) { 9778be712e4SBarry Smith const PetscScalar *sbuff; 978*c376f201SBarry Smith 979*c376f201SBarry Smith for (PetscInt kk = 0, gid = Istart; kk < nloc; kk++, gid++) { 9808be712e4SBarry Smith PetscScalar vval = lid_matched[kk] ? 1.0 : 0.0; 981*c376f201SBarry Smith 9828be712e4SBarry Smith PetscCall(VecSetValues(locMaxEdge, 1, &gid, &vval, INSERT_VALUES)); /* set with GID */ 9838be712e4SBarry Smith } 9848be712e4SBarry Smith PetscCall(VecAssemblyBegin(locMaxEdge)); 9858be712e4SBarry Smith PetscCall(VecAssemblyEnd(locMaxEdge)); 9868be712e4SBarry Smith PetscCall(VecScatterBegin(mpimat->Mvctx, locMaxEdge, ghostMaxEdge, INSERT_VALUES, SCATTER_FORWARD)); 9878be712e4SBarry Smith PetscCall(VecScatterEnd(mpimat->Mvctx, locMaxEdge, ghostMaxEdge, INSERT_VALUES, SCATTER_FORWARD)); 9888be712e4SBarry Smith PetscCall(VecGetArrayRead(ghostMaxEdge, &sbuff)); 989*c376f201SBarry Smith for (PetscInt kk = 0; kk < num_ghosts; kk++) { lghost_matched[kk] = (PetscBool)(PetscRealPart(sbuff[kk]) != 0.0); } 9908be712e4SBarry Smith PetscCall(VecRestoreArrayRead(ghostMaxEdge, &sbuff)); 9918be712e4SBarry Smith } 9928be712e4SBarry Smith /* compute 'locMaxEdge' inside sub iteration b/c max weight can drop as neighbors are matched */ 993*c376f201SBarry Smith for (PetscInt kk = 0, gid = Istart; kk < nloc; kk++, gid++) { 9948be712e4SBarry Smith PetscReal max_e = 0., tt; 9958be712e4SBarry Smith PetscScalar vval; 9968be712e4SBarry Smith const PetscInt lid = kk; 997*c376f201SBarry Smith PetscMPIInt max_pe = rank, pe, n; 998*c376f201SBarry Smith 9998be712e4SBarry Smith ii = matA->i; 10008be712e4SBarry Smith n = ii[lid + 1] - ii[lid]; 10018e3a54c0SPierre Jolivet aj = PetscSafePointerPlusOffset(matA->j, ii[lid]); 10028e3a54c0SPierre Jolivet ap = PetscSafePointerPlusOffset(matA->a, ii[lid]); 1003*c376f201SBarry Smith for (PetscInt jj = 0; jj < n; jj++) { 10048be712e4SBarry Smith PetscInt lidj = aj[jj]; 1005*c376f201SBarry Smith 10068be712e4SBarry Smith if (lid_matched[lidj]) continue; /* this is new - can change local max */ 10078be712e4SBarry Smith if (lidj != lid && PetscRealPart(ap[jj]) > max_e) max_e = PetscRealPart(ap[jj]); 10088be712e4SBarry Smith } 10098be712e4SBarry Smith if (lid_cprowID && (ix = lid_cprowID[lid]) != -1) { /* if I have any ghost neighbors */ 10108be712e4SBarry Smith ii = matB->compressedrow.i; 10118be712e4SBarry Smith n = ii[ix + 1] - ii[ix]; 10128be712e4SBarry Smith ap = matB->a + ii[ix]; 10138be712e4SBarry Smith aj = matB->j + ii[ix]; 1014*c376f201SBarry Smith for (PetscInt jj = 0; jj < n; jj++) { 10158be712e4SBarry Smith PetscInt lidj = aj[jj]; 1016*c376f201SBarry Smith 10178be712e4SBarry Smith if (lghost_matched[lidj]) continue; 10188be712e4SBarry Smith if ((tt = PetscRealPart(ap[jj])) > max_e) max_e = tt; 10198be712e4SBarry Smith } 10208be712e4SBarry Smith } 10218be712e4SBarry Smith vval = (PetscScalar)max_e; 10228be712e4SBarry Smith PetscCall(VecSetValues(locMaxEdge, 1, &gid, &vval, INSERT_VALUES)); /* set with GID */ 10238be712e4SBarry Smith // max PE with max edge 10248be712e4SBarry Smith if (lid_cprowID && (ix = lid_cprowID[lid]) != -1) { /* if I have any ghost neighbors */ 10258be712e4SBarry Smith ii = matB->compressedrow.i; 10268be712e4SBarry Smith n = ii[ix + 1] - ii[ix]; 10278be712e4SBarry Smith ap = matB->a + ii[ix]; 10288be712e4SBarry Smith aj = matB->j + ii[ix]; 1029*c376f201SBarry Smith for (PetscInt jj = 0; jj < n; jj++) { 10308be712e4SBarry Smith PetscInt lidj = aj[jj]; 1031*c376f201SBarry Smith 10328be712e4SBarry Smith if (lghost_matched[lidj]) continue; 10338be712e4SBarry Smith if ((pe = lghost_pe[aj[jj]]) > max_pe && PetscRealPart(ap[jj]) >= max_e - MY_MEPS) { max_pe = pe; } 10348be712e4SBarry Smith } 10358be712e4SBarry Smith } 10368be712e4SBarry Smith vval = (PetscScalar)max_pe; 10378be712e4SBarry Smith PetscCall(VecSetValues(locMaxPE, 1, &gid, &vval, INSERT_VALUES)); 10388be712e4SBarry Smith } 10398be712e4SBarry Smith PetscCall(VecAssemblyBegin(locMaxEdge)); 10408be712e4SBarry Smith PetscCall(VecAssemblyEnd(locMaxEdge)); 10418be712e4SBarry Smith PetscCall(VecAssemblyBegin(locMaxPE)); 10428be712e4SBarry Smith PetscCall(VecAssemblyEnd(locMaxPE)); 10438be712e4SBarry Smith /* compute 'lghost_max_ew' and 'lghost_max_pe' to get ready for next iteration*/ 10448be712e4SBarry Smith if (isMPI) { 10458be712e4SBarry Smith const PetscScalar *buf; 1046*c376f201SBarry Smith 10478be712e4SBarry Smith PetscCall(VecScatterBegin(mpimat->Mvctx, locMaxEdge, ghostMaxEdge, INSERT_VALUES, SCATTER_FORWARD)); 10488be712e4SBarry Smith PetscCall(VecScatterEnd(mpimat->Mvctx, locMaxEdge, ghostMaxEdge, INSERT_VALUES, SCATTER_FORWARD)); 10498be712e4SBarry Smith PetscCall(VecScatterBegin(mpimat->Mvctx, locMaxPE, ghostMaxPE, INSERT_VALUES, SCATTER_FORWARD)); 10508be712e4SBarry Smith PetscCall(VecScatterEnd(mpimat->Mvctx, locMaxPE, ghostMaxPE, INSERT_VALUES, SCATTER_FORWARD)); 10518be712e4SBarry Smith PetscCall(VecGetArrayRead(ghostMaxPE, &buf)); 1052*c376f201SBarry Smith for (PetscInt kk = 0; kk < num_ghosts; kk++) { 10538be712e4SBarry Smith lghost_max_pe[kk] = (PetscMPIInt)PetscRealPart(buf[kk]); // the MAX proc of the ghost now 10548be712e4SBarry Smith } 10558be712e4SBarry Smith PetscCall(VecRestoreArrayRead(ghostMaxPE, &buf)); 10568be712e4SBarry Smith } 10578be712e4SBarry Smith // if no active edges, stop 10588be712e4SBarry Smith if (gn_act_n[0] < 1) break; 10598be712e4SBarry Smith // inc and check (self stopping iteration 1060*c376f201SBarry Smith PetscCheck(old_num_edge != gn_act_n[0], PETSC_COMM_SELF, PETSC_ERR_SUP, "HEM stalled step %" PetscInt_FMT "/%" PetscInt_FMT, sub_it + 1, n_sub_its); 10618be712e4SBarry Smith sub_it++; 1062*c376f201SBarry Smith PetscCheck(sub_it < n_sub_its, PETSC_COMM_SELF, PETSC_ERR_SUP, "failed to finish HEM step %" PetscInt_FMT "/%" PetscInt_FMT, sub_it + 1, n_sub_its); 1063452acf8bSMark Adams old_num_edge = gn_act_n[0]; 10648be712e4SBarry Smith } /* sub_it loop */ 10658be712e4SBarry Smith /* clean up iteration */ 10668be712e4SBarry Smith PetscCall(PetscFree(Edges)); 1067*c376f201SBarry Smith if (isMPI) { // can be hoisted 10688be712e4SBarry Smith PetscCall(VecRestoreArrayRead(ghostMaxEdge, &lghost_max_ew)); 10698be712e4SBarry Smith PetscCall(VecDestroy(&ghostMaxEdge)); 10708be712e4SBarry Smith PetscCall(VecDestroy(&ghostMaxPE)); 1071*c376f201SBarry Smith PetscCall(PetscFree4(lghost_matched, lghost_pe, lghost_gid, lghost_max_pe)); 10728be712e4SBarry Smith } 10738be712e4SBarry Smith PetscCall(VecDestroy(&locMaxEdge)); 10748be712e4SBarry Smith PetscCall(VecDestroy(&locMaxPE)); 10758be712e4SBarry Smith /* create next graph */ 10768be712e4SBarry Smith { 10778be712e4SBarry Smith Vec diag; 1078*c376f201SBarry Smith 10798be712e4SBarry Smith /* add identity for unmatched vertices so they stay alive */ 1080*c376f201SBarry Smith for (PetscInt kk = 0, gid1, gid = Istart; kk < nloc; kk++, gid++) { 10818be712e4SBarry Smith if (!lid_matched[kk]) { 10828be712e4SBarry Smith const PetscInt lid = kk; 10838be712e4SBarry Smith PetscCDIntNd *pos; 1084*c376f201SBarry Smith 10858be712e4SBarry Smith PetscCall(PetscCDGetHeadPos(agg_llists, lid, &pos)); 1086*c376f201SBarry Smith PetscCheck(pos, PETSC_COMM_SELF, PETSC_ERR_PLIB, "empty list in singleton: %" PetscInt_FMT, gid); 10878be712e4SBarry Smith PetscCall(PetscCDIntNdGetID(pos, &gid1)); 1088*c376f201SBarry Smith PetscCheck(gid1 == gid, PETSC_COMM_SELF, PETSC_ERR_PLIB, "first in list (%" PetscInt_FMT ") in singleton not %" PetscInt_FMT, gid1, gid); 10898be712e4SBarry Smith PetscCall(MatSetValues(P, 1, &gid, 1, &gid, &one, INSERT_VALUES)); 10908be712e4SBarry Smith } 10918be712e4SBarry Smith } 10928be712e4SBarry Smith PetscCall(MatAssemblyBegin(P, MAT_FINAL_ASSEMBLY)); 10938be712e4SBarry Smith PetscCall(MatAssemblyEnd(P, MAT_FINAL_ASSEMBLY)); 10948be712e4SBarry Smith 10958be712e4SBarry Smith /* project to make new graph with collapsed edges */ 10968be712e4SBarry Smith PetscCall(MatPtAP(cMat, P, MAT_INITIAL_MATRIX, 1.0, &tMat)); 10978be712e4SBarry Smith PetscCall(MatDestroy(&P)); 10988be712e4SBarry Smith PetscCall(MatDestroy(&cMat)); 10998be712e4SBarry Smith cMat = tMat; 11008be712e4SBarry Smith PetscCall(MatCreateVecs(cMat, &diag, NULL)); 110153673ba5SSatish Balay PetscCall(MatGetDiagonal(cMat, diag)); 11028be712e4SBarry Smith PetscCall(VecReciprocal(diag)); 11038be712e4SBarry Smith PetscCall(VecSqrtAbs(diag)); 11048be712e4SBarry Smith PetscCall(MatDiagonalScale(cMat, diag, diag)); 11058be712e4SBarry Smith PetscCall(VecDestroy(&diag)); 11068be712e4SBarry Smith } 11078be712e4SBarry Smith } /* coarsen iterator */ 11088be712e4SBarry Smith 11098be712e4SBarry Smith /* make fake matrix with Mat->B only for smoothed agg QR. Need this if we make an aux graph (ie, PtAP) with k > 1 */ 11108be712e4SBarry Smith if (size > 1) { 11118be712e4SBarry Smith Mat mat; 11128be712e4SBarry Smith PetscCDIntNd *pos; 11138be712e4SBarry Smith PetscInt NN, MM, jj = 0, mxsz = 0; 11148be712e4SBarry Smith 1115*c376f201SBarry Smith for (PetscInt kk = 0; kk < nloc; kk++) { 11168be712e4SBarry Smith PetscCall(PetscCDCountAt(agg_llists, kk, &jj)); 11178be712e4SBarry Smith if (jj > mxsz) mxsz = jj; 11188be712e4SBarry Smith } 11198be712e4SBarry Smith PetscCall(MatGetSize(a_Gmat, &MM, &NN)); 11208be712e4SBarry Smith if (mxsz > MM - nloc) mxsz = MM - nloc; 11218be712e4SBarry Smith /* matrix of ghost adj for square graph */ 11228be712e4SBarry Smith PetscCall(MatCreateAIJ(comm, nloc, nloc, PETSC_DETERMINE, PETSC_DETERMINE, 0, NULL, mxsz, NULL, &mat)); 1123*c376f201SBarry Smith for (PetscInt lid = 0, gid = Istart; lid < nloc; lid++, gid++) { 11248be712e4SBarry Smith PetscCall(PetscCDGetHeadPos(agg_llists, lid, &pos)); 11258be712e4SBarry Smith while (pos) { 11268be712e4SBarry Smith PetscInt gid1; 1127*c376f201SBarry Smith 11288be712e4SBarry Smith PetscCall(PetscCDIntNdGetID(pos, &gid1)); 11298be712e4SBarry Smith PetscCall(PetscCDGetNextPos(agg_llists, lid, &pos)); 1130*c376f201SBarry Smith if (gid1 < Istart || gid1 >= Istart + nloc) PetscCall(MatSetValues(mat, 1, &gid, 1, &gid1, &one, ADD_VALUES)); 11318be712e4SBarry Smith } 11328be712e4SBarry Smith } 11338be712e4SBarry Smith PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY)); 11348be712e4SBarry Smith PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY)); 11358be712e4SBarry Smith PetscCall(PetscCDSetMat(agg_llists, mat)); 11368be712e4SBarry Smith PetscCall(PetscCDDestroy(ghost_deleted_list)); 11378be712e4SBarry Smith if (rbuff_sz) PetscCall(PetscFree(rbuff)); // always true 11388be712e4SBarry Smith } 11398be712e4SBarry Smith // move BCs into some node 11408be712e4SBarry Smith if (bc_list) { 11418be712e4SBarry Smith PetscCDIntNd *pos; 1142*c376f201SBarry Smith 11438be712e4SBarry Smith PetscCall(PetscCDGetHeadPos(bc_list, 0, &pos)); 11448be712e4SBarry Smith while (pos) { 11458be712e4SBarry Smith PetscInt gid1; 1146*c376f201SBarry Smith 11478be712e4SBarry Smith PetscCall(PetscCDIntNdGetID(pos, &gid1)); 11488be712e4SBarry Smith PetscCall(PetscCDGetNextPos(bc_list, 0, &pos)); 11498be712e4SBarry Smith PetscCall(PetscCDAppendID(agg_llists, bc_agg, gid1)); 11508be712e4SBarry Smith } 11518be712e4SBarry Smith PetscCall(PetscCDRemoveAllAt(bc_list, 0)); 11528be712e4SBarry Smith PetscCall(PetscCDDestroy(bc_list)); 11538be712e4SBarry Smith } 11548be712e4SBarry Smith { 11558be712e4SBarry Smith // check sizes -- all vertices must get in graph 11568be712e4SBarry Smith PetscInt sz, globalsz, MM; 1157*c376f201SBarry Smith 11588be712e4SBarry Smith PetscCall(MatGetSize(a_Gmat, &MM, NULL)); 11598be712e4SBarry Smith PetscCall(PetscCDCount(agg_llists, &sz)); 11608be712e4SBarry Smith PetscCall(MPIU_Allreduce(&sz, &globalsz, 1, MPIU_INT, MPI_SUM, comm)); 1161*c376f201SBarry Smith PetscCheck(MM == globalsz, comm, PETSC_ERR_SUP, "lost %" PetscInt_FMT " equations ?", (MM - globalsz)); 11628be712e4SBarry Smith } 11638be712e4SBarry Smith // cleanup 11648be712e4SBarry Smith PetscCall(MatDestroy(&cMat)); 1165*c376f201SBarry Smith PetscCall(PetscFree3(lid_matched, lid_cprowID, lid_max_pe)); 11668be712e4SBarry Smith PetscCall(ISDestroy(&info_is)); 11678be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 11688be712e4SBarry Smith } 11698be712e4SBarry Smith 11708be712e4SBarry Smith /* 11718be712e4SBarry Smith HEM coarsen, simple greedy. 11728be712e4SBarry Smith */ 11738be712e4SBarry Smith static PetscErrorCode MatCoarsenApply_HEM(MatCoarsen coarse) 11748be712e4SBarry Smith { 11758be712e4SBarry Smith Mat mat = coarse->graph; 11768be712e4SBarry Smith 11778be712e4SBarry Smith PetscFunctionBegin; 11788be712e4SBarry Smith PetscCall(MatCoarsenApply_HEM_private(mat, coarse->max_it, coarse->threshold, &coarse->agg_lists)); 11798be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 11808be712e4SBarry Smith } 11818be712e4SBarry Smith 11828be712e4SBarry Smith static PetscErrorCode MatCoarsenView_HEM(MatCoarsen coarse, PetscViewer viewer) 11838be712e4SBarry Smith { 11848be712e4SBarry Smith PetscMPIInt rank; 11858be712e4SBarry Smith PetscBool iascii; 11868be712e4SBarry Smith 11878be712e4SBarry Smith PetscFunctionBegin; 11888be712e4SBarry Smith PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)coarse), &rank)); 11898be712e4SBarry Smith PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 11908be712e4SBarry Smith if (iascii) { 11918be712e4SBarry Smith PetscCDIntNd *pos, *pos2; 1192e0b7e82fSBarry Smith PetscViewerFormat format; 1193e0b7e82fSBarry Smith 1194*c376f201SBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " matching steps with threshold = %g\n", coarse->max_it, (double)coarse->threshold)); 1195e0b7e82fSBarry Smith PetscCall(PetscViewerGetFormat(viewer, &format)); 1196e0b7e82fSBarry Smith if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 1197bfa5bdfcSBarry Smith if (coarse->agg_lists) { 11988be712e4SBarry Smith PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 11998be712e4SBarry Smith for (PetscInt kk = 0; kk < coarse->agg_lists->size; kk++) { 12008be712e4SBarry Smith PetscCall(PetscCDGetHeadPos(coarse->agg_lists, kk, &pos)); 1201*c376f201SBarry Smith if ((pos2 = pos)) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "selected local %" PetscInt_FMT ": ", kk)); 12028be712e4SBarry Smith while (pos) { 12038be712e4SBarry Smith PetscInt gid1; 1204*c376f201SBarry Smith 12058be712e4SBarry Smith PetscCall(PetscCDIntNdGetID(pos, &gid1)); 12068be712e4SBarry Smith PetscCall(PetscCDGetNextPos(coarse->agg_lists, kk, &pos)); 1207*c376f201SBarry Smith PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " %" PetscInt_FMT " ", gid1)); 12088be712e4SBarry Smith } 12098be712e4SBarry Smith if (pos2) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n")); 12108be712e4SBarry Smith } 12118be712e4SBarry Smith PetscCall(PetscViewerFlush(viewer)); 12128be712e4SBarry Smith PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 1213bfa5bdfcSBarry Smith } else { 1214bfa5bdfcSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, " HEM aggregator lists are not available\n")); 1215bfa5bdfcSBarry Smith } 12168be712e4SBarry Smith } 1217e0b7e82fSBarry Smith } 12188be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 12198be712e4SBarry Smith } 12208be712e4SBarry Smith 1221*c376f201SBarry Smith /*MC 1222*c376f201SBarry Smith MATCOARSENHEM - A coarsener that uses HEM a simple greedy coarsener 1223*c376f201SBarry Smith 1224*c376f201SBarry Smith Level: beginner 1225*c376f201SBarry Smith 1226*c376f201SBarry Smith .seealso: `MatCoarsen`, `MatCoarsenMISKSetDistance()`, `MatCoarsenApply()`, `MatCoarsenSetType()`, `MatCoarsenType`, `MatCoarsenCreate()`, `MATCOARSENMISK`, `MATCOARSENMIS` 1227*c376f201SBarry Smith M*/ 1228*c376f201SBarry Smith 12298be712e4SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenCreate_HEM(MatCoarsen coarse) 12308be712e4SBarry Smith { 12318be712e4SBarry Smith PetscFunctionBegin; 12328be712e4SBarry Smith coarse->ops->apply = MatCoarsenApply_HEM; 12338be712e4SBarry Smith coarse->ops->view = MatCoarsenView_HEM; 12348be712e4SBarry Smith coarse->max_it = 4; 12358be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 12368be712e4SBarry Smith } 1237