xref: /petsc/src/vec/is/sf/tests/ex19.c (revision d9acb416d05abeed0a33bde3a81aeb2ea0364f6a)
1 static char help[] = "Test leaf sorting in PetscSFSetGraph()\n\n";
2 
3 #include <petscsf.h>
4 
5 typedef struct {
6   MPI_Comm      comm;
7   PetscMPIInt   rank, size;
8   PetscInt      leaveStep, nLeavesPerRank;
9   PetscBool     contiguousLeaves;
10   PetscCopyMode localmode, remotemode;
11   PetscInt     *ilocal;
12   PetscSFNode  *iremote;
13 } AppCtx;
14 
15 static PetscErrorCode GetOptions(MPI_Comm comm, AppCtx *ctx)
16 {
17   PetscFunctionBegin;
18   ctx->comm             = comm;
19   ctx->nLeavesPerRank   = 4;
20   ctx->leaveStep        = 1;
21   ctx->contiguousLeaves = PETSC_FALSE;
22   ctx->localmode        = PETSC_OWN_POINTER;
23   ctx->remotemode       = PETSC_OWN_POINTER;
24   ctx->ilocal           = NULL;
25   ctx->iremote          = NULL;
26   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n_leaves_per_rank", &ctx->nLeavesPerRank, NULL));
27   PetscCall(PetscOptionsGetInt(NULL, NULL, "-leave_step", &ctx->leaveStep, NULL));
28   PetscCall(PetscOptionsGetEnum(NULL, NULL, "-localmode", PetscCopyModes, (PetscEnum *)&ctx->localmode, NULL));
29   PetscCall(PetscOptionsGetEnum(NULL, NULL, "-remotemode", PetscCopyModes, (PetscEnum *)&ctx->remotemode, NULL));
30   ctx->contiguousLeaves = (PetscBool)(ctx->leaveStep == 1);
31   PetscCallMPI(MPI_Comm_size(comm, &ctx->size));
32   PetscCallMPI(MPI_Comm_rank(comm, &ctx->rank));
33   PetscFunctionReturn(PETSC_SUCCESS);
34 }
35 
36 static PetscErrorCode PetscSFCheckEqual_Private(PetscSF sf0, PetscSF sf1)
37 {
38   PetscInt  nRoot, nLeave;
39   Vec       vecRoot0, vecLeave0, vecRoot1, vecLeave1;
40   MPI_Comm  comm;
41   PetscBool flg;
42 
43   PetscFunctionBegin;
44   PetscCall(PetscObjectGetComm((PetscObject)sf0, &comm));
45   PetscCall(PetscSFGetGraph(sf0, &nRoot, NULL, NULL, NULL));
46   PetscCall(PetscSFGetLeafRange(sf0, NULL, &nLeave));
47   nLeave++;
48   PetscCall(VecCreateFromOptions(comm, NULL, 1, nRoot, PETSC_DECIDE, &vecRoot0));
49   PetscCall(VecCreateFromOptions(comm, NULL, 1, nLeave, PETSC_DECIDE, &vecLeave0));
50   PetscCall(VecDuplicate(vecRoot0, &vecRoot1));
51   PetscCall(VecDuplicate(vecLeave0, &vecLeave1));
52   {
53     PetscRandom rand;
54 
55     PetscCall(PetscRandomCreate(comm, &rand));
56     PetscCall(PetscRandomSetFromOptions(rand));
57     PetscCall(VecSetRandom(vecRoot0, rand));
58     PetscCall(VecSetRandom(vecLeave0, rand));
59     PetscCall(VecCopy(vecRoot0, vecRoot1));
60     PetscCall(VecCopy(vecLeave0, vecLeave1));
61     PetscCall(PetscRandomDestroy(&rand));
62   }
63 
64   PetscCall(VecScatterBegin(sf0, vecRoot0, vecLeave0, ADD_VALUES, SCATTER_FORWARD));
65   PetscCall(VecScatterEnd(sf0, vecRoot0, vecLeave0, ADD_VALUES, SCATTER_FORWARD));
66   PetscCall(VecScatterBegin(sf1, vecRoot1, vecLeave1, ADD_VALUES, SCATTER_FORWARD));
67   PetscCall(VecScatterEnd(sf1, vecRoot1, vecLeave1, ADD_VALUES, SCATTER_FORWARD));
68   PetscCall(VecEqual(vecLeave0, vecLeave1, &flg));
69   PetscCheck(flg, comm, PETSC_ERR_PLIB, "leave vectors differ");
70 
71   PetscCall(VecScatterBegin(sf0, vecLeave0, vecRoot0, ADD_VALUES, SCATTER_REVERSE));
72   PetscCall(VecScatterEnd(sf0, vecLeave0, vecRoot0, ADD_VALUES, SCATTER_REVERSE));
73   PetscCall(VecScatterBegin(sf1, vecLeave1, vecRoot1, ADD_VALUES, SCATTER_REVERSE));
74   PetscCall(VecScatterEnd(sf1, vecLeave1, vecRoot1, ADD_VALUES, SCATTER_REVERSE));
75   PetscCall(VecEqual(vecRoot0, vecRoot1, &flg));
76   PetscCheck(flg, comm, PETSC_ERR_PLIB, "root vectors differ");
77 
78   PetscCall(VecDestroy(&vecRoot0));
79   PetscCall(VecDestroy(&vecRoot1));
80   PetscCall(VecDestroy(&vecLeave0));
81   PetscCall(VecDestroy(&vecLeave1));
82   PetscFunctionReturn(PETSC_SUCCESS);
83 }
84 
85 PetscErrorCode CreateSF0(AppCtx *ctx, PetscSF *sf0)
86 {
87   PetscInt     j, k, r;
88   PetscInt     nLeaves = ctx->nLeavesPerRank * ctx->size;
89   PetscInt     nroots  = ctx->nLeavesPerRank;
90   PetscSF      sf;
91   PetscInt    *ilocal;
92   PetscSFNode *iremote;
93 
94   PetscFunctionBegin;
95   PetscCall(PetscMalloc1(nLeaves + 1, &ctx->ilocal));
96   PetscCall(PetscMalloc1(nLeaves, &ctx->iremote));
97   ilocal          = ctx->ilocal;
98   iremote         = ctx->iremote;
99   ilocal[nLeaves] = -ctx->leaveStep;
100   PetscCall(PetscSFCreate(ctx->comm, &sf));
101   for (r = 0, j = nLeaves - 1; r < ctx->size; r++) {
102     for (k = 0; k < ctx->nLeavesPerRank; k++, j--) {
103       ilocal[j]        = ilocal[j + 1] + ctx->leaveStep;
104       iremote[j].rank  = r;
105       iremote[j].index = k;
106     }
107   }
108   PetscCall(PetscSFSetGraph(sf, nroots, nLeaves, ilocal, ctx->localmode, iremote, ctx->remotemode));
109   {
110     const PetscInt *tlocal;
111     PetscBool       sorted;
112 
113     PetscCall(PetscSFGetGraph(sf, NULL, NULL, &tlocal, NULL));
114     PetscCheck(!ctx->contiguousLeaves || !tlocal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ilocal=NULL expected for contiguous case");
115     if (tlocal) {
116       PetscCall(PetscSortedInt(nLeaves, tlocal, &sorted));
117       PetscCheck(sorted, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ilocal expected to be sorted");
118     }
119   }
120   *sf0 = sf;
121   PetscFunctionReturn(PETSC_SUCCESS);
122 }
123 
124 PetscErrorCode CreateSF1(AppCtx *ctx, PetscSF *sf1)
125 {
126   PetscInt     j, k, r;
127   PetscInt    *ilocal = NULL;
128   PetscSFNode *iremote;
129   PetscInt     nLeaves = ctx->nLeavesPerRank * ctx->size;
130   PetscInt     nroots  = ctx->nLeavesPerRank;
131   PetscSF      sf;
132 
133   PetscFunctionBegin;
134   ilocal = NULL;
135   if (!ctx->contiguousLeaves) PetscCall(PetscCalloc1(nLeaves + 1, &ilocal));
136   PetscCall(PetscMalloc1(nLeaves, &iremote));
137   PetscCall(PetscSFCreate(ctx->comm, &sf));
138   for (r = 0, j = 0; r < ctx->size; r++) {
139     for (k = 0; k < ctx->nLeavesPerRank; k++, j++) {
140       if (!ctx->contiguousLeaves) ilocal[j + 1] = ilocal[j] + ctx->leaveStep;
141       iremote[j].rank  = r;
142       iremote[j].index = k;
143     }
144   }
145   PetscCheck(j == nLeaves, PETSC_COMM_SELF, PETSC_ERR_PLIB, "j != nLeaves");
146   PetscCall(PetscSFSetGraph(sf, nroots, nLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
147   if (ctx->contiguousLeaves) {
148     const PetscInt *tlocal;
149 
150     PetscCall(PetscSFGetGraph(sf, NULL, NULL, &tlocal, NULL));
151     PetscCheck(!tlocal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ilocal=NULL expected for contiguous case");
152   }
153   *sf1 = sf;
154   PetscFunctionReturn(PETSC_SUCCESS);
155 }
156 
157 int main(int argc, char **argv)
158 {
159   AppCtx   ctx;
160   PetscSF  sf0, sf1;
161   MPI_Comm comm;
162 
163   PetscFunctionBeginUser;
164   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
165   comm = PETSC_COMM_WORLD;
166   PetscCall(GetOptions(comm, &ctx));
167 
168   PetscCall(CreateSF0(&ctx, &sf0));
169   PetscCall(CreateSF1(&ctx, &sf1));
170   PetscCall(PetscSFViewFromOptions(sf0, NULL, "-sf0_view"));
171   PetscCall(PetscSFViewFromOptions(sf1, NULL, "-sf1_view"));
172   PetscCall(PetscSFCheckEqual_Private(sf0, sf1));
173 
174   if (ctx.localmode != PETSC_OWN_POINTER) PetscCall(PetscFree(ctx.ilocal));
175   if (ctx.remotemode != PETSC_OWN_POINTER) PetscCall(PetscFree(ctx.iremote));
176   PetscCall(PetscSFDestroy(&sf0));
177   PetscCall(PetscSFDestroy(&sf1));
178   PetscCall(PetscFinalize());
179   return 0;
180 }
181 
182 /*TEST
183   testset:
184     suffix: 1
185     nsize: {{1 3}}
186     args: -n_leaves_per_rank {{0 5}} -leave_step {{1 3}}
187     test:
188       suffix: a
189       args: -localmode {{COPY_VALUES OWN_POINTER}} -remotemode {{COPY_VALUES OWN_POINTER}}
190     test:
191       suffix: b
192       args: -localmode USE_POINTER -remotemode {{COPY_VALUES OWN_POINTER USE_POINTER}}
193 TEST*/
194