Skip to Content.
Sympa Menu

cgal-discuss - Re: [cgal-discuss] Slicing/splitting/cutting polyhedral meshes

Subject: CGAL users discussion list

List archive

Re: [cgal-discuss] Slicing/splitting/cutting polyhedral meshes


Chronological Thread 
  • From: Vinicius Azevedo <>
  • To: cgal-discuss <>
  • Subject: Re: [cgal-discuss] Slicing/splitting/cutting polyhedral meshes
  • Date: Wed, 29 Apr 2015 14:48:22 -0400

OK, so I think that I found a potential problem. Inside the polyhedron_cut_plane_3 function we have:

Halfedge_handle start = h;
do {
   h->vertex()->point() = construct_point( plane, h->next()->opposite()->facet()->plane(), h->opposite()->facet()->plane());
   h = h->next();
} while ( h != start);

This function basically computes the intersection point from 3 planes: the cutting plane and the two planes that pass through the current edge. The function construct_point is:

 operator()( const Plane_3& p, const Plane_3& q, const Plane_3& r) const { 
   Object obj = intersection( p, q);
   Line_3 line;
   if ( assign( line, obj)) {
        obj = intersection( r, line);
        Point_3 pt;
        if ( assign( pt, obj)) {
            return pt;
        }
   }
   std::cerr << "ERROR: coplanar planes used for computing " "intersecting point." << std::endl;
   std::cerr << "       Return ORIGIN. Don't trust result." << std::endl;
   return ORIGIN;
}

This function tries to find the intersection point by successively computing intersections between the planes. But this algorithm will fail if r and q planes are colinear. That means that if an edge has two colinear planes (e.g. a simple plane), that algorithm will not work. So perhaps a function that better computes the intersection between a plane and an edge has to be used to solved this problem

Any ideas? Is what I said correct?

Vinicius.
 

Vinicius C. Azevedo.

On Wed, Apr 29, 2015 at 11:15 AM, Vinicius Azevedo <> wrote:
Sure, here's the code snippet: http://pastebin.com/C892FnCu

The interesting thing though is that the halfedge handle that I'm passing is not one that directly crosses the polyline intersection.If I try to pass an edge that is being actually intersected by the polyline I get an assertion error.

Vinicius C. Azevedo.

On Wed, Apr 29, 2015 at 11:01 AM, Sebastien Loriot (GeometryFactory) <> wrote:
On 04/29/2015 04:56 PM, Vinicius Azevedo wrote:
I was confused by the cut/slice examples on the Polyhedron demo. I
thought you was mentioning that, my bad. That's what I was looking for,
thanks! I have some questions though on how to use this function:

  * What should I pass to the function as starting halfedge? Is it an
    edge returned by any_intersection?
I think so.

  * Does this function add vertices from the polyline that is generated
    from the intersection of the polyhedron with the plane? That's not
    what I'm getting on my results
    (https://dl.dropboxusercontent.com/u/5919734/bunnyCut1.png &&
    https://dl.dropboxusercontent.com/u/5919734/bunnyCut2.png). It only
    cuts the mesh by following pre-existent vertices. Also I'm having
    the error "coplanar planes used for computing the intersection
    plane, Return ORIGIN, Don't trust the result".

I never used the function. Could you share your code so that I can
give it a try?

Sebastien.

Thanks for your help so far!

Vinicius C. Azevedo.

On Tue, Apr 28, 2015 at 10:24 AM, Sebastien Loriot (GeometryFactory)
< <mailto:>> wrote:

    On 04/28/2015 02:15 PM, Vinicius Azevedo wrote:

        Hi Sebastian,

        Thanks for your answer. I didn't expressed myself so well. What
        I want
        is to actually cut the mesh into two different independent parts
        with a
        given cut plane. In this video
        https://www.youtube.com/watch?v=IOnk7PBv1AU and paper
        http://www.math.ucla.edu/~jteran/papers/WJST14.pdf, they present
        some
        techniques to do that robustly. I saw the examples that you
        mentioned,
        and it seems that they only compute the intersections of planes and
        Polyhedrons, but don't actually split it into different parts.


    I believe this is what polyhedron_cut_plane_3.h is doing. The comment
    for this function is:
    // Cuts the polyhedron `poly' at plane `plane' starting at halfedge `h'.
    // Traces the intersection curve of `plane' with `poly' starting at `h',
    // cuts `poly' along that intersection curve and deletes the (connected)
    // component on the negative side of the plane. The hole created along
    // the intersection curve is filled with a new face containing the plane
    // and the points in the vertices are computed.

    Sebastien.


        Vinicius C. Azevedo.

        On Tue, Apr 28, 2015 at 1:48 AM, Sebastien Loriot (GeometryFactory)
        < <mailto:>
        <mailto: <mailto:>>> wrote:

             On 04/27/2015 04:31 PM, Vinicius Azevedo wrote:

                 Hi all,

                 I would like to know if anyone has a code snippet for
                 slicing/splitting/cutting polyhedral meshes in two
        different
                 parts with
                 arbitrary planes. I've seen some old threads about it

        (http://cgal-discuss.949826.n4.nabble.com/Slices-in-a-mesh-td951695.html),
                 but didn't find any code examples so far.

                 Thanks!
                 Vinicius C. Azevedo.



             The thread you mention is to compute the intersection of a
        polyhedron
             with a plane which is one or several polygons. There is a
        non documented
             functor that is doing this in CGAL/Polyhedron_slicer_3.h
             A better and documented version will be included in CGAL 4.8

             In your case, it seems you want the intersection of a
        polyhedron with
             a halfspace.
             There is another documented function
        'CGAL::polyhedron_cut_plane_3'
             I never used that is defined in CGAL/polyhedron_cut_plane_3.h
             that should be doing what you want.

             Sebastien.

             --
             You are currently subscribed to cgal-discuss.
             To unsubscribe or access the archives, go to
        https://sympa.inria.fr/sympa/info/cgal-discuss





    --
    You are currently subscribed to cgal-discuss.
    To unsubscribe or access the archives, go to
    https://sympa.inria.fr/sympa/info/cgal-discuss





--
You are currently subscribed to cgal-discuss.
To unsubscribe or access the archives, go to
https://sympa.inria.fr/sympa/info/cgal-discuss







Archive powered by MHonArc 2.6.18.

Top of Page