[vtkusers] Export mesh of cutting plane

Fernando Nellmeldin f.nellmeldin at open-engineering.com
Fri Sep 29 07:23:05 EDT 2017


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

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


More information about the vtkusers mailing list