xref: /petsc/src/mat/graphops/coarsen/impls/hem/hem.c (revision e0b7e82fd3cf27fce84cc3e37e8d70a5c36a2d4e)
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;
398be712e4SBarry Smith     n                  = n->next;
408be712e4SBarry Smith     PetscCall(PetscFree(lstn));
418be712e4SBarry Smith   }
428be712e4SBarry Smith   if (ail->pool_list.array) PetscCall(PetscFree(ail->pool_list.array));
438be712e4SBarry Smith   PetscCall(PetscFree(ail->array));
448be712e4SBarry Smith   if (ail->mat) PetscCall(MatDestroy(&ail->mat));
458be712e4SBarry Smith   /* delete this (+agg+pool array) */
468be712e4SBarry Smith   PetscCall(PetscFree(ail));
478be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
488be712e4SBarry Smith }
498be712e4SBarry Smith 
508be712e4SBarry Smith /* PetscCDSetChunkSize
518be712e4SBarry Smith  */
528be712e4SBarry Smith PetscErrorCode PetscCDSetChunkSize(PetscCoarsenData *ail, PetscInt a_sz)
538be712e4SBarry Smith {
548be712e4SBarry Smith   PetscFunctionBegin;
558be712e4SBarry Smith   ail->chk_sz = a_sz;
568be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
578be712e4SBarry Smith }
588be712e4SBarry Smith 
598be712e4SBarry Smith /*  PetscCDGetNewNode
608be712e4SBarry Smith  */
618be712e4SBarry Smith static PetscErrorCode PetscCDGetNewNode(PetscCoarsenData *ail, PetscCDIntNd **a_out, PetscInt a_id)
628be712e4SBarry Smith {
638be712e4SBarry Smith   PetscFunctionBegin;
648be712e4SBarry Smith   *a_out = NULL; /* squelch -Wmaybe-uninitialized */
658be712e4SBarry Smith   if (ail->extra_nodes) {
668be712e4SBarry Smith     PetscCDIntNd *node = ail->extra_nodes;
678be712e4SBarry Smith     ail->extra_nodes   = node->next;
688be712e4SBarry Smith     node->gid          = a_id;
698be712e4SBarry Smith     node->next         = NULL;
708be712e4SBarry Smith     *a_out             = node;
718be712e4SBarry Smith   } else {
728be712e4SBarry Smith     if (!ail->pool_list.array) {
738be712e4SBarry Smith       if (!ail->chk_sz) ail->chk_sz = 10; /* use a chuck size of ail->size? */
748be712e4SBarry Smith       PetscCall(PetscMalloc1(ail->chk_sz, &ail->pool_list.array));
758be712e4SBarry Smith       ail->new_node       = ail->pool_list.array;
768be712e4SBarry Smith       ail->new_left       = ail->chk_sz;
778be712e4SBarry Smith       ail->new_node->next = NULL;
788be712e4SBarry Smith     } else if (!ail->new_left) {
798be712e4SBarry Smith       PetscCDArrNd *node;
808be712e4SBarry Smith       PetscCall(PetscMalloc(ail->chk_sz * sizeof(PetscCDIntNd) + sizeof(PetscCDArrNd), &node));
818be712e4SBarry Smith       node->array         = (PetscCDIntNd *)(node + 1);
828be712e4SBarry Smith       node->next          = ail->pool_list.next;
838be712e4SBarry Smith       ail->pool_list.next = node;
848be712e4SBarry Smith       ail->new_left       = ail->chk_sz;
858be712e4SBarry Smith       ail->new_node       = node->array;
868be712e4SBarry Smith     }
878be712e4SBarry Smith     ail->new_node->gid  = a_id;
888be712e4SBarry Smith     ail->new_node->next = NULL;
898be712e4SBarry Smith     *a_out              = ail->new_node++;
908be712e4SBarry Smith     ail->new_left--;
918be712e4SBarry Smith   }
928be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
938be712e4SBarry Smith }
948be712e4SBarry Smith 
958be712e4SBarry Smith /* PetscCDIntNdSetID
968be712e4SBarry Smith  */
978be712e4SBarry Smith PetscErrorCode PetscCDIntNdSetID(PetscCDIntNd *a_this, PetscInt a_id)
988be712e4SBarry Smith {
998be712e4SBarry Smith   PetscFunctionBegin;
1008be712e4SBarry Smith   a_this->gid = a_id;
1018be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
1028be712e4SBarry Smith }
1038be712e4SBarry Smith 
1048be712e4SBarry Smith /* PetscCDIntNdGetID
1058be712e4SBarry Smith  */
1068be712e4SBarry Smith PetscErrorCode PetscCDIntNdGetID(const PetscCDIntNd *a_this, PetscInt *a_gid)
1078be712e4SBarry Smith {
1088be712e4SBarry Smith   PetscFunctionBegin;
1098be712e4SBarry Smith   *a_gid = a_this->gid;
1108be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
1118be712e4SBarry Smith }
1128be712e4SBarry Smith 
1138be712e4SBarry Smith /* PetscCDGetHeadPos
1148be712e4SBarry Smith  */
1158be712e4SBarry Smith PetscErrorCode PetscCDGetHeadPos(const PetscCoarsenData *ail, PetscInt a_idx, PetscCDIntNd **pos)
1168be712e4SBarry Smith {
1178be712e4SBarry Smith   PetscFunctionBegin;
1188be712e4SBarry Smith   PetscCheck(a_idx < ail->size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "a_idx >= ail->size: a_idx=%" PetscInt_FMT ".", a_idx);
1198be712e4SBarry Smith   *pos = ail->array[a_idx];
1208be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
1218be712e4SBarry Smith }
1228be712e4SBarry Smith 
1238be712e4SBarry Smith /* PetscCDGetNextPos
1248be712e4SBarry Smith  */
1258be712e4SBarry Smith PetscErrorCode PetscCDGetNextPos(const PetscCoarsenData *ail, PetscInt l_idx, PetscCDIntNd **pos)
1268be712e4SBarry Smith {
1278be712e4SBarry Smith   PetscFunctionBegin;
128f4f49eeaSPierre Jolivet   PetscCheck(*pos, PETSC_COMM_SELF, PETSC_ERR_PLIB, "NULL input position.");
1298be712e4SBarry Smith   *pos = (*pos)->next;
1308be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
1318be712e4SBarry Smith }
1328be712e4SBarry Smith 
1338be712e4SBarry Smith /* PetscCDAppendID
1348be712e4SBarry Smith  */
1358be712e4SBarry Smith PetscErrorCode PetscCDAppendID(PetscCoarsenData *ail, PetscInt a_idx, PetscInt a_id)
1368be712e4SBarry Smith {
1378be712e4SBarry Smith   PetscCDIntNd *n, *n2;
1388be712e4SBarry Smith 
1398be712e4SBarry Smith   PetscFunctionBegin;
1408be712e4SBarry Smith   PetscCall(PetscCDGetNewNode(ail, &n, a_id));
1418be712e4SBarry Smith   PetscCheck(a_idx < ail->size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Index %" PetscInt_FMT " out of range.", a_idx);
1428be712e4SBarry Smith   if (!(n2 = ail->array[a_idx])) ail->array[a_idx] = n;
1438be712e4SBarry Smith   else {
1448be712e4SBarry Smith     do {
1458be712e4SBarry Smith       if (!n2->next) {
1468be712e4SBarry Smith         n2->next = n;
1478be712e4SBarry Smith         PetscCheck(!n->next, PETSC_COMM_SELF, PETSC_ERR_PLIB, "n should not have a next");
1488be712e4SBarry Smith         break;
1498be712e4SBarry Smith       }
1508be712e4SBarry Smith       n2 = n2->next;
1518be712e4SBarry Smith     } while (n2);
1528be712e4SBarry Smith     PetscCheck(n2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "n2 should be non-null");
1538be712e4SBarry Smith   }
1548be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
1558be712e4SBarry Smith }
1568be712e4SBarry Smith 
1578be712e4SBarry Smith /* PetscCDAppendNode
1588be712e4SBarry Smith  */
1598be712e4SBarry Smith PetscErrorCode PetscCDAppendNode(PetscCoarsenData *ail, PetscInt a_idx, PetscCDIntNd *a_n)
1608be712e4SBarry Smith {
1618be712e4SBarry Smith   PetscCDIntNd *n2;
1628be712e4SBarry Smith 
1638be712e4SBarry Smith   PetscFunctionBegin;
1648be712e4SBarry Smith   PetscCheck(a_idx < ail->size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Index %" PetscInt_FMT " out of range.", a_idx);
1658be712e4SBarry Smith   if (!(n2 = ail->array[a_idx])) ail->array[a_idx] = a_n;
1668be712e4SBarry Smith   else {
1678be712e4SBarry Smith     do {
1688be712e4SBarry Smith       if (!n2->next) {
1698be712e4SBarry Smith         n2->next  = a_n;
1708be712e4SBarry Smith         a_n->next = NULL;
1718be712e4SBarry Smith         break;
1728be712e4SBarry Smith       }
1738be712e4SBarry Smith       n2 = n2->next;
1748be712e4SBarry Smith     } while (n2);
1758be712e4SBarry Smith     PetscCheck(n2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "n2 should be non-null");
1768be712e4SBarry Smith   }
1778be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
1788be712e4SBarry Smith }
1798be712e4SBarry Smith 
1808be712e4SBarry Smith /* PetscCDRemoveNextNode: a_last->next, this exposes single linked list structure to API (not used)
1818be712e4SBarry Smith  */
1828be712e4SBarry Smith PetscErrorCode PetscCDRemoveNextNode(PetscCoarsenData *ail, PetscInt a_idx, PetscCDIntNd *a_last)
1838be712e4SBarry Smith {
1848be712e4SBarry Smith   PetscCDIntNd *del;
1858be712e4SBarry Smith 
1868be712e4SBarry Smith   PetscFunctionBegin;
1878be712e4SBarry Smith   PetscCheck(a_idx < ail->size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Index %" PetscInt_FMT " out of range.", a_idx);
1888be712e4SBarry Smith   PetscCheck(a_last->next, PETSC_COMM_SELF, PETSC_ERR_PLIB, "a_last should have a next");
1898be712e4SBarry Smith   del          = a_last->next;
1908be712e4SBarry Smith   a_last->next = del->next;
1918be712e4SBarry Smith   /* del->next = NULL; -- this still used in a iterator so keep it intact -- need to fix this with a double linked list */
1928be712e4SBarry Smith   /* could reuse n2 but PetscCDAppendNode sometimes uses it */
1938be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
1948be712e4SBarry Smith }
1958be712e4SBarry Smith 
1968be712e4SBarry Smith /* PetscCDPrint
1978be712e4SBarry Smith  */
1988be712e4SBarry Smith PetscErrorCode PetscCDPrint(const PetscCoarsenData *ail, PetscInt my0, MPI_Comm comm)
1998be712e4SBarry Smith {
2008be712e4SBarry Smith   PetscCDIntNd *n, *n2;
2018be712e4SBarry Smith   PetscInt      ii;
2028be712e4SBarry Smith 
2038be712e4SBarry Smith   PetscFunctionBegin;
2048be712e4SBarry Smith   for (ii = 0; ii < ail->size; ii++) {
2058be712e4SBarry Smith     n2 = n = ail->array[ii];
2068be712e4SBarry Smith     if (n) PetscCall(PetscSynchronizedPrintf(comm, "list %" PetscInt_FMT ":", ii + my0));
2078be712e4SBarry Smith     while (n) {
2088be712e4SBarry Smith       PetscCall(PetscSynchronizedPrintf(comm, " %" PetscInt_FMT, n->gid));
2098be712e4SBarry Smith       n = n->next;
2108be712e4SBarry Smith     }
2118be712e4SBarry Smith     if (n2) PetscCall(PetscSynchronizedPrintf(comm, "\n"));
2128be712e4SBarry Smith   }
2138be712e4SBarry Smith   PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
2148be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
2158be712e4SBarry Smith }
2168be712e4SBarry Smith 
2178be712e4SBarry Smith /* PetscCDMoveAppend - take list in a_srcidx and appends to destidx
2188be712e4SBarry Smith  */
2198be712e4SBarry Smith PetscErrorCode PetscCDMoveAppend(PetscCoarsenData *ail, PetscInt a_destidx, PetscInt a_srcidx)
2208be712e4SBarry Smith {
2218be712e4SBarry Smith   PetscCDIntNd *n;
2228be712e4SBarry Smith 
2238be712e4SBarry Smith   PetscFunctionBegin;
2248be712e4SBarry Smith   PetscCheck(a_srcidx < ail->size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Index %" PetscInt_FMT " out of range.", a_srcidx);
2258be712e4SBarry Smith   PetscCheck(a_destidx < ail->size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Index %" PetscInt_FMT " out of range.", a_destidx);
2268be712e4SBarry Smith   PetscCheck(a_destidx != a_srcidx, PETSC_COMM_SELF, PETSC_ERR_PLIB, "a_destidx==a_srcidx %" PetscInt_FMT ".", a_destidx);
2278be712e4SBarry Smith   n = ail->array[a_destidx];
2288be712e4SBarry Smith   if (!n) ail->array[a_destidx] = ail->array[a_srcidx];
2298be712e4SBarry Smith   else {
2308be712e4SBarry Smith     do {
2318be712e4SBarry Smith       if (!n->next) {
2328be712e4SBarry Smith         n->next = ail->array[a_srcidx]; // append
2338be712e4SBarry Smith         break;
2348be712e4SBarry Smith       }
2358be712e4SBarry Smith       n = n->next;
2368be712e4SBarry Smith     } while (1);
2378be712e4SBarry Smith   }
2388be712e4SBarry Smith   ail->array[a_srcidx] = NULL; // empty
2398be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
2408be712e4SBarry Smith }
2418be712e4SBarry Smith 
2428be712e4SBarry Smith /* PetscCDRemoveAllAt - empty one list and move data to cache
2438be712e4SBarry Smith  */
2448be712e4SBarry Smith PetscErrorCode PetscCDRemoveAllAt(PetscCoarsenData *ail, PetscInt a_idx)
2458be712e4SBarry Smith {
2468be712e4SBarry Smith   PetscCDIntNd *rem, *n1;
2478be712e4SBarry Smith 
2488be712e4SBarry Smith   PetscFunctionBegin;
2498be712e4SBarry Smith   PetscCheck(a_idx < ail->size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Index %" PetscInt_FMT " out of range.", a_idx);
2508be712e4SBarry Smith   rem               = ail->array[a_idx];
2518be712e4SBarry Smith   ail->array[a_idx] = NULL;
2528be712e4SBarry Smith   if (!(n1 = ail->extra_nodes)) ail->extra_nodes = rem;
2538be712e4SBarry Smith   else {
2548be712e4SBarry Smith     while (n1->next) n1 = n1->next;
2558be712e4SBarry Smith     n1->next = rem;
2568be712e4SBarry Smith   }
2578be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
2588be712e4SBarry Smith }
2598be712e4SBarry Smith 
2608be712e4SBarry Smith /* PetscCDCountAt
2618be712e4SBarry Smith  */
2628be712e4SBarry Smith PetscErrorCode PetscCDCountAt(const PetscCoarsenData *ail, PetscInt a_idx, PetscInt *a_sz)
2638be712e4SBarry Smith {
2648be712e4SBarry Smith   PetscCDIntNd *n1;
2658be712e4SBarry Smith   PetscInt      sz = 0;
2668be712e4SBarry Smith 
2678be712e4SBarry Smith   PetscFunctionBegin;
2688be712e4SBarry Smith   PetscCheck(a_idx < ail->size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Index %" PetscInt_FMT " out of range.", a_idx);
2698be712e4SBarry Smith   n1 = ail->array[a_idx];
2708be712e4SBarry Smith   while (n1) {
2718be712e4SBarry Smith     n1 = n1->next;
2728be712e4SBarry Smith     sz++;
2738be712e4SBarry Smith   }
2748be712e4SBarry Smith   *a_sz = sz;
2758be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
2768be712e4SBarry Smith }
2778be712e4SBarry Smith 
2788be712e4SBarry Smith /* PetscCDSize
2798be712e4SBarry Smith  */
2808be712e4SBarry Smith PetscErrorCode PetscCDCount(const PetscCoarsenData *ail, PetscInt *a_sz)
2818be712e4SBarry Smith {
2828be712e4SBarry Smith   PetscInt sz = 0;
2838be712e4SBarry Smith 
2848be712e4SBarry Smith   PetscFunctionBegin;
2858be712e4SBarry Smith   for (int ii = 0; ii < ail->size; ii++) {
2868be712e4SBarry Smith     PetscCDIntNd *n1 = ail->array[ii];
2878be712e4SBarry Smith     while (n1) {
2888be712e4SBarry Smith       n1 = n1->next;
2898be712e4SBarry Smith       sz++;
2908be712e4SBarry Smith     }
2918be712e4SBarry Smith   }
2928be712e4SBarry Smith   *a_sz = sz;
2938be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
2948be712e4SBarry Smith }
2958be712e4SBarry Smith 
2968be712e4SBarry Smith /* PetscCDIsEmptyAt - Is the list empty? (not used)
2978be712e4SBarry Smith  */
2988be712e4SBarry Smith PetscErrorCode PetscCDIsEmptyAt(const PetscCoarsenData *ail, PetscInt a_idx, PetscBool *a_e)
2998be712e4SBarry Smith {
3008be712e4SBarry Smith   PetscFunctionBegin;
3018be712e4SBarry Smith   PetscCheck(a_idx < ail->size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Index %" PetscInt_FMT " out of range.", a_idx);
3028be712e4SBarry Smith   *a_e = (PetscBool)(ail->array[a_idx] == NULL);
3038be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
3048be712e4SBarry Smith }
3058be712e4SBarry Smith 
3068be712e4SBarry Smith /* PetscCDGetNonemptyIS - used for C-F methods
3078be712e4SBarry Smith  */
3088be712e4SBarry Smith PetscErrorCode PetscCDGetNonemptyIS(PetscCoarsenData *ail, IS *a_mis)
3098be712e4SBarry Smith {
3108be712e4SBarry Smith   PetscCDIntNd *n;
3118be712e4SBarry Smith   PetscInt      ii, kk;
3128be712e4SBarry Smith   PetscInt     *permute;
3138be712e4SBarry Smith 
3148be712e4SBarry Smith   PetscFunctionBegin;
3158be712e4SBarry Smith   for (ii = kk = 0; ii < ail->size; ii++) {
3168be712e4SBarry Smith     n = ail->array[ii];
3178be712e4SBarry Smith     if (n) kk++;
3188be712e4SBarry Smith   }
3198be712e4SBarry Smith   PetscCall(PetscMalloc1(kk, &permute));
3208be712e4SBarry Smith   for (ii = kk = 0; ii < ail->size; ii++) {
3218be712e4SBarry Smith     n = ail->array[ii];
3228be712e4SBarry Smith     if (n) permute[kk++] = ii;
3238be712e4SBarry Smith   }
3248be712e4SBarry Smith   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, kk, permute, PETSC_OWN_POINTER, a_mis));
3258be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
3268be712e4SBarry Smith }
3278be712e4SBarry Smith 
3288be712e4SBarry Smith /* PetscCDGetMat
3298be712e4SBarry Smith  */
3308be712e4SBarry Smith PetscErrorCode PetscCDGetMat(PetscCoarsenData *ail, Mat *a_mat)
3318be712e4SBarry Smith {
3328be712e4SBarry Smith   PetscFunctionBegin;
3338be712e4SBarry Smith   *a_mat = ail->mat;
3348be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
3358be712e4SBarry Smith }
3368be712e4SBarry Smith 
3378be712e4SBarry Smith /* PetscCDSetMat
3388be712e4SBarry Smith  */
3398be712e4SBarry Smith PetscErrorCode PetscCDSetMat(PetscCoarsenData *ail, Mat a_mat)
3408be712e4SBarry Smith {
3418be712e4SBarry Smith   PetscFunctionBegin;
3428be712e4SBarry Smith   if (ail->mat) {
3438be712e4SBarry Smith     PetscCall(MatDestroy(&ail->mat)); //should not happen
3448be712e4SBarry Smith   }
3458be712e4SBarry Smith   ail->mat = a_mat;
3468be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
3478be712e4SBarry Smith }
3488be712e4SBarry Smith 
3498926f930SMark Adams /* PetscCDClearMat
3508926f930SMark Adams  */
3518926f930SMark Adams PetscErrorCode PetscCDClearMat(PetscCoarsenData *ail)
3528926f930SMark Adams {
3538926f930SMark Adams   PetscFunctionBegin;
3548926f930SMark Adams   ail->mat = NULL;
3558926f930SMark Adams   PetscFunctionReturn(PETSC_SUCCESS);
3568926f930SMark Adams }
3578926f930SMark Adams 
3588be712e4SBarry Smith /* PetscCDGetASMBlocks - get IS of aggregates for ASM smoothers
3598be712e4SBarry Smith  */
3608be712e4SBarry Smith PetscErrorCode PetscCDGetASMBlocks(const PetscCoarsenData *ail, const PetscInt a_bs, PetscInt *a_sz, IS **a_local_is)
3618be712e4SBarry Smith {
3628be712e4SBarry Smith   PetscCDIntNd *n;
3638be712e4SBarry Smith   PetscInt      lsz, ii, kk, *idxs, jj, gid;
3648be712e4SBarry Smith   IS           *is_loc = NULL;
3658be712e4SBarry Smith 
3668be712e4SBarry Smith   PetscFunctionBegin;
3678be712e4SBarry Smith   for (ii = kk = 0; ii < ail->size; ii++) {
3688be712e4SBarry Smith     if (ail->array[ii]) kk++;
3698be712e4SBarry Smith   }
3708be712e4SBarry Smith   *a_sz = kk;
3718be712e4SBarry Smith   PetscCall(PetscMalloc1(kk, &is_loc));
3728be712e4SBarry Smith   for (ii = kk = 0; ii < ail->size; ii++) {
3738be712e4SBarry Smith     for (lsz = 0, n = ail->array[ii]; n; lsz++, n = n->next) /* void */
3748be712e4SBarry Smith       ;
3758be712e4SBarry Smith     if (lsz) {
3768be712e4SBarry Smith       PetscCall(PetscMalloc1(a_bs * lsz, &idxs));
3778be712e4SBarry Smith       for (lsz = 0, n = ail->array[ii]; n; n = n->next) {
3788be712e4SBarry Smith         PetscCall(PetscCDIntNdGetID(n, &gid));
3798be712e4SBarry Smith         for (jj = 0; jj < a_bs; lsz++, jj++) idxs[lsz] = a_bs * gid + jj;
3808be712e4SBarry Smith       }
3818be712e4SBarry Smith       PetscCall(ISCreateGeneral(PETSC_COMM_SELF, lsz, idxs, PETSC_OWN_POINTER, &is_loc[kk++]));
3828be712e4SBarry Smith     }
3838be712e4SBarry Smith   }
3848be712e4SBarry Smith   PetscCheck(*a_sz == kk, PETSC_COMM_SELF, PETSC_ERR_PLIB, "*a_sz %" PetscInt_FMT " != kk %" PetscInt_FMT, *a_sz, kk);
3858be712e4SBarry Smith   *a_local_is = is_loc; /* out */
3868be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
3878be712e4SBarry Smith }
3888be712e4SBarry Smith 
3898be712e4SBarry Smith /* edge for priority queue */
3908be712e4SBarry Smith typedef struct edge_tag {
3918be712e4SBarry Smith   PetscReal weight;
3928be712e4SBarry Smith   PetscInt  lid0, gid1, ghost1_idx;
3938be712e4SBarry Smith } Edge;
3948be712e4SBarry Smith 
3958be712e4SBarry Smith #define MY_MEPS (PETSC_MACHINE_EPSILON * 100)
3968be712e4SBarry Smith static int gamg_hem_compare(const void *a, const void *b)
3978be712e4SBarry Smith {
3988be712e4SBarry Smith   PetscReal va = ((Edge *)a)->weight, vb = ((Edge *)b)->weight;
3998be712e4SBarry Smith   return (va <= vb - MY_MEPS) ? 1 : (va > vb + MY_MEPS) ? -1 : 0; /* 0 for equal */
4008be712e4SBarry Smith }
4018be712e4SBarry Smith 
4028be712e4SBarry Smith /*
4038be712e4SBarry Smith   MatCoarsenApply_HEM_private - parallel heavy edge matching
4048be712e4SBarry Smith 
4058be712e4SBarry Smith   Input Parameter:
4068be712e4SBarry Smith    . a_Gmat - global matrix of the graph
4078be712e4SBarry Smith    . n_iter - number of matching iterations
4088be712e4SBarry Smith    . threshold - threshold for filtering graphs
4098be712e4SBarry Smith 
4108be712e4SBarry Smith   Output Parameter:
4118be712e4SBarry Smith    . a_locals_llist - array of list of local nodes rooted at local node
4128be712e4SBarry Smith */
4138be712e4SBarry Smith static PetscErrorCode MatCoarsenApply_HEM_private(Mat a_Gmat, const PetscInt n_iter, const PetscReal threshold, PetscCoarsenData **a_locals_llist)
4148be712e4SBarry Smith {
4158be712e4SBarry Smith #define REQ_BF_SIZE 100
4168be712e4SBarry Smith   PetscBool         isMPI;
4178be712e4SBarry Smith   MPI_Comm          comm;
4188be712e4SBarry Smith   PetscInt          ix, *ii, *aj, Iend, my0, ncomm_procs, bc_agg = -1, *rbuff = NULL, rbuff_sz = 0;
4198be712e4SBarry Smith   PetscMPIInt       rank, size, comm_procs[REQ_BF_SIZE], *lid_max_pe;
4208be712e4SBarry Smith   const PetscInt    nloc = a_Gmat->rmap->n, request_size = PetscCeilReal((PetscReal)sizeof(MPI_Request) / (PetscReal)sizeof(PetscInt));
4218be712e4SBarry Smith   PetscInt         *lid_cprowID;
4228be712e4SBarry Smith   PetscBool        *lid_matched;
4238be712e4SBarry Smith   Mat_SeqAIJ       *matA, *matB = NULL;
4248be712e4SBarry Smith   Mat_MPIAIJ       *mpimat     = NULL;
4258be712e4SBarry Smith   PetscScalar       one        = 1.;
4268be712e4SBarry Smith   PetscCoarsenData *agg_llists = NULL, *ghost_deleted_list = NULL, *bc_list = NULL;
4278be712e4SBarry Smith   Mat               cMat, tMat, P;
4288be712e4SBarry Smith   MatScalar        *ap;
4298be712e4SBarry Smith   IS                info_is;
4308be712e4SBarry Smith 
4318be712e4SBarry Smith   PetscFunctionBegin;
4328be712e4SBarry Smith   PetscCall(PetscObjectGetComm((PetscObject)a_Gmat, &comm));
4338be712e4SBarry Smith   PetscCallMPI(MPI_Comm_rank(comm, &rank));
4348be712e4SBarry Smith   PetscCallMPI(MPI_Comm_size(comm, &size));
4358be712e4SBarry Smith   PetscCall(MatGetOwnershipRange(a_Gmat, &my0, &Iend));
4368be712e4SBarry Smith   PetscCall(ISCreate(comm, &info_is));
4378926f930SMark Adams   PetscCall(PetscInfo(info_is, "Start %" PetscInt_FMT " iterations of HEM.\n", n_iter));
4388be712e4SBarry Smith 
4398be712e4SBarry Smith   PetscCall(PetscMalloc1(nloc, &lid_matched));
4408be712e4SBarry Smith   PetscCall(PetscMalloc1(nloc, &lid_cprowID));
4418be712e4SBarry Smith   PetscCall(PetscMalloc1(nloc, &lid_max_pe));
4428be712e4SBarry Smith 
4438be712e4SBarry Smith   PetscCall(PetscCDCreate(nloc, &agg_llists));
4448be712e4SBarry Smith   PetscCall(PetscCDSetChunkSize(agg_llists, nloc + 1));
4458be712e4SBarry Smith   *a_locals_llist = agg_llists;
4468be712e4SBarry Smith   /* add self to all lists */
4478be712e4SBarry Smith   for (int kk = 0; kk < nloc; kk++) PetscCall(PetscCDAppendID(agg_llists, kk, my0 + kk));
4488be712e4SBarry Smith   /* make a copy of the graph, this gets destroyed in iterates */
4498be712e4SBarry Smith   PetscCall(MatDuplicate(a_Gmat, MAT_COPY_VALUES, &cMat));
4508be712e4SBarry Smith   PetscCall(MatConvert(cMat, MATAIJ, MAT_INPLACE_MATRIX, &cMat));
4518be712e4SBarry Smith   isMPI = (PetscBool)(size > 1);
4528be712e4SBarry Smith   if (isMPI) {
4538be712e4SBarry Smith     /* list of deleted ghosts, should compress this */
4548be712e4SBarry Smith     PetscCall(PetscCDCreate(size, &ghost_deleted_list));
4558be712e4SBarry Smith     PetscCall(PetscCDSetChunkSize(ghost_deleted_list, 100));
4568be712e4SBarry Smith   }
4578be712e4SBarry Smith   for (int iter = 0; iter < n_iter; iter++) {
4588be712e4SBarry Smith     const PetscScalar *lghost_max_ew, *lid_max_ew;
4598be712e4SBarry Smith     PetscBool         *lghost_matched;
4608be712e4SBarry Smith     PetscMPIInt       *lghost_pe, *lghost_max_pe;
4618be712e4SBarry Smith     Vec                locMaxEdge, ghostMaxEdge, ghostMaxPE, locMaxPE;
4628be712e4SBarry Smith     PetscInt          *lghost_gid, nEdges, nEdges0, num_ghosts = 0;
4638be712e4SBarry Smith     Edge              *Edges;
464452acf8bSMark Adams     const int          n_sub_its = 1000; // in case of a bug, stop at some point
4658be712e4SBarry Smith     /* get submatrices of cMat */
4668be712e4SBarry Smith     for (int 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 */
4848be712e4SBarry Smith     for (int 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;
4978be712e4SBarry Smith       PetscCall(MatCreateVecs(cMat, &vec, NULL));
4988be712e4SBarry Smith       PetscCall(VecGetLocalSize(mpimat->lvec, &num_ghosts));
4998be712e4SBarry Smith       /* lghost_matched */
5008be712e4SBarry Smith       for (PetscInt kk = 0, gid = my0; kk < nloc; kk++, gid++) {
5018be712e4SBarry Smith         PetscScalar vval = lid_matched[kk] ? 1.0 : 0.0;
5028be712e4SBarry Smith         PetscCall(VecSetValues(vec, 1, &gid, &vval, INSERT_VALUES));
5038be712e4SBarry Smith       }
5048be712e4SBarry Smith       PetscCall(VecAssemblyBegin(vec));
5058be712e4SBarry Smith       PetscCall(VecAssemblyEnd(vec));
5068be712e4SBarry Smith       PetscCall(VecScatterBegin(mpimat->Mvctx, vec, mpimat->lvec, INSERT_VALUES, SCATTER_FORWARD));
5078be712e4SBarry Smith       PetscCall(VecScatterEnd(mpimat->Mvctx, vec, mpimat->lvec, INSERT_VALUES, SCATTER_FORWARD));
5088be712e4SBarry Smith       PetscCall(VecGetArrayRead(mpimat->lvec, &buf)); /* get proc ID in 'buf' */
5098be712e4SBarry Smith       PetscCall(PetscMalloc1(num_ghosts, &lghost_matched));
5108be712e4SBarry Smith       for (int kk = 0; kk < num_ghosts; kk++) {
5118be712e4SBarry Smith         lghost_matched[kk] = (PetscBool)(PetscRealPart(buf[kk]) != 0); // the proc of the ghost for now
5128be712e4SBarry Smith       }
5138be712e4SBarry Smith       PetscCall(VecRestoreArrayRead(mpimat->lvec, &buf));
5148be712e4SBarry Smith       /* lghost_pe */
5158be712e4SBarry Smith       vval = (PetscScalar)(rank);
5168be712e4SBarry Smith       for (PetscInt kk = 0, gid = my0; kk < nloc; kk++, gid++) PetscCall(VecSetValues(vec, 1, &gid, &vval, INSERT_VALUES)); /* set with GID */
5178be712e4SBarry Smith       PetscCall(VecAssemblyBegin(vec));
5188be712e4SBarry Smith       PetscCall(VecAssemblyEnd(vec));
5198be712e4SBarry Smith       PetscCall(VecScatterBegin(mpimat->Mvctx, vec, mpimat->lvec, INSERT_VALUES, SCATTER_FORWARD));
5208be712e4SBarry Smith       PetscCall(VecScatterEnd(mpimat->Mvctx, vec, mpimat->lvec, INSERT_VALUES, SCATTER_FORWARD));
5218be712e4SBarry Smith       PetscCall(VecGetArrayRead(mpimat->lvec, &buf)); /* get proc ID in 'buf' */
5228be712e4SBarry Smith       PetscCall(PetscMalloc1(num_ghosts, &lghost_pe));
5238be712e4SBarry Smith       for (int 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 */
5268be712e4SBarry Smith       for (PetscInt kk = 0, gid = my0; kk < nloc; kk++, gid++) {
5278be712e4SBarry Smith         vval = (PetscScalar)(gid);
5288be712e4SBarry Smith         PetscCall(VecSetValues(vec, 1, &gid, &vval, INSERT_VALUES)); /* set with GID */
5298be712e4SBarry Smith       }
5308be712e4SBarry Smith       PetscCall(VecAssemblyBegin(vec));
5318be712e4SBarry Smith       PetscCall(VecAssemblyEnd(vec));
5328be712e4SBarry Smith       PetscCall(VecScatterBegin(mpimat->Mvctx, vec, mpimat->lvec, INSERT_VALUES, SCATTER_FORWARD));
5338be712e4SBarry Smith       PetscCall(VecScatterEnd(mpimat->Mvctx, vec, mpimat->lvec, INSERT_VALUES, SCATTER_FORWARD));
5348be712e4SBarry Smith       PetscCall(VecDestroy(&vec));
5358be712e4SBarry Smith       PetscCall(VecGetArrayRead(mpimat->lvec, &buf)); /* get proc ID in 'lghost_gid' */
5368be712e4SBarry Smith       PetscCall(PetscMalloc1(num_ghosts, &lghost_gid));
5378be712e4SBarry Smith       for (int 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)
5418be712e4SBarry Smith     for (int 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;
5448be712e4SBarry Smith       for (int k = 0; k < ncomm_procs && idx == -1; k++)
5458be712e4SBarry Smith         if (comm_procs[k] == proc) idx = k;
5468be712e4SBarry Smith       if (idx == -1) { comm_procs[ncomm_procs++] = proc; }
5478be712e4SBarry Smith       PetscCheck(ncomm_procs != REQ_BF_SIZE, PETSC_COMM_SELF, PETSC_ERR_SUP, "Receive request array too small: %d", (int)ncomm_procs);
5488be712e4SBarry Smith     }
5498be712e4SBarry Smith     /* count edges, compute initial 'locMaxEdge', 'locMaxPE' */
5508be712e4SBarry Smith     nEdges0 = 0;
5518be712e4SBarry Smith     for (PetscInt kk = 0, gid = my0; kk < nloc; kk++, gid++) {
5528be712e4SBarry Smith       PetscReal   max_e = 0., tt;
5538be712e4SBarry Smith       PetscScalar vval;
5548be712e4SBarry Smith       PetscInt    lid = kk, max_pe = rank, pe, n;
5558be712e4SBarry Smith       ii = matA->i;
5568be712e4SBarry Smith       n  = ii[lid + 1] - ii[lid];
5578e3a54c0SPierre Jolivet       aj = PetscSafePointerPlusOffset(matA->j, ii[lid]);
5588e3a54c0SPierre Jolivet       ap = PetscSafePointerPlusOffset(matA->a, ii[lid]);
5598be712e4SBarry Smith       for (int jj = 0; jj < n; jj++) {
5608be712e4SBarry Smith         PetscInt lidj = aj[jj];
5618be712e4SBarry Smith         if ((tt = PetscRealPart(ap[jj])) > threshold && lidj != lid) {
5628be712e4SBarry Smith           if (tt > max_e) max_e = tt;
5638be712e4SBarry Smith           if (lidj > lid) nEdges0++;
5648be712e4SBarry Smith         }
5658be712e4SBarry Smith       }
5668be712e4SBarry Smith       if ((ix = lid_cprowID[lid]) != -1) { /* if I have any ghost neighbors */
5678be712e4SBarry Smith         ii = matB->compressedrow.i;
5688be712e4SBarry Smith         n  = ii[ix + 1] - ii[ix];
5698be712e4SBarry Smith         ap = matB->a + ii[ix];
5708be712e4SBarry Smith         aj = matB->j + ii[ix];
5718be712e4SBarry Smith         for (int jj = 0; jj < n; jj++) {
5728be712e4SBarry Smith           if ((tt = PetscRealPart(ap[jj])) > threshold) {
5738be712e4SBarry Smith             if (tt > max_e) max_e = tt;
5748be712e4SBarry Smith             nEdges0++;
5758be712e4SBarry Smith             if ((pe = lghost_pe[aj[jj]]) > max_pe) max_pe = pe;
5768be712e4SBarry Smith           }
5778be712e4SBarry Smith         }
5788be712e4SBarry Smith       }
5798be712e4SBarry Smith       vval = max_e;
5808be712e4SBarry Smith       PetscCall(VecSetValues(locMaxEdge, 1, &gid, &vval, INSERT_VALUES));
5818be712e4SBarry Smith       vval = (PetscScalar)max_pe;
5828be712e4SBarry Smith       PetscCall(VecSetValues(locMaxPE, 1, &gid, &vval, INSERT_VALUES));
5838be712e4SBarry Smith       if (iter == 0 && max_e <= MY_MEPS) { // add BCs to fake aggregate
5848be712e4SBarry Smith         lid_matched[lid] = PETSC_TRUE;
5858be712e4SBarry Smith         if (bc_agg == -1) {
5868be712e4SBarry Smith           bc_agg = lid;
5878be712e4SBarry Smith           PetscCall(PetscCDCreate(1, &bc_list));
5888be712e4SBarry Smith         }
5898be712e4SBarry Smith         PetscCall(PetscCDRemoveAllAt(agg_llists, lid));
5908be712e4SBarry Smith         PetscCall(PetscCDAppendID(bc_list, 0, my0 + lid));
5918be712e4SBarry Smith       }
5928be712e4SBarry Smith     }
5938be712e4SBarry Smith     PetscCall(VecAssemblyBegin(locMaxEdge));
5948be712e4SBarry Smith     PetscCall(VecAssemblyEnd(locMaxEdge));
5958be712e4SBarry Smith     PetscCall(VecAssemblyBegin(locMaxPE));
5968be712e4SBarry Smith     PetscCall(VecAssemblyEnd(locMaxPE));
5978be712e4SBarry Smith     /* make 'ghostMaxEdge_max_ew', 'lghost_max_pe' */
5988be712e4SBarry Smith     if (mpimat) {
5998be712e4SBarry Smith       const PetscScalar *buf;
6008be712e4SBarry Smith       PetscCall(VecDuplicate(mpimat->lvec, &ghostMaxEdge));
6018be712e4SBarry Smith       PetscCall(VecScatterBegin(mpimat->Mvctx, locMaxEdge, ghostMaxEdge, INSERT_VALUES, SCATTER_FORWARD));
6028be712e4SBarry Smith       PetscCall(VecScatterEnd(mpimat->Mvctx, locMaxEdge, ghostMaxEdge, INSERT_VALUES, SCATTER_FORWARD));
6038be712e4SBarry Smith 
6048be712e4SBarry Smith       PetscCall(VecDuplicate(mpimat->lvec, &ghostMaxPE));
6058be712e4SBarry Smith       PetscCall(VecScatterBegin(mpimat->Mvctx, locMaxPE, ghostMaxPE, INSERT_VALUES, SCATTER_FORWARD));
6068be712e4SBarry Smith       PetscCall(VecScatterEnd(mpimat->Mvctx, locMaxPE, ghostMaxPE, INSERT_VALUES, SCATTER_FORWARD));
6078be712e4SBarry Smith       PetscCall(VecGetArrayRead(ghostMaxPE, &buf));
6088be712e4SBarry Smith       PetscCall(PetscMalloc1(num_ghosts, &lghost_max_pe));
6098be712e4SBarry Smith       for (int kk = 0; kk < num_ghosts; kk++) lghost_max_pe[kk] = (PetscMPIInt)PetscRealPart(buf[kk]); // the MAX proc of the ghost now
6108be712e4SBarry Smith       PetscCall(VecRestoreArrayRead(ghostMaxPE, &buf));
6118be712e4SBarry Smith     }
6128be712e4SBarry Smith     { // make lid_max_pe
6138be712e4SBarry Smith       const PetscScalar *buf;
6148be712e4SBarry Smith       PetscCall(VecGetArrayRead(locMaxPE, &buf));
6158be712e4SBarry Smith       for (int kk = 0; kk < nloc; kk++) lid_max_pe[kk] = (PetscMPIInt)PetscRealPart(buf[kk]); // the MAX proc of the ghost now
6168be712e4SBarry Smith       PetscCall(VecRestoreArrayRead(locMaxPE, &buf));
6178be712e4SBarry Smith     }
6188be712e4SBarry Smith     /* setup sorted list of edges, and make 'Edges' */
6198be712e4SBarry Smith     PetscCall(PetscMalloc1(nEdges0, &Edges));
6208be712e4SBarry Smith     nEdges = 0;
6218be712e4SBarry Smith     for (int kk = 0, n; kk < nloc; kk++) {
6228be712e4SBarry Smith       const PetscInt lid = kk;
6238be712e4SBarry Smith       PetscReal      tt;
6248be712e4SBarry Smith       ii = matA->i;
6258be712e4SBarry Smith       n  = ii[lid + 1] - ii[lid];
6268e3a54c0SPierre Jolivet       aj = PetscSafePointerPlusOffset(matA->j, ii[lid]);
6278e3a54c0SPierre Jolivet       ap = PetscSafePointerPlusOffset(matA->a, ii[lid]);
6288be712e4SBarry Smith       for (int jj = 0; jj < n; jj++) {
6298be712e4SBarry Smith         PetscInt lidj = aj[jj];
6308be712e4SBarry Smith         if ((tt = PetscRealPart(ap[jj])) > threshold && lidj != lid) {
6318be712e4SBarry Smith           if (lidj > lid) {
6328be712e4SBarry Smith             Edges[nEdges].lid0       = lid;
6338be712e4SBarry Smith             Edges[nEdges].gid1       = lidj + my0;
6348be712e4SBarry Smith             Edges[nEdges].ghost1_idx = -1;
6358be712e4SBarry Smith             Edges[nEdges].weight     = tt;
6368be712e4SBarry Smith             nEdges++;
6378be712e4SBarry Smith           }
6388be712e4SBarry Smith         }
6398be712e4SBarry Smith       }
6408be712e4SBarry Smith       if ((ix = lid_cprowID[lid]) != -1) { /* if I have any ghost neighbor */
6418be712e4SBarry Smith         ii = matB->compressedrow.i;
6428be712e4SBarry Smith         n  = ii[ix + 1] - ii[ix];
6438be712e4SBarry Smith         ap = matB->a + ii[ix];
6448be712e4SBarry Smith         aj = matB->j + ii[ix];
6458be712e4SBarry Smith         for (int jj = 0; jj < n; jj++) {
6468be712e4SBarry Smith           if ((tt = PetscRealPart(ap[jj])) > threshold) {
6478be712e4SBarry Smith             Edges[nEdges].lid0       = lid;
6488be712e4SBarry Smith             Edges[nEdges].gid1       = lghost_gid[aj[jj]];
6498be712e4SBarry Smith             Edges[nEdges].ghost1_idx = aj[jj];
6508be712e4SBarry Smith             Edges[nEdges].weight     = tt;
6518be712e4SBarry Smith             nEdges++;
6528be712e4SBarry Smith           }
6538be712e4SBarry Smith         }
6548be712e4SBarry Smith       }
6558be712e4SBarry Smith     }
6568be712e4SBarry Smith     PetscCheck(nEdges == nEdges0, PETSC_COMM_SELF, PETSC_ERR_SUP, "nEdges != nEdges0: %d %d", (int)nEdges0, (int)nEdges);
657810441c8SPierre Jolivet     if (Edges) qsort(Edges, nEdges, sizeof(Edge), gamg_hem_compare);
6588be712e4SBarry Smith 
6598926f930SMark Adams     PetscCall(PetscInfo(info_is, "[%d] HEM iteration %d with %d edges\n", rank, iter, (int)nEdges));
6608be712e4SBarry Smith 
6618be712e4SBarry Smith     /* projection matrix */
6628be712e4SBarry Smith     PetscCall(MatCreate(comm, &P));
6638be712e4SBarry Smith     PetscCall(MatSetType(P, MATAIJ));
6648be712e4SBarry Smith     PetscCall(MatSetSizes(P, nloc, nloc, PETSC_DETERMINE, PETSC_DETERMINE));
6658be712e4SBarry Smith     PetscCall(MatMPIAIJSetPreallocation(P, 1, NULL, 1, NULL));
6668be712e4SBarry Smith     PetscCall(MatSeqAIJSetPreallocation(P, 1, NULL));
6678be712e4SBarry Smith     PetscCall(MatSetUp(P));
6688be712e4SBarry Smith     /* process - communicate - process */
669452acf8bSMark Adams     for (int sub_it = 0, old_num_edge = 0; /* sub_it < n_sub_its */; /* sub_it++ */) {
6708be712e4SBarry Smith       PetscInt    nactive_edges = 0, n_act_n[3], gn_act_n[3];
6718be712e4SBarry Smith       PetscMPIInt tag1, tag2;
6728be712e4SBarry Smith       PetscCall(VecGetArrayRead(locMaxEdge, &lid_max_ew));
6738be712e4SBarry Smith       if (isMPI) {
6748be712e4SBarry Smith         PetscCall(VecGetArrayRead(ghostMaxEdge, &lghost_max_ew));
6758be712e4SBarry Smith         PetscCall(PetscCommGetNewTag(comm, &tag1));
6768be712e4SBarry Smith         PetscCall(PetscCommGetNewTag(comm, &tag2));
6778be712e4SBarry Smith       }
6788be712e4SBarry Smith       for (int kk = 0; kk < nEdges; kk++) {
6798be712e4SBarry Smith         /* HEM */
6808be712e4SBarry Smith         const Edge    *e    = &Edges[kk];
6818be712e4SBarry Smith         const PetscInt lid0 = e->lid0, gid1 = e->gid1, ghost1_idx = e->ghost1_idx, gid0 = lid0 + my0, lid1 = gid1 - my0;
6828be712e4SBarry Smith         PetscBool      isOK = PETSC_TRUE, print = PETSC_FALSE;
683452acf8bSMark Adams         if (print) PetscCall(PetscSynchronizedPrintf(comm, "\t[%d] edge (%d %d), %d %d %d\n", rank, (int)gid0, (int)gid1, lid_matched[lid0], (ghost1_idx != -1 && lghost_matched[ghost1_idx]), (ghost1_idx == -1 && lid_matched[lid1])));
6848be712e4SBarry Smith         /* skip if either vertex is matched already */
6858be712e4SBarry Smith         if (lid_matched[lid0] || (ghost1_idx != -1 && lghost_matched[ghost1_idx]) || (ghost1_idx == -1 && lid_matched[lid1])) continue;
6868be712e4SBarry Smith 
6878be712e4SBarry Smith         nactive_edges++;
6888be712e4SBarry 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]));
689452acf8bSMark Adams         if (print) PetscCall(PetscSynchronizedPrintf(comm, "\t[%d] active edge (%d %d), diff0 = %10.4e\n", rank, (int)gid0, (int)gid1, (double)(PetscRealPart(lid_max_ew[lid0]) - (double)e->weight)));
6908be712e4SBarry Smith         // smaller edge, lid_max_ew get updated - e0
6918be712e4SBarry Smith         if (PetscRealPart(lid_max_ew[lid0]) > e->weight + MY_MEPS) {
6928be712e4SBarry Smith           if (print)
693452acf8bSMark Adams             PetscCall(PetscSynchronizedPrintf(comm, "\t\t[%d] 1) e0 SKIPPING small edge %20.14e edge (%d %d), diff = %10.4e to proc %d. max = %20.14e, w = %20.14e\n", rank, (double)e->weight, (int)gid0, (int)gid1, (double)(PetscRealPart(lid_max_ew[lid0]) - e->weight), ghost1_idx != -1 ? (int)lghost_pe[ghost1_idx] : rank, (double)PetscRealPart(lid_max_ew[lid0]),
6948be712e4SBarry Smith                                               (double)e->weight));
6958be712e4SBarry Smith           continue; // we are basically filter edges here
6968be712e4SBarry Smith         }
6978be712e4SBarry Smith         // e1 - local
6988be712e4SBarry Smith         if (ghost1_idx == -1) {
6998be712e4SBarry Smith           if (PetscRealPart(lid_max_ew[lid1]) > e->weight + MY_MEPS) {
7008be712e4SBarry Smith             if (print)
701452acf8bSMark Adams               PetscCall(PetscSynchronizedPrintf(comm, "\t\t%c[%d] 2) e1 SKIPPING small local edge %20.14e edge (%d %d), diff = %10.4e\n", ghost1_idx != -1 ? '\t' : ' ', rank, (double)e->weight, (int)gid0, (int)gid1, (double)(PetscRealPart(lid_max_ew[lid1]) - e->weight)));
7028be712e4SBarry Smith             continue; // we are basically filter edges here
7038be712e4SBarry Smith           }
7048be712e4SBarry Smith         } else { // e1 - ghost
7058be712e4SBarry Smith           /* see if edge might get matched on other proc */
7068be712e4SBarry Smith           PetscReal g_max_e1 = PetscRealPart(lghost_max_ew[ghost1_idx]);
7078be712e4SBarry Smith           if (print)
708452acf8bSMark Adams             PetscCall(PetscSynchronizedPrintf(comm, "\t\t\t[%d] CHECK GHOST e1, edge (%d %d), 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, (int)gid0, (int)gid1, (double)PetscRealPart(lid_max_ew[lid0]),
7098be712e4SBarry Smith                                               (double)e->weight, (double)(PetscRealPart(lghost_max_ew[ghost1_idx]) - e->weight), (int)lghost_pe[ghost1_idx], lid_max_pe[lid0], lghost_max_pe[ghost1_idx]));
7108be712e4SBarry Smith           if (g_max_e1 > e->weight + MY_MEPS) {
711452acf8bSMark Adams             /* 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, (int)gid0, (int)gid1, g_max_e1 - e->weight, (int)lghost_pe[ghost1_idx], lghost_max_pe[ghost1_idx], g_max_e1, e->weight )); */
7128be712e4SBarry Smith             continue;
713452acf8bSMark Adams           } else if (g_max_e1 >= e->weight - MY_MEPS && lghost_pe[ghost1_idx] > rank) { // is 'lghost_max_pe[ghost1_idx] > rank' needed?
7148be712e4SBarry Smith             /* check for max_ea == to this edge and larger processor that will deal with this */
7158be712e4SBarry Smith             if (print)
716452acf8bSMark Adams               PetscCall(PetscSynchronizedPrintf(comm, "\t\t\t[%d] ghost e1 SKIPPING EQUAL (%d %d), diff = %10.4e from larger proc %d with max pe %d. max = %20.14e, w = %20.14e\n", rank, (int)gid0, (int)gid1,
7178be712e4SBarry Smith                                                 (double)(PetscRealPart(lid_max_ew[lid0]) - (double)e->weight), (int)lghost_pe[ghost1_idx], (int)lghost_max_pe[ghost1_idx], (double)g_max_e1, (double)e->weight));
7188be712e4SBarry Smith             isOK = PETSC_FALSE; // this guy could delete me
7198be712e4SBarry Smith             continue;
7208be712e4SBarry Smith           } else {
721452acf8bSMark Adams             /* 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, (int)gid0, (int)gid1, g_max_e1 - e->weight, (int)lghost_pe[ghost1_idx], lghost_max_pe[ghost1_idx], g_max_e1, e->weight )); */
7228be712e4SBarry Smith           }
7238be712e4SBarry Smith         }
7248be712e4SBarry Smith         /* check ghost for v0 */
7258be712e4SBarry Smith         if (isOK) {
7268be712e4SBarry Smith           PetscReal max_e, ew;
7278be712e4SBarry Smith           if ((ix = lid_cprowID[lid0]) != -1) { /* if I have any ghost neighbors */
7288be712e4SBarry Smith             int n;
7298be712e4SBarry Smith             ii = matB->compressedrow.i;
7308be712e4SBarry Smith             n  = ii[ix + 1] - ii[ix];
7318be712e4SBarry Smith             ap = matB->a + ii[ix];
7328be712e4SBarry Smith             aj = matB->j + ii[ix];
7338be712e4SBarry Smith             for (int jj = 0; jj < n && isOK; jj++) {
7348be712e4SBarry Smith               PetscInt lidj = aj[jj];
7358be712e4SBarry Smith               if (lghost_matched[lidj]) continue;
7368be712e4SBarry Smith               ew = PetscRealPart(ap[jj]);
7378be712e4SBarry Smith               if (ew <= threshold) continue;
7388be712e4SBarry Smith               max_e = PetscRealPart(lghost_max_ew[lidj]);
7398be712e4SBarry Smith               /* check for max_e == to this edge and larger processor that will deal with this */
740452acf8bSMark Adams               if (ew >= PetscRealPart(lid_max_ew[lid0]) - MY_MEPS && lghost_max_pe[lidj] > rank) isOK = PETSC_FALSE;
7418926f930SMark Adams               PetscCheck(ew <= max_e + MY_MEPS, PETSC_COMM_SELF, PETSC_ERR_SUP, "edge weight %e > max %e. ncols = %d, gid0 = %d, gid1 = %d", (double)PetscRealPart(ew), (double)PetscRealPart(max_e), (int)n, (int)(lid0 + my0), (int)lghost_gid[lidj]);
7428be712e4SBarry Smith               if (print)
743452acf8bSMark Adams                 PetscCall(PetscSynchronizedPrintf(comm, "\t\t\t\t[%d] e0: looked at ghost adj (%d %d), 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, (int)gid0, (int)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));
7448be712e4SBarry Smith             }
745452acf8bSMark Adams             if (!isOK && print) PetscCall(PetscSynchronizedPrintf(comm, "\t\t[%d] skip edge (%d %d) from ghost inspection\n", rank, (int)gid0, (int)gid1));
7468be712e4SBarry Smith           }
7478be712e4SBarry Smith           /* check local v1 */
7488be712e4SBarry Smith           if (ghost1_idx == -1) {
7498be712e4SBarry Smith             if ((ix = lid_cprowID[lid1]) != -1) { /* if I have any ghost neighbors */
7508be712e4SBarry Smith               int n;
7518be712e4SBarry Smith               ii = matB->compressedrow.i;
7528be712e4SBarry Smith               n  = ii[ix + 1] - ii[ix];
7538be712e4SBarry Smith               ap = matB->a + ii[ix];
7548be712e4SBarry Smith               aj = matB->j + ii[ix];
7558be712e4SBarry Smith               for (int jj = 0; jj < n && isOK; jj++) {
7568be712e4SBarry Smith                 PetscInt lidj = aj[jj];
7578be712e4SBarry Smith                 if (lghost_matched[lidj]) continue;
7588be712e4SBarry Smith                 ew = PetscRealPart(ap[jj]);
7598be712e4SBarry Smith                 if (ew <= threshold) continue;
7608be712e4SBarry Smith                 max_e = PetscRealPart(lghost_max_ew[lidj]);
7618be712e4SBarry Smith                 /* check for max_e == to this edge and larger processor that will deal with this */
762452acf8bSMark Adams                 if (ew >= PetscRealPart(lid_max_ew[lid1]) - MY_MEPS && lghost_max_pe[lidj] > rank) isOK = PETSC_FALSE;
7638be712e4SBarry 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));
7648be712e4SBarry Smith                 if (print)
765452acf8bSMark Adams                   PetscCall(PetscSynchronizedPrintf(comm, "\t\t\t\t\t[%d] e1: looked at ghost adj (%d %d), diff = %10.4e, ghost on proc %d (max %d)\n", rank, (int)gid0, (int)lghost_gid[lidj], (double)(max_e - ew), lghost_pe[lidj], lghost_max_pe[lidj]));
7668be712e4SBarry Smith               }
7678be712e4SBarry Smith             }
768452acf8bSMark Adams             if (!isOK && print) PetscCall(PetscSynchronizedPrintf(comm, "\t\t[%d] skip edge (%d %d) from ghost inspection\n", rank, (int)gid0, (int)gid1));
7698be712e4SBarry Smith           }
7708be712e4SBarry Smith         }
7718be712e4SBarry Smith         PetscReal e1_max_w = (ghost1_idx == -1 ? PetscRealPart(lid_max_ew[lid0]) : PetscRealPart(lghost_max_ew[ghost1_idx]));
7728be712e4SBarry Smith         if (print)
773d8b4a066SPierre Jolivet           PetscCall(PetscSynchronizedPrintf(comm, "\t[%d] MATCHING (%d %d) e1 max weight = %e, e1 weight diff %e, %s. isOK = %d\n", rank, (int)gid0, (int)gid1, (double)e1_max_w, (double)(e1_max_w - e->weight), ghost1_idx == -1 ? "local" : "ghost", isOK));
7748be712e4SBarry Smith         /* do it */
7758be712e4SBarry Smith         if (isOK) {
7768be712e4SBarry Smith           if (ghost1_idx == -1) {
7778be712e4SBarry Smith             PetscCheck(!lid_matched[lid1], PETSC_COMM_SELF, PETSC_ERR_SUP, "local %d is matched", (int)gid1);
7788be712e4SBarry Smith             lid_matched[lid1] = PETSC_TRUE;                       /* keep track of what we've done this round */
7798be712e4SBarry Smith             PetscCall(PetscCDMoveAppend(agg_llists, lid0, lid1)); // takes lid1's list and appends to lid0's
7808be712e4SBarry Smith           } else {
7818be712e4SBarry Smith             /* add gid1 to list of ghost deleted by me -- I need their children */
7828be712e4SBarry Smith             PetscMPIInt proc = lghost_pe[ghost1_idx];
7838be712e4SBarry Smith             PetscCheck(!lghost_matched[ghost1_idx], PETSC_COMM_SELF, PETSC_ERR_SUP, "ghost %d is matched", (int)lghost_gid[ghost1_idx]);
7848be712e4SBarry Smith             lghost_matched[ghost1_idx] = PETSC_TRUE;
7858be712e4SBarry Smith             PetscCall(PetscCDAppendID(ghost_deleted_list, proc, ghost1_idx)); /* cache to send messages */
7868be712e4SBarry Smith             PetscCall(PetscCDAppendID(ghost_deleted_list, proc, lid0));
7878be712e4SBarry Smith           }
7888be712e4SBarry Smith           lid_matched[lid0] = PETSC_TRUE; /* keep track of what we've done this round */
7898be712e4SBarry Smith           /* set projection */
7908be712e4SBarry Smith           PetscCall(MatSetValues(P, 1, &gid0, 1, &gid0, &one, INSERT_VALUES));
7918be712e4SBarry Smith           PetscCall(MatSetValues(P, 1, &gid1, 1, &gid0, &one, INSERT_VALUES));
792452acf8bSMark Adams           //PetscCall(PetscPrintf(comm,"\t %d.%d) match active EDGE %d : (%d %d)\n",iter,sub_it, (int)nactive_edges, (int)gid0, (int)gid1));
7938be712e4SBarry Smith         } /* matched */
7948be712e4SBarry Smith       } /* edge loop */
795452acf8bSMark Adams       PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
7968be712e4SBarry Smith       if (isMPI) PetscCall(VecRestoreArrayRead(ghostMaxEdge, &lghost_max_ew));
7978be712e4SBarry Smith       PetscCall(VecRestoreArrayRead(locMaxEdge, &lid_max_ew));
7988be712e4SBarry Smith       // count active for test, latter, update deleted ghosts
7998be712e4SBarry Smith       n_act_n[0] = nactive_edges;
8008be712e4SBarry Smith       if (ghost_deleted_list) PetscCall(PetscCDCount(ghost_deleted_list, &n_act_n[2]));
8018be712e4SBarry Smith       else n_act_n[2] = 0;
8028be712e4SBarry Smith       PetscCall(PetscCDCount(agg_llists, &n_act_n[1]));
8038be712e4SBarry Smith       PetscCall(MPIU_Allreduce(n_act_n, gn_act_n, 3, MPIU_INT, MPI_SUM, comm));
8048be712e4SBarry Smith       PetscCall(PetscInfo(info_is, "[%d] %d.%d) nactive edges=%" PetscInt_FMT ", ncomm_procs=%d, nEdges=%d, %" PetscInt_FMT " deleted ghosts, N=%" PetscInt_FMT "\n", rank, iter, sub_it, gn_act_n[0], (int)ncomm_procs, (int)nEdges, gn_act_n[2], gn_act_n[1]));
8058be712e4SBarry Smith       /* deal with deleted ghost */
8068be712e4SBarry Smith       if (isMPI) {
8078be712e4SBarry Smith         PetscCDIntNd *pos;
8088be712e4SBarry Smith         PetscInt     *sbuffs1[REQ_BF_SIZE], ndel;
8098be712e4SBarry Smith         PetscInt     *sbuffs2[REQ_BF_SIZE];
8108be712e4SBarry Smith         MPI_Status    status;
8118be712e4SBarry Smith 
8128be712e4SBarry Smith         /* send deleted ghosts */
8138be712e4SBarry Smith         for (int proc_idx = 0; proc_idx < ncomm_procs; proc_idx++) {
8148be712e4SBarry Smith           const PetscMPIInt proc = comm_procs[proc_idx];
8158be712e4SBarry Smith           PetscInt         *sbuff, *pt, scount;
8168be712e4SBarry Smith           MPI_Request      *request;
8178be712e4SBarry Smith           /* count ghosts */
8188be712e4SBarry Smith           PetscCall(PetscCDCountAt(ghost_deleted_list, proc, &ndel));
8198be712e4SBarry Smith           ndel /= 2; // two entries for each proc
8208be712e4SBarry Smith           scount = 2 + 2 * ndel;
8218be712e4SBarry Smith           PetscCall(PetscMalloc1(scount + request_size, &sbuff));
8228be712e4SBarry Smith           /* save requests */
8238be712e4SBarry Smith           sbuffs1[proc_idx] = sbuff;
8248be712e4SBarry Smith           request           = (MPI_Request *)sbuff;
8258be712e4SBarry Smith           sbuff = pt = (PetscInt *)(sbuff + request_size);
8268be712e4SBarry Smith           /* write [ndel, proc, n*[gid1,gid0] */
8278be712e4SBarry Smith           *pt++ = ndel; // number of deleted to send
8288be712e4SBarry Smith           *pt++ = rank; // proc (not used)
8298be712e4SBarry Smith           PetscCall(PetscCDGetHeadPos(ghost_deleted_list, proc, &pos));
8308be712e4SBarry Smith           while (pos) {
8318be712e4SBarry Smith             PetscInt lid0, ghost_idx, gid1;
8328be712e4SBarry Smith             PetscCall(PetscCDIntNdGetID(pos, &ghost_idx));
8338be712e4SBarry Smith             gid1 = lghost_gid[ghost_idx];
8348be712e4SBarry Smith             PetscCall(PetscCDGetNextPos(ghost_deleted_list, proc, &pos));
8358be712e4SBarry Smith             PetscCall(PetscCDIntNdGetID(pos, &lid0));
8368be712e4SBarry Smith             PetscCall(PetscCDGetNextPos(ghost_deleted_list, proc, &pos));
8378be712e4SBarry Smith             *pt++ = gid1;
8388be712e4SBarry Smith             *pt++ = lid0 + my0; // gid0
8398be712e4SBarry Smith           }
8408be712e4SBarry Smith           PetscCheck(pt - sbuff == scount, PETSC_COMM_SELF, PETSC_ERR_SUP, "sbuff-pt != scount: %d", (int)(pt - sbuff));
8418be712e4SBarry Smith           /* MPI_Isend:  tag1 [ndel, proc, n*[gid1,gid0] ] */
8428be712e4SBarry Smith           PetscCallMPI(MPI_Isend(sbuff, scount, MPIU_INT, proc, tag1, comm, request));
8438be712e4SBarry Smith           PetscCall(PetscCDRemoveAllAt(ghost_deleted_list, proc)); // done with this list
8448be712e4SBarry Smith         }
8458be712e4SBarry Smith         /* receive deleted, send back partial aggregates, clear lists */
8468be712e4SBarry Smith         for (int proc_idx = 0; proc_idx < ncomm_procs; proc_idx++) {
8478be712e4SBarry Smith           PetscCallMPI(MPI_Probe(comm_procs[proc_idx] /* MPI_ANY_SOURCE */, tag1, comm, &status));
8488be712e4SBarry Smith           {
8498be712e4SBarry Smith             PetscInt         *pt, *pt2, *pt3, *sbuff, tmp;
8508be712e4SBarry Smith             MPI_Request      *request;
8518be712e4SBarry Smith             int               rcount, scount, ndel;
8528be712e4SBarry Smith             const PetscMPIInt proc = status.MPI_SOURCE;
8538be712e4SBarry Smith             PetscCallMPI(MPI_Get_count(&status, MPIU_INT, &rcount));
8548be712e4SBarry Smith             if (rcount > rbuff_sz) {
8558be712e4SBarry Smith               if (rbuff) PetscCall(PetscFree(rbuff));
8568be712e4SBarry Smith               PetscCall(PetscMalloc1(rcount, &rbuff));
8578be712e4SBarry Smith               rbuff_sz = rcount;
8588be712e4SBarry Smith             }
8598be712e4SBarry Smith             /* MPI_Recv: tag1 [ndel, proc, ndel*[gid1,gid0] ] */
8608be712e4SBarry Smith             PetscCallMPI(MPI_Recv(rbuff, rcount, MPIU_INT, proc, tag1, comm, &status));
8618be712e4SBarry Smith             /* read and count sends *[lid0, n, n*[gid] ] */
8628be712e4SBarry Smith             pt     = rbuff;
8638be712e4SBarry Smith             scount = 0;
8648be712e4SBarry Smith             ndel   = *pt++; // number of deleted to recv
8658be712e4SBarry Smith             tmp    = *pt++; // proc (not used)
8668be712e4SBarry Smith             while (ndel--) {
8678be712e4SBarry Smith               PetscInt gid1 = *pt++, lid1 = gid1 - my0;
8688be712e4SBarry Smith               int      gh_gid0 = *pt++; // gid on other proc (not used here to count)
8698be712e4SBarry Smith               PetscCheck(lid1 >= 0 && lid1 < nloc, PETSC_COMM_SELF, PETSC_ERR_SUP, "received ghost deleted %d", (int)gid1);
8708be712e4SBarry Smith               PetscCheck(!lid_matched[lid1], PETSC_COMM_SELF, PETSC_ERR_PLIB, "%d) received matched local gid %" PetscInt_FMT ",%d, with ghost (lid) %" PetscInt_FMT " from proc %d", sub_it, gid1, gh_gid0, tmp, proc);
8718be712e4SBarry Smith               lid_matched[lid1] = PETSC_TRUE;                    /* keep track of what we've done this round */
8728be712e4SBarry Smith               PetscCall(PetscCDCountAt(agg_llists, lid1, &tmp)); // n
8738be712e4SBarry Smith               /* PetscCheck(tmp == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "sending %d (!= 1) size aggregate. gid-0 %d, from %d (gid-1 %d)", (int)tmp, (int) gid, proc, gh_gid0); */
8748be712e4SBarry Smith               scount += tmp + 2; // lid0, n, n*[gid]
8758be712e4SBarry Smith             }
8768be712e4SBarry Smith             PetscCheck((pt - rbuff) == rcount, PETSC_COMM_SELF, PETSC_ERR_SUP, "receive buffer size != num read: %d; rcount: %d", (int)(pt - rbuff), rcount);
8778be712e4SBarry Smith             /* send tag2: *[gid0, n, n*[gid] ] */
8788be712e4SBarry Smith             PetscCall(PetscMalloc1(scount + request_size, &sbuff));
8798be712e4SBarry Smith             sbuffs2[proc_idx] = sbuff; /* cache request */
8808be712e4SBarry Smith             request           = (MPI_Request *)sbuff;
8818be712e4SBarry Smith             pt2 = sbuff = sbuff + request_size;
8828be712e4SBarry Smith             // read again: n, proc, n*[gid1,gid0]
8838be712e4SBarry Smith             pt   = rbuff;
8848be712e4SBarry Smith             ndel = *pt++;
8858be712e4SBarry Smith             tmp  = *pt++; // proc (not used)
8868be712e4SBarry Smith             while (ndel--) {
8878be712e4SBarry Smith               PetscInt gid1 = *pt++, lid1 = gid1 - my0, gh_gid0 = *pt++;
8888be712e4SBarry Smith               /* write [gid0, aggSz, aggSz[gid] ] */
8898be712e4SBarry Smith               *pt2++ = gh_gid0;
8908be712e4SBarry Smith               pt3    = pt2++; /* save pointer for later */
8918be712e4SBarry Smith               PetscCall(PetscCDGetHeadPos(agg_llists, lid1, &pos));
8928be712e4SBarry Smith               while (pos) {
8938be712e4SBarry Smith                 PetscInt gid;
8948be712e4SBarry Smith                 PetscCall(PetscCDIntNdGetID(pos, &gid));
8958be712e4SBarry Smith                 PetscCall(PetscCDGetNextPos(agg_llists, lid1, &pos));
8968be712e4SBarry Smith                 *pt2++ = gid;
8978be712e4SBarry Smith               }
8988be712e4SBarry Smith               *pt3 = (pt2 - pt3) - 1;
8998be712e4SBarry Smith               /* clear list */
9008be712e4SBarry Smith               PetscCall(PetscCDRemoveAllAt(agg_llists, lid1));
9018be712e4SBarry Smith             }
9028be712e4SBarry Smith             PetscCheck((pt2 - sbuff) == scount, PETSC_COMM_SELF, PETSC_ERR_SUP, "buffer size != num write: %d %d", (int)(pt2 - sbuff), (int)scount);
9038be712e4SBarry Smith             /* MPI_Isend: requested data tag2 *[lid0, n, n*[gid1] ] */
9048be712e4SBarry Smith             PetscCallMPI(MPI_Isend(sbuff, scount, MPIU_INT, proc, tag2, comm, request));
9058be712e4SBarry Smith           }
9068be712e4SBarry Smith         } // proc_idx
9078be712e4SBarry Smith         /* receive tag2 *[gid0, n, n*[gid] ] */
9088be712e4SBarry Smith         for (int proc_idx = 0; proc_idx < ncomm_procs; proc_idx++) {
9098be712e4SBarry Smith           PetscMPIInt proc;
9108be712e4SBarry Smith           PetscInt   *pt;
9118be712e4SBarry Smith           int         rcount;
9128be712e4SBarry Smith           PetscCallMPI(MPI_Probe(comm_procs[proc_idx] /* MPI_ANY_SOURCE */, tag2, comm, &status));
9138be712e4SBarry Smith           PetscCallMPI(MPI_Get_count(&status, MPIU_INT, &rcount));
9148be712e4SBarry Smith           if (rcount > rbuff_sz) {
9158be712e4SBarry Smith             if (rbuff) PetscCall(PetscFree(rbuff));
9168be712e4SBarry Smith             PetscCall(PetscMalloc1(rcount, &rbuff));
9178be712e4SBarry Smith             rbuff_sz = rcount;
9188be712e4SBarry Smith           }
9198be712e4SBarry Smith           proc = status.MPI_SOURCE;
9208be712e4SBarry Smith           /* MPI_Recv:  tag1 [n, proc, n*[gid1,lid0] ] */
9218be712e4SBarry Smith           PetscCallMPI(MPI_Recv(rbuff, rcount, MPIU_INT, proc, tag2, comm, &status));
9228be712e4SBarry Smith           pt = rbuff;
9238be712e4SBarry Smith           while (pt - rbuff < rcount) {
9248be712e4SBarry Smith             PetscInt gid0 = *pt++, n = *pt++;
9258be712e4SBarry Smith             while (n--) {
9268be712e4SBarry Smith               PetscInt gid1 = *pt++;
9278be712e4SBarry Smith               PetscCall(PetscCDAppendID(agg_llists, gid0 - my0, gid1));
9288be712e4SBarry Smith             }
9298be712e4SBarry Smith           }
9308be712e4SBarry Smith           PetscCheck((pt - rbuff) == rcount, PETSC_COMM_SELF, PETSC_ERR_SUP, "recv buffer size != num read: %d %d", (int)(pt - rbuff), (int)rcount);
9318be712e4SBarry Smith         }
9328be712e4SBarry Smith         /* wait for tag1 isends */
9338be712e4SBarry Smith         for (int proc_idx = 0; proc_idx < ncomm_procs; proc_idx++) {
9348be712e4SBarry Smith           MPI_Request *request;
9358be712e4SBarry Smith           request = (MPI_Request *)sbuffs1[proc_idx];
9368be712e4SBarry Smith           PetscCallMPI(MPI_Wait(request, &status));
9378be712e4SBarry Smith           PetscCall(PetscFree(sbuffs1[proc_idx]));
9388be712e4SBarry Smith         }
9398be712e4SBarry Smith         /* wait for tag2 isends */
9408be712e4SBarry Smith         for (int proc_idx = 0; proc_idx < ncomm_procs; proc_idx++) {
9418be712e4SBarry Smith           MPI_Request *request = (MPI_Request *)sbuffs2[proc_idx];
9428be712e4SBarry Smith           PetscCallMPI(MPI_Wait(request, &status));
9438be712e4SBarry Smith           PetscCall(PetscFree(sbuffs2[proc_idx]));
9448be712e4SBarry Smith         }
9458be712e4SBarry Smith       } /* MPI */
9468be712e4SBarry Smith       /* set 'lghost_matched' - use locMaxEdge, ghostMaxEdge (recomputed next) */
9478be712e4SBarry Smith       if (isMPI) {
9488be712e4SBarry Smith         const PetscScalar *sbuff;
9498be712e4SBarry Smith         for (PetscInt kk = 0, gid = my0; kk < nloc; kk++, gid++) {
9508be712e4SBarry Smith           PetscScalar vval = lid_matched[kk] ? 1.0 : 0.0;
9518be712e4SBarry Smith           PetscCall(VecSetValues(locMaxEdge, 1, &gid, &vval, INSERT_VALUES)); /* set with GID */
9528be712e4SBarry Smith         }
9538be712e4SBarry Smith         PetscCall(VecAssemblyBegin(locMaxEdge));
9548be712e4SBarry Smith         PetscCall(VecAssemblyEnd(locMaxEdge));
9558be712e4SBarry Smith         PetscCall(VecScatterBegin(mpimat->Mvctx, locMaxEdge, ghostMaxEdge, INSERT_VALUES, SCATTER_FORWARD));
9568be712e4SBarry Smith         PetscCall(VecScatterEnd(mpimat->Mvctx, locMaxEdge, ghostMaxEdge, INSERT_VALUES, SCATTER_FORWARD));
9578be712e4SBarry Smith         PetscCall(VecGetArrayRead(ghostMaxEdge, &sbuff));
9588be712e4SBarry Smith         for (int kk = 0; kk < num_ghosts; kk++) { lghost_matched[kk] = (PetscBool)(PetscRealPart(sbuff[kk]) != 0.0); }
9598be712e4SBarry Smith         PetscCall(VecRestoreArrayRead(ghostMaxEdge, &sbuff));
9608be712e4SBarry Smith       }
9618be712e4SBarry Smith       /* compute 'locMaxEdge' inside sub iteration b/c max weight can drop as neighbors are matched */
9628be712e4SBarry Smith       for (PetscInt kk = 0, gid = my0; kk < nloc; kk++, gid++) {
9638be712e4SBarry Smith         PetscReal      max_e = 0., tt;
9648be712e4SBarry Smith         PetscScalar    vval;
9658be712e4SBarry Smith         const PetscInt lid    = kk;
9668be712e4SBarry Smith         int            max_pe = rank, pe, n;
9678be712e4SBarry Smith         ii                    = matA->i;
9688be712e4SBarry Smith         n                     = ii[lid + 1] - ii[lid];
9698e3a54c0SPierre Jolivet         aj                    = PetscSafePointerPlusOffset(matA->j, ii[lid]);
9708e3a54c0SPierre Jolivet         ap                    = PetscSafePointerPlusOffset(matA->a, ii[lid]);
9718be712e4SBarry Smith         for (int jj = 0; jj < n; jj++) {
9728be712e4SBarry Smith           PetscInt lidj = aj[jj];
9738be712e4SBarry Smith           if (lid_matched[lidj]) continue; /* this is new - can change local max */
9748be712e4SBarry Smith           if (lidj != lid && PetscRealPart(ap[jj]) > max_e) max_e = PetscRealPart(ap[jj]);
9758be712e4SBarry Smith         }
9768be712e4SBarry Smith         if (lid_cprowID && (ix = lid_cprowID[lid]) != -1) { /* if I have any ghost neighbors */
9778be712e4SBarry Smith           ii = matB->compressedrow.i;
9788be712e4SBarry Smith           n  = ii[ix + 1] - ii[ix];
9798be712e4SBarry Smith           ap = matB->a + ii[ix];
9808be712e4SBarry Smith           aj = matB->j + ii[ix];
9818be712e4SBarry Smith           for (int jj = 0; jj < n; jj++) {
9828be712e4SBarry Smith             PetscInt lidj = aj[jj];
9838be712e4SBarry Smith             if (lghost_matched[lidj]) continue;
9848be712e4SBarry Smith             if ((tt = PetscRealPart(ap[jj])) > max_e) max_e = tt;
9858be712e4SBarry Smith           }
9868be712e4SBarry Smith         }
9878be712e4SBarry Smith         vval = (PetscScalar)max_e;
9888be712e4SBarry Smith         PetscCall(VecSetValues(locMaxEdge, 1, &gid, &vval, INSERT_VALUES)); /* set with GID */
9898be712e4SBarry Smith         // max PE with max edge
9908be712e4SBarry Smith         if (lid_cprowID && (ix = lid_cprowID[lid]) != -1) { /* if I have any ghost neighbors */
9918be712e4SBarry Smith           ii = matB->compressedrow.i;
9928be712e4SBarry Smith           n  = ii[ix + 1] - ii[ix];
9938be712e4SBarry Smith           ap = matB->a + ii[ix];
9948be712e4SBarry Smith           aj = matB->j + ii[ix];
9958be712e4SBarry Smith           for (int jj = 0; jj < n; jj++) {
9968be712e4SBarry Smith             PetscInt lidj = aj[jj];
9978be712e4SBarry Smith             if (lghost_matched[lidj]) continue;
9988be712e4SBarry Smith             if ((pe = lghost_pe[aj[jj]]) > max_pe && PetscRealPart(ap[jj]) >= max_e - MY_MEPS) { max_pe = pe; }
9998be712e4SBarry Smith           }
10008be712e4SBarry Smith         }
10018be712e4SBarry Smith         vval = (PetscScalar)max_pe;
10028be712e4SBarry Smith         PetscCall(VecSetValues(locMaxPE, 1, &gid, &vval, INSERT_VALUES));
10038be712e4SBarry Smith       }
10048be712e4SBarry Smith       PetscCall(VecAssemblyBegin(locMaxEdge));
10058be712e4SBarry Smith       PetscCall(VecAssemblyEnd(locMaxEdge));
10068be712e4SBarry Smith       PetscCall(VecAssemblyBegin(locMaxPE));
10078be712e4SBarry Smith       PetscCall(VecAssemblyEnd(locMaxPE));
10088be712e4SBarry Smith       /* compute 'lghost_max_ew' and 'lghost_max_pe' to get ready for next iteration*/
10098be712e4SBarry Smith       if (isMPI) {
10108be712e4SBarry Smith         const PetscScalar *buf;
10118be712e4SBarry Smith         PetscCall(VecScatterBegin(mpimat->Mvctx, locMaxEdge, ghostMaxEdge, INSERT_VALUES, SCATTER_FORWARD));
10128be712e4SBarry Smith         PetscCall(VecScatterEnd(mpimat->Mvctx, locMaxEdge, ghostMaxEdge, INSERT_VALUES, SCATTER_FORWARD));
10138be712e4SBarry Smith         PetscCall(VecScatterBegin(mpimat->Mvctx, locMaxPE, ghostMaxPE, INSERT_VALUES, SCATTER_FORWARD));
10148be712e4SBarry Smith         PetscCall(VecScatterEnd(mpimat->Mvctx, locMaxPE, ghostMaxPE, INSERT_VALUES, SCATTER_FORWARD));
10158be712e4SBarry Smith         PetscCall(VecGetArrayRead(ghostMaxPE, &buf));
10168be712e4SBarry Smith         for (int kk = 0; kk < num_ghosts; kk++) {
10178be712e4SBarry Smith           lghost_max_pe[kk] = (PetscMPIInt)PetscRealPart(buf[kk]); // the MAX proc of the ghost now
10188be712e4SBarry Smith         }
10198be712e4SBarry Smith         PetscCall(VecRestoreArrayRead(ghostMaxPE, &buf));
10208be712e4SBarry Smith       }
10218be712e4SBarry Smith       // if no active edges, stop
10228be712e4SBarry Smith       if (gn_act_n[0] < 1) break;
10238be712e4SBarry Smith       // inc and check (self stopping iteration
1024452acf8bSMark Adams       PetscCheck(old_num_edge != gn_act_n[0], PETSC_COMM_SELF, PETSC_ERR_SUP, "HEM stalled step %d/%d", sub_it + 1, n_sub_its);
10258be712e4SBarry Smith       sub_it++;
10268be712e4SBarry Smith       PetscCheck(sub_it < n_sub_its, PETSC_COMM_SELF, PETSC_ERR_SUP, "failed to finish HEM step %d/%d", sub_it + 1, n_sub_its);
1027452acf8bSMark Adams       old_num_edge = gn_act_n[0];
10288be712e4SBarry Smith     } /* sub_it loop */
10298be712e4SBarry Smith     /* clean up iteration */
10308be712e4SBarry Smith     PetscCall(PetscFree(Edges));
10318be712e4SBarry Smith     if (mpimat) { // can be hoisted
10328be712e4SBarry Smith       PetscCall(VecRestoreArrayRead(ghostMaxEdge, &lghost_max_ew));
10338be712e4SBarry Smith       PetscCall(VecDestroy(&ghostMaxEdge));
10348be712e4SBarry Smith       PetscCall(VecDestroy(&ghostMaxPE));
10358be712e4SBarry Smith       PetscCall(PetscFree(lghost_pe));
10368be712e4SBarry Smith       PetscCall(PetscFree(lghost_gid));
10378be712e4SBarry Smith       PetscCall(PetscFree(lghost_matched));
10388be712e4SBarry Smith       PetscCall(PetscFree(lghost_max_pe));
10398be712e4SBarry Smith     }
10408be712e4SBarry Smith     PetscCall(VecDestroy(&locMaxEdge));
10418be712e4SBarry Smith     PetscCall(VecDestroy(&locMaxPE));
10428be712e4SBarry Smith     /* create next graph */
10438be712e4SBarry Smith     {
10448be712e4SBarry Smith       Vec diag;
10458be712e4SBarry Smith       /* add identity for unmatched vertices so they stay alive */
10468be712e4SBarry Smith       for (PetscInt kk = 0, gid1, gid = my0; kk < nloc; kk++, gid++) {
10478be712e4SBarry Smith         if (!lid_matched[kk]) {
10488be712e4SBarry Smith           const PetscInt lid = kk;
10498be712e4SBarry Smith           PetscCDIntNd  *pos;
10508be712e4SBarry Smith           PetscCall(PetscCDGetHeadPos(agg_llists, lid, &pos));
10518be712e4SBarry Smith           PetscCheck(pos, PETSC_COMM_SELF, PETSC_ERR_PLIB, "empty list in singleton: %d", (int)gid);
10528be712e4SBarry Smith           PetscCall(PetscCDIntNdGetID(pos, &gid1));
10538be712e4SBarry Smith           PetscCheck(gid1 == gid, PETSC_COMM_SELF, PETSC_ERR_PLIB, "first in list (%d) in singleton not %d", (int)gid1, (int)gid);
10548be712e4SBarry Smith           PetscCall(MatSetValues(P, 1, &gid, 1, &gid, &one, INSERT_VALUES));
10558be712e4SBarry Smith         }
10568be712e4SBarry Smith       }
10578be712e4SBarry Smith       PetscCall(MatAssemblyBegin(P, MAT_FINAL_ASSEMBLY));
10588be712e4SBarry Smith       PetscCall(MatAssemblyEnd(P, MAT_FINAL_ASSEMBLY));
10598be712e4SBarry Smith 
10608be712e4SBarry Smith       /* project to make new graph with collapsed edges */
10618be712e4SBarry Smith       PetscCall(MatPtAP(cMat, P, MAT_INITIAL_MATRIX, 1.0, &tMat));
10628be712e4SBarry Smith       PetscCall(MatDestroy(&P));
10638be712e4SBarry Smith       PetscCall(MatDestroy(&cMat));
10648be712e4SBarry Smith       cMat = tMat;
10658be712e4SBarry Smith       PetscCall(MatCreateVecs(cMat, &diag, NULL));
106653673ba5SSatish Balay       PetscCall(MatGetDiagonal(cMat, diag));
10678be712e4SBarry Smith       PetscCall(VecReciprocal(diag));
10688be712e4SBarry Smith       PetscCall(VecSqrtAbs(diag));
10698be712e4SBarry Smith       PetscCall(MatDiagonalScale(cMat, diag, diag));
10708be712e4SBarry Smith       PetscCall(VecDestroy(&diag));
10718be712e4SBarry Smith     }
10728be712e4SBarry Smith   } /* coarsen iterator */
10738be712e4SBarry Smith 
10748be712e4SBarry 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 */
10758be712e4SBarry Smith   if (size > 1) {
10768be712e4SBarry Smith     Mat           mat;
10778be712e4SBarry Smith     PetscCDIntNd *pos;
10788be712e4SBarry Smith     PetscInt      NN, MM, jj = 0, mxsz = 0;
10798be712e4SBarry Smith 
10808be712e4SBarry Smith     for (int kk = 0; kk < nloc; kk++) {
10818be712e4SBarry Smith       PetscCall(PetscCDCountAt(agg_llists, kk, &jj));
10828be712e4SBarry Smith       if (jj > mxsz) mxsz = jj;
10838be712e4SBarry Smith     }
10848be712e4SBarry Smith     PetscCall(MatGetSize(a_Gmat, &MM, &NN));
10858be712e4SBarry Smith     if (mxsz > MM - nloc) mxsz = MM - nloc;
10868be712e4SBarry Smith     /* matrix of ghost adj for square graph */
10878be712e4SBarry Smith     PetscCall(MatCreateAIJ(comm, nloc, nloc, PETSC_DETERMINE, PETSC_DETERMINE, 0, NULL, mxsz, NULL, &mat));
10888be712e4SBarry Smith     for (PetscInt lid = 0, gid = my0; lid < nloc; lid++, gid++) {
10898be712e4SBarry Smith       PetscCall(PetscCDGetHeadPos(agg_llists, lid, &pos));
10908be712e4SBarry Smith       while (pos) {
10918be712e4SBarry Smith         PetscInt gid1;
10928be712e4SBarry Smith         PetscCall(PetscCDIntNdGetID(pos, &gid1));
10938be712e4SBarry Smith         PetscCall(PetscCDGetNextPos(agg_llists, lid, &pos));
10948be712e4SBarry Smith         if (gid1 < my0 || gid1 >= my0 + nloc) PetscCall(MatSetValues(mat, 1, &gid, 1, &gid1, &one, ADD_VALUES));
10958be712e4SBarry Smith       }
10968be712e4SBarry Smith     }
10978be712e4SBarry Smith     PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
10988be712e4SBarry Smith     PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
10998be712e4SBarry Smith     PetscCall(PetscCDSetMat(agg_llists, mat));
11008be712e4SBarry Smith     PetscCall(PetscCDDestroy(ghost_deleted_list));
11018be712e4SBarry Smith     if (rbuff_sz) PetscCall(PetscFree(rbuff)); // always true
11028be712e4SBarry Smith   }
11038be712e4SBarry Smith   // move BCs into some node
11048be712e4SBarry Smith   if (bc_list) {
11058be712e4SBarry Smith     PetscCDIntNd *pos;
11068be712e4SBarry Smith     PetscCall(PetscCDGetHeadPos(bc_list, 0, &pos));
11078be712e4SBarry Smith     while (pos) {
11088be712e4SBarry Smith       PetscInt gid1;
11098be712e4SBarry Smith       PetscCall(PetscCDIntNdGetID(pos, &gid1));
11108be712e4SBarry Smith       PetscCall(PetscCDGetNextPos(bc_list, 0, &pos));
11118be712e4SBarry Smith       PetscCall(PetscCDAppendID(agg_llists, bc_agg, gid1));
11128be712e4SBarry Smith     }
11138be712e4SBarry Smith     PetscCall(PetscCDRemoveAllAt(bc_list, 0));
11148be712e4SBarry Smith     PetscCall(PetscCDDestroy(bc_list));
11158be712e4SBarry Smith   }
11168be712e4SBarry Smith   {
11178be712e4SBarry Smith     // check sizes -- all vertices must get in graph
11188be712e4SBarry Smith     PetscInt sz, globalsz, MM;
11198be712e4SBarry Smith     PetscCall(MatGetSize(a_Gmat, &MM, NULL));
11208be712e4SBarry Smith     PetscCall(PetscCDCount(agg_llists, &sz));
11218be712e4SBarry Smith     PetscCall(MPIU_Allreduce(&sz, &globalsz, 1, MPIU_INT, MPI_SUM, comm));
11228be712e4SBarry Smith     PetscCheck(MM == globalsz, comm, PETSC_ERR_SUP, "lost %d equations ?", (int)(MM - globalsz));
11238be712e4SBarry Smith   }
11248be712e4SBarry Smith   // cleanup
11258be712e4SBarry Smith   PetscCall(MatDestroy(&cMat));
11268be712e4SBarry Smith   PetscCall(PetscFree(lid_cprowID));
11278be712e4SBarry Smith   PetscCall(PetscFree(lid_max_pe));
11288be712e4SBarry Smith   PetscCall(PetscFree(lid_matched));
11298be712e4SBarry Smith   PetscCall(ISDestroy(&info_is));
11308be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
11318be712e4SBarry Smith }
11328be712e4SBarry Smith 
11338be712e4SBarry Smith /*
11348be712e4SBarry Smith    HEM coarsen, simple greedy.
11358be712e4SBarry Smith */
11368be712e4SBarry Smith static PetscErrorCode MatCoarsenApply_HEM(MatCoarsen coarse)
11378be712e4SBarry Smith {
11388be712e4SBarry Smith   Mat mat = coarse->graph;
11398be712e4SBarry Smith 
11408be712e4SBarry Smith   PetscFunctionBegin;
11418be712e4SBarry Smith   PetscCall(MatCoarsenApply_HEM_private(mat, coarse->max_it, coarse->threshold, &coarse->agg_lists));
11428be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
11438be712e4SBarry Smith }
11448be712e4SBarry Smith 
11458be712e4SBarry Smith static PetscErrorCode MatCoarsenView_HEM(MatCoarsen coarse, PetscViewer viewer)
11468be712e4SBarry Smith {
11478be712e4SBarry Smith   PetscMPIInt rank;
11488be712e4SBarry Smith   PetscBool   iascii;
11498be712e4SBarry Smith 
11508be712e4SBarry Smith   PetscFunctionBegin;
11518be712e4SBarry Smith   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)coarse), &rank));
11528be712e4SBarry Smith   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
11538be712e4SBarry Smith   if (iascii) {
11548be712e4SBarry Smith     PetscCDIntNd     *pos, *pos2;
1155*e0b7e82fSBarry Smith     PetscViewerFormat format;
1156*e0b7e82fSBarry Smith 
11578be712e4SBarry Smith     PetscCall(PetscViewerASCIIPrintf(viewer, "%d matching steps with threshold = %g\n", (int)coarse->max_it, (double)coarse->threshold));
1158*e0b7e82fSBarry Smith     PetscCall(PetscViewerGetFormat(viewer, &format));
1159*e0b7e82fSBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
11608be712e4SBarry Smith       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
11618be712e4SBarry Smith       for (PetscInt kk = 0; kk < coarse->agg_lists->size; kk++) {
11628be712e4SBarry Smith         PetscCall(PetscCDGetHeadPos(coarse->agg_lists, kk, &pos));
11638be712e4SBarry Smith         if ((pos2 = pos)) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "selected local %d: ", (int)kk));
11648be712e4SBarry Smith         while (pos) {
11658be712e4SBarry Smith           PetscInt gid1;
11668be712e4SBarry Smith           PetscCall(PetscCDIntNdGetID(pos, &gid1));
11678be712e4SBarry Smith           PetscCall(PetscCDGetNextPos(coarse->agg_lists, kk, &pos));
11688be712e4SBarry Smith           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " %d ", (int)gid1));
11698be712e4SBarry Smith         }
11708be712e4SBarry Smith         if (pos2) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n"));
11718be712e4SBarry Smith       }
11728be712e4SBarry Smith       PetscCall(PetscViewerFlush(viewer));
11738be712e4SBarry Smith       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
11748be712e4SBarry Smith     }
1175*e0b7e82fSBarry Smith   }
11768be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
11778be712e4SBarry Smith }
11788be712e4SBarry Smith 
11798be712e4SBarry Smith /*
11808be712e4SBarry Smith   MatCoarsenCreate_HEM - A coarsener that uses HEM a simple greedy coarsener
11818be712e4SBarry Smith */
11828be712e4SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenCreate_HEM(MatCoarsen coarse)
11838be712e4SBarry Smith {
11848be712e4SBarry Smith   PetscFunctionBegin;
11858be712e4SBarry Smith   coarse->ops->apply = MatCoarsenApply_HEM;
11868be712e4SBarry Smith   coarse->ops->view  = MatCoarsenView_HEM;
11878be712e4SBarry Smith   coarse->max_it     = 4;
11888be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
11898be712e4SBarry Smith }
1190