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