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