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