xref: /petsc/src/sys/tests/ex21.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,NULL,help));
195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscProcessTree(n,mask,parentId,&Nlevels,&Level,&Levelcnt,&Idbylevel,&Column));
20c4762a1bSJed Brown   for (i=0; i<n; i++) {
21c4762a1bSJed Brown     if (!mask[i]) {
225f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," %" PetscInt_FMT " ",Level[i]));
23c4762a1bSJed Brown     }
24c4762a1bSJed Brown   }
255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nNumber of levels %" PetscInt_FMT "\n",Nlevels));
26c4762a1bSJed Brown   for (i=0; i<Nlevels; i++) {
275f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nLevel %" PetscInt_FMT " ",i));
28c4762a1bSJed Brown     for (j=0; j<Levelcnt[i]; j++) {
295f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%" PetscInt_FMT " ",Idbylevel[cnt++]));
30c4762a1bSJed Brown     }
31c4762a1bSJed Brown   }
325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nColumn of each node"));
33c4762a1bSJed Brown   for (i=0; i<n; i++) {
34c4762a1bSJed Brown     if (!mask[i]) {
355f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," %" PetscInt_FMT " ",Column[i]));
36c4762a1bSJed Brown     }
37c4762a1bSJed Brown   }
385f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n"));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(Level));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(Levelcnt));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(Idbylevel));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(Column));
43*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
44*b122ec5aSJacob Faibussowitsch   return 0;
45c4762a1bSJed Brown }
46c4762a1bSJed Brown 
47c4762a1bSJed Brown /*TEST
48c4762a1bSJed Brown 
49c4762a1bSJed Brown    test:
50c4762a1bSJed Brown 
51c4762a1bSJed Brown TEST*/
52