[petsc-users] How to get cell number from node number?
Stefano Zampini
stefano.zampini at gmail.com
Fri Nov 30 08:40:11 CST 2018
Just to add to what Matt wrote: it seems your question is to get the cell numbers given a node (=vertex) numbering
In general, there are many cells that share the same vertex. If you have a vertex id (vertex), you can traverse the cells that share that vertex by still using DMPlexGetTransitiveClosure, with a PETSC_FALSE argument to get the so-called “star” of the point
ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
ierr = DMPlexGetTransitiveClosure(dm, vertex, PETSC_FALSE, &clSize, &closure);CHKERRQ(ierr);
for (cl = 0; cl < clSize*2; cl += 2) {
const PetscInt point = closure[cl];
if (point >= cStart && point < cEnd) {
<this is a cell>
}
}
ierr = DMPlexRestoreTransitiveClosure(dm, vertex, PETSC_FALSE, &clSize, &closure);CHKERRQ(ierr);
> On Nov 30, 2018, at 4:28 PM, Matthew Knepley <knepley at gmail.com> wrote:
>
> On Fri, Nov 30, 2018 at 8:19 AM barry via petsc-users <petsc-users at mcs.anl.gov <mailto:petsc-users at mcs.anl.gov>> wrote:
> Sorry for the imprecisely description.
> I have an unstructured grid (DMPLEX) with triangle mesh.
>
> cell number (2 dim)={0}; line number (1 dim)={1, 2, 3}; node number (0 dim)={4, 5, 6}
>
>
> So, you have a 2D interpolated simplex mesh in DMPlex, which means it has vertices, edges, and cells.
> If you want the vertices on a given cell, you would get the closure and then filter out any points which are
> not vertices. For example,
>
> PetscInt *closure = NULL;
> PetscInt clSize, cl;
>
> ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
> ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr);
> for (cl = 0; cl < clSize*2; cl += 2) {
> const PetscInt point = closure[cl];
>
> if (point >= vStart && point < vEnd) {
> <this is a vertex>
> }
> }
> ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr);
>
> If you do not use edges, then its much easier. Just call DMPlexGetCone(). If you do use edges, but also
> get the vertices all the time, you can make an index for this query, which will speed it up greatly at the
> cost of some memory.
>
> Thanks,
>
> Matt
> Thank you,
>
> Barry
>
> On 11/30/18 8:33 PM, Stefano Zampini wrote:
>>
>>
>> Il giorno Ven 30 Nov 2018, 15:05 Tsung-Hsing Chen via petsc-users <petsc-users at mcs.anl.gov <mailto:petsc-users at mcs.anl.gov>> ha scritto:
>> Hi,
>> Is there any function that can get cell number from the given node number exist already?
>>
>> Cell number of what? DMPLEX? And what is "node number"?
>>
>>
>> Thank you,
>> Barry
>
>
> --
> 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/20181130/1db08f29/attachment.html>
More information about the petsc-users
mailing list