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 5dd7309edSJed Brown static PetscErrorCode ViewOffsets(DM dm, Vec X) 6dd7309edSJed Brown { 7dd7309edSJed Brown PetscInt num_elem, elem_size, num_comp, num_dof; 8dd7309edSJed Brown PetscInt *elem_restr_offsets; 9dd7309edSJed Brown const PetscScalar *x = NULL; 10dd7309edSJed Brown const char *name; 11dd7309edSJed Brown 12dd7309edSJed Brown PetscFunctionBegin; 139566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)dm, &name)); 149566063dSJacob Faibussowitsch PetscCall(DMPlexGetLocalOffsets(dm, NULL, 0, 0, 0, &num_elem, &elem_size, &num_comp, &num_dof, &elem_restr_offsets)); 15*d0609cedSBarry Smith PetscCall(PetscPrintf(PETSC_COMM_SELF,"DM %s offsets: num_elem %" PetscInt_FMT ", size %" PetscInt_FMT 16*d0609cedSBarry Smith ", comp %" PetscInt_FMT ", dof %" PetscInt_FMT "\n",name, num_elem, elem_size, num_comp, num_dof)); 179566063dSJacob Faibussowitsch if (X) PetscCall(VecGetArrayRead(X, &x)); 18dd7309edSJed Brown for (PetscInt c=0; c<num_elem; c++) { 199566063dSJacob Faibussowitsch PetscCall(PetscIntView(elem_size, &elem_restr_offsets[c*elem_size], PETSC_VIEWER_STDOUT_SELF)); 20dd7309edSJed Brown if (x) { 21dd7309edSJed Brown for (PetscInt i=0; i<elem_size; i++) { 229566063dSJacob Faibussowitsch PetscCall(PetscScalarView(num_comp, &x[elem_restr_offsets[c*elem_size+i]], PETSC_VIEWER_STDERR_SELF)); 23dd7309edSJed Brown } 24dd7309edSJed Brown } 25dd7309edSJed Brown } 269566063dSJacob Faibussowitsch if (X) PetscCall(VecRestoreArrayRead(X, &x)); 279566063dSJacob Faibussowitsch PetscCall(PetscFree(elem_restr_offsets)); 28dd7309edSJed Brown PetscFunctionReturn(0); 29dd7309edSJed Brown } 30dd7309edSJed Brown 31c4762a1bSJed Brown int main(int argc, char **argv) 32c4762a1bSJed Brown { 33c4762a1bSJed Brown DM dm; 34c4762a1bSJed Brown PetscSection section; 35c4762a1bSJed Brown PetscFE fe; 3630602db0SMatthew G. Knepley PetscInt dim,c,cStart,cEnd; 37dd7309edSJed Brown PetscBool view_coord = PETSC_FALSE, tensor = PETSC_TRUE; 38c4762a1bSJed Brown 399566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 40*d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Tensor closure restrictions", "DMPLEX"); 419566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-closure_tensor", "Apply DMPlexSetClosurePermutationTensor", "ex8.c", tensor, &tensor, NULL)); 429566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-view_coord", "View coordinates of element closures", "ex8.c", view_coord, &view_coord, NULL)); 43*d0609cedSBarry Smith PetscOptionsEnd(); 4446181b2aSJed Brown 459566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &dm)); 469566063dSJacob Faibussowitsch PetscCall(DMSetType(dm, DMPLEX)); 479566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 489566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 499566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 50c4762a1bSJed Brown 519566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF,dim,1,PETSC_FALSE,NULL,PETSC_DETERMINE,&fe)); 529566063dSJacob Faibussowitsch PetscCall(DMAddField(dm,NULL,(PetscObject)fe)); 539566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 549566063dSJacob Faibussowitsch if (tensor) PetscCall(DMPlexSetClosurePermutationTensor(dm,PETSC_DETERMINE,NULL)); 559566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm,§ion)); 569566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm,0,&cStart,&cEnd)); 57c4762a1bSJed Brown for (c=cStart; c<cEnd; c++) { 58c4762a1bSJed Brown PetscInt numindices,*indices; 599566063dSJacob Faibussowitsch PetscCall(DMPlexGetClosureIndices(dm,section,section,c,PETSC_TRUE,&numindices,&indices,NULL,NULL)); 609566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Element #%" PetscInt_FMT "\n",c-cStart)); 619566063dSJacob Faibussowitsch PetscCall(PetscIntView(numindices,indices,PETSC_VIEWER_STDOUT_SELF)); 629566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreClosureIndices(dm,section,section,c,PETSC_TRUE,&numindices,&indices,NULL,NULL)); 63c4762a1bSJed Brown } 6446181b2aSJed Brown if (view_coord) { 6546181b2aSJed Brown DM cdm; 6646181b2aSJed Brown Vec X; 6746181b2aSJed Brown PetscInt cdim; 689566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm,&cdm)); 699566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)cdm, "coords")); 709566063dSJacob Faibussowitsch if (tensor) PetscCall(DMPlexSetClosurePermutationTensor(cdm,PETSC_DETERMINE,NULL)); 719566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &X)); 729566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &cdim)); 7346181b2aSJed Brown for (c=cStart; c<cEnd; c++) { 74dd7309edSJed Brown PetscScalar *x = NULL; 7546181b2aSJed Brown PetscInt ndof; 769566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(cdm, NULL, X, c, &ndof, &x)); 7746181b2aSJed Brown PetscCheck(ndof % cdim == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "ndof not divisible by cdim"); 789566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Element #%" PetscInt_FMT " coordinates\n",c-cStart)); 7946181b2aSJed Brown for (PetscInt i=0; i<ndof; i+= cdim) { 809566063dSJacob Faibussowitsch PetscCall(PetscScalarView(cdim, &x[i], PETSC_VIEWER_STDOUT_SELF)); 8146181b2aSJed Brown } 829566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(cdm, NULL, X, c, &ndof, &x)); 8346181b2aSJed Brown } 849566063dSJacob Faibussowitsch PetscCall(ViewOffsets(dm, NULL)); 859566063dSJacob Faibussowitsch PetscCall(ViewOffsets(cdm, X)); 8646181b2aSJed Brown } 879566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 889566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 899566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 90b122ec5aSJacob Faibussowitsch return 0; 91c4762a1bSJed Brown } 92c4762a1bSJed Brown 93c4762a1bSJed Brown /*TEST 94c4762a1bSJed Brown 95c4762a1bSJed Brown test: 96c4762a1bSJed Brown suffix: 1d_q2 9730602db0SMatthew G. Knepley args: -dm_plex_dim 1 -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_box_faces 2 98c4762a1bSJed Brown test: 99c4762a1bSJed Brown suffix: 2d_q1 10030602db0SMatthew G. Knepley args: -dm_plex_dim 2 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 2,2 101c4762a1bSJed Brown test: 102c4762a1bSJed Brown suffix: 2d_q2 10330602db0SMatthew G. Knepley args: -dm_plex_dim 2 -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_box_faces 2,2 104c4762a1bSJed Brown test: 105c4762a1bSJed Brown suffix: 2d_q3 10630602db0SMatthew G. Knepley args: -dm_plex_dim 2 -petscspace_degree 3 -dm_plex_simplex 0 -dm_plex_box_faces 1,1 107c4762a1bSJed Brown test: 108c4762a1bSJed Brown suffix: 3d_q1 10930602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1 11046181b2aSJed Brown test: 11146181b2aSJed Brown suffix: 1d_q1_periodic 112a05c9aa3SJed Brown requires: !complex 113dd7309edSJed 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 11446181b2aSJed Brown test: 11546181b2aSJed Brown suffix: 2d_q1_periodic 116a05c9aa3SJed Brown requires: !complex 117dd7309edSJed 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 11846181b2aSJed Brown test: 119a05c9aa3SJed Brown suffix: 3d_q1_periodic 120a05c9aa3SJed Brown requires: !complex 121a05c9aa3SJed Brown args: -dm_plex_dim 3 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 3,2,1 -dm_plex_box_bd periodic,none,none -dm_view -view_coord 122a05c9aa3SJed Brown test: 123a05c9aa3SJed Brown suffix: 3d_q2_periodic # not actually periodic because only 2 cells 12446181b2aSJed 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 125c4762a1bSJed Brown 126c4762a1bSJed Brown TEST*/ 127