cgal/Straight_skeleton_2/include/CGAL/Polygon_offset_builder_2.h

205 lines
7.0 KiB
C++

// Copyright (c) 2005-2008 Fernando Luis Cacciola Carballal. All rights reserved.
//
// This file is part of CGAL (www.cgal.org).
//
// $URL$
// $Id$
// SPDX-License-Identifier: GPL-3.0-or-later OR LicenseRef-Commercial
//
// Author(s) : Fernando Cacciola <fernando_cacciola@ciudad.com.ar>
//
#ifndef CGAL_POLYGON_OFFSET_BUILDER_2_H
#define CGAL_POLYGON_OFFSET_BUILDER_2_H 1
#include <CGAL/license/Straight_skeleton_2.h>
#include <CGAL/Straight_skeleton_2/Straight_skeleton_aux.h>
#include <CGAL/Polygon_offset_builder_traits_2.h>
#include <optional>
#include <memory>
#include <vector>
namespace CGAL {
template<class Traits_, class SSkel_>
struct Default_polygon_offset_builder_2_visitor
{
typedef Traits_ Traits ;
typedef SSkel_ SSkel ;
typedef typename SSkel::Halfedge_const_handle Halfedge_const_handle ;
typedef typename SSkel::Vertex_const_handle Vertex_const_handle ;
typedef typename Traits::FT FT ;
typedef typename Traits::Point_2 Point_2 ;
void on_construction_started ( FT ) const {}
void on_offset_contour_started() const {}
void on_offset_point ( Point_2 const&, Halfedge_const_handle ) const {}
Point_2 on_offset_point_overflowed( Halfedge_const_handle aBisector ) const
{
::CGAL::warning_fail( "", __FILE__, __LINE__, "! Unable to compute polygon offset point due to overflow !" ) ;
return aBisector->vertex()->point() ;
}
void on_offset_contour_finished ( bool /*is_complete*/ ) const {}
void on_construction_finished () const {}
void on_error( char const* ) const {}
} ;
template<class Ss_, class Traits_, class Container_, class Visitor_ = Default_polygon_offset_builder_2_visitor<Traits_,Ss_> >
class Polygon_offset_builder_2
{
public :
typedef Ss_ Ss ;
typedef Traits_ Traits ;
typedef Container_ Container ;
typedef Visitor_ Visitor ;
typedef typename Traits::Kernel K ;
typedef typename Traits::FT FT ;
typedef typename Traits::Point_2 Point_2 ;
typedef std::optional<Point_2> OptionalPoint_2 ;
typedef std::shared_ptr<Container> ContainerPtr ;
Polygon_offset_builder_2( Ss const& aSs, Traits const& aTraits = Traits(), Visitor const& aVisitor = Visitor() ) ;
template<class OutputIterator>
OutputIterator construct_offset_contours( FT aTime, OutputIterator aOut ) ;
struct Bisector_data
{
Bisector_data() : IsVisited(false), IsUsedSeed(false) {}
bool IsVisited ;
bool IsUsedSeed ;
} ;
private:
typedef typename Ss::Vertex Vertex ;
typedef typename Ss::Halfedge_handle Halfedge_handle ;
typedef typename Ss::Halfedge_const_handle Halfedge_const_handle ;
typedef typename Ss::Vertex_const_handle Vertex_const_handle ;
typedef std::vector<Halfedge_const_handle> Halfedge_vector ;
typedef typename Traits::Segment_2 Segment_2 ;
typedef typename Traits::Trisegment_2 Trisegment_2 ;
typedef typename Traits::Trisegment_2_ptr Trisegment_2_ptr ;
typedef CGAL_SS_i::Triedge<Halfedge_handle> Triedge ;
enum Hook_position { SOURCE, TARGET, INSIDE } ;
static char const* Hook_position2Str( Hook_position aPos )
{
return aPos == SOURCE ? "source" : ( aPos == TARGET ? "target" : "inside" ) ;
}
Halfedge_const_handle LocateHook( FT aTime, Halfedge_const_handle aBisector, bool aIncludeLastBisector, Hook_position& rPos ) ;
void AddOffsetVertex( FT aTime, Halfedge_const_handle aHook, ContainerPtr aPoly ) ;
template<class OutputIterator>
OutputIterator TraceOffsetPolygon( FT aTime, Halfedge_const_handle aHook, OutputIterator aOut ) ;
Halfedge_const_handle LocateSeed( FT aTime, Halfedge_const_handle aBorder ) ;
Halfedge_const_handle LocateSeed( FT aTime ) ;
Bisector_data const& GetBisectorData ( Halfedge_const_handle aBisector ) const { return mBisectorData[aBisector->id()] ; }
Bisector_data& GetBisectorData ( Halfedge_const_handle aBisector ) { return mBisectorData[aBisector->id()] ; }
bool IsVisited( Halfedge_const_handle aBisector ) { return GetBisectorData(aBisector).IsVisited ; }
bool IsUsedSeed( Halfedge_const_handle aBisector ) { return GetBisectorData(aBisector).IsUsedSeed ; }
void Visit( Halfedge_const_handle aBisector ) { GetBisectorData(aBisector).IsVisited = true ; }
void SetIsUsedSeed( Halfedge_const_handle aBisector ) { GetBisectorData(aBisector).IsUsedSeed = true ; }
inline Segment_2 CreateSegment ( Halfedge_const_handle aH ) const
{
Point_2 s = aH->opposite()->vertex()->point() ;
Point_2 t = aH->vertex()->point() ;
return K().construct_segment_2_object()(s,t);
}
Trisegment_2_ptr GetTrisegment ( Vertex_const_handle aNode ) const ;
public:
Comparison_result Compare_offset_against_event_time( FT aT, Vertex_const_handle aNode ) const
{
CGAL_precondition( aNode->is_skeleton() ) ;
Comparison_result r = aNode->has_infinite_time() ? SMALLER
: static_cast<Comparison_result>(mTraits.compare_offset_against_event_time_2_object()(aT,GetTrisegment(aNode)));
return r ;
}
std::optional<Point_2> Construct_offset_point( FT aT, Halfedge_const_handle aBisector ) const
{
CGAL_assertion(aBisector->is_bisector());
CGAL_assertion(handle_assigned(aBisector->opposite()));
CGAL_assertion(aBisector->opposite()->is_bisector());
Halfedge_const_handle lBorderA = aBisector->defining_contour_edge();
Halfedge_const_handle lBorderB = aBisector->opposite()->defining_contour_edge();
Vertex_const_handle lNodeS = aBisector->opposite()->vertex();
Vertex_const_handle lNodeT = aBisector->vertex();
// If aBisector is not a border bisector the offset point construction needs to get to seed event
Trisegment_2_ptr lSeedEvent ;
if ( aBisector->is_inner_bisector() )
{
CGAL_assertion ( lNodeS->is_skeleton() ) ;
CGAL_assertion ( lNodeT->is_skeleton() ) ;
Vertex_const_handle lSeedNode = aBisector->slope() == POSITIVE ? lNodeS : lNodeT ;
lSeedEvent = GetTrisegment(lSeedNode) ;
CGAL_POLYOFFSET_TRACE(3,"Seed node for " << e2str(*aBisector) << " is " << v2str(*lSeedNode) << " event=" << lSeedEvent ) ;
}
OptionalPoint_2 p = mTraits.construct_offset_point_2_object()( aT,
CreateSegment(lBorderA),
lBorderA->weight(),
CreateSegment(lBorderB),
lBorderB->weight(),
lSeedEvent );
return p ;
}
private:
void ResetBisectorData();
Traits const& mTraits ;
Visitor const& mVisitor ;
Halfedge_vector mBorders ;
std::vector<Bisector_data> mBisectorData;
OptionalPoint_2 mLastPoint ;
CGAL_POLYOFFSET_DEBUG_CODE( int mStepID ; )
};
} // end namespace CGAL
#include <CGAL/Straight_skeleton_2/Polygon_offset_builder_2_impl.h>
#endif // CGAL_POLYGON_OFFSET_BUILDER_2_H //
// EOF //