1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2d9deefdfSMatthew G. Knepley 3d9deefdfSMatthew G. Knepley /*@C 496ca5757SLisandro Dalcin DMPlexInvertCell - Flips cell orientations since Plex stores some of them internally with outward normals. 5d9deefdfSMatthew G. Knepley 6d9deefdfSMatthew G. Knepley Input Parameters: 796ca5757SLisandro Dalcin + cellType - The cell type 8d9deefdfSMatthew G. Knepley - cone - The incoming cone 9d9deefdfSMatthew G. Knepley 10d9deefdfSMatthew G. Knepley Output Parameter: 11d9deefdfSMatthew G. Knepley . cone - The inverted cone (in-place) 12d9deefdfSMatthew G. Knepley 13d9deefdfSMatthew G. Knepley Level: developer 14d9deefdfSMatthew G. Knepley 15d9deefdfSMatthew G. Knepley .seealso: DMPlexGenerate() 16d9deefdfSMatthew G. Knepley @*/ 1796ca5757SLisandro Dalcin PetscErrorCode DMPlexInvertCell(DMPolytopeType cellType, PetscInt cone[]) 18d9deefdfSMatthew G. Knepley { 1996ca5757SLisandro Dalcin #define SWAPCONE(cone,i,j) \ 2096ca5757SLisandro Dalcin do { \ 2196ca5757SLisandro Dalcin PetscInt _cone_tmp; \ 2296ca5757SLisandro Dalcin _cone_tmp = (cone)[i]; \ 2396ca5757SLisandro Dalcin (cone)[i] = (cone)[j]; \ 2496ca5757SLisandro Dalcin (cone)[j] = _cone_tmp; \ 2596ca5757SLisandro Dalcin } while (0) 26d9deefdfSMatthew G. Knepley 27d9deefdfSMatthew G. Knepley PetscFunctionBegin; 2896ca5757SLisandro Dalcin switch (cellType) { 2996ca5757SLisandro Dalcin case DM_POLYTOPE_POINT: break; 3096ca5757SLisandro Dalcin case DM_POLYTOPE_SEGMENT: break; 3196ca5757SLisandro Dalcin case DM_POLYTOPE_POINT_PRISM_TENSOR: break; 3296ca5757SLisandro Dalcin case DM_POLYTOPE_TRIANGLE: break; 3396ca5757SLisandro Dalcin case DM_POLYTOPE_QUADRILATERAL: break; 3496ca5757SLisandro Dalcin case DM_POLYTOPE_SEG_PRISM_TENSOR: SWAPCONE(cone,2,3); break; 3596ca5757SLisandro Dalcin case DM_POLYTOPE_TETRAHEDRON: SWAPCONE(cone,0,1); break; 3696ca5757SLisandro Dalcin case DM_POLYTOPE_HEXAHEDRON: SWAPCONE(cone,1,3); break; 3796ca5757SLisandro Dalcin case DM_POLYTOPE_TRI_PRISM: SWAPCONE(cone,1,2); break; 3896ca5757SLisandro Dalcin case DM_POLYTOPE_TRI_PRISM_TENSOR: break; 3996ca5757SLisandro Dalcin case DM_POLYTOPE_QUAD_PRISM_TENSOR: break; 40d9deefdfSMatthew G. Knepley default: break; 41d9deefdfSMatthew G. Knepley } 42d9deefdfSMatthew G. Knepley PetscFunctionReturn(0); 4396ca5757SLisandro Dalcin #undef SWAPCONE 4496ca5757SLisandro Dalcin } 4596ca5757SLisandro Dalcin 4696ca5757SLisandro Dalcin /*@C 4796ca5757SLisandro Dalcin DMPlexReorderCell - Flips cell orientations since Plex stores some of them internally with outward normals. 4896ca5757SLisandro Dalcin 4996ca5757SLisandro Dalcin Input Parameters: 5096ca5757SLisandro Dalcin + dm - The DMPlex object 5196ca5757SLisandro Dalcin . cell - The cell 5296ca5757SLisandro Dalcin - cone - The incoming cone 5396ca5757SLisandro Dalcin 5496ca5757SLisandro Dalcin Output Parameter: 5596ca5757SLisandro Dalcin . cone - The reordered cone (in-place) 5696ca5757SLisandro Dalcin 5796ca5757SLisandro Dalcin Level: developer 5896ca5757SLisandro Dalcin 5996ca5757SLisandro Dalcin .seealso: DMPlexGenerate() 6096ca5757SLisandro Dalcin @*/ 6196ca5757SLisandro Dalcin PetscErrorCode DMPlexReorderCell(DM dm, PetscInt cell, PetscInt cone[]) 6296ca5757SLisandro Dalcin { 6396ca5757SLisandro Dalcin DMPolytopeType cellType; 6496ca5757SLisandro Dalcin 6596ca5757SLisandro Dalcin PetscFunctionBegin; 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCellType(dm, cell, &cellType)); 67*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexInvertCell(cellType, cone)); 6896ca5757SLisandro Dalcin PetscFunctionReturn(0); 69d9deefdfSMatthew G. Knepley } 70d9deefdfSMatthew G. Knepley 7194ef8ddeSSatish Balay /*@C 72d9deefdfSMatthew G. Knepley DMPlexTriangleSetOptions - Set the options used for the Triangle mesh generator 73d9deefdfSMatthew G. Knepley 74d9deefdfSMatthew G. Knepley Not Collective 75d9deefdfSMatthew G. Knepley 76d9deefdfSMatthew G. Knepley Inputs Parameters: 77d9deefdfSMatthew G. Knepley + dm - The DMPlex object 78d9deefdfSMatthew G. Knepley - opts - The command line options 79d9deefdfSMatthew G. Knepley 80d9deefdfSMatthew G. Knepley Level: developer 81d9deefdfSMatthew G. Knepley 82d9deefdfSMatthew G. Knepley .seealso: DMPlexTetgenSetOptions(), DMPlexGenerate() 83d9deefdfSMatthew G. Knepley @*/ 84d9deefdfSMatthew G. Knepley PetscErrorCode DMPlexTriangleSetOptions(DM dm, const char *opts) 85d9deefdfSMatthew G. Knepley { 86d9deefdfSMatthew G. Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 87d9deefdfSMatthew G. Knepley 88d9deefdfSMatthew G. Knepley PetscFunctionBegin; 89d9deefdfSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 90d9deefdfSMatthew G. Knepley PetscValidPointer(opts, 2); 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(mesh->triangleOpts)); 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrallocpy(opts, &mesh->triangleOpts)); 93d9deefdfSMatthew G. Knepley PetscFunctionReturn(0); 94d9deefdfSMatthew G. Knepley } 95d9deefdfSMatthew G. Knepley 9694ef8ddeSSatish Balay /*@C 97d9deefdfSMatthew G. Knepley DMPlexTetgenSetOptions - Set the options used for the Tetgen mesh generator 98d9deefdfSMatthew G. Knepley 99d9deefdfSMatthew G. Knepley Not Collective 100d9deefdfSMatthew G. Knepley 101d9deefdfSMatthew G. Knepley Inputs Parameters: 102d9deefdfSMatthew G. Knepley + dm - The DMPlex object 103d9deefdfSMatthew G. Knepley - opts - The command line options 104d9deefdfSMatthew G. Knepley 105d9deefdfSMatthew G. Knepley Level: developer 106d9deefdfSMatthew G. Knepley 107d9deefdfSMatthew G. Knepley .seealso: DMPlexTriangleSetOptions(), DMPlexGenerate() 108d9deefdfSMatthew G. Knepley @*/ 109d9deefdfSMatthew G. Knepley PetscErrorCode DMPlexTetgenSetOptions(DM dm, const char *opts) 110d9deefdfSMatthew G. Knepley { 111d9deefdfSMatthew G. Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 112d9deefdfSMatthew G. Knepley 113d9deefdfSMatthew G. Knepley PetscFunctionBegin; 114d9deefdfSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 115d9deefdfSMatthew G. Knepley PetscValidPointer(opts, 2); 116*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(mesh->tetgenOpts)); 117*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrallocpy(opts, &mesh->tetgenOpts)); 118d9deefdfSMatthew G. Knepley PetscFunctionReturn(0); 119d9deefdfSMatthew G. Knepley } 120d9deefdfSMatthew G. Knepley 121d9deefdfSMatthew G. Knepley /*@C 122d9deefdfSMatthew G. Knepley DMPlexGenerate - Generates a mesh. 123d9deefdfSMatthew G. Knepley 124d9deefdfSMatthew G. Knepley Not Collective 125d9deefdfSMatthew G. Knepley 126d9deefdfSMatthew G. Knepley Input Parameters: 127d9deefdfSMatthew G. Knepley + boundary - The DMPlex boundary object 128d9deefdfSMatthew G. Knepley . name - The mesh generation package name 129d9deefdfSMatthew G. Knepley - interpolate - Flag to create intermediate mesh elements 130d9deefdfSMatthew G. Knepley 131d9deefdfSMatthew G. Knepley Output Parameter: 132d9deefdfSMatthew G. Knepley . mesh - The DMPlex object 133d9deefdfSMatthew G. Knepley 1342c4dacf2SBarry Smith Options Database: 1352c4dacf2SBarry Smith + -dm_plex_generate <name> - package to generate mesh, for example, triangle, ctetgen or tetgen 136c0517cd5SMatthew G. Knepley - -dm_generator <name> - package to generate mesh, for example, triangle, ctetgen or tetgen 1372c4dacf2SBarry Smith 138d9deefdfSMatthew G. Knepley Level: intermediate 139d9deefdfSMatthew G. Knepley 140d9deefdfSMatthew G. Knepley .seealso: DMPlexCreate(), DMRefine() 141d9deefdfSMatthew G. Knepley @*/ 142d9deefdfSMatthew G. Knepley PetscErrorCode DMPlexGenerate(DM boundary, const char name[], PetscBool interpolate, DM *mesh) 143d9deefdfSMatthew G. Knepley { 144c0517cd5SMatthew G. Knepley DMGeneratorFunctionList fl; 1453f2a96e3SMatthew G. Knepley char genname[PETSC_MAX_PATH_LEN]; 1463f2a96e3SMatthew G. Knepley const char *suggestions; 147d9deefdfSMatthew G. Knepley PetscInt dim; 1483a074057SBarry Smith PetscBool flg; 149d9deefdfSMatthew G. Knepley 150d9deefdfSMatthew G. Knepley PetscFunctionBegin; 151d9deefdfSMatthew G. Knepley PetscValidHeaderSpecific(boundary, DM_CLASSID, 1); 152064a246eSJacob Faibussowitsch PetscValidLogicalCollectiveBool(boundary, interpolate, 3); 153*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(boundary, &dim)); 154*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(((PetscObject) boundary)->options,((PetscObject) boundary)->prefix, "-dm_generator", genname, sizeof(genname), &flg)); 155d9deefdfSMatthew G. Knepley if (flg) name = genname; 1562c4dacf2SBarry Smith else { 157*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(((PetscObject) boundary)->options,((PetscObject) boundary)->prefix, "-dm_plex_generate", genname, sizeof(genname), &flg)); 1582c4dacf2SBarry Smith if (flg) name = genname; 1592c4dacf2SBarry Smith } 1603a074057SBarry Smith 161c0517cd5SMatthew G. Knepley fl = DMGenerateList; 162d9deefdfSMatthew G. Knepley if (name) { 1633a074057SBarry Smith while (fl) { 164*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcmp(fl->name,name,&flg)); 1653a074057SBarry Smith if (flg) { 166*5f80ce2aSJacob Faibussowitsch CHKERRQ((*fl->generate)(boundary,interpolate,mesh)); 1673a074057SBarry Smith PetscFunctionReturn(0); 168d9deefdfSMatthew G. Knepley } 1693a074057SBarry Smith fl = fl->next; 170d9deefdfSMatthew G. Knepley } 17198921bdaSJacob Faibussowitsch SETERRQ(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); 1723a074057SBarry Smith } else { 1733a074057SBarry Smith while (fl) { 1743a074057SBarry Smith if (boundary->dim == fl->dim) { 175*5f80ce2aSJacob Faibussowitsch CHKERRQ((*fl->generate)(boundary,interpolate,mesh)); 1763a074057SBarry Smith PetscFunctionReturn(0); 1773a074057SBarry Smith } 1783a074057SBarry Smith fl = fl->next; 1793a074057SBarry Smith } 1802c4dacf2SBarry Smith suggestions = ""; 1812c4dacf2SBarry Smith if (boundary->dim+1 == 2) suggestions = " You may need to add --download-triangle to your ./configure options"; 1822c4dacf2SBarry Smith else if (boundary->dim+1 == 3) suggestions = " You may need to add --download-ctetgen or --download-tetgen in your ./configure options"; 18398921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No grid generator of dimension %D registered%s",boundary->dim+1,suggestions); 1843a074057SBarry Smith } 1853a074057SBarry Smith } 186