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