1c4762a1bSJed Brown static char help[] = "Tests for creation of submeshes\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmplex.h> 4c4762a1bSJed Brown 5*30602db0SMatthew G. Knepley PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm) 6c4762a1bSJed Brown { 7c4762a1bSJed Brown PetscErrorCode ierr; 8c4762a1bSJed Brown 9c4762a1bSJed Brown PetscFunctionBegin; 10*30602db0SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 11*30602db0SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 12c4762a1bSJed Brown ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 13c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 14c4762a1bSJed Brown PetscFunctionReturn(0); 15c4762a1bSJed Brown } 16c4762a1bSJed Brown 17c4762a1bSJed Brown PetscErrorCode CreateSubmesh(DM dm, PetscBool start, DM *subdm) 18c4762a1bSJed Brown { 19c4762a1bSJed Brown DMLabel label, map; 20c4762a1bSJed Brown PetscInt cStart, cEnd, cStartSub, cEndSub, c; 21c4762a1bSJed Brown PetscErrorCode ierr; 22c4762a1bSJed Brown 23c4762a1bSJed Brown PetscFunctionBegin; 24c4762a1bSJed Brown ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 25c4762a1bSJed Brown ierr = DMCreateLabel(dm, "cells");CHKERRQ(ierr); 26c4762a1bSJed Brown ierr = DMGetLabel(dm, "cells", &label);CHKERRQ(ierr); 27c4762a1bSJed Brown ierr = DMLabelClearStratum(label, 1);CHKERRQ(ierr); 28c4762a1bSJed Brown if (start) {cStartSub = cStart; cEndSub = cEnd/2;} 29c4762a1bSJed Brown else {cStartSub = cEnd/2; cEndSub = cEnd;} 30c4762a1bSJed Brown for (c = cStartSub; c < cEndSub; ++c) {ierr = DMLabelSetValue(label, c, 1);CHKERRQ(ierr);} 31c4762a1bSJed Brown ierr = DMPlexFilter(dm, label, 1, subdm);CHKERRQ(ierr); 32c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) *subdm, "Submesh");CHKERRQ(ierr); 33c4762a1bSJed Brown ierr = DMViewFromOptions(*subdm, NULL, "-dm_view");CHKERRQ(ierr); 34c4762a1bSJed Brown ierr = DMPlexGetSubpointMap(*subdm, &map);CHKERRQ(ierr); 35c4762a1bSJed Brown ierr = DMLabelView(map, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 36c4762a1bSJed Brown PetscFunctionReturn(0); 37c4762a1bSJed Brown } 38c4762a1bSJed Brown 39c4762a1bSJed Brown int main(int argc, char **argv) 40c4762a1bSJed Brown { 41c4762a1bSJed Brown DM dm, subdm; 42c4762a1bSJed Brown PetscErrorCode ierr; 43c4762a1bSJed Brown 44c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 45*30602db0SMatthew G. Knepley ierr = CreateMesh(PETSC_COMM_WORLD, &dm);CHKERRQ(ierr); 46c4762a1bSJed Brown ierr = CreateSubmesh(dm, PETSC_TRUE, &subdm);CHKERRQ(ierr); 47c4762a1bSJed Brown ierr = DMSetFromOptions(subdm);CHKERRQ(ierr); 48c4762a1bSJed Brown ierr = DMDestroy(&subdm);CHKERRQ(ierr); 49c4762a1bSJed Brown ierr = CreateSubmesh(dm, PETSC_FALSE, &subdm);CHKERRQ(ierr); 50c4762a1bSJed Brown ierr = DMSetFromOptions(subdm);CHKERRQ(ierr); 51c4762a1bSJed Brown ierr = DMDestroy(&subdm);CHKERRQ(ierr); 52c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 53c4762a1bSJed Brown ierr = PetscFinalize(); 54c4762a1bSJed Brown return ierr; 55c4762a1bSJed Brown } 56c4762a1bSJed Brown 57c4762a1bSJed Brown /*TEST 58c4762a1bSJed Brown 59c4762a1bSJed Brown test: 60c4762a1bSJed Brown suffix: 0 61c4762a1bSJed Brown requires: triangle 62*30602db0SMatthew G. Knepley args: -dm_coord_space 0 -dm_view ascii::ascii_info_detail -dm_plex_check_all 63c4762a1bSJed Brown 64c4762a1bSJed Brown TEST*/ 65