xref: /petsc/src/dm/impls/plex/tests/ex20.c (revision 6a98f8dc3f2c9149905a87dc2e9d0fedaf64e09a)
1 
2 const char help[] = "Test DMPlex implementation of DMAdaptLabel().\n\n";
3 
4 #include <petscdm.h>
5 #include <petscdmplex.h>
6 
7 int main(int argc, char **argv)
8 {
9   DM             dm, dmAdapt;
10   DMLabel        adaptLabel;
11   PetscInt       dim, nfaces, faces[3], cStart, cEnd;
12   PetscBool      interpolate;
13   PetscErrorCode ierr;
14 
15   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
16   dim         = 2;
17   nfaces      = 3;
18   interpolate = PETSC_TRUE;
19   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ex20",NULL);CHKERRQ(ierr);
20   ierr = PetscOptionsRangeInt("-dim","domain dimension",NULL,dim,&dim,NULL,1,3);CHKERRQ(ierr);
21   ierr = PetscOptionsBoundedInt("-nfaces","number of faces per dimension",NULL,nfaces,&nfaces,NULL,0);CHKERRQ(ierr);
22   ierr = PetscOptionsEnd();CHKERRQ(ierr);
23   faces[0] = faces[1] = faces[2] = nfaces;
24   ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD,dim,PETSC_TRUE,faces,NULL,NULL,NULL,interpolate,&dm);CHKERRQ(ierr);
25   ierr = PetscObjectSetName((PetscObject)dm,"Pre Adaptation Mesh");CHKERRQ(ierr);
26   ierr = DMViewFromOptions(dm,NULL,"-pre_adapt_dm_view");CHKERRQ(ierr);
27   ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr);
28   ierr = DMLabelCreate(PETSC_COMM_SELF,"adapt",&adaptLabel);CHKERRQ(ierr);
29   ierr = DMLabelSetDefaultValue(adaptLabel,DM_ADAPT_COARSEN);CHKERRQ(ierr);
30   if (cEnd > cStart) {ierr = DMLabelSetValue(adaptLabel,cStart,DM_ADAPT_REFINE);CHKERRQ(ierr);}
31   ierr = DMAdaptLabel(dm,adaptLabel,&dmAdapt);CHKERRQ(ierr);
32   ierr = PetscObjectSetName((PetscObject)dmAdapt,"Post Adaptation Mesh");CHKERRQ(ierr);
33   ierr = DMViewFromOptions(dmAdapt,NULL,"-post_adapt_dm_view");CHKERRQ(ierr);
34   ierr = DMDestroy(&dmAdapt);CHKERRQ(ierr);
35   ierr = DMLabelDestroy(&adaptLabel);CHKERRQ(ierr);
36   ierr = DMDestroy(&dm);CHKERRQ(ierr);
37   ierr = PetscFinalize();
38   return ierr;
39 }
40 
41 /*TEST
42 
43   test:
44     suffix: 2d
45     requires: triangle !single
46     args: -dim 2 -pre_adapt_dm_view ascii::ascii_info_detail -post_adapt_dm_view ascii::ascii_info_detail
47   test:
48     suffix: 3d_tetgen
49     requires: tetgen complex
50     args: -dim 3 -pre_adapt_dm_view ascii::ascii_info_detail -post_adapt_dm_view ascii::ascii_info_detail
51   test:
52     suffix: 3d_ctetgen
53     requires: ctetgen !complex !single
54     args: -dim 3 -pre_adapt_dm_view ascii::ascii_info_detail -post_adapt_dm_view ascii::ascii_info_detail
55 
56 TEST*/
57