G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsCurveLoop< T > Class Template Reference

Detailed Description

template<class T>
class gismo::gsCurveLoop< T >

A closed loop given by a collection of curves.

The loop may be oriented clockwise (CW) or counterclockwise (CCW).

Template Parameters
Tcoefficient type

Public Types

typedef memory::shared_ptr< gsCurveLoopPtr
 Shared pointer for gsCurveLoop.
 
typedef memory::unique_ptr< gsCurveLoopuPtr
 Unique pointer for gsCurveLoop.
 

Public Member Functions

 gsCurveLoop ()
 Default empty constructor.
 
 gsCurveLoop (const std::vector< gsVector3d< T > * > points3D, const std::vector< bool > isConvex, T margin, gsVector3d< T > *outNormal)
 
 gsCurveLoop (const std::vector< T > angles3D, const std::vector< T > lengths3D, const std::vector< bool > isConvex, T margin, const bool unitSquare=false)
 
 gsCurveLoop (gsCurve< T > *cc)
 Constructor by one curve.
 
 gsCurveLoop (std::vector< gsCurve< T > * > const &curves)
 Constructor by a list of curves forming a loop.
 
isOn

returns true if the points in

Parameters
uare ON the curve ideally it should store in
paramResultthe parametric values associated to those points.
is_ccw

returns true if the curve is oriented counterclockwise

reverse

changes the orientation of the curve

gsCurve< T >::uPtr singleCurve () const
 Return the curve-loop as a single new B-Spline curve.
 
std::ostream & print (std::ostream &os) const
 Prints the object as a string.
 
gsMatrix< T > normal (int const &c, gsMatrix< T > const &u)
 
uPtr split (int startIndex, int endIndex, gsCurve< T > *newCurveThisFace, gsCurve< T > *newCurveNewFace)
 
std::vector< T > lineIntersections (int const &direction, T const &abscissa)
 
int numCurves () const
 get Number of curves
 
gsMatrix< T > sample (int npoints=50, int numEndPoints=2) const
 Sample npoints uniformly distributed (in parameter domain) points on the loop.
 
std::vector< T > domainSizes () const
 return a vector containing whose nth entry is the length of the domin of the nth curve.
 
gsMatrix< T > splitCurve (size_t curveId, T lengthRatio=.5)
 
gsVector3d< T > initFrom3DPlaneFit (const std::vector< gsVector3d< T > * > points3D, T margin)
 Initialize a curve loop from some 3d vertices by projecting them onto the plane of best fit and constructing a polygon.
 
bool initFrom3DByAngles (const std::vector< gsVector3d< T > * > points3D, const std::vector< bool > isConvex, T margin)
 
void initFromIsConvex (const std::vector< bool > isConvex, T margin)
 Initialize a curve loop from a list of bools indicating whether the corners are convex or not.
 
void flip1 (T minu=0, T maxu=1)
 flip a curve loop in the u direction
 
static bool approximatingPolygon (const std::vector< T > &signedAngles, const std::vector< T > &lengths, T margin, gsMatrix< T > &result)
 
static void adjustPolygonToUnitSquare (gsMatrix< T > &corners, T const margin)
 

Constructor & Destructor Documentation

◆ gsCurveLoop() [1/2]

template<class T >
gsCurveLoop ( const std::vector< gsVector3d< T > * >  points3D,
const std::vector< bool >  isConvex,
margin,
gsVector3d< T > *  outNormal 
)

Constructor from a sequence of points in 3D. Produces a polygon in 2D whose angles are in proportion to the corresponding angles in space. You need to pass in a sequence of bools indicating which angles are convex.

◆ gsCurveLoop() [2/2]

template<class T >
gsCurveLoop ( const std::vector< T >  angles3D,
const std::vector< T >  lengths3D,
const std::vector< bool >  isConvex,
margin,
const bool  unitSquare = false 
)

Constructor from a sequence of angles and lengths (measured in 3D). If unitSquare==true, when there are four vertices and isConvex are all true, then return the spline curve loop of degree 1 of the unit square with each edge being one spline curve.

Member Function Documentation

◆ adjustPolygonToUnitSquare()

template<class T >
void adjustPolygonToUnitSquare ( gsMatrix< T > &  corners,
T const  margin 
)
static

Scale and shift a polygon (represented by a matrix whose rows are its corners so that it will be contained in the square [0, 1]x[0, 1]. Modifies the polygon in place. Specifying margin > 0 will make the polygon fit in the interior of the square.

◆ approximatingPolygon()

template<class T >
bool approximatingPolygon ( const std::vector< T > &  signedAngles,
const std::vector< T > &  lengths,
margin,
gsMatrix< T > &  result 
)
static

Computes the corners of a polygon that matches the specified angles and relative side lengths as well as possible. Fits the resulting polygon inside the square [0, 1] x [0, 1]. Fills a gsmatrix whose rows are the corners of the polygon. Returns true on success.

◆ initFrom3DByAngles()

template<class T >
bool initFrom3DByAngles ( const std::vector< gsVector3d< T > * >  points3D,
const std::vector< bool >  isConvex,
margin 
)

Initialize a curve loop from some 3d vertices by trying to match the angles as well as possible. Return true on success.

◆ lineIntersections()

template<class T >
std::vector< T > lineIntersections ( int const &  direction,
T const &  abscissa 
)

Computes the intersections with the axis-aligned line with constant value abscissa in direction.

◆ normal()

template<class T >
gsMatrix< T > normal ( int const &  c,
gsMatrix< T > const &  u 
)

Computes the normal vector of curve i at parameter u (is outer in case of CCW orientation)

◆ split()

template<class T >
gsCurveLoop< T >::uPtr split ( int  startIndex,
int  endIndex,
gsCurve< T > *  newCurveThisFace,
gsCurve< T > *  newCurveNewFace 
)

Split a second curve loop off from this one by inserting two new curves. Curve startIndex will be included in the new loop, while curve endIndex will be kept in the original loop.

◆ splitCurve()

template<class T >
gsMatrix< T > splitCurve ( size_t  curveId,
lengthRatio = .5 
)

split the curveId^th curve in the loop into two curves. Return the point where it splits.

Parameters
curveId
lengthRatiothe ratio between the lengths of the first new curve and that of the original curve