1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[]= "Tests ISSetBlockSize() on ISBlock().\n\n"; 3*c4762a1bSJed Brown 4*c4762a1bSJed Brown #include <petscis.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 PetscErrorCode ierr; 10*c4762a1bSJed Brown PetscInt bs = 2,n = 3,ix[3] = {1,7,9}; 11*c4762a1bSJed Brown const PetscInt *indices; 12*c4762a1bSJed Brown IS is; 13*c4762a1bSJed Brown PetscBool broken = PETSC_FALSE; 14*c4762a1bSJed Brown 15*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 16*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-broken",&broken,NULL);CHKERRQ(ierr); 17*c4762a1bSJed Brown ierr = ISCreateBlock(PETSC_COMM_SELF,bs,n,ix,PETSC_COPY_VALUES,&is);CHKERRQ(ierr); 18*c4762a1bSJed Brown ierr = ISGetIndices(is,&indices);CHKERRQ(ierr); 19*c4762a1bSJed Brown ierr = PetscIntView(bs*3,indices,NULL);CHKERRQ(ierr); 20*c4762a1bSJed Brown ierr = ISRestoreIndices(is,&indices);CHKERRQ(ierr); 21*c4762a1bSJed Brown if (broken) { 22*c4762a1bSJed Brown ierr = ISSetBlockSize(is,3);CHKERRQ(ierr); 23*c4762a1bSJed Brown ierr = ISGetIndices(is,&indices);CHKERRQ(ierr); 24*c4762a1bSJed Brown ierr = PetscIntView(bs*3,indices,NULL);CHKERRQ(ierr); 25*c4762a1bSJed Brown ierr = ISRestoreIndices(is,&indices);CHKERRQ(ierr); 26*c4762a1bSJed Brown } 27*c4762a1bSJed Brown ierr = ISDestroy(&is);CHKERRQ(ierr); 28*c4762a1bSJed Brown 29*c4762a1bSJed Brown ierr = PetscFinalize(); 30*c4762a1bSJed Brown return ierr; 31*c4762a1bSJed Brown } 32*c4762a1bSJed Brown 33*c4762a1bSJed Brown 34*c4762a1bSJed Brown 35*c4762a1bSJed Brown /*TEST 36*c4762a1bSJed Brown 37*c4762a1bSJed Brown test: 38*c4762a1bSJed Brown 39*c4762a1bSJed Brown test: 40*c4762a1bSJed Brown suffix: 2 41*c4762a1bSJed Brown args: -broken 42*c4762a1bSJed Brown filter: Error: grep -o "[0]PETSC ERROR: Object is in wrong state" 43*c4762a1bSJed Brown 44*c4762a1bSJed Brown TEST*/ 45