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