xref: /petsc/src/dm/impls/plex/tutorials/ex8.c (revision dd7309eda20e1908cb4ae8f4a64e23c6acefadd5)
1c4762a1bSJed Brown static char help[] = "Element closure restrictions in tensor/lexicographic/spectral-element ordering using DMPlex\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdmplex.h>
4c4762a1bSJed Brown 
5*dd7309edSJed Brown static PetscErrorCode ViewOffsets(DM dm, Vec X)
6*dd7309edSJed Brown {
7*dd7309edSJed Brown   PetscInt num_elem, elem_size, num_comp, num_dof;
8*dd7309edSJed Brown   PetscInt *elem_restr_offsets;
9*dd7309edSJed Brown   const PetscScalar *x = NULL;
10*dd7309edSJed Brown   const char *name;
11*dd7309edSJed Brown   PetscErrorCode ierr;
12*dd7309edSJed Brown 
13*dd7309edSJed Brown   PetscFunctionBegin;
14*dd7309edSJed Brown   ierr = PetscObjectGetName((PetscObject)dm, &name);CHKERRQ(ierr);
15*dd7309edSJed Brown   ierr = DMPlexGetLocalOffsets(dm, NULL, 0, 0, 0, &num_elem, &elem_size, &num_comp, &num_dof, &elem_restr_offsets);CHKERRQ(ierr);
16*dd7309edSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF,"DM %s offsets: num_elem %" PetscInt_FMT ", size %" PetscInt_FMT
17*dd7309edSJed Brown                      ", comp %" PetscInt_FMT ", dof %" PetscInt_FMT "\n",
18*dd7309edSJed Brown                      name, num_elem, elem_size, num_comp, num_dof);CHKERRQ(ierr);
19*dd7309edSJed Brown   if (X) {ierr = VecGetArrayRead(X, &x);CHKERRQ(ierr);}
20*dd7309edSJed Brown   for (PetscInt c=0; c<num_elem; c++) {
21*dd7309edSJed Brown     ierr = PetscIntView(elem_size, &elem_restr_offsets[c*elem_size], PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
22*dd7309edSJed Brown     if (x) {
23*dd7309edSJed Brown       for (PetscInt i=0; i<elem_size; i++) {
24*dd7309edSJed Brown         ierr = PetscScalarView(num_comp, &x[elem_restr_offsets[c*elem_size+i]], PETSC_VIEWER_STDERR_SELF);CHKERRQ(ierr);
25*dd7309edSJed Brown       }
26*dd7309edSJed Brown     }
27*dd7309edSJed Brown   }
28*dd7309edSJed Brown   if (X) {ierr = VecRestoreArrayRead(X, &x);CHKERRQ(ierr);}
29*dd7309edSJed Brown   ierr = PetscFree(elem_restr_offsets);CHKERRQ(ierr);
30*dd7309edSJed Brown   PetscFunctionReturn(0);
31*dd7309edSJed Brown }
32*dd7309edSJed Brown 
33c4762a1bSJed Brown int main(int argc, char **argv)
34c4762a1bSJed Brown {
35c4762a1bSJed Brown   DM             dm;
36c4762a1bSJed Brown   PetscSection   section;
37c4762a1bSJed Brown   PetscFE        fe;
3830602db0SMatthew G. Knepley   PetscInt       dim,c,cStart,cEnd;
39*dd7309edSJed Brown   PetscBool      view_coord = PETSC_FALSE, tensor = PETSC_TRUE;
40c4762a1bSJed Brown   PetscErrorCode ierr;
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
4346181b2aSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Tensor closure restrictions", "DMPLEX");CHKERRQ(ierr);
44*dd7309edSJed Brown   ierr = PetscOptionsBool("-closure_tensor", "Apply DMPlexSetClosurePermutationTensor", "ex8.c", tensor, &tensor, NULL);CHKERRQ(ierr);
4546181b2aSJed Brown   ierr = PetscOptionsBool("-view_coord", "View coordinates of element closures", "ex8.c", view_coord, &view_coord, NULL);CHKERRQ(ierr);
4646181b2aSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
4746181b2aSJed Brown 
4830602db0SMatthew G. Knepley   ierr = DMCreate(PETSC_COMM_WORLD, &dm);CHKERRQ(ierr);
4930602db0SMatthew G. Knepley   ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr);
5030602db0SMatthew G. Knepley   ierr = DMSetFromOptions(dm);CHKERRQ(ierr);
5130602db0SMatthew G. Knepley   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
5230602db0SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   ierr = PetscFECreateDefault(PETSC_COMM_SELF,dim,1,PETSC_FALSE,NULL,PETSC_DETERMINE,&fe);CHKERRQ(ierr);
55c4762a1bSJed Brown   ierr = DMAddField(dm,NULL,(PetscObject)fe);CHKERRQ(ierr);
56c4762a1bSJed Brown   ierr = DMCreateDS(dm);CHKERRQ(ierr);
57*dd7309edSJed Brown   if (tensor) {ierr = DMPlexSetClosurePermutationTensor(dm,PETSC_DETERMINE,NULL);CHKERRQ(ierr);}
58c4762a1bSJed Brown   ierr = DMGetLocalSection(dm,&section);CHKERRQ(ierr);
59c4762a1bSJed Brown   ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr);
60c4762a1bSJed Brown   for (c=cStart; c<cEnd; c++) {
61c4762a1bSJed Brown     PetscInt numindices,*indices;
6271f0bbf9SMatthew G. Knepley     ierr = DMPlexGetClosureIndices(dm,section,section,c,PETSC_TRUE,&numindices,&indices,NULL,NULL);CHKERRQ(ierr);
6346181b2aSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF,"Element #%" PetscInt_FMT "\n",c-cStart);CHKERRQ(ierr);
64c4762a1bSJed Brown     ierr = PetscIntView(numindices,indices,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
6571f0bbf9SMatthew G. Knepley     ierr = DMPlexRestoreClosureIndices(dm,section,section,c,PETSC_TRUE,&numindices,&indices,NULL,NULL);CHKERRQ(ierr);
66c4762a1bSJed Brown   }
6746181b2aSJed Brown   if (view_coord) {
6846181b2aSJed Brown     DM cdm;
6946181b2aSJed Brown     Vec X;
7046181b2aSJed Brown     PetscInt cdim;
7146181b2aSJed Brown     ierr = DMGetCoordinateDM(dm,&cdm);CHKERRQ(ierr);
72*dd7309edSJed Brown     ierr = PetscObjectSetName((PetscObject)cdm, "coords");CHKERRQ(ierr);
73*dd7309edSJed Brown     if (tensor) {ierr = DMPlexSetClosurePermutationTensor(cdm,PETSC_DETERMINE,NULL);CHKERRQ(ierr);}
7446181b2aSJed Brown     ierr = DMGetCoordinatesLocal(dm, &X);CHKERRQ(ierr);
7546181b2aSJed Brown     ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr);
7646181b2aSJed Brown     for (c=cStart; c<cEnd; c++) {
77*dd7309edSJed Brown       PetscScalar *x = NULL;
7846181b2aSJed Brown       PetscInt ndof;
7946181b2aSJed Brown       ierr = DMPlexVecGetClosure(cdm, NULL, X, c, &ndof, &x);CHKERRQ(ierr);
8046181b2aSJed Brown       PetscCheck(ndof % cdim == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "ndof not divisible by cdim");
8146181b2aSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF,"Element #%" PetscInt_FMT " coordinates\n",c-cStart);CHKERRQ(ierr);
8246181b2aSJed Brown       for (PetscInt i=0; i<ndof; i+= cdim) {
8346181b2aSJed Brown         ierr = PetscScalarView(cdim, &x[i], PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
8446181b2aSJed Brown       }
8546181b2aSJed Brown       ierr = DMPlexVecRestoreClosure(cdm, NULL, X, c, &ndof, &x);CHKERRQ(ierr);
8646181b2aSJed Brown     }
87*dd7309edSJed Brown     ierr = ViewOffsets(dm, NULL);CHKERRQ(ierr);
88*dd7309edSJed Brown     ierr = ViewOffsets(cdm, X);CHKERRQ(ierr);
8946181b2aSJed Brown   }
90c4762a1bSJed Brown   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
91c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
92c4762a1bSJed Brown   ierr = PetscFinalize();
93c4762a1bSJed Brown   return ierr;
94c4762a1bSJed Brown }
95c4762a1bSJed Brown 
96c4762a1bSJed Brown /*TEST
97c4762a1bSJed Brown 
98c4762a1bSJed Brown   test:
99c4762a1bSJed Brown     suffix: 1d_q2
10030602db0SMatthew G. Knepley     args: -dm_plex_dim 1 -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_box_faces 2
101c4762a1bSJed Brown   test:
102c4762a1bSJed Brown     suffix: 2d_q1
10330602db0SMatthew G. Knepley     args: -dm_plex_dim 2 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 2,2
104c4762a1bSJed Brown   test:
105c4762a1bSJed Brown     suffix: 2d_q2
10630602db0SMatthew G. Knepley     args: -dm_plex_dim 2 -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_box_faces 2,2
107c4762a1bSJed Brown   test:
108c4762a1bSJed Brown     suffix: 2d_q3
10930602db0SMatthew G. Knepley     args: -dm_plex_dim 2 -petscspace_degree 3 -dm_plex_simplex 0 -dm_plex_box_faces 1,1
110c4762a1bSJed Brown   test:
111c4762a1bSJed Brown     suffix: 3d_q1
11230602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1
11346181b2aSJed Brown   test:
11446181b2aSJed Brown     suffix: 1d_q1_periodic
115*dd7309edSJed Brown     args: -dm_plex_dim 1 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 3 -dm_plex_box_bd periodic -dm_view -view_coord
11646181b2aSJed Brown   test:
11746181b2aSJed Brown     suffix: 2d_q1_periodic
118*dd7309edSJed Brown     args: -dm_plex_dim 2 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 3,2 -dm_plex_box_bd periodic,none -dm_view -view_coord
11946181b2aSJed Brown   test:
12046181b2aSJed Brown     suffix: 3d_q2_periodic
12146181b2aSJed Brown     args: -dm_plex_dim 3 -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_box_faces 2,2,2 -dm_plex_box_bd periodic,none,periodic -dm_view
122c4762a1bSJed Brown 
123c4762a1bSJed Brown TEST*/
124