1*c4762a1bSJed Brown static char help[]= "Tests ISView() and ISLoad() \n\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petscis.h> 4*c4762a1bSJed Brown #include <petscviewer.h> 5*c4762a1bSJed Brown 6*c4762a1bSJed Brown int main(int argc,char **argv) 7*c4762a1bSJed Brown { 8*c4762a1bSJed Brown PetscErrorCode ierr; 9*c4762a1bSJed Brown PetscInt n = 3, izero[3] = {0,0,0}, j, i; 10*c4762a1bSJed Brown PetscInt ix[3][3][3] = {{{3,5,4},{1,7,9},{0,2,8}}, 11*c4762a1bSJed Brown {{0,2,8},{3,5,4},{1,7,9}}, 12*c4762a1bSJed Brown {{1,7,9},{0,2,8},{3,5,4}}}; 13*c4762a1bSJed Brown IS isx[3],il; 14*c4762a1bSJed Brown PetscMPIInt size,rank; 15*c4762a1bSJed Brown PetscViewer vx,vl; 16*c4762a1bSJed Brown PetscBool equal; 17*c4762a1bSJed Brown 18*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 19*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 20*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 21*c4762a1bSJed Brown if (size > 3) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_SIZ,"Example only works with up to three processes"); 22*c4762a1bSJed Brown 23*c4762a1bSJed Brown for (i = 0; i < 3; i++) { 24*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_WORLD,n,ix[i][rank],PETSC_COPY_VALUES,&isx[i]);CHKERRQ(ierr); 25*c4762a1bSJed Brown } 26*c4762a1bSJed Brown 27*c4762a1bSJed Brown for (j = 0; j < 3; j++) { 28*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile",FILE_MODE_WRITE,&vx);CHKERRQ(ierr); 29*c4762a1bSJed Brown ierr = ISView(isx[0],vx);CHKERRQ(ierr); 30*c4762a1bSJed Brown ierr = PetscViewerDestroy(&vx);CHKERRQ(ierr); 31*c4762a1bSJed Brown 32*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile",FILE_MODE_READ,&vl);CHKERRQ(ierr); 33*c4762a1bSJed Brown ierr = ISCreate(PETSC_COMM_WORLD,&il);CHKERRQ(ierr); 34*c4762a1bSJed Brown ierr = ISLoad(il,vl);CHKERRQ(ierr); 35*c4762a1bSJed Brown ierr = ISEqual(il,isx[0],&equal);CHKERRQ(ierr); 36*c4762a1bSJed Brown if (!equal) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Iteration %D - Index set loaded from file does not match",j); 37*c4762a1bSJed Brown ierr = ISDestroy(&il);CHKERRQ(ierr); 38*c4762a1bSJed Brown ierr = PetscViewerDestroy(&vl);CHKERRQ(ierr); 39*c4762a1bSJed Brown 40*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile",FILE_MODE_APPEND,&vx);CHKERRQ(ierr); 41*c4762a1bSJed Brown ierr = ISView(isx[1],vx);CHKERRQ(ierr); 42*c4762a1bSJed Brown ierr = ISView(isx[2],vx);CHKERRQ(ierr); 43*c4762a1bSJed Brown ierr = PetscViewerDestroy(&vx);CHKERRQ(ierr); 44*c4762a1bSJed Brown 45*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile",FILE_MODE_READ,&vl);CHKERRQ(ierr); 46*c4762a1bSJed Brown for (i = 0; i < 3; i++) { 47*c4762a1bSJed Brown ierr = ISCreate(PETSC_COMM_WORLD,&il);CHKERRQ(ierr); 48*c4762a1bSJed Brown ierr = ISLoad(il,vl);CHKERRQ(ierr); 49*c4762a1bSJed Brown ierr = ISEqual(il,isx[i],&equal);CHKERRQ(ierr); 50*c4762a1bSJed Brown if (!equal) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Iteration %D - Index set %D loaded from file does not match",j,i); 51*c4762a1bSJed Brown ierr = ISDestroy(&il);CHKERRQ(ierr); 52*c4762a1bSJed Brown } 53*c4762a1bSJed Brown ierr = PetscViewerDestroy(&vl);CHKERRQ(ierr); 54*c4762a1bSJed Brown 55*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile",FILE_MODE_READ,&vl);CHKERRQ(ierr); 56*c4762a1bSJed Brown for (i = 0; i < 3; i++) { 57*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_WORLD,n,izero,PETSC_COPY_VALUES,&il);CHKERRQ(ierr); 58*c4762a1bSJed Brown ierr = ISLoad(il,vl);CHKERRQ(ierr); 59*c4762a1bSJed Brown ierr = ISEqual(il,isx[i],&equal);CHKERRQ(ierr); 60*c4762a1bSJed Brown if (!equal) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Iteration %D - Index set %D loaded from file does not match",j,i); 61*c4762a1bSJed Brown ierr = ISDestroy(&il);CHKERRQ(ierr); 62*c4762a1bSJed Brown } 63*c4762a1bSJed Brown ierr = PetscViewerDestroy(&vl);CHKERRQ(ierr); 64*c4762a1bSJed Brown } 65*c4762a1bSJed Brown 66*c4762a1bSJed Brown for (j = 0; j < 3; j++) { 67*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile_noheader",FILE_MODE_WRITE,&vx);CHKERRQ(ierr); 68*c4762a1bSJed Brown ierr = PetscViewerBinarySetSkipHeader(vx,PETSC_TRUE);CHKERRQ(ierr); 69*c4762a1bSJed Brown for (i = 0; i < 3; i++) { 70*c4762a1bSJed Brown ierr = ISView(isx[i],vx);CHKERRQ(ierr); 71*c4762a1bSJed Brown } 72*c4762a1bSJed Brown ierr = PetscViewerDestroy(&vx);CHKERRQ(ierr); 73*c4762a1bSJed Brown 74*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile_noheader",FILE_MODE_READ,&vl);CHKERRQ(ierr); 75*c4762a1bSJed Brown ierr = PetscViewerBinarySetSkipHeader(vl,PETSC_TRUE);CHKERRQ(ierr); 76*c4762a1bSJed Brown for (i = 0; i < 3; i++) { 77*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_WORLD,n,izero,PETSC_COPY_VALUES,&il);CHKERRQ(ierr); 78*c4762a1bSJed Brown ierr = ISLoad(il,vl);CHKERRQ(ierr); 79*c4762a1bSJed Brown ierr = ISEqual(il,isx[i],&equal);CHKERRQ(ierr); 80*c4762a1bSJed Brown if (!equal) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Iteration %D - Index set %D loaded from file does not match",j,i); 81*c4762a1bSJed Brown ierr = ISDestroy(&il);CHKERRQ(ierr); 82*c4762a1bSJed Brown } 83*c4762a1bSJed Brown ierr = PetscViewerDestroy(&vl);CHKERRQ(ierr); 84*c4762a1bSJed Brown } 85*c4762a1bSJed Brown 86*c4762a1bSJed Brown for (i = 0; i < 3; i++) { 87*c4762a1bSJed Brown ierr = ISDestroy(&isx[i]);CHKERRQ(ierr); 88*c4762a1bSJed Brown } 89*c4762a1bSJed Brown 90*c4762a1bSJed Brown 91*c4762a1bSJed Brown for (j = 0; j < 2; j++) { 92*c4762a1bSJed Brown const char *filename = (j == 0) ? "testfile_isstride" : "testfile_isblock"; 93*c4762a1bSJed Brown PetscInt blocksize = (j == 0) ? 1 : size; 94*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_WRITE,&vx);CHKERRQ(ierr); 95*c4762a1bSJed Brown for (i = 0; i < 3; i++) { 96*c4762a1bSJed Brown if (j == 0) { 97*c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,n,rank,rank+1,&isx[i]);CHKERRQ(ierr); 98*c4762a1bSJed Brown } else { 99*c4762a1bSJed Brown ierr = ISCreateBlock(PETSC_COMM_WORLD,blocksize,n,ix[i][rank],PETSC_COPY_VALUES,&isx[i]);CHKERRQ(ierr); 100*c4762a1bSJed Brown } 101*c4762a1bSJed Brown ierr = ISView(isx[i],vx);CHKERRQ(ierr); 102*c4762a1bSJed Brown ierr = ISToGeneral(isx[i]);CHKERRQ(ierr); 103*c4762a1bSJed Brown } 104*c4762a1bSJed Brown ierr = PetscViewerDestroy(&vx);CHKERRQ(ierr); 105*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&vl);CHKERRQ(ierr); 106*c4762a1bSJed Brown for (i = 0; i < 3; i++) { 107*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_WORLD,blocksize*n,izero,PETSC_COPY_VALUES,&il);CHKERRQ(ierr); 108*c4762a1bSJed Brown ierr = ISLoad(il,vl);CHKERRQ(ierr); 109*c4762a1bSJed Brown ierr = ISEqual(il,isx[i],&equal);CHKERRQ(ierr); 110*c4762a1bSJed Brown if (!equal) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Iteration %D - Index set %D loaded from file does not match",j,i); 111*c4762a1bSJed Brown ierr = ISDestroy(&il);CHKERRQ(ierr); 112*c4762a1bSJed Brown } 113*c4762a1bSJed Brown ierr = PetscViewerDestroy(&vl);CHKERRQ(ierr); 114*c4762a1bSJed Brown for (i = 0; i < 3; i++) { 115*c4762a1bSJed Brown ierr = ISDestroy(&isx[i]);CHKERRQ(ierr); 116*c4762a1bSJed Brown } 117*c4762a1bSJed Brown } 118*c4762a1bSJed Brown 119*c4762a1bSJed Brown ierr = PetscFinalize(); 120*c4762a1bSJed Brown return ierr; 121*c4762a1bSJed Brown } 122*c4762a1bSJed Brown 123*c4762a1bSJed Brown /*TEST 124*c4762a1bSJed Brown 125*c4762a1bSJed Brown testset: 126*c4762a1bSJed Brown args: -viewer_binary_mpiio 0 127*c4762a1bSJed Brown output_file: output/ex2_1.out 128*c4762a1bSJed Brown test: 129*c4762a1bSJed Brown suffix: stdio_1 130*c4762a1bSJed Brown nsize: 1 131*c4762a1bSJed Brown test: 132*c4762a1bSJed Brown suffix: stdio_2 133*c4762a1bSJed Brown nsize: 2 134*c4762a1bSJed Brown test: 135*c4762a1bSJed Brown suffix: stdio_3 136*c4762a1bSJed Brown nsize: 3 137*c4762a1bSJed Brown 138*c4762a1bSJed Brown testset: 139*c4762a1bSJed Brown requires: mpiio 140*c4762a1bSJed Brown args: -viewer_binary_mpiio 1 141*c4762a1bSJed Brown output_file: output/ex2_1.out 142*c4762a1bSJed Brown test: 143*c4762a1bSJed Brown suffix: mpiio_1 144*c4762a1bSJed Brown nsize: 1 145*c4762a1bSJed Brown test: 146*c4762a1bSJed Brown suffix: mpiio_2 147*c4762a1bSJed Brown nsize: 2 148*c4762a1bSJed Brown test: 149*c4762a1bSJed Brown suffix: mpiio_3 150*c4762a1bSJed Brown nsize: 3 151*c4762a1bSJed Brown 152*c4762a1bSJed Brown TEST*/ 153