1f210f596SVaclav Hapla static char help[] = "Test PetscSFConcatenate()\n\n";
2f210f596SVaclav Hapla
3f210f596SVaclav Hapla #include <petscsf.h>
4f210f596SVaclav Hapla
5f210f596SVaclav Hapla typedef struct {
6f210f596SVaclav Hapla MPI_Comm comm;
7f210f596SVaclav Hapla PetscMPIInt rank, size;
8c1730cc1SVaclav Hapla PetscInt leaveStep, nsfs, n;
9c1730cc1SVaclav Hapla PetscBool sparseLeaves;
10c1730cc1SVaclav Hapla PetscBool compare;
11c1730cc1SVaclav Hapla PetscBool irregular;
12c1730cc1SVaclav Hapla PetscSFConcatenateRootMode rootMode;
13c1730cc1SVaclav Hapla PetscViewer viewer;
14f210f596SVaclav Hapla } AppCtx;
15f210f596SVaclav Hapla
GetOptions(MPI_Comm comm,AppCtx * ctx)16d71ae5a4SJacob Faibussowitsch static PetscErrorCode GetOptions(MPI_Comm comm, AppCtx *ctx)
17d71ae5a4SJacob Faibussowitsch {
18c1730cc1SVaclav Hapla PetscViewerFormat format;
19c1730cc1SVaclav Hapla
20f210f596SVaclav Hapla PetscFunctionBegin;
21f210f596SVaclav Hapla ctx->comm = comm;
22f210f596SVaclav Hapla ctx->nsfs = 3;
23c1730cc1SVaclav Hapla ctx->n = 1;
24f210f596SVaclav Hapla ctx->leaveStep = 1;
25f210f596SVaclav Hapla ctx->sparseLeaves = PETSC_FALSE;
26c1730cc1SVaclav Hapla ctx->compare = PETSC_FALSE;
27c1730cc1SVaclav Hapla ctx->irregular = PETSC_FALSE;
28c1730cc1SVaclav Hapla ctx->rootMode = PETSCSF_CONCATENATE_ROOTMODE_LOCAL;
29c1730cc1SVaclav Hapla ctx->viewer = NULL;
309566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nsfs", &ctx->nsfs, NULL));
31c1730cc1SVaclav Hapla PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &ctx->n, NULL));
329566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-leave_step", &ctx->leaveStep, NULL));
33c1730cc1SVaclav Hapla PetscCall(PetscOptionsGetBool(NULL, NULL, "-irregular", &ctx->irregular, NULL));
34c1730cc1SVaclav Hapla PetscCall(PetscOptionsGetBool(NULL, NULL, "-compare_to_reference", &ctx->compare, NULL));
35c1730cc1SVaclav Hapla PetscCall(PetscOptionsGetEnum(NULL, NULL, "-root_mode", PetscSFConcatenateRootModes, (PetscEnum *)&ctx->rootMode, NULL));
36648c30bcSBarry Smith PetscCall(PetscOptionsCreateViewer(comm, NULL, NULL, "-sf_view", &ctx->viewer, &format, NULL));
37c1730cc1SVaclav Hapla if (ctx->viewer) PetscCall(PetscViewerPushFormat(ctx->viewer, format));
38f210f596SVaclav Hapla ctx->sparseLeaves = (PetscBool)(ctx->leaveStep != 1);
399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &ctx->size));
409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &ctx->rank));
413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
42f210f596SVaclav Hapla }
43f210f596SVaclav Hapla
PetscSFCheckEqual_Private(PetscSF sf0,PetscSF sf1)44d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSFCheckEqual_Private(PetscSF sf0, PetscSF sf1)
45d71ae5a4SJacob Faibussowitsch {
46f210f596SVaclav Hapla PetscInt nRoot, nLeave;
47f210f596SVaclav Hapla Vec vecRoot0, vecLeave0, vecRoot1, vecLeave1;
48f210f596SVaclav Hapla MPI_Comm comm;
49f210f596SVaclav Hapla PetscBool flg;
50f210f596SVaclav Hapla
51f210f596SVaclav Hapla PetscFunctionBegin;
529566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sf0, &comm));
539566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf0, &nRoot, NULL, NULL, NULL));
549566063dSJacob Faibussowitsch PetscCall(PetscSFGetLeafRange(sf0, NULL, &nLeave));
55f210f596SVaclav Hapla nLeave++;
5677433607SBarry Smith PetscCall(VecCreateFromOptions(comm, NULL, 1, nRoot, PETSC_DECIDE, &vecRoot0));
5777433607SBarry Smith PetscCall(VecCreateFromOptions(comm, NULL, 1, nLeave, PETSC_DECIDE, &vecLeave0));
589566063dSJacob Faibussowitsch PetscCall(VecDuplicate(vecRoot0, &vecRoot1));
599566063dSJacob Faibussowitsch PetscCall(VecDuplicate(vecLeave0, &vecLeave1));
60f210f596SVaclav Hapla {
61f210f596SVaclav Hapla PetscRandom rand;
62f210f596SVaclav Hapla
639566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(comm, &rand));
649566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand));
659566063dSJacob Faibussowitsch PetscCall(VecSetRandom(vecRoot0, rand));
669566063dSJacob Faibussowitsch PetscCall(VecSetRandom(vecLeave0, rand));
679566063dSJacob Faibussowitsch PetscCall(VecCopy(vecRoot0, vecRoot1));
689566063dSJacob Faibussowitsch PetscCall(VecCopy(vecLeave0, vecLeave1));
699566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand));
70f210f596SVaclav Hapla }
71f210f596SVaclav Hapla
729566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(sf0, vecRoot0, vecLeave0, ADD_VALUES, SCATTER_FORWARD));
739566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(sf0, vecRoot0, vecLeave0, ADD_VALUES, SCATTER_FORWARD));
749566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(sf1, vecRoot1, vecLeave1, ADD_VALUES, SCATTER_FORWARD));
759566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(sf1, vecRoot1, vecLeave1, ADD_VALUES, SCATTER_FORWARD));
769566063dSJacob Faibussowitsch PetscCall(VecEqual(vecLeave0, vecLeave1, &flg));
77f210f596SVaclav Hapla PetscCheck(flg, comm, PETSC_ERR_PLIB, "leave vectors differ");
78f210f596SVaclav Hapla
799566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(sf0, vecLeave0, vecRoot0, ADD_VALUES, SCATTER_REVERSE));
809566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(sf0, vecLeave0, vecRoot0, ADD_VALUES, SCATTER_REVERSE));
819566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(sf1, vecLeave1, vecRoot1, ADD_VALUES, SCATTER_REVERSE));
829566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(sf1, vecLeave1, vecRoot1, ADD_VALUES, SCATTER_REVERSE));
839566063dSJacob Faibussowitsch PetscCall(VecEqual(vecRoot0, vecRoot1, &flg));
84f210f596SVaclav Hapla PetscCheck(flg, comm, PETSC_ERR_PLIB, "root vectors differ");
85f210f596SVaclav Hapla
869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vecRoot0));
879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vecRoot1));
889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vecLeave0));
899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vecLeave1));
903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
91f210f596SVaclav Hapla }
92f210f596SVaclav Hapla
PetscSFViewCustom(PetscSF sf,PetscViewer viewer)93c1730cc1SVaclav Hapla static PetscErrorCode PetscSFViewCustom(PetscSF sf, PetscViewer viewer)
94d71ae5a4SJacob Faibussowitsch {
956497c311SBarry Smith PetscMPIInt rank, nranks;
966497c311SBarry Smith PetscInt i, nroots, nleaves;
97c1730cc1SVaclav Hapla const PetscInt *ilocal;
98c1730cc1SVaclav Hapla const PetscSFNode *iremote;
99c1730cc1SVaclav Hapla PetscLayout rootLayout;
100c1730cc1SVaclav Hapla PetscInt *gremote;
101c1730cc1SVaclav Hapla
102c1730cc1SVaclav Hapla PetscFunctionBegin;
103c1730cc1SVaclav Hapla PetscCall(PetscSFSetUp(sf));
104c1730cc1SVaclav Hapla PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, &iremote));
105c1730cc1SVaclav Hapla PetscCall(PetscSFGetRootRanks(sf, &nranks, NULL, NULL, NULL, NULL));
106c1730cc1SVaclav Hapla PetscCall(PetscSFGetGraphLayout(sf, &rootLayout, NULL, NULL, &gremote));
107c1730cc1SVaclav Hapla PetscCheck(nroots == rootLayout->n, PetscObjectComm((PetscObject)sf), PETSC_ERR_PLIB, "Assertion failed: nroots == rootLayout->n");
108c1730cc1SVaclav Hapla PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)sf, viewer));
109c1730cc1SVaclav Hapla PetscCall(PetscViewerASCIIPushTab(viewer));
110c1730cc1SVaclav Hapla PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf), &rank));
111c1730cc1SVaclav Hapla PetscCall(PetscViewerASCIIPushSynchronized(viewer));
112c1730cc1SVaclav Hapla if (rank == 0) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "rank #leaves #roots\n"));
113c1730cc1SVaclav Hapla PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%2d] %7" PetscInt_FMT " %6" PetscInt_FMT "\n", rank, nleaves, nroots));
114c1730cc1SVaclav Hapla PetscCall(PetscViewerFlush(viewer));
115c1730cc1SVaclav Hapla if (rank == 0) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "leaves roots roots in global numbering\n"));
116c1730cc1SVaclav Hapla for (i = 0; i < nleaves; i++)
1176497c311SBarry Smith PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "(%2d, %2d) <- (%2d, %2" PetscInt_FMT ") = %3" PetscInt_FMT "\n", rank, (PetscMPIInt)(ilocal ? ilocal[i] : i), (PetscMPIInt)iremote[i].rank, iremote[i].index, gremote[i]));
118c1730cc1SVaclav Hapla PetscCall(PetscViewerFlush(viewer));
119c1730cc1SVaclav Hapla PetscCall(PetscViewerASCIIPopSynchronized(viewer));
120c1730cc1SVaclav Hapla PetscCall(PetscViewerASCIIPopTab(viewer));
121c1730cc1SVaclav Hapla PetscCall(PetscLayoutDestroy(&rootLayout));
122c1730cc1SVaclav Hapla PetscCall(PetscFree(gremote));
1233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
124c1730cc1SVaclav Hapla }
125c1730cc1SVaclav Hapla
CreateReferenceSF_Regular(AppCtx * ctx,PetscSF * refSF)126c1730cc1SVaclav Hapla PetscErrorCode CreateReferenceSF_Regular(AppCtx *ctx, PetscSF *refSF)
127c1730cc1SVaclav Hapla {
128c1730cc1SVaclav Hapla PetscInt j;
129f210f596SVaclav Hapla PetscInt *ilocal = NULL;
130c1730cc1SVaclav Hapla PetscInt nLeaves = ctx->nsfs * ctx->n * ctx->size;
131c1730cc1SVaclav Hapla PetscInt nroots = ctx->n * ctx->nsfs;
132f210f596SVaclav Hapla PetscSF sf;
133f210f596SVaclav Hapla
134f210f596SVaclav Hapla PetscFunctionBegin;
135f210f596SVaclav Hapla ilocal = NULL;
13648a46eb9SPierre Jolivet if (ctx->sparseLeaves) PetscCall(PetscCalloc1(nLeaves + 1, &ilocal));
1379566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(ctx->comm, &sf));
138c1730cc1SVaclav Hapla for (j = 0; j < nLeaves; j++) {
139c1730cc1SVaclav Hapla if (ctx->sparseLeaves) ilocal[j + 1] = ilocal[j] + ctx->leaveStep;
140c1730cc1SVaclav Hapla }
141c1730cc1SVaclav Hapla switch (ctx->rootMode) {
142c1730cc1SVaclav Hapla case PETSCSF_CONCATENATE_ROOTMODE_SHARED:
143c1730cc1SVaclav Hapla case PETSCSF_CONCATENATE_ROOTMODE_LOCAL: {
144c1730cc1SVaclav Hapla PetscInt i, k;
145c1730cc1SVaclav Hapla PetscMPIInt r;
146c1730cc1SVaclav Hapla PetscSFNode *iremote;
147c1730cc1SVaclav Hapla
148c1730cc1SVaclav Hapla PetscCall(PetscCalloc1(nLeaves, &iremote));
149f210f596SVaclav Hapla for (i = 0, j = 0; i < ctx->nsfs; i++) {
150f210f596SVaclav Hapla for (r = 0; r < ctx->size; r++) {
151c1730cc1SVaclav Hapla for (k = 0; k < ctx->n; k++, j++) {
152f210f596SVaclav Hapla iremote[j].rank = r;
153c1730cc1SVaclav Hapla iremote[j].index = k + i * ctx->n;
154f210f596SVaclav Hapla }
155f210f596SVaclav Hapla }
156f210f596SVaclav Hapla }
1579566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, nroots, nLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
158c1730cc1SVaclav Hapla } break;
159c1730cc1SVaclav Hapla case PETSCSF_CONCATENATE_ROOTMODE_GLOBAL: {
160c1730cc1SVaclav Hapla PetscLayout map = NULL;
161c1730cc1SVaclav Hapla PetscInt *gremote;
162c1730cc1SVaclav Hapla
163c1730cc1SVaclav Hapla PetscCall(PetscLayoutCreateFromSizes(ctx->comm, nroots, PETSC_DECIDE, 1, &map));
164c1730cc1SVaclav Hapla PetscCall(PetscMalloc1(nLeaves, &gremote));
165c1730cc1SVaclav Hapla for (j = 0; j < nLeaves; j++) gremote[j] = j;
166c1730cc1SVaclav Hapla PetscCall(PetscSFSetGraphLayout(sf, map, nLeaves, ilocal, PETSC_OWN_POINTER, gremote));
167c1730cc1SVaclav Hapla PetscCall(PetscFree(gremote));
168c1730cc1SVaclav Hapla PetscCall(PetscLayoutDestroy(&map));
169c1730cc1SVaclav Hapla } break;
170c1730cc1SVaclav Hapla default:
171c1730cc1SVaclav Hapla SETERRQ(ctx->comm, PETSC_ERR_SUP, "unsupported rootmode %d", ctx->rootMode);
172c1730cc1SVaclav Hapla }
173c1730cc1SVaclav Hapla PetscCall(PetscObjectSetName((PetscObject)sf, "reference_sf"));
174c1730cc1SVaclav Hapla if (ctx->viewer) PetscCall(PetscSFViewCustom(sf, ctx->viewer));
175f210f596SVaclav Hapla *refSF = sf;
1763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
177f210f596SVaclav Hapla }
178f210f596SVaclav Hapla
CreateSFs_Irregular(AppCtx * ctx,PetscSF * newSFs[],PetscInt * leafOffsets[])179c1730cc1SVaclav Hapla PetscErrorCode CreateSFs_Irregular(AppCtx *ctx, PetscSF *newSFs[], PetscInt *leafOffsets[])
180d71ae5a4SJacob Faibussowitsch {
181f210f596SVaclav Hapla PetscInt i;
182f210f596SVaclav Hapla PetscInt *lOffsets = NULL;
183f210f596SVaclav Hapla PetscSF *sfs;
184c1730cc1SVaclav Hapla PetscInt nLeaves = ctx->n * ctx->size + (ctx->size - 1) * ctx->size / 2;
185c1730cc1SVaclav Hapla PetscInt nroots = ctx->n + ctx->rank + ctx->nsfs - 1 + ctx->size - 1;
186f210f596SVaclav Hapla
187f210f596SVaclav Hapla PetscFunctionBegin;
18848a46eb9SPierre Jolivet if (ctx->sparseLeaves) PetscCall(PetscCalloc1(ctx->nsfs + 1, &lOffsets));
1899566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ctx->nsfs, &sfs));
190f210f596SVaclav Hapla for (i = 0; i < ctx->nsfs; i++) {
191c1730cc1SVaclav Hapla PetscSF sf;
192f210f596SVaclav Hapla PetscInt j, k;
193f210f596SVaclav Hapla PetscMPIInt r;
194f210f596SVaclav Hapla PetscInt *ilocal = NULL;
195f210f596SVaclav Hapla PetscSFNode *iremote;
196c1730cc1SVaclav Hapla char name[32];
197f210f596SVaclav Hapla
19848a46eb9SPierre Jolivet if (ctx->sparseLeaves) PetscCall(PetscCalloc1(nLeaves + 1, &ilocal));
1999566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nLeaves, &iremote));
200c1730cc1SVaclav Hapla for (r = ctx->size - 1, j = 0; r >= 0; r--) {
201c1730cc1SVaclav Hapla for (k = 0; k < ctx->n + r; k++, j++) {
202ad540459SPierre Jolivet if (ctx->sparseLeaves) ilocal[j + 1] = ilocal[j] + ctx->leaveStep;
203f210f596SVaclav Hapla iremote[j].rank = r;
204c1730cc1SVaclav Hapla iremote[j].index = k + i + ctx->rank;
205f210f596SVaclav Hapla }
206f210f596SVaclav Hapla }
207f210f596SVaclav Hapla if (ctx->sparseLeaves) lOffsets[i + 1] = lOffsets[i] + ilocal[j];
208f210f596SVaclav Hapla
209c1730cc1SVaclav Hapla PetscCall(PetscSFCreate(ctx->comm, &sf));
210c1730cc1SVaclav Hapla PetscCall(PetscSFSetGraph(sf, nroots, nLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
211c1730cc1SVaclav Hapla PetscCall(PetscSNPrintf(name, sizeof(name), "sf_%" PetscInt_FMT, i));
212c1730cc1SVaclav Hapla PetscCall(PetscObjectSetName((PetscObject)sf, name));
213c1730cc1SVaclav Hapla if (ctx->viewer) PetscCall(PetscSFViewCustom(sf, ctx->viewer));
214c1730cc1SVaclav Hapla sfs[i] = sf;
215c1730cc1SVaclav Hapla }
216c1730cc1SVaclav Hapla *newSFs = sfs;
217c1730cc1SVaclav Hapla *leafOffsets = lOffsets;
2183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
219c1730cc1SVaclav Hapla }
220c1730cc1SVaclav Hapla
CreateSFs_Regular(AppCtx * ctx,PetscSF * newSFs[],PetscInt * leafOffsets[])221c1730cc1SVaclav Hapla PetscErrorCode CreateSFs_Regular(AppCtx *ctx, PetscSF *newSFs[], PetscInt *leafOffsets[])
222c1730cc1SVaclav Hapla {
223c1730cc1SVaclav Hapla PetscInt i;
224c1730cc1SVaclav Hapla PetscInt *lOffsets = NULL;
225c1730cc1SVaclav Hapla PetscInt nLeaves = ctx->n * ctx->size;
226c1730cc1SVaclav Hapla PetscSF *sfs;
227c1730cc1SVaclav Hapla PetscSFConcatenateRootMode mode = ctx->compare ? ctx->rootMode : PETSCSF_CONCATENATE_ROOTMODE_LOCAL;
228c1730cc1SVaclav Hapla
229c1730cc1SVaclav Hapla PetscFunctionBegin;
230c1730cc1SVaclav Hapla if (ctx->sparseLeaves) PetscCall(PetscCalloc1(ctx->nsfs + 1, &lOffsets));
231c1730cc1SVaclav Hapla PetscCall(PetscCalloc1(ctx->nsfs, &sfs));
232c1730cc1SVaclav Hapla for (i = 0; i < ctx->nsfs; i++) {
233c1730cc1SVaclav Hapla PetscSF sf;
234c1730cc1SVaclav Hapla PetscInt j;
235c1730cc1SVaclav Hapla PetscInt *ilocal = NULL;
236c1730cc1SVaclav Hapla char name[32];
237c1730cc1SVaclav Hapla
238c1730cc1SVaclav Hapla PetscCall(PetscSFCreate(ctx->comm, &sf));
239c1730cc1SVaclav Hapla if (ctx->sparseLeaves) {
240c1730cc1SVaclav Hapla PetscCall(PetscCalloc1(nLeaves + 1, &ilocal));
241aa624791SPierre Jolivet for (j = 0; j < nLeaves; j++) ilocal[j + 1] = ilocal[j] + ctx->leaveStep;
242c1730cc1SVaclav Hapla lOffsets[i + 1] = lOffsets[i] + ilocal[nLeaves];
243c1730cc1SVaclav Hapla }
244c1730cc1SVaclav Hapla switch (mode) {
245c1730cc1SVaclav Hapla case PETSCSF_CONCATENATE_ROOTMODE_LOCAL: {
246c1730cc1SVaclav Hapla PetscInt k, nroots = ctx->n;
247c1730cc1SVaclav Hapla PetscMPIInt r;
248c1730cc1SVaclav Hapla PetscSFNode *iremote;
249c1730cc1SVaclav Hapla
250c1730cc1SVaclav Hapla PetscCall(PetscMalloc1(nLeaves, &iremote));
251c1730cc1SVaclav Hapla for (r = 0, j = 0; r < ctx->size; r++) {
252c1730cc1SVaclav Hapla for (k = 0; k < ctx->n; k++, j++) {
253c1730cc1SVaclav Hapla iremote[j].rank = r;
254c1730cc1SVaclav Hapla iremote[j].index = k;
255c1730cc1SVaclav Hapla }
256c1730cc1SVaclav Hapla }
257c1730cc1SVaclav Hapla PetscCall(PetscSFSetGraph(sf, nroots, nLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
258c1730cc1SVaclav Hapla } break;
259c1730cc1SVaclav Hapla case PETSCSF_CONCATENATE_ROOTMODE_SHARED: {
260c1730cc1SVaclav Hapla PetscInt k, nroots = ctx->n * ctx->nsfs;
261c1730cc1SVaclav Hapla PetscMPIInt r;
262c1730cc1SVaclav Hapla PetscSFNode *iremote;
263c1730cc1SVaclav Hapla
264c1730cc1SVaclav Hapla PetscCall(PetscMalloc1(nLeaves, &iremote));
265c1730cc1SVaclav Hapla for (r = 0, j = 0; r < ctx->size; r++) {
266c1730cc1SVaclav Hapla for (k = 0; k < ctx->n; k++, j++) {
267c1730cc1SVaclav Hapla iremote[j].rank = r;
268c1730cc1SVaclav Hapla iremote[j].index = k + i * ctx->n;
269c1730cc1SVaclav Hapla }
270c1730cc1SVaclav Hapla }
271c1730cc1SVaclav Hapla PetscCall(PetscSFSetGraph(sf, nroots, nLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
272c1730cc1SVaclav Hapla } break;
273c1730cc1SVaclav Hapla case PETSCSF_CONCATENATE_ROOTMODE_GLOBAL: {
274c1730cc1SVaclav Hapla PetscInt nroots = ctx->n;
275c1730cc1SVaclav Hapla PetscLayout map = NULL;
276c1730cc1SVaclav Hapla PetscInt *gremote;
277c1730cc1SVaclav Hapla
278c1730cc1SVaclav Hapla PetscCall(PetscLayoutCreateFromSizes(ctx->comm, nroots, PETSC_DECIDE, 1, &map));
279c1730cc1SVaclav Hapla PetscCall(PetscMalloc1(nLeaves, &gremote));
280c1730cc1SVaclav Hapla for (j = 0; j < nLeaves; j++) gremote[j] = j;
281c1730cc1SVaclav Hapla PetscCall(PetscSFSetGraphLayout(sf, map, nLeaves, ilocal, PETSC_OWN_POINTER, gremote));
282c1730cc1SVaclav Hapla PetscCall(PetscFree(gremote));
283c1730cc1SVaclav Hapla PetscCall(PetscLayoutDestroy(&map));
284c1730cc1SVaclav Hapla } break;
285c1730cc1SVaclav Hapla default:
286c1730cc1SVaclav Hapla SETERRQ(ctx->comm, PETSC_ERR_SUP, "unsupported rootmode %d", ctx->rootMode);
287c1730cc1SVaclav Hapla }
288c1730cc1SVaclav Hapla PetscCall(PetscSNPrintf(name, sizeof(name), "sf_%" PetscInt_FMT, i));
289c1730cc1SVaclav Hapla PetscCall(PetscObjectSetName((PetscObject)sf, name));
290c1730cc1SVaclav Hapla if (ctx->viewer) PetscCall(PetscSFViewCustom(sf, ctx->viewer));
291c1730cc1SVaclav Hapla sfs[i] = sf;
292f210f596SVaclav Hapla }
293f210f596SVaclav Hapla *newSFs = sfs;
294f210f596SVaclav Hapla *leafOffsets = lOffsets;
2953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
296f210f596SVaclav Hapla }
297f210f596SVaclav Hapla
DestroySFs(AppCtx * ctx,PetscSF * sfs[])298d71ae5a4SJacob Faibussowitsch PetscErrorCode DestroySFs(AppCtx *ctx, PetscSF *sfs[])
299d71ae5a4SJacob Faibussowitsch {
300f210f596SVaclav Hapla PetscInt i;
301f210f596SVaclav Hapla
302f210f596SVaclav Hapla PetscFunctionBegin;
30348a46eb9SPierre Jolivet for (i = 0; i < ctx->nsfs; i++) PetscCall(PetscSFDestroy(&(*sfs)[i]));
3049566063dSJacob Faibussowitsch PetscCall(PetscFree(*sfs));
3053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
306f210f596SVaclav Hapla }
307f210f596SVaclav Hapla
main(int argc,char ** argv)308d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
309d71ae5a4SJacob Faibussowitsch {
310c1730cc1SVaclav Hapla AppCtx ctx_;
311c1730cc1SVaclav Hapla AppCtx *ctx = &ctx_;
312c1730cc1SVaclav Hapla PetscSF sf;
313c1730cc1SVaclav Hapla PetscSF *sfs = NULL;
314c1730cc1SVaclav Hapla PetscInt *leafOffsets = NULL;
315f210f596SVaclav Hapla MPI_Comm comm;
316f210f596SVaclav Hapla
317327415f7SBarry Smith PetscFunctionBeginUser;
3189566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
319f210f596SVaclav Hapla comm = PETSC_COMM_WORLD;
320c1730cc1SVaclav Hapla PetscCall(GetOptions(comm, ctx));
321f210f596SVaclav Hapla
322c1730cc1SVaclav Hapla if (ctx->irregular) {
323c1730cc1SVaclav Hapla PetscCall(CreateSFs_Irregular(ctx, &sfs, &leafOffsets));
324c1730cc1SVaclav Hapla } else {
325c1730cc1SVaclav Hapla PetscCall(CreateSFs_Regular(ctx, &sfs, &leafOffsets));
326c1730cc1SVaclav Hapla }
327c1730cc1SVaclav Hapla PetscCall(PetscSFConcatenate(comm, ctx->nsfs, sfs, ctx->rootMode, leafOffsets, &sf));
328c1730cc1SVaclav Hapla PetscCall(PetscObjectSetName((PetscObject)sf, "result_sf"));
329c1730cc1SVaclav Hapla if (ctx->viewer) {
330c1730cc1SVaclav Hapla PetscCall(PetscPrintf(comm, "rootMode = %s:\n", PetscSFConcatenateRootModes[ctx->rootMode]));
331c1730cc1SVaclav Hapla PetscCall(PetscSFViewCustom(sf, ctx->viewer));
332c1730cc1SVaclav Hapla }
333c1730cc1SVaclav Hapla if (ctx->compare) {
334c1730cc1SVaclav Hapla PetscSF sfRef;
335c1730cc1SVaclav Hapla
336c1730cc1SVaclav Hapla PetscAssert(!ctx->irregular, comm, PETSC_ERR_SUP, "Combination -compare_to_reference true -irregular true not implemented");
337c1730cc1SVaclav Hapla PetscCall(CreateReferenceSF_Regular(ctx, &sfRef));
3389566063dSJacob Faibussowitsch PetscCall(PetscSFCheckEqual_Private(sf, sfRef));
339c1730cc1SVaclav Hapla PetscCall(PetscSFDestroy(&sfRef));
340c1730cc1SVaclav Hapla }
341c1730cc1SVaclav Hapla PetscCall(DestroySFs(ctx, &sfs));
3429566063dSJacob Faibussowitsch PetscCall(PetscFree(leafOffsets));
3439566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf));
344c1730cc1SVaclav Hapla if (ctx->viewer) {
345c1730cc1SVaclav Hapla PetscCall(PetscViewerPopFormat(ctx->viewer));
346648c30bcSBarry Smith PetscCall(PetscViewerDestroy(&ctx->viewer));
347c1730cc1SVaclav Hapla }
3489566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
349b122ec5aSJacob Faibussowitsch return 0;
350f210f596SVaclav Hapla }
351f210f596SVaclav Hapla
352f210f596SVaclav Hapla /*TEST
353f210f596SVaclav Hapla test:
354f210f596SVaclav Hapla nsize: {{1 3}}
355c1730cc1SVaclav Hapla args: -compare_to_reference -nsfs {{1 3}} -n {{0 1 5}} -leave_step {{1 3}} -root_mode {{local shared global}}
356*3886731fSPierre Jolivet output_file: output/empty.out
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