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