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