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 PetscErrorCode ierr; 8c4762a1bSJed Brown PetscSF sfA, sfB, sfBA; 9c4762a1bSJed Brown PetscInt nrootsA, nleavesA, nrootsB, nleavesB; 10c4762a1bSJed Brown PetscInt *ilocalA, *ilocalB; 11c4762a1bSJed Brown PetscSFNode *iremoteA, *iremoteB; 12c4762a1bSJed Brown Vec a, b, ba; 13c4762a1bSJed Brown const PetscScalar *arrayR; 14c4762a1bSJed Brown PetscScalar *arrayW; 15c4762a1bSJed Brown PetscMPIInt size; 16c4762a1bSJed Brown PetscInt i; 17c4762a1bSJed Brown PetscInt maxleafB; 18c4762a1bSJed Brown PetscBool flag = PETSC_FALSE; 19c4762a1bSJed Brown 20c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 21*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 224643bc0bSVaclav Hapla PetscCheck(size == 1,PETSC_COMM_WORLD, PETSC_ERR_USER, "Only coded for one MPI process"); 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 29*5f80ce2aSJacob Faibussowitsch CHKERRQ(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 } 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 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 67*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFSetGraph(sfA, nrootsA, nleavesA, ilocalA, PETSC_OWN_POINTER, iremoteA, PETSC_OWN_POINTER)); 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFSetGraph(sfB, nrootsB, nleavesB, ilocalB, PETSC_OWN_POINTER, iremoteB, PETSC_OWN_POINTER)); 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFSetUp(sfA)); 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFSetUp(sfB)); 71*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)sfA, "sfA")); 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)sfB, "sfB")); 73c4762a1bSJed Brown 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_WORLD, nrootsA, &a)); 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_WORLD, nleavesA, &b)); 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFGetLeafRange(sfB, NULL, &maxleafB)); 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_WORLD, maxleafB+1, &ba)); 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(a, &arrayW)); 79c4762a1bSJed Brown for (i = 0; i < nrootsA; i++) { 80c4762a1bSJed Brown arrayW[i] = (PetscScalar)i; 81c4762a1bSJed Brown } 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(a, &arrayW)); 83c4762a1bSJed Brown 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Initial Vec A\n")); 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(a, NULL)); 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(a, &arrayR)); 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(b, &arrayW)); 88c4762a1bSJed Brown 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFBcastBegin(sfA, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE)); 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFBcastEnd(sfA, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE)); 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(b, &arrayW)); 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(a, &arrayR)); 93*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "\nBroadcast A->B over sfA\n")); 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(b, NULL)); 95c4762a1bSJed Brown 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(b, &arrayR)); 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(ba, &arrayW)); 98c4762a1bSJed Brown arrayW[0] = 10.0; /* Not touched by bcast */ 99*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFBcastBegin(sfB, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE)); 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFBcastEnd(sfB, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE)); 101*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(b, &arrayR)); 102*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(ba, &arrayW)); 103c4762a1bSJed Brown 104*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "\nBroadcast B->BA over sfB\n")); 105*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(ba, NULL)); 106c4762a1bSJed Brown 107*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFCompose(sfA, sfB, &sfBA)); 108*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFSetFromOptions(sfBA)); 109*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)sfBA, "(sfB o sfA)")); 110*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(a, &arrayR)); 111*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(ba, &arrayW)); 112c4762a1bSJed Brown arrayW[0] = 11.0; /* Not touched by bcast */ 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFBcastBegin(sfBA, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE)); 114*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFBcastEnd(sfBA, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE)); 115*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(ba, &arrayW)); 116*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(a, &arrayR)); 117*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "\nBroadcast A->BA over sfBA (sfB o sfA)\n")); 118*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(ba, NULL)); 119c4762a1bSJed Brown 120*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&ba)); 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b)); 122*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&a)); 123c4762a1bSJed Brown 124*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFView(sfA, NULL)); 125*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFView(sfB, NULL)); 126*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFView(sfBA, NULL)); 127*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFDestroy(&sfA)); 128*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFDestroy(&sfB)); 129*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFDestroy(&sfBA)); 130c4762a1bSJed Brown 131c4762a1bSJed Brown ierr = PetscFinalize(); 132c4762a1bSJed Brown return ierr; 133c4762a1bSJed Brown } 134c4762a1bSJed Brown 135c4762a1bSJed Brown /*TEST 136c4762a1bSJed Brown 137c4762a1bSJed Brown test: 138c4762a1bSJed Brown suffix: 1 139c4762a1bSJed Brown 140c4762a1bSJed Brown test: 141c4762a1bSJed Brown suffix: 2 142c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 143c4762a1bSJed Brown args: -sparse_sfB 144c4762a1bSJed Brown 145c4762a1bSJed Brown test: 146c4762a1bSJed Brown suffix: 2_window 147c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 148c4762a1bSJed Brown output_file: output/ex4_2.out 149c4762a1bSJed Brown args: -sparse_sfB -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 150dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 151c4762a1bSJed Brown 152c4762a1bSJed Brown # The nightly test suite with MPICH uses ch3:sock, which is broken when winsize == 0 in some of the processes 153c4762a1bSJed Brown test: 154c4762a1bSJed Brown suffix: 2_window_shared 155c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 156c4762a1bSJed Brown output_file: output/ex4_2.out 157c4762a1bSJed Brown args: -sparse_sfB -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor shared 158dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) !defined(PETSC_HAVE_MPICH_NUMVERSION) defined(PETSC_HAVE_MPI_ONE_SIDED) 159c4762a1bSJed Brown 160c4762a1bSJed Brown TEST*/ 161