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