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