1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Demonstrates constructing an application ordering.\n\n"; 3*c4762a1bSJed Brown 4*c4762a1bSJed Brown #include <petscsys.h> 5*c4762a1bSJed Brown #include <petscao.h> 6*c4762a1bSJed Brown #include <petscviewer.h> 7*c4762a1bSJed Brown 8*c4762a1bSJed Brown int main(int argc,char **argv) 9*c4762a1bSJed Brown { 10*c4762a1bSJed Brown PetscErrorCode ierr; 11*c4762a1bSJed Brown PetscInt i,n = 5; 12*c4762a1bSJed Brown PetscInt getpetsc[] = {0,3,4},getapp[] = {2,1,9,7}; 13*c4762a1bSJed Brown PetscInt getpetsc1[] = {0,3,4},getapp1[] = {2,1,9,7}; 14*c4762a1bSJed Brown PetscInt getpetsc2[] = {0,3,4},getapp2[] = {2,1,9,7}; 15*c4762a1bSJed Brown PetscInt getpetsc3[] = {0,3,4},getapp3[] = {2,1,9,7}; 16*c4762a1bSJed Brown PetscInt getpetsc4[] = {0,3,4},getapp4[] = {2,1,9,7}; 17*c4762a1bSJed Brown PetscMPIInt rank,size; 18*c4762a1bSJed Brown IS ispetsc,isapp; 19*c4762a1bSJed Brown AO ao; 20*c4762a1bSJed Brown const PetscInt *app; 21*c4762a1bSJed Brown 22*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 23*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr); 24*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 25*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 26*c4762a1bSJed Brown 27*c4762a1bSJed Brown /* create the index sets */ 28*c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,n,rank,size,&isapp);CHKERRQ(ierr); 29*c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,n,n*rank,1,&ispetsc);CHKERRQ(ierr); /* natural numbering */ 30*c4762a1bSJed Brown 31*c4762a1bSJed Brown /* create the application ordering */ 32*c4762a1bSJed Brown ierr = AOCreateBasicIS(isapp,ispetsc,&ao);CHKERRQ(ierr); 33*c4762a1bSJed Brown ierr = AOView(ao,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 34*c4762a1bSJed Brown 35*c4762a1bSJed Brown ierr = AOPetscToApplication(ao,4,getapp);CHKERRQ(ierr); 36*c4762a1bSJed Brown ierr = AOApplicationToPetsc(ao,3,getpetsc);CHKERRQ(ierr); 37*c4762a1bSJed Brown 38*c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d] 2,1,9,7 PetscToApplication %D %D %D %D\n",rank,getapp[0],getapp[1],getapp[2],getapp[3]);CHKERRQ(ierr); 39*c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d] 0,3,4 ApplicationToPetsc %D %D %D\n",rank,getpetsc[0],getpetsc[1],getpetsc[2]);CHKERRQ(ierr); 40*c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr); 41*c4762a1bSJed Brown ierr = AODestroy(&ao);CHKERRQ(ierr); 42*c4762a1bSJed Brown 43*c4762a1bSJed Brown /* test MemoryScalable ao */ 44*c4762a1bSJed Brown /*-------------------------*/ 45*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\nTest AOCreateMemoryScalable: \n");CHKERRQ(ierr); 46*c4762a1bSJed Brown ierr = AOCreateMemoryScalableIS(isapp,ispetsc,&ao);CHKERRQ(ierr); 47*c4762a1bSJed Brown ierr = AOView(ao,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 48*c4762a1bSJed Brown 49*c4762a1bSJed Brown ierr = AOPetscToApplication(ao,4,getapp1);CHKERRQ(ierr); 50*c4762a1bSJed Brown ierr = AOApplicationToPetsc(ao,3,getpetsc1);CHKERRQ(ierr); 51*c4762a1bSJed Brown 52*c4762a1bSJed Brown /* Check accuracy */; 53*c4762a1bSJed Brown for (i=0; i<4; i++) { 54*c4762a1bSJed Brown if (getapp1[i] != getapp[i]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"getapp1 %d != getapp %d",getapp1[i],getapp[i]); 55*c4762a1bSJed Brown } 56*c4762a1bSJed Brown for (i=0; i<3; i++) { 57*c4762a1bSJed Brown if (getpetsc1[i] != getpetsc[i]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"getpetsc1 %d != getpetsc %d",getpetsc1[i],getpetsc[i]); 58*c4762a1bSJed Brown } 59*c4762a1bSJed Brown 60*c4762a1bSJed Brown ierr = AODestroy(&ao);CHKERRQ(ierr); 61*c4762a1bSJed Brown 62*c4762a1bSJed Brown /* test MemoryScalable ao: ispetsc = NULL */ 63*c4762a1bSJed Brown /*-----------------------------------------------*/ 64*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\nTest AOCreateMemoryScalable with ispetsc=NULL:\n");CHKERRQ(ierr); 65*c4762a1bSJed Brown ierr = AOCreateMemoryScalableIS(isapp,NULL,&ao);CHKERRQ(ierr); 66*c4762a1bSJed Brown 67*c4762a1bSJed Brown ierr = AOView(ao,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 68*c4762a1bSJed Brown 69*c4762a1bSJed Brown ierr = AOPetscToApplication(ao,4,getapp2);CHKERRQ(ierr); 70*c4762a1bSJed Brown ierr = AOApplicationToPetsc(ao,3,getpetsc2);CHKERRQ(ierr); 71*c4762a1bSJed Brown 72*c4762a1bSJed Brown /* Check accuracy */; 73*c4762a1bSJed Brown for (i=0; i<4; i++) { 74*c4762a1bSJed Brown if (getapp2[i] != getapp[i]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"getapp2 %d != getapp %d",getapp2[i],getapp[i]); 75*c4762a1bSJed Brown } 76*c4762a1bSJed Brown for (i=0; i<3; i++) { 77*c4762a1bSJed Brown if (getpetsc2[i] != getpetsc[i]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"getpetsc2 %d != getpetsc %d",getpetsc2[i],getpetsc[i]); 78*c4762a1bSJed Brown } 79*c4762a1bSJed Brown ierr = AODestroy(&ao);CHKERRQ(ierr); 80*c4762a1bSJed Brown 81*c4762a1bSJed Brown /* test AOCreateMemoryScalable() ao: */ 82*c4762a1bSJed Brown ierr = ISGetIndices(isapp,&app);CHKERRQ(ierr); 83*c4762a1bSJed Brown ierr = AOCreateMemoryScalable(PETSC_COMM_WORLD,n,app,NULL,&ao);CHKERRQ(ierr); 84*c4762a1bSJed Brown ierr = ISRestoreIndices(isapp,&app);CHKERRQ(ierr); 85*c4762a1bSJed Brown 86*c4762a1bSJed Brown ierr = AOPetscToApplication(ao,4,getapp4);CHKERRQ(ierr); 87*c4762a1bSJed Brown ierr = AOApplicationToPetsc(ao,3,getpetsc4);CHKERRQ(ierr); 88*c4762a1bSJed Brown 89*c4762a1bSJed Brown /* Check accuracy */; 90*c4762a1bSJed Brown for (i=0; i<4; i++) { 91*c4762a1bSJed Brown if (getapp4[i] != getapp[i]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"getapp4 %d != getapp %d",getapp4[i],getapp[i]); 92*c4762a1bSJed Brown } 93*c4762a1bSJed Brown for (i=0; i<3; i++) { 94*c4762a1bSJed Brown if (getpetsc4[i] != getpetsc[i]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"getpetsc4 %d != getpetsc %d",getpetsc4[i],getpetsc[i]); 95*c4762a1bSJed Brown } 96*c4762a1bSJed Brown ierr = AODestroy(&ao);CHKERRQ(ierr); 97*c4762a1bSJed Brown 98*c4762a1bSJed Brown /* test general API */ 99*c4762a1bSJed Brown /*------------------*/ 100*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\nTest general API: \n");CHKERRQ(ierr); 101*c4762a1bSJed Brown ierr = AOCreate(PETSC_COMM_WORLD,&ao);CHKERRQ(ierr); 102*c4762a1bSJed Brown ierr = AOSetIS(ao,isapp,ispetsc);CHKERRQ(ierr); 103*c4762a1bSJed Brown ierr = AOSetType(ao,AOMEMORYSCALABLE);CHKERRQ(ierr); 104*c4762a1bSJed Brown ierr = AOSetFromOptions(ao);CHKERRQ(ierr); 105*c4762a1bSJed Brown 106*c4762a1bSJed Brown /* ispetsc and isapp are nolonger used. */ 107*c4762a1bSJed Brown ierr = ISDestroy(&ispetsc);CHKERRQ(ierr); 108*c4762a1bSJed Brown ierr = ISDestroy(&isapp);CHKERRQ(ierr); 109*c4762a1bSJed Brown 110*c4762a1bSJed Brown ierr = AOPetscToApplication(ao,4,getapp3);CHKERRQ(ierr); 111*c4762a1bSJed Brown ierr = AOApplicationToPetsc(ao,3,getpetsc3);CHKERRQ(ierr); 112*c4762a1bSJed Brown 113*c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d] 2,1,9,7 PetscToApplication %D %D %D %D\n",rank,getapp3[0],getapp3[1],getapp3[2],getapp3[3]);CHKERRQ(ierr); 114*c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d] 0,3,4 ApplicationToPetsc %D %D %D\n",rank,getpetsc3[0],getpetsc3[1],getpetsc3[2]);CHKERRQ(ierr); 115*c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr); 116*c4762a1bSJed Brown 117*c4762a1bSJed Brown /* Check accuracy */; 118*c4762a1bSJed Brown for (i=0; i<4; i++) { 119*c4762a1bSJed Brown if (getapp3[i] != getapp[i]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"getapp3 %d != getapp %d",getapp3[i],getapp[i]); 120*c4762a1bSJed Brown } 121*c4762a1bSJed Brown for (i=0; i<3; i++) { 122*c4762a1bSJed Brown if (getpetsc3[i] != getpetsc[i]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"getpetsc3 %d != getpetsc %d",getpetsc3[i],getpetsc[i]); 123*c4762a1bSJed Brown } 124*c4762a1bSJed Brown 125*c4762a1bSJed Brown ierr = AODestroy(&ao);CHKERRQ(ierr); 126*c4762a1bSJed Brown ierr = PetscFinalize(); 127*c4762a1bSJed Brown return ierr; 128*c4762a1bSJed Brown } 129*c4762a1bSJed Brown 130*c4762a1bSJed Brown 131*c4762a1bSJed Brown 132*c4762a1bSJed Brown 133*c4762a1bSJed Brown 134*c4762a1bSJed Brown /*TEST 135*c4762a1bSJed Brown 136*c4762a1bSJed Brown test: 137*c4762a1bSJed Brown 138*c4762a1bSJed Brown test: 139*c4762a1bSJed Brown suffix: 2 140*c4762a1bSJed Brown nsize: 2 141*c4762a1bSJed Brown 142*c4762a1bSJed Brown test: 143*c4762a1bSJed Brown suffix: 3 144*c4762a1bSJed Brown nsize: 3 145*c4762a1bSJed Brown 146*c4762a1bSJed Brown test: 147*c4762a1bSJed Brown suffix: 4 148*c4762a1bSJed Brown nsize: 3 149*c4762a1bSJed Brown args: -ao_type basic 150*c4762a1bSJed Brown output_file: output/ex1_3.out 151*c4762a1bSJed Brown 152*c4762a1bSJed Brown TEST*/ 153