xref: /petsc/src/vec/is/ao/tests/ex1.c (revision 08401ef684002a709c6d3db98a0c9f54a8bcf1ec)
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 
219566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
249566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   /* create the index sets */
279566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_WORLD,n,rank,size,&isapp));
289566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_WORLD,n,n*rank,1,&ispetsc)); /* natural numbering */
29c4762a1bSJed Brown 
30c4762a1bSJed Brown   /* create the application ordering */
319566063dSJacob Faibussowitsch   PetscCall(AOCreateBasicIS(isapp,ispetsc,&ao));
329566063dSJacob Faibussowitsch   PetscCall(AOView(ao,PETSC_VIEWER_STDOUT_WORLD));
33c4762a1bSJed Brown 
349566063dSJacob Faibussowitsch   PetscCall(AOPetscToApplication(ao,4,getapp));
359566063dSJacob Faibussowitsch   PetscCall(AOApplicationToPetsc(ao,3,getpetsc));
36c4762a1bSJed Brown 
379566063dSJacob 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]));
389566063dSJacob 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]));
399566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
409566063dSJacob Faibussowitsch   PetscCall(AODestroy(&ao));
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   /* test MemoryScalable ao */
43c4762a1bSJed Brown   /*-------------------------*/
449566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nTest AOCreateMemoryScalable: \n"));
459566063dSJacob Faibussowitsch   PetscCall(AOCreateMemoryScalableIS(isapp,ispetsc,&ao));
469566063dSJacob Faibussowitsch   PetscCall(AOView(ao,PETSC_VIEWER_STDOUT_WORLD));
47c4762a1bSJed Brown 
489566063dSJacob Faibussowitsch   PetscCall(AOPetscToApplication(ao,4,getapp1));
499566063dSJacob Faibussowitsch   PetscCall(AOApplicationToPetsc(ao,3,getpetsc1));
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   /* Check accuracy */;
52c4762a1bSJed Brown   for (i=0; i<4; i++) {
53*08401ef6SPierre Jolivet     PetscCheck(getapp1[i] == getapp[i],PETSC_COMM_SELF,PETSC_ERR_USER,"getapp1 %" PetscInt_FMT " != getapp %" PetscInt_FMT,getapp1[i],getapp[i]);
54c4762a1bSJed Brown   }
55c4762a1bSJed Brown   for (i=0; i<3; i++) {
56*08401ef6SPierre Jolivet     PetscCheck(getpetsc1[i] == getpetsc[i],PETSC_COMM_SELF,PETSC_ERR_USER,"getpetsc1 %" PetscInt_FMT " != getpetsc %" PetscInt_FMT,getpetsc1[i],getpetsc[i]);
57c4762a1bSJed Brown   }
58c4762a1bSJed Brown 
599566063dSJacob Faibussowitsch   PetscCall(AODestroy(&ao));
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   /* test MemoryScalable ao: ispetsc = NULL */
62c4762a1bSJed Brown   /*-----------------------------------------------*/
639566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nTest AOCreateMemoryScalable with ispetsc=NULL:\n"));
649566063dSJacob Faibussowitsch   PetscCall(AOCreateMemoryScalableIS(isapp,NULL,&ao));
65c4762a1bSJed Brown 
669566063dSJacob Faibussowitsch   PetscCall(AOView(ao,PETSC_VIEWER_STDOUT_WORLD));
67c4762a1bSJed Brown 
689566063dSJacob Faibussowitsch   PetscCall(AOPetscToApplication(ao,4,getapp2));
699566063dSJacob Faibussowitsch   PetscCall(AOApplicationToPetsc(ao,3,getpetsc2));
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   /* Check accuracy */;
72c4762a1bSJed Brown   for (i=0; i<4; i++) {
73*08401ef6SPierre Jolivet     PetscCheck(getapp2[i] == getapp[i],PETSC_COMM_SELF,PETSC_ERR_USER,"getapp2 %" PetscInt_FMT " != getapp %" PetscInt_FMT,getapp2[i],getapp[i]);
74c4762a1bSJed Brown   }
75c4762a1bSJed Brown   for (i=0; i<3; i++) {
76*08401ef6SPierre Jolivet     PetscCheck(getpetsc2[i] == getpetsc[i],PETSC_COMM_SELF,PETSC_ERR_USER,"getpetsc2 %" PetscInt_FMT " != getpetsc %" PetscInt_FMT,getpetsc2[i],getpetsc[i]);
77c4762a1bSJed Brown   }
789566063dSJacob Faibussowitsch   PetscCall(AODestroy(&ao));
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   /* test AOCreateMemoryScalable() ao: */
819566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isapp,&app));
829566063dSJacob Faibussowitsch   PetscCall(AOCreateMemoryScalable(PETSC_COMM_WORLD,n,app,NULL,&ao));
839566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isapp,&app));
84c4762a1bSJed Brown 
859566063dSJacob Faibussowitsch   PetscCall(AOPetscToApplication(ao,4,getapp4));
869566063dSJacob Faibussowitsch   PetscCall(AOApplicationToPetsc(ao,3,getpetsc4));
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   /* Check accuracy */;
89c4762a1bSJed Brown   for (i=0; i<4; i++) {
90*08401ef6SPierre Jolivet     PetscCheck(getapp4[i] == getapp[i],PETSC_COMM_SELF,PETSC_ERR_USER,"getapp4 %" PetscInt_FMT " != getapp %" PetscInt_FMT,getapp4[i],getapp[i]);
91c4762a1bSJed Brown   }
92c4762a1bSJed Brown   for (i=0; i<3; i++) {
93*08401ef6SPierre Jolivet     PetscCheck(getpetsc4[i] == getpetsc[i],PETSC_COMM_SELF,PETSC_ERR_USER,"getpetsc4 %" PetscInt_FMT " != getpetsc %" PetscInt_FMT,getpetsc4[i],getpetsc[i]);
94c4762a1bSJed Brown   }
959566063dSJacob Faibussowitsch   PetscCall(AODestroy(&ao));
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   /* test general API */
98c4762a1bSJed Brown   /*------------------*/
999566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nTest general API: \n"));
1009566063dSJacob Faibussowitsch   PetscCall(AOCreate(PETSC_COMM_WORLD,&ao));
1019566063dSJacob Faibussowitsch   PetscCall(AOSetIS(ao,isapp,ispetsc));
1029566063dSJacob Faibussowitsch   PetscCall(AOSetType(ao,AOMEMORYSCALABLE));
1039566063dSJacob Faibussowitsch   PetscCall(AOSetFromOptions(ao));
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   /* ispetsc and isapp are nolonger used. */
1069566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ispetsc));
1079566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isapp));
108c4762a1bSJed Brown 
1099566063dSJacob Faibussowitsch   PetscCall(AOPetscToApplication(ao,4,getapp3));
1109566063dSJacob Faibussowitsch   PetscCall(AOApplicationToPetsc(ao,3,getpetsc3));
111c4762a1bSJed Brown 
1129566063dSJacob 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]));
1139566063dSJacob 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]));
1149566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   /* Check accuracy */;
117c4762a1bSJed Brown   for (i=0; i<4; i++) {
118*08401ef6SPierre Jolivet     PetscCheck(getapp3[i] == getapp[i],PETSC_COMM_SELF,PETSC_ERR_USER,"getapp3 %" PetscInt_FMT " != getapp %" PetscInt_FMT,getapp3[i],getapp[i]);
119c4762a1bSJed Brown   }
120c4762a1bSJed Brown   for (i=0; i<3; i++) {
121*08401ef6SPierre Jolivet     PetscCheck(getpetsc3[i] == getpetsc[i],PETSC_COMM_SELF,PETSC_ERR_USER,"getpetsc3 %" PetscInt_FMT " != getpetsc %" PetscInt_FMT,getpetsc3[i],getpetsc[i]);
122c4762a1bSJed Brown   }
123c4762a1bSJed Brown 
1249566063dSJacob Faibussowitsch   PetscCall(AODestroy(&ao));
1259566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
126b122ec5aSJacob Faibussowitsch   return 0;
127c4762a1bSJed Brown }
128c4762a1bSJed Brown 
129c4762a1bSJed Brown /*TEST
130c4762a1bSJed Brown 
131c4762a1bSJed Brown    test:
132c4762a1bSJed Brown 
133c4762a1bSJed Brown    test:
134c4762a1bSJed Brown       suffix: 2
135c4762a1bSJed Brown       nsize: 2
136c4762a1bSJed Brown 
137c4762a1bSJed Brown    test:
138c4762a1bSJed Brown       suffix: 3
139c4762a1bSJed Brown       nsize: 3
140c4762a1bSJed Brown 
141c4762a1bSJed Brown    test:
142c4762a1bSJed Brown       suffix: 4
143c4762a1bSJed Brown       nsize: 3
144c4762a1bSJed Brown       args: -ao_type basic
145c4762a1bSJed Brown       output_file: output/ex1_3.out
146c4762a1bSJed Brown 
147c4762a1bSJed Brown TEST*/
148