[vtk-developers] vtkCellArray memory layout
Berk Geveci
berk.geveci at kitware.com
Wed Mar 20 19:35:54 EDT 2013
Actually, now I remember having a discussion with Will where he explained
the reasoning behind this design. In vtkPolyData the locations array is
created only as needed. In rendering, all a mapper needs to do is to
iterate over triangles or polygons to render them. The vtkCellArray
structure is a good fit for this. Without needing locations.
vtkUnstructuredGrid adds the locations array to provide random access to
vtkCellArray.
Given this particular use case, maybe we should consider splitting the
implementation of cell arrays for polydata and unstructured grids? This
would allow for separate optimized implementations to address each data
object's needs.
PS: I get the feeling that this discussion is going to get out of hand
pretty quickly since redesigning these data structures requires a lot of
thought. As much as I enjoy a good discussion on data structures, I am not
sure that we are prepared to go this deep in changing core data structures
in VTK right now.
On Wed, Mar 20, 2013 at 5:22 PM, David Lonie <david.lonie at kitware.com>wrote:
> Hi list,
>
> I have some questions about the vtkCellArray implementation, as well as
> some concerns about it after seeing some widespread abuses of its
> "internal" access methods.
>
> First, some background:
>
> The vtkCellArray currently stores cell information in the form
>
> [ cellSize1 id id ... cellSize2 id id ... cellSizeN id id ... ]
>
> where the cell sizes tell the number of points in each cell, and the ids
> are offsets into an array of points. So a real vtkCellArray may look
> something like:
>
> [(2) 0 1 (4) 3 2 0 4 (3) 1 3 5 (3) 3 5 7]
>
> where the numbers in parenthesis are the sizes of the cells.
>
> The docs do point out the current implementation is "totally inadequate"
> for random access, and I'm trying to understand what benefit it offers at
> the expense of random cell lookup? There is a vtkCellTypes class that will
> build up an offset array like the one I propose above, but this approach
> effectively doubles the amount of memory needed to store the cell sizes
> (per vtkCellTypes instance, which is often used by parent containers as
> well as filters), and time must be spent generating these supplemental
> structures each time the data changes.
>
> It seems to me that another implementation would be more useful. In
> particular, splitting the array into two, one with cell offsets and another
> with just point ids, transforming the above list into:
>
> [0 2 6 9 12] (offsets, last entry points to end of ids)
> [(0 1) (3 2 0 4) (1 3 5) (3 5 7)] (ids)
>
> The parentheses in ids just group cells for readability.
>
> This offers random access by cell id for free, and eliminates the
> time/memory spent building the supplemental vtkCellTypes objects. There is
> a slightly higher penalty for sequential access than before, due to having
> to iterate through two arrays, but I doubt that this would have serious
> performance impacts.
>
> I've also noticed that there are a number of "random access" methods like
> vtkCellArray::GetCell(vtkIdType loc, ...) that take an offset into the
> "mixed" internal array. This breaks encapsulation and is rather unintuitive
> -- I think it's safe to say that most developers would expect the "loc"
> argument to be a cell id. Since these methods depend on implementation
> details they are labeled as internal in the docs. However, they have public
> visibility and are used commonly throughout various filters, which often
> build up a structure similar to vtkCellTypes in order to perform their
> lookups.
>
> I'm asking at this time because I'm looking at some changes that may
> modify the vtkCellArray API/implementation anyway, and it'd be nice to fix
> up this mix of sequential access, unintuitive random-access, and auxiliary
> lookup structures into a more friendly and intuitive API. These changes
> would not be effected until after VTK 6.0 is out.
>
> Thoughts? Discussions? Alternate perspectives? All are welcome :-)
>
> Thanks,
> Dave
>
> _______________________________________________
> Powered by www.kitware.com
>
> Visit other Kitware open-source projects at
> http://www.kitware.com/opensource/opensource.html
>
> Follow this link to subscribe/unsubscribe:
> http://www.vtk.org/mailman/listinfo/vtk-developers
>
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://public.kitware.com/pipermail/vtk-developers/attachments/20130320/7e13aca6/attachment.html>
More information about the vtk-developers
mailing list