xref: /petsc/src/vec/is/ao/tests/ex5.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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 
11c4762a1bSJed Brown int main(int argc, char *argv[])
12c4762a1bSJed Brown {
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 
25*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv, (char *) 0, help));
26c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
275f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm, &rank));
28c4762a1bSJed Brown 
295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutCreate(comm, &layout));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutSetSize(layout, n_global));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutSetLocalSize(layout, PETSC_DECIDE));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutSetUp(layout));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutGetLocalSize(layout, &local_size));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutGetRange(layout, &start, &end));
35c4762a1bSJed Brown 
365f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(local_size,&app_indices));
375f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(local_size,&petsc_indices));
38c4762a1bSJed Brown   /*  Add values for local indices for usual states */
39c4762a1bSJed Brown   for (i = 0; i < local_size; ++i) {
40c4762a1bSJed Brown     app_indices[i] = start + i;
41c4762a1bSJed Brown     petsc_indices[i] = end -1 - i;
42c4762a1bSJed Brown   }
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   /* Create the AO object that maps from lexicographic ordering to Petsc Vec ordering */
455f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(comm, local_size, &app_indices[0], PETSC_COPY_VALUES, &app_is));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(comm, local_size, &petsc_indices[0], PETSC_COPY_VALUES, &petsc_is));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(AOCreate(comm, &app2petsc));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(AOSetIS(app2petsc, app_is, petsc_is));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(AOSetType(app2petsc, AOMEMORYSCALABLE));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(AOSetFromOptions(app2petsc));
515f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&app_is));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&petsc_is));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(AOView(app2petsc, PETSC_VIEWER_STDOUT_WORLD));
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   /* Test AOApplicationToPetsc */
565f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(n_loc,&ia));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(n_loc,&ia0));
58dd400576SPatrick Sanan   if (rank == 0) {
59c4762a1bSJed Brown     ia[0] = 0;
60c4762a1bSJed Brown     ia[1] = -1;
61c4762a1bSJed Brown     ia[2] = 1;
62c4762a1bSJed Brown     ia[3] = 2;
63c4762a1bSJed Brown     ia[4] = -1;
64c4762a1bSJed Brown     ia[5] = 4;
65c4762a1bSJed Brown     ia[6] = 5;
66c4762a1bSJed Brown     ia[7] = 6;
67c4762a1bSJed Brown   } else {
68c4762a1bSJed Brown     ia[0] = -1;
69c4762a1bSJed Brown     ia[1] = 8;
70c4762a1bSJed Brown     ia[2] = 9;
71c4762a1bSJed Brown     ia[3] = 10;
72c4762a1bSJed Brown     ia[4] = -1;
73c4762a1bSJed Brown     ia[5] = 12;
74c4762a1bSJed Brown     ia[6] = 13;
75c4762a1bSJed Brown     ia[7] = 14;
76c4762a1bSJed Brown   }
775f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscArraycpy(ia0,ia,n_loc));
78c4762a1bSJed Brown 
795f80ce2aSJacob Faibussowitsch   CHKERRQ(AOApplicationToPetsc(app2petsc, n_loc, ia));
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   for (i=0; i<n_loc; ++i) {
825f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"proc = %d : %" PetscInt_FMT " -> %" PetscInt_FMT " \n", rank, ia0[i], ia[i]));
83c4762a1bSJed Brown   }
845f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(AODestroy(&app2petsc));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutDestroy(&layout));
875f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(app_indices));
885f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(petsc_indices));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(ia));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(ia0));
91*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
92*b122ec5aSJacob Faibussowitsch   return 0;
93c4762a1bSJed Brown }
94c4762a1bSJed Brown 
95c4762a1bSJed Brown /*TEST
96c4762a1bSJed Brown 
97c4762a1bSJed Brown    test:
98c4762a1bSJed Brown      nsize: 2
99c4762a1bSJed Brown 
100c4762a1bSJed Brown TEST*/
101