1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 3 /*@C 4 DMPlexInvertCell - Flips cell orientations since Plex stores some of them internally with outward normals. 5 6 Input Parameters: 7 + cellType - The cell type 8 - cone - The incoming cone 9 10 Output Parameter: 11 . cone - The inverted cone (in-place) 12 13 Level: developer 14 15 .seealso: DMPlexGenerate() 16 @*/ 17 PetscErrorCode DMPlexInvertCell(DMPolytopeType cellType, PetscInt cone[]) 18 { 19 #define SWAPCONE(cone,i,j) \ 20 do { \ 21 PetscInt _cone_tmp; \ 22 _cone_tmp = (cone)[i]; \ 23 (cone)[i] = (cone)[j]; \ 24 (cone)[j] = _cone_tmp; \ 25 } while (0) 26 27 PetscFunctionBegin; 28 switch (cellType) { 29 case DM_POLYTOPE_POINT: break; 30 case DM_POLYTOPE_SEGMENT: break; 31 case DM_POLYTOPE_POINT_PRISM_TENSOR: break; 32 case DM_POLYTOPE_TRIANGLE: break; 33 case DM_POLYTOPE_QUADRILATERAL: break; 34 case DM_POLYTOPE_SEG_PRISM_TENSOR: SWAPCONE(cone,2,3); break; 35 case DM_POLYTOPE_TETRAHEDRON: SWAPCONE(cone,0,1); break; 36 case DM_POLYTOPE_HEXAHEDRON: SWAPCONE(cone,1,3); break; 37 case DM_POLYTOPE_TRI_PRISM: SWAPCONE(cone,1,2); break; 38 case DM_POLYTOPE_TRI_PRISM_TENSOR: break; 39 case DM_POLYTOPE_QUAD_PRISM_TENSOR: break; 40 default: break; 41 } 42 PetscFunctionReturn(0); 43 #undef SWAPCONE 44 } 45 46 /*@C 47 DMPlexReorderCell - Flips cell orientations since Plex stores some of them internally with outward normals. 48 49 Input Parameters: 50 + dm - The DMPlex object 51 . cell - The cell 52 - cone - The incoming cone 53 54 Output Parameter: 55 . cone - The reordered cone (in-place) 56 57 Level: developer 58 59 .seealso: DMPlexGenerate() 60 @*/ 61 PetscErrorCode DMPlexReorderCell(DM dm, PetscInt cell, PetscInt cone[]) 62 { 63 DMPolytopeType cellType; 64 PetscErrorCode ierr; 65 66 PetscFunctionBegin; 67 ierr = DMPlexGetCellType(dm, cell, &cellType);CHKERRQ(ierr); 68 ierr = DMPlexInvertCell(cellType, cone);CHKERRQ(ierr); 69 PetscFunctionReturn(0); 70 } 71 72 /*@C 73 DMPlexTriangleSetOptions - Set the options used for the Triangle mesh generator 74 75 Not Collective 76 77 Inputs Parameters: 78 + dm - The DMPlex object 79 - opts - The command line options 80 81 Level: developer 82 83 .seealso: DMPlexTetgenSetOptions(), DMPlexGenerate() 84 @*/ 85 PetscErrorCode DMPlexTriangleSetOptions(DM dm, const char *opts) 86 { 87 DM_Plex *mesh = (DM_Plex*) dm->data; 88 PetscErrorCode ierr; 89 90 PetscFunctionBegin; 91 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 92 PetscValidPointer(opts, 2); 93 ierr = PetscFree(mesh->triangleOpts);CHKERRQ(ierr); 94 ierr = PetscStrallocpy(opts, &mesh->triangleOpts);CHKERRQ(ierr); 95 PetscFunctionReturn(0); 96 } 97 98 /*@C 99 DMPlexTetgenSetOptions - Set the options used for the Tetgen mesh generator 100 101 Not Collective 102 103 Inputs Parameters: 104 + dm - The DMPlex object 105 - opts - The command line options 106 107 Level: developer 108 109 .seealso: DMPlexTriangleSetOptions(), DMPlexGenerate() 110 @*/ 111 PetscErrorCode DMPlexTetgenSetOptions(DM dm, const char *opts) 112 { 113 DM_Plex *mesh = (DM_Plex*) dm->data; 114 PetscErrorCode ierr; 115 116 PetscFunctionBegin; 117 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 118 PetscValidPointer(opts, 2); 119 ierr = PetscFree(mesh->tetgenOpts);CHKERRQ(ierr); 120 ierr = PetscStrallocpy(opts, &mesh->tetgenOpts);CHKERRQ(ierr); 121 PetscFunctionReturn(0); 122 } 123 124 /*@C 125 DMPlexGenerate - Generates a mesh. 126 127 Not Collective 128 129 Input Parameters: 130 + boundary - The DMPlex boundary object 131 . name - The mesh generation package name 132 - interpolate - Flag to create intermediate mesh elements 133 134 Output Parameter: 135 . mesh - The DMPlex object 136 137 Options Database: 138 + -dm_plex_generate <name> - package to generate mesh, for example, triangle, ctetgen or tetgen 139 - -dm_generator <name> - package to generate mesh, for example, triangle, ctetgen or tetgen 140 141 Level: intermediate 142 143 .seealso: DMPlexCreate(), DMRefine() 144 @*/ 145 PetscErrorCode DMPlexGenerate(DM boundary, const char name[], PetscBool interpolate, DM *mesh) 146 { 147 DMGeneratorFunctionList fl; 148 char genname[PETSC_MAX_PATH_LEN]; 149 const char *suggestions; 150 PetscInt dim; 151 PetscBool flg; 152 PetscErrorCode ierr; 153 154 PetscFunctionBegin; 155 PetscValidHeaderSpecific(boundary, DM_CLASSID, 1); 156 PetscValidLogicalCollectiveBool(boundary, interpolate, 3); 157 ierr = DMGetDimension(boundary, &dim);CHKERRQ(ierr); 158 ierr = PetscOptionsGetString(((PetscObject) boundary)->options,((PetscObject) boundary)->prefix, "-dm_generator", genname, sizeof(genname), &flg);CHKERRQ(ierr); 159 if (flg) name = genname; 160 else { 161 ierr = PetscOptionsGetString(((PetscObject) boundary)->options,((PetscObject) boundary)->prefix, "-dm_plex_generate", genname, sizeof(genname), &flg);CHKERRQ(ierr); 162 if (flg) name = genname; 163 } 164 165 fl = DMGenerateList; 166 if (name) { 167 while (fl) { 168 ierr = PetscStrcmp(fl->name,name,&flg);CHKERRQ(ierr); 169 if (flg) { 170 ierr = (*fl->generate)(boundary,interpolate,mesh);CHKERRQ(ierr); 171 PetscFunctionReturn(0); 172 } 173 fl = fl->next; 174 } 175 SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Grid generator %s not registered; you may need to add --download-%s to your ./configure options",name,name); 176 } else { 177 while (fl) { 178 if (boundary->dim == fl->dim) { 179 ierr = (*fl->generate)(boundary,interpolate,mesh);CHKERRQ(ierr); 180 PetscFunctionReturn(0); 181 } 182 fl = fl->next; 183 } 184 suggestions = ""; 185 if (boundary->dim+1 == 2) suggestions = " You may need to add --download-triangle to your ./configure options"; 186 else if (boundary->dim+1 == 3) suggestions = " You may need to add --download-ctetgen or --download-tetgen in your ./configure options"; 187 SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No grid generator of dimension %D registered%s",boundary->dim+1,suggestions); 188 } 189 } 190