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