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*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 20*9566063dSJacob Faibussowitsch PetscCallMPI(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 23*9566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(PETSC_COMM_WORLD, &sfA)); 24*9566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(PETSC_COMM_WORLD, &sfB)); 25*9566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(sfA)); 26*9566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(sfB)); 27c4762a1bSJed Brown 28*9566063dSJacob Faibussowitsch PetscCall(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 } 40*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleavesA, &ilocalA)); 41*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleavesA, &iremoteA)); 42*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleavesB, &ilocalB)); 43*9566063dSJacob Faibussowitsch PetscCall(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 66*9566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sfA, nrootsA, nleavesA, ilocalA, PETSC_OWN_POINTER, iremoteA, PETSC_OWN_POINTER)); 67*9566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sfB, nrootsB, nleavesB, ilocalB, PETSC_OWN_POINTER, iremoteB, PETSC_OWN_POINTER)); 68*9566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(sfA)); 69*9566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(sfB)); 70*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)sfA, "sfA")); 71*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)sfB, "sfB")); 72c4762a1bSJed Brown 73*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_WORLD, nrootsA, &a)); 74*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_WORLD, nleavesA, &b)); 75*9566063dSJacob Faibussowitsch PetscCall(PetscSFGetLeafRange(sfB, NULL, &maxleafB)); 76*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_WORLD, maxleafB+1, &ba)); 77*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(a, &arrayW)); 78c4762a1bSJed Brown for (i = 0; i < nrootsA; i++) { 79c4762a1bSJed Brown arrayW[i] = (PetscScalar)i; 80c4762a1bSJed Brown } 81*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(a, &arrayW)); 82c4762a1bSJed Brown 83*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial Vec A\n")); 84*9566063dSJacob Faibussowitsch PetscCall(VecView(a, NULL)); 85*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(a, &arrayR)); 86*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(b, &arrayW)); 87c4762a1bSJed Brown 88*9566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sfA, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE)); 89*9566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sfA, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE)); 90*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(b, &arrayW)); 91*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(a, &arrayR)); 92*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nBroadcast A->B over sfA\n")); 93*9566063dSJacob Faibussowitsch PetscCall(VecView(b, NULL)); 94c4762a1bSJed Brown 95*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(b, &arrayR)); 96*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(ba, &arrayW)); 97c4762a1bSJed Brown arrayW[0] = 10.0; /* Not touched by bcast */ 98*9566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sfB, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE)); 99*9566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sfB, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE)); 100*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(b, &arrayR)); 101*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ba, &arrayW)); 102c4762a1bSJed Brown 103*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nBroadcast B->BA over sfB\n")); 104*9566063dSJacob Faibussowitsch PetscCall(VecView(ba, NULL)); 105c4762a1bSJed Brown 106*9566063dSJacob Faibussowitsch PetscCall(PetscSFCompose(sfA, sfB, &sfBA)); 107*9566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(sfBA)); 108*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)sfBA, "(sfB o sfA)")); 109*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(a, &arrayR)); 110*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(ba, &arrayW)); 111c4762a1bSJed Brown arrayW[0] = 11.0; /* Not touched by bcast */ 112*9566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sfBA, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE)); 113*9566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sfBA, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE)); 114*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ba, &arrayW)); 115*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(a, &arrayR)); 116*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nBroadcast A->BA over sfBA (sfB o sfA)\n")); 117*9566063dSJacob Faibussowitsch PetscCall(VecView(ba, NULL)); 118c4762a1bSJed Brown 119*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ba)); 120*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 121*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&a)); 122c4762a1bSJed Brown 123*9566063dSJacob Faibussowitsch PetscCall(PetscSFView(sfA, NULL)); 124*9566063dSJacob Faibussowitsch PetscCall(PetscSFView(sfB, NULL)); 125*9566063dSJacob Faibussowitsch PetscCall(PetscSFView(sfBA, NULL)); 126*9566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfA)); 127*9566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfB)); 128*9566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfBA)); 129c4762a1bSJed Brown 130*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 131b122ec5aSJacob 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