[vtkusers] Export mesh of cutting plane

Fernando Nellmeldin f.nellmeldin at open-engineering.com
Thu Oct 19 04:20:00 EDT 2017


Hello. Thank you for your answer. The method 1 worked !
I don't like method 2 because it depends a lot on the input data, which can
change... (that angle thing is a little too much for our case).
We'll use method 1 for the moment!
My problem was that the input of vtkCutter wasn't the output of
vtkTableBasedClipDataSet but the original model itself, then the mesh in
the cut got calculated twice and using different methods.

Thanks a lot!

On 19 October 2017 at 08:04, kenichiro yoshimi <rccm.kyoshimi at gmail.com>
wrote:

> Hello,
>
> Firstly the cutting algorithms are different between
> vtkTableBasedClipDataSet and vtkCutter, thus in general the generated
> cross sections are quite different.
>
> I have no idea how to extract a section surface directly from output
> of vtkTableBasedClipDataSet, however here are two insufficient ways:
> 1) vtkTableBasedClipDataSet -> vtkCutter (Cutter.tar.gz)
> 2) vtkTableBasedClipDataSet -> vtkGeometryFilter -> vtkPolyDataNormals
> -> vtkConnectivityFilter (Cutter2.tar.gz)
>
> ---
> 1) This method conducts double cutting using the same plane in turn.
> But in order to avoid the broken polygons as output of vtkCutter, the
> position of a cutting plain need to be slightly shifted from the used
> in vtkTableBasedClipDataSet. This means good accuracy of results are
> not be retained.
> 2) This method divides the surface by feature edges (vtkPolyDataNormals)
> and then selects the desired surface among them based on polygonal
> connectivity (vtkConnectivityFilter). Here, the origin of a cutting
> plane need to be known and be only on the cross section of objects so
> as to extract the cross section according to the criteria "extract the
> region closest to the origin". In addition, in vtkPolyDataNormals the
> feature angle to determine feature edges need to be adjusted in
> accordance with each dataset.
>
> I'm sorry for my poor English.
>
> 2017-10-18 23:00 GMT+09:00 Fernando Nellmeldin
> <f.nellmeldin at open-engineering.com>:
> > Hello, sorry for coming back to this topic, but I have found some cases
> > where this doesn't work.
> >
> > I attached a screenshot.
> > On top is the 3d mesh cut by a plane parallel to the screen. We can see
> that
> > are some elements (hexahedrons) that got cut and they were split in sets
> of
> > tetrahedrons.
> > At the bottom, we see the result of calling vtkCutter ->
> > vtkCleanUnstructuredGrid -> vtkXMLUnstructuredGridWriter. The result are
> all
> > quads and all the generated nodes by the cut are gone, and we see no
> > triangles.
> > I would expect to have the same result that we see on the top, but with
> 2d
> > elements instead of 3d.
> >
> > Here's my code:
> >
> > void cutModel(vtkSmartPointer<vtkUnstructuredGrid> model, double
> normal[3])
> > {
> > vtkSmartPointer<vtkPlane> plane = vtkSmartPointer<vtkPlane>::New();
> > plane->SetOrigin(0.0, 0.0, 0.0);
> > plane->SetNormal(normal);
> >
> > vtkSmartPointer<vtkCutter> cutter = vtkSmartPointer<vtkCutter>::New();
> > cutter->SetCutFunction(plane);
> > cutter->SetInputData(model);
> > cutter->GenerateTrianglesOff();
> > cutter->Update();
> >
> > vtkSmartPointer<vtkCleanUnstructuredGrid> cleaner =
> > vtkSmartPointer<vtkCleanUnstructuredGrid>::New();
> > cleaner->SetInputConnection(cutter->GetOutputPort());
> > cleaner->Update();
> > vtkSmartPointer<vtkXMLUnstructuredGridWriter> writer =
> > vtkSmartPointer<vtkXMLUnstructuredGridWriter>::New();
> > writer->SetFileName("planecut.vtu");
> > writer->SetInputConnection(cleaner->GetOutputPort());
> > writer->SetDataModeToAscii();
> > writer->Write();
> > }
> >
> >
> >
> > Thank you.
> >
> >
> > On 29 September 2017 at 13:23, Fernando Nellmeldin
> > <f.nellmeldin at open-engineering.com> wrote:
> >>
> >> Hello, that worked as a charm, thank you very much for your time!!
> >>
> >> On 29 September 2017 at 11:31, kenichiro yoshimi <
> rccm.kyoshimi at gmail.com>
> >> wrote:
> >>>
> >>> Hi Fernando,
> >>>
> >>> This will be accomplished by using vtkCleanUnstructuredGrid class
> >>> found in ParaView without triangulation. I have attached a simple
> >>> code.
> >>>
> >>> Thanks
> >>>
> >>> 2017-09-29 16:37 GMT+09:00 Fernando Nellmeldin
> >>> <f.nellmeldin at open-engineering.com>:
> >>> > Hello. That did the job! (although I get a triangle mesh while the
> >>> > initial
> >>> > surface was full of sliced hexahedron (quads).
> >>> >
> >>> > Thank you!
> >>> >
> >>> > On 29 September 2017 at 04:05, kenichiro yoshimi
> >>> > <rccm.kyoshimi at gmail.com>
> >>> > wrote:
> >>> >>
> >>> >> Hello,
> >>> >>
> >>> >> Using vtkDataSetTriangleFilter is a simple way which allows you to
> >>> >> convert polyData to unstructuredGrid. Something like this
> (VTK-8.0.0):
> >>> >> ---
> >>> >>   // Create cutter
> >>> >>   vtkSmartPointer<vtkCutter> cutter =
> >>> >>     vtkSmartPointer<vtkCutter>::New();
> >>> >>   cutter->SetCutFunction(plane);
> >>> >>   cutter->SetInputConnection(model->GetOutputPort());
> >>> >>
> >>> >>   // Convert polyData to unstructuredGrid
> >>> >>   vtkSmartPointer<vtkDataSetTriangleFilter> triFilter =
> >>> >>     vtkSmartPointer<vtkDataSetTriangleFilter>::New();
> >>> >>   triFilter->SetInputConnection(cutter->GetOutputPort());
> >>> >>   triFilter->Update();
> >>> >>
> >>> >>   vtkSmartPointer<vtkXMLUnstructuredGridWriter> writer =
> >>> >>     vtkSmartPointer<vtkXMLUnstructuredGridWriter>::New();
> >>> >>   writer->SetFileName("unstructuredGrid.vtu");
> >>> >>   writer->SetInputConnection(triFilter->GetOutputPort());
> >>> >>   writer->Write();
> >>> >> ---
> >>> >>
> >>> >> Regards
> >>> >>
> >>> >> 2017-09-28 23:43 GMT+09:00 Fernando Nellmeldin
> >>> >> <f.nellmeldin at open-engineering.com>:
> >>> >> > Hello.
> >>> >> > I have a method to cut a volume in a vtkUnstructuredGrid given a
> >>> >> > vtkPlane.
> >>> >> > I'm using vtkTableBasedClipDataSet to extract half of the model,
> >>> >> > also as
> >>> >> > vtkUnstructuredGrid. This is working.
> >>> >> >
> >>> >> > Now, I would like to obtain ONLY the surface of the cut, than
> means:
> >>> >> > I
> >>> >> > want
> >>> >> > a vtkUnstructuredGrid with the nodes and faces being on the
> surface
> >>> >> > of
> >>> >> > the
> >>> >> > cut.
> >>> >> > How can this be achieved?
> >>> >> >
> >>> >> > I tried using vtkCutter with vtkStripper but the mesh I get
> doesn't
> >>> >> > have
> >>> >> > any
> >>> >> > cells, only the points and their associated data.
> >>> >> >
> >>> >> > Here's how my code looks like.
> >>> >> >
> >>> >> > Thank You!
> >>> >> >
> >>> >> > vtkSmartPointer<vtkUnstructuredGrid> model =
> >>> >> > vtkSmartPointer<vtkUnstructuredGrid>::New();
> >>> >> > // fill the ugrid ...
> >>> >> >
> >>> >> > // Create The Plane to cut
> >>> >> > vtkSmartPointer<vtkPlane> plane = vtkSmartPointer<vtkPlane>::
> New();
> >>> >> > plane->SetOrigin(0.0, 0.0, 0.0);
> >>> >> > plane->SetNormal(1.0, 0, 0);
> >>> >> >
> >>> >> > // Code to cut the model in half
> >>> >> > vtkSmartPointer<vtkTableBasedClipDataSet> clipDataSet =
> >>> >> > vtkSmartPointer<vtkTableBasedClipDataSet>::New();
> >>> >> > clipDataSet->SetClipFunction(plane);
> >>> >> > clipDataSet->InsideOutOn();
> >>> >> > clipDataSet->GenerateClippedOutputOn();
> >>> >> > clipDataSet->SetInputConnection(model->GetProducerPort());
> >>> >> > clipDataSet->Update();
> >>> >> >
> >>> >> > // This has half the model, this is Ok but not what I want
> >>> >> > vtkSmartPointer<vtkUnstructuredGrid> halfModel =
> >>> >> > clipDataSet->GetOutput();
> >>> >> >
> >>> >> > ////****////
> >>> >> > // BEGINS code not working to extract the surface of the cut
> >>> >> > vtkSmartPointer<vtkCutter> cutter =
> >>> >> > vtkSmartPointer<vtkCutter>::New();
> >>> >> > cutter->SetCutFunction(plane);
> >>> >> > cutter->SetInputConnection(model->GetProducerPort());
> >>> >> > cutter->Update();
> >>> >> >
> >>> >> > vtkSmartPointer<vtkStripper> stripper =
> >>> >> > vtkSmartPointer<vtkStripper>::New();
> >>> >> > stripper->SetInputConnection(cutter->GetOutputPort());
> >>> >> > stripper->Update();
> >>> >> >
> >>> >> > vtkSmartPointer<vtkPolyData> pd =
> >>> >> > vtkSmartPointer<vtkPolyData>::New();
> >>> >> > pd->SetPoints(stripper->GetOutput()->GetPoints());
> >>> >> > pd->SetPolys(stripper->GetOutput()->GetLines());
> >>> >> >
> >>> >> > vtkSmartPointer<vtkUnstructuredGrid> clipped =
> >>> >> > vtkSmartPointer<vtkUnstructuredGrid>::New();
> >>> >> > clipped->ShallowCopy(pd);
> >>> >> >
> >>> >> > // Here I see that there is no Cells
> >>> >> > clipped->Print(std::cout);
> >>> >> >
> >>> >> > // ENDS code to extract the surface of the cut
> >>> >> > ////****////
> >>> >> >
> >>> >> > // Write the unstructured grid
> >>> >> > vtkSmartPointer<vtkXMLUnstructuredGridWriter> writer =
> >>> >> > vtkSmartPointer<vtkXMLUnstructuredGridWriter>::New();
> >>> >> > writer->SetFileName("d:/unstructuredGrid.vtu");
> >>> >> > writer->SetInput(clipped);
> >>> >> > writer->SetDataModeToAscii();
> >>> >> > writer->Write();
> >>> >> > // file doesnt have any cell, only points and pointdata
> >>> >> >
> >>> >> > _______________________________________________
> >>> >> > Powered by www.kitware.com
> >>> >> >
> >>> >> > Visit other Kitware open-source projects at
> >>> >> > http://www.kitware.com/opensource/opensource.html
> >>> >> >
> >>> >> > Please keep messages on-topic and check the VTK FAQ at:
> >>> >> > http://www.vtk.org/Wiki/VTK_FAQ
> >>> >> >
> >>> >> > Search the list archives at: http://markmail.org/search/?q=
> vtkusers
> >>> >> >
> >>> >> > Follow this link to subscribe/unsubscribe:
> >>> >> > http://public.kitware.com/mailman/listinfo/vtkusers
> >>> >> >
> >>> >
> >>> >
> >>> >
> >>> >
> >>> > --
> >>> > Fernando NELLMELDIN
> >>> > Software Engineer
> >>> > _______________________________________________________________
> >>> >
> >>> > Open Engineering s.a.
> >>> >
> >>> > Rue Bois Saint-Jean 15/1
> >>> > B-4102 Seraing (Belgium)
> >>> > Tel: +32.4.353.30.34
> >>> >
> >>> > http://www.open-engineering.com
> >>> > https://www.linkedin.com/company/open-engineering?trk=
> biz-companies-cym
> >>> >
> >>> > ____________________________________________________________
> _____________
> >>
> >>
> >>
> >>
> >> --
> >> Fernando NELLMELDIN
> >> Software Engineer
> >> _______________________________________________________________
> >>
> >> Open Engineering s.a.
> >>
> >> Rue Bois Saint-Jean 15/1
> >> B-4102 Seraing (Belgium)
> >> Tel: +32.4.353.30.34
> >>
> >> http://www.open-engineering.com
> >> https://www.linkedin.com/company/open-engineering?trk=biz-companies-cym
> >> ____________________________________________________________
> _____________
> >
> >
> >
> >
> > --
> > Fernando NELLMELDIN
> > Software Engineer
> > _______________________________________________________________
> >
> > Open Engineering s.a.
> >
> > Rue Bois Saint-Jean 15/1
> > B-4102 Seraing (Belgium)
> > Tel: +32.4.353.30.34
> >
> > http://www.open-engineering.com
> > https://www.linkedin.com/company/open-engineering?trk=biz-companies-cym
> > ____________________________________________________________
> _____________
>



-- 
*Fernando NELLMELDIN*
Software Engineer
*_______________________________________________________________*

*Open Engineering s.a.*

Rue Bois Saint-Jean 15/1
B-4102 Seraing (Belgium)
Tel: +32.4.353.30.34

http://www.open-engineering.com
https://www.linkedin.com/company/open-engineering?trk=biz-companies-cym

*_________________________________________________________________________*
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://public.kitware.com/pipermail/vtkusers/attachments/20171019/9c4479f5/attachment.html>


More information about the vtkusers mailing list