xref: /petsc/src/vec/is/ao/tests/ex7.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Demonstrates constructing an application ordering.\n\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscao.h>
5*c4762a1bSJed Brown #include <petscviewer.h>
6*c4762a1bSJed Brown 
7*c4762a1bSJed Brown int main(int argc,char **argv)
8*c4762a1bSJed Brown {
9*c4762a1bSJed Brown   PetscInt       n = 5;
10*c4762a1bSJed Brown   PetscErrorCode ierr;
11*c4762a1bSJed Brown   PetscMPIInt    rank,size;
12*c4762a1bSJed Brown   IS             ispetsc,isapp;
13*c4762a1bSJed Brown   AO             ao;
14*c4762a1bSJed Brown 
15*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
16*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr);
17*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
18*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
19*c4762a1bSJed Brown 
20*c4762a1bSJed Brown   /* create the index sets */
21*c4762a1bSJed Brown   ierr = ISCreateStride(PETSC_COMM_WORLD,n,rank,size,&ispetsc);CHKERRQ(ierr);
22*c4762a1bSJed Brown   ierr = ISCreateStride(PETSC_COMM_WORLD,n,n*rank,1,&isapp);CHKERRQ(ierr);
23*c4762a1bSJed Brown 
24*c4762a1bSJed Brown   /* create the application ordering */
25*c4762a1bSJed Brown   ierr = AOCreateBasicIS(isapp,ispetsc,&ao);CHKERRQ(ierr);
26*c4762a1bSJed Brown 
27*c4762a1bSJed Brown 
28*c4762a1bSJed Brown   ierr = AOView(ao,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
29*c4762a1bSJed Brown 
30*c4762a1bSJed Brown   ierr = ISView(ispetsc,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
31*c4762a1bSJed Brown   ierr = ISView(isapp,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
32*c4762a1bSJed Brown   ierr = AOPetscToApplicationIS(ao,ispetsc);CHKERRQ(ierr);
33*c4762a1bSJed Brown   ierr = ISView(isapp,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
34*c4762a1bSJed Brown   ierr = ISView(ispetsc,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
35*c4762a1bSJed Brown 
36*c4762a1bSJed Brown 
37*c4762a1bSJed Brown   ierr = ISDestroy(&ispetsc);CHKERRQ(ierr);
38*c4762a1bSJed Brown   ierr = ISDestroy(&isapp);CHKERRQ(ierr);
39*c4762a1bSJed Brown 
40*c4762a1bSJed Brown   ierr = AODestroy(&ao);CHKERRQ(ierr);
41*c4762a1bSJed Brown   ierr = PetscFinalize();
42*c4762a1bSJed Brown   return ierr;
43*c4762a1bSJed Brown }
44*c4762a1bSJed Brown 
45*c4762a1bSJed Brown /*TEST
46*c4762a1bSJed Brown 
47*c4762a1bSJed Brown    test:
48*c4762a1bSJed Brown       nsize: 2
49*c4762a1bSJed Brown 
50*c4762a1bSJed Brown TEST*/
51