1 static const char help[] = "Tests for Plex transforms, including regular refinement"; 2 3 #include <petscdmplex.h> 4 #include <petscsf.h> 5 6 #include <petsc/private/dmpleximpl.h> 7 8 static PetscErrorCode LabelPoints(DM dm) 9 { 10 DMLabel label; 11 PetscInt pStart, pEnd, p; 12 PetscBool flg = PETSC_FALSE; 13 PetscErrorCode ierr; 14 15 PetscFunctionBegin; 16 ierr = PetscOptionsGetBool(NULL, NULL, "-label_mesh", &flg, NULL);CHKERRQ(ierr); 17 if (!flg) PetscFunctionReturn(0); 18 ierr = DMCreateLabel(dm, "test");CHKERRQ(ierr); 19 ierr = DMGetLabel(dm, "test", &label);CHKERRQ(ierr); 20 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 21 for (p = pStart; p < pEnd; ++p) { 22 ierr = DMLabelSetValue(label, p, p);CHKERRQ(ierr); 23 } 24 PetscFunctionReturn(0); 25 } 26 27 static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm) 28 { 29 PetscErrorCode ierr; 30 31 PetscFunctionBegin; 32 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 33 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 34 ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 35 ierr = LabelPoints(*dm);CHKERRQ(ierr); 36 ierr = PetscObjectSetOptionsPrefix((PetscObject) *dm, "post_label_");CHKERRQ(ierr); 37 ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 38 ierr = PetscObjectSetOptionsPrefix((PetscObject) *dm, NULL);CHKERRQ(ierr); 39 ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 40 PetscFunctionReturn(0); 41 } 42 43 int main(int argc, char **argv) 44 { 45 DM dm; 46 PetscErrorCode ierr; 47 48 ierr = PetscInitialize(&argc, &argv, NULL, help); if (ierr) return ierr; 49 ierr = CreateMesh(PETSC_COMM_WORLD, &dm);CHKERRQ(ierr); 50 ierr = DMDestroy(&dm);CHKERRQ(ierr); 51 ierr = PetscFinalize(); 52 return ierr; 53 } 54 55 /*TEST 56 test: 57 suffix: ref_seg 58 args: -dm_plex_reference_cell_domain -dm_plex_cell segment -dm_refine 1 -dm_plex_check_all 59 60 test: 61 suffix: ref_tri 62 args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_refine 2 -dm_plex_check_all 63 64 test: 65 suffix: box_tri 66 requires: triangle 67 nsize: {{1 3 5}} 68 args: -dm_distribute -dm_plex_box_faces 3,3 -dm_refine 2 -dm_plex_check_all 69 70 test: 71 suffix: ref_quad 72 args: -dm_plex_reference_cell_domain -dm_plex_cell quadrilateral -dm_refine 2 -dm_plex_check_all 73 74 test: 75 suffix: box_quad 76 nsize: {{1 3 5}} 77 args: -dm_distribute -dm_plex_box_faces 3,3 -dm_plex_simplex 0 -dm_refine 2 -dm_plex_check_all 78 79 test: 80 suffix: ref_tet 81 args: -dm_plex_reference_cell_domain -dm_plex_cell tetrahedron -dm_refine 2 -dm_plex_check_all 82 83 test: 84 suffix: box_tet 85 requires: ctetgen 86 nsize: {{1 3 5}} 87 args: -dm_distribute -dm_plex_dim 3 -dm_plex_box_faces 3,3,3 -dm_refine 2 -dm_plex_check_all 88 89 test: 90 suffix: ref_hex 91 args: -dm_plex_reference_cell_domain -dm_plex_cell hexahedron -dm_refine 2 -dm_plex_check_all 92 93 test: 94 suffix: box_hex 95 nsize: {{1 3 5}} 96 args: -dm_distribute -dm_plex_dim 3 -dm_plex_box_faces 3,3,3 -dm_plex_simplex 0 -dm_refine 2 -dm_plex_check_all 97 98 test: 99 suffix: ref_trip 100 args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism -dm_refine 2 -dm_plex_check_all 101 102 test: 103 suffix: ref_tquad 104 args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quad -dm_refine 2 -dm_plex_check_all 105 106 test: 107 suffix: ref_ttrip 108 args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_triangular_prism -dm_refine 2 -dm_plex_check_all 109 110 test: 111 suffix: ref_tquadp 112 args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quadrilateral_prism -dm_refine 2 -dm_plex_check_all 113 114 test: 115 suffix: ref_pyramid 116 args: -dm_plex_reference_cell_domain -dm_plex_cell pyramid -dm_refine 2 -dm_plex_check_all 117 118 testset: 119 args: -dm_coord_space 0 -dm_plex_transform_type refine_tobox -dm_plex_check_all 120 121 test: 122 suffix: ref_tri_tobox 123 args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_refine 2 124 125 test: 126 suffix: box_tri_tobox 127 requires: triangle 128 nsize: {{1 3 5}} 129 args: -dm_distribute -dm_plex_box_faces 3,3 -dm_refine 2 130 131 test: 132 suffix: ref_tet_tobox 133 args: -dm_plex_reference_cell_domain -dm_plex_cell tetrahedron -dm_refine 2 134 135 test: 136 suffix: box_tet_tobox 137 requires: ctetgen 138 nsize: {{1 3 5}} 139 args: -dm_distribute -dm_plex_dim 3 -dm_plex_box_faces 3,3,3 -dm_refine 2 140 141 test: 142 suffix: ref_trip_tobox 143 args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism -dm_refine 2 144 145 test: 146 suffix: ref_ttrip_tobox 147 args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_triangular_prism -dm_refine 2 148 149 test: 150 suffix: ref_tquadp_tobox 151 args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quadrilateral_prism -dm_refine 2 152 153 testset: 154 args: -dm_coord_space 0 -label_mesh -post_label_dm_extrude 2 -post_label_dm_plex_check_all -dm_view ::ascii_info_detail 155 156 test: 157 suffix: extrude_quad 158 args: -dm_plex_simplex 0 159 160 TEST*/ 161