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