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