xref: /petsc/src/sys/tests/ex21.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
1 
2 static char help[] = "Tests PetscTreeProcess()";
3 
4 #include <petscsys.h>
5 
6 /*
7                           2              6
8                     1         4
9                     5
10 */
11 int main(int argc,char **argv)
12 {
13   PetscErrorCode ierr;
14   PetscInt       n          = 7,cnt = 0,i,j;
15   PetscBool      mask[]     = {PETSC_TRUE,PETSC_FALSE,PETSC_FALSE,PETSC_TRUE,PETSC_FALSE,PETSC_FALSE,PETSC_FALSE};
16   PetscInt       parentId[] = {-1,         2,         0,         -1,         2,         1,         0};
17   PetscInt       Nlevels,*Level,*Levelcnt,*Idbylevel,*Column;
18 
19   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
20   CHKERRQ(PetscProcessTree(n,mask,parentId,&Nlevels,&Level,&Levelcnt,&Idbylevel,&Column));
21   for (i=0; i<n; i++) {
22     if (!mask[i]) {
23       CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," %" PetscInt_FMT " ",Level[i]));
24     }
25   }
26   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nNumber of levels %" PetscInt_FMT "\n",Nlevels));
27   for (i=0; i<Nlevels; i++) {
28     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nLevel %" PetscInt_FMT " ",i));
29     for (j=0; j<Levelcnt[i]; j++) {
30       CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%" PetscInt_FMT " ",Idbylevel[cnt++]));
31     }
32   }
33   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nColumn of each node"));
34   for (i=0; i<n; i++) {
35     if (!mask[i]) {
36       CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," %" PetscInt_FMT " ",Column[i]));
37     }
38   }
39   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n"));
40   CHKERRQ(PetscFree(Level));
41   CHKERRQ(PetscFree(Levelcnt));
42   CHKERRQ(PetscFree(Idbylevel));
43   CHKERRQ(PetscFree(Column));
44   ierr = PetscFinalize();
45   return ierr;
46 }
47 
48 /*TEST
49 
50    test:
51 
52 TEST*/
53