xref: /petsc/src/vec/is/sf/tests/ex4.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown static char help[]= "Test PetscSFFCompose when the ilocal array is not the identity\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscsf.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown int main(int argc, char **argv)
6c4762a1bSJed Brown {
7c4762a1bSJed Brown   PetscSF            sfA, sfB, sfBA;
8c4762a1bSJed Brown   PetscInt           nrootsA, nleavesA, nrootsB, nleavesB;
9c4762a1bSJed Brown   PetscInt          *ilocalA, *ilocalB;
10c4762a1bSJed Brown   PetscSFNode       *iremoteA, *iremoteB;
11c4762a1bSJed Brown   Vec                a, b, ba;
12c4762a1bSJed Brown   const PetscScalar *arrayR;
13c4762a1bSJed Brown   PetscScalar       *arrayW;
14c4762a1bSJed Brown   PetscMPIInt        size;
15c4762a1bSJed Brown   PetscInt           i;
16c4762a1bSJed Brown   PetscInt           maxleafB;
17c4762a1bSJed Brown   PetscBool          flag = PETSC_FALSE;
18c4762a1bSJed Brown 
19*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,NULL,help));
205f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
214643bc0bSVaclav Hapla   PetscCheck(size == 1,PETSC_COMM_WORLD, PETSC_ERR_USER, "Only coded for one MPI process");
22c4762a1bSJed Brown 
235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFCreate(PETSC_COMM_WORLD, &sfA));
245f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFCreate(PETSC_COMM_WORLD, &sfB));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetFromOptions(sfA));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetFromOptions(sfB));
27c4762a1bSJed Brown 
285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-sparse_sfB",&flag,NULL));
29c4762a1bSJed Brown 
30c4762a1bSJed Brown   if (flag) {
31c4762a1bSJed Brown     /* sfA permutes indices, sfB has sparse leaf space. */
32c4762a1bSJed Brown     nrootsA = 3;
33c4762a1bSJed Brown     nleavesA = 3;
34c4762a1bSJed Brown     nrootsB = 3;
35c4762a1bSJed Brown     nleavesB = 2;
36c4762a1bSJed Brown   } else {
37c4762a1bSJed Brown     /* sfA reverses indices, sfB is identity */
38c4762a1bSJed Brown     nrootsA = nrootsB = nleavesA = nleavesB = 4;
39c4762a1bSJed Brown   }
405f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nleavesA, &ilocalA));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nleavesA, &iremoteA));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nleavesB, &ilocalB));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nleavesB, &iremoteB));
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   for (i = 0; i < nleavesA; i++) {
46c4762a1bSJed Brown     iremoteA[i].rank = 0;
47c4762a1bSJed Brown     iremoteA[i].index = i;
48c4762a1bSJed Brown     if (flag) {
49c4762a1bSJed Brown       ilocalA[i] = (i + 1) % nleavesA;
50c4762a1bSJed Brown     } else {
51c4762a1bSJed Brown       ilocalA[i] = nleavesA - i - 1;
52c4762a1bSJed Brown     }
53c4762a1bSJed Brown   }
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   for (i = 0; i < nleavesB; i++) {
56c4762a1bSJed Brown     iremoteB[i].rank = 0;
57c4762a1bSJed Brown     if (flag) {
58c4762a1bSJed Brown       ilocalB[i] = nleavesB - i;
59c4762a1bSJed Brown       iremoteB[i].index = nleavesB - i - 1;
60c4762a1bSJed Brown     } else {
61c4762a1bSJed Brown       ilocalB[i] = i;
62c4762a1bSJed Brown       iremoteB[i].index = i;
63c4762a1bSJed Brown     }
64c4762a1bSJed Brown   }
65c4762a1bSJed Brown 
665f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetGraph(sfA, nrootsA, nleavesA, ilocalA, PETSC_OWN_POINTER, iremoteA, PETSC_OWN_POINTER));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetGraph(sfB, nrootsB, nleavesB, ilocalB, PETSC_OWN_POINTER, iremoteB, PETSC_OWN_POINTER));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetUp(sfA));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetUp(sfB));
705f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)sfA, "sfA"));
715f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)sfB, "sfB"));
72c4762a1bSJed Brown 
735f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_WORLD, nrootsA, &a));
745f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_WORLD, nleavesA, &b));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFGetLeafRange(sfB, NULL, &maxleafB));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_WORLD, maxleafB+1, &ba));
775f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(a, &arrayW));
78c4762a1bSJed Brown   for (i = 0; i < nrootsA; i++) {
79c4762a1bSJed Brown     arrayW[i] = (PetscScalar)i;
80c4762a1bSJed Brown   }
815f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(a, &arrayW));
82c4762a1bSJed Brown 
835f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Initial Vec A\n"));
845f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(a, NULL));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(a, &arrayR));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(b, &arrayW));
87c4762a1bSJed Brown 
885f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFBcastBegin(sfA, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFBcastEnd(sfA, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(b, &arrayW));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(a, &arrayR));
925f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "\nBroadcast A->B over sfA\n"));
935f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(b, NULL));
94c4762a1bSJed Brown 
955f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(b, &arrayR));
965f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(ba, &arrayW));
97c4762a1bSJed Brown   arrayW[0] = 10.0;             /* Not touched by bcast */
985f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFBcastBegin(sfB, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE));
995f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFBcastEnd(sfB, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE));
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(b, &arrayR));
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(ba, &arrayW));
102c4762a1bSJed Brown 
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "\nBroadcast B->BA over sfB\n"));
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(ba, NULL));
105c4762a1bSJed Brown 
1065f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFCompose(sfA, sfB, &sfBA));
1075f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetFromOptions(sfBA));
1085f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)sfBA, "(sfB o sfA)"));
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(a, &arrayR));
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(ba, &arrayW));
111c4762a1bSJed Brown   arrayW[0] = 11.0;             /* Not touched by bcast */
1125f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFBcastBegin(sfBA, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE));
1135f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFBcastEnd(sfBA, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE));
1145f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(ba, &arrayW));
1155f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(a, &arrayR));
1165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "\nBroadcast A->BA over sfBA (sfB o sfA)\n"));
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(ba, NULL));
118c4762a1bSJed Brown 
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&ba));
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&b));
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&a));
122c4762a1bSJed Brown 
1235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFView(sfA, NULL));
1245f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFView(sfB, NULL));
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFView(sfBA, NULL));
1265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFDestroy(&sfA));
1275f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFDestroy(&sfB));
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFDestroy(&sfBA));
129c4762a1bSJed Brown 
130*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
131*b122ec5aSJacob Faibussowitsch   return 0;
132c4762a1bSJed Brown }
133c4762a1bSJed Brown 
134c4762a1bSJed Brown /*TEST
135c4762a1bSJed Brown 
136c4762a1bSJed Brown    test:
137c4762a1bSJed Brown      suffix: 1
138c4762a1bSJed Brown 
139c4762a1bSJed Brown    test:
140c4762a1bSJed Brown      suffix: 2
141c4762a1bSJed Brown      filter: grep -v "type" | grep -v "sort"
142c4762a1bSJed Brown      args: -sparse_sfB
143c4762a1bSJed Brown 
144c4762a1bSJed Brown    test:
145c4762a1bSJed Brown      suffix: 2_window
146c4762a1bSJed Brown      filter: grep -v "type" | grep -v "sort"
147c4762a1bSJed Brown      output_file: output/ex4_2.out
148c4762a1bSJed Brown      args: -sparse_sfB -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
149dfd57a17SPierre Jolivet      requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
150c4762a1bSJed Brown 
151c4762a1bSJed Brown    # The nightly test suite with MPICH uses ch3:sock, which is broken when winsize == 0 in some of the processes
152c4762a1bSJed Brown    test:
153c4762a1bSJed Brown      suffix: 2_window_shared
154c4762a1bSJed Brown      filter: grep -v "type" | grep -v "sort"
155c4762a1bSJed Brown      output_file: output/ex4_2.out
156c4762a1bSJed Brown      args: -sparse_sfB -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor shared
157dfd57a17SPierre Jolivet      requires: defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) !defined(PETSC_HAVE_MPICH_NUMVERSION) defined(PETSC_HAVE_MPI_ONE_SIDED)
158c4762a1bSJed Brown 
159c4762a1bSJed Brown TEST*/
160