G+Smo  24.08.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gsTrimSurface< T > Class Template Reference

Detailed Description

template<class T>
class gismo::gsTrimSurface< T >

Class for a trim surface.

Public Member Functions

gsCurveLoop< T > & boundaryLoop ()
 Get a boundary curve loop.
 
const gsCurveLoop< T > & boundaryLoop () const
 Look at (non-const) boundaryLoop.
 
gsTrimSurface< T > * clone () const
 Clone function. Used to make a copy of the (derived) geometry.
 
gsVector3d< T > cornerNormal (int const &sourceID) const
 Compute the surface normal at a corner.
 
void cuttingAngles (int const &sourceID, int const &targetID, T *angle, T *angle1, T *angle2, bool const &isCCWviewFromNormal=true) const
 Compute angle discrepancy between the angles defined by the edge source-target with the sides of the angle wrt source
 
gsBSpline< T > cuttingCurve (int const &sourceID, int const &targetID) const
 Define a spline curve connecting source-target
 
gsMatrix< T > derivatives (int sourceID) const
 
gsMatrix< T > evalCurve (int loopNumber, int curveNumber, const gsMatrix< T > &u) const
 Look at evalCurve_into.
 
void evalCurve_into (int loopNumber, int curveNumber, const gsMatrix< T > &u, gsMatrix< T > &result) const
 
gsMatrix< T > evalSurface (const gsMatrix< T > &u) const
 
void evalSurface_into (const gsMatrix< T > &u, gsMatrix< T > &result) const
 
gsCurve< T > & getCurve (int loopNumber, int curveNumber) const
 Returns curveNumber-th curve in loopNumber-th loop.
 
getLengthOfCurve (const int loopNumber, const int curveNumber, const T eps=1e-6, const int nmbSegments=50) const
 Computes length of curveNumber-th curve in loopNumber-th loop.
 
void getPhysicalyUniformCurveParameters (const int loopNumber, const int curveNumber, const int nmbParams, gsVector< T > &parameters, const T eps=1e-6)
 
 gsTrimSurface ()
 Default empty constructor.
 
 gsTrimSurface (typename gsSurface< T >::Ptr g, gsPlanarDomain< T > *d)
 Constructor by a shared pointer.
 
 gsTrimSurface (gsMatrix< T > const &corner, int patchDeg1, int patchDeg2, int curveDeg)
 
nearestPoint (int loopNumber, int curveNumber, int nTrialPoints, int nIterations, const gsMatrix< T > &spacePoint)
 find the parameter of the nearest point on a curve to a given point in space
 
std::ostream & print (std::ostream &os) const
 Prints the object as a string.
 
gsMatrix< T > sampleNormal (int loopNumber, int curveNumber, size_t npoints) const
 sample standard unit normals along a trimming curve
 
gsMatrix< T > splitCurve (size_t loopId, size_t curveId, T lengthRatio=.5)
 
gsVector< T > TangentCoefs_bisect (int const &sourceID) const
 Compute the (scale of) tangent of the curve emanateing from a vertex in the parameter whose image bisects the corresponding angle.
 
gsVector< T > TangentCoefs_bisect (int const &sourceID, gsVector3d< T > normal) const
 
gsMatrix< T > TangentCoefs_next (int const &sourceID) const
 Return the coefficients of the representation of the unit tangent of the edge ENAMATING from vertex sourceID in terms of the standard basis of the tangent space.
 
gsMatrix< T > TangentCoefs_prev (int const &sourceID) const
 Return the coefficients of the representation of the unit tangent of the edge COMING to vertex sourceID in terms of the standard basis of the tangent space.
 
memory::unique_ptr< gsMesh< T > > toMesh (int npoints=50) const
 Return a triangulation of the trimmed surface.
 
gsMatrix< T > trimCurTangents (int loopN, int curveN, size_t npoints) const
 Return the tangent vectors of the trimming curve curveNumber in trimming loop loopNumber.
 
gsMatrix< T > unitNormal (gsMatrix< T > point) const
 Compute the unit normal vector of the trimmed surface at a point in the parameter domain.
 
gsMatrix< T > UnitTangentCoefs_next (int const &sourceID, gsMatrix< T > const &corJacobian) const
 Return the coefficients of the representation of the unit tangent of the edge ENAMATING from vertex sourceID in terms of the standard basis of the tangent space.
 
gsMatrix< T > UnitTangentCoefs_prev (int const &sourceID, gsMatrix< T > const &corJacobian) const
 Return the coefficients of the representation of the unit tangent of the edge COMING to vertex sourceID in terms of the standard basis of the tangent space.
 

Private Member Functions

arcLength (const gsCurve< T > &curve, const T a, const T b, const int quadPoints=4) const
 
void evalCurve_into (const gsCurve< T > &curve, const gsMatrix< T > &u, gsMatrix< T > &result) const
 
