xref: /petsc/src/vec/is/ao/tests/ex1.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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   PetscErrorCode ierr;
11c4762a1bSJed Brown   PetscInt       i,n = 5;
12c4762a1bSJed Brown   PetscInt       getpetsc[]  = {0,3,4},getapp[]  = {2,1,9,7};
13c4762a1bSJed Brown   PetscInt       getpetsc1[] = {0,3,4},getapp1[] = {2,1,9,7};
14c4762a1bSJed Brown   PetscInt       getpetsc2[] = {0,3,4},getapp2[] = {2,1,9,7};
15c4762a1bSJed Brown   PetscInt       getpetsc3[] = {0,3,4},getapp3[] = {2,1,9,7};
16c4762a1bSJed Brown   PetscInt       getpetsc4[] = {0,3,4},getapp4[] = {2,1,9,7};
17c4762a1bSJed Brown   PetscMPIInt    rank,size;
18c4762a1bSJed Brown   IS             ispetsc,isapp;
19c4762a1bSJed Brown   AO             ao;
20c4762a1bSJed Brown   const PetscInt *app;
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
23c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr);
24ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
25ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
26c4762a1bSJed Brown 
27c4762a1bSJed Brown   /* create the index sets */
28c4762a1bSJed Brown   ierr = ISCreateStride(PETSC_COMM_WORLD,n,rank,size,&isapp);CHKERRQ(ierr);
29c4762a1bSJed Brown   ierr = ISCreateStride(PETSC_COMM_WORLD,n,n*rank,1,&ispetsc);CHKERRQ(ierr); /* natural numbering */
30c4762a1bSJed Brown 
31c4762a1bSJed Brown   /* create the application ordering */
32c4762a1bSJed Brown   ierr = AOCreateBasicIS(isapp,ispetsc,&ao);CHKERRQ(ierr);
33c4762a1bSJed Brown   ierr = AOView(ao,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   ierr = AOPetscToApplication(ao,4,getapp);CHKERRQ(ierr);
36c4762a1bSJed Brown   ierr = AOApplicationToPetsc(ao,3,getpetsc);CHKERRQ(ierr);
37c4762a1bSJed Brown 
3805d37114SJacob Faibussowitsch   ierr = 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]);CHKERRQ(ierr);
3905d37114SJacob Faibussowitsch   ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d] 0,3,4 ApplicationToPetsc %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n",rank,getpetsc[0],getpetsc[1],getpetsc[2]);CHKERRQ(ierr);
40c4762a1bSJed Brown   ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr);
41c4762a1bSJed Brown   ierr = AODestroy(&ao);CHKERRQ(ierr);
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   /* test MemoryScalable ao */
44c4762a1bSJed Brown   /*-------------------------*/
45c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"\nTest AOCreateMemoryScalable: \n");CHKERRQ(ierr);
46c4762a1bSJed Brown   ierr = AOCreateMemoryScalableIS(isapp,ispetsc,&ao);CHKERRQ(ierr);
47c4762a1bSJed Brown   ierr = AOView(ao,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   ierr = AOPetscToApplication(ao,4,getapp1);CHKERRQ(ierr);
50c4762a1bSJed Brown   ierr = AOApplicationToPetsc(ao,3,getpetsc1);CHKERRQ(ierr);
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   /* Check accuracy */;
53c4762a1bSJed Brown   for (i=0; i<4; i++) {
54*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(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++) {
57*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(getpetsc1[i] != getpetsc[i],PETSC_COMM_SELF,PETSC_ERR_USER,"getpetsc1 %" PetscInt_FMT " != getpetsc %" PetscInt_FMT,getpetsc1[i],getpetsc[i]);
58c4762a1bSJed Brown   }
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   ierr = AODestroy(&ao);CHKERRQ(ierr);
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   /* test MemoryScalable ao: ispetsc = NULL */
63c4762a1bSJed Brown   /*-----------------------------------------------*/
64c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"\nTest AOCreateMemoryScalable with ispetsc=NULL:\n");CHKERRQ(ierr);
65c4762a1bSJed Brown   ierr = AOCreateMemoryScalableIS(isapp,NULL,&ao);CHKERRQ(ierr);
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   ierr = AOView(ao,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   ierr = AOPetscToApplication(ao,4,getapp2);CHKERRQ(ierr);
70c4762a1bSJed Brown   ierr = AOApplicationToPetsc(ao,3,getpetsc2);CHKERRQ(ierr);
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   /* Check accuracy */;
73c4762a1bSJed Brown   for (i=0; i<4; i++) {
74*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(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++) {
77*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(getpetsc2[i] != getpetsc[i],PETSC_COMM_SELF,PETSC_ERR_USER,"getpetsc2 %" PetscInt_FMT " != getpetsc %" PetscInt_FMT,getpetsc2[i],getpetsc[i]);
78c4762a1bSJed Brown   }
79c4762a1bSJed Brown   ierr = AODestroy(&ao);CHKERRQ(ierr);
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   /* test AOCreateMemoryScalable() ao: */
82c4762a1bSJed Brown   ierr = ISGetIndices(isapp,&app);CHKERRQ(ierr);
83c4762a1bSJed Brown   ierr = AOCreateMemoryScalable(PETSC_COMM_WORLD,n,app,NULL,&ao);CHKERRQ(ierr);
84c4762a1bSJed Brown   ierr = ISRestoreIndices(isapp,&app);CHKERRQ(ierr);
85c4762a1bSJed Brown 
86c4762a1bSJed Brown   ierr = AOPetscToApplication(ao,4,getapp4);CHKERRQ(ierr);
87c4762a1bSJed Brown   ierr = AOApplicationToPetsc(ao,3,getpetsc4);CHKERRQ(ierr);
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   /* Check accuracy */;
90c4762a1bSJed Brown   for (i=0; i<4; i++) {
91*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(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++) {
94*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(getpetsc4[i] != getpetsc[i],PETSC_COMM_SELF,PETSC_ERR_USER,"getpetsc4 %" PetscInt_FMT " != getpetsc %" PetscInt_FMT,getpetsc4[i],getpetsc[i]);
95c4762a1bSJed Brown   }
96c4762a1bSJed Brown   ierr = AODestroy(&ao);CHKERRQ(ierr);
97c4762a1bSJed Brown 
98c4762a1bSJed Brown   /* test general API */
99c4762a1bSJed Brown   /*------------------*/
100c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"\nTest general API: \n");CHKERRQ(ierr);
101c4762a1bSJed Brown   ierr = AOCreate(PETSC_COMM_WORLD,&ao);CHKERRQ(ierr);
102c4762a1bSJed Brown   ierr = AOSetIS(ao,isapp,ispetsc);CHKERRQ(ierr);
103c4762a1bSJed Brown   ierr = AOSetType(ao,AOMEMORYSCALABLE);CHKERRQ(ierr);
104c4762a1bSJed Brown   ierr = AOSetFromOptions(ao);CHKERRQ(ierr);
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   /* ispetsc and isapp are nolonger used. */
107c4762a1bSJed Brown   ierr = ISDestroy(&ispetsc);CHKERRQ(ierr);
108c4762a1bSJed Brown   ierr = ISDestroy(&isapp);CHKERRQ(ierr);
109c4762a1bSJed Brown 
110c4762a1bSJed Brown   ierr = AOPetscToApplication(ao,4,getapp3);CHKERRQ(ierr);
111c4762a1bSJed Brown   ierr = AOApplicationToPetsc(ao,3,getpetsc3);CHKERRQ(ierr);
112c4762a1bSJed Brown 
11305d37114SJacob Faibussowitsch   ierr = 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]);CHKERRQ(ierr);
11405d37114SJacob Faibussowitsch   ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d] 0,3,4 ApplicationToPetsc %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n",rank,getpetsc3[0],getpetsc3[1],getpetsc3[2]);CHKERRQ(ierr);
115c4762a1bSJed Brown   ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr);
116c4762a1bSJed Brown 
117c4762a1bSJed Brown   /* Check accuracy */;
118c4762a1bSJed Brown   for (i=0; i<4; i++) {
119*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(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++) {
122*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(getpetsc3[i] != getpetsc[i],PETSC_COMM_SELF,PETSC_ERR_USER,"getpetsc3 %" PetscInt_FMT " != getpetsc %" PetscInt_FMT,getpetsc3[i],getpetsc[i]);
123c4762a1bSJed Brown   }
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   ierr = AODestroy(&ao);CHKERRQ(ierr);
126c4762a1bSJed Brown   ierr = PetscFinalize();
127c4762a1bSJed Brown   return ierr;
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