xref: /petsc/src/vec/is/sf/tests/ex18.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1f210f596SVaclav Hapla 
2f210f596SVaclav Hapla static char help[] = "Test PetscSFConcatenate()\n\n";
3f210f596SVaclav Hapla 
4f210f596SVaclav Hapla #include <petscsf.h>
5f210f596SVaclav Hapla 
6f210f596SVaclav Hapla typedef struct {
7f210f596SVaclav Hapla   MPI_Comm    comm;
8f210f596SVaclav Hapla   PetscMPIInt rank, size;
9f210f596SVaclav Hapla   PetscInt    leaveStep, nsfs, nLeavesPerRank;
10f210f596SVaclav Hapla   PetscBool   shareRoots, sparseLeaves;
11f210f596SVaclav Hapla } AppCtx;
12f210f596SVaclav Hapla 
139371c9d4SSatish Balay static PetscErrorCode GetOptions(MPI_Comm comm, AppCtx *ctx) {
14f210f596SVaclav Hapla   PetscFunctionBegin;
15f210f596SVaclav Hapla   ctx->comm           = comm;
16f210f596SVaclav Hapla   ctx->nsfs           = 3;
17f210f596SVaclav Hapla   ctx->nLeavesPerRank = 4;
18f210f596SVaclav Hapla   ctx->leaveStep      = 1;
19f210f596SVaclav Hapla   ctx->shareRoots     = PETSC_FALSE;
20f210f596SVaclav Hapla   ctx->sparseLeaves   = PETSC_FALSE;
219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nsfs", &ctx->nsfs, NULL));
229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n_leaves_per_rank", &ctx->nLeavesPerRank, NULL));
239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-leave_step", &ctx->leaveStep, NULL));
249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-share_roots", &ctx->shareRoots, NULL));
25f210f596SVaclav Hapla   ctx->sparseLeaves = (PetscBool)(ctx->leaveStep != 1);
269566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &ctx->size));
279566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &ctx->rank));
28f210f596SVaclav Hapla   PetscFunctionReturn(0);
29f210f596SVaclav Hapla }
30f210f596SVaclav Hapla 
319371c9d4SSatish Balay static PetscErrorCode PetscSFCheckEqual_Private(PetscSF sf0, PetscSF sf1) {
32f210f596SVaclav Hapla   PetscInt  nRoot, nLeave;
33f210f596SVaclav Hapla   Vec       vecRoot0, vecLeave0, vecRoot1, vecLeave1;
34f210f596SVaclav Hapla   MPI_Comm  comm;
35f210f596SVaclav Hapla   PetscBool flg;
36f210f596SVaclav Hapla 
37f210f596SVaclav Hapla   PetscFunctionBegin;
389566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)sf0, &comm));
399566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf0, &nRoot, NULL, NULL, NULL));
409566063dSJacob Faibussowitsch   PetscCall(PetscSFGetLeafRange(sf0, NULL, &nLeave));
41f210f596SVaclav Hapla   nLeave++;
429566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(comm, nRoot, PETSC_DECIDE, &vecRoot0));
439566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(comm, nLeave, PETSC_DECIDE, &vecLeave0));
449566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(vecRoot0, &vecRoot1));
459566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(vecLeave0, &vecLeave1));
46f210f596SVaclav Hapla   {
47f210f596SVaclav Hapla     PetscRandom rand;
48f210f596SVaclav Hapla 
499566063dSJacob Faibussowitsch     PetscCall(PetscRandomCreate(comm, &rand));
509566063dSJacob Faibussowitsch     PetscCall(PetscRandomSetFromOptions(rand));
519566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(vecRoot0, rand));
529566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(vecLeave0, rand));
539566063dSJacob Faibussowitsch     PetscCall(VecCopy(vecRoot0, vecRoot1));
549566063dSJacob Faibussowitsch     PetscCall(VecCopy(vecLeave0, vecLeave1));
559566063dSJacob Faibussowitsch     PetscCall(PetscRandomDestroy(&rand));
56f210f596SVaclav Hapla   }
57f210f596SVaclav Hapla 
589566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(sf0, vecRoot0, vecLeave0, ADD_VALUES, SCATTER_FORWARD));
599566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(sf0, vecRoot0, vecLeave0, ADD_VALUES, SCATTER_FORWARD));
609566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(sf1, vecRoot1, vecLeave1, ADD_VALUES, SCATTER_FORWARD));
619566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(sf1, vecRoot1, vecLeave1, ADD_VALUES, SCATTER_FORWARD));
629566063dSJacob Faibussowitsch   PetscCall(VecEqual(vecLeave0, vecLeave1, &flg));
63f210f596SVaclav Hapla   PetscCheck(flg, comm, PETSC_ERR_PLIB, "leave vectors differ");
64f210f596SVaclav Hapla 
659566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(sf0, vecLeave0, vecRoot0, ADD_VALUES, SCATTER_REVERSE));
669566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(sf0, vecLeave0, vecRoot0, ADD_VALUES, SCATTER_REVERSE));
679566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(sf1, vecLeave1, vecRoot1, ADD_VALUES, SCATTER_REVERSE));
689566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(sf1, vecLeave1, vecRoot1, ADD_VALUES, SCATTER_REVERSE));
699566063dSJacob Faibussowitsch   PetscCall(VecEqual(vecRoot0, vecRoot1, &flg));
70f210f596SVaclav Hapla   PetscCheck(flg, comm, PETSC_ERR_PLIB, "root vectors differ");
71f210f596SVaclav Hapla 
729566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vecRoot0));
739566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vecRoot1));
749566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vecLeave0));
759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vecLeave1));
76f210f596SVaclav Hapla   PetscFunctionReturn(0);
77f210f596SVaclav Hapla }
78f210f596SVaclav Hapla 
799371c9d4SSatish Balay PetscErrorCode CreateReferenceSF(AppCtx *ctx, PetscSF *refSF) {
80f210f596SVaclav Hapla   PetscInt     i, j, k, r;
81f210f596SVaclav Hapla   PetscInt    *ilocal = NULL;
82f210f596SVaclav Hapla   PetscSFNode *iremote;
83f210f596SVaclav Hapla   PetscInt     nLeaves = ctx->nsfs * ctx->nLeavesPerRank * ctx->size;
84f210f596SVaclav Hapla   PetscInt     nroots  = ctx->nLeavesPerRank * ctx->nsfs;
85f210f596SVaclav Hapla   PetscSF      sf;
86f210f596SVaclav Hapla 
87f210f596SVaclav Hapla   PetscFunctionBegin;
88f210f596SVaclav Hapla   ilocal = NULL;
89*48a46eb9SPierre Jolivet   if (ctx->sparseLeaves) PetscCall(PetscCalloc1(nLeaves + 1, &ilocal));
909566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nLeaves, &iremote));
919566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(ctx->comm, &sf));
92f210f596SVaclav Hapla   for (i = 0, j = 0; i < ctx->nsfs; i++) {
93f210f596SVaclav Hapla     for (r = 0; r < ctx->size; r++) {
94f210f596SVaclav Hapla       for (k = 0; k < ctx->nLeavesPerRank; k++, j++) {
959371c9d4SSatish Balay         if (ctx->sparseLeaves) { ilocal[j + 1] = ilocal[j] + ctx->leaveStep; }
96f210f596SVaclav Hapla         iremote[j].rank  = r;
97f210f596SVaclav Hapla         iremote[j].index = k + i * ctx->nLeavesPerRank;
98f210f596SVaclav Hapla       }
99f210f596SVaclav Hapla     }
100f210f596SVaclav Hapla   }
1019566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(sf, nroots, nLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
102f210f596SVaclav Hapla   *refSF = sf;
103f210f596SVaclav Hapla   PetscFunctionReturn(0);
104f210f596SVaclav Hapla }
105f210f596SVaclav Hapla 
1069371c9d4SSatish Balay PetscErrorCode CreateSFs(AppCtx *ctx, PetscSF *newSFs[], PetscInt *leafOffsets[]) {
107f210f596SVaclav Hapla   PetscInt  i;
108f210f596SVaclav Hapla   PetscInt *lOffsets = NULL;
109f210f596SVaclav Hapla   PetscSF  *sfs;
110f210f596SVaclav Hapla   PetscInt  nLeaves = ctx->nLeavesPerRank * ctx->size;
111f210f596SVaclav Hapla   PetscInt  nroots  = ctx->shareRoots ? ctx->nLeavesPerRank * ctx->nsfs : ctx->nLeavesPerRank;
112f210f596SVaclav Hapla 
113f210f596SVaclav Hapla   PetscFunctionBegin;
114*48a46eb9SPierre Jolivet   if (ctx->sparseLeaves) PetscCall(PetscCalloc1(ctx->nsfs + 1, &lOffsets));
1159566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ctx->nsfs, &sfs));
116f210f596SVaclav Hapla   for (i = 0; i < ctx->nsfs; i++) {
117f210f596SVaclav Hapla     PetscInt     j, k;
118f210f596SVaclav Hapla     PetscMPIInt  r;
119f210f596SVaclav Hapla     PetscInt    *ilocal = NULL;
120f210f596SVaclav Hapla     PetscSFNode *iremote;
121f210f596SVaclav Hapla 
122*48a46eb9SPierre Jolivet     if (ctx->sparseLeaves) PetscCall(PetscCalloc1(nLeaves + 1, &ilocal));
1239566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nLeaves, &iremote));
124f210f596SVaclav Hapla     for (r = 0, j = 0; r < ctx->size; r++) {
125f210f596SVaclav Hapla       for (k = 0; k < ctx->nLeavesPerRank; k++, j++) {
1269371c9d4SSatish Balay         if (ctx->sparseLeaves) { ilocal[j + 1] = ilocal[j] + ctx->leaveStep; }
127f210f596SVaclav Hapla         iremote[j].rank  = r;
128f210f596SVaclav Hapla         iremote[j].index = ctx->shareRoots ? k + i * ctx->nLeavesPerRank : k;
129f210f596SVaclav Hapla       }
130f210f596SVaclav Hapla     }
131f210f596SVaclav Hapla     if (ctx->sparseLeaves) lOffsets[i + 1] = lOffsets[i] + ilocal[j];
132f210f596SVaclav Hapla 
1339566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(ctx->comm, &sfs[i]));
1349566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(sfs[i], nroots, nLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
135f210f596SVaclav Hapla   }
136f210f596SVaclav Hapla   *newSFs      = sfs;
137f210f596SVaclav Hapla   *leafOffsets = lOffsets;
138f210f596SVaclav Hapla   PetscFunctionReturn(0);
139f210f596SVaclav Hapla }
140f210f596SVaclav Hapla 
1419371c9d4SSatish Balay PetscErrorCode DestroySFs(AppCtx *ctx, PetscSF *sfs[]) {
142f210f596SVaclav Hapla   PetscInt i;
143f210f596SVaclav Hapla 
144f210f596SVaclav Hapla   PetscFunctionBegin;
145*48a46eb9SPierre Jolivet   for (i = 0; i < ctx->nsfs; i++) PetscCall(PetscSFDestroy(&(*sfs)[i]));
1469566063dSJacob Faibussowitsch   PetscCall(PetscFree(*sfs));
147f210f596SVaclav Hapla   PetscFunctionReturn(0);
148f210f596SVaclav Hapla }
149f210f596SVaclav Hapla 
1509371c9d4SSatish Balay int main(int argc, char **argv) {
151f210f596SVaclav Hapla   AppCtx    ctx;
152f210f596SVaclav Hapla   PetscSF   sf, sfRef;
153f210f596SVaclav Hapla   PetscSF  *sfs;
154f210f596SVaclav Hapla   PetscInt *leafOffsets;
155f210f596SVaclav Hapla   MPI_Comm  comm;
156f210f596SVaclav Hapla 
157327415f7SBarry Smith   PetscFunctionBeginUser;
1589566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
159f210f596SVaclav Hapla   comm = PETSC_COMM_WORLD;
1609566063dSJacob Faibussowitsch   PetscCall(GetOptions(comm, &ctx));
161f210f596SVaclav Hapla 
1629566063dSJacob Faibussowitsch   PetscCall(CreateSFs(&ctx, &sfs, &leafOffsets));
1639566063dSJacob Faibussowitsch   PetscCall(PetscSFConcatenate(comm, ctx.nsfs, sfs, ctx.shareRoots, leafOffsets, &sf));
1649566063dSJacob Faibussowitsch   PetscCall(CreateReferenceSF(&ctx, &sfRef));
1659566063dSJacob Faibussowitsch   PetscCall(PetscSFCheckEqual_Private(sf, sfRef));
166f210f596SVaclav Hapla 
1679566063dSJacob Faibussowitsch   PetscCall(DestroySFs(&ctx, &sfs));
1689566063dSJacob Faibussowitsch   PetscCall(PetscFree(leafOffsets));
1699566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sf));
1709566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfRef));
1719566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
172b122ec5aSJacob Faibussowitsch   return 0;
173f210f596SVaclav Hapla }
174f210f596SVaclav Hapla 
175f210f596SVaclav Hapla /*TEST
176f210f596SVaclav Hapla   test:
177f210f596SVaclav Hapla     nsize: {{1 3}}
178f210f596SVaclav Hapla     args: -nsfs {{1 3}} -n_leaves_per_rank {{0 1 5}} -leave_step {{1 3}} -share_roots {{true false}}
179f210f596SVaclav Hapla TEST*/
180