findParameter (const gsCurve< T > &curve, const T arc, T curArc, T lowParam, T uppParam, const T eps)
 
void fromArcsToParams (const int loopNumber, const int curveNumber, const gsVector< T > &arcs, gsVector< T > &parameters, const T eps)
 
getLengthOfCurve (const gsCurve< T > &curve, gsMatrix< T > &params, bool linear=false) const
 

Constructor & Destructor Documentation

gsTrimSurface ( gsMatrix< T > const &  corner,
int  patchDeg1,
int  patchDeg2,
int  curveDeg 
)

Construct a trivial trimmed surface from 4 vertices Inputs: corner (the 4 corner vertices), patchDeg1 and patchDeg2 (spline degree in dim 1 and 2), curveDeg (spline degree of all trimming curves)

Member Function Documentation

T arcLength ( const gsCurve< T > &  curve,
const T  a,
const T  b,
const int  quadPoints = 4 
) const
private

Computes the arc length between surface( curve(a) ) and surface( curve(b) ).

Integral is evaluated with the number of quadPoints.

gsMatrix< T > derivatives ( int  sourceID) const

Compute the partial derivatives of the surface parametrization at a corner point of a trimmed patch specified by the index of the trimming curve sourceID emanating from the corner point

void evalCurve_into ( int  loopNumber,
int  curveNumber,
const gsMatrix< T > &  u,
gsMatrix< T > &  result 
) const
inline

Evaluates curveNumber-th curve from loopNumber-th loop.

Curve is evaluated at parameter u. Parameter u is in curve domain and resulting points are in physical space ( result = surf(curve(u)) ).

void evalCurve_into ( const gsCurve< T > &  curve,
const gsMatrix< T > &  u,
gsMatrix< T > &  result 
) const
private

Evaluates curve (in physical space) at parameter u.

result = surface( curve(u) )

gsMatrix<T> evalSurface ( const gsMatrix< T > &  u) const
inline

Look at evalSurface_into

Parameters
uparameters
Returns
points on surface at parameter u
void evalSurface_into ( const gsMatrix< T > &  u,
gsMatrix< T > &  result 
) const
inline

Evaluates surface.

Note that this is NOT trimmed surface. (You can evaluate at parameter which is trimmed away.)

Parameters
uparameters
resultpoints on surface at parameter u
T findParameter ( const gsCurve< T > &  curve,
const T  arc,
curArc,
lowParam,
uppParam,
const T  eps 
)
private

Let be C( t ) = surface( curve( t ) ). Function finds a parameter t0 of the curve C, where C( t0 ) lies arc away from the C( 0 ).

Let P1 lies curArc away from C( 0 ) on the curve C. Let P2 lies less than arc away from C( 0 ) on the curve.

uppParam corresponds to P1 lowParam corresponts to P2

void fromArcsToParams ( const int  loopNumber,
const int  curveNumber,
const gsVector< T > &  arcs,
gsVector< T > &  parameters,
const T  eps 
)
private

Let curve be a curveNumber-th curve in loopNumber-th loop. Let be C( t ) = surface( curve( t ) ), t in [parameter domain] Function returns parameters t_i such that C( t_i ) = P_i and L( C( 0 ), C( t_i ) ) = arcs( i ) where L( x, y ) is an arc length of a curve C between x and y.

We return t_i in a vector parameters.

T getLengthOfCurve ( const gsCurve< T > &  curve,
gsMatrix< T > &  params,
bool  linear = false 
) const
private

Computes the length of curve: surface( curve(t) ).

Parameters
curve
paramsprecomputed parameters (we evaluate length of the curve from the begining till the end of params values)
lineartrue - we do linear approximation false - we compute integral of derivative surface( curve(t) )'
void getPhysicalyUniformCurveParameters ( const int  loopNumber,
const int  curveNumber,
const int  nmbParams,
gsVector< T > &  parameters,
const T  eps = 1e-6 
)
inline

Returns parameters t_i, i = 1, 2, ..., nmbParams, such that arc between surface( curve (t_i) ) and surface( curve (t_{i + 1}) ) is equal for all i.

With other words: Parameters are maped into points that are uniform in physical space.

curve is the curveNumber-th curve in loopNumber-th loop.

gsMatrix<T> splitCurve ( size_t  loopId,
size_t  curveId,
lengthRatio = .5 
)
inline

split the curveId^th curve in the loopId^th loop of the planar domain into two curves

Parameters
loopIdspecifies the loop
curveIdspecifies the curve in the loop
lengthRatiothe ratio of the lengths of the first new curve and of the original curve
gsVector< T > TangentCoefs_bisect ( int const &  sourceID,
gsVector3d< T >  normal 
) const

Compute the (scale of) tangent of the curve emanateing from a vertex in the parameter whose image bisects the corresponding angle In addition, accounting for the bisecting vector defined according to the normal