<div dir="ltr"><div dir="ltr">On Tue, Mar 12, 2019 at 9:25 AM <<a href="mailto:finnkochinski@keemail.me">finnkochinski@keemail.me</a>> wrote:<br></div><div class="gmail_quote"><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex">
<div>
<div>Thank you,<br></div><div>I added DMPlexInterpolate(). Now I get a different error in DMPlexInterpolate():<br></div><div>...<br></div></div></blockquote><div><br></div><div>Okay, you did not need interpolate. You already specified all the levels. However, your orientations are wrong.</div><div>Get rid of that code.</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div><div></div><div>[0]PETSC ERROR: Dimension 0 not supported<br></div><div>...<br></div><div><br></div><div>The source code and full output are attached. Anybody able to fix this?<br></div></div></blockquote><div><br></div><div>The problem is a numbering convention in the library. Plex will accept any consistent DAG. However, if you</div><div>want to use other things in the library, like geometry routines, then there is a convention on numbering (which</div><div>makes many things simpler). We require that you number contiguously:</div><div><br></div><div> Cells: [0, Nc)</div><div> Vertices: [Nc, Nc+Nv)</div><div> Edges: [Nc+Nv, Nc+Nv+Ne)</div><div> Faces: [Nc+Nv+Ne, Nc+Nv+Ne+Nf)</div><div><br></div><div>You numbered the faces in the vertices slot, so the geometry routines got confused.</div><div><br></div><div> Thanks,</div><div><br></div><div> Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div><div></div><div>regards<br></div><div>Chris<br></div><div><br></div><div>-- <br></div><div> Securely sent with Tutanota. Get your own encrypted, ad-free mailbox: <br></div><div> <a href="https://tutanota.com" target="_blank">https://tutanota.com</a><br></div><div><br></div><div><br></div><div>Mar 12, 2019, 1:08 PM by <a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>:<br></div><blockquote class="gmail-m_1585930363849037245tutanota_quote" style="border-left:1px solid rgb(147,163,184);padding-left:10px;margin-left:5px"><div dir="ltr"><div dir="ltr">On Tue, Mar 12, 2019 at 9:03 AM Chris Finn via petsc-users <<a href="mailto:petsc-users@mcs.anl.gov" rel="noopener noreferrer" target="_blank">petsc-users@mcs.anl.gov</a>> wrote:<br></div><div><blockquote style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div><div>Hello,<br></div><div>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:<br></div></div></blockquote><div><br></div><div>All the geometry stuff requires that you interpolate the mesh I think, Just use DMPlexInterpolate().<br></div><div><br></div><div> Thanks,<br></div><div><br></div><div> Matt<br></div><div> <br></div><blockquote style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div><div><br></div><div>...<br></div><div>[0]PETSC ERROR: --------------------- Error Message -----------------<br></div><div>[0]PETSC ERROR: Argument out of range<br></div><div>[0]PETSC ERROR: Cannot handle faces with 1 vertices<br></div><div>...<br></div><div>(full output is attached).<br></div><div><br></div><div>What am I doing wrong?<br></div><div>regards<br></div><div>Chris<br></div><div><br></div><div>Here is the code:<br></div><div>static char help[] = "No help \n";<br></div><div><br></div><div>#include <petscdmplex.h><br></div><div><br></div><div>#undef __FUNCT__<br></div><div>#define __FUNCT__ "main"<br></div><div>int main(int argc,char **args)<br></div><div>{<br></div><div> PetscErrorCode ierr;<br></div><div> PetscInitialize(&argc,&args,(char*)0,help);<br></div><div><br></div><div> DM dm;<br></div><div> int cStart,cEnd;<br></div><div> int fStart,fEnd;<br></div><div><br></div><div> int depth = 3;<br></div><div> int dim = 3;<br></div><div><br></div><div> PetscInt numPoints[4] = {1,4,6,4};<br></div><div> PetscInt coneSize[15] = {4,3,3,3,3,2,2,2,2,2,2,0,0,0,0};<br></div><div> 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};<br></div><div> PetscInt coneOrientations[28] = {0 };<br></div><div> PetscScalar vertexCoords[12] = {0,0,0, 1,0,0, 0,1,0, 0,0,1};<br></div><div><br></div><div> DMCreate(PETSC_COMM_WORLD, &dm);<br></div><div> DMSetType(dm, DMPLEX);<br></div><div> DMSetDimension(dm,dim);<br></div><div> DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords);<br></div><div><br></div><div> DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);<br></div><div> DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);<br></div><div><br></div><div> for (int k =cStart;k<cEnd;k++){<br></div><div> double vol;<br></div><div> double centroid[3];<br></div><div> double normal[3];<br></div><div> ierr = DMPlexComputeCellGeometryFVM(dm, k, &vol, centroid,NULL);CHKERRQ(ierr);<br></div><div> printf("FVM: V=%f c=(%f %f %f) n=(%f %f %f)\n",vol,centroid[0],centroid[1],centroid[2],<br></div><div> normal[0],normal[1],normal[2]);<br></div><div> }<br></div><div> ierr = PetscFinalize();<br></div><div> return 0;<br></div><div>}<br></div><div><br></div><div><br></div></div></blockquote></div><div><br></div><div><br></div><div>-- <br></div><div dir="ltr"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><div><div>What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br></div><div>-- Norbert Wiener<br></div></div><div><br></div><div><a href="http://www.cse.buffalo.edu/~knepley/" rel="noopener noreferrer" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div></blockquote><div><br></div> </div>
</blockquote></div><br clear="all"><div><br></div>-- <br><div dir="ltr" class="gmail_signature"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><div>What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>-- Norbert Wiener</div><div><br></div><div><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div>