1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Demonstrates constructing an application ordering.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscsys.h> 5c4762a1bSJed Brown #include <petscao.h> 6c4762a1bSJed Brown #include <petscviewer.h> 7c4762a1bSJed Brown 8c4762a1bSJed Brown int main(int argc,char **argv) 9c4762a1bSJed Brown { 10c4762a1bSJed Brown PetscInt i,n = 5; 11c4762a1bSJed Brown PetscInt getpetsc[] = {0,3,4},getapp[] = {2,1,9,7}; 12c4762a1bSJed Brown PetscInt getpetsc1[] = {0,3,4},getapp1[] = {2,1,9,7}; 13c4762a1bSJed Brown PetscInt getpetsc2[] = {0,3,4},getapp2[] = {2,1,9,7}; 14c4762a1bSJed Brown PetscInt getpetsc3[] = {0,3,4},getapp3[] = {2,1,9,7}; 15c4762a1bSJed Brown PetscInt getpetsc4[] = {0,3,4},getapp4[] = {2,1,9,7}; 16c4762a1bSJed Brown PetscMPIInt rank,size; 17c4762a1bSJed Brown IS ispetsc,isapp; 18c4762a1bSJed Brown AO ao; 19c4762a1bSJed Brown const PetscInt *app; 20c4762a1bSJed Brown 21*327415f7SBarry Smith PetscFunctionBeginUser; 229566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 239566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 26c4762a1bSJed Brown 27c4762a1bSJed Brown /* create the index sets */ 289566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,n,rank,size,&isapp)); 299566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,n,n*rank,1,&ispetsc)); /* natural numbering */ 30c4762a1bSJed Brown 31c4762a1bSJed Brown /* create the application ordering */ 329566063dSJacob Faibussowitsch PetscCall(AOCreateBasicIS(isapp,ispetsc,&ao)); 339566063dSJacob Faibussowitsch PetscCall(AOView(ao,PETSC_VIEWER_STDOUT_WORLD)); 34c4762a1bSJed Brown 359566063dSJacob Faibussowitsch PetscCall(AOPetscToApplication(ao,4,getapp)); 369566063dSJacob Faibussowitsch PetscCall(AOApplicationToPetsc(ao,3,getpetsc)); 37c4762a1bSJed Brown 389566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d] 2,1,9,7 PetscToApplication %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n",rank,getapp[0],getapp[1],getapp[2],getapp[3])); 399566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d] 0,3,4 ApplicationToPetsc %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n",rank,getpetsc[0],getpetsc[1],getpetsc[2])); 409566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 419566063dSJacob Faibussowitsch PetscCall(AODestroy(&ao)); 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* test MemoryScalable ao */ 44c4762a1bSJed Brown /*-------------------------*/ 459566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nTest AOCreateMemoryScalable: \n")); 469566063dSJacob Faibussowitsch PetscCall(AOCreateMemoryScalableIS(isapp,ispetsc,&ao)); 479566063dSJacob Faibussowitsch PetscCall(AOView(ao,PETSC_VIEWER_STDOUT_WORLD)); 48c4762a1bSJed Brown 499566063dSJacob Faibussowitsch PetscCall(AOPetscToApplication(ao,4,getapp1)); 509566063dSJacob Faibussowitsch PetscCall(AOApplicationToPetsc(ao,3,getpetsc1)); 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* Check accuracy */; 53c4762a1bSJed Brown for (i=0; i<4; i++) { 5408401ef6SPierre Jolivet PetscCheck(getapp1[i] == getapp[i],PETSC_COMM_SELF,PETSC_ERR_USER,"getapp1 %" PetscInt_FMT " != getapp %" PetscInt_FMT,getapp1[i],getapp[i]); 55c4762a1bSJed Brown } 56c4762a1bSJed Brown for (i=0; i<3; i++) { 5708401ef6SPierre Jolivet PetscCheck(getpetsc1[i] == getpetsc[i],PETSC_COMM_SELF,PETSC_ERR_USER,"getpetsc1 %" PetscInt_FMT " != getpetsc %" PetscInt_FMT,getpetsc1[i],getpetsc[i]); 58c4762a1bSJed Brown } 59c4762a1bSJed Brown 609566063dSJacob Faibussowitsch PetscCall(AODestroy(&ao)); 61c4762a1bSJed Brown 62c4762a1bSJed Brown /* test MemoryScalable ao: ispetsc = NULL */ 63c4762a1bSJed Brown /*-----------------------------------------------*/ 649566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nTest AOCreateMemoryScalable with ispetsc=NULL:\n")); 659566063dSJacob Faibussowitsch PetscCall(AOCreateMemoryScalableIS(isapp,NULL,&ao)); 66c4762a1bSJed Brown 679566063dSJacob Faibussowitsch PetscCall(AOView(ao,PETSC_VIEWER_STDOUT_WORLD)); 68c4762a1bSJed Brown 699566063dSJacob Faibussowitsch PetscCall(AOPetscToApplication(ao,4,getapp2)); 709566063dSJacob Faibussowitsch PetscCall(AOApplicationToPetsc(ao,3,getpetsc2)); 71c4762a1bSJed Brown 72c4762a1bSJed Brown /* Check accuracy */; 73c4762a1bSJed Brown for (i=0; i<4; i++) { 7408401ef6SPierre Jolivet PetscCheck(getapp2[i] == getapp[i],PETSC_COMM_SELF,PETSC_ERR_USER,"getapp2 %" PetscInt_FMT " != getapp %" PetscInt_FMT,getapp2[i],getapp[i]); 75c4762a1bSJed Brown } 76c4762a1bSJed Brown for (i=0; i<3; i++) { 7708401ef6SPierre Jolivet PetscCheck(getpetsc2[i] == getpetsc[i],PETSC_COMM_SELF,PETSC_ERR_USER,"getpetsc2 %" PetscInt_FMT " != getpetsc %" PetscInt_FMT,getpetsc2[i],getpetsc[i]); 78c4762a1bSJed Brown } 799566063dSJacob Faibussowitsch PetscCall(AODestroy(&ao)); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* test AOCreateMemoryScalable() ao: */ 829566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isapp,&app)); 839566063dSJacob Faibussowitsch PetscCall(AOCreateMemoryScalable(PETSC_COMM_WORLD,n,app,NULL,&ao)); 849566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isapp,&app)); 85c4762a1bSJed Brown 869566063dSJacob Faibussowitsch PetscCall(AOPetscToApplication(ao,4,getapp4)); 879566063dSJacob Faibussowitsch PetscCall(AOApplicationToPetsc(ao,3,getpetsc4)); 88c4762a1bSJed Brown 89c4762a1bSJed Brown /* Check accuracy */; 90c4762a1bSJed Brown for (i=0; i<4; i++) { 9108401ef6SPierre Jolivet PetscCheck(getapp4[i] == getapp[i],PETSC_COMM_SELF,PETSC_ERR_USER,"getapp4 %" PetscInt_FMT " != getapp %" PetscInt_FMT,getapp4[i],getapp[i]); 92c4762a1bSJed Brown } 93c4762a1bSJed Brown for (i=0; i<3; i++) { 9408401ef6SPierre Jolivet PetscCheck(getpetsc4[i] == getpetsc[i],PETSC_COMM_SELF,PETSC_ERR_USER,"getpetsc4 %" PetscInt_FMT " != getpetsc %" PetscInt_FMT,getpetsc4[i],getpetsc[i]); 95c4762a1bSJed Brown } 969566063dSJacob Faibussowitsch PetscCall(AODestroy(&ao)); 97c4762a1bSJed Brown 98c4762a1bSJed Brown /* test general API */ 99c4762a1bSJed Brown /*------------------*/ 1009566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nTest general API: \n")); 1019566063dSJacob Faibussowitsch PetscCall(AOCreate(PETSC_COMM_WORLD,&ao)); 1029566063dSJacob Faibussowitsch PetscCall(AOSetIS(ao,isapp,ispetsc)); 1039566063dSJacob Faibussowitsch PetscCall(AOSetType(ao,AOMEMORYSCALABLE)); 1049566063dSJacob Faibussowitsch PetscCall(AOSetFromOptions(ao)); 105c4762a1bSJed Brown 106c4762a1bSJed Brown /* ispetsc and isapp are nolonger used. */ 1079566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ispetsc)); 1089566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isapp)); 109c4762a1bSJed Brown 1109566063dSJacob Faibussowitsch PetscCall(AOPetscToApplication(ao,4,getapp3)); 1119566063dSJacob Faibussowitsch PetscCall(AOApplicationToPetsc(ao,3,getpetsc3)); 112c4762a1bSJed Brown 1139566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d] 2,1,9,7 PetscToApplication %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n",rank,getapp3[0],getapp3[1],getapp3[2],getapp3[3])); 1149566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d] 0,3,4 ApplicationToPetsc %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n",rank,getpetsc3[0],getpetsc3[1],getpetsc3[2])); 1159566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 116c4762a1bSJed Brown 117c4762a1bSJed Brown /* Check accuracy */; 118c4762a1bSJed Brown for (i=0; i<4; i++) { 11908401ef6SPierre Jolivet PetscCheck(getapp3[i] == getapp[i],PETSC_COMM_SELF,PETSC_ERR_USER,"getapp3 %" PetscInt_FMT " != getapp %" PetscInt_FMT,getapp3[i],getapp[i]); 120c4762a1bSJed Brown } 121c4762a1bSJed Brown for (i=0; i<3; i++) { 12208401ef6SPierre Jolivet PetscCheck(getpetsc3[i] == getpetsc[i],PETSC_COMM_SELF,PETSC_ERR_USER,"getpetsc3 %" PetscInt_FMT " != getpetsc %" PetscInt_FMT,getpetsc3[i],getpetsc[i]); 123c4762a1bSJed Brown } 124c4762a1bSJed Brown 1259566063dSJacob Faibussowitsch PetscCall(AODestroy(&ao)); 1269566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 127b122ec5aSJacob Faibussowitsch return 0; 128c4762a1bSJed Brown } 129c4762a1bSJed Brown 130c4762a1bSJed Brown /*TEST 131c4762a1bSJed Brown 132c4762a1bSJed Brown test: 133c4762a1bSJed Brown 134c4762a1bSJed Brown test: 135c4762a1bSJed Brown suffix: 2 136c4762a1bSJed Brown nsize: 2 137c4762a1bSJed Brown 138c4762a1bSJed Brown test: 139c4762a1bSJed Brown suffix: 3 140c4762a1bSJed Brown nsize: 3 141c4762a1bSJed Brown 142c4762a1bSJed Brown test: 143c4762a1bSJed Brown suffix: 4 144c4762a1bSJed Brown nsize: 3 145c4762a1bSJed Brown args: -ao_type basic 146c4762a1bSJed Brown output_file: output/ex1_3.out 147c4762a1bSJed Brown 148c4762a1bSJed Brown TEST*/ 149