xref: /petsc/src/vec/is/ao/tests/ex5.c (revision 09b68a49ed2854d1e4985cc2aa6af33c7c4e69b3)
1c4762a1bSJed Brown /*
2c4762a1bSJed Brown  Created by Huy Vo on 12/3/18.
3c4762a1bSJed Brown */
4c4762a1bSJed Brown static char help[] = "Test memory scalable AO.\n\n";
5c4762a1bSJed Brown 
6c4762a1bSJed Brown #include <petsc.h>
7c4762a1bSJed Brown #include <petscvec.h>
8c4762a1bSJed Brown #include <petscmat.h>
9c4762a1bSJed Brown #include <petscao.h>
10c4762a1bSJed Brown 
main(int argc,char * argv[])11d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[])
12d71ae5a4SJacob Faibussowitsch {
13c4762a1bSJed Brown   PetscInt       n_global = 16;
14c4762a1bSJed Brown   MPI_Comm       comm;
15c4762a1bSJed Brown   PetscLayout    layout;
16c4762a1bSJed Brown   PetscInt       local_size;
17c4762a1bSJed Brown   PetscInt       start, end;
18c4762a1bSJed Brown   PetscMPIInt    rank;
19c4762a1bSJed Brown   PetscInt      *app_indices, *petsc_indices, *ia, *ia0;
20c4762a1bSJed Brown   PetscInt       i;
21c4762a1bSJed Brown   AO             app2petsc;
22c4762a1bSJed Brown   IS             app_is, petsc_is;
23c4762a1bSJed Brown   const PetscInt n_loc = 8;
24c4762a1bSJed Brown 
25327415f7SBarry Smith   PetscFunctionBeginUser;
26c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
27c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
289566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
29c4762a1bSJed Brown 
309566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &layout));
319566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetSize(layout, n_global));
329566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(layout, PETSC_DECIDE));
339566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(layout));
349566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetLocalSize(layout, &local_size));
359566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(layout, &start, &end));
36c4762a1bSJed Brown 
379566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(local_size, &app_indices));
389566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(local_size, &petsc_indices));
39c4762a1bSJed Brown   /*  Add values for local indices for usual states */
40c4762a1bSJed Brown   for (i = 0; i < local_size; ++i) {
41c4762a1bSJed Brown     app_indices[i]   = start + i;
42c4762a1bSJed Brown     petsc_indices[i] = end - 1 - i;
43c4762a1bSJed Brown   }
44c4762a1bSJed Brown 
45*f0b74427SPierre Jolivet   /* Create the AO object that maps from lexicographic ordering to PETSc Vec ordering */
469566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(comm, local_size, &app_indices[0], PETSC_COPY_VALUES, &app_is));
479566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(comm, local_size, &petsc_indices[0], PETSC_COPY_VALUES, &petsc_is));
489566063dSJacob Faibussowitsch   PetscCall(AOCreate(comm, &app2petsc));
499566063dSJacob Faibussowitsch   PetscCall(AOSetIS(app2petsc, app_is, petsc_is));
509566063dSJacob Faibussowitsch   PetscCall(AOSetType(app2petsc, AOMEMORYSCALABLE));
519566063dSJacob Faibussowitsch   PetscCall(AOSetFromOptions(app2petsc));
529566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&app_is));
539566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&petsc_is));
549566063dSJacob Faibussowitsch   PetscCall(AOView(app2petsc, PETSC_VIEWER_STDOUT_WORLD));
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   /* Test AOApplicationToPetsc */
579566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n_loc, &ia));
589566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n_loc, &ia0));
59dd400576SPatrick Sanan   if (rank == 0) {
60c4762a1bSJed Brown     ia[0] = 0;
61c4762a1bSJed Brown     ia[1] = -1;
62c4762a1bSJed Brown     ia[2] = 1;
63c4762a1bSJed Brown     ia[3] = 2;
64c4762a1bSJed Brown     ia[4] = -1;
65c4762a1bSJed Brown     ia[5] = 4;
66c4762a1bSJed Brown     ia[6] = 5;
67c4762a1bSJed Brown     ia[7] = 6;
68c4762a1bSJed Brown   } else {
69c4762a1bSJed Brown     ia[0] = -1;
70c4762a1bSJed Brown     ia[1] = 8;
71c4762a1bSJed Brown     ia[2] = 9;
72c4762a1bSJed Brown     ia[3] = 10;
73c4762a1bSJed Brown     ia[4] = -1;
74c4762a1bSJed Brown     ia[5] = 12;
75c4762a1bSJed Brown     ia[6] = 13;
76c4762a1bSJed Brown     ia[7] = 14;
77c4762a1bSJed Brown   }
789566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(ia0, ia, n_loc));
79c4762a1bSJed Brown 
809566063dSJacob Faibussowitsch   PetscCall(AOApplicationToPetsc(app2petsc, n_loc, ia));
81c4762a1bSJed Brown 
8248a46eb9SPierre Jolivet   for (i = 0; i < n_loc; ++i) PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "proc = %d : %" PetscInt_FMT " -> %" PetscInt_FMT " \n", rank, ia0[i], ia[i]));
839566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
849566063dSJacob Faibussowitsch   PetscCall(AODestroy(&app2petsc));
859566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&layout));
869566063dSJacob Faibussowitsch   PetscCall(PetscFree(app_indices));
879566063dSJacob Faibussowitsch   PetscCall(PetscFree(petsc_indices));
889566063dSJacob Faibussowitsch   PetscCall(PetscFree(ia));
899566063dSJacob Faibussowitsch   PetscCall(PetscFree(ia0));
909566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
91b122ec5aSJacob Faibussowitsch   return 0;
92c4762a1bSJed Brown }
93c4762a1bSJed Brown 
94c4762a1bSJed Brown /*TEST
95c4762a1bSJed Brown 
96c4762a1bSJed Brown    test:
97c4762a1bSJed Brown      nsize: 2
98c4762a1bSJed Brown 
99c4762a1bSJed Brown TEST*/
100