xref: /petsc/src/vec/is/sf/tests/ex5.c (revision 9566063d113dddea24716c546802770db7481bc0)
1c4762a1bSJed Brown static char help[]= "Test PetscSFFCompose when the ilocal arrays are not identity nor dense\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petsc.h>
4c4762a1bSJed Brown #include <petscsf.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc, char **argv)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   PetscSF      sfA, sfB, sfBA, sfAAm, sfBBm, sfAm, sfBm;
9c4762a1bSJed Brown   PetscInt     nrootsA, nleavesA, nrootsB, nleavesB;
10c4762a1bSJed Brown   PetscInt    *ilocalA, *ilocalB;
11c4762a1bSJed Brown   PetscSFNode *iremoteA, *iremoteB;
12c4762a1bSJed Brown   PetscMPIInt  rank,size;
13c4762a1bSJed Brown   PetscInt     i,m,n,k,nl = 2,mA,mB,nldataA,nldataB;
14c4762a1bSJed Brown   PetscInt    *rdA,*rdB,*ldA,*ldB;
15c4762a1bSJed Brown   PetscBool    inverse    = PETSC_FALSE;
16c4762a1bSJed Brown 
17*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
18*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-nl",&nl,NULL));
19*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-explicit_inverse",&inverse,NULL));
20*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
21*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
22c4762a1bSJed Brown 
23*9566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(PETSC_COMM_WORLD, &sfA));
24*9566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(PETSC_COMM_WORLD, &sfB));
25*9566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sfA));
26*9566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sfB));
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   n = 4*nl*size;
29c4762a1bSJed Brown   m = 2*nl;
30c4762a1bSJed Brown   k = nl;
31c4762a1bSJed Brown 
32dd400576SPatrick Sanan   nldataA = rank == 0 ? n : 0;
33c4762a1bSJed Brown   nldataB = 3*nl;
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   nrootsA  = m;
36dd400576SPatrick Sanan   nleavesA = rank == 0 ? size*m : 0;
37dd400576SPatrick Sanan   nrootsB  = rank == 0 ? n : 0;
38c4762a1bSJed Brown   nleavesB = k;
39c4762a1bSJed Brown 
40*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleavesA, &ilocalA));
41*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleavesA, &iremoteA));
42*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleavesB, &ilocalB));
43*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleavesB, &iremoteB));
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   /* sf A bcast is equivalent to a sparse gather on process 0
46c4762a1bSJed Brown      process 0 receives data in the middle [nl,3*nl] of the leaf data array for A */
47c4762a1bSJed Brown   for (i = 0; i < nleavesA; i++) {
48c4762a1bSJed Brown     iremoteA[i].rank = i/m;
49c4762a1bSJed Brown     iremoteA[i].index = i%m;
50c4762a1bSJed Brown     ilocalA[i] = nl + i/m * 4*nl + i%m;
51c4762a1bSJed Brown   }
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   /* sf B bcast is equivalent to a sparse scatter from process 0
54c4762a1bSJed Brown      process 0 sends data from [nl,2*nl] of the leaf data array for A
55c4762a1bSJed Brown      each process receives, in reverse order, in the middle [nl,2*nl] of the leaf data array for B */
56c4762a1bSJed Brown   for (i = 0; i < nleavesB; i++) {
57c4762a1bSJed Brown     iremoteB[i].rank = 0;
58c4762a1bSJed Brown     iremoteB[i].index = rank * 4*nl + nl + i%m;
59c4762a1bSJed Brown     ilocalB[i] = 2*nl - i - 1;
60c4762a1bSJed Brown   }
61*9566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(sfA, nrootsA, nleavesA, ilocalA, PETSC_OWN_POINTER, iremoteA, PETSC_OWN_POINTER));
62*9566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(sfB, nrootsB, nleavesB, ilocalB, PETSC_OWN_POINTER, iremoteB, PETSC_OWN_POINTER));
63*9566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sfA));
64*9566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sfB));
65*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)sfA, "sfA"));
66*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)sfB, "sfB"));
67*9566063dSJacob Faibussowitsch   PetscCall(PetscSFViewFromOptions(sfA, NULL, "-view"));
68*9566063dSJacob Faibussowitsch   PetscCall(PetscSFViewFromOptions(sfB, NULL, "-view"));
69c4762a1bSJed Brown 
70*9566063dSJacob Faibussowitsch   PetscCall(PetscSFGetLeafRange(sfA, NULL, &mA));
71*9566063dSJacob Faibussowitsch   PetscCall(PetscSFGetLeafRange(sfB, NULL, &mB));
72*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nrootsA, &rdA, nldataA, &ldA));
73*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nrootsB, &rdB, nldataB, &ldB));
74c4762a1bSJed Brown   for (i = 0; i < nrootsA; i++) rdA[i] = m*rank + i;
75c4762a1bSJed Brown   for (i = 0; i < nldataA; i++) ldA[i] = -1;
76c4762a1bSJed Brown   for (i = 0; i < nldataB; i++) ldB[i] = -1;
77c4762a1bSJed Brown 
78*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "BcastB(BcastA)\n"));
79*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "A: root data\n"));
80*9566063dSJacob Faibussowitsch   PetscCall(PetscIntView(nrootsA, rdA, PETSC_VIEWER_STDOUT_WORLD));
81*9566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sfA, MPIU_INT, rdA, ldA,MPI_REPLACE));
82*9566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(sfA, MPIU_INT, rdA, ldA,MPI_REPLACE));
83*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "A: leaf data (all)\n"));
84*9566063dSJacob Faibussowitsch   PetscCall(PetscIntView(nldataA, ldA, PETSC_VIEWER_STDOUT_WORLD));
85*9566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sfB, MPIU_INT, ldA, ldB,MPI_REPLACE));
86*9566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(sfB, MPIU_INT, ldA, ldB,MPI_REPLACE));
87*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "B: leaf data (all)\n"));
88*9566063dSJacob Faibussowitsch   PetscCall(PetscIntView(nldataB, ldB, PETSC_VIEWER_STDOUT_WORLD));
89c4762a1bSJed Brown 
90*9566063dSJacob Faibussowitsch   PetscCall(PetscSFCompose(sfA, sfB, &sfBA));
91*9566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sfBA));
92*9566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sfBA));
93*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)sfBA, "sfBA"));
94*9566063dSJacob Faibussowitsch   PetscCall(PetscSFViewFromOptions(sfBA, NULL, "-view"));
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   for (i = 0; i < nldataB; i++) ldB[i] = -1;
97*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "BcastBA\n"));
98*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "BA: root data\n"));
99*9566063dSJacob Faibussowitsch   PetscCall(PetscIntView(nrootsA, rdA, PETSC_VIEWER_STDOUT_WORLD));
100*9566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sfBA, MPIU_INT, rdA, ldB,MPI_REPLACE));
101*9566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(sfBA, MPIU_INT, rdA, ldB,MPI_REPLACE));
102*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "BA: leaf data (all)\n"));
103*9566063dSJacob Faibussowitsch   PetscCall(PetscIntView(nldataB, ldB, PETSC_VIEWER_STDOUT_WORLD));
104c4762a1bSJed Brown 
105*9566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateInverseSF(sfA, &sfAm));
106*9566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sfAm));
107*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)sfAm, "sfAm"));
108*9566063dSJacob Faibussowitsch   PetscCall(PetscSFViewFromOptions(sfAm, NULL, "-view"));
109c4762a1bSJed Brown 
110c4762a1bSJed Brown   if (!inverse) {
111*9566063dSJacob Faibussowitsch     PetscCall(PetscSFComposeInverse(sfA, sfA, &sfAAm));
112c4762a1bSJed Brown   } else {
113*9566063dSJacob Faibussowitsch     PetscCall(PetscSFCompose(sfA, sfAm, &sfAAm));
114c4762a1bSJed Brown   }
115*9566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sfAAm));
116*9566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sfAAm));
117*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)sfAAm, "sfAAm"));
118*9566063dSJacob Faibussowitsch   PetscCall(PetscSFViewFromOptions(sfAAm, NULL, "-view"));
119c4762a1bSJed Brown 
120*9566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateInverseSF(sfB, &sfBm));
121*9566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sfBm));
122*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)sfBm, "sfBm"));
123*9566063dSJacob Faibussowitsch   PetscCall(PetscSFViewFromOptions(sfBm, NULL, "-view"));
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   if (!inverse) {
126*9566063dSJacob Faibussowitsch     PetscCall(PetscSFComposeInverse(sfB, sfB, &sfBBm));
127c4762a1bSJed Brown   } else {
128*9566063dSJacob Faibussowitsch     PetscCall(PetscSFCompose(sfB, sfBm, &sfBBm));
129c4762a1bSJed Brown   }
130*9566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sfBBm));
131*9566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sfBBm));
132*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)sfBBm, "sfBBm"));
133*9566063dSJacob Faibussowitsch   PetscCall(PetscSFViewFromOptions(sfBBm, NULL, "-view"));
134c4762a1bSJed Brown 
135*9566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rdA, ldA));
136*9566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rdB, ldB));
137c4762a1bSJed Brown 
138*9566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfA));
139*9566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfB));
140*9566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfBA));
141*9566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfAm));
142*9566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfBm));
143*9566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfAAm));
144*9566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfBBm));
145c4762a1bSJed Brown 
146*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
147b122ec5aSJacob Faibussowitsch   return 0;
148c4762a1bSJed Brown }
149c4762a1bSJed Brown 
150c4762a1bSJed Brown /*TEST
151c4762a1bSJed Brown 
152c4762a1bSJed Brown    test:
153c4762a1bSJed Brown      suffix: 1
154c4762a1bSJed Brown      args: -view -explicit_inverse {{0 1}}
155c4762a1bSJed Brown 
156c4762a1bSJed Brown    test:
157c4762a1bSJed Brown      nsize: 7
158c4762a1bSJed Brown      filter: grep -v "type" | grep -v "sort"
159c4762a1bSJed Brown      suffix: 2
160c4762a1bSJed Brown      args: -view -nl 5 -explicit_inverse {{0 1}}
161c4762a1bSJed Brown 
162c4762a1bSJed Brown    # we cannot test for -sf_window_flavor dynamic because SFCompose with sparse leaves may change the root data pointer only locally, and this is not supported by the dynamic case
163c4762a1bSJed Brown    test:
164c4762a1bSJed Brown      nsize: 7
165c4762a1bSJed Brown      suffix: 2_window
166c4762a1bSJed Brown      filter: grep -v "type" | grep -v "sort"
167c4762a1bSJed Brown      output_file: output/ex5_2.out
168c4762a1bSJed Brown      args: -view -nl 5 -explicit_inverse {{0 1}} -sf_type window -sf_window_sync {{fence lock active}} -sf_window_flavor {{create allocate}}
169dfd57a17SPierre Jolivet      requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
170c4762a1bSJed Brown 
171c4762a1bSJed Brown    # The nightly test suite with MPICH uses ch3:sock, which is broken when winsize == 0 in some of the processes
172c4762a1bSJed Brown    test:
173c4762a1bSJed Brown      nsize: 7
174c4762a1bSJed Brown      suffix: 2_window_shared
175c4762a1bSJed Brown      filter: grep -v "type" | grep -v "sort"
176c4762a1bSJed Brown      output_file: output/ex5_2.out
177c4762a1bSJed Brown      args: -view -nl 5 -explicit_inverse {{0 1}} -sf_type window -sf_window_sync {{fence lock active}} -sf_window_flavor shared
178dfd57a17SPierre Jolivet      requires: defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) !defined(PETSC_HAVE_MPICH_NUMVERSION) defined(PETSC_HAVE_MPI_ONE_SIDED)
179c4762a1bSJed Brown 
180c4762a1bSJed Brown TEST*/
181