xref: /petsc/src/vec/is/sf/tests/ex5.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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;
19*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nl",&nl,NULL));
20*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-explicit_inverse",&inverse,NULL));
21*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
22*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
23c4762a1bSJed Brown 
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFCreate(PETSC_COMM_WORLD, &sfA));
25*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFCreate(PETSC_COMM_WORLD, &sfB));
26*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetFromOptions(sfA));
27*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetFromOptions(sfB));
28c4762a1bSJed Brown 
29c4762a1bSJed Brown   n = 4*nl*size;
30c4762a1bSJed Brown   m = 2*nl;
31c4762a1bSJed Brown   k = nl;
32c4762a1bSJed Brown 
33dd400576SPatrick Sanan   nldataA = rank == 0 ? n : 0;
34c4762a1bSJed Brown   nldataB = 3*nl;
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   nrootsA  = m;
37dd400576SPatrick Sanan   nleavesA = rank == 0 ? size*m : 0;
38dd400576SPatrick Sanan   nrootsB  = rank == 0 ? n : 0;
39c4762a1bSJed Brown   nleavesB = k;
40c4762a1bSJed Brown 
41*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nleavesA, &ilocalA));
42*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nleavesA, &iremoteA));
43*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nleavesB, &ilocalB));
44*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nleavesB, &iremoteB));
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   }
62*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetGraph(sfA, nrootsA, nleavesA, ilocalA, PETSC_OWN_POINTER, iremoteA, PETSC_OWN_POINTER));
63*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetGraph(sfB, nrootsB, nleavesB, ilocalB, PETSC_OWN_POINTER, iremoteB, PETSC_OWN_POINTER));
64*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetUp(sfA));
65*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetUp(sfB));
66*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)sfA, "sfA"));
67*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)sfB, "sfB"));
68*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFViewFromOptions(sfA, NULL, "-view"));
69*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFViewFromOptions(sfB, NULL, "-view"));
70c4762a1bSJed Brown 
71*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFGetLeafRange(sfA, NULL, &mA));
72*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFGetLeafRange(sfB, NULL, &mB));
73*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(nrootsA, &rdA, nldataA, &ldA));
74*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(nrootsB, &rdB, nldataB, &ldB));
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 
79*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "BcastB(BcastA)\n"));
80*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "A: root data\n"));
81*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscIntView(nrootsA, rdA, PETSC_VIEWER_STDOUT_WORLD));
82*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFBcastBegin(sfA, MPIU_INT, rdA, ldA,MPI_REPLACE));
83*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFBcastEnd(sfA, MPIU_INT, rdA, ldA,MPI_REPLACE));
84*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "A: leaf data (all)\n"));
85*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscIntView(nldataA, ldA, PETSC_VIEWER_STDOUT_WORLD));
86*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFBcastBegin(sfB, MPIU_INT, ldA, ldB,MPI_REPLACE));
87*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFBcastEnd(sfB, MPIU_INT, ldA, ldB,MPI_REPLACE));
88*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "B: leaf data (all)\n"));
89*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscIntView(nldataB, ldB, PETSC_VIEWER_STDOUT_WORLD));
90c4762a1bSJed Brown 
91*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFCompose(sfA, sfB, &sfBA));
92*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetFromOptions(sfBA));
93*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetUp(sfBA));
94*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)sfBA, "sfBA"));
95*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFViewFromOptions(sfBA, NULL, "-view"));
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   for (i = 0; i < nldataB; i++) ldB[i] = -1;
98*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "BcastBA\n"));
99*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "BA: root data\n"));
100*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscIntView(nrootsA, rdA, PETSC_VIEWER_STDOUT_WORLD));
101*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFBcastBegin(sfBA, MPIU_INT, rdA, ldB,MPI_REPLACE));
102*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFBcastEnd(sfBA, MPIU_INT, rdA, ldB,MPI_REPLACE));
103*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "BA: leaf data (all)\n"));
104*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscIntView(nldataB, ldB, PETSC_VIEWER_STDOUT_WORLD));
105c4762a1bSJed Brown 
106*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFCreateInverseSF(sfA, &sfAm));
107*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetFromOptions(sfAm));
108*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)sfAm, "sfAm"));
109*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFViewFromOptions(sfAm, NULL, "-view"));
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   if (!inverse) {
112*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFComposeInverse(sfA, sfA, &sfAAm));
113c4762a1bSJed Brown   } else {
114*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFCompose(sfA, sfAm, &sfAAm));
115c4762a1bSJed Brown   }
116*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetFromOptions(sfAAm));
117*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetUp(sfAAm));
118*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)sfAAm, "sfAAm"));
119*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFViewFromOptions(sfAAm, NULL, "-view"));
120c4762a1bSJed Brown 
121*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFCreateInverseSF(sfB, &sfBm));
122*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetFromOptions(sfBm));
123*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)sfBm, "sfBm"));
124*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFViewFromOptions(sfBm, NULL, "-view"));
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   if (!inverse) {
127*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFComposeInverse(sfB, sfB, &sfBBm));
128c4762a1bSJed Brown   } else {
129*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFCompose(sfB, sfBm, &sfBBm));
130c4762a1bSJed Brown   }
131*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetFromOptions(sfBBm));
132*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetUp(sfBBm));
133*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)sfBBm, "sfBBm"));
134*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFViewFromOptions(sfBBm, NULL, "-view"));
135c4762a1bSJed Brown 
136*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(rdA, ldA));
137*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(rdB, ldB));
138c4762a1bSJed Brown 
139*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFDestroy(&sfA));
140*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFDestroy(&sfB));
141*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFDestroy(&sfBA));
142*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFDestroy(&sfAm));
143*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFDestroy(&sfBm));
144*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFDestroy(&sfAAm));
145*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFDestroy(&sfBBm));
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}}
171dfd57a17SPierre 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
180dfd57a17SPierre 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