[petsc-users] DMPlexComputeCellGeometryFVM: "Cannot handle faces with 1 vertices"

Matthew Knepley knepley at gmail.com
Tue Mar 12 08:08:31 CDT 2019


On Tue, Mar 12, 2019 at 9:03 AM Chris Finn via petsc-users <
petsc-users at mcs.anl.gov> wrote:

> Hello,
> with the code below, I create a tetrahedron using DMPlexCreateFromDAG,
> then I try to run DMPlexComputeCellGeometryFVM on this cell. The latter
> call fails with output:
>

All the geometry stuff requires that you interpolate the mesh I think, Just
use DMPlexInterpolate().

  Thanks,

    Matt


> ...
> [0]PETSC ERROR: --------------------- Error Message -----------------
> [0]PETSC ERROR: Argument out of range
> [0]PETSC ERROR: Cannot handle faces with 1 vertices
> ...
> (full output is attached).
>
> What am I doing wrong?
> regards
> Chris
>
> Here is the code:
> static char help[] = "No help \n";
>
> #include <petscdmplex.h>
>
> #undef __FUNCT__
> #define __FUNCT__ "main"
> int main(int argc,char **args)
> {
>   PetscErrorCode ierr;
>   PetscInitialize(&argc,&args,(char*)0,help);
>
>   DM                dm;
>   int cStart,cEnd;
>   int fStart,fEnd;
>
>   int depth = 3;
>   int dim = 3;
>
>   PetscInt    numPoints[4]        = {1,4,6,4};
>   PetscInt    coneSize[15]         = {4,3,3,3,3,2,2,2,2,2,2,0,0,0,0};
>   PetscInt    cones[28]            = {1,2,3,4, 5,9,8, 9,6,10, 10,8,7,
> 5,6,7, 11,12, 12,13, 13,11, 11,14, 12,14, 13,14};
>   PetscInt    coneOrientations[28] = {0 };
>   PetscScalar vertexCoords[12]     = {0,0,0, 1,0,0, 0,1,0, 0,0,1};
>
>   DMCreate(PETSC_COMM_WORLD, &dm);
>   DMSetType(dm, DMPLEX);
>   DMSetDimension(dm,dim);
>   DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones,
> coneOrientations, vertexCoords);
>
>   DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
>   DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);
>
>   for (int k =cStart;k<cEnd;k++){
>     double vol;
>     double centroid[3];
>     double normal[3];
>     ierr = DMPlexComputeCellGeometryFVM(dm, k, &vol,
> centroid,NULL);CHKERRQ(ierr);
>     printf("FVM: V=%f c=(%f %f %f) n=(%f %f
> %f)\n",vol,centroid[0],centroid[1],centroid[2],
>       normal[0],normal[1],normal[2]);
>   }
>   ierr = PetscFinalize();
>   return 0;
> }
>
>
>

-- 
What most experimenters take for granted before they begin their
experiments is infinitely more interesting than any results to which their
experiments lead.
-- Norbert Wiener

https://www.cse.buffalo.edu/~knepley/ <http://www.cse.buffalo.edu/~knepley/>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20190312/0a4462de/attachment.html>


More information about the petsc-users mailing list