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