1*c4762a1bSJed Brown static const char help[] = "Test freeing of MPI types in PetscSF\n\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petscvec.h> 4*c4762a1bSJed Brown #include <petscsf.h> 5*c4762a1bSJed Brown #include <petscviewer.h> 6*c4762a1bSJed Brown 7*c4762a1bSJed Brown int main(int argc, char **argv) 8*c4762a1bSJed Brown { 9*c4762a1bSJed Brown PetscInt ierr; 10*c4762a1bSJed Brown PetscSF sf; 11*c4762a1bSJed Brown Vec A,Aout; 12*c4762a1bSJed Brown PetscScalar *bufA; 13*c4762a1bSJed Brown PetscScalar *bufAout; 14*c4762a1bSJed Brown PetscMPIInt rank, size; 15*c4762a1bSJed Brown PetscInt nroots, nleaves; 16*c4762a1bSJed Brown PetscInt i; 17*c4762a1bSJed Brown PetscInt *ilocal; 18*c4762a1bSJed Brown PetscSFNode *iremote; 19*c4762a1bSJed Brown PetscBool test_dupped_type; 20*c4762a1bSJed Brown MPI_Datatype contig; 21*c4762a1bSJed Brown 22*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 23*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 24*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 25*c4762a1bSJed Brown 26*c4762a1bSJed Brown if (size != 1) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_USER, "Only coded for one MPI process"); 27*c4762a1bSJed Brown 28*c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,"","PetscSF type freeing options","none");CHKERRQ(ierr); 29*c4762a1bSJed Brown test_dupped_type = PETSC_FALSE; 30*c4762a1bSJed Brown ierr = PetscOptionsBool("-test_dupped_type", "Test dupped input type","",test_dupped_type,&test_dupped_type,NULL);CHKERRQ(ierr); 31*c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 32*c4762a1bSJed Brown 33*c4762a1bSJed Brown ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); 34*c4762a1bSJed Brown ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); 35*c4762a1bSJed Brown 36*c4762a1bSJed Brown nleaves = 1; 37*c4762a1bSJed Brown nroots = 1; 38*c4762a1bSJed Brown ierr = PetscMalloc1(nleaves,&ilocal);CHKERRQ(ierr); 39*c4762a1bSJed Brown 40*c4762a1bSJed Brown for (i = 0; i<nleaves; i++) { 41*c4762a1bSJed Brown ilocal[i] = i; 42*c4762a1bSJed Brown } 43*c4762a1bSJed Brown 44*c4762a1bSJed Brown ierr = PetscMalloc1(nleaves,&iremote);CHKERRQ(ierr); 45*c4762a1bSJed Brown iremote[0].rank = 0; 46*c4762a1bSJed Brown iremote[0].index = 0; 47*c4762a1bSJed Brown ierr = PetscSFSetGraph(sf,nroots,nleaves,ilocal,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 48*c4762a1bSJed Brown ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 49*c4762a1bSJed Brown ierr = PetscSFView(sf,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 50*c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 51*c4762a1bSJed Brown ierr = VecSetSizes(A,4,PETSC_DETERMINE);CHKERRQ(ierr); 52*c4762a1bSJed Brown ierr = VecSetFromOptions(A);CHKERRQ(ierr); 53*c4762a1bSJed Brown ierr = VecSetUp(A);CHKERRQ(ierr); 54*c4762a1bSJed Brown 55*c4762a1bSJed Brown ierr = VecDuplicate(A,&Aout);CHKERRQ(ierr); 56*c4762a1bSJed Brown ierr = VecGetArray(A,&bufA);CHKERRQ(ierr); 57*c4762a1bSJed Brown for (i=0; i<4; i++) { 58*c4762a1bSJed Brown bufA[i] = (PetscScalar)i; 59*c4762a1bSJed Brown } 60*c4762a1bSJed Brown ierr = VecRestoreArray(A,&bufA);CHKERRQ(ierr); 61*c4762a1bSJed Brown 62*c4762a1bSJed Brown ierr = VecGetArrayRead(A,(const PetscScalar**)&bufA);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = VecGetArray(Aout,&bufAout);CHKERRQ(ierr); 64*c4762a1bSJed Brown 65*c4762a1bSJed Brown ierr = MPI_Type_contiguous(4, MPIU_SCALAR, &contig);CHKERRQ(ierr); 66*c4762a1bSJed Brown ierr = MPI_Type_commit(&contig);CHKERRQ(ierr); 67*c4762a1bSJed Brown 68*c4762a1bSJed Brown if (test_dupped_type) { 69*c4762a1bSJed Brown MPI_Datatype tmp; 70*c4762a1bSJed Brown ierr = MPI_Type_dup(contig, &tmp);CHKERRQ(ierr); 71*c4762a1bSJed Brown ierr = MPI_Type_free(&contig);CHKERRQ(ierr); 72*c4762a1bSJed Brown contig = tmp; 73*c4762a1bSJed Brown } 74*c4762a1bSJed Brown for (i=0;i<10000;i++) { 75*c4762a1bSJed Brown ierr = PetscSFBcastBegin(sf,contig,bufA,bufAout);CHKERRQ(ierr); 76*c4762a1bSJed Brown ierr = PetscSFBcastEnd(sf,contig,bufA,bufAout);CHKERRQ(ierr); 77*c4762a1bSJed Brown } 78*c4762a1bSJed Brown ierr = VecRestoreArrayRead(A,(const PetscScalar**)&bufA);CHKERRQ(ierr); 79*c4762a1bSJed Brown ierr = VecRestoreArray(Aout,&bufAout);CHKERRQ(ierr); 80*c4762a1bSJed Brown 81*c4762a1bSJed Brown ierr = VecView(Aout,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 82*c4762a1bSJed Brown ierr = VecDestroy(&A);CHKERRQ(ierr); 83*c4762a1bSJed Brown ierr = VecDestroy(&Aout);CHKERRQ(ierr); 84*c4762a1bSJed Brown ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 85*c4762a1bSJed Brown ierr = MPI_Type_free(&contig);CHKERRQ(ierr); 86*c4762a1bSJed Brown ierr = PetscFinalize(); 87*c4762a1bSJed Brown return ierr; 88*c4762a1bSJed Brown } 89*c4762a1bSJed Brown 90*c4762a1bSJed Brown /*TEST 91*c4762a1bSJed Brown 92*c4762a1bSJed Brown test: 93*c4762a1bSJed Brown suffix: basic 94*c4762a1bSJed Brown args: -sf_type basic 95*c4762a1bSJed Brown 96*c4762a1bSJed Brown test: 97*c4762a1bSJed Brown suffix: basic_dupped 98*c4762a1bSJed Brown args: -test_dupped_type -sf_type basic 99*c4762a1bSJed Brown 100*c4762a1bSJed Brown test: 101*c4762a1bSJed Brown suffix: window 102*c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 103*c4762a1bSJed Brown args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create allocate dynamic}} 104*c4762a1bSJed Brown requires: define(PETSC_HAVE_MPI_ONE_SIDED) 105*c4762a1bSJed Brown 106*c4762a1bSJed Brown test: 107*c4762a1bSJed Brown suffix: window_dupped 108*c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 109*c4762a1bSJed Brown args: -test_dupped_type -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create allocate dynamic}} 110*c4762a1bSJed Brown requires: define(PETSC_HAVE_MPI_ONE_SIDED) 111*c4762a1bSJed Brown 112*c4762a1bSJed Brown test: 113*c4762a1bSJed Brown suffix: window_shared 114*c4762a1bSJed Brown output_file: output/ex3_window.out 115*c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 116*c4762a1bSJed Brown args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor shared 117*c4762a1bSJed Brown requires: define(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) define(PETSC_HAVE_MPI_ONE_SIDED) 118*c4762a1bSJed Brown 119*c4762a1bSJed Brown test: 120*c4762a1bSJed Brown suffix: window_dupped_shared 121*c4762a1bSJed Brown output_file: output/ex3_window_dupped.out 122*c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 123*c4762a1bSJed Brown args: -test_dupped_type -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor shared 124*c4762a1bSJed Brown requires: define(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) define(PETSC_HAVE_MPI_ONE_SIDED) 125*c4762a1bSJed Brown 126*c4762a1bSJed Brown TEST*/ 127