1*14357523SPierre Jolivet static char help[] = "Test PetscSFCompose() when the ilocal array is not the identity\n\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown #include <petscsf.h>
4c4762a1bSJed Brown
main(int argc,char ** argv)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
6d71ae5a4SJacob Faibussowitsch {
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
19327415f7SBarry Smith PetscFunctionBeginUser;
209566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
22bd158744SPierre Jolivet PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only coded for one MPI process");
23c4762a1bSJed Brown
249566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(PETSC_COMM_WORLD, &sfA));
259566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(PETSC_COMM_WORLD, &sfB));
269566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(sfA));
279566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(sfB));
28c4762a1bSJed Brown
299566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-sparse_sfB", &flag, NULL));
30c4762a1bSJed Brown
31c4762a1bSJed Brown if (flag) {
32c4762a1bSJed Brown /* sfA permutes indices, sfB has sparse leaf space. */
33c4762a1bSJed Brown nrootsA = 3;
34c4762a1bSJed Brown nleavesA = 3;
35c4762a1bSJed Brown nrootsB = 3;
36c4762a1bSJed Brown nleavesB = 2;
37c4762a1bSJed Brown } else {
38c4762a1bSJed Brown /* sfA reverses indices, sfB is identity */
39c4762a1bSJed Brown nrootsA = nrootsB = nleavesA = nleavesB = 4;
40c4762a1bSJed Brown }
419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleavesA, &ilocalA));
429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleavesA, &iremoteA));
439566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleavesB, &ilocalB));
449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleavesB, &iremoteB));
45c4762a1bSJed Brown
46c4762a1bSJed Brown for (i = 0; i < nleavesA; i++) {
47c4762a1bSJed Brown iremoteA[i].rank = 0;
48c4762a1bSJed Brown iremoteA[i].index = i;
49c4762a1bSJed Brown if (flag) {
50c4762a1bSJed Brown ilocalA[i] = (i + 1) % nleavesA;
51c4762a1bSJed Brown } else {
52c4762a1bSJed Brown ilocalA[i] = nleavesA - i - 1;
53c4762a1bSJed Brown }
54c4762a1bSJed Brown }
55c4762a1bSJed Brown
56c4762a1bSJed Brown for (i = 0; i < nleavesB; i++) {
57c4762a1bSJed Brown iremoteB[i].rank = 0;
58c4762a1bSJed Brown if (flag) {
59c4762a1bSJed Brown ilocalB[i] = nleavesB - i;
60c4762a1bSJed Brown iremoteB[i].index = nleavesB - i - 1;
61c4762a1bSJed Brown } else {
62c4762a1bSJed Brown ilocalB[i] = i;
63c4762a1bSJed Brown iremoteB[i].index = i;
64c4762a1bSJed Brown }
65c4762a1bSJed Brown }
66c4762a1bSJed Brown
679566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sfA, nrootsA, nleavesA, ilocalA, PETSC_OWN_POINTER, iremoteA, PETSC_OWN_POINTER));
689566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sfB, nrootsB, nleavesB, ilocalB, PETSC_OWN_POINTER, iremoteB, PETSC_OWN_POINTER));
699566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(sfA));
709566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(sfB));
719566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)sfA, "sfA"));
729566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)sfB, "sfB"));
73c4762a1bSJed Brown
749566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_WORLD, nrootsA, &a));
759566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_WORLD, nleavesA, &b));
769566063dSJacob Faibussowitsch PetscCall(PetscSFGetLeafRange(sfB, NULL, &maxleafB));
779566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_WORLD, maxleafB + 1, &ba));
789566063dSJacob Faibussowitsch PetscCall(VecGetArray(a, &arrayW));
79ad540459SPierre Jolivet for (i = 0; i < nrootsA; i++) arrayW[i] = (PetscScalar)i;
809566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(a, &arrayW));
81c4762a1bSJed Brown
829566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial Vec A\n"));
839566063dSJacob Faibussowitsch PetscCall(VecView(a, NULL));
849566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(a, &arrayR));
859566063dSJacob Faibussowitsch PetscCall(VecGetArray(b, &arrayW));
86c4762a1bSJed Brown
879566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sfA, MPIU_SCALAR, arrayR, arrayW, MPI_REPLACE));
889566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sfA, MPIU_SCALAR, arrayR, arrayW, MPI_REPLACE));
899566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(b, &arrayW));
909566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(a, &arrayR));
919566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nBroadcast A->B over sfA\n"));
929566063dSJacob Faibussowitsch PetscCall(VecView(b, NULL));
93c4762a1bSJed Brown
949566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(b, &arrayR));
959566063dSJacob Faibussowitsch PetscCall(VecGetArray(ba, &arrayW));
96c4762a1bSJed Brown arrayW[0] = 10.0; /* Not touched by bcast */
979566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sfB, MPIU_SCALAR, arrayR, arrayW, MPI_REPLACE));
989566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sfB, MPIU_SCALAR, arrayR, arrayW, MPI_REPLACE));
999566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(b, &arrayR));
1009566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ba, &arrayW));
101c4762a1bSJed Brown
1029566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nBroadcast B->BA over sfB\n"));
1039566063dSJacob Faibussowitsch PetscCall(VecView(ba, NULL));
104c4762a1bSJed Brown
1059566063dSJacob Faibussowitsch PetscCall(PetscSFCompose(sfA, sfB, &sfBA));
1069566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(sfBA));
1079566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)sfBA, "(sfB o sfA)"));
1089566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(a, &arrayR));
1099566063dSJacob Faibussowitsch PetscCall(VecGetArray(ba, &arrayW));
110c4762a1bSJed Brown arrayW[0] = 11.0; /* Not touched by bcast */
1119566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sfBA, MPIU_SCALAR, arrayR, arrayW, MPI_REPLACE));
1129566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sfBA, MPIU_SCALAR, arrayR, arrayW, MPI_REPLACE));
1139566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ba, &arrayW));
1149566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(a, &arrayR));
1159566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nBroadcast A->BA over sfBA (sfB o sfA)\n"));
1169566063dSJacob Faibussowitsch PetscCall(VecView(ba, NULL));
117c4762a1bSJed Brown
1189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ba));
1199566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b));
1209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&a));
121c4762a1bSJed Brown
1229566063dSJacob Faibussowitsch PetscCall(PetscSFView(sfA, NULL));
1239566063dSJacob Faibussowitsch PetscCall(PetscSFView(sfB, NULL));
1249566063dSJacob Faibussowitsch PetscCall(PetscSFView(sfBA, NULL));
1259566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfA));
1269566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfB));
1279566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfBA));
128c4762a1bSJed Brown
1299566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
130b122ec5aSJacob Faibussowitsch return 0;
131c4762a1bSJed Brown }
132c4762a1bSJed Brown
133c4762a1bSJed Brown /*TEST
134c4762a1bSJed Brown
135c4762a1bSJed Brown test:
136c4762a1bSJed Brown suffix: 1
137c4762a1bSJed Brown
138c4762a1bSJed Brown test:
139c4762a1bSJed Brown suffix: 2
140c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort"
141c4762a1bSJed Brown args: -sparse_sfB
142c4762a1bSJed Brown
143c4762a1bSJed Brown test:
144c4762a1bSJed Brown suffix: 2_window
145c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort"
146c4762a1bSJed Brown output_file: output/ex4_2.out
147c4762a1bSJed Brown args: -sparse_sfB -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
148dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
149c4762a1bSJed Brown
150c4762a1bSJed Brown # The nightly test suite with MPICH uses ch3:sock, which is broken when winsize == 0 in some of the processes
151c4762a1bSJed Brown test:
152c4762a1bSJed Brown suffix: 2_window_shared
153c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort"
154c4762a1bSJed Brown output_file: output/ex4_2.out
155c4762a1bSJed Brown args: -sparse_sfB -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor shared
156e3c15826SJunchao Zhang requires: defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) !defined(PETSC_HAVE_MPICH) defined(PETSC_HAVE_MPI_ONE_SIDED)
157c4762a1bSJed Brown
158c4762a1bSJed Brown TEST*/
159