xref: /petsc/src/vec/is/ao/tests/ex5.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
1 /*
2  Created by Huy Vo on 12/3/18.
3 */
4 static char help[] = "Test memory scalable AO.\n\n";
5 
6 #include<petsc.h>
7 #include<petscvec.h>
8 #include<petscmat.h>
9 #include<petscao.h>
10 
11 int main(int argc, char *argv[])
12 {
13   PetscInt              n_global = 16;
14   MPI_Comm              comm;
15   PetscLayout           layout;
16   PetscInt              local_size;
17   PetscInt              start, end;
18   PetscMPIInt           rank;
19   PetscInt              *app_indices,*petsc_indices,*ia,*ia0;
20   PetscInt              i;
21   AO                    app2petsc;
22   IS                    app_is, petsc_is;
23   const PetscInt        n_loc = 8;
24 
25   PetscFunctionBeginUser;
26   PetscCall(PetscInitialize(&argc, &argv, (char *) 0, help));
27   comm = PETSC_COMM_WORLD;
28   PetscCallMPI(MPI_Comm_rank(comm, &rank));
29 
30   PetscCall(PetscLayoutCreate(comm, &layout));
31   PetscCall(PetscLayoutSetSize(layout, n_global));
32   PetscCall(PetscLayoutSetLocalSize(layout, PETSC_DECIDE));
33   PetscCall(PetscLayoutSetUp(layout));
34   PetscCall(PetscLayoutGetLocalSize(layout, &local_size));
35   PetscCall(PetscLayoutGetRange(layout, &start, &end));
36 
37   PetscCall(PetscMalloc1(local_size,&app_indices));
38   PetscCall(PetscMalloc1(local_size,&petsc_indices));
39   /*  Add values for local indices for usual states */
40   for (i = 0; i < local_size; ++i) {
41     app_indices[i] = start + i;
42     petsc_indices[i] = end -1 - i;
43   }
44 
45   /* Create the AO object that maps from lexicographic ordering to Petsc Vec ordering */
46   PetscCall(ISCreateGeneral(comm, local_size, &app_indices[0], PETSC_COPY_VALUES, &app_is));
47   PetscCall(ISCreateGeneral(comm, local_size, &petsc_indices[0], PETSC_COPY_VALUES, &petsc_is));
48   PetscCall(AOCreate(comm, &app2petsc));
49   PetscCall(AOSetIS(app2petsc, app_is, petsc_is));
50   PetscCall(AOSetType(app2petsc, AOMEMORYSCALABLE));
51   PetscCall(AOSetFromOptions(app2petsc));
52   PetscCall(ISDestroy(&app_is));
53   PetscCall(ISDestroy(&petsc_is));
54   PetscCall(AOView(app2petsc, PETSC_VIEWER_STDOUT_WORLD));
55 
56   /* Test AOApplicationToPetsc */
57   PetscCall(PetscMalloc1(n_loc,&ia));
58   PetscCall(PetscMalloc1(n_loc,&ia0));
59   if (rank == 0) {
60     ia[0] = 0;
61     ia[1] = -1;
62     ia[2] = 1;
63     ia[3] = 2;
64     ia[4] = -1;
65     ia[5] = 4;
66     ia[6] = 5;
67     ia[7] = 6;
68   } else {
69     ia[0] = -1;
70     ia[1] = 8;
71     ia[2] = 9;
72     ia[3] = 10;
73     ia[4] = -1;
74     ia[5] = 12;
75     ia[6] = 13;
76     ia[7] = 14;
77   }
78   PetscCall(PetscArraycpy(ia0,ia,n_loc));
79 
80   PetscCall(AOApplicationToPetsc(app2petsc, n_loc, ia));
81 
82   for (i=0; i<n_loc; ++i) {
83     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"proc = %d : %" PetscInt_FMT " -> %" PetscInt_FMT " \n", rank, ia0[i], ia[i]));
84   }
85   PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
86   PetscCall(AODestroy(&app2petsc));
87   PetscCall(PetscLayoutDestroy(&layout));
88   PetscCall(PetscFree(app_indices));
89   PetscCall(PetscFree(petsc_indices));
90   PetscCall(PetscFree(ia));
91   PetscCall(PetscFree(ia0));
92   PetscCall(PetscFinalize());
93   return 0;
94 }
95 
96 /*TEST
97 
98    test:
99      nsize: 2
100 
101 TEST*/
102