1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[]= "Tests ISSetBlockSize() on ISBlock().\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscis.h> 5c4762a1bSJed Brown #include <petscviewer.h> 6c4762a1bSJed Brown 7c4762a1bSJed Brown int main(int argc,char **argv) 8c4762a1bSJed Brown { 9c4762a1bSJed Brown PetscInt bs = 2,n = 3,ix[3] = {1,7,9}; 10c4762a1bSJed Brown const PetscInt *indices; 11c4762a1bSJed Brown IS is; 12c4762a1bSJed Brown PetscBool broken = PETSC_FALSE; 13c4762a1bSJed Brown 14*327415f7SBarry Smith PetscFunctionBeginUser; 159566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 169566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-broken",&broken,NULL)); 179566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,n,ix,PETSC_COPY_VALUES,&is)); 189566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is,&indices)); 199566063dSJacob Faibussowitsch PetscCall(PetscIntView(bs*3,indices,NULL)); 209566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is,&indices)); 21c4762a1bSJed Brown if (broken) { 229566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize(is,3)); 239566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is,&indices)); 249566063dSJacob Faibussowitsch PetscCall(PetscIntView(bs*3,indices,NULL)); 259566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is,&indices)); 26c4762a1bSJed Brown } 279566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 28c4762a1bSJed Brown 299566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 30b122ec5aSJacob Faibussowitsch return 0; 31c4762a1bSJed Brown } 32c4762a1bSJed Brown 33c4762a1bSJed Brown /*TEST 34c4762a1bSJed Brown 35c4762a1bSJed Brown test: 36c4762a1bSJed Brown 37c4762a1bSJed Brown test: 38c4762a1bSJed Brown suffix: 2 39c4762a1bSJed Brown args: -broken 40c4762a1bSJed Brown filter: Error: grep -o "[0]PETSC ERROR: Object is in wrong state" 41c4762a1bSJed Brown 42c4762a1bSJed Brown TEST*/ 43