xref: /petsc/src/vec/is/sf/tests/ex18.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
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;
9c1730cc1SVaclav Hapla   PetscInt                   leaveStep, nsfs, n;
10c1730cc1SVaclav Hapla   PetscBool                  sparseLeaves;
11c1730cc1SVaclav Hapla   PetscBool                  compare;
12c1730cc1SVaclav Hapla   PetscBool                  irregular;
13c1730cc1SVaclav Hapla   PetscSFConcatenateRootMode rootMode;
14c1730cc1SVaclav Hapla   PetscViewer                viewer;
15f210f596SVaclav Hapla } AppCtx;
16f210f596SVaclav Hapla 
17d71ae5a4SJacob Faibussowitsch static PetscErrorCode GetOptions(MPI_Comm comm, AppCtx *ctx)
18d71ae5a4SJacob Faibussowitsch {
19c1730cc1SVaclav Hapla   PetscViewerFormat format;
20c1730cc1SVaclav Hapla 
21f210f596SVaclav Hapla   PetscFunctionBegin;
22f210f596SVaclav Hapla   ctx->comm         = comm;
23f210f596SVaclav Hapla   ctx->nsfs         = 3;
24c1730cc1SVaclav Hapla   ctx->n            = 1;
25f210f596SVaclav Hapla   ctx->leaveStep    = 1;
26f210f596SVaclav Hapla   ctx->sparseLeaves = PETSC_FALSE;
27c1730cc1SVaclav Hapla   ctx->compare      = PETSC_FALSE;
28c1730cc1SVaclav Hapla   ctx->irregular    = PETSC_FALSE;
29c1730cc1SVaclav Hapla   ctx->rootMode     = PETSCSF_CONCATENATE_ROOTMODE_LOCAL;
30c1730cc1SVaclav Hapla   ctx->viewer       = NULL;
319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nsfs", &ctx->nsfs, NULL));
32c1730cc1SVaclav Hapla   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &ctx->n, NULL));
339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-leave_step", &ctx->leaveStep, NULL));
34c1730cc1SVaclav Hapla   PetscCall(PetscOptionsGetBool(NULL, NULL, "-irregular", &ctx->irregular, NULL));
35c1730cc1SVaclav Hapla   PetscCall(PetscOptionsGetBool(NULL, NULL, "-compare_to_reference", &ctx->compare, NULL));
36c1730cc1SVaclav Hapla   PetscCall(PetscOptionsGetEnum(NULL, NULL, "-root_mode", PetscSFConcatenateRootModes, (PetscEnum *)&ctx->rootMode, NULL));
37c1730cc1SVaclav Hapla   PetscCall(PetscOptionsGetViewer(comm, NULL, NULL, "-sf_view", &ctx->viewer, &format, NULL));
38c1730cc1SVaclav Hapla   if (ctx->viewer) PetscCall(PetscViewerPushFormat(ctx->viewer, format));
39f210f596SVaclav Hapla   ctx->sparseLeaves = (PetscBool)(ctx->leaveStep != 1);
409566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &ctx->size));
419566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &ctx->rank));
42*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
43f210f596SVaclav Hapla }
44f210f596SVaclav Hapla 
45d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSFCheckEqual_Private(PetscSF sf0, PetscSF sf1)
46d71ae5a4SJacob Faibussowitsch {
47f210f596SVaclav Hapla   PetscInt  nRoot, nLeave;
48f210f596SVaclav Hapla   Vec       vecRoot0, vecLeave0, vecRoot1, vecLeave1;
49f210f596SVaclav Hapla   MPI_Comm  comm;
50f210f596SVaclav Hapla   PetscBool flg;
51f210f596SVaclav Hapla 
52f210f596SVaclav Hapla   PetscFunctionBegin;
539566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)sf0, &comm));
549566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf0, &nRoot, NULL, NULL, NULL));
559566063dSJacob Faibussowitsch   PetscCall(PetscSFGetLeafRange(sf0, NULL, &nLeave));
56f210f596SVaclav Hapla   nLeave++;
579566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(comm, nRoot, PETSC_DECIDE, &vecRoot0));
589566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(comm, nLeave, PETSC_DECIDE, &vecLeave0));
599566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(vecRoot0, &vecRoot1));
609566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(vecLeave0, &vecLeave1));
61f210f596SVaclav Hapla   {
62f210f596SVaclav Hapla     PetscRandom rand;
63f210f596SVaclav Hapla 
649566063dSJacob Faibussowitsch     PetscCall(PetscRandomCreate(comm, &rand));
659566063dSJacob Faibussowitsch     PetscCall(PetscRandomSetFromOptions(rand));
669566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(vecRoot0, rand));
679566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(vecLeave0, rand));
689566063dSJacob Faibussowitsch     PetscCall(VecCopy(vecRoot0, vecRoot1));
699566063dSJacob Faibussowitsch     PetscCall(VecCopy(vecLeave0, vecLeave1));
709566063dSJacob Faibussowitsch     PetscCall(PetscRandomDestroy(&rand));
71f210f596SVaclav Hapla   }
72f210f596SVaclav Hapla 
739566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(sf0, vecRoot0, vecLeave0, ADD_VALUES, SCATTER_FORWARD));
749566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(sf0, vecRoot0, vecLeave0, ADD_VALUES, SCATTER_FORWARD));
759566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(sf1, vecRoot1, vecLeave1, ADD_VALUES, SCATTER_FORWARD));
769566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(sf1, vecRoot1, vecLeave1, ADD_VALUES, SCATTER_FORWARD));
779566063dSJacob Faibussowitsch   PetscCall(VecEqual(vecLeave0, vecLeave1, &flg));
78f210f596SVaclav Hapla   PetscCheck(flg, comm, PETSC_ERR_PLIB, "leave vectors differ");
79f210f596SVaclav Hapla 
809566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(sf0, vecLeave0, vecRoot0, ADD_VALUES, SCATTER_REVERSE));
819566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(sf0, vecLeave0, vecRoot0, ADD_VALUES, SCATTER_REVERSE));
829566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(sf1, vecLeave1, vecRoot1, ADD_VALUES, SCATTER_REVERSE));
839566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(sf1, vecLeave1, vecRoot1, ADD_VALUES, SCATTER_REVERSE));
849566063dSJacob Faibussowitsch   PetscCall(VecEqual(vecRoot0, vecRoot1, &flg));
85f210f596SVaclav Hapla   PetscCheck(flg, comm, PETSC_ERR_PLIB, "root vectors differ");
86f210f596SVaclav Hapla 
879566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vecRoot0));
889566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vecRoot1));
899566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vecLeave0));
909566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vecLeave1));
91*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
92f210f596SVaclav Hapla }
93f210f596SVaclav Hapla 
94c1730cc1SVaclav Hapla static PetscErrorCode PetscSFViewCustom(PetscSF sf, PetscViewer viewer)
95d71ae5a4SJacob Faibussowitsch {
96c1730cc1SVaclav Hapla   PetscMPIInt        rank;
97c1730cc1SVaclav Hapla   PetscInt           i, nroots, nleaves, nranks;
98c1730cc1SVaclav Hapla   const PetscInt    *ilocal;
99c1730cc1SVaclav Hapla   const PetscSFNode *iremote;
100c1730cc1SVaclav Hapla   PetscLayout        rootLayout;
101c1730cc1SVaclav Hapla   PetscInt          *gremote;
102c1730cc1SVaclav Hapla 
103c1730cc1SVaclav Hapla   PetscFunctionBegin;
104c1730cc1SVaclav Hapla   PetscCall(PetscSFSetUp(sf));
105c1730cc1SVaclav Hapla   PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, &iremote));
106c1730cc1SVaclav Hapla   PetscCall(PetscSFGetRootRanks(sf, &nranks, NULL, NULL, NULL, NULL));
107c1730cc1SVaclav Hapla   PetscCall(PetscSFGetGraphLayout(sf, &rootLayout, NULL, NULL, &gremote));
108c1730cc1SVaclav Hapla   PetscCheck(nroots == rootLayout->n, PetscObjectComm((PetscObject)sf), PETSC_ERR_PLIB, "Assertion failed: nroots == rootLayout->n");
109c1730cc1SVaclav Hapla   PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)sf, viewer));
110c1730cc1SVaclav Hapla   PetscCall(PetscViewerASCIIPushTab(viewer));
111c1730cc1SVaclav Hapla   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf), &rank));
112c1730cc1SVaclav Hapla   PetscCall(PetscViewerASCIIPushSynchronized(viewer));
113c1730cc1SVaclav Hapla   if (rank == 0) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "rank #leaves #roots\n"));
114c1730cc1SVaclav Hapla   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%2d] %7" PetscInt_FMT " %6" PetscInt_FMT "\n", rank, nleaves, nroots));
115c1730cc1SVaclav Hapla   PetscCall(PetscViewerFlush(viewer));
116c1730cc1SVaclav Hapla   if (rank == 0) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "leaves      roots       roots in global numbering\n"));
117c1730cc1SVaclav Hapla   for (i = 0; i < nleaves; i++)
118c1730cc1SVaclav Hapla     PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "(%2d, %2" PetscInt_FMT ") <- (%2" PetscInt_FMT ", %2" PetscInt_FMT ")  = %3" PetscInt_FMT "\n", rank, ilocal ? ilocal[i] : i, iremote[i].rank, iremote[i].index, gremote[i]));
119c1730cc1SVaclav Hapla   PetscCall(PetscViewerFlush(viewer));
120c1730cc1SVaclav Hapla   PetscCall(PetscViewerASCIIPopSynchronized(viewer));
121c1730cc1SVaclav Hapla   PetscCall(PetscViewerASCIIPopTab(viewer));
122c1730cc1SVaclav Hapla   PetscCall(PetscLayoutDestroy(&rootLayout));
123c1730cc1SVaclav Hapla   PetscCall(PetscFree(gremote));
124*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
125c1730cc1SVaclav Hapla }
126c1730cc1SVaclav Hapla 
127c1730cc1SVaclav Hapla PetscErrorCode CreateReferenceSF_Regular(AppCtx *ctx, PetscSF *refSF)
128c1730cc1SVaclav Hapla {
129c1730cc1SVaclav Hapla   PetscInt  j;
130f210f596SVaclav Hapla   PetscInt *ilocal  = NULL;
131c1730cc1SVaclav Hapla   PetscInt  nLeaves = ctx->nsfs * ctx->n * ctx->size;
132c1730cc1SVaclav Hapla   PetscInt  nroots  = ctx->n * ctx->nsfs;
133f210f596SVaclav Hapla   PetscSF   sf;
134f210f596SVaclav Hapla 
135f210f596SVaclav Hapla   PetscFunctionBegin;
136f210f596SVaclav Hapla   ilocal = NULL;
13748a46eb9SPierre Jolivet   if (ctx->sparseLeaves) PetscCall(PetscCalloc1(nLeaves + 1, &ilocal));
1389566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(ctx->comm, &sf));
139c1730cc1SVaclav Hapla   for (j = 0; j < nLeaves; j++) {
140c1730cc1SVaclav Hapla     if (ctx->sparseLeaves) ilocal[j + 1] = ilocal[j] + ctx->leaveStep;
141c1730cc1SVaclav Hapla   }
142c1730cc1SVaclav Hapla   switch (ctx->rootMode) {
143c1730cc1SVaclav Hapla   case PETSCSF_CONCATENATE_ROOTMODE_SHARED:
144c1730cc1SVaclav Hapla   case PETSCSF_CONCATENATE_ROOTMODE_LOCAL: {
145c1730cc1SVaclav Hapla     PetscInt     i, k;
146c1730cc1SVaclav Hapla     PetscMPIInt  r;
147c1730cc1SVaclav Hapla     PetscSFNode *iremote;
148c1730cc1SVaclav Hapla 
149c1730cc1SVaclav Hapla     PetscCall(PetscCalloc1(nLeaves, &iremote));
150f210f596SVaclav Hapla     for (i = 0, j = 0; i < ctx->nsfs; i++) {
151f210f596SVaclav Hapla       for (r = 0; r < ctx->size; r++) {
152c1730cc1SVaclav Hapla         for (k = 0; k < ctx->n; k++, j++) {
153f210f596SVaclav Hapla           iremote[j].rank  = r;
154c1730cc1SVaclav Hapla           iremote[j].index = k + i * ctx->n;
155f210f596SVaclav Hapla         }
156f210f596SVaclav Hapla       }
157f210f596SVaclav Hapla     }
1589566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(sf, nroots, nLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
159c1730cc1SVaclav Hapla   } break;
160c1730cc1SVaclav Hapla   case PETSCSF_CONCATENATE_ROOTMODE_GLOBAL: {
161c1730cc1SVaclav Hapla     PetscLayout map = NULL;
162c1730cc1SVaclav Hapla     PetscInt   *gremote;
163c1730cc1SVaclav Hapla 
164c1730cc1SVaclav Hapla     PetscCall(PetscLayoutCreateFromSizes(ctx->comm, nroots, PETSC_DECIDE, 1, &map));
165c1730cc1SVaclav Hapla     PetscCall(PetscMalloc1(nLeaves, &gremote));
166c1730cc1SVaclav Hapla     for (j = 0; j < nLeaves; j++) gremote[j] = j;
167c1730cc1SVaclav Hapla     PetscCall(PetscSFSetGraphLayout(sf, map, nLeaves, ilocal, PETSC_OWN_POINTER, gremote));
168c1730cc1SVaclav Hapla     PetscCall(PetscFree(gremote));
169c1730cc1SVaclav Hapla     PetscCall(PetscLayoutDestroy(&map));
170c1730cc1SVaclav Hapla   } break;
171c1730cc1SVaclav Hapla   default:
172c1730cc1SVaclav Hapla     SETERRQ(ctx->comm, PETSC_ERR_SUP, "unsupported rootmode %d", ctx->rootMode);
173c1730cc1SVaclav Hapla   }
174c1730cc1SVaclav Hapla   PetscCall(PetscObjectSetName((PetscObject)sf, "reference_sf"));
175c1730cc1SVaclav Hapla   if (ctx->viewer) PetscCall(PetscSFViewCustom(sf, ctx->viewer));
176f210f596SVaclav Hapla   *refSF = sf;
177*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
178f210f596SVaclav Hapla }
179f210f596SVaclav Hapla 
180c1730cc1SVaclav Hapla PetscErrorCode CreateSFs_Irregular(AppCtx *ctx, PetscSF *newSFs[], PetscInt *leafOffsets[])
181d71ae5a4SJacob Faibussowitsch {
182f210f596SVaclav Hapla   PetscInt  i;
183f210f596SVaclav Hapla   PetscInt *lOffsets = NULL;
184f210f596SVaclav Hapla   PetscSF  *sfs;
185c1730cc1SVaclav Hapla   PetscInt  nLeaves = ctx->n * ctx->size + (ctx->size - 1) * ctx->size / 2;
186c1730cc1SVaclav Hapla   PetscInt  nroots  = ctx->n + ctx->rank + ctx->nsfs - 1 + ctx->size - 1;
187f210f596SVaclav Hapla 
188f210f596SVaclav Hapla   PetscFunctionBegin;
18948a46eb9SPierre Jolivet   if (ctx->sparseLeaves) PetscCall(PetscCalloc1(ctx->nsfs + 1, &lOffsets));
1909566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ctx->nsfs, &sfs));
191f210f596SVaclav Hapla   for (i = 0; i < ctx->nsfs; i++) {
192c1730cc1SVaclav Hapla     PetscSF      sf;
193f210f596SVaclav Hapla     PetscInt     j, k;
194f210f596SVaclav Hapla     PetscMPIInt  r;
195f210f596SVaclav Hapla     PetscInt    *ilocal = NULL;
196f210f596SVaclav Hapla     PetscSFNode *iremote;
197c1730cc1SVaclav Hapla     char         name[32];
198f210f596SVaclav Hapla 
19948a46eb9SPierre Jolivet     if (ctx->sparseLeaves) PetscCall(PetscCalloc1(nLeaves + 1, &ilocal));
2009566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nLeaves, &iremote));
201c1730cc1SVaclav Hapla     for (r = ctx->size - 1, j = 0; r >= 0; r--) {
202c1730cc1SVaclav Hapla       for (k = 0; k < ctx->n + r; k++, j++) {
203ad540459SPierre Jolivet         if (ctx->sparseLeaves) ilocal[j + 1] = ilocal[j] + ctx->leaveStep;
204f210f596SVaclav Hapla         iremote[j].rank  = r;
205c1730cc1SVaclav Hapla         iremote[j].index = k + i + ctx->rank;
206f210f596SVaclav Hapla       }
207f210f596SVaclav Hapla     }
208f210f596SVaclav Hapla     if (ctx->sparseLeaves) lOffsets[i + 1] = lOffsets[i] + ilocal[j];
209f210f596SVaclav Hapla 
210c1730cc1SVaclav Hapla     PetscCall(PetscSFCreate(ctx->comm, &sf));
211c1730cc1SVaclav Hapla     PetscCall(PetscSFSetGraph(sf, nroots, nLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
212c1730cc1SVaclav Hapla     PetscCall(PetscSNPrintf(name, sizeof(name), "sf_%" PetscInt_FMT, i));
213c1730cc1SVaclav Hapla     PetscCall(PetscObjectSetName((PetscObject)sf, name));
214c1730cc1SVaclav Hapla     if (ctx->viewer) PetscCall(PetscSFViewCustom(sf, ctx->viewer));
215c1730cc1SVaclav Hapla     sfs[i] = sf;
216c1730cc1SVaclav Hapla   }
217c1730cc1SVaclav Hapla   *newSFs      = sfs;
218c1730cc1SVaclav Hapla   *leafOffsets = lOffsets;
219*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
220c1730cc1SVaclav Hapla }
221c1730cc1SVaclav Hapla 
222c1730cc1SVaclav Hapla PetscErrorCode CreateSFs_Regular(AppCtx *ctx, PetscSF *newSFs[], PetscInt *leafOffsets[])
223c1730cc1SVaclav Hapla {
224c1730cc1SVaclav Hapla   PetscInt                   i;
225c1730cc1SVaclav Hapla   PetscInt                  *lOffsets = NULL;
226c1730cc1SVaclav Hapla   PetscInt                   nLeaves  = ctx->n * ctx->size;
227c1730cc1SVaclav Hapla   PetscSF                   *sfs;
228c1730cc1SVaclav Hapla   PetscSFConcatenateRootMode mode = ctx->compare ? ctx->rootMode : PETSCSF_CONCATENATE_ROOTMODE_LOCAL;
229c1730cc1SVaclav Hapla 
230c1730cc1SVaclav Hapla   PetscFunctionBegin;
231c1730cc1SVaclav Hapla   if (ctx->sparseLeaves) PetscCall(PetscCalloc1(ctx->nsfs + 1, &lOffsets));
232c1730cc1SVaclav Hapla   PetscCall(PetscCalloc1(ctx->nsfs, &sfs));
233c1730cc1SVaclav Hapla   for (i = 0; i < ctx->nsfs; i++) {
234c1730cc1SVaclav Hapla     PetscSF   sf;
235c1730cc1SVaclav Hapla     PetscInt  j;
236c1730cc1SVaclav Hapla     PetscInt *ilocal = NULL;
237c1730cc1SVaclav Hapla     char      name[32];
238c1730cc1SVaclav Hapla 
239c1730cc1SVaclav Hapla     PetscCall(PetscSFCreate(ctx->comm, &sf));
240c1730cc1SVaclav Hapla     if (ctx->sparseLeaves) {
241c1730cc1SVaclav Hapla       PetscCall(PetscCalloc1(nLeaves + 1, &ilocal));
242c1730cc1SVaclav Hapla       for (j = 0; j < nLeaves; j++) { ilocal[j + 1] = ilocal[j] + ctx->leaveStep; }
243c1730cc1SVaclav Hapla       lOffsets[i + 1] = lOffsets[i] + ilocal[nLeaves];
244c1730cc1SVaclav Hapla     }
245c1730cc1SVaclav Hapla     switch (mode) {
246c1730cc1SVaclav Hapla     case PETSCSF_CONCATENATE_ROOTMODE_LOCAL: {
247c1730cc1SVaclav Hapla       PetscInt     k, nroots = ctx->n;
248c1730cc1SVaclav Hapla       PetscMPIInt  r;
249c1730cc1SVaclav Hapla       PetscSFNode *iremote;
250c1730cc1SVaclav Hapla 
251c1730cc1SVaclav Hapla       PetscCall(PetscMalloc1(nLeaves, &iremote));
252c1730cc1SVaclav Hapla       for (r = 0, j = 0; r < ctx->size; r++) {
253c1730cc1SVaclav Hapla         for (k = 0; k < ctx->n; k++, j++) {
254c1730cc1SVaclav Hapla           iremote[j].rank  = r;
255c1730cc1SVaclav Hapla           iremote[j].index = k;
256c1730cc1SVaclav Hapla         }
257c1730cc1SVaclav Hapla       }
258c1730cc1SVaclav Hapla       PetscCall(PetscSFSetGraph(sf, nroots, nLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
259c1730cc1SVaclav Hapla     } break;
260c1730cc1SVaclav Hapla     case PETSCSF_CONCATENATE_ROOTMODE_SHARED: {
261c1730cc1SVaclav Hapla       PetscInt     k, nroots = ctx->n * ctx->nsfs;
262c1730cc1SVaclav Hapla       PetscMPIInt  r;
263c1730cc1SVaclav Hapla       PetscSFNode *iremote;
264c1730cc1SVaclav Hapla 
265c1730cc1SVaclav Hapla       PetscCall(PetscMalloc1(nLeaves, &iremote));
266c1730cc1SVaclav Hapla       for (r = 0, j = 0; r < ctx->size; r++) {
267c1730cc1SVaclav Hapla         for (k = 0; k < ctx->n; k++, j++) {
268c1730cc1SVaclav Hapla           iremote[j].rank  = r;
269c1730cc1SVaclav Hapla           iremote[j].index = k + i * ctx->n;
270c1730cc1SVaclav Hapla         }
271c1730cc1SVaclav Hapla       }
272c1730cc1SVaclav Hapla       PetscCall(PetscSFSetGraph(sf, nroots, nLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
273c1730cc1SVaclav Hapla     } break;
274c1730cc1SVaclav Hapla     case PETSCSF_CONCATENATE_ROOTMODE_GLOBAL: {
275c1730cc1SVaclav Hapla       PetscInt    nroots = ctx->n;
276c1730cc1SVaclav Hapla       PetscLayout map    = NULL;
277c1730cc1SVaclav Hapla       PetscInt   *gremote;
278c1730cc1SVaclav Hapla 
279c1730cc1SVaclav Hapla       PetscCall(PetscLayoutCreateFromSizes(ctx->comm, nroots, PETSC_DECIDE, 1, &map));
280c1730cc1SVaclav Hapla       PetscCall(PetscMalloc1(nLeaves, &gremote));
281c1730cc1SVaclav Hapla       for (j = 0; j < nLeaves; j++) gremote[j] = j;
282c1730cc1SVaclav Hapla       PetscCall(PetscSFSetGraphLayout(sf, map, nLeaves, ilocal, PETSC_OWN_POINTER, gremote));
283c1730cc1SVaclav Hapla       PetscCall(PetscFree(gremote));
284c1730cc1SVaclav Hapla       PetscCall(PetscLayoutDestroy(&map));
285c1730cc1SVaclav Hapla     } break;
286c1730cc1SVaclav Hapla     default:
287c1730cc1SVaclav Hapla       SETERRQ(ctx->comm, PETSC_ERR_SUP, "unsupported rootmode %d", ctx->rootMode);
288c1730cc1SVaclav Hapla     }
289c1730cc1SVaclav Hapla     PetscCall(PetscSNPrintf(name, sizeof(name), "sf_%" PetscInt_FMT, i));
290c1730cc1SVaclav Hapla     PetscCall(PetscObjectSetName((PetscObject)sf, name));
291c1730cc1SVaclav Hapla     if (ctx->viewer) PetscCall(PetscSFViewCustom(sf, ctx->viewer));
292c1730cc1SVaclav Hapla     sfs[i] = sf;
293f210f596SVaclav Hapla   }
294f210f596SVaclav Hapla   *newSFs      = sfs;
295f210f596SVaclav Hapla   *leafOffsets = lOffsets;
296*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
297f210f596SVaclav Hapla }
298f210f596SVaclav Hapla 
299d71ae5a4SJacob Faibussowitsch PetscErrorCode DestroySFs(AppCtx *ctx, PetscSF *sfs[])
300d71ae5a4SJacob Faibussowitsch {
301f210f596SVaclav Hapla   PetscInt i;
302f210f596SVaclav Hapla 
303f210f596SVaclav Hapla   PetscFunctionBegin;
30448a46eb9SPierre Jolivet   for (i = 0; i < ctx->nsfs; i++) PetscCall(PetscSFDestroy(&(*sfs)[i]));
3059566063dSJacob Faibussowitsch   PetscCall(PetscFree(*sfs));
306*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
307f210f596SVaclav Hapla }
308f210f596SVaclav Hapla 
309d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
310d71ae5a4SJacob Faibussowitsch {
311c1730cc1SVaclav Hapla   AppCtx    ctx_;
312c1730cc1SVaclav Hapla   AppCtx   *ctx = &ctx_;
313c1730cc1SVaclav Hapla   PetscSF   sf;
314c1730cc1SVaclav Hapla   PetscSF  *sfs         = NULL;
315c1730cc1SVaclav Hapla   PetscInt *leafOffsets = NULL;
316f210f596SVaclav Hapla   MPI_Comm  comm;
317f210f596SVaclav Hapla 
318327415f7SBarry Smith   PetscFunctionBeginUser;
3199566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
320f210f596SVaclav Hapla   comm = PETSC_COMM_WORLD;
321c1730cc1SVaclav Hapla   PetscCall(GetOptions(comm, ctx));
322f210f596SVaclav Hapla 
323c1730cc1SVaclav Hapla   if (ctx->irregular) {
324c1730cc1SVaclav Hapla     PetscCall(CreateSFs_Irregular(ctx, &sfs, &leafOffsets));
325c1730cc1SVaclav Hapla   } else {
326c1730cc1SVaclav Hapla     PetscCall(CreateSFs_Regular(ctx, &sfs, &leafOffsets));
327c1730cc1SVaclav Hapla   }
328c1730cc1SVaclav Hapla   PetscCall(PetscSFConcatenate(comm, ctx->nsfs, sfs, ctx->rootMode, leafOffsets, &sf));
329c1730cc1SVaclav Hapla   PetscCall(PetscObjectSetName((PetscObject)sf, "result_sf"));
330c1730cc1SVaclav Hapla   if (ctx->viewer) {
331c1730cc1SVaclav Hapla     PetscCall(PetscPrintf(comm, "rootMode = %s:\n", PetscSFConcatenateRootModes[ctx->rootMode]));
332c1730cc1SVaclav Hapla     PetscCall(PetscSFViewCustom(sf, ctx->viewer));
333c1730cc1SVaclav Hapla   }
334c1730cc1SVaclav Hapla   if (ctx->compare) {
335c1730cc1SVaclav Hapla     PetscSF sfRef;
336c1730cc1SVaclav Hapla 
337c1730cc1SVaclav Hapla     PetscAssert(!ctx->irregular, comm, PETSC_ERR_SUP, "Combination  -compare_to_reference true -irregular true  not implemented");
338c1730cc1SVaclav Hapla     PetscCall(CreateReferenceSF_Regular(ctx, &sfRef));
3399566063dSJacob Faibussowitsch     PetscCall(PetscSFCheckEqual_Private(sf, sfRef));
340c1730cc1SVaclav Hapla     PetscCall(PetscSFDestroy(&sfRef));
341c1730cc1SVaclav Hapla   }
342c1730cc1SVaclav Hapla   PetscCall(DestroySFs(ctx, &sfs));
3439566063dSJacob Faibussowitsch   PetscCall(PetscFree(leafOffsets));
3449566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sf));
345c1730cc1SVaclav Hapla   if (ctx->viewer) {
346c1730cc1SVaclav Hapla     PetscCall(PetscViewerPopFormat(ctx->viewer));
347c1730cc1SVaclav Hapla     PetscCall(PetscViewerDestroy(&ctx->viewer));
348c1730cc1SVaclav Hapla   }
3499566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
350b122ec5aSJacob Faibussowitsch   return 0;
351f210f596SVaclav Hapla }
352f210f596SVaclav Hapla 
353f210f596SVaclav Hapla /*TEST
354f210f596SVaclav Hapla   test:
355f210f596SVaclav Hapla     nsize: {{1 3}}
356c1730cc1SVaclav Hapla     args: -compare_to_reference -nsfs {{1 3}} -n {{0 1 5}} -leave_step {{1 3}} -root_mode {{local shared global}}
357c1730cc1SVaclav Hapla 
358c1730cc1SVaclav Hapla   test:
359c1730cc1SVaclav Hapla     suffix: 2
360c1730cc1SVaclav Hapla     nsize: 2
361c1730cc1SVaclav Hapla     args: -irregular {{false true}separate output} -sf_view -nsfs 3 -n 1 -leave_step {{1 3}separate output} -root_mode {{local shared global}separate output}
362f210f596SVaclav Hapla TEST*/
363