1c4762a1bSJed Brown static char help[]= "Test PetscSFFCompose when the ilocal array is not the identity\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; 10c4762a1bSJed Brown PetscInt nrootsA, nleavesA, nrootsB, nleavesB; 11c4762a1bSJed Brown PetscInt *ilocalA, *ilocalB; 12c4762a1bSJed Brown PetscSFNode *iremoteA, *iremoteB; 13c4762a1bSJed Brown Vec a, b, ba; 14c4762a1bSJed Brown const PetscScalar *arrayR; 15c4762a1bSJed Brown PetscScalar *arrayW; 16c4762a1bSJed Brown PetscMPIInt size; 17c4762a1bSJed Brown PetscInt i; 18c4762a1bSJed Brown PetscInt maxleafB; 19c4762a1bSJed Brown PetscBool flag = PETSC_FALSE; 20c4762a1bSJed Brown 21c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 22ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 23c4762a1bSJed Brown 24c4762a1bSJed Brown if (size != 1) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_USER, "Only coded for one MPI process"); 25c4762a1bSJed Brown 26c4762a1bSJed Brown ierr = PetscSFCreate(PETSC_COMM_WORLD, &sfA);CHKERRQ(ierr); 27c4762a1bSJed Brown ierr = PetscSFCreate(PETSC_COMM_WORLD, &sfB);CHKERRQ(ierr); 28c4762a1bSJed Brown ierr = PetscSFSetFromOptions(sfA);CHKERRQ(ierr); 29c4762a1bSJed Brown ierr = PetscSFSetFromOptions(sfB);CHKERRQ(ierr); 30c4762a1bSJed Brown 31c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-sparse_sfB",&flag,NULL);CHKERRQ(ierr); 32c4762a1bSJed Brown 33c4762a1bSJed Brown if (flag) { 34c4762a1bSJed Brown /* sfA permutes indices, sfB has sparse leaf space. */ 35c4762a1bSJed Brown nrootsA = 3; 36c4762a1bSJed Brown nleavesA = 3; 37c4762a1bSJed Brown nrootsB = 3; 38c4762a1bSJed Brown nleavesB = 2; 39c4762a1bSJed Brown } else { 40c4762a1bSJed Brown /* sfA reverses indices, sfB is identity */ 41c4762a1bSJed Brown nrootsA = nrootsB = nleavesA = nleavesB = 4; 42c4762a1bSJed Brown } 43c4762a1bSJed Brown ierr = PetscMalloc1(nleavesA, &ilocalA);CHKERRQ(ierr); 44c4762a1bSJed Brown ierr = PetscMalloc1(nleavesA, &iremoteA);CHKERRQ(ierr); 45c4762a1bSJed Brown ierr = PetscMalloc1(nleavesB, &ilocalB);CHKERRQ(ierr); 46c4762a1bSJed Brown ierr = PetscMalloc1(nleavesB, &iremoteB);CHKERRQ(ierr); 47c4762a1bSJed Brown 48c4762a1bSJed Brown for (i = 0; i < nleavesA; i++) { 49c4762a1bSJed Brown iremoteA[i].rank = 0; 50c4762a1bSJed Brown iremoteA[i].index = i; 51c4762a1bSJed Brown if (flag) { 52c4762a1bSJed Brown ilocalA[i] = (i + 1) % nleavesA; 53c4762a1bSJed Brown } else { 54c4762a1bSJed Brown ilocalA[i] = nleavesA - i - 1; 55c4762a1bSJed Brown } 56c4762a1bSJed Brown } 57c4762a1bSJed Brown 58c4762a1bSJed Brown for (i = 0; i < nleavesB; i++) { 59c4762a1bSJed Brown iremoteB[i].rank = 0; 60c4762a1bSJed Brown if (flag) { 61c4762a1bSJed Brown ilocalB[i] = nleavesB - i; 62c4762a1bSJed Brown iremoteB[i].index = nleavesB - i - 1; 63c4762a1bSJed Brown } else { 64c4762a1bSJed Brown ilocalB[i] = i; 65c4762a1bSJed Brown iremoteB[i].index = i; 66c4762a1bSJed Brown } 67c4762a1bSJed Brown } 68c4762a1bSJed Brown 69c4762a1bSJed Brown ierr = PetscSFSetGraph(sfA, nrootsA, nleavesA, ilocalA, PETSC_OWN_POINTER, iremoteA, PETSC_OWN_POINTER);CHKERRQ(ierr); 70c4762a1bSJed Brown ierr = PetscSFSetGraph(sfB, nrootsB, nleavesB, ilocalB, PETSC_OWN_POINTER, iremoteB, PETSC_OWN_POINTER);CHKERRQ(ierr); 71c4762a1bSJed Brown ierr = PetscSFSetUp(sfA);CHKERRQ(ierr); 72c4762a1bSJed Brown ierr = PetscSFSetUp(sfB);CHKERRQ(ierr); 73c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)sfA, "sfA");CHKERRQ(ierr); 74c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)sfB, "sfB");CHKERRQ(ierr); 75c4762a1bSJed Brown 76c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_WORLD, nrootsA, &a);CHKERRQ(ierr); 77c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_WORLD, nleavesA, &b);CHKERRQ(ierr); 78c4762a1bSJed Brown ierr = PetscSFGetLeafRange(sfB, NULL, &maxleafB);CHKERRQ(ierr); 79c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_WORLD, maxleafB+1, &ba);CHKERRQ(ierr); 80c4762a1bSJed Brown ierr = VecGetArray(a, &arrayW);CHKERRQ(ierr); 81c4762a1bSJed Brown for (i = 0; i < nrootsA; i++) { 82c4762a1bSJed Brown arrayW[i] = (PetscScalar)i; 83c4762a1bSJed Brown } 84c4762a1bSJed Brown ierr = VecRestoreArray(a, &arrayW);CHKERRQ(ierr); 85c4762a1bSJed Brown 86c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial Vec A\n");CHKERRQ(ierr); 87c4762a1bSJed Brown ierr = VecView(a, NULL);CHKERRQ(ierr); 88c4762a1bSJed Brown ierr = VecGetArrayRead(a, &arrayR);CHKERRQ(ierr); 89c4762a1bSJed Brown ierr = VecGetArray(b, &arrayW);CHKERRQ(ierr); 90c4762a1bSJed Brown 91ad227feaSJunchao Zhang ierr = PetscSFBcastBegin(sfA, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE);CHKERRQ(ierr); 92ad227feaSJunchao Zhang ierr = PetscSFBcastEnd(sfA, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE);CHKERRQ(ierr); 93c4762a1bSJed Brown ierr = VecRestoreArray(b, &arrayW);CHKERRQ(ierr); 94c4762a1bSJed Brown ierr = VecRestoreArrayRead(a, &arrayR);CHKERRQ(ierr); 95c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "\nBroadcast A->B over sfA\n");CHKERRQ(ierr); 96c4762a1bSJed Brown ierr = VecView(b, NULL);CHKERRQ(ierr); 97c4762a1bSJed Brown 98c4762a1bSJed Brown ierr = VecGetArrayRead(b, &arrayR);CHKERRQ(ierr); 99c4762a1bSJed Brown ierr = VecGetArray(ba, &arrayW);CHKERRQ(ierr); 100c4762a1bSJed Brown arrayW[0] = 10.0; /* Not touched by bcast */ 101ad227feaSJunchao Zhang ierr = PetscSFBcastBegin(sfB, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE);CHKERRQ(ierr); 102ad227feaSJunchao Zhang ierr = PetscSFBcastEnd(sfB, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE);CHKERRQ(ierr); 103c4762a1bSJed Brown ierr = VecRestoreArrayRead(b, &arrayR);CHKERRQ(ierr); 104c4762a1bSJed Brown ierr = VecRestoreArray(ba, &arrayW);CHKERRQ(ierr); 105c4762a1bSJed Brown 106c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "\nBroadcast B->BA over sfB\n");CHKERRQ(ierr); 107c4762a1bSJed Brown ierr = VecView(ba, NULL);CHKERRQ(ierr); 108c4762a1bSJed Brown 109c4762a1bSJed Brown ierr = PetscSFCompose(sfA, sfB, &sfBA);CHKERRQ(ierr); 110c4762a1bSJed Brown ierr = PetscSFSetFromOptions(sfBA);CHKERRQ(ierr); 111c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)sfBA, "(sfB o sfA)");CHKERRQ(ierr); 112c4762a1bSJed Brown ierr = VecGetArrayRead(a, &arrayR);CHKERRQ(ierr); 113c4762a1bSJed Brown ierr = VecGetArray(ba, &arrayW);CHKERRQ(ierr); 114c4762a1bSJed Brown arrayW[0] = 11.0; /* Not touched by bcast */ 115ad227feaSJunchao Zhang ierr = PetscSFBcastBegin(sfBA, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE);CHKERRQ(ierr); 116ad227feaSJunchao Zhang ierr = PetscSFBcastEnd(sfBA, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE);CHKERRQ(ierr); 117c4762a1bSJed Brown ierr = VecRestoreArray(ba, &arrayW);CHKERRQ(ierr); 118c4762a1bSJed Brown ierr = VecRestoreArrayRead(a, &arrayR);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "\nBroadcast A->BA over sfBA (sfB o sfA)\n");CHKERRQ(ierr); 120c4762a1bSJed Brown ierr = VecView(ba, NULL);CHKERRQ(ierr); 121c4762a1bSJed Brown 122c4762a1bSJed Brown ierr = VecDestroy(&ba);CHKERRQ(ierr); 123c4762a1bSJed Brown ierr = VecDestroy(&b);CHKERRQ(ierr); 124c4762a1bSJed Brown ierr = VecDestroy(&a);CHKERRQ(ierr); 125c4762a1bSJed Brown 126c4762a1bSJed Brown ierr = PetscSFView(sfA, NULL);CHKERRQ(ierr); 127c4762a1bSJed Brown ierr = PetscSFView(sfB, NULL);CHKERRQ(ierr); 128c4762a1bSJed Brown ierr = PetscSFView(sfBA, NULL);CHKERRQ(ierr); 129c4762a1bSJed Brown ierr = PetscSFDestroy(&sfA);CHKERRQ(ierr); 130c4762a1bSJed Brown ierr = PetscSFDestroy(&sfB);CHKERRQ(ierr); 131c4762a1bSJed Brown ierr = PetscSFDestroy(&sfBA);CHKERRQ(ierr); 132c4762a1bSJed Brown 133c4762a1bSJed Brown ierr = PetscFinalize(); 134c4762a1bSJed Brown 135c4762a1bSJed Brown return ierr; 136c4762a1bSJed Brown } 137c4762a1bSJed Brown 138c4762a1bSJed Brown /*TEST 139c4762a1bSJed Brown 140c4762a1bSJed Brown test: 141c4762a1bSJed Brown suffix: 1 142c4762a1bSJed Brown 143c4762a1bSJed Brown test: 144c4762a1bSJed Brown suffix: 2 145c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 146c4762a1bSJed Brown args: -sparse_sfB 147c4762a1bSJed Brown 148c4762a1bSJed Brown test: 149c4762a1bSJed Brown suffix: 2_window 150c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 151c4762a1bSJed Brown output_file: output/ex4_2.out 152c4762a1bSJed Brown args: -sparse_sfB -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 153*dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 154c4762a1bSJed Brown 155c4762a1bSJed Brown # The nightly test suite with MPICH uses ch3:sock, which is broken when winsize == 0 in some of the processes 156c4762a1bSJed Brown test: 157c4762a1bSJed Brown suffix: 2_window_shared 158c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 159c4762a1bSJed Brown output_file: output/ex4_2.out 160c4762a1bSJed Brown args: -sparse_sfB -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor shared 161*dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) !defined(PETSC_HAVE_MPICH_NUMVERSION) defined(PETSC_HAVE_MPI_ONE_SIDED) 162c4762a1bSJed Brown 163c4762a1bSJed Brown TEST*/ 164