[vtkusers] Solid Voxelization with VTK
Berti Krüger
berti_krueger at hotmail.com
Fri Sep 28 22:57:42 EDT 2018
Hello everyone.
One last update to this already long thread:
Meanwhile i implemented a majority voting approach using the voxelization function from the last post because i had problems voxelizing meshes with degenerated normals or polygons with wrong orientation:
[cid:35c9abf9-ea5c-48e4-98d5-db27ec909ccd]
Since the original algorithm from the last post more ore less simply "rasterizes" the polygon data in 3D along the z-direction, no voxelization artifacts from wrong normals result in this direction.
So to get rid of any artifacts in any direction, i voxelize the mesh three times in three different directions.
With the voxelization function given in the last post, i voxelize the mesh as it was given and save the result in an vtkImageData object. Then i simply rotate the original mesh 90° around the y-Axis, voxelize it again and save the result in another vtkImageData object. After that i rotate the original mesh (unrotated version) again but this time 90° around the x-Axis, voxelize it and save the result again in another vtkImageData object.
At the end i have three vtkImageData objects with results of the voxelization of the mesh in three different directions (along the x-axis, y-axis, z-axis).
The resulting vtkImageData objects are of course all rotated 90° in different directions too, so a direction comparision of the voxelization results is not possible. But by using rotations by 90° all that changes for a direct comparision is that the x, y, z components of the coordinates of the different results simply switch (e.g. x <-> z or x <-> y). So comparisions between the results can be done without any additional computations (no rotations of the vtkImageData or so).
I then compare each voxelization result with another. If at least two or all three voxelization results have a voxel set in a specific cell then the final result is that this cell is set. Otherwise it is unset.
The result:
[cid:784b647e-ef6e-461d-a16c-17b1974642a7]
It is even quite fast (only around three times slower than the simple voxelization step which is clear since it does three times the work) and much faster than the raycasting approach.
Btw: Does one of the filters (vtkPolyDataToImageStencil, vtkImageStencilToImage) use SMP?
I think it could be easily exploited because of the parallel nature of it.
Cheers,
Berti
________________________________
Von: vtkusers <vtkusers-bounces at public.kitware.com> im Auftrag von Berti Krüger <berti_krueger at hotmail.com>
Gesendet: Mittwoch, 26. September 2018 07:04
An: David Gobbi
Cc: VTK Users
Betreff: Re: [vtkusers] Solid Voxelization with VTK
Hi David,
thanks again for your help.
I added the method calls on the normalFilter as you have suggested but it didn't change the wrong voxelization of the Eiffel tower model.
It did indead change the result to something different:
[cid:ea753225-8561-4c0f-a7b0-2bb764dff664] [cid:a3cf6ebe-b37d-42ec-b0aa-fc502c405607] [cid:4439e9e1-f381-4580-8b22-005d383ff492]
left image: result with wrong orientation and without normalsFilter->AutoOrientNormalsOn()
center image: mesh with wrong orientation and normalsFilter->AutoOrientNormalsOn()
right image: meshlab corrected version
I attached both Eiffel tower models (with the wrong normals / polygon orientation and the corrected meshlab version) if you want to try it yourself.
The vtkStripper is indead unecessary. I tried various meshes with and without using the vtkStripper filter step inbetween and it didn't make any difference. Since i took it from the Slicer3D source code Csaba pointed me to, i don't have any idea why the Slicer Guys use it there. Maybe they have some special kind of polydata where it is needed. For my general purpose meshes it does not seem to make a difference. The processing is, of course, faster without this additional step.
The complete source code for voxelization and visualization and the CMakeList.txt is attached.
The resulting code for the voxelization function:
bool ConvertPolyDataToBinaryLabelMap(vtkSmartPointer<vtkPolyData> closedSurfacePolyData,
vtkSmartPointer<vtkImageData> binaryLabelMap)
{
// Check for consistency
if (closedSurfacePolyData->GetNumberOfPoints() < 2 || closedSurfacePolyData->GetNumberOfCells() < 2)
{
std::cout << "Convert: Cannot create binary labelmap from surface with number of points: "
<< closedSurfacePolyData->GetNumberOfPoints() << " and number of cells: "
<< closedSurfacePolyData->GetNumberOfCells() << std::endl;
return false;
}
// Compute polydata normals
//
// The purpose of the vtkPolyDataNormals filter here is to enforce a
// consistent orientation of the polygons (via the method ConsistencyOn)
// and to get an automatic determination of the correct orientation of the
// normals on the polydata surface (via the method AutoOrientNormalsOn).
//
// To increase the performance of the filter we turn off the splitting of
// sharp edges (via the method SplittingOff) since it is not necessary for
// getting a correct result here.
//
// The application of the filter on the polydata is purely optional if all
// the polygons are already correctly oriented.
vtkNew<vtkPolyDataNormals> normalFilter;
normalFilter->SetInputData(closedSurfacePolyData);
normalFilter->ConsistencyOn();
normalFilter->AutoOrientNormalsOn();
normalFilter->SplittingOff();
// Make sure that we have a clean triangle polydata using the
// vtkTriangleFilter which generates triangles from input polygons.
// Only needed if the input polygons are not flat.
vtkNew<vtkTriangleFilter> triangle;
triangle->SetInputConnection(normalFilter->GetOutputPort());
// Convert polydata to stencil
vtkNew<vtkPolyDataToImageStencil> polyDataToImageStencil;
polyDataToImageStencil->SetInputConnection(triangle->GetOutputPort());
polyDataToImageStencil->SetOutputSpacing(binaryLabelMap->GetSpacing());
polyDataToImageStencil->SetOutputOrigin(binaryLabelMap->GetOrigin());
polyDataToImageStencil->SetOutputWholeExtent(binaryLabelMap->GetExtent());
polyDataToImageStencil->Update();
// Convert stencil to image
vtkNew<vtkImageStencilToImage> imageStencilToImage;
imageStencilToImage->SetInputConnection(polyDataToImageStencil->GetOutputPort());
imageStencilToImage->SetOutsideValue(0);
imageStencilToImage->SetInsideValue(1);
imageStencilToImage->SetOutput(binaryLabelMap);
imageStencilToImage->Update();
return true;
}
Thanks again for everything.
Cheers,
Berti
________________________________
Von: David Gobbi <david.gobbi at gmail.com>
Gesendet: Dienstag, 25. September 2018 10:20
An: Berti Krüger
Cc: Csaba Pinter; VTK Users; Bill Lorensen
Betreff: Re: [vtkusers] Solid Voxelization with VTK
Hi Berti,
Just a few extra comments about the three "preprocessing" steps that you apply in your sample code:
For vtkPolyDataNormals, you might want to add normalsFilter->AutoOrientNormalsOn(), since this option is specifically designed to fix inside-out shapes like the Eiffel tower model. And normalsFilter->SplittingOff() will make this filter work faster. It would also be good if the comment mentioned that the purpose of the filter is to enforce consistent orientation of the polygons, and that the filter is optional if all the polygons are already correctly oriented.
The vtkStripper should be removed, unless you actually saw that it provides some benefit. The vtkPolyDataToImageStencil class should give exactly the same results for triangles as for strips. And vtkTriangleFilter is only needed if the input polygons are not flat.
Thanks again for the code,
- David
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <https://public.kitware.com/pipermail/vtkusers/attachments/20180929/62ef5c89/attachment-0001.html>
-------------- next part --------------
A non-text attachment was scrubbed...
Name: eiffel_tower_wrong_orientation.jpg
Type: image/jpeg
Size: 88782 bytes
Desc: eiffel_tower_wrong_orientation.jpg
URL: <https://public.kitware.com/pipermail/vtkusers/attachments/20180929/62ef5c89/attachment-0004.jpg>
-------------- next part --------------
A non-text attachment was scrubbed...
Name: eiffel_tower_right_orientation.jpg
Type: image/jpeg
Size: 37648 bytes
Desc: eiffel_tower_right_orientation.jpg
URL: <https://public.kitware.com/pipermail/vtkusers/attachments/20180929/62ef5c89/attachment-0005.jpg>
-------------- next part --------------
A non-text attachment was scrubbed...
Name: eiffel_tower_false.png
Type: image/png
Size: 54477 bytes
Desc: eiffel_tower_false.png
URL: <https://public.kitware.com/pipermail/vtkusers/attachments/20180929/62ef5c89/attachment-0001.png>
-------------- next part --------------
A non-text attachment was scrubbed...
Name: alienqueen_voxelization.jpg
Type: image/jpeg
Size: 65454 bytes
Desc: alienqueen_voxelization.jpg
URL: <https://public.kitware.com/pipermail/vtkusers/attachments/20180929/62ef5c89/attachment-0006.jpg>
-------------- next part --------------
A non-text attachment was scrubbed...
Name: alienqueen_majority_voting.jpg
Type: image/jpeg
Size: 101543 bytes
Desc: alienqueen_majority_voting.jpg
URL: <https://public.kitware.com/pipermail/vtkusers/attachments/20180929/62ef5c89/attachment-0007.jpg>
More information about the vtkusers
mailing list