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 */
PetscCDCreate(PetscInt a_size,PetscCoarsenData ** a_out)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 */
PetscCDDestroy(PetscCoarsenData * ail)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;
39c376f201SBarry 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 */
PetscCDSetChunkSize(PetscCoarsenData * ail,PetscInt a_sz)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 */
PetscCDGetNewNode(PetscCoarsenData * ail,PetscCDIntNd ** a_out,PetscInt a_id)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;
68c376f201SBarry 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;
82c376f201SBarry 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 */
PetscCDIntNdSetID(PetscCDIntNd * a_this,PetscInt a_id)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 */
PetscCDIntNdGetID(const PetscCDIntNd * a_this,PetscInt * a_gid)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 */
PetscCDGetHeadPos(const PetscCoarsenData * ail,PetscInt a_idx,PetscCDIntNd ** pos)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 */
PetscCDGetNextPos(const PetscCoarsenData * ail,PetscInt l_idx,PetscCDIntNd ** pos)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 */
PetscCDAppendID(PetscCoarsenData * ail,PetscInt a_idx,PetscInt a_id)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 */
PetscCDAppendNode(PetscCoarsenData * ail,PetscInt a_idx,PetscCDIntNd * a_n)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 */
PetscCDRemoveNextNode(PetscCoarsenData * ail,PetscInt a_idx,PetscCDIntNd * a_last)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 */
PetscCDPrint(const PetscCoarsenData * ail,PetscInt Istart,MPI_Comm comm)201c376f201SBarry 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];
209c376f201SBarry 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 */
PetscCDMoveAppend(PetscCoarsenData * ail,PetscInt a_destidx,PetscInt a_srcidx)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 */
PetscCDRemoveAllAt(PetscCoarsenData * ail,PetscInt a_idx)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 */
PetscCDCountAt(const PetscCoarsenData * ail,PetscInt a_idx,PetscInt * a_sz)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 */
PetscCDCount(const PetscCoarsenData * ail,PetscInt * a_sz)2838be712e4SBarry Smith PetscErrorCode PetscCDCount(const PetscCoarsenData *ail, PetscInt *a_sz)
2848be712e4SBarry Smith {
2858be712e4SBarry Smith PetscInt sz = 0;
2868be712e4SBarry Smith
2878be712e4SBarry Smith PetscFunctionBegin;
288c376f201SBarry Smith for (PetscInt ii = 0; ii < ail->size; ii++) {
2898be712e4SBarry Smith PetscCDIntNd *n1 = ail->array[ii];
290c376f201SBarry 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 */
PetscCDIsEmptyAt(const PetscCoarsenData * ail,PetscInt a_idx,PetscBool * a_e)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 */
PetscCDGetNonemptyIS(PetscCoarsenData * ail,IS * a_mis)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 */
PetscCDGetMat(PetscCoarsenData * ail,Mat * a_mat)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 */
PetscCDSetMat(PetscCoarsenData * ail,Mat a_mat)3438be712e4SBarry Smith PetscErrorCode PetscCDSetMat(PetscCoarsenData *ail, Mat a_mat)
3448be712e4SBarry Smith {
3458be712e4SBarry Smith PetscFunctionBegin;
3463a7d0413SPierre Jolivet if (ail->mat) PetscCall(MatDestroy(&ail->mat)); //should not happen
3478be712e4SBarry Smith ail->mat = a_mat;
3488be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
3498be712e4SBarry Smith }
3508be712e4SBarry Smith
3518926f930SMark Adams /* PetscCDClearMat
3528926f930SMark Adams */
PetscCDClearMat(PetscCoarsenData * ail)3538926f930SMark Adams PetscErrorCode PetscCDClearMat(PetscCoarsenData *ail)
3548926f930SMark Adams {
3558926f930SMark Adams PetscFunctionBegin;
3568926f930SMark Adams ail->mat = NULL;
3578926f930SMark Adams PetscFunctionReturn(PETSC_SUCCESS);
3588926f930SMark Adams }
3598926f930SMark Adams
3608be712e4SBarry Smith /* PetscCDGetASMBlocks - get IS of aggregates for ASM smoothers
3618be712e4SBarry Smith */
PetscCDGetASMBlocks(const PetscCoarsenData * ail,const PetscInt a_bs,PetscInt * a_sz,IS ** a_local_is)3628be712e4SBarry Smith PetscErrorCode PetscCDGetASMBlocks(const PetscCoarsenData *ail, const PetscInt a_bs, PetscInt *a_sz, IS **a_local_is)
3638be712e4SBarry Smith {
3648be712e4SBarry Smith PetscCDIntNd *n;
3658be712e4SBarry Smith PetscInt lsz, ii, kk, *idxs, jj, gid;
3668be712e4SBarry Smith IS *is_loc = NULL;
3678be712e4SBarry Smith
3688be712e4SBarry Smith PetscFunctionBegin;
3698be712e4SBarry Smith for (ii = kk = 0; ii < ail->size; ii++) {
3708be712e4SBarry Smith if (ail->array[ii]) kk++;
3718be712e4SBarry Smith }
3728be712e4SBarry Smith *a_sz = kk;
3738be712e4SBarry Smith PetscCall(PetscMalloc1(kk, &is_loc));
3748be712e4SBarry Smith for (ii = kk = 0; ii < ail->size; ii++) {
3758be712e4SBarry Smith for (lsz = 0, n = ail->array[ii]; n; lsz++, n = n->next) /* void */
3768be712e4SBarry Smith ;
3778be712e4SBarry Smith if (lsz) {
3788be712e4SBarry Smith PetscCall(PetscMalloc1(a_bs * lsz, &idxs));
3798be712e4SBarry Smith for (lsz = 0, n = ail->array[ii]; n; n = n->next) {
3808be712e4SBarry Smith PetscCall(PetscCDIntNdGetID(n, &gid));
3818be712e4SBarry Smith for (jj = 0; jj < a_bs; lsz++, jj++) idxs[lsz] = a_bs * gid + jj;
3828be712e4SBarry Smith }
3838be712e4SBarry Smith PetscCall(ISCreateGeneral(PETSC_COMM_SELF, lsz, idxs, PETSC_OWN_POINTER, &is_loc[kk++]));
3848be712e4SBarry Smith }
3858be712e4SBarry Smith }
3868be712e4SBarry Smith PetscCheck(*a_sz == kk, PETSC_COMM_SELF, PETSC_ERR_PLIB, "*a_sz %" PetscInt_FMT " != kk %" PetscInt_FMT, *a_sz, kk);
3878be712e4SBarry Smith *a_local_is = is_loc; /* out */
3888be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
3898be712e4SBarry Smith }
3908be712e4SBarry Smith
3918be712e4SBarry Smith /* edge for priority queue */
3928be712e4SBarry Smith typedef struct edge_tag {
3938be712e4SBarry Smith PetscReal weight;
3948be712e4SBarry Smith PetscInt lid0, gid1, ghost1_idx;
3958be712e4SBarry Smith } Edge;
3968be712e4SBarry Smith
3978be712e4SBarry Smith #define MY_MEPS (PETSC_MACHINE_EPSILON * 100)
gamg_hem_compare(const void * a,const void * b)3988be712e4SBarry Smith static int gamg_hem_compare(const void *a, const void *b)
3998be712e4SBarry Smith {
4008be712e4SBarry Smith PetscReal va = ((Edge *)a)->weight, vb = ((Edge *)b)->weight;
4018be712e4SBarry Smith return (va <= vb - MY_MEPS) ? 1 : (va > vb + MY_MEPS) ? -1 : 0; /* 0 for equal */
4028be712e4SBarry Smith }
4038be712e4SBarry Smith
4048be712e4SBarry Smith /*
4058be712e4SBarry Smith MatCoarsenApply_HEM_private - parallel heavy edge matching
4068be712e4SBarry Smith
4078be712e4SBarry Smith Input Parameter:
4088be712e4SBarry Smith . a_Gmat - global matrix of the graph
4098be712e4SBarry Smith . n_iter - number of matching iterations
4108be712e4SBarry Smith . threshold - threshold for filtering graphs
4118be712e4SBarry Smith
4128be712e4SBarry Smith Output Parameter:
4138be712e4SBarry Smith . a_locals_llist - array of list of local nodes rooted at local node
4148be712e4SBarry Smith */
MatCoarsenApply_HEM_private(Mat a_Gmat,const PetscInt n_iter,const PetscReal threshold,PetscCoarsenData ** a_locals_llist)4158be712e4SBarry Smith static PetscErrorCode MatCoarsenApply_HEM_private(Mat a_Gmat, const PetscInt n_iter, const PetscReal threshold, PetscCoarsenData **a_locals_llist)
4168be712e4SBarry Smith {
4178be712e4SBarry Smith #define REQ_BF_SIZE 100
4188be712e4SBarry Smith PetscBool isMPI;
4198be712e4SBarry Smith MPI_Comm comm;
420c376f201SBarry Smith PetscInt ix, *ii, *aj, Istart, bc_agg = -1, *rbuff = NULL, rbuff_sz = 0;
421c376f201SBarry Smith PetscMPIInt rank, size, comm_procs[REQ_BF_SIZE], ncomm_procs, *lid_max_pe;
422c376f201SBarry Smith const PetscInt nloc = a_Gmat->rmap->n, request_size = PetscCeilInt((int)sizeof(MPI_Request), (int)sizeof(PetscInt));
4238be712e4SBarry Smith PetscInt *lid_cprowID;
4248be712e4SBarry Smith PetscBool *lid_matched;
4258be712e4SBarry Smith Mat_SeqAIJ *matA, *matB = NULL;
4268be712e4SBarry Smith Mat_MPIAIJ *mpimat = NULL;
4278be712e4SBarry Smith PetscScalar one = 1.;
4288be712e4SBarry Smith PetscCoarsenData *agg_llists = NULL, *ghost_deleted_list = NULL, *bc_list = NULL;
4298be712e4SBarry Smith Mat cMat, tMat, P;
4308be712e4SBarry Smith MatScalar *ap;
4318be712e4SBarry Smith IS info_is;
4328be712e4SBarry Smith
4338be712e4SBarry Smith PetscFunctionBegin;
4348be712e4SBarry Smith PetscCall(PetscObjectGetComm((PetscObject)a_Gmat, &comm));
4358be712e4SBarry Smith PetscCallMPI(MPI_Comm_rank(comm, &rank));
4368be712e4SBarry Smith PetscCallMPI(MPI_Comm_size(comm, &size));
437c376f201SBarry Smith PetscCall(MatGetOwnershipRange(a_Gmat, &Istart, NULL));
4388be712e4SBarry Smith PetscCall(ISCreate(comm, &info_is));
4398926f930SMark Adams PetscCall(PetscInfo(info_is, "Start %" PetscInt_FMT " iterations of HEM.\n", n_iter));
4408be712e4SBarry Smith
441c376f201SBarry Smith PetscCall(PetscMalloc3(nloc, &lid_matched, nloc, &lid_cprowID, nloc, &lid_max_pe));
4428be712e4SBarry Smith PetscCall(PetscCDCreate(nloc, &agg_llists));
4438be712e4SBarry Smith PetscCall(PetscCDSetChunkSize(agg_llists, nloc + 1));
4448be712e4SBarry Smith *a_locals_llist = agg_llists;
4458be712e4SBarry Smith /* add self to all lists */
446c376f201SBarry Smith for (PetscInt kk = 0; kk < nloc; kk++) PetscCall(PetscCDAppendID(agg_llists, kk, Istart + kk));
4478be712e4SBarry Smith /* make a copy of the graph, this gets destroyed in iterates */
4488be712e4SBarry Smith PetscCall(MatDuplicate(a_Gmat, MAT_COPY_VALUES, &cMat));
4498be712e4SBarry Smith PetscCall(MatConvert(cMat, MATAIJ, MAT_INPLACE_MATRIX, &cMat));
4508be712e4SBarry Smith isMPI = (PetscBool)(size > 1);
4518be712e4SBarry Smith if (isMPI) {
4528be712e4SBarry Smith /* list of deleted ghosts, should compress this */
4538be712e4SBarry Smith PetscCall(PetscCDCreate(size, &ghost_deleted_list));
4548be712e4SBarry Smith PetscCall(PetscCDSetChunkSize(ghost_deleted_list, 100));
4558be712e4SBarry Smith }
456c376f201SBarry Smith for (PetscInt iter = 0; iter < n_iter; iter++) {
4578be712e4SBarry Smith const PetscScalar *lghost_max_ew, *lid_max_ew;
4588be712e4SBarry Smith PetscBool *lghost_matched;
4598be712e4SBarry Smith PetscMPIInt *lghost_pe, *lghost_max_pe;
4608be712e4SBarry Smith Vec locMaxEdge, ghostMaxEdge, ghostMaxPE, locMaxPE;
4618be712e4SBarry Smith PetscInt *lghost_gid, nEdges, nEdges0, num_ghosts = 0;
4628be712e4SBarry Smith Edge *Edges;
463c376f201SBarry Smith const PetscInt n_sub_its = 1000; // in case of a bug, stop at some point
464c376f201SBarry Smith
4658be712e4SBarry Smith /* get submatrices of cMat */
466c376f201SBarry Smith for (PetscInt kk = 0; kk < nloc; kk++) lid_cprowID[kk] = -1;
4678be712e4SBarry Smith if (isMPI) {
4688be712e4SBarry Smith mpimat = (Mat_MPIAIJ *)cMat->data;
4698be712e4SBarry Smith matA = (Mat_SeqAIJ *)mpimat->A->data;
4708be712e4SBarry Smith matB = (Mat_SeqAIJ *)mpimat->B->data;
4718be712e4SBarry Smith if (!matB->compressedrow.use) {
4728be712e4SBarry Smith /* force construction of compressed row data structure since code below requires it */
4738be712e4SBarry Smith PetscCall(MatCheckCompressedRow(mpimat->B, matB->nonzerorowcnt, &matB->compressedrow, matB->i, mpimat->B->rmap->n, -1.0));
4748be712e4SBarry Smith }
4758be712e4SBarry Smith /* set index into compressed row 'lid_cprowID' */
4768be712e4SBarry Smith for (ix = 0; ix < matB->compressedrow.nrows; ix++) {
4778be712e4SBarry Smith PetscInt *ridx = matB->compressedrow.rindex, lid = ridx[ix];
4788be712e4SBarry Smith if (ridx[ix] >= 0) lid_cprowID[lid] = ix;
4798be712e4SBarry Smith }
4808be712e4SBarry Smith } else {
4818be712e4SBarry Smith matA = (Mat_SeqAIJ *)cMat->data;
4828be712e4SBarry Smith }
4838be712e4SBarry Smith /* set matched flags: true for empty list */
484c376f201SBarry Smith for (PetscInt kk = 0; kk < nloc; kk++) {
4858be712e4SBarry Smith PetscCall(PetscCDCountAt(agg_llists, kk, &ix));
4868be712e4SBarry Smith if (ix > 0) lid_matched[kk] = PETSC_FALSE;
4878be712e4SBarry Smith else lid_matched[kk] = PETSC_TRUE; // call deleted gids as matched
4888be712e4SBarry Smith }
4898be712e4SBarry Smith /* max edge and pe vecs */
4908be712e4SBarry Smith PetscCall(MatCreateVecs(cMat, &locMaxEdge, NULL));
4918be712e4SBarry Smith PetscCall(MatCreateVecs(cMat, &locMaxPE, NULL));
4928be712e4SBarry Smith /* get 'lghost_pe' & 'lghost_gid' & init. 'lghost_matched' using 'mpimat->lvec' */
4938be712e4SBarry Smith if (isMPI) {
4948be712e4SBarry Smith Vec vec;
4958be712e4SBarry Smith PetscScalar vval;
4968be712e4SBarry Smith const PetscScalar *buf;
497c376f201SBarry Smith
4988be712e4SBarry Smith PetscCall(MatCreateVecs(cMat, &vec, NULL));
4998be712e4SBarry Smith PetscCall(VecGetLocalSize(mpimat->lvec, &num_ghosts));
5008be712e4SBarry Smith /* lghost_matched */
501c376f201SBarry Smith for (PetscInt kk = 0, gid = Istart; kk < nloc; kk++, gid++) {
5028be712e4SBarry Smith PetscScalar vval = lid_matched[kk] ? 1.0 : 0.0;
503c376f201SBarry Smith
5048be712e4SBarry Smith PetscCall(VecSetValues(vec, 1, &gid, &vval, INSERT_VALUES));
5058be712e4SBarry Smith }
5068be712e4SBarry Smith PetscCall(VecAssemblyBegin(vec));
5078be712e4SBarry Smith PetscCall(VecAssemblyEnd(vec));
5088be712e4SBarry Smith PetscCall(VecScatterBegin(mpimat->Mvctx, vec, mpimat->lvec, INSERT_VALUES, SCATTER_FORWARD));
5098be712e4SBarry Smith PetscCall(VecScatterEnd(mpimat->Mvctx, vec, mpimat->lvec, INSERT_VALUES, SCATTER_FORWARD));
5108be712e4SBarry Smith PetscCall(VecGetArrayRead(mpimat->lvec, &buf)); /* get proc ID in 'buf' */
511c376f201SBarry Smith PetscCall(PetscMalloc4(num_ghosts, &lghost_matched, num_ghosts, &lghost_pe, num_ghosts, &lghost_gid, num_ghosts, &lghost_max_pe));
512c376f201SBarry Smith
5133a7d0413SPierre Jolivet for (PetscInt kk = 0; kk < num_ghosts; kk++) lghost_matched[kk] = (PetscBool)(PetscRealPart(buf[kk]) != 0); // the proc of the ghost for now
5148be712e4SBarry Smith PetscCall(VecRestoreArrayRead(mpimat->lvec, &buf));
5158be712e4SBarry Smith /* lghost_pe */
51657508eceSPierre Jolivet vval = (PetscScalar)rank;
517c376f201SBarry Smith for (PetscInt kk = 0, gid = Istart; kk < nloc; kk++, gid++) PetscCall(VecSetValues(vec, 1, &gid, &vval, INSERT_VALUES)); /* set with GID */
5188be712e4SBarry Smith PetscCall(VecAssemblyBegin(vec));
5198be712e4SBarry Smith PetscCall(VecAssemblyEnd(vec));
5208be712e4SBarry Smith PetscCall(VecScatterBegin(mpimat->Mvctx, vec, mpimat->lvec, INSERT_VALUES, SCATTER_FORWARD));
5218be712e4SBarry Smith PetscCall(VecScatterEnd(mpimat->Mvctx, vec, mpimat->lvec, INSERT_VALUES, SCATTER_FORWARD));
5228be712e4SBarry Smith PetscCall(VecGetArrayRead(mpimat->lvec, &buf)); /* get proc ID in 'buf' */
523c376f201SBarry Smith for (PetscInt kk = 0; kk < num_ghosts; kk++) lghost_pe[kk] = (PetscMPIInt)PetscRealPart(buf[kk]); // the proc of the ghost for now
5248be712e4SBarry Smith PetscCall(VecRestoreArrayRead(mpimat->lvec, &buf));
5258be712e4SBarry Smith /* lghost_gid */
526c376f201SBarry Smith for (PetscInt kk = 0, gid = Istart; kk < nloc; kk++, gid++) {
52757508eceSPierre Jolivet vval = (PetscScalar)gid;
528c376f201SBarry Smith
5298be712e4SBarry Smith PetscCall(VecSetValues(vec, 1, &gid, &vval, INSERT_VALUES)); /* set with GID */
5308be712e4SBarry Smith }
5318be712e4SBarry Smith PetscCall(VecAssemblyBegin(vec));
5328be712e4SBarry Smith PetscCall(VecAssemblyEnd(vec));
5338be712e4SBarry Smith PetscCall(VecScatterBegin(mpimat->Mvctx, vec, mpimat->lvec, INSERT_VALUES, SCATTER_FORWARD));
5348be712e4SBarry Smith PetscCall(VecScatterEnd(mpimat->Mvctx, vec, mpimat->lvec, INSERT_VALUES, SCATTER_FORWARD));
5358be712e4SBarry Smith PetscCall(VecDestroy(&vec));
5368be712e4SBarry Smith PetscCall(VecGetArrayRead(mpimat->lvec, &buf)); /* get proc ID in 'lghost_gid' */
537c376f201SBarry Smith for (PetscInt kk = 0; kk < num_ghosts; kk++) lghost_gid[kk] = (PetscInt)PetscRealPart(buf[kk]);
5388be712e4SBarry Smith PetscCall(VecRestoreArrayRead(mpimat->lvec, &buf));
5398be712e4SBarry Smith }
5408be712e4SBarry Smith // get 'comm_procs' (could hoist)
541c376f201SBarry Smith for (PetscInt kk = 0; kk < REQ_BF_SIZE; kk++) comm_procs[kk] = -1;
5428be712e4SBarry Smith for (ix = 0, ncomm_procs = 0; ix < num_ghosts; ix++) {
5438be712e4SBarry Smith PetscMPIInt proc = lghost_pe[ix], idx = -1;
544c376f201SBarry Smith
5456497c311SBarry Smith for (PetscMPIInt k = 0; k < ncomm_procs && idx == -1; k++)
5468be712e4SBarry Smith if (comm_procs[k] == proc) idx = k;
5476497c311SBarry Smith if (idx == -1) comm_procs[ncomm_procs++] = proc;
548c376f201SBarry Smith PetscCheck(ncomm_procs != REQ_BF_SIZE, PETSC_COMM_SELF, PETSC_ERR_SUP, "Receive request array too small: %d", ncomm_procs);
5498be712e4SBarry Smith }
5508be712e4SBarry Smith /* count edges, compute initial 'locMaxEdge', 'locMaxPE' */
5518be712e4SBarry Smith nEdges0 = 0;
552c376f201SBarry Smith for (PetscInt kk = 0, gid = Istart; kk < nloc; kk++, gid++) {
5538be712e4SBarry Smith PetscReal max_e = 0., tt;
5548be712e4SBarry Smith PetscScalar vval;
5558be712e4SBarry Smith PetscInt lid = kk, max_pe = rank, pe, n;
556c376f201SBarry Smith
5578be712e4SBarry Smith ii = matA->i;
5588be712e4SBarry Smith n = ii[lid + 1] - ii[lid];
5598e3a54c0SPierre Jolivet aj = PetscSafePointerPlusOffset(matA->j, ii[lid]);
5608e3a54c0SPierre Jolivet ap = PetscSafePointerPlusOffset(matA->a, ii[lid]);
561c376f201SBarry Smith for (PetscInt jj = 0; jj < n; jj++) {
5628be712e4SBarry Smith PetscInt lidj = aj[jj];
563c376f201SBarry Smith
5648be712e4SBarry Smith if ((tt = PetscRealPart(ap[jj])) > threshold && lidj != lid) {
5658be712e4SBarry Smith if (tt > max_e) max_e = tt;
5668be712e4SBarry Smith if (lidj > lid) nEdges0++;
5678be712e4SBarry Smith }
5688be712e4SBarry Smith }
5698be712e4SBarry Smith if ((ix = lid_cprowID[lid]) != -1) { /* if I have any ghost neighbors */
5708be712e4SBarry Smith ii = matB->compressedrow.i;
5718be712e4SBarry Smith n = ii[ix + 1] - ii[ix];
5728be712e4SBarry Smith ap = matB->a + ii[ix];
5738be712e4SBarry Smith aj = matB->j + ii[ix];
574c376f201SBarry Smith for (PetscInt jj = 0; jj < n; jj++) {
5758be712e4SBarry Smith if ((tt = PetscRealPart(ap[jj])) > threshold) {
5768be712e4SBarry Smith if (tt > max_e) max_e = tt;
5778be712e4SBarry Smith nEdges0++;
5788be712e4SBarry Smith if ((pe = lghost_pe[aj[jj]]) > max_pe) max_pe = pe;
5798be712e4SBarry Smith }
5808be712e4SBarry Smith }
5818be712e4SBarry Smith }
5828be712e4SBarry Smith vval = max_e;
5838be712e4SBarry Smith PetscCall(VecSetValues(locMaxEdge, 1, &gid, &vval, INSERT_VALUES));
5848be712e4SBarry Smith vval = (PetscScalar)max_pe;
5858be712e4SBarry Smith PetscCall(VecSetValues(locMaxPE, 1, &gid, &vval, INSERT_VALUES));
5868be712e4SBarry Smith if (iter == 0 && max_e <= MY_MEPS) { // add BCs to fake aggregate
5878be712e4SBarry Smith lid_matched[lid] = PETSC_TRUE;
5888be712e4SBarry Smith if (bc_agg == -1) {
5898be712e4SBarry Smith bc_agg = lid;
5908be712e4SBarry Smith PetscCall(PetscCDCreate(1, &bc_list));
5918be712e4SBarry Smith }
5928be712e4SBarry Smith PetscCall(PetscCDRemoveAllAt(agg_llists, lid));
593c376f201SBarry Smith PetscCall(PetscCDAppendID(bc_list, 0, Istart + lid));
5948be712e4SBarry Smith }
5958be712e4SBarry Smith }
5968be712e4SBarry Smith PetscCall(VecAssemblyBegin(locMaxEdge));
5978be712e4SBarry Smith PetscCall(VecAssemblyEnd(locMaxEdge));
5988be712e4SBarry Smith PetscCall(VecAssemblyBegin(locMaxPE));
5998be712e4SBarry Smith PetscCall(VecAssemblyEnd(locMaxPE));
6008be712e4SBarry Smith /* make 'ghostMaxEdge_max_ew', 'lghost_max_pe' */
601c376f201SBarry Smith if (isMPI) {
6028be712e4SBarry Smith const PetscScalar *buf;
603c376f201SBarry Smith
6048be712e4SBarry Smith PetscCall(VecDuplicate(mpimat->lvec, &ghostMaxEdge));
6058be712e4SBarry Smith PetscCall(VecScatterBegin(mpimat->Mvctx, locMaxEdge, ghostMaxEdge, INSERT_VALUES, SCATTER_FORWARD));
6068be712e4SBarry Smith PetscCall(VecScatterEnd(mpimat->Mvctx, locMaxEdge, ghostMaxEdge, INSERT_VALUES, SCATTER_FORWARD));
6078be712e4SBarry Smith
6088be712e4SBarry Smith PetscCall(VecDuplicate(mpimat->lvec, &ghostMaxPE));
6098be712e4SBarry Smith PetscCall(VecScatterBegin(mpimat->Mvctx, locMaxPE, ghostMaxPE, INSERT_VALUES, SCATTER_FORWARD));
6108be712e4SBarry Smith PetscCall(VecScatterEnd(mpimat->Mvctx, locMaxPE, ghostMaxPE, INSERT_VALUES, SCATTER_FORWARD));
6118be712e4SBarry Smith PetscCall(VecGetArrayRead(ghostMaxPE, &buf));
612c376f201SBarry Smith for (PetscInt kk = 0; kk < num_ghosts; kk++) lghost_max_pe[kk] = (PetscMPIInt)PetscRealPart(buf[kk]); // the MAX proc of the ghost now
6138be712e4SBarry Smith PetscCall(VecRestoreArrayRead(ghostMaxPE, &buf));
6148be712e4SBarry Smith }
6158be712e4SBarry Smith { // make lid_max_pe
6168be712e4SBarry Smith const PetscScalar *buf;
617c376f201SBarry Smith
6188be712e4SBarry Smith PetscCall(VecGetArrayRead(locMaxPE, &buf));
619c376f201SBarry Smith for (PetscInt kk = 0; kk < nloc; kk++) lid_max_pe[kk] = (PetscMPIInt)PetscRealPart(buf[kk]); // the MAX proc of the ghost now
6208be712e4SBarry Smith PetscCall(VecRestoreArrayRead(locMaxPE, &buf));
6218be712e4SBarry Smith }
6228be712e4SBarry Smith /* setup sorted list of edges, and make 'Edges' */
6238be712e4SBarry Smith PetscCall(PetscMalloc1(nEdges0, &Edges));
6248be712e4SBarry Smith nEdges = 0;
625c376f201SBarry Smith for (PetscInt kk = 0, n; kk < nloc; kk++) {
6268be712e4SBarry Smith const PetscInt lid = kk;
6278be712e4SBarry Smith PetscReal tt;
628c376f201SBarry Smith
6298be712e4SBarry Smith ii = matA->i;
6308be712e4SBarry Smith n = ii[lid + 1] - ii[lid];
6318e3a54c0SPierre Jolivet aj = PetscSafePointerPlusOffset(matA->j, ii[lid]);
6328e3a54c0SPierre Jolivet ap = PetscSafePointerPlusOffset(matA->a, ii[lid]);
633c376f201SBarry Smith for (PetscInt jj = 0; jj < n; jj++) {
6348be712e4SBarry Smith PetscInt lidj = aj[jj];
635c376f201SBarry Smith
6368be712e4SBarry Smith if ((tt = PetscRealPart(ap[jj])) > threshold && lidj != lid) {
6378be712e4SBarry Smith if (lidj > lid) {
6388be712e4SBarry Smith Edges[nEdges].lid0 = lid;
639c376f201SBarry Smith Edges[nEdges].gid1 = lidj + Istart;
6408be712e4SBarry Smith Edges[nEdges].ghost1_idx = -1;
6418be712e4SBarry Smith Edges[nEdges].weight = tt;
6428be712e4SBarry Smith nEdges++;
6438be712e4SBarry Smith }
6448be712e4SBarry Smith }
6458be712e4SBarry Smith }
6468be712e4SBarry Smith if ((ix = lid_cprowID[lid]) != -1) { /* if I have any ghost neighbor */
6478be712e4SBarry Smith ii = matB->compressedrow.i;
6488be712e4SBarry Smith n = ii[ix + 1] - ii[ix];
6498be712e4SBarry Smith ap = matB->a + ii[ix];
6508be712e4SBarry Smith aj = matB->j + ii[ix];
651c376f201SBarry Smith for (PetscInt jj = 0; jj < n; jj++) {
6528be712e4SBarry Smith if ((tt = PetscRealPart(ap[jj])) > threshold) {
6538be712e4SBarry Smith Edges[nEdges].lid0 = lid;
6548be712e4SBarry Smith Edges[nEdges].gid1 = lghost_gid[aj[jj]];
6558be712e4SBarry Smith Edges[nEdges].ghost1_idx = aj[jj];
6568be712e4SBarry Smith Edges[nEdges].weight = tt;
6578be712e4SBarry Smith nEdges++;
6588be712e4SBarry Smith }
6598be712e4SBarry Smith }
6608be712e4SBarry Smith }
6618be712e4SBarry Smith }
662c376f201SBarry Smith PetscCheck(nEdges == nEdges0, PETSC_COMM_SELF, PETSC_ERR_SUP, "nEdges != nEdges0: %" PetscInt_FMT " %" PetscInt_FMT, nEdges0, nEdges);
663810441c8SPierre Jolivet if (Edges) qsort(Edges, nEdges, sizeof(Edge), gamg_hem_compare);
6648be712e4SBarry Smith
665c376f201SBarry Smith PetscCall(PetscInfo(info_is, "[%d] HEM iteration %" PetscInt_FMT " with %" PetscInt_FMT " edges\n", rank, iter, nEdges));
6668be712e4SBarry Smith
6678be712e4SBarry Smith /* projection matrix */
6688be712e4SBarry Smith PetscCall(MatCreate(comm, &P));
6698be712e4SBarry Smith PetscCall(MatSetType(P, MATAIJ));
6708be712e4SBarry Smith PetscCall(MatSetSizes(P, nloc, nloc, PETSC_DETERMINE, PETSC_DETERMINE));
6718be712e4SBarry Smith PetscCall(MatMPIAIJSetPreallocation(P, 1, NULL, 1, NULL));
6728be712e4SBarry Smith PetscCall(MatSeqAIJSetPreallocation(P, 1, NULL));
6738be712e4SBarry Smith PetscCall(MatSetUp(P));
6748be712e4SBarry Smith /* process - communicate - process */
675c376f201SBarry Smith for (PetscInt sub_it = 0, old_num_edge = 0; /* sub_it < n_sub_its */; /* sub_it++ */) {
6768be712e4SBarry Smith PetscInt nactive_edges = 0, n_act_n[3], gn_act_n[3];
6778be712e4SBarry Smith PetscMPIInt tag1, tag2;
678c376f201SBarry Smith
6798be712e4SBarry Smith PetscCall(VecGetArrayRead(locMaxEdge, &lid_max_ew));
6808be712e4SBarry Smith if (isMPI) {
6818be712e4SBarry Smith PetscCall(VecGetArrayRead(ghostMaxEdge, &lghost_max_ew));
6828be712e4SBarry Smith PetscCall(PetscCommGetNewTag(comm, &tag1));
6838be712e4SBarry Smith PetscCall(PetscCommGetNewTag(comm, &tag2));
6848be712e4SBarry Smith }
685c376f201SBarry Smith for (PetscInt kk = 0; kk < nEdges; kk++) {
6868be712e4SBarry Smith const Edge *e = &Edges[kk];
687c376f201SBarry Smith const PetscInt lid0 = e->lid0, gid1 = e->gid1, ghost1_idx = e->ghost1_idx, gid0 = lid0 + Istart, lid1 = gid1 - Istart;
6888be712e4SBarry Smith PetscBool isOK = PETSC_TRUE, print = PETSC_FALSE;
689c376f201SBarry Smith
690c376f201SBarry Smith if (print)
691c376f201SBarry 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"));
6928be712e4SBarry Smith /* skip if either vertex is matched already */
6938be712e4SBarry Smith if (lid_matched[lid0] || (ghost1_idx != -1 && lghost_matched[ghost1_idx]) || (ghost1_idx == -1 && lid_matched[lid1])) continue;
6948be712e4SBarry Smith
6958be712e4SBarry Smith nactive_edges++;
6968be712e4SBarry 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]));
697835f2295SStefano Zampini 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]) - e->weight)));
6988be712e4SBarry Smith // smaller edge, lid_max_ew get updated - e0
6998be712e4SBarry Smith if (PetscRealPart(lid_max_ew[lid0]) > e->weight + MY_MEPS) {
7008be712e4SBarry Smith if (print)
701c376f201SBarry 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]),
7028be712e4SBarry Smith (double)e->weight));
7038be712e4SBarry Smith continue; // we are basically filter edges here
7048be712e4SBarry Smith }
7058be712e4SBarry Smith // e1 - local
7068be712e4SBarry Smith if (ghost1_idx == -1) {
7078be712e4SBarry Smith if (PetscRealPart(lid_max_ew[lid1]) > e->weight + MY_MEPS) {
7088be712e4SBarry Smith if (print)
709c376f201SBarry 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)));
7108be712e4SBarry Smith continue; // we are basically filter edges here
7118be712e4SBarry Smith }
7128be712e4SBarry Smith } else { // e1 - ghost
7138be712e4SBarry Smith /* see if edge might get matched on other proc */
7148be712e4SBarry Smith PetscReal g_max_e1 = PetscRealPart(lghost_max_ew[ghost1_idx]);
715c376f201SBarry Smith
7168be712e4SBarry Smith if (print)
717c376f201SBarry 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]),
718c376f201SBarry 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]));
7198be712e4SBarry Smith if (g_max_e1 > e->weight + MY_MEPS) {
720c376f201SBarry 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 )); */
7218be712e4SBarry Smith continue;
722452acf8bSMark Adams } else if (g_max_e1 >= e->weight - MY_MEPS && lghost_pe[ghost1_idx] > rank) { // is 'lghost_max_pe[ghost1_idx] > rank' needed?
7238be712e4SBarry Smith /* check for max_ea == to this edge and larger processor that will deal with this */
7248be712e4SBarry Smith if (print)
725835f2295SStefano Zampini 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, (double)(PetscRealPart(lid_max_ew[lid0]) - e->weight), lghost_pe[ghost1_idx], lghost_max_pe[ghost1_idx], (double)g_max_e1,
726835f2295SStefano Zampini (double)e->weight));
7278be712e4SBarry Smith continue;
7288be712e4SBarry Smith } else {
729c376f201SBarry 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 )); */
7308be712e4SBarry Smith }
7318be712e4SBarry Smith }
7328be712e4SBarry Smith /* check ghost for v0 */
7338be712e4SBarry Smith if (isOK) {
7348be712e4SBarry Smith PetscReal max_e, ew;
735c376f201SBarry Smith
7368be712e4SBarry Smith if ((ix = lid_cprowID[lid0]) != -1) { /* if I have any ghost neighbors */
737c376f201SBarry Smith PetscInt n;
738c376f201SBarry Smith
7398be712e4SBarry Smith ii = matB->compressedrow.i;
7408be712e4SBarry Smith n = ii[ix + 1] - ii[ix];
7418be712e4SBarry Smith ap = matB->a + ii[ix];
7428be712e4SBarry Smith aj = matB->j + ii[ix];
743c376f201SBarry Smith for (PetscInt jj = 0; jj < n && isOK; jj++) {
7448be712e4SBarry Smith PetscInt lidj = aj[jj];
745c376f201SBarry Smith
7468be712e4SBarry Smith if (lghost_matched[lidj]) continue;
7478be712e4SBarry Smith ew = PetscRealPart(ap[jj]);
7488be712e4SBarry Smith if (ew <= threshold) continue;
7498be712e4SBarry Smith max_e = PetscRealPart(lghost_max_ew[lidj]);
750c376f201SBarry Smith
7518be712e4SBarry Smith /* check for max_e == to this edge and larger processor that will deal with this */
752452acf8bSMark Adams if (ew >= PetscRealPart(lid_max_ew[lid0]) - MY_MEPS && lghost_max_pe[lidj] > rank) isOK = PETSC_FALSE;
753c376f201SBarry 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]);
7548be712e4SBarry Smith if (print)
75557508eceSPierre Jolivet 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));
7568be712e4SBarry Smith }
757c376f201SBarry Smith if (!isOK && print) PetscCall(PetscSynchronizedPrintf(comm, "\t\t[%d] skip edge (%" PetscInt_FMT " %" PetscInt_FMT ") from ghost inspection\n", rank, gid0, gid1));
7588be712e4SBarry Smith }
7598be712e4SBarry Smith /* check local v1 */
7608be712e4SBarry Smith if (ghost1_idx == -1) {
7618be712e4SBarry Smith if ((ix = lid_cprowID[lid1]) != -1) { /* if I have any ghost neighbors */
762c376f201SBarry Smith PetscInt n;
763c376f201SBarry Smith
7648be712e4SBarry Smith ii = matB->compressedrow.i;
7658be712e4SBarry Smith n = ii[ix + 1] - ii[ix];
7668be712e4SBarry Smith ap = matB->a + ii[ix];
7678be712e4SBarry Smith aj = matB->j + ii[ix];
768c376f201SBarry Smith for (PetscInt jj = 0; jj < n && isOK; jj++) {
7698be712e4SBarry Smith PetscInt lidj = aj[jj];
770c376f201SBarry Smith
7718be712e4SBarry Smith if (lghost_matched[lidj]) continue;
7728be712e4SBarry Smith ew = PetscRealPart(ap[jj]);
7738be712e4SBarry Smith if (ew <= threshold) continue;
7748be712e4SBarry Smith max_e = PetscRealPart(lghost_max_ew[lidj]);
7758be712e4SBarry Smith /* check for max_e == to this edge and larger processor that will deal with this */
776452acf8bSMark Adams if (ew >= PetscRealPart(lid_max_ew[lid1]) - MY_MEPS && lghost_max_pe[lidj] > rank) isOK = PETSC_FALSE;
7778be712e4SBarry 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));
7788be712e4SBarry Smith if (print)
779c376f201SBarry 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]));
7808be712e4SBarry Smith }
7818be712e4SBarry Smith }
782c376f201SBarry Smith if (!isOK && print) PetscCall(PetscSynchronizedPrintf(comm, "\t\t[%d] skip edge (%" PetscInt_FMT " %" PetscInt_FMT ") from ghost inspection\n", rank, gid0, gid1));
7838be712e4SBarry Smith }
7848be712e4SBarry Smith }
7858be712e4SBarry Smith PetscReal e1_max_w = (ghost1_idx == -1 ? PetscRealPart(lid_max_ew[lid0]) : PetscRealPart(lghost_max_ew[ghost1_idx]));
7868be712e4SBarry Smith if (print)
787c376f201SBarry 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));
7888be712e4SBarry Smith /* do it */
7898be712e4SBarry Smith if (isOK) {
7908be712e4SBarry Smith if (ghost1_idx == -1) {
791c376f201SBarry Smith PetscCheck(!lid_matched[lid1], PETSC_COMM_SELF, PETSC_ERR_SUP, "local %" PetscInt_FMT " is matched", gid1);
7928be712e4SBarry Smith lid_matched[lid1] = PETSC_TRUE; /* keep track of what we've done this round */
7938be712e4SBarry Smith PetscCall(PetscCDMoveAppend(agg_llists, lid0, lid1)); // takes lid1's list and appends to lid0's
7948be712e4SBarry Smith } else {
7958be712e4SBarry Smith /* add gid1 to list of ghost deleted by me -- I need their children */
7968be712e4SBarry Smith PetscMPIInt proc = lghost_pe[ghost1_idx];
797c376f201SBarry Smith PetscCheck(!lghost_matched[ghost1_idx], PETSC_COMM_SELF, PETSC_ERR_SUP, "ghost %" PetscInt_FMT " is matched", lghost_gid[ghost1_idx]);
7988be712e4SBarry Smith lghost_matched[ghost1_idx] = PETSC_TRUE;
7998be712e4SBarry Smith PetscCall(PetscCDAppendID(ghost_deleted_list, proc, ghost1_idx)); /* cache to send messages */
8008be712e4SBarry Smith PetscCall(PetscCDAppendID(ghost_deleted_list, proc, lid0));
8018be712e4SBarry Smith }
8028be712e4SBarry Smith lid_matched[lid0] = PETSC_TRUE; /* keep track of what we've done this round */
8038be712e4SBarry Smith /* set projection */
8048be712e4SBarry Smith PetscCall(MatSetValues(P, 1, &gid0, 1, &gid0, &one, INSERT_VALUES));
8058be712e4SBarry Smith PetscCall(MatSetValues(P, 1, &gid1, 1, &gid0, &one, INSERT_VALUES));
806c376f201SBarry 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));
8078be712e4SBarry Smith } /* matched */
8088be712e4SBarry Smith } /* edge loop */
809452acf8bSMark Adams PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
8108be712e4SBarry Smith if (isMPI) PetscCall(VecRestoreArrayRead(ghostMaxEdge, &lghost_max_ew));
8118be712e4SBarry Smith PetscCall(VecRestoreArrayRead(locMaxEdge, &lid_max_ew));
8128be712e4SBarry Smith // count active for test, latter, update deleted ghosts
8138be712e4SBarry Smith n_act_n[0] = nactive_edges;
8148be712e4SBarry Smith if (ghost_deleted_list) PetscCall(PetscCDCount(ghost_deleted_list, &n_act_n[2]));
8158be712e4SBarry Smith else n_act_n[2] = 0;
8168be712e4SBarry Smith PetscCall(PetscCDCount(agg_llists, &n_act_n[1]));
817462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(n_act_n, gn_act_n, 3, MPIU_INT, MPI_SUM, comm));
818c376f201SBarry 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]));
8198be712e4SBarry Smith /* deal with deleted ghost */
8208be712e4SBarry Smith if (isMPI) {
8218be712e4SBarry Smith PetscCDIntNd *pos;
8228be712e4SBarry Smith PetscInt *sbuffs1[REQ_BF_SIZE], ndel;
8238be712e4SBarry Smith PetscInt *sbuffs2[REQ_BF_SIZE];
8248be712e4SBarry Smith MPI_Status status;
8258be712e4SBarry Smith
8268be712e4SBarry Smith /* send deleted ghosts */
827c376f201SBarry Smith for (PetscInt proc_idx = 0; proc_idx < ncomm_procs; proc_idx++) {
8288be712e4SBarry Smith const PetscMPIInt proc = comm_procs[proc_idx];
8298be712e4SBarry Smith PetscInt *sbuff, *pt, scount;
8308be712e4SBarry Smith MPI_Request *request;
831c376f201SBarry Smith
8328be712e4SBarry Smith /* count ghosts */
8338be712e4SBarry Smith PetscCall(PetscCDCountAt(ghost_deleted_list, proc, &ndel));
8348be712e4SBarry Smith ndel /= 2; // two entries for each proc
8358be712e4SBarry Smith scount = 2 + 2 * ndel;
8368be712e4SBarry Smith PetscCall(PetscMalloc1(scount + request_size, &sbuff));
8378be712e4SBarry Smith /* save requests */
8388be712e4SBarry Smith sbuffs1[proc_idx] = sbuff;
8398be712e4SBarry Smith request = (MPI_Request *)sbuff;
840835f2295SStefano Zampini sbuff = pt = sbuff + request_size;
8418be712e4SBarry Smith /* write [ndel, proc, n*[gid1,gid0] */
8428be712e4SBarry Smith *pt++ = ndel; // number of deleted to send
8438be712e4SBarry Smith *pt++ = rank; // proc (not used)
8448be712e4SBarry Smith PetscCall(PetscCDGetHeadPos(ghost_deleted_list, proc, &pos));
8458be712e4SBarry Smith while (pos) {
8468be712e4SBarry Smith PetscInt lid0, ghost_idx, gid1;
847c376f201SBarry Smith
8488be712e4SBarry Smith PetscCall(PetscCDIntNdGetID(pos, &ghost_idx));
8498be712e4SBarry Smith gid1 = lghost_gid[ghost_idx];
8508be712e4SBarry Smith PetscCall(PetscCDGetNextPos(ghost_deleted_list, proc, &pos));
8518be712e4SBarry Smith PetscCall(PetscCDIntNdGetID(pos, &lid0));
8528be712e4SBarry Smith PetscCall(PetscCDGetNextPos(ghost_deleted_list, proc, &pos));
8538be712e4SBarry Smith *pt++ = gid1;
854c376f201SBarry Smith *pt++ = lid0 + Istart; // gid0
8558be712e4SBarry Smith }
856835f2295SStefano Zampini PetscCheck(pt - sbuff == scount, PETSC_COMM_SELF, PETSC_ERR_SUP, "sbuff-pt != scount: %zu", pt - sbuff);
8576497c311SBarry Smith /* MPIU_Isend: tag1 [ndel, proc, n*[gid1,gid0] ] */
8586497c311SBarry Smith PetscCallMPI(MPIU_Isend(sbuff, scount, MPIU_INT, proc, tag1, comm, request));
8598be712e4SBarry Smith PetscCall(PetscCDRemoveAllAt(ghost_deleted_list, proc)); // done with this list
8608be712e4SBarry Smith }
8618be712e4SBarry Smith /* receive deleted, send back partial aggregates, clear lists */
862c376f201SBarry Smith for (PetscInt proc_idx = 0; proc_idx < ncomm_procs; proc_idx++) {
8638be712e4SBarry Smith PetscCallMPI(MPI_Probe(comm_procs[proc_idx] /* MPI_ANY_SOURCE */, tag1, comm, &status));
8648be712e4SBarry Smith {
8658be712e4SBarry Smith PetscInt *pt, *pt2, *pt3, *sbuff, tmp;
8668be712e4SBarry Smith MPI_Request *request;
867c376f201SBarry Smith PetscMPIInt rcount, scount;
8688be712e4SBarry Smith const PetscMPIInt proc = status.MPI_SOURCE;
869c376f201SBarry Smith
8708be712e4SBarry Smith PetscCallMPI(MPI_Get_count(&status, MPIU_INT, &rcount));
8718be712e4SBarry Smith if (rcount > rbuff_sz) {
8728be712e4SBarry Smith if (rbuff) PetscCall(PetscFree(rbuff));
8738be712e4SBarry Smith PetscCall(PetscMalloc1(rcount, &rbuff));
8748be712e4SBarry Smith rbuff_sz = rcount;
8758be712e4SBarry Smith }
8768be712e4SBarry Smith /* MPI_Recv: tag1 [ndel, proc, ndel*[gid1,gid0] ] */
8778be712e4SBarry Smith PetscCallMPI(MPI_Recv(rbuff, rcount, MPIU_INT, proc, tag1, comm, &status));
8788be712e4SBarry Smith /* read and count sends *[lid0, n, n*[gid] ] */
8798be712e4SBarry Smith pt = rbuff;
8808be712e4SBarry Smith scount = 0;
8818be712e4SBarry Smith ndel = *pt++; // number of deleted to recv
8828be712e4SBarry Smith tmp = *pt++; // proc (not used)
8838be712e4SBarry Smith while (ndel--) {
884c376f201SBarry Smith PetscInt gid1 = *pt++, lid1 = gid1 - Istart;
8856497c311SBarry Smith PetscInt gh_gid0 = *pt++; // gid on other proc (not used here to count)
886c376f201SBarry Smith
887c376f201SBarry Smith PetscCheck(lid1 >= 0 && lid1 < nloc, PETSC_COMM_SELF, PETSC_ERR_SUP, "received ghost deleted %" PetscInt_FMT, gid1);
8886497c311SBarry Smith PetscCheck(!lid_matched[lid1], PETSC_COMM_SELF, PETSC_ERR_PLIB, "%" PetscInt_FMT ") received matched local gid %" PetscInt_FMT ",%" PetscInt_FMT ", with ghost (lid) %" PetscInt_FMT " from proc %d", sub_it, gid1, gh_gid0, tmp, proc);
8898be712e4SBarry Smith lid_matched[lid1] = PETSC_TRUE; /* keep track of what we've done this round */
8908be712e4SBarry Smith PetscCall(PetscCDCountAt(agg_llists, lid1, &tmp)); // n
8918be712e4SBarry Smith scount += tmp + 2; // lid0, n, n*[gid]
8928be712e4SBarry Smith }
893c376f201SBarry Smith PetscCheck((pt - rbuff) == (ptrdiff_t)rcount, PETSC_COMM_SELF, PETSC_ERR_SUP, "receive buffer size != num read: %zu; rcount: %d", pt - rbuff, rcount);
8948be712e4SBarry Smith /* send tag2: *[gid0, n, n*[gid] ] */
8958be712e4SBarry Smith PetscCall(PetscMalloc1(scount + request_size, &sbuff));
8968be712e4SBarry Smith sbuffs2[proc_idx] = sbuff; /* cache request */
8978be712e4SBarry Smith request = (MPI_Request *)sbuff;
8988be712e4SBarry Smith pt2 = sbuff = sbuff + request_size;
8998be712e4SBarry Smith // read again: n, proc, n*[gid1,gid0]
9008be712e4SBarry Smith pt = rbuff;
9018be712e4SBarry Smith ndel = *pt++;
9028be712e4SBarry Smith tmp = *pt++; // proc (not used)
9038be712e4SBarry Smith while (ndel--) {
904c376f201SBarry Smith PetscInt gid1 = *pt++, lid1 = gid1 - Istart, gh_gid0 = *pt++;
905c376f201SBarry Smith
9068be712e4SBarry Smith /* write [gid0, aggSz, aggSz[gid] ] */
9078be712e4SBarry Smith *pt2++ = gh_gid0;
9088be712e4SBarry Smith pt3 = pt2++; /* save pointer for later */
9098be712e4SBarry Smith PetscCall(PetscCDGetHeadPos(agg_llists, lid1, &pos));
9108be712e4SBarry Smith while (pos) {
9118be712e4SBarry Smith PetscInt gid;
912c376f201SBarry Smith
9138be712e4SBarry Smith PetscCall(PetscCDIntNdGetID(pos, &gid));
9148be712e4SBarry Smith PetscCall(PetscCDGetNextPos(agg_llists, lid1, &pos));
9158be712e4SBarry Smith *pt2++ = gid;
9168be712e4SBarry Smith }
917835f2295SStefano Zampini PetscCall(PetscIntCast(pt2 - pt3 - 1, pt3));
9188be712e4SBarry Smith /* clear list */
9198be712e4SBarry Smith PetscCall(PetscCDRemoveAllAt(agg_llists, lid1));
9208be712e4SBarry Smith }
921c376f201SBarry Smith PetscCheck((pt2 - sbuff) == (ptrdiff_t)scount, PETSC_COMM_SELF, PETSC_ERR_SUP, "buffer size != num write: %zu %d", pt2 - sbuff, scount);
9226497c311SBarry Smith /* MPIU_Isend: requested data tag2 *[lid0, n, n*[gid1] ] */
9236497c311SBarry Smith PetscCallMPI(MPIU_Isend(sbuff, scount, MPIU_INT, proc, tag2, comm, request));
9248be712e4SBarry Smith }
9258be712e4SBarry Smith } // proc_idx
9268be712e4SBarry Smith /* receive tag2 *[gid0, n, n*[gid] ] */
927c376f201SBarry Smith for (PetscMPIInt proc_idx = 0; proc_idx < ncomm_procs; proc_idx++) {
9288be712e4SBarry Smith PetscMPIInt proc;
9298be712e4SBarry Smith PetscInt *pt;
9308be712e4SBarry Smith int rcount;
931c376f201SBarry Smith
9328be712e4SBarry Smith PetscCallMPI(MPI_Probe(comm_procs[proc_idx] /* MPI_ANY_SOURCE */, tag2, comm, &status));
9338be712e4SBarry Smith PetscCallMPI(MPI_Get_count(&status, MPIU_INT, &rcount));
9348be712e4SBarry Smith if (rcount > rbuff_sz) {
9358be712e4SBarry Smith if (rbuff) PetscCall(PetscFree(rbuff));
9368be712e4SBarry Smith PetscCall(PetscMalloc1(rcount, &rbuff));
9378be712e4SBarry Smith rbuff_sz = rcount;
9388be712e4SBarry Smith }
9398be712e4SBarry Smith proc = status.MPI_SOURCE;
9408be712e4SBarry Smith /* MPI_Recv: tag1 [n, proc, n*[gid1,lid0] ] */
9418be712e4SBarry Smith PetscCallMPI(MPI_Recv(rbuff, rcount, MPIU_INT, proc, tag2, comm, &status));
9428be712e4SBarry Smith pt = rbuff;
9438be712e4SBarry Smith while (pt - rbuff < rcount) {
9448be712e4SBarry Smith PetscInt gid0 = *pt++, n = *pt++;
945c376f201SBarry Smith
9468be712e4SBarry Smith while (n--) {
9478be712e4SBarry Smith PetscInt gid1 = *pt++;
948c376f201SBarry Smith
949c376f201SBarry Smith PetscCall(PetscCDAppendID(agg_llists, gid0 - Istart, gid1));
9508be712e4SBarry Smith }
9518be712e4SBarry Smith }
952c376f201SBarry Smith PetscCheck((pt - rbuff) == (ptrdiff_t)rcount, PETSC_COMM_SELF, PETSC_ERR_SUP, "recv buffer size != num read: %zu %d", pt - rbuff, rcount);
9538be712e4SBarry Smith }
9548be712e4SBarry Smith /* wait for tag1 isends */
955c376f201SBarry Smith for (PetscMPIInt proc_idx = 0; proc_idx < ncomm_procs; proc_idx++) {
956c376f201SBarry Smith MPI_Request *request = (MPI_Request *)sbuffs1[proc_idx];
957c376f201SBarry Smith
9588be712e4SBarry Smith PetscCallMPI(MPI_Wait(request, &status));
9598be712e4SBarry Smith PetscCall(PetscFree(sbuffs1[proc_idx]));
9608be712e4SBarry Smith }
9618be712e4SBarry Smith /* wait for tag2 isends */
962c376f201SBarry Smith for (PetscMPIInt proc_idx = 0; proc_idx < ncomm_procs; proc_idx++) {
9638be712e4SBarry Smith MPI_Request *request = (MPI_Request *)sbuffs2[proc_idx];
964c376f201SBarry Smith
9658be712e4SBarry Smith PetscCallMPI(MPI_Wait(request, &status));
9668be712e4SBarry Smith PetscCall(PetscFree(sbuffs2[proc_idx]));
9678be712e4SBarry Smith }
9688be712e4SBarry Smith } /* MPI */
9698be712e4SBarry Smith /* set 'lghost_matched' - use locMaxEdge, ghostMaxEdge (recomputed next) */
9708be712e4SBarry Smith if (isMPI) {
9718be712e4SBarry Smith const PetscScalar *sbuff;
972c376f201SBarry Smith
973c376f201SBarry Smith for (PetscInt kk = 0, gid = Istart; kk < nloc; kk++, gid++) {
9748be712e4SBarry Smith PetscScalar vval = lid_matched[kk] ? 1.0 : 0.0;
975c376f201SBarry Smith
9768be712e4SBarry Smith PetscCall(VecSetValues(locMaxEdge, 1, &gid, &vval, INSERT_VALUES)); /* set with GID */
9778be712e4SBarry Smith }
9788be712e4SBarry Smith PetscCall(VecAssemblyBegin(locMaxEdge));
9798be712e4SBarry Smith PetscCall(VecAssemblyEnd(locMaxEdge));
9808be712e4SBarry Smith PetscCall(VecScatterBegin(mpimat->Mvctx, locMaxEdge, ghostMaxEdge, INSERT_VALUES, SCATTER_FORWARD));
9818be712e4SBarry Smith PetscCall(VecScatterEnd(mpimat->Mvctx, locMaxEdge, ghostMaxEdge, INSERT_VALUES, SCATTER_FORWARD));
9828be712e4SBarry Smith PetscCall(VecGetArrayRead(ghostMaxEdge, &sbuff));
983*ac530a7eSPierre Jolivet for (PetscInt kk = 0; kk < num_ghosts; kk++) lghost_matched[kk] = (PetscBool)(PetscRealPart(sbuff[kk]) != 0.0);
9848be712e4SBarry Smith PetscCall(VecRestoreArrayRead(ghostMaxEdge, &sbuff));
9858be712e4SBarry Smith }
9868be712e4SBarry Smith /* compute 'locMaxEdge' inside sub iteration b/c max weight can drop as neighbors are matched */
987c376f201SBarry Smith for (PetscInt kk = 0, gid = Istart; kk < nloc; kk++, gid++) {
9888be712e4SBarry Smith PetscReal max_e = 0., tt;
9898be712e4SBarry Smith PetscScalar vval;
9908be712e4SBarry Smith const PetscInt lid = kk;
991c376f201SBarry Smith PetscMPIInt max_pe = rank, pe, n;
992c376f201SBarry Smith
9938be712e4SBarry Smith ii = matA->i;
9946497c311SBarry Smith PetscCall(PetscMPIIntCast(ii[lid + 1] - ii[lid], &n));
9958e3a54c0SPierre Jolivet aj = PetscSafePointerPlusOffset(matA->j, ii[lid]);
9968e3a54c0SPierre Jolivet ap = PetscSafePointerPlusOffset(matA->a, ii[lid]);
9976497c311SBarry Smith for (PetscMPIInt jj = 0; jj < n; jj++) {
9988be712e4SBarry Smith PetscInt lidj = aj[jj];
999c376f201SBarry Smith
10008be712e4SBarry Smith if (lid_matched[lidj]) continue; /* this is new - can change local max */
10018be712e4SBarry Smith if (lidj != lid && PetscRealPart(ap[jj]) > max_e) max_e = PetscRealPart(ap[jj]);
10028be712e4SBarry Smith }
10038be712e4SBarry Smith if (lid_cprowID && (ix = lid_cprowID[lid]) != -1) { /* if I have any ghost neighbors */
10048be712e4SBarry Smith ii = matB->compressedrow.i;
10056497c311SBarry Smith PetscCall(PetscMPIIntCast(ii[ix + 1] - ii[ix], &n));
10068be712e4SBarry Smith ap = matB->a + ii[ix];
10078be712e4SBarry Smith aj = matB->j + ii[ix];
10086497c311SBarry Smith for (PetscMPIInt jj = 0; jj < n; jj++) {
10098be712e4SBarry Smith PetscInt lidj = aj[jj];
1010c376f201SBarry Smith
10118be712e4SBarry Smith if (lghost_matched[lidj]) continue;
10128be712e4SBarry Smith if ((tt = PetscRealPart(ap[jj])) > max_e) max_e = tt;
10138be712e4SBarry Smith }
10148be712e4SBarry Smith }
1015835f2295SStefano Zampini vval = max_e;
10168be712e4SBarry Smith PetscCall(VecSetValues(locMaxEdge, 1, &gid, &vval, INSERT_VALUES)); /* set with GID */
10178be712e4SBarry Smith // max PE with max edge
10188be712e4SBarry Smith if (lid_cprowID && (ix = lid_cprowID[lid]) != -1) { /* if I have any ghost neighbors */
10198be712e4SBarry Smith ii = matB->compressedrow.i;
10206497c311SBarry Smith PetscCall(PetscMPIIntCast(ii[ix + 1] - ii[ix], &n));
10218be712e4SBarry Smith ap = matB->a + ii[ix];
10228be712e4SBarry Smith aj = matB->j + ii[ix];
1023c376f201SBarry Smith for (PetscInt jj = 0; jj < n; jj++) {
10248be712e4SBarry Smith PetscInt lidj = aj[jj];
1025c376f201SBarry Smith
10268be712e4SBarry Smith if (lghost_matched[lidj]) continue;
1027*ac530a7eSPierre Jolivet if ((pe = lghost_pe[aj[jj]]) > max_pe && PetscRealPart(ap[jj]) >= max_e - MY_MEPS) max_pe = pe;
10288be712e4SBarry Smith }
10298be712e4SBarry Smith }
1030835f2295SStefano Zampini vval = max_pe;
10318be712e4SBarry Smith PetscCall(VecSetValues(locMaxPE, 1, &gid, &vval, INSERT_VALUES));
10328be712e4SBarry Smith }
10338be712e4SBarry Smith PetscCall(VecAssemblyBegin(locMaxEdge));
10348be712e4SBarry Smith PetscCall(VecAssemblyEnd(locMaxEdge));
10358be712e4SBarry Smith PetscCall(VecAssemblyBegin(locMaxPE));
10368be712e4SBarry Smith PetscCall(VecAssemblyEnd(locMaxPE));
10378be712e4SBarry Smith /* compute 'lghost_max_ew' and 'lghost_max_pe' to get ready for next iteration*/
10388be712e4SBarry Smith if (isMPI) {
10398be712e4SBarry Smith const PetscScalar *buf;
1040c376f201SBarry Smith
10418be712e4SBarry Smith PetscCall(VecScatterBegin(mpimat->Mvctx, locMaxEdge, ghostMaxEdge, INSERT_VALUES, SCATTER_FORWARD));
10428be712e4SBarry Smith PetscCall(VecScatterEnd(mpimat->Mvctx, locMaxEdge, ghostMaxEdge, INSERT_VALUES, SCATTER_FORWARD));
10438be712e4SBarry Smith PetscCall(VecScatterBegin(mpimat->Mvctx, locMaxPE, ghostMaxPE, INSERT_VALUES, SCATTER_FORWARD));
10448be712e4SBarry Smith PetscCall(VecScatterEnd(mpimat->Mvctx, locMaxPE, ghostMaxPE, INSERT_VALUES, SCATTER_FORWARD));
10458be712e4SBarry Smith PetscCall(VecGetArrayRead(ghostMaxPE, &buf));
10463a7d0413SPierre Jolivet for (PetscInt kk = 0; kk < num_ghosts; kk++) lghost_max_pe[kk] = (PetscMPIInt)PetscRealPart(buf[kk]); // the MAX proc of the ghost now
10478be712e4SBarry Smith PetscCall(VecRestoreArrayRead(ghostMaxPE, &buf));
10488be712e4SBarry Smith }
10498be712e4SBarry Smith // if no active edges, stop
10508be712e4SBarry Smith if (gn_act_n[0] < 1) break;
10518be712e4SBarry Smith // inc and check (self stopping iteration
1052c376f201SBarry 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);
10538be712e4SBarry Smith sub_it++;
1054c376f201SBarry 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);
1055452acf8bSMark Adams old_num_edge = gn_act_n[0];
10568be712e4SBarry Smith } /* sub_it loop */
10578be712e4SBarry Smith /* clean up iteration */
10588be712e4SBarry Smith PetscCall(PetscFree(Edges));
1059c376f201SBarry Smith if (isMPI) { // can be hoisted
10608be712e4SBarry Smith PetscCall(VecRestoreArrayRead(ghostMaxEdge, &lghost_max_ew));
10618be712e4SBarry Smith PetscCall(VecDestroy(&ghostMaxEdge));
10628be712e4SBarry Smith PetscCall(VecDestroy(&ghostMaxPE));
1063c376f201SBarry Smith PetscCall(PetscFree4(lghost_matched, lghost_pe, lghost_gid, lghost_max_pe));
10648be712e4SBarry Smith }
10658be712e4SBarry Smith PetscCall(VecDestroy(&locMaxEdge));
10668be712e4SBarry Smith PetscCall(VecDestroy(&locMaxPE));
10678be712e4SBarry Smith /* create next graph */
10688be712e4SBarry Smith {
10698be712e4SBarry Smith Vec diag;
1070c376f201SBarry Smith
10718be712e4SBarry Smith /* add identity for unmatched vertices so they stay alive */
1072c376f201SBarry Smith for (PetscInt kk = 0, gid1, gid = Istart; kk < nloc; kk++, gid++) {
10738be712e4SBarry Smith if (!lid_matched[kk]) {
10748be712e4SBarry Smith const PetscInt lid = kk;
10758be712e4SBarry Smith PetscCDIntNd *pos;
1076c376f201SBarry Smith
10778be712e4SBarry Smith PetscCall(PetscCDGetHeadPos(agg_llists, lid, &pos));
1078c376f201SBarry Smith PetscCheck(pos, PETSC_COMM_SELF, PETSC_ERR_PLIB, "empty list in singleton: %" PetscInt_FMT, gid);
10798be712e4SBarry Smith PetscCall(PetscCDIntNdGetID(pos, &gid1));
1080c376f201SBarry Smith PetscCheck(gid1 == gid, PETSC_COMM_SELF, PETSC_ERR_PLIB, "first in list (%" PetscInt_FMT ") in singleton not %" PetscInt_FMT, gid1, gid);
10818be712e4SBarry Smith PetscCall(MatSetValues(P, 1, &gid, 1, &gid, &one, INSERT_VALUES));
10828be712e4SBarry Smith }
10838be712e4SBarry Smith }
10848be712e4SBarry Smith PetscCall(MatAssemblyBegin(P, MAT_FINAL_ASSEMBLY));
10858be712e4SBarry Smith PetscCall(MatAssemblyEnd(P, MAT_FINAL_ASSEMBLY));
10868be712e4SBarry Smith
10878be712e4SBarry Smith /* project to make new graph with collapsed edges */
10888be712e4SBarry Smith PetscCall(MatPtAP(cMat, P, MAT_INITIAL_MATRIX, 1.0, &tMat));
10898be712e4SBarry Smith PetscCall(MatDestroy(&P));
10908be712e4SBarry Smith PetscCall(MatDestroy(&cMat));
10918be712e4SBarry Smith cMat = tMat;
10928be712e4SBarry Smith PetscCall(MatCreateVecs(cMat, &diag, NULL));
109353673ba5SSatish Balay PetscCall(MatGetDiagonal(cMat, diag));
10948be712e4SBarry Smith PetscCall(VecReciprocal(diag));
10958be712e4SBarry Smith PetscCall(VecSqrtAbs(diag));
10968be712e4SBarry Smith PetscCall(MatDiagonalScale(cMat, diag, diag));
10978be712e4SBarry Smith PetscCall(VecDestroy(&diag));
10988be712e4SBarry Smith }
10998be712e4SBarry Smith } /* coarsen iterator */
11008be712e4SBarry Smith
11018be712e4SBarry 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 */
11028be712e4SBarry Smith if (size > 1) {
11038be712e4SBarry Smith Mat mat;
11048be712e4SBarry Smith PetscCDIntNd *pos;
11058be712e4SBarry Smith PetscInt NN, MM, jj = 0, mxsz = 0;
11068be712e4SBarry Smith
1107c376f201SBarry Smith for (PetscInt kk = 0; kk < nloc; kk++) {
11088be712e4SBarry Smith PetscCall(PetscCDCountAt(agg_llists, kk, &jj));
11098be712e4SBarry Smith if (jj > mxsz) mxsz = jj;
11108be712e4SBarry Smith }
11118be712e4SBarry Smith PetscCall(MatGetSize(a_Gmat, &MM, &NN));
11128be712e4SBarry Smith if (mxsz > MM - nloc) mxsz = MM - nloc;
11138be712e4SBarry Smith /* matrix of ghost adj for square graph */
11148be712e4SBarry Smith PetscCall(MatCreateAIJ(comm, nloc, nloc, PETSC_DETERMINE, PETSC_DETERMINE, 0, NULL, mxsz, NULL, &mat));
1115c376f201SBarry Smith for (PetscInt lid = 0, gid = Istart; lid < nloc; lid++, gid++) {
11168be712e4SBarry Smith PetscCall(PetscCDGetHeadPos(agg_llists, lid, &pos));
11178be712e4SBarry Smith while (pos) {
11188be712e4SBarry Smith PetscInt gid1;
1119c376f201SBarry Smith
11208be712e4SBarry Smith PetscCall(PetscCDIntNdGetID(pos, &gid1));
11218be712e4SBarry Smith PetscCall(PetscCDGetNextPos(agg_llists, lid, &pos));
1122c376f201SBarry Smith if (gid1 < Istart || gid1 >= Istart + nloc) PetscCall(MatSetValues(mat, 1, &gid, 1, &gid1, &one, ADD_VALUES));
11238be712e4SBarry Smith }
11248be712e4SBarry Smith }
11258be712e4SBarry Smith PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
11268be712e4SBarry Smith PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
11278be712e4SBarry Smith PetscCall(PetscCDSetMat(agg_llists, mat));
11288be712e4SBarry Smith PetscCall(PetscCDDestroy(ghost_deleted_list));
11298be712e4SBarry Smith if (rbuff_sz) PetscCall(PetscFree(rbuff)); // always true
11308be712e4SBarry Smith }
11318be712e4SBarry Smith // move BCs into some node
11328be712e4SBarry Smith if (bc_list) {
11338be712e4SBarry Smith PetscCDIntNd *pos;
1134c376f201SBarry Smith
11358be712e4SBarry Smith PetscCall(PetscCDGetHeadPos(bc_list, 0, &pos));
11368be712e4SBarry Smith while (pos) {
11378be712e4SBarry Smith PetscInt gid1;
1138c376f201SBarry Smith
11398be712e4SBarry Smith PetscCall(PetscCDIntNdGetID(pos, &gid1));
11408be712e4SBarry Smith PetscCall(PetscCDGetNextPos(bc_list, 0, &pos));
11418be712e4SBarry Smith PetscCall(PetscCDAppendID(agg_llists, bc_agg, gid1));
11428be712e4SBarry Smith }
11438be712e4SBarry Smith PetscCall(PetscCDRemoveAllAt(bc_list, 0));
11448be712e4SBarry Smith PetscCall(PetscCDDestroy(bc_list));
11458be712e4SBarry Smith }
11468be712e4SBarry Smith {
11478be712e4SBarry Smith // check sizes -- all vertices must get in graph
11488be712e4SBarry Smith PetscInt sz, globalsz, MM;
1149c376f201SBarry Smith
11508be712e4SBarry Smith PetscCall(MatGetSize(a_Gmat, &MM, NULL));
11518be712e4SBarry Smith PetscCall(PetscCDCount(agg_llists, &sz));
1152462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&sz, &globalsz, 1, MPIU_INT, MPI_SUM, comm));
115357508eceSPierre Jolivet PetscCheck(MM == globalsz, comm, PETSC_ERR_SUP, "lost %" PetscInt_FMT " equations ?", MM - globalsz);
11548be712e4SBarry Smith }
11558be712e4SBarry Smith // cleanup
11568be712e4SBarry Smith PetscCall(MatDestroy(&cMat));
1157c376f201SBarry Smith PetscCall(PetscFree3(lid_matched, lid_cprowID, lid_max_pe));
11588be712e4SBarry Smith PetscCall(ISDestroy(&info_is));
11598be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
11608be712e4SBarry Smith }
11618be712e4SBarry Smith
11628be712e4SBarry Smith /*
11638be712e4SBarry Smith HEM coarsen, simple greedy.
11648be712e4SBarry Smith */
MatCoarsenApply_HEM(MatCoarsen coarse)11658be712e4SBarry Smith static PetscErrorCode MatCoarsenApply_HEM(MatCoarsen coarse)
11668be712e4SBarry Smith {
11678be712e4SBarry Smith Mat mat = coarse->graph;
11688be712e4SBarry Smith
11698be712e4SBarry Smith PetscFunctionBegin;
11708be712e4SBarry Smith PetscCall(MatCoarsenApply_HEM_private(mat, coarse->max_it, coarse->threshold, &coarse->agg_lists));
11718be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
11728be712e4SBarry Smith }
11738be712e4SBarry Smith
MatCoarsenView_HEM(MatCoarsen coarse,PetscViewer viewer)11748be712e4SBarry Smith static PetscErrorCode MatCoarsenView_HEM(MatCoarsen coarse, PetscViewer viewer)
11758be712e4SBarry Smith {
11768be712e4SBarry Smith PetscMPIInt rank;
11779f196a02SMartin Diehl PetscBool isascii;
11788be712e4SBarry Smith
11798be712e4SBarry Smith PetscFunctionBegin;
11808be712e4SBarry Smith PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)coarse), &rank));
11819f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
11829f196a02SMartin Diehl if (isascii) {
11838be712e4SBarry Smith PetscCDIntNd *pos, *pos2;
1184e0b7e82fSBarry Smith PetscViewerFormat format;
1185e0b7e82fSBarry Smith
1186c376f201SBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " matching steps with threshold = %g\n", coarse->max_it, (double)coarse->threshold));
1187e0b7e82fSBarry Smith PetscCall(PetscViewerGetFormat(viewer, &format));
1188e0b7e82fSBarry Smith if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1189bfa5bdfcSBarry Smith if (coarse->agg_lists) {
11908be712e4SBarry Smith PetscCall(PetscViewerASCIIPushSynchronized(viewer));
11918be712e4SBarry Smith for (PetscInt kk = 0; kk < coarse->agg_lists->size; kk++) {
11928be712e4SBarry Smith PetscCall(PetscCDGetHeadPos(coarse->agg_lists, kk, &pos));
1193c376f201SBarry Smith if ((pos2 = pos)) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "selected local %" PetscInt_FMT ": ", kk));
11948be712e4SBarry Smith while (pos) {
11958be712e4SBarry Smith PetscInt gid1;
1196c376f201SBarry Smith
11978be712e4SBarry Smith PetscCall(PetscCDIntNdGetID(pos, &gid1));
11988be712e4SBarry Smith PetscCall(PetscCDGetNextPos(coarse->agg_lists, kk, &pos));
1199c376f201SBarry Smith PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " %" PetscInt_FMT " ", gid1));
12008be712e4SBarry Smith }
12018be712e4SBarry Smith if (pos2) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n"));
12028be712e4SBarry Smith }
12038be712e4SBarry Smith PetscCall(PetscViewerFlush(viewer));
12048be712e4SBarry Smith PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1205bfa5bdfcSBarry Smith } else {
1206bfa5bdfcSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, " HEM aggregator lists are not available\n"));
1207bfa5bdfcSBarry Smith }
12088be712e4SBarry Smith }
1209e0b7e82fSBarry Smith }
12108be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
12118be712e4SBarry Smith }
12128be712e4SBarry Smith
1213c376f201SBarry Smith /*MC
1214c376f201SBarry Smith MATCOARSENHEM - A coarsener that uses HEM a simple greedy coarsener
1215c376f201SBarry Smith
1216c376f201SBarry Smith Level: beginner
1217c376f201SBarry Smith
1218c376f201SBarry Smith .seealso: `MatCoarsen`, `MatCoarsenMISKSetDistance()`, `MatCoarsenApply()`, `MatCoarsenSetType()`, `MatCoarsenType`, `MatCoarsenCreate()`, `MATCOARSENMISK`, `MATCOARSENMIS`
1219c376f201SBarry Smith M*/
1220c376f201SBarry Smith
MatCoarsenCreate_HEM(MatCoarsen coarse)12218be712e4SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenCreate_HEM(MatCoarsen coarse)
12228be712e4SBarry Smith {
12238be712e4SBarry Smith PetscFunctionBegin;
12248be712e4SBarry Smith coarse->ops->apply = MatCoarsenApply_HEM;
12258be712e4SBarry Smith coarse->ops->view = MatCoarsenView_HEM;
12268be712e4SBarry Smith coarse->max_it = 4;
12278be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
12288be712e4SBarry Smith }
1229