1*3c429eb5SStefano Zampini static char help[] = "Tests ISRenumber.\n\n"; 2*3c429eb5SStefano Zampini 3*3c429eb5SStefano Zampini #include <petscis.h> 4*3c429eb5SStefano Zampini 5*3c429eb5SStefano Zampini PetscErrorCode TestRenumber(IS is, IS mult) 6*3c429eb5SStefano Zampini { 7*3c429eb5SStefano Zampini IS nis; 8*3c429eb5SStefano Zampini PetscInt N; 9*3c429eb5SStefano Zampini PetscErrorCode ierr; 10*3c429eb5SStefano Zampini 11*3c429eb5SStefano Zampini PetscFunctionBegin; 12*3c429eb5SStefano Zampini ierr = PetscPrintf(PetscObjectComm((PetscObject)is),"\n-----------------\n");CHKERRQ(ierr); 13*3c429eb5SStefano Zampini ierr = PetscPrintf(PetscObjectComm((PetscObject)is),"\nInitial\n");CHKERRQ(ierr); 14*3c429eb5SStefano Zampini ierr = ISView(is,NULL);CHKERRQ(ierr); 15*3c429eb5SStefano Zampini if (mult) { 16*3c429eb5SStefano Zampini ierr = PetscPrintf(PetscObjectComm((PetscObject)is),"\nMult\n");CHKERRQ(ierr); 17*3c429eb5SStefano Zampini ierr = ISView(mult,NULL);CHKERRQ(ierr); 18*3c429eb5SStefano Zampini } 19*3c429eb5SStefano Zampini ierr = ISRenumber(is,mult,&N,NULL);CHKERRQ(ierr); 20*3c429eb5SStefano Zampini ierr = PetscPrintf(PetscObjectComm((PetscObject)is),"\nRenumbered, unique entries %" PetscInt_FMT "\n",N);CHKERRQ(ierr); 21*3c429eb5SStefano Zampini ierr = ISRenumber(is,mult,NULL,&nis);CHKERRQ(ierr); 22*3c429eb5SStefano Zampini ierr = ISView(nis,NULL);CHKERRQ(ierr); 23*3c429eb5SStefano Zampini ierr = ISDestroy(&nis);CHKERRQ(ierr); 24*3c429eb5SStefano Zampini PetscFunctionReturn(0); 25*3c429eb5SStefano Zampini } 26*3c429eb5SStefano Zampini 27*3c429eb5SStefano Zampini int main(int argc, char **argv) 28*3c429eb5SStefano Zampini { 29*3c429eb5SStefano Zampini IS is; 30*3c429eb5SStefano Zampini PetscErrorCode ierr; 31*3c429eb5SStefano Zampini PetscMPIInt size, rank; 32*3c429eb5SStefano Zampini 33*3c429eb5SStefano Zampini ierr = PetscInitialize(&argc, &argv, NULL, help); if (ierr) return ierr; 34*3c429eb5SStefano Zampini ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 35*3c429eb5SStefano Zampini ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 36*3c429eb5SStefano Zampini 37*3c429eb5SStefano Zampini for (PetscInt c = 0; c < 3; c++) { 38*3c429eb5SStefano Zampini IS mult = NULL; 39*3c429eb5SStefano Zampini 40*3c429eb5SStefano Zampini ierr = ISCreateStride(PETSC_COMM_WORLD,0,0,0,&is);CHKERRQ(ierr); 41*3c429eb5SStefano Zampini if (c) { 42*3c429eb5SStefano Zampini PetscInt n; 43*3c429eb5SStefano Zampini ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr); 44*3c429eb5SStefano Zampini ierr = ISCreateStride(PETSC_COMM_WORLD,n,c-2,0,&mult);CHKERRQ(ierr); 45*3c429eb5SStefano Zampini } 46*3c429eb5SStefano Zampini ierr = TestRenumber(is,mult);CHKERRQ(ierr); 47*3c429eb5SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 48*3c429eb5SStefano Zampini ierr = ISDestroy(&mult);CHKERRQ(ierr); 49*3c429eb5SStefano Zampini 50*3c429eb5SStefano Zampini ierr = ISCreateStride(PETSC_COMM_WORLD,2,-rank-1,-4,&is);CHKERRQ(ierr); 51*3c429eb5SStefano Zampini if (c) { 52*3c429eb5SStefano Zampini PetscInt n; 53*3c429eb5SStefano Zampini ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr); 54*3c429eb5SStefano Zampini ierr = ISCreateStride(PETSC_COMM_WORLD,n,c-2,0,&mult);CHKERRQ(ierr); 55*3c429eb5SStefano Zampini } 56*3c429eb5SStefano Zampini ierr = TestRenumber(is,mult);CHKERRQ(ierr); 57*3c429eb5SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 58*3c429eb5SStefano Zampini ierr = ISDestroy(&mult);CHKERRQ(ierr); 59*3c429eb5SStefano Zampini 60*3c429eb5SStefano Zampini ierr = ISCreateStride(PETSC_COMM_WORLD,10,4+rank,2,&is);CHKERRQ(ierr); 61*3c429eb5SStefano Zampini if (c) { 62*3c429eb5SStefano Zampini PetscInt n; 63*3c429eb5SStefano Zampini ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr); 64*3c429eb5SStefano Zampini ierr = ISCreateStride(PETSC_COMM_WORLD,n,c-2,1,&mult);CHKERRQ(ierr); 65*3c429eb5SStefano Zampini } 66*3c429eb5SStefano Zampini ierr = TestRenumber(is,mult);CHKERRQ(ierr); 67*3c429eb5SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 68*3c429eb5SStefano Zampini ierr = ISDestroy(&mult);CHKERRQ(ierr); 69*3c429eb5SStefano Zampini 70*3c429eb5SStefano Zampini ierr = ISCreateStride(PETSC_COMM_WORLD,10,-rank-1,2,&is);CHKERRQ(ierr); 71*3c429eb5SStefano Zampini if (c) { 72*3c429eb5SStefano Zampini PetscInt n; 73*3c429eb5SStefano Zampini ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr); 74*3c429eb5SStefano Zampini ierr = ISCreateStride(PETSC_COMM_WORLD,n,c-2,1,&mult);CHKERRQ(ierr); 75*3c429eb5SStefano Zampini } 76*3c429eb5SStefano Zampini ierr = TestRenumber(is,mult);CHKERRQ(ierr); 77*3c429eb5SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 78*3c429eb5SStefano Zampini ierr = ISDestroy(&mult);CHKERRQ(ierr); 79*3c429eb5SStefano Zampini } 80*3c429eb5SStefano Zampini /* Finalize */ 81*3c429eb5SStefano Zampini ierr = PetscFinalize(); 82*3c429eb5SStefano Zampini return ierr; 83*3c429eb5SStefano Zampini } 84*3c429eb5SStefano Zampini 85*3c429eb5SStefano Zampini /*TEST 86*3c429eb5SStefano Zampini 87*3c429eb5SStefano Zampini test: 88*3c429eb5SStefano Zampini suffix: 1 89*3c429eb5SStefano Zampini nsize: {{1 2}separate output} 90*3c429eb5SStefano Zampini 91*3c429eb5SStefano Zampini TEST*/ 92