1 static char help[] = "Tests for periodic mesh output\n\n"; 2 3 #include <petscdmplex.h> 4 5 PetscErrorCode CheckMesh(DM dm) 6 { 7 PetscReal detJ, J[9]; 8 PetscReal vol; 9 PetscInt dim, depth, cStart, cEnd, c; 10 11 PetscFunctionBegin; 12 PetscCall(DMGetDimension(dm, &dim)); 13 PetscCall(DMPlexGetDepth(dm, &depth)); 14 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 15 for (c = cStart; c < cEnd; ++c) { 16 PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, NULL, J, NULL, &detJ)); 17 PetscCheck(detJ > 0.0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh cell %" PetscInt_FMT " is inverted, |J| = %g", c, (double)detJ); 18 if (depth > 1) { 19 PetscCall(DMPlexComputeCellGeometryFVM(dm, c, &vol, NULL, NULL)); 20 PetscCheck(vol > 0.0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh cell %" PetscInt_FMT " is inverted, vol = %g", c, (double)vol); 21 } 22 } 23 PetscFunctionReturn(0); 24 } 25 26 PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm) 27 { 28 PetscFunctionBegin; 29 PetscCall(DMCreate(comm, dm)); 30 PetscCall(DMSetType(*dm, DMPLEX)); 31 PetscCall(DMSetFromOptions(*dm)); 32 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 33 PetscFunctionReturn(0); 34 } 35 36 int main(int argc, char **argv) 37 { 38 DM dm; 39 40 PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 41 PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm)); 42 PetscCall(CheckMesh(dm)); 43 PetscCall(DMDestroy(&dm)); 44 PetscCall(PetscFinalize()); 45 return 0; 46 } 47 48 /*TEST 49 50 test: 51 suffix: 0 52 args: -dm_plex_simplex 0 -dm_plex_box_faces 3,1,0 -dm_plex_box_bd periodic,none -dm_view ::ascii_info_detail 53 test: 54 suffix: 1 55 nsize: 2 56 args: -dm_plex_simplex 0 -dm_plex_box_faces 3,1,0 -dm_plex_box_bd periodic,none -petscpartitioner_type simple -dm_view ::ascii_info_detail 57 test: 58 suffix: 2 59 nsize: 2 60 args: -dm_plex_simplex 0 -dm_plex_box_faces 6,2,0 -dm_plex_box_bd periodic,none -petscpartitioner_type simple -dm_view ::ascii_info_detail 61 test: 62 suffix: 3 63 nsize: 4 64 args: -dm_plex_simplex 0 -dm_plex_box_faces 6,2,0 -dm_plex_box_bd periodic,none -petscpartitioner_type simple -dm_view ::ascii_info_detail 65 test: 66 suffix: 4 67 nsize: 2 68 args: -dm_plex_simplex 0 -dm_plex_box_faces 3,1,0 -dm_plex_box_bd periodic,none -dm_plex_periodic_cut -petscpartitioner_type simple -dm_view ::ascii_info_detail 69 test: 70 suffix: 5 71 nsize: 2 72 args: -dm_plex_simplex 0 -dm_plex_box_faces 6,2,0 -dm_plex_box_bd periodic,none -dm_plex_periodic_cut -petscpartitioner_type simple -dm_view ::ascii_info_detail 73 test: 74 suffix: 6 75 nsize: 4 76 args: -dm_plex_simplex 0 -dm_plex_box_faces 6,2,0 -dm_plex_box_bd periodic,none -dm_plex_periodic_cut -petscpartitioner_type simple -dm_view ::ascii_info_detail 77 78 TEST*/ 79