xref: /petsc/src/sys/tests/ex21.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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   PetscInt       n          = 7,cnt = 0,i,j;
14c4762a1bSJed Brown   PetscBool      mask[]     = {PETSC_TRUE,PETSC_FALSE,PETSC_FALSE,PETSC_TRUE,PETSC_FALSE,PETSC_FALSE,PETSC_FALSE};
15c4762a1bSJed Brown   PetscInt       parentId[] = {-1,         2,         0,         -1,         2,         1,         0};
16c4762a1bSJed Brown   PetscInt       Nlevels,*Level,*Levelcnt,*Idbylevel,*Column;
17c4762a1bSJed Brown 
18*327415f7SBarry Smith   PetscFunctionBeginUser;
199566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
209566063dSJacob Faibussowitsch   PetscCall(PetscProcessTree(n,mask,parentId,&Nlevels,&Level,&Levelcnt,&Idbylevel,&Column));
21c4762a1bSJed Brown   for (i=0; i<n; i++) {
22c4762a1bSJed Brown     if (!mask[i]) {
239566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD," %" PetscInt_FMT " ",Level[i]));
24c4762a1bSJed Brown     }
25c4762a1bSJed Brown   }
269566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nNumber of levels %" PetscInt_FMT "\n",Nlevels));
27c4762a1bSJed Brown   for (i=0; i<Nlevels; i++) {
289566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nLevel %" PetscInt_FMT " ",i));
29c4762a1bSJed Brown     for (j=0; j<Levelcnt[i]; j++) {
309566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%" PetscInt_FMT " ",Idbylevel[cnt++]));
31c4762a1bSJed Brown     }
32c4762a1bSJed Brown   }
339566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nColumn of each node"));
34c4762a1bSJed Brown   for (i=0; i<n; i++) {
35c4762a1bSJed Brown     if (!mask[i]) {
369566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD," %" PetscInt_FMT " ",Column[i]));
37c4762a1bSJed Brown     }
38c4762a1bSJed Brown   }
399566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n"));
409566063dSJacob Faibussowitsch   PetscCall(PetscFree(Level));
419566063dSJacob Faibussowitsch   PetscCall(PetscFree(Levelcnt));
429566063dSJacob Faibussowitsch   PetscCall(PetscFree(Idbylevel));
439566063dSJacob Faibussowitsch   PetscCall(PetscFree(Column));
449566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
45b122ec5aSJacob Faibussowitsch   return 0;
46c4762a1bSJed Brown }
47c4762a1bSJed Brown 
48c4762a1bSJed Brown /*TEST
49c4762a1bSJed Brown 
50c4762a1bSJed Brown    test:
51c4762a1bSJed Brown 
52c4762a1bSJed Brown TEST*/
53