TrimAffineTransTrimCurves

(trim_aux.c:1919)

Prototype:

  void TrimAffineTransTrimCurves(TrimCrvStruct *TrimCrvList,
                                 CagdRType OldUMin,
                                 CagdRType OldUMax,
                                 CagdRType OldVMin,
                                 CagdRType OldVMax,
                                 CagdRType NewUMin,
                                 CagdRType NewUMax,
                                 CagdRType NewVMin,
                                 CagdRType NewVMax)


Description:

Map the given trimming curves into a new domain, in place.

Parameters:

TrimCrvList: Trimming curves to affinely map.
OldUMin, OldUMax, OldVMin, OldVMax: Domain to map trimming curves from.
NewUMin, NewUMax, NewVMin, NewVMax: Domain to map trimming curves to.


Returned Value:

void


See Also:

TrimAffineTransTrimSrf

Keywords:

trimming curves


TrimAffineTransTrimSrf

(trim_aux.c:1971)

Prototype:

  TrimSrfStruct *TrimAffineTransTrimSrf(const TrimSrfStruct *CTrimSrf,
                                        CagdRType NewUMin,
                                        CagdRType NewUMax,
                                        CagdRType NewVMin,
                                        CagdRType NewVMax)


Description:

Maps the given trimmed surface into a new domain.

Parameters:

CTrimSrf: Trimmed surface to affinely map its parametric domain.
NewUMin, NewUMax, NewVMin, NewVMax: New parametric domain to map to.


Returned Value:

TrimSrfStruct *: A trimmed surface that is geometrically identical but with new different parametric domain.


See Also:

BspKnotAffineTransOrder2 TrimAffineTransTrimCurves

Keywords:

trimming curves


TrimAllPrisaSrfs

(tr_prisa.c:59)

Prototype:

  TrimSrfStruct *TrimAllPrisaSrfs(const TrimSrfStruct *TSrfs,
                                  int SamplesPerCurve,
                                  CagdRType Epsilon,
                                  CagdSrfDirType Dir,
                                  CagdVType Space)


Description:

Computes a piecewise ruled surface approximation to a given set of trimmed surfaces with given Epsilon, and lay them out "nicely" onto the XY plane, by approximating each ruled surface as a developable surface with SamplesPerCurve samples. Dir controls the direction of ruled approximation, SpaceScale and Offset controls the placement of the different planar pieces. Prisa is the hebrew word for the process of flattening out a three dimensional surface. I have still to find an english word for it.

Parameters:

TSrfs: To approximate and flatten out.
SamplesPerCurve: uring the approximation of a ruled surface as a developable surface.
Epsilon: Accuracy control for the piecewise ruled surface approximation. if Epsilon us positive, the surfaces are laid down on the plane, otherwise they are return as 3-space ruled surfaces and form a piecewise ruled-surface approximation to Srfs.
Dir: Direction of ruled/developable surface approximation. Either U or V.
Space: A vector in the XY plane to denote the amount of translation from one flattened out surface to the next.


Returned Value:

TrimSrfStruct *: A list of planar trimmed surfaces denoting the layout (prisa) of the given TSrfs to the accuracy requested.


See Also:

TrimPiecewiseRuledSrfApprox TrimPrisaRuledSrf SymbAllPrisaSrfs

Keywords:

layout prisa


TrimClassifyTrimCrvsOrientation

(trim_aux.c:1429)

Prototype:

  void TrimClassifyTrimCrvsOrientation(TrimCrvStruct *TCrvs, CagdRType Tol)


Description:

Go over all given loops and classify the individual trimming curves as "reversed" if indeed the curve is reversed in the loop. Assumes the trimming curves were already linked into loops.

Parameters:

TCrvs: Trimming curve to classify their orientation, in place.
Tol: Tolerance of end points comparisons.


Returned Value:

void


See Also:

TrimLinkTrimmingCurves2Loops TrimLinkTrimmingCurves2Loops1 TrimLinkTrimmingCurves2Loops2 TrimMergeTrimmingCurves2Loops TrimClassifyTrimLoopOrient

Keywords:

trimming curves


TrimClassifyTrimCurveOrient

(trim2ply.c:2069)

Prototype:

  CagdBType TrimClassifyTrimCurveOrient(const CagdCrvStruct *UVCrv)


Description:

Given a closed, piecewise linear trimming curve, returns TRUE if the curve is clockwise, FALSE if counter clockwise. Orientation is determined by computing the signed area of the polygon and examining the sign of the result.

Parameters:

UVCrv: Trimming curve to examine its orientation.


Returned Value:

CagdBType: TRUE if the curve is clockwise, FALSE if counter clockwise.


See Also:

CagdCrvAreaPoly

Keywords:




TrimClassifyTrimLoopOrient

(trim2ply.c:2094)

Prototype:

  CagdBType TrimClassifyTrimLoopOrient(const TrimCrvSegStruct *TSegs)


Description:

Given a closed loop consisting of several, piecewise linear, trimming curves, returns TRUE if loop is clockwise, FALSE if counter clockwise. Orientation is determined by computing the signed area of the polygon and examining the sign of the result.

Parameters:

TSegs: Trimming loop to examine its orientation.


Returned Value:

CagdBType: TRUE if the curve is clockwise, FALSE if counter clockwise.


See Also:

TrimClassifyTrimmingLoops TrimClassifyTrimCurveOrient TrimClassifyTrimCrvsOrientation

Keywords:




TrimClassifyTrimmingLoops

(trim2ply.c:1849)

Prototype:

  int TrimClassifyTrimmingLoops(TrimCrvStruct **TrimLoops)


Description:

Classify the given trimming curve loops into a hierarchy. All curves with even nesting levels are considered outside loops and are oriented clockwise. All curves with odd nesting level are considered islands and are oriented counterclockwise. An island Ci of outside loop Cj will be placed in an "_subTrims" attribute under Ci.

Parameters:

TrimLoops: Input loops, reorganized in place.


Returned Value:

int: TRUE if successful, FALSE otherwise.


See Also:

TrimCrvsHierarchy2Polys TrimClassifyTrimCurveOrient TrimCrvFreeWithSubTrims TrimCrvFreeListWithSubTrims

Keywords:




TrimClipSrfToTrimCrvs

(trim_aux.c:2805)

Prototype:

  TrimSrfStruct *TrimClipSrfToTrimCrvs(TrimSrfStruct *TrimSrf)


Description:

Extract the minimal region of the tensor product surface that contains the domain that is prescribed by the trimming curves. The return surface represents the exact same geometry as the input surface but possible with a small size.

Parameters:

TrimSrf: Input trimmed surface to extract a minimal valid region


Returned Value:

TrimSrfStruct *: Returns a possible smaller surface with similar geometry.


See Also:

TrimSrfTrimCrvSquareDomain TrimSrfTrimCrvAllDomain

Keywords:




TrimCnvrtBsp2BzrSrf

(trim_sub.c:1126)

Prototype:

  TrimSrfStruct *TrimCnvrtBsp2BzrSrf(const TrimSrfStruct *TrimSrf)


Description:

Subdivides a trimmed surface into patches in which all the surfaces are Bezier-like (that is, have no internal knots in either U or V direction). The trimming curves are also obviously subdivided. The result is returned as a list of trimmed Bezier-like surfaces.

Parameters:

TrimSrf: A trimmed Bezier or B-spline surface to subdivide into Bezier-like patches.


Returned Value:

TrimSrfStruct *: The subdivided surfaces, represented as a linked list.


Keywords:

subdivision


TrimCoerceTrimUVCrv2Plane

(trim_aux.c:1062)

Prototype:

  int TrimCoerceTrimUVCrv2Plane(TrimCrvSegStruct *TSeg)


Description:

Make sure the UV curve in the trimming curve segment is planar (Z == 0). If not, the UV is projected to the XY plane, in place.

Parameters:

TSeg: Trimming curve segment to coerce to teh XY plane, in place.


Returned Value:

int: TRUE if curve segment was modified.


See Also:



Keywords:

trimming curves


TrimCrv2Polyline

(trim_aux.c:2094)

Prototype:

  CagdPolylineStruct *TrimCrv2Polyline(const CagdCrvStruct *TrimCrv,
                                       CagdRType TolSamples,
                                       SymbCrvApproxMethodType Method,
                                       CagdBType OptiLin)


Description:

Routine to approx. a single curve as a polyline with TolSamples samples/tolerance. Polyline is always E3 CagdPolylineStruct type. NULL is returned in case of an error, otherwise CagdPolylineStruct.

Parameters:

TrimCrv: To approximate as a polyline.
TolSamples: Tolerance of approximation error (Method = 2) or Number of samples to compute on polyline (Method = 0, 1).
Method: 0 - TolSamples are set uniformly in parametric space, 1 - TolSamples are set optimally, considering the isocurve's curvature. 2 - TolSamples sets the maximum error allowed between the piecewise linear approximation and original curve.
OptiLin: If TRUE, optimize linear curves.


Returned Value:

CagdPolylineStruct *: A polyline representing the piecewise linear approximation from, or NULL in case of an error.


See Also:

BspCrv2Polyline BzrCrv2Polyline IritCurve2Polylines SymbCrv2Polyline TrimCrvs2Polylines

Keywords:

piecewise linear approximation polyline


TrimCrvAgainstTrimCrvs

(trim_iso.c:594)

Prototype:

  CagdCrvStruct *TrimCrvAgainstTrimCrvs(CagdCrvStruct *UVCrv,
                                        const TrimSrfStruct *TrimSrf,
                                        CagdRType Eps)


Description:

Trim a given curve in UV space of a trimmed surface to the valid domain only. Returned is a list of segments of UV curve that is inside TrimSrf.

Parameters:

UVCrv: A curve in the parametric space to trim. Used in place and Freed at the end. Can be a non isoparametric curve.
TrimSrf: A trimmed surface.
Eps: Tolerance of approximation.


Returned Value:

CagdCrvStruct *: A list of trimmed segments of UVCrv inside TrimSrf.


See Also:

TrimIntersectTrimCrvIsoVals

Keywords:




TrimCrvBBox

(trim_aux.c:140)

Prototype:

  int TrimCrvBBox(const TrimCrvStruct *TCrv, int UV, CagdBBoxStruct *BBox)


Description:

Computes a bounding box for a freeform trimmed curve.

Parameters:

TCrv: To compute a bounding box for.
UV: TRUE to process the UV curve, FALSE for the Euclidean curve.
BBox: Where bounding information is to be saved.


Returned Value:

int: TRUE if successful, FALSE otherwise


See Also:

TrimSrfBBox TrimCrvListBBox CagdCrvBBox

Keywords:

bbox bounding box


TrimCrvCopy

(trim_gen.c:295)

Prototype:

  TrimCrvStruct *TrimCrvCopy(const TrimCrvStruct *TrimCrv)


Description:

Duplicates a trimming curve structure.

Parameters:

TrimCrv: A trimming curve to duplicate.


Returned Value:

TrimCrvStruct *: A trimming curve structure.


Keywords:

allocation


TrimCrvCopyList

(trim_gen.c:323)

Prototype:

  TrimCrvStruct *TrimCrvCopyList(const TrimCrvStruct *TrimCrvList)


Description:

Allocates and copies a list of trimming curve structures.

Parameters:

TrimCrvList: To be copied.


Returned Value:

TrimCrvStruct *: A duplicated list of trimming curves.


Keywords:

copy


TrimCrvFree

(trim_gen.c:352)

Prototype:

  void TrimCrvFree(TrimCrvStruct *TrimCrv)


Description:

Deallocates a trimming curve structure.

Parameters:

TrimCrv: trimming curve to free.


Returned Value:

void


Keywords:

allocation


TrimCrvFreeList

(trim_gen.c:372)

Prototype:

  void TrimCrvFreeList(TrimCrvStruct *TrimCrvList)


Description:

Deallocates a list of trimming curve structures.

Parameters:

TrimCrvList: list of trimming curve to free.


Returned Value:

void


Keywords:

allocation


TrimCrvFreeListWithSubTrims

(trim2ply.c:2039)

Prototype:

  void TrimCrvFreeListWithSubTrims(TrimCrvStruct *TrimCrv)


Description:

Frees a list of trimming curve loops, which may or may not contain internal trimming curve loops. Each trimming curve in the list is freed by calling TrimCrvFreeRecursive.

Parameters:

TrimCrv: A list of trimming curves to free.


Returned Value:

void


See Also:

TrimCrvFreeWithSubTrims TrimClassifyTrimmingLoops

Keywords:




TrimCrvFreeWithSubTrims

(trim2ply.c:2007)

Prototype:

  void TrimCrvFreeWithSubTrims(TrimCrvStruct *TrimCrv)


Description:

Frees a trimming curve and all its internal trimming curve loops, recursively (by TrimCrvFreeListRecursive). The list of internal trimming curve loops are assumed to be in the _subTrims attribute.

Parameters:

TrimCrv: A trimming curve to free.


Returned Value:

void


See Also:

TrimCrvFreeListWithSubTrims TrimClassifyTrimmingLoops

Keywords:




TrimCrvListBBox

(trim_aux.c:163)

Prototype:

  int TrimCrvListBBox(const TrimCrvStruct *TCrvs, int UV, CagdBBoxStruct *BBox)


Description:

Computes a bounding box for a list of freeform trimmed curves.

Parameters:

TCrvs: To compute a bounding box for.
UV: TRUE to process the UV curve, FALSE for the Euclidean curve.
BBox: Where bounding information is to be saved.


Returned Value:

int: TRUE if successful, FALSE otherwise


See Also:

TrimCrvBBox CagdCrvBBox GMBBSetBBoxPrecise

Keywords:

bbox bounding box


TrimCrvNew

(trim_gen.c:269)

Prototype:

  TrimCrvStruct *TrimCrvNew(TrimCrvSegStruct *TrimCrvSegList)


Description:

Allocates a trimming curve structure.

Parameters:

TrimCrvSegList: List of trimming curve segments forming the trimming curve.


Returned Value:

TrimCrvStruct *: A trimmig curve.


Keywords:

allocation


TrimCrvSegBBox

(trim_aux.c:69)

Prototype:

  int TrimCrvSegBBox(const TrimCrvSegStruct *TCrvSeg,
                     int UV,
                     CagdBBoxStruct *BBox)


Description:

Computes a bounding box for a freeform trimmed curve.

Parameters:

TCrvSeg: To compute a bounding box for.
UV: TRUE to process the UV curve, FALSE for the Euclidean curve.
BBox: Where bounding information is to be saved.


Returned Value:

int: TRUE if successful, FALSE otherwise


See Also:

TrimSrfBBox TrimCrvListBBox CagdCrvBBox TrimCrvBBox TrimCrvListBBox

Keywords:

bbox bounding box


TrimCrvSegCopy

(trim_gen.c:164)

Prototype:

  TrimCrvSegStruct *TrimCrvSegCopy(const TrimCrvSegStruct *TrimCrvSeg)


Description:

Duplicates a trimming curve segment structure.

Parameters:

TrimCrvSeg: trimming curve segment to duplicate.


Returned Value:

TrimCrvSegStruct *: A trimming curve segment structure.


Keywords:

allocation


TrimCrvSegCopyList

(trim_gen.c:194)

Prototype:

  TrimCrvSegStruct *TrimCrvSegCopyList(const TrimCrvSegStruct *TrimCrvSegList)


Description:

Allocates and copies a list of trimming curve segment structures.

Parameters:

TrimCrvSegList: To be copied.


Returned Value:

TrimCrvSegStruct *: A duplicated list of trimming curve segments.


Keywords:

copy


TrimCrvSegFree

(trim_gen.c:223)

Prototype:

  void TrimCrvSegFree(TrimCrvSegStruct *TrimCrvSeg)


Description:

Deallocates a trimming curve segment structure.

Parameters:

TrimCrvSeg: trimming curve segment to free.


Returned Value:

void


Keywords:

allocation


TrimCrvSegFreeList

(trim_gen.c:244)

Prototype:

  void TrimCrvSegFreeList(TrimCrvSegStruct *TrimCrvSegList)


Description:

Deallocates a list of trimming curve segment structures.

Parameters:

TrimCrvSegList: list of trimming curve segments to free.


Returned Value:

void


Keywords:

allocation


TrimCrvSegListBBox

(trim_aux.c:106)

Prototype:

  int TrimCrvSegListBBox(const TrimCrvSegStruct *TCrvSegs,
                         int UV,
                         CagdBBoxStruct *BBox)


Description:

Computes a bounding box for a list of freeform trimmed curves.

Parameters:

TCrvSegs: To compute a bounding box for.
UV: TRUE to process the UV curve, FALSE for the Euclidean curve.
BBox: Where bounding information is to be saved.


Returned Value:

int: TRUE if successful, FALSE otherwise


See Also:

TrimCrvSegBBox CagdCrvBBox GMBBSetBBoxPrecise TrimCrvBBox TrimCrvListBBox

Keywords:

bbox bounding box


TrimCrvSegListReverse

(trim2ply.c:1690)

Prototype:

  TrimCrvSegStruct *TrimCrvSegListReverse(TrimCrvSegStruct *TSegs)


Description:

Reverse all the curves in the given list of trimming curve segments, in place.

Parameters:

TSegs: Crv segments to reverse, in place.


Returned Value:

TrimCrvSegStruct *: Reversed list.


See Also:

TrimCrvSegReverse

Keywords:




TrimCrvSegNew

(trim_gen.c:55)

Prototype:

  TrimCrvSegStruct *TrimCrvSegNew(CagdCrvStruct *UVCrv, CagdCrvStruct *EucCrv)


Description:

Allocates a trimming curve segment structure. Allows periodic and float end conditions - converts them to open end. Input curves are used in place.

Parameters:

UVCrv: A UV curve. Only the E2/P2 potion of the curve is considered.
EucCrv: Optional Euclidean curve. Must be an E3 curve.


Returned Value:

TrimCrvSegStruct *: A trimming curve segment structure.


Keywords:

allocation


TrimCrvSegNewList

(trim_gen.c:125)

Prototype:

  TrimCrvSegStruct *TrimCrvSegNewList(CagdCrvStruct *UVCrvs,
                                      CagdCrvStruct *EucCrvs)


Description:

Allocates a list of trimming curve segment structure. Input curves are used in place.

Parameters:

UVCrvs: A list of UV curves. Only the E2/P2 potion of the curve is considered.
EucCrvs: Optional list of Euclidean curves. Must be an E3 curve. If provided, list length must be the same as UVCrvs length.


Returned Value:

TrimCrvSegStruct *: A list of trimming curve segment structures.


Keywords:

allocation


TrimCrvSegReverse

(trim2ply.c:1658)

Prototype:

  void TrimCrvSegReverse(TrimCrvSegStruct *TSeg)


Description:

Reverse the given trimming curve segment, in place.

Parameters:

TSeg: Crv segment to reverse (its UV, and Euclidean curve, if exists).


Returned Value:

void


See Also:

TrimCrvSegListReverse

Keywords:




TrimCrvTrimParamList

(trim_iso.c:234)

Prototype:

  CagdCrvStruct *TrimCrvTrimParamList(CagdCrvStruct *Crv,
                                      TrimIsoInterStruct *InterList)


Description:

Trim Crv at the domains prescribed in the intersection list InterList Both Crv and InterList are FREED in this routine.

Parameters:

Crv: To trim out according to the prescribed intersections.
InterList: List of intersections, as parameters into Crv.


Returned Value:

CagdCrvStruct *: List of trimmed curves. May be empty (NULL).


Keywords:

curves isoparametric curves


TrimCrvs2Polylines

(trim_aux.c:2036)

Prototype:

  CagdPolylineStruct *TrimCrvs2Polylines(TrimSrfStruct *TrimSrf,
                                         CagdBType ParamSpace,
                                         CagdRType TolSamples,
                                         SymbCrvApproxMethodType Method)


Description:

Routine to convert the trimming curves of a trimmed surface to polylines. Polyline are always E3 of CagdPolylineStruct type. NULL is returned in case of an error, otherwise list of CagdPolylineStruct.

Parameters:

TrimSrf: To extract isoparametric curves from.
ParamSpace: TRUE for curves in parametric space, FALSE of 3D Euclidean space.
TolSamples: Tolerance of approximation error (Method = 2) or Number of samples to compute on polyline (Method = 0, 1).
Method: 0 - TolSamples are set uniformly in parametric space, 1 - TolSamples are set optimally, considering the isocurve's curvature. 2 - TolSamples sets the maximum error allowed between the piecewise linear approximation and original curve.


Returned Value:

CagdPolylineStruct *: List of polylines representing a piecewise linear approximation of the extracted isoparamteric curves or NULL is case of an error.


See Also:

TrimCrv2Polyline TrimEvalTrimCrvToEuclid

Keywords:

trimming curves


TrimCrvsHierarchy2Polys

(trim2ply.c:911)

Prototype:

  IPPolygonStruct *TrimCrvsHierarchy2Polys(TrimCrvStruct *TrimLoops)


Description:

Converts a hierarchy of trimming curves/loops into closed, simple, polygons. The input trimming curves are destroyed by this function. A trimming curve inside another trimming curve are chained into one by adding a bidirection line segment between the two curves. A trimming loop will have its contained trimming loops in an attribute "_subTrims".

Parameters:

TrimLoops: A linked list of trimming loops hierarchy (a 'forrest').


Returned Value:

IPPolygonStruct *: Piecewise linear polylines approximating the given trimming curves hierarchy.


See Also:

TrimClassifyTrimmingLoops

Keywords:

polygonization surface approximation


TrimDbg

(trim_dbg.c:32)

Prototype:

  void TrimDbg(const void *Obj)


Description:

Prints trimmed surfaces to stderr. Should be linked to programs for debugging purposes, so trimmed surfaces may be inspected from a debugger.

Parameters:

Obj: A trimmed surface - to be printed to stderr.


Returned Value:

void


Keywords:

debugging


TrimDbgTCrvSegs

(trim_dbg.c:103)

Prototype:

  void TrimDbgTCrvSegs(const TrimCrvSegStruct *TrimSegs)


Description:

Prints the trimming curves of a given trimmed surface..

Parameters:

TrimSegs: Trimming curve segments to print.


Returned Value:

void


See Also:

TrimDbg

Keywords:

debugging


TrimDbgTCrvs

(trim_dbg.c:63)

Prototype:

  void TrimDbgTCrvs(const TrimCrvStruct *TrimCrv)


Description:

Prints the trimming curves of a given trimmed surface..

Parameters:

TrimCrv: Trimming curves to print.


Returned Value:

void


See Also:

TrimDbg

Keywords:

debugging


TrimDbgVerifyTSeg

(trim_dbg.c:131)

Prototype:

  int TrimDbgVerifyTSeg(const TrimCrvSegStruct *TrimSeg)


Description:

Verify the correctness of a trimming curve segment.

Parameters:

TrimSeg: Trimming curve segment to verify.


Returned Value:

int: TRUE if verified, FALSE otherwise.


See Also:

TrimDbgVerifyUVCrv

Keywords:




TrimDbgVerifyUVCrv

(trim_dbg.c:157)

Prototype:

  int TrimDbgVerifyUVCrv(const CagdCrvStruct *TrimCrv)


Description:

Verify the correctness of a trimming curve.

Parameters:

TrimCrv: Trimming curve to verify.


Returned Value:

int: TRUE if verified, FALSE otherwise.


See Also:

TrimDbgVerifyUVCrv

Keywords:




TrimDescribeError

(trim_err.c:59)

Prototype:

  const char *TrimDescribeError(TrimFatalErrorType ErrorNum)


Description:

Returns a string describing a the given error. Errors can be raised by any member of this trim library as well as other users. Raised error will cause an invokation of TrimFatalError function which decides how to handle this error. TrimFatalError can for example, invoke this routine with the error type, print the appropriate message and quit the program.

Parameters:

ErrorNum: Type of the error that was raised.


Returned Value:

const char *: A string describing the error type.


Keywords:

error handling


TrimEnsureNoSingleTrimCrvLoops

(trim2ply.c:1763)

Prototype:

  int TrimEnsureNoSingleTrimCrvLoops(TrimCrvStruct **TrimLoops)
  


Description:

Divides closed trimming loops consisting of a single curve, into two. In other words, after this functions, all trimming loops will have at least two trimming curves.

Parameters:

TrimLoops: Input loops, to possibly subdivided, in place.


Returned Value:

int: TRUE if successful, FALSE otherwise.


Keywords:




TrimEvalTrimCrvToEuclid

(trim_aux.c:2189)

Prototype:

  CagdCrvStruct *TrimEvalTrimCrvToEuclid(const CagdSrfStruct *Srf,
                                         const CagdCrvStruct *UVCrv)


Description:

Computes the composed Euclidean curve of Srf(UVCrv). The resulting curve is either computed using a piecewise linear approximation or by symbolically composing it onto the surface. See TrimSetEuclidComposedFromUV for a way to control this computation.

Parameters:

Srf: To compute the Euclidean UVCrv for.
UVCrv: A curve in the parametric space of Srf.


Returned Value:

CagdCrvStruct *: A Euclidean curve in Srf, following UVCrv.


See Also:

TrimCrvs2Polylines TrimSetEuclidComposedFromUV TrimEvalTrimCrvToEuclid2

Keywords:




TrimEvalTrimCrvToEuclid2

(trim_aux.c:2218)

Prototype:

  CagdCrvStruct *TrimEvalTrimCrvToEuclid2(const CagdSrfStruct *Srf,
                                          const CagdCrvStruct *UVCrv,
                                          CagdCrvStruct **UVCrvLinear)


Description:

Computes the composed Euclidean curve of TrimSrf(UVCrv). The resulting curve is either computed using a piecewise linear approximation or by symbolically composing it onto the surface. See TrimSetEuclidComposedFromUV for a way to control this computation.

Parameters:

Srf: To compute the Euclidean UVCrv for.
UVCrv: A curve in the parametric space of TrimSrf.
UVCrvLinear: f not NULL, updated with the piecewise linear UV space curve, of the same length as the returned Euclidean curve.


Returned Value:

CagdCrvStruct *: A Euclidean curve in TrimSrf, following UVCrv.


See Also:

TrimCrvs2Polylines TrimSetEuclidComposedFromUV TrimEvalTrimCrvToEuclid

Keywords:




TrimExtendTrimmingDomain

(trim_gen.c:1893)

Prototype:

  CagdCrvStruct *TrimExtendTrimmingDomain(const CagdSrfStruct *Srf,
                                          CagdCrvStruct *TrimCrvs,
                                          CagdRType Extnt,
                                          CagdBType MergeCrvs)


Description:

Given a trimming domain, extend it in all UV directions by Extnt amount interior to domain curves are not affected. Trimming curves are assuming not-merged. That is, the UMin and VMin boundary trimming curves are not merged into one curve, etc.

Parameters:

Srf: Surfaces that owns these trimming curves.
TrimCrvs: Trimming curves to extend, in place.
Extnt: Extension amount.
MergeCrvs: If TRUE, aims to merge the given curves into loops, etc.


Returned Value:

CagdCrvStruct *: Extended trimming curves of extended domain.


See Also:



Keywords:




TrimFatalError

(trim_ftl.c:56)

Prototype:

  void TrimFatalError(TrimFatalErrorType ErrID)


Description:

Trap Trim_lib errors right here. Provides a default error handler for the trim library. Gets an error description using TrimDescribeError, prints it and exit the program using exit.

Parameters:

ErrID: Error type that was raised.


Returned Value:

void


Keywords:

error handling


TrimFindClosestTrimCurve2UV

(trim_gen.c:1746)

Prototype:

  CagdRType TrimFindClosestTrimCurve2UV(TrimCrvStruct *TCrvs,
                                        const CagdUVType UV,
                                        TrimCrvSegStruct **ClosestTSeg)


Description:

Find the closest location on the trimming curves TCrvs, to the given UV locations, in the parametric space.

Parameters:

TCrvs: Trimming curves to search for the closest location.
UV: The UV location to find a closest location on TCrvs.
ClosestTSeg: Will be updated with the closest TCrv in TCrvs.


Returned Value:

CagdRType: Parameter value of the closest curve, ClosestTCrv.


See Also:



Keywords:




TrimGetFullDomainTrimCrv

(trim_gen.c:1696)

Prototype:

  const TrimCrvSegStruct *TrimGetFullDomainTrimCrv(const TrimSrfStruct *TSrf)


Description:

Identifies the outer loop of the trimmed surfaces and returns a reference to it if it covers all the tensor product surface domain. A NULL is returned if outer loops is not entire tensor product domain.

Parameters:

TSrf: Trimmed surface to look for its outer full loop.


Returned Value:

const TrimCrvSegStruct *: Return outer full loop if exists.


See Also:

TrimGetLargestTrimmedSrf TrimGetOuterTrimCrv

Keywords:




TrimGetLargestTrimmedSrf

(trim_gen.c:1596)

Prototype:

  TrimSrfStruct *TrimGetLargestTrimmedSrf(TrimSrfStruct **TSrfs, int Extract)


Description:

Find the largest trimmed surface in the given list. While there notion of 'largest' is not really well define, we simply use here the heuristics of largest == trimmed surface with longest trimming curves.

Parameters:

TSrfs: List of trimmed usrface to find the 'largest'.
Extract: TRUE to extract and return the 'largest' trimmed surface from TSrfs, FALSE to only find it.


Returned Value:

TrimSrfStruct *: List of trimmed surface to search for 'largest'.


See Also:

TrimGetOuterTrimCrv TrimGetFullDomainTrimCrv

Keywords:




TrimGetOuterTrimCrv

(trim_gen.c:1649)

Prototype:

  const TrimCrvSegStruct *TrimGetOuterTrimCrv(const TrimSrfStruct *TSrf)


Description:

Returns the outer trimming curve found in TSrf. If their is more than one outer curve, one of then is returned.

Parameters:

TSrf: Trimmed surface to look for its outer loop.


Returned Value:

const TrimCrvSegStruct *: Outer loop.


See Also:

TrimGetLargestTrimmedSrf TrimGetFullDomainTrimCrv

Keywords:




TrimGetTrimCrvLinearApprox

(trim_iso.c:910)

Prototype:

  CagdRType TrimGetTrimCrvLinearApprox(void)


Description:

Get the current tolerance used when approximating higher order trimming curves using piecewise linear approximation.

Parameters:

None


Returned Value:

CagdRType: Sought tolerance.


See Also:

TrimSetTrimCrvLinearApprox

Keywords:




TrimGetTrimmingCurves

(trim_aux.c:787)

Prototype:

  CagdCrvStruct *TrimGetTrimmingCurves(const TrimSrfStruct *TrimSrf,
                                       CagdBType ParamSpace,
                                       CagdBType EvalEuclid)


Description:

Extracts the trimming curves of the given trimmed surface.

Parameters:

TrimSrf: Trimmed surface to extract trimming curves from.
ParamSpace: TRUE for curves in parametric space, FALSE for 3D Euclidean space.
EvalEuclid: If TRUE and ParamSpace is FALSE, evaluate Euclidean curve even if one exists.


Returned Value:

CagdCrvStruct *: List of trimming curves of TrimSrf.


See Also:

TrimManageTrimmingCurvesDegrees TrimGetTrimmingCurves2

Keywords:

trimming curves


TrimGetTrimmingCurves2

(trim_aux.c:819)

Prototype:

  CagdCrvStruct *TrimGetTrimmingCurves2(const TrimCrvStruct *TrimCrvList,
                                        const TrimSrfStruct *TrimSrf,
                                        CagdBType ParamSpace,
                                        CagdBType EvalEuclid)


Description:

Extracts the trimming curves as curves from the given trimming curves.

Parameters:

TrimCrvList: Trimming curves to extract trimming curves as curves.
TrimSrf: Trimmed surface to extract trimming curves from. This parameter is optional and used only if EvalEuclid and/or !ParamSpace.
ParamSpace: TRUE for curves in parametric space, FALSE of 3D Euclidean space.
EvalEuclid: If TRUE and ParamSpace is FALSE, evaluate Euclidean curve even if one exists.


Returned Value:

CagdCrvStruct *: List of trimming curves as curves.


See Also:

TrimManageTrimmingCurvesDegrees TrimGetTrimmingCurves

Keywords:

trimming curves


TrimIntersectCrvsIsoVals

(trim_iso.c:506)

Prototype:

  TrimIsoInterStruct **TrimIntersectCrvsIsoVals(const CagdCrvStruct *UVCrvs,
                                                int Dir,
                                                CagdRType *IsoParams,
                                                int NumOfIsocurves)


Description:

Computes the intersections of given UV curves with the ordered isoparametric values prescribed by IsoParams, in axis Axis.

Parameters:

UVCrvs: UV curves to intersect. Must be piecewise linear.
Dir: Either U or V.
IsoParams: Vector of isoparametric values in direction Dir.
NumOfIsocurves: ize of vector IsoParams.


Returned Value:

TrimIsoInterStruct **: A vector of size NumOfIsocurves, each slot contains a list of intersection parameter values.


See Also:

TrimIntersectTrimCrvIsoVals

Keywords:

isoparametric curves


TrimIntersectTrimCrvIsoVals

(trim_iso.c:303)

Prototype:

  TrimIsoInterStruct **TrimIntersectTrimCrvIsoVals(const TrimSrfStruct *TrimSrf,
                                                   int Dir,
                                                   CagdRType *OrigIsoParams,
                                                   int NumOfIsocurves,
                                                   CagdBType Perturb)


Description:

Computes the intersections of the trimming curves of TrimSrf with the ordered isoparametric values prescribed by OrigIsoParams, in axis Axis.

Parameters:

TrimSrf: Trimmed surface to consider.
Dir: Either U or V.
OrigIsoParams: Vectors of isoparametric values in direction Dir.
NumOfIsocurves: ize of vector IsoParams.
Perturb: TRUE to epsilon-perturb the iso-param values.


Returned Value:

TrimIsoInterStruct **: A vector of size NumOfIsocurves, each contains a list of intersection parameter values.


See Also:

TrimSrf2Polylines TrimCrvAgainstTrimCrvs

Keywords:

isoparametric curves


TrimIsPointInsideTrimCrvs

(trim_aux.c:2518)

Prototype:

  CagdBType TrimIsPointInsideTrimCrvs(const TrimCrvStruct *TrimCrvs,
                                      CagdUVType UV)


Description:

Returns TRUE if the given UV value is inside the domain prescribed by the trimming curves.

Parameters:

TrimCrvs: Trimming curves to consider.
UV: Parametric location.


Returned Value:

CagdBType: TRUE if inside, FALSE otherwise.


See Also:

TrimIsPointInsideTrimUVCrv TrimIsPointInsideTrimSrf MdlIsPointInsideTrimSrf

Keywords:

inclusion


TrimIsPointInsideTrimSrf

(trim_aux.c:2493)

Prototype:

  CagdBType TrimIsPointInsideTrimSrf(const TrimSrfStruct *TrimSrf, CagdUVType UV)


Description:

Returns TRUE if the given UV value is inside the trimmed surface's parametric domain.

Parameters:

TrimSrf: Trimmed surface to consider.
UV: Parametric location.


Returned Value:

CagdBType: TRUE if inside, FALSE otherwise.


See Also:

TrimIsPointInsideTrimUVCrv TrimIsPointInsideTrimCrvs

Keywords:

inclusion


TrimIsPointInsideTrimUVCrv

(trim_aux.c:2587)

Prototype:

  int TrimIsPointInsideTrimUVCrv(const CagdCrvStruct *UVCrv, CagdUVType UV)


Description:

Returns the number of times a ray in the -V direction from UV crosses UVCrv.

Parameters:

UVCrv: Trimming curve to consider.
UV: Parametric location.


Returned Value:

int: Number of crossings of UVCrv by a ray from UV in -V dir.


See Also:

TrimIsPointInsideTrimCrvs TrimIsPointInsideTrimSrf MdlIsPointInsideTrimSrf TrimIsPointInsideTrimUVCrvs

Keywords:

inclusion


TrimIsPointInsideTrimUVCrvs

(trim_aux.c:2557)

Prototype:

  int TrimIsPointInsideTrimUVCrvs(const CagdCrvStruct *UVCrvs, CagdUVType UV)


Description:

Returns the number of times a ray in the -V direction from UV crosses list of curves UVCrvs.

Parameters:

UVCrvs: Trimming curves to consider.
UV: Parametric location.


Returned Value:

int: Number of crossings of UVCrvs by a ray from UV in -V dir.


See Also:

TrimIsPointInsideTrimCrvs TrimIsPointInsideTrimSrf MdlIsPointInsideTrimSrf TrimIsPointInsideTrimUVCrv

Keywords:

inclusion


TrimLinkTrimmingCurves2Loops

(trim_aux.c:1157)

Prototype:

  TrimCrvStruct *TrimLinkTrimmingCurves2Loops(const TrimCrvStruct *TCrvs,
                                              CagdRType MaxTol,
                                              CagdBType *ClosedLoops)


Description:

Link all given curves into loops. Loops can be still open as they can start/end on the boundary, a curve someone else will have to complete... Aims to link the curves into loops from a very tight tolerance and until MaxTol or the loop was detected as closed.

Parameters:

TCrvs: Trimming curves to link into loops. Each, assumed to hold one or more trimming segments.
MaxTol: Maximal tolerance of end points comparisons.
ClosedLoops: Optional return parameter - if not NULL, Will be set to TRUE only if all linked loops are indeed closed to within tolerance.


Returned Value:

TrimCrvStruct *: Trimming curve segments linked into loops.


See Also:

TrimLinkTrimmingCurves2Loops2 TrimLinkTrimmingCurves2Loops1 TrimMergeTrimmingCurves2Loops

Keywords:

trimming curves


TrimLinkTrimmingCurves2Loops1

(trim_aux.c:1107)

Prototype:

  TrimCrvStruct *TrimLinkTrimmingCurves2Loops1(const TrimCrvSegStruct *TSegs,
                                               CagdRType MaxTol,
                                               CagdBType *ClosedLoops)


Description:

Link all given curve segments into loops. Loops can be still open as they can start/end on the boundary, a curve someone else will have to complete... Aims to link the curves into loops from a very tight tolerance and until MaxTol or the loop was detected as closed.

Parameters:

TSegs: Trimming curve segments to link into loops. Each, must hold one trimming segements.
MaxTol: Maximal tolerance of end points comparisons.
ClosedLoops: Optional return parameter - if not NULL, wWill be set to TRUE only if all linked loops are indeed closed to within tolerance.


Returned Value:

TrimCrvStruct *: Trimming curve segments linked into loops.


See Also:

TrimLinkTrimmingCurves2Loops2 TrimMergeTrimmingCurves2Loops

Keywords:

trimming curves


TrimLinkTrimmingCurves2Loops2

(trim_aux.c:1282)

Prototype:

  TrimCrvStruct *TrimLinkTrimmingCurves2Loops2(TrimCrvStruct *TCrvs,
                                               CagdRType Tol,
                                               CagdBType *ClosedLoops)


Description:

Link all given curves into loops. Loops can be still open as they can start/end on the boundary, a curve someone else will have to complete...

Parameters:

TCrvs: Trimming curves to link into loops. Each can hold a list of trimming segments.
Tol: Tolerance of end points comparisons and link.
ClosedLoops: Will be set to TRUE only if all linked loops are indeed closed to within tolerance.


Returned Value:

TrimCrvStruct *: Trimming curve segments linked into loops.


See Also:

TrimLinkTrimmingCurves2Loops TrimMergeTrimmingCurves2Loops TrimLinkTrimmingCurves2Loops2Inc

Keywords:

trimming curves


TrimLoopUV2Weight

(trimcntr.c:1072)

Prototype:

  CagdRType TrimLoopUV2Weight(const IrtRType *UV,
                              IrtRType *BndryUV,
                              CagdRType UMin,
                              CagdRType UMax,
                              CagdRType VMin,
                              CagdRType VMax,
                              CagdBType Last)


Description:

Computes a unique weight according to the location where the vertex meets the boundary: VMin is (0-1), UMax is (1-2), VMax is (2-3), UMin is (3-4).

Parameters:

UV: UV boundary location To assign a weight.
BndryUV: If not NULL, updated with the closest point on the boundary to UV.
UMin, UMax, VMin, VMax: Domain of surface.
Last: Weight 0 and 4 are distinguished by UV being Last (Weight is 4) or not (Weight is 0).


Returned Value:

CagdRType: Weight of UV.


See Also:

TrimSrfsFromContours TrimSrfsFromTrimPlsHierarchy TrimLoopWeight2UVToData TrimLoopWeightRelationInside

Keywords:




TrimLoopWeight2UVToData

(trimcntr.c:1140)

Prototype:

  CagdRType *TrimLoopWeight2UVToData(IrtRType Wgt,
                                     CagdRType UMin,
                                     CagdRType UMax,
                                     CagdRType VMin,
                                     CagdRType VMax,
                                     CagdUVType UV)


Description:

Computes the UV value associated with the given weight on the given boundary. Weight W must be in the range [0, 4]: VMin is (0-1), UMax is (1-2), VMax is (2-3), UMin is (3-4).

Parameters:

Wgt: Weight to convert to a UV boundary location.
UMin, UMax, VMin, VMax: Domain of surface.
UV: UV of Wgt.


Returned Value:

CagdRType *: UV of Wgt.


See Also:

TrimSrfsFromContours TrimSrfsFromTrimPlsHierarchy TrimLoopUV2Weight TrimLoopWeightRelationInside

Keywords:




TrimLoopWeightRelationInside

(trimcntr.c:1032)

Prototype:

  int TrimLoopWeightRelationInside(CagdRType V1,
                                   CagdRType V2,
                                   CagdRType V)


Description:

Returns TRUE iff V is between V1 and V2, going counter clockwise along the boundary of the parametric domain.

Parameters:

V1, V2: The ordered end points to verify that V is inside.
V: The vertex to examine if between V1 and V2.


Returned Value:

int: TRUE if inside (between V1 and V2).


See Also:

TrimSrfsFromContours TrimSrfsFromTrimPlsHierarchy TrimLoopWeight2UVToData TrimLoopUV2Weight

Keywords:




TrimManageTrimmingCurvesDegrees

(trim_aux.c:885)

Prototype:

  TrimSrfStruct *TrimManageTrimmingCurvesDegrees(TrimSrfStruct *TrimSrf,
                                                 int FitOrder,
                                                 CagdBType EvalEuclid)


Description:

Manage all trimming curves of a given surface as follows: 1. If FitOrder == 2, a piecewise linear (approximation) is computed for any higher order trimming curve using globals _TrimUVCrvApproxMethod and _TrimUVCrvApproxTolSamples approximation specs. 2. If FitOrder > 2, a single polynomial of FitOrder order is fitted for all trimming curves. 3. If FitOrder < 2, nothing is done.

Parameters:

TrimSrf: Trimmed surface to extract trimming curves from, in place.
FitOrder: The order of the newly fitted trimming curves.
EvalEuclid: If TRUE reevaluate Euclidean curve as well.


Returned Value:

TrimSrfStruct *: The trimmed surface, modified in place, that holds the newly created trimming curves.


See Also:

TrimGetTrimmingCurves TrimSetTrimCrvLinearApprox TrimCrv2Polyline

Keywords:

trimming curves


TrimMatch2ndCrvLenSpeedAs1stCrv

(trim2ply.c:1230)

Prototype:

  void TrimMatch2ndCrvLenSpeedAs1stCrv(CagdCrvStruct **Crv1,
                                       CagdCrvStruct **Crv2,
                                       const CagdSrfStruct *Srf1,
                                       const CagdSrfStruct *Srf2)


Description:

Change the length of the given 2nd curve to be the same as the 1st Crv, by refinement. 1st Crv, Crv1, is assumed to be longer (or equal length).

Parameters:

Crv1: Base curve to force Crv2 to match its length. If Crv1 is an isoparametric curve it is forced to a uniform speed. Crv1 can also be refined in this function.
Crv2: 2nd curve to change its length to follow Crv1, in place.
Srf1, Srf2: f not NULL, Srf1 to compose with Crv1 (that is assumed a UV curve in Srf1), and we project the points of Crv1 onto the E3 curve of Crv2 or (Srf2(Crv2)). This, so we establish a precise match of the parametrization.


Returned Value:

void


See Also:

MvarProjUVCrvOnE3CrvMatchSpeed

Keywords:




TrimMergePolylines

(trim2ply.c:1438)

Prototype:

  IPPolygonStruct *TrimMergePolylines(IPPolygonStruct *Polys, IrtRType Eps)


Description:

Merges separated polylines into longer ones, in place, as possible. Given a list of polylines, matches end points and merged as possible polylines with common end points, in place.

Parameters:

Polys: Polylines to merge, in place.
Eps: Epsilon of similarity to merge points at.


Returned Value:

IPPolygonStruct *: Merged as possible polylines.


See Also:

GMMergeGeometry MvarPolyMergePolylines GMMergePolylines GMMergePolylines2

Keywords:

merge polyline


TrimMergeTrimmingCurves2Loops

(trim_aux.c:1708)

Prototype:

  TrimCrvStruct *TrimMergeTrimmingCurves2Loops(const TrimCrvStruct *TrimCrvs)


Description:

Merges all given trimming curves into closed loops, in place. Only UV curves are merged and Euclidean trimming curves are purged.

Parameters:

TrimCrvs: Trimming curves to merge into loops.


Returned Value:

TrimCrvStruct *: Trimming curves merged into loops.


See Also:

TrimMergeTrimmingCurves2Loops2 TrimLinkTrimmingCurves2Loops

Keywords:

trimming curves


TrimMergeTrimmingCurves2Loops2

(trim_aux.c:1793)

Prototype:

  CagdCrvStruct *TrimMergeTrimmingCurves2Loops2(CagdCrvStruct *UVCrvs,
                                                CagdRType Tol)


Description:

Merges all given curves into closed loops, in place.

Parameters:

UVCrvs: Curves to merge into loops, in place.
Tol: Tolerance of end points comparisons.


Returned Value:

CagdCrvStruct *: Curves merged into loops.


See Also:

TrimMergeTrimmingCurves2Loops TrimLinkTrimmingCurves2Loops

Keywords:

trimming curves


TrimOrderTrimCrvSegsInLoop

(trim2ply.c:1717)

Prototype:

  TrimCrvSegStruct *TrimOrderTrimCrvSegsInLoop(TrimCrvSegStruct *TSegs)


Description:

Make sure the given loop is ordered so the end of the i'th segment is the beginning of the i+1'th segment by comapring end points.

Parameters:

TSegs: List of curve segments to reorder, in place.


Returned Value:

TrimCrvSegStruct *: Ordered list.


See Also:



Keywords:




TrimOrientTrimingCrvs

(trim2ply.c:1815)

Prototype:

  CagdBType TrimOrientTrimingCrvs(TrimSrfStruct *TSrf)


Description:

Properly orient the trimming curves of the trimmed surface so walking along them, the valid portion of the trimmed surface will be on the same side, in the parametric domain.

Parameters:

TSrf: Trimmed surface to properly orient its trimming curves.


Returned Value:

CagdBType: TRUE if successful, FALSE otherwise.


Keywords:




TrimPiecewiseRuledSrfApprox

(tr_prisa.c:130)

Prototype:

  TrimSrfStruct *TrimPiecewiseRuledSrfApprox(const TrimSrfStruct *CTSrf,
                                             CagdBType ConsistentDir,
                                             CagdRType Epsilon,
                                             CagdSrfDirType Dir)


Description:

Constructs a piecewise ruled surface approximation to the given trimmed surface, TSrf, in the given direction, Dir, that is close to the surface to within Epsilon. If ConsitentDir then ruled surface parametrization is set to be the same as original surface TSrf. Otherwise, ruling dir is always CAGD_CONST_V_DIR. Surface is assumed to have point types E3 or P3 only.

Parameters:

CTSrf: To approximate using piecewise ruled surfaces.
ConsistentDir: o we want parametrization to be the same as TSrf?
Epsilon: Accuracy of piecewise ruled surface approximation.
Dir: Direction of piecewise ruled surface approximation. Either U or V.


Returned Value:

TrimSrfStruct *: A list of trimmed ruled surfaces approximating TSrf to within Epsilon in direction Dir.


See Also:

TrimAllPrisaSrfs TrimPrisaRuledSrf SymbAllPrisaSrfs

Keywords:

layout prisa ruled surface approximation


TrimPointInsideTrimmedCrvsToData

(trim_aux.c:547)

Prototype:

  CagdRType *TrimPointInsideTrimmedCrvsToData(TrimCrvStruct *TrimCrvList,
                                              const TrimSrfStruct *TSrf,
                                              CagdUVType UVRetVal)


Description:

Finds a point inside a set of trimmed crvs. Returned is a UV location

Parameters:

TrimCrvList: To find a location inside it.
TSrf: If provided, will attempt to find a point inside the trimmed curve from the surface boundary. If NULL, an interior point to the trimming curves will be selected.
UVRetVal: A location in the parametric space of the surface that is part of the valid trimmed surface domain.


Returned Value:

CagdRType *: A location in the parametric space of the surface that is part of the valid trimmed surface domain. NULl is returned if failed to compute.


Keywords:




TrimPolylines2LinTrimCrvs

(trimcntr.c:637)

Prototype:

  TrimCrvStruct *TrimPolylines2LinTrimCrvs(const IPPolygonStruct *Polys)


Description:

Returns a list of linear Bspline curves constructed from given polylines.

Parameters:

Polys: Input polylines to convert to lienar bspline curves.


Returned Value:

TrimCrvStruct *: Linear Bspline curves representing Polys.


Keywords:




TrimPrisaRuledSrf

(tr_prisa.c:381)

Prototype:

  TrimSrfStruct *TrimPrisaRuledSrf(const TrimSrfStruct *TSrf,
                                   int SamplesPerCurve,
                                   CagdRType Space,
                                   CagdVType Offset,
                                   CagdSrfDirType Dir)


Description:

Layout a single trimmed ruled surface, by approximating it as a set of polygons. The given trimmed ruled surface might be non-developable, in which case approximation will be of a surface with no twist. The trimmed ruled surface is assumed to be constructed using CagdRuledSrf and that the ruled direction is consistent and is always CAGD_CONST_V_DIR.

Parameters:

TSrf: A trimmed ruled surface to layout flat on the XY plane.
SamplesPerCurve: uring the approximation of a ruled surface as a developable surface.
Space: Increment on Y on the offset vector, after this surface was placed in the XY plane.
Offset: A vector in the XY plane to denote the amount of translation for the flatten surface in the XY plane.
Dir: Direction of piecewise ruled surface approximation. Either U or V.


Returned Value:

TrimSrfStruct *: A planar trimmed surface in the XY plane approximating the falttening process of TSrf.


See Also:

TrimPiecewiseRuledSrfApprox TrimAllPrisaSrfs SymbAllPrisaSrfs

Keywords:

layout prisa


TrimRemovEucTrimCrvs

(trim_aux.c:2456)

Prototype:

  void TrimRemovEucTrimCrvs(TrimSrfStruct *TSrf)


Description:

Remove the Euclidean curves from the trimming curves in the given data.

Parameters:

TSrf: A trimmed surface to remove all trimming curves from.


Returned Value:

void


See Also:

MdlRemovEucTrimCrvs

Keywords:




TrimRemoveCrvSegTrimCrvSegs

(trim_sub.c:854)

Prototype:

  int TrimRemoveCrvSegTrimCrvSegs(TrimCrvSegStruct *TrimCrvSeg,
                                  TrimCrvSegStruct **TrimCrvSegs)


Description:

Removes but not delete the given trimming crv segment from the list of trimming curve segments pointed by TrimCrvSegs.

Parameters:

TrimCrvSeg: Segment to delete.
TrimCrvSegs: List of trimming curve segments to delete TrimCrvSeg from.


Returned Value:

int: TRUE if found and removed, FALSE otherwise.


See Also:

TrimRemoveCrvSegTrimCrvs

Keywords:




TrimRemoveCrvSegTrimCrvs

(trim_sub.c:801)

Prototype:

  int TrimRemoveCrvSegTrimCrvs(TrimCrvSegStruct *TrimCrvSeg,
                               TrimCrvStruct **TrimCrvs)


Description:

Removes but not delete the given trimming crv segment from the list of trimming curves point by TrimCrvs.

Parameters:

TrimCrvSeg: Segment to delete.
TrimCrvs: List of trimming curves to delete TrimCrvSeg from.


Returned Value:

int: TRUE if found and removed, FALSE otherwise.


See Also:

TrimRemoveCrvSegTrimCrvSegs

Keywords:




TrimSetEuclidComposedFromUV

(trim_aux.c:2396)

Prototype:

  int TrimSetEuclidComposedFromUV(int EuclidComposedFromUV)


Description:

Sets the way Euclidean trimming curves are computed from parametric trimming curves. Either by symbolic composition (TRUE) or by piecewise linear approximation of trimming curves (FALSE).

Parameters:

EuclidComposedFromUV: Do we want symbolic composition for Euclidean curves, or should we piecewise linear sample the UV trimming curves.


Returned Value:

int: Old value of way of Euclidean curve's computation


See Also:

TrimCrvs2Polylines TrimSetEuclidLinearFromUV

Keywords:




TrimSetEuclidLinearFromUV

(trim_aux.c:2427)

Prototype:

  int TrimSetEuclidLinearFromUV(int EuclidLinearFromUV)


Description:

Sets the way Euclidean trimming curves are computed from parametric trimming curves. TRUE to ensure output is piecewise linear.

Parameters:

EuclidLinearFromUV: TRUE to ensure Euclidean curves evaluated are made piecewise linear.


Returned Value:

int: Old value of piecewise linear requirement.


See Also:

TrimCrvs2Polylines TrimSetEuclidComposedFromUV

Keywords:




TrimSetFatalErrorFunc

(trim_ftl.c:28)

Prototype:

  TrimSetErrorFuncType TrimSetFatalErrorFunc(TrimSetErrorFuncType ErrorFunc)


Description:

Sets the error function to be used by Trim_lib.

Parameters:

ErrorFunc: New error function to use.


Returned Value:

TrimSetErrorFuncType: Old error function reference.


Keywords:

error handling


TrimSetNumTrimVrtcsInCell

(trim2pl2.c:340)

Prototype:

  int TrimSetNumTrimVrtcsInCell(int NumTrimVrtcsInCell)


Description:

Sets the way trimming curves contributes to each cell in the domain by setting the number of vertices trimming curves can contribute to each such cell.

Parameters:

NumTrimVrtcsInCell: Number of requested trimming vertices in a cell.


Returned Value:

int: Old value of way of num of cells.


Keywords:




TrimSetTrimCrvLinearApprox

(trim_iso.c:882)

Prototype:

  SymbCrvApproxMethodType TrimSetTrimCrvLinearApprox(CagdRType UVTolSamples,
                                                SymbCrvApproxMethodType UVMethod)


Description:

Sets the tolerances to use when approximating higher order trimming curves using piecewise linear approximation, for intersection computation.

Parameters:

UVTolSamples: Piecewise linear approximation of high order trimming curves - number of samples per curve or tolerance.
UVMethod: Method of sampling.


Returned Value:

SymbCrvApproxMethodType: Old method of curve sampling.


See Also:

TrimGetTrimCrvLinearApprox SymbCrv2Polyline TrimCrv2Polyline

Keywords:




TrimSrf2Curves

(trim_iso.c:114)

Prototype:

  CagdCrvStruct *TrimSrf2Curves(TrimSrfStruct *TrimSrf,
                                int NumOfIsocurves[2])


Description:

Routine to extract from a trimmed surface NumOfIsoline isocurve list in each param. direction. Iso parametric curves are sampled equally spaced in parametric space. NULL is returned in case of an error, otherwise list of CagdCrvStruct. As the isoparametric curves are trimmed according to the trimming curves the resulting number of curves is arbitrary.

Parameters:

TrimSrf: To extract isoparametric curves from.
NumOfIsocurves: In each (U or V) direction.


Returned Value:

CagdCrvStruct *: List of extracted isoparametric curves. These curves inherit the order and continuity of the original Srf. NULL is returned in case of an error.


Keywords:

curves isoparametric curves


TrimSrf2KnotCurves

(trim_iso.c:686)

Prototype:

  CagdCrvStruct *TrimSrf2KnotCurves(TrimSrfStruct *TrimSrf)


Description:

Extracts a list of knot curves along all knots in U and V of TrimSrf.

Parameters:

TrimSrf: To extract a knot curves from.


Returned Value:

CagdCrvStruct *: Extracted knot curves.


See Also:

CagdSrf2CtrlMesh CagdSrf2KnotLines CagdSrf2KnotPolylines CagdSrf2KnotCurves

Keywords:

control mesh knots knot lines knot curves


TrimSrf2Polygons2

(trim2pl2.c:178)

Prototype:

  IPPolygonStruct *TrimSrf2Polygons2(const TrimSrfStruct *CTrimSrf,
                                     CagdSrf2PlsInfoStrct *TessInfo)


Description:

Routine to convert a single trimmed surface to set of triangles approximating it. FineNess is a fineness control on result and the larger is more triangles may result. A value of 10 is a good start value. NULL is returned in case of an error, otherwise list of IPPolygonStruct. This routine looks for C1 discontinuities in the surface and splits it into C1 continuous patches to invoke TrimC1Srf2Polygons to gen. polygons.

Parameters:

CTrimSrf: To approximate into triangles.
TessInfo: All auxiliary information/state to tessellate.


Returned Value:

IPPolygonStruct *: A list of polygons with optional normal and/or UV parametric information. NULL is returned in case of an error.


See Also:

BspSrf2PolygonSetErrFunc BzrSrf2Polygons IritSurface2Polygons IritTrimSrf2Polygons CagdSrf2Polygons BspSrf2Polygons BspC1Srf2Polygons TrimSetNumTrimVrtcsInCell

Keywords:

polygonization surface approximation


TrimSrf2Polylines

(trim_iso.c:75)

Prototype:

  CagdPolylineStruct *TrimSrf2Polylines(TrimSrfStruct *TrimSrf,
                                        int NumOfIsocurves[2],
                                        CagdRType TolSamples,
                                        SymbCrvApproxMethodType Method,
                                        CagdSrf2PlsInfoStrct *TessInfo)


Description:

Routine to convert a single trimmed surface to NumOfIsolines polylines in each parametric direction with TolSamples samples/tolerance in each isoparametric curve. Polyline are always E3 of CagdPolylineStruct type. NULL is returned in case of an error, otherwise list of CagdPolylineStruct. Attempt is made to extract isolines along C1 discontinuities first.

Parameters:

TrimSrf: To extract isoparametric curves from.
NumOfIsocurves: n each (U or V) direction.
TolSamples: Tolerance of approximation error (Method = 2) or Number of samples to compute on polyline (Method = 0, 1).
Method: 0 - TolSamples are set uniformly in parametric space, 1 - TolSamples are set optimally, considering the isocurve's curvature. 2 - TolSamples sets the maximum error allowed between the piecewise linear approximation and original curve.
TessInfo: Auxiliary struct which hold the parameters for the creation of the polys. Can be NULL for default params.


Returned Value:

CagdPolylineStruct *: List of polylines representing a piecewise linear approximation of the extracted isoparametric curves or NULL is case of an error.


See Also:

TrimCrv2Polyline

Keywords:

isoparametric curves


TrimSrfAdap2Polygons

(trim2ply.c:143)

Prototype:

  IPPolygonStruct *TrimSrfAdap2Polygons(const TrimSrfStruct *TrimSrf,
                                        CagdSrf2PlsInfoStrct *TessInfo)


Description:

Routine to convert a single trimmed surface to set of polygons approximating it. Tolerance is a tolerance control on result, typically related to the the accuracy of the approximation. A value of 0.1 is a good rough start. NULL is returned in case of an error or use of call back function to get a hold over the created polygons, otherwise list of IPPolygonStruct. This routine looks for C1 discontinuities in the surface and splits it into C1 continuous patches first.

Parameters:

TrimSrf: To approximate into triangles.
TessInfo: All auxiliary information/state to tessellate.


Returned Value:

IPPolygonStruct *: A list of polygons with optional normal and/or UV parametric information. NULL is returned in case of an error or if use of call back function to collect the polygons.


See Also:

CagdSrf2PolygonSetErrFunc CagdSrfAdap2PolyDefErrFunc CagdSrf2PolyAdapSetErrFunc CagdSrf2PolyAdapSetAuxDataFunc CagdSrf2PolyAdapSetPolyGenFunc CagdSrfAdap2Polygons TrimSrf2Polygons

Keywords:

polygonization surface approximation


TrimSrfBBox

(trim_aux.c:283)

Prototype:

  CagdBBoxStruct *TrimSrfBBox(const TrimSrfStruct *TSrf, CagdBBoxStruct *BBox)


Description:

Computes a bounding box for a freeform trimmed surface.

Parameters:

TSrf: To compute a bounding box for.
BBox: Where bounding information is to be saved.


Returned Value:

CagdBBoxStruct *: BBox.


See Also:

TrimSrfListBBox CagdSrfBBox GMBBSetBBoxPrecise

Keywords:

bbox bounding box


TrimSrfCnvrt2BzrRglrSrf

(untrim.c:192)

Prototype:

  CagdSrfStruct *TrimSrfCnvrt2BzrRglrSrf(const TrimSrfStruct *TrimSrf)


Description:

Given a trimmed surface - subdivides it into trimmed Bezier surfaces and convert the trimmed surface into regular Bezier/Bspline patches. Returns a list of non-trimmed Bezier surfaces, that together, are identical geometrically to the input trimmed surface.

Parameters:

TrimSrf: A trimmed Bezier or a Bspline surface to convert to non-trimmed Bezier surfaces.


Returned Value:

CagdSrfStruct *: The non trimmed regular Bezier surfaces.


See Also:

TrimSrfSubdivAtParam TrimSrfCnvrt2BzrTrimSrf

Keywords:

subdivision


TrimSrfCnvrt2BzrRglrSrf2

(untrim.c:293)

Prototype:

  CagdSrfStruct *TrimSrfCnvrt2BzrRglrSrf2(const TrimSrfStruct *TSrf,
                                          int ComposeE3,
                                          int OnlyBzrSrfs,
                                          CagdRType Eps)


Description:

Divides the given trimmed surface TSrf to a list of tensor products. The result is tiling the original trimmed surface to within machine precision. The given TSrf surface is recursively divided until the trimming curves are simple, in which case the trimmed surface is converted into a regular tensor product surface. A trimmed surface is considered simple if it has a single trimming loop that is double monotone with U or V (That is from the minimal point to the maximal point in U or V we have monotone progress in two separated paths).

Parameters:

TSrf: Trimmed surface to decompose into a list of tensor product (and no trimming) surfaces.
ComposeE3: TRUE to compose the tiles into TSrf, FALSE to return the surface tiles in the parametric domain of TSrf.
OnlyBzrSrfs: TRUE to force only Bezier tensor products in result.
Eps: Tolerance of the decomposition.


Returned Value:

CagdSrfStruct *: A list of regular tensor product surfaces that represents the same region as TSrf, to within machine precision.


See Also:

TrimCrvIsParameterizableDomain Trim2DSrfFromDoubleMonotoneTrimCrv TrimSrfCnvrt2BzrRglrSrf

Keywords:




TrimSrfCnvrt2BzrTrimSrf

(untrim.c:62)

Prototype:

  TrimSrfStruct *TrimSrfCnvrt2BzrTrimSrf(const TrimSrfStruct *TrimSrf)


Description:

Given a trimmed surface - subdivide it into trimmed Bezier surfaces (each spanning domain [0, 1]^2). Returns a list of trimmed Bezier surfaces, that together, are identical geometrically to the input trimmed surface.

Parameters:

TrimSrf: A Bezier or a B-spline trimmed surface to convert to Bezier.


Returned Value:

TrimSrfStruct *: The subdivided Bezier trimmed surfaces.


See Also:

TrimSrfSubdivAtParam

Keywords:

subdivision


TrimSrfCnvrt2TensorProdSrf

(untrim.c:358)

Prototype:

  CagdSrfStruct *TrimSrfCnvrt2TensorProdSrf(const TrimSrfStruct *TSrf,
                                            int ComposeE3,
                                            CagdRType Eps)


Description:

Divides the given trimmed surface TSrf to a list of tensor product srfs. The result is tiling the original trimmed surface to within machine precision. The given TSrf surface is recursively divided until the trimming curves are simple, in which case the trimmed surface is converted into a regular tensor product surface. A trimmed surface is considered simple if it has a single trimming loop that is double monotone with U or V (That is from the minimal point to the maximal point in U or V we have monotone progress in two separated paths).

Parameters:

TSrf: Trimmed surface to decompose into a list of tensor product (and no trimming) surfaces.
ComposeE3: TRUE to compose the tiles into TSrf, FALSE to return the surface tiles in the parametric domain of TSrf.
Eps: Tolerance of the decomposition.


Returned Value:

CagdSrfStruct *: A list of regular tensor product surfaces that represents the same region as TSrf, to within machine precision.


See Also:

TrimCrvIsParameterizableDomain Trim2DSrfFromDoubleMonotoneTrimCrv TrimSrfCnvrt2BzrRglrSrf TrimSrfCnvrt2BzrRglrSrf2

Keywords:




TrimSrfCopy

(trim_gen.c:1240)

Prototype:

  TrimSrfStruct *TrimSrfCopy(const TrimSrfStruct *TrimSrf)


Description:

Duplicates a trimming surface structure.

Parameters:

TrimSrf: A trimming surface to duplicate.


Returned Value:

TrimSrfStruct *: A trimming surface structure.


Keywords:

allocation


TrimSrfCopyList

(trim_gen.c:1268)

Prototype:

  TrimSrfStruct *TrimSrfCopyList(const TrimSrfStruct *TrimSrfList)


Description:

Allocates and copies a list of trimming surface structures.

Parameters:

TrimSrfList: To be copied.


Returned Value:

TrimSrfStruct *: A duplicated list of trimming surfaces.


Keywords:

copy


TrimSrfDegreeRaise

(trim_aux.c:403)

Prototype:

  TrimSrfStruct *TrimSrfDegreeRaise(const TrimSrfStruct *TrimSrf,
                                    CagdSrfDirType Dir)


Description:

Returns a new trimmed surface representing the same surface as TrimSrf but with its degree raised by one.

Parameters:

TrimSrf: To raise its degree.
Dir: Direction of degree raising. Either U or V.


Returned Value:

TrimSrfStruct *: A surface with same geometry as Srf but with one degree higher.


See Also:

CagdSrfDegreeRaise

Keywords:

degree raising


TrimSrfDomain

(trim_aux.c:261)

Prototype:

  void TrimSrfDomain(const TrimSrfStruct *TrimSrf,
                     CagdRType *UMin,
                     CagdRType *UMax,
                     CagdRType *VMin,
                     CagdRType *VMax)


Description:

Returns the parametric domain of a trimmed surface.

Parameters:

TrimSrf: To get its parametric domain.
UMin: Where to put the minimal U domain's boundary.
UMax: Where to put the maximal U domain's boundary.
VMin: Where to put the minimal V domain's boundary.
VMax: Where to put the maximal V domain's boundary.


Returned Value:

void


See Also:

CagdSrfDomain

Keywords:

domain parametric domain


TrimSrfEvalMalloc

(trim_aux.c:345)

Prototype:

  CagdRType *TrimSrfEvalMalloc(const TrimSrfStruct *TrimSrf,
                               CagdRType u,
                               CagdRType v)


Description:

Given a trimmed surface and parameter values u, v, evaluate the surface at (u, v). No test is made to make sure (u, v) is in the untrimmed domain.

Parameters:

TrimSrf: o evaluate at the given parametric location (u, v).
u, v: The parameter values at which TrimSrf is to be evaluated.


Returned Value:

CagdRType *: A vector holding all the coefficients of all components of surface TrimSrf's point type. If, for example, TrimSrf's point type is P2, the W, X, and Y will be saved in the first three locations of the returned vector. The first location (index 0) of the returned vector is reserved for the rational coefficient W and XYZ always starts at second location of the returned vector (index 1). In dynamic memory.


See Also:

CagdSrfEval TrimSrfEvalToData

Keywords:

evaluation


TrimSrfEvalToData

(trim_aux.c:378)

Prototype:

  void TrimSrfEvalToData(const TrimSrfStruct *TrimSrf,
                         CagdRType u,
                         CagdRType v,
                         CagdRType *Pt)


Description:

Given a trimmed surface and parameter values u, v, evaluate the surface at (u, v). No test is made to make sure (u, v) is in the untrimmed domain.

Parameters:

TrimSrf: o evaluate at the given parametric location (u, v).
u, v: The parameter values at which TrimSrf is to be evaluated.
Pt: A vector holding all the coefficients of all components of surface TrimSrf's point type. If, for example, TrimSrf's point type is P2, the W, X, and Y will be saved in the first three locations of the returned vector. The first location (index 0) of the returned vector is reserved for the rational coefficient W and XYZ always starts at second location of the returned vector (index 1).


Returned Value:

void


See Also:

CagdSrfEval TrimSrfEvalMalloc

Keywords:

evaluation


TrimSrfFree

(trim_gen.c:1297)

Prototype:

  void TrimSrfFree(TrimSrfStruct *TrimSrf)


Description:

Deallocates a trimmed surface structure.

Parameters:

TrimSrf: trimmed surface to free.


Returned Value:

void


Keywords:

allocation


TrimSrfFreeEuclideanTrimCrvs

(trim_aux.c:424)

Prototype:

  void TrimSrfFreeEuclideanTrimCrvs(TrimSrfStruct *TrimSrf)


Description:

Frees all Euclidean trimming curves of trimmed surface TrimSrf.

Parameters:

TrimSrf: To free all its Euclidean trimming curves. UV curves are left as is.


Returned Value:

void


Keywords:




TrimSrfFreeList

(trim_gen.c:1321)

Prototype:

  void TrimSrfFreeList(TrimSrfStruct *TrimSrfList)


Description:

Deallocates a list of trimmed surface structures.

Parameters:

TrimSrfList: list of trimmed surface to free.


Returned Value:

void


Keywords:

allocation


TrimSrfFromE3TrimmingCurves

(trim_gen.c:570)

Prototype:

  TrimSrfStruct *TrimSrfFromE3TrimmingCurves(TrimCrvStruct *TCrvs,
                                             const IrtPlnType Plane)


Description:

Constructs a planar trimmed surface, large enough to accomodate the TCrvs as rimming curves. The trimming curves are assumed to hold E3 curves only and are used in place.

Parameters:

TCrvs: Trimming curves to use, with only E3 information. Assumed first trimming curve is the outer largest one. Used in place and UV curves are added on the fly.
Plane: Optional place to consider. If NULL computed from TCrvs.


Returned Value:

TrimSrfStruct *: Constructed planar trimmed surface.


See Also:

TrimSrfNew CagdPrimPlaneFromE3Crv MvarInverseCrvOnSrfProj

Keywords:




TrimSrfFromSrf

(trimcntr.c:85)

Prototype:

  TrimSrfStruct *TrimSrfFromSrf(CagdSrfStruct *Srf, int SingleTCrv)


Description:

Build a trimmed surface from the given tensor product surface Srf, that spans all the domain of Srf.

Parameters:

Srf: Surface to convert into a trimmed surface, in place.
SingleTCrv: TRUE a single boundary trimming e loop. or FALSE to have four boundary trimming curves (for UMin/Max, VMin/Max). * *


Returned Value:

TrimSrfStruct *: Constructed trimmed surface.


See Also:

TrimSrfNew

Keywords:




TrimSrfListBBox

(trim_aux.c:305)

Prototype:

  CagdBBoxStruct *TrimSrfListBBox(const TrimSrfStruct *TSrfs, CagdBBoxStruct *BBox)


Description:

Computes a bounding box for a list of freeform trimmed surfaces.

Parameters:

TSrfs: To compute a bounding box for.
BBox: Where bounding information is to be saved.


Returned Value:

CagdBBoxStruct *: BBox.


See Also:

TrimSrfBBox CagdSrfBBox GMBBSetBBoxPrecise

Keywords:

bbox bounding box


TrimSrfListMatTransform

(trim_gen.c:1492)

Prototype:

  TrimSrfStruct *TrimSrfListMatTransform(const TrimSrfStruct *TrimSrfs,
                                         CagdMType Mat)


Description:

Transforms, in place, the given trimmed surface as specified by a homogeneous matrix Mat. Applies an homogeneous transformation, to the given list of surfaces Srfs as specified by homogeneous transformation Mat.

Parameters:

TrimSrfs: To be transformed.
Mat: Defining the transformation.


Returned Value:

TrimSrfStruct *: Returned transformed surfaces.


See Also:

TrimSrfMatTransform TrimSrfMatTransform2

Keywords:

scaling rotation translation transformations


TrimSrfMatTransform

(trim_gen.c:1438)

Prototype:

  TrimSrfStruct *TrimSrfMatTransform(const TrimSrfStruct *TrimSrf,
                                     CagdMType Mat)


Description:

Transforms the given trimmed surface as specified by a homogeneous matrix Mat.

Parameters:

TrimSrf: Trimmed surface to transform.
Mat: Homogeneous transformation to apply to trimmed surface.


Returned Value:

TrimSrfStruct *: Transformed trimmed surface.


See Also:

TrimSrfListMatTransform TrimSrfMatTransform2

Keywords:

Trimmed surface


TrimSrfMatTransform2

(trim_gen.c:1392)

Prototype:

  void TrimSrfMatTransform2(TrimSrfStruct *TrimSrf, CagdMType Mat)


Description:

Transforms, in place, the given trimmed surface as specified by a homogeneous matrix Mat.

Parameters:

TrimSrf: Trimmed surface to transform.
Mat: Homogeneous transformation to apply to trimmed surface.


Returned Value:

void


See Also:

TrimSrfListMatTransform TrimSrfMatTransform

Keywords:

Trimmed surface


TrimSrfNew

(trim_gen.c:403)

Prototype:

  TrimSrfStruct *TrimSrfNew(CagdSrfStruct *Srf,
                            TrimCrvStruct *TrimCrvList,
                            CagdBType HasTopLvlTrim)


Description:

Constructor for a trimmed surface.

Parameters:

Srf: Surface to make into a trimmed surface. Used in place.
TrimCrvList: A list of trimming curves, used in place.
HasTopLvlTrim: Do we have a top level outer most trimming curve?


Returned Value:

TrimSrfStruct *: The trimmed surface.


See Also:

TrimSrfNew2 TrimSrfNew3

Keywords:

allocation


TrimSrfNew2

(trim_gen.c:472)

Prototype:

  TrimSrfStruct *TrimSrfNew2(CagdSrfStruct *Srf,
                             CagdCrvStruct *TrimCrvList,
                             CagdBType HasTopLvlTrim)


Description:

Constructor for a trimmed surface.

Parameters:

Srf: Surface to make into a trimmed surface. Used in place.
TrimCrvList: A list of trimming curves, as regular curves, used in place.
HasTopLvlTrim: Do we have a top level outer most trimming curve?


Returned Value:

TrimSrfStruct *: The trimmed surface.


See Also:

TrimSrfNew TrimSrfNew3

Keywords:

allocation


TrimSrfNew3

(trim_gen.c:512)

Prototype:

  TrimSrfStruct *TrimSrfNew3(CagdSrfStruct *Srf,
                             CagdCrvStruct *TrimCrvList,
                             CagdBType HasTopLvlTrim)


Description:

Constructor for a trimmed surface. Same as TrimSrfNew2 but also verify and forces all trimming curves to be in the domain of Srf.

Parameters:

Srf: Surface to make into a trimmed surface. Used in place.
TrimCrvList: A list of trimming curves, used in place.
HasTopLvlTrim: Do we have a top level outer most trimming curve?


Returned Value:

TrimSrfStruct *: The trimmed surface.


See Also:

TrimSrfNew TrimSrfNew2

Keywords:

allocation


TrimSrfNumOfTrimCrvSegs

(trim_aux.c:222)

Prototype:

  int TrimSrfNumOfTrimCrvSegs(const TrimSrfStruct *TSrf)


Description:

Counts how many trimming curves segments this trimmed surface has.

Parameters:

TSrf: The trimmed surface to consider.


Returned Value:

int: Number of trimming curve segments.


See Also:

TrimSrfNumOfTrimLoops

Keywords:




TrimSrfNumOfTrimLoops

(trim_aux.c:195)

Prototype:

  int TrimSrfNumOfTrimLoops(const TrimSrfStruct *TSrf)


Description:

Counts how many trimming loops this trimmed surface has.

Parameters:

TSrf: The trimmed surface to consider.


Returned Value:

int: Number of trimming loops.


See Also:

TrimSrfNumOfTrimCrvSegs

Keywords:




TrimSrfRefineAtParams

(trim_aux.c:665)

Prototype:

  TrimSrfStruct *TrimSrfRefineAtParams(const TrimSrfStruct *TrimSrf,
                                       CagdSrfDirType Dir,
                                       CagdBType Replace,
                                       CagdRType *t,
                                       int n)


Description:

Given a trimmed surface - refines it at the given n knots as defined by vector t. If Replace is TRUE, the values in t replaces current knot vector. Returns pointer to refined surface (Note a Bezier surface will be converted into a Bspline surface).

Parameters:

TrimSrf: To refine.
Dir: Direction of refinement. Either U or V.
Replace: If TRUE, t holds knots in exactly the same length as the length of the knot vector of Srf and t simply replaces the knot vector.
t: Vector of knots with length of n.
n: Length of vector t.


Returned Value:

TrimSrfStruct *: A refined surface of TrimSrf after insertion of all the knots as specified by vector t of length n.


See Also:

CagdSrfRefineAtParams

Keywords:

refinement subdivision


TrimSrfRegionFromTrimSrf

(trim_aux.c:462)

Prototype:

  TrimSrfStruct *TrimSrfRegionFromTrimSrf(TrimSrfStruct *TrimSrf,
                                          CagdRType t1,
                                          CagdRType t2,
                                          CagdSrfDirType Dir)


Description:

Given a trimmed surface - extracts a sub-region within the domain specified by t1 and t2, in the direction Dir.

Parameters:

TrimSrf: To extract a sub-region from.
t1, t2: Parametric domain boundaries of sub-region.
Dir: Direction of region extraction. Either U or V.


Returned Value:

TrimSrfStruct *: Sub-region extracted from TrimSrf from t1 to t2.


See Also:

TrimSrfSubdivAtParam CagdSrfRegionFromSrf

Keywords:

regions subdivision


TrimSrfReverse

(trim_aux.c:690)

Prototype:

  TrimSrfStruct *TrimSrfReverse(const TrimSrfStruct *TrimSrf)


Description:

Returns a new trimmed surface that is the reversed surface of TrimSrf by reversing the control mesh and the knot vector (if Bspline surface) of TrimSrf in the U direction, as well as its trimming curves. See also CagdSrfReverse and BspKnotReverse.

Parameters:

TrimSrf: To be reversed.


Returned Value:

TrimSrfStruct *: Reversed surface of TrimSrf.


See Also:

TrimSrfReverse2 CagdSrfReverse

Keywords:

reverse


TrimSrfReverse2

(trim_aux.c:737)

Prototype:

  TrimSrfStruct *TrimSrfReverse2(const TrimSrfStruct *TrimSrf)


Description:

Returns a new trimmed surface that is the reversed surface of Srf by flipping the U and the V directions of the surface, as well as flipping them in the trimming curves. See also BspKnotReverse.

Parameters:

TrimSrf: To be reversed.


Returned Value:

TrimSrfStruct *: Reversed surface of TrimSrf.


See Also:

TrimSrfReverse CagdSrfReverse2

Keywords:

reverse


TrimSrfSetStateTrimCrvsManagement

(trim_sub.c:76)

Prototype:

  int TrimSrfSetStateTrimCrvsManagement(int TrimmingFitOrder)


Description:

Controls the way future trimmed surfaces subdivisions are made: If <2, trimming curves are left as is. If 2, all trimming curves are made piecewise linear as a side effect. If 2>, each trimming curve is fitted using a single polynomial of that degree.

Parameters:

TrimmingFitOrder: Order to fit trimming curves as defined above.


Returned Value:

int: Old state of this flag.


See Also:



Keywords:




TrimSrfSubdivAtInnerLoops

(trim_sub.c:1186)

Prototype:

  TrimSrfStruct *TrimSrfSubdivAtInnerLoops(TrimSrfStruct *TSrf)


Description:

Subdivides a trimmed surface into sub-surfaces that have no interal trimming loops. This function works by trying to detect a single internal trimming loop, and subdivide the entire trimmed surface at the middle of the internal loop's bounding box. This step typically removes at least one internal trimming loop. The function then calls itself recursively, until none of the sub-surfaces have any internal trimming loops.

Parameters:

TSrf: The trimmed surface to subdivide.


Returned Value:

TrimSrfStruct *: A list of trimmed surfaces with no internal trimming loops.


See Also:

TrimRemoveCrvSegTrimCrvs TrimSrfSubdivValAtInnerLoop

Keywords:




TrimSrfSubdivAtParam

(trim_sub.c:114)

Prototype:

  TrimSrfStruct *TrimSrfSubdivAtParam(const TrimSrfStruct *TSrf,
                                      CagdRType t,
                                      CagdSrfDirType Dir)


Description:

Given a trimmed surface - subdivides it into two sub-surfaces at given parametric value t in the given direction Dir. Returns pointer to a list of two trimmed surfaces, at most. It can very well may happen that the subdivided surface is completely trimmed out and hence nothing is returned for it.

Parameters:

TSrf: To subdivide at the prescribed parameter value t.
t: The parameter to subdivide the curve Crv at.
Dir: Direction of subdivision. Either U or V.


Returned Value:

TrimSrfStruct *: The subdivided surfaces. Usually two, but can have only one, if other is totally trimmed away.


See Also:

TrimSrfSubdivTrimmingCrvs

Keywords:

subdivision


TrimSrfSubdivTrimCrvsAtInnerLoops

(trim_sub.c:1243)

Prototype:

  TrimCrvStruct *TrimSrfSubdivTrimCrvsAtInnerLoops(const TrimCrvStruct *TCrvs)


Description:

Subdivides the given trimmed curves into sub-curves that have no interal trimming loops. This function works by detecting internal trimming loops, and subdividing the given trimming curves at the middle of the internal loop. This step removes at least one internal trimming loop. The function then calls itself recursively, until none of the trimmings have any internal trimming loops.

Parameters:

TCrvs: The trimmed curves to subdivide.


Returned Value:

TrimCrvStruct *: A list of trimmed curves with no internal trimming loops.


See Also:

TrimSrfSubdivInnerLoops TrimSrfSubdivValAtInnerLoop

Keywords:




TrimSrfSubdivTrimmingCrvs

(trim_sub.c:229)

Prototype:

  int TrimSrfSubdivTrimmingCrvs(const TrimCrvStruct *TrimCrvs,
                                CagdRType t,
                                CagdSrfDirType Dir,
                                TrimCrvStruct **TrimCrvs1,
                                TrimCrvStruct **TrimCrvs2)


Description:

Given a set of trimming curves - subdivides them into two groups below and above the subdividing line in direction Dir at parameter t.

Parameters:

TrimCrvs: To subdivide at the prescribed parameter value t.
t: The parameter to subdivide the curve Crv at.
Dir: Direction of subdivision. Either U or V.
TrimCrvs1: eturned first half of trimming curves, < t. Could be NULL.
TrimCrvs2: eturned second half of trimming curves, > t. Could be NULL.


Returned Value:

int: TRUE if successful and have two halves. FALSE if failed or have only one half.


See Also:

TrimSrfSubdivAtParam

Keywords:

subdivision


TrimSrfSubdivValAtInnerLoop

(trim_sub.c:1294)

Prototype:

  CagdSrfDirType TrimSrfSubdivValAtInnerLoop(const TrimCrvStruct *TCrvs,
                                             CagdRType *SubdivVal)


Description:

Given the trimmed curves, find an internal trimming loop, if any, and compute a location to divide the trimming curves (surface) so as to eliminate that internal loop.

Parameters:

TCrvs: The trimmed curves to subdivide.
SubdivVal: Will be updated with the subdivision value if found one.


Returned Value:

CagdSrfDirType: Subdivision direction to divide at or CAGD_NO_DIR if no internal loop was found.


See Also:

TrimSrfSubdivInnerLoops TrimSrfSubdivTrimCrvsAtInnerLoops

Keywords:




TrimSrfTransform

(trim_gen.c:1350)

Prototype:

  void TrimSrfTransform(TrimSrfStruct *TrimSrf,
                        const CagdRType *Translate,
                        CagdRType Scale)


Description:

Linearly transforms, in place, given trimmed surface as specified by Translate and Scale.

Parameters:

TrimSrf: Trimmed surface to transform.
Translate: Translation factor. Can be NULL for non.
Scale: Scaling factor.


Returned Value:

void


Keywords:




TrimSrfTrimCrvAllDomain

(trim_aux.c:2766)

Prototype:

  CagdBType TrimSrfTrimCrvAllDomain(const TrimSrfStruct *TrimSrf)


Description:

Examine the trimming curves of the given trimmed surface and returns TRUE iff the valid domain of trimming equals the entire surface domain.

Parameters:

TrimSrf: Trimmed surface to examine.


Returned Value:

CagdBType: TRUE if entire surface domain, FALSE otherwise.


See Also:

TrimSrfTrimCrvSquareDomain TrimClipSrfToTrimCrvs

Keywords:




TrimSrfTrimCrvSquareDomain

(trim_aux.c:2667)

Prototype:

  CagdBType TrimSrfTrimCrvSquareDomain(const TrimCrvStruct *TrimCrvList,
                                       CagdRType *UMin,
                                       CagdRType *UMax,
                                       CagdRType *VMin,
                                       CagdRType *VMax)


Description:

Examine the trimming curves of the given trimmed surface and returns TRUE iff the trimmed domain is a sub isoparametric square. In such a case the U/VMin/Max are set to the domain of the square.

Parameters:

TrimCrvList: Trimming curves to examine.
UMin, UMax, VMin, VMax: Domain of square, if return TRUE


Returned Value:

CagdBType: TRUE if a isoparametric square domain, FALSE otherwise.


See Also:

TrimSrfTrimCrvAllDomain TrimClipSrfToTrimCrvs

Keywords:




TrimSrfVerifyTrimCrvsValidity

(trim_gen.c:664)

Prototype:

  int TrimSrfVerifyTrimCrvsValidity(TrimSrfStruct *TrimSrf)


Description:

Verify that all trimming curves are indeed in the parametric domain of the surface and that all of them matches neighboring curves.

Parameters:

TrimSrf: To verify the validity of the trimming curves. This includes the verification of the continuity of the trimming loops and the inclusion in the domain of the trimming curves.


Returned Value:

int: TRUE if valid, FALSE if cannot correct the trimming curves.


Keywords:




TrimSrfsFromContours

(trimcntr.c:169)

Prototype:

  TrimSrfStruct *TrimSrfsFromContours(const CagdSrfStruct *Srf,
                                      const IPPolygonStruct *CCntrs)


Description:

Creates a set of trimmed surfaces as defined by given set of contours that can contain either closed or open contours. Open contours must terminate at the boundary of the parametric domain of the surface. Closed contours must be completely contained in the parametric domain with last point equals first.

Parameters:

Srf: To trim into pieces.
CCntrs: Polylines to use as separating edges.


Returned Value:

TrimSrfStruct *: List of trimmed surface pieces.


See Also:

TrimSrfsFromContours2 TrimSrfsFromTrimPlsHierarchy UserDivideSrfAtInterCrvs

Keywords:




TrimSrfsFromContours2

(trimcntr.c:395)

Prototype:

  TrimSrfStruct *TrimSrfsFromContours2(const CagdSrfStruct *Srf,
                                       const CagdCrvStruct *CCntrs)


Description:

Same as TrimSrfsFromContours after converting the curves to polylines.

Parameters:

Srf: To trim into pieces.
CCntrs: Curves to use as separating trimming curves.


Returned Value:

TrimSrfStruct *: List of trimmed surface pieces.


See Also:

TrimSrfsFromContours TrimSrfsFromTrimPlsHierarchy

Keywords:




TrimSrfsFromTrimPlsHierarchy

(trimcntr.c:443)

Prototype:

  TrimSrfStruct *TrimSrfsFromTrimPlsHierarchy(IPPolygonStruct *TopLevel,
                                              IPPolygonStruct *TrimPls,
                                              const CagdSrfStruct *Srf)


Description:

Construct trimmed surface from the given hierarchy of trimming polylines. If TopLevel is provided, it serves as the top level outer loop and both Odd and Even nested trimmed surfaces are extracted. If TopLevel is NULL only Odd nested trimmed surfaces are extracted.

Parameters:

TopLevel: The top level outer loop or NULL if none.
TrimPls: Hierarchy of trimming polylines.
Srf: Surface to trim out.


Returned Value:

TrimSrfStruct *: List of trimmed surface out of the given counters.


See Also:

TrimSrfsFromContours

Keywords:




TrimSrfsSame

(trim_gen.c:1526)

Prototype:

  CagdBType TrimSrfsSame(const TrimSrfStruct *TSrf1,
                         const TrimSrfStruct *TSrf2,
                         CagdRType Eps)


Description:

Compare the two trimmed surfaces for similarity.

Parameters:

TSrf1, TSrf2: The two trimmed surfaces to compare.
Eps: Tolerance of equality.


Returned Value:

CagdBType: TRUE if trimmed surfaces are the same, FALSE otehrwise.


See Also:

CagdSrfsSame

Keywords:




TrimUntrimSetLineSweepOutputCrvPairs

(untrim.c:1170)

Prototype:

  CagdBType TrimUntrimSetLineSweepOutputCrvPairs(CagdBType NewValue)


Description:

Configures whether the line-sweep algorithm should output surfaces or pairs of curves.

Parameters:

NewValue: TRUE to output pairs of curves, FALSE for surfaces.


Returned Value:

CagdBType: The normal at the middle of the surface patch.


See Also:

CagdQuadCrvToQuadsLineSweep

Keywords:

untrimming.


TrimUntrimmingResultFree

(untrim.c:1057)

Prototype:

  void TrimUntrimmingResultFree(TrimUntrimResultStruct *Untrim)


Description:

Deallocates an untrimming result structure.

Parameters:

Untrim: The untrimming result structure.


Returned Value:

void


Keywords:

untrimming


TrimUntrimmingResultFreeList

(untrim.c:1078)

Prototype:

  void TrimUntrimmingResultFreeList(TrimUntrimResultStruct *Untrim)


Description:

Deallocates a list of untrimming result structures.

Parameters:

Untrim: The list of untrimming result structures.


Returned Value:

void


Keywords:

untrimming


TrimUntrimmingResultToObj

(untrim.c:1112)

Prototype:

  IPObjectStruct *TrimUntrimmingResultToObj(
                                        const TrimUntrimResultStruct *Untrimmed)


Description:

Converts a list of untrimming result structures to Irit objects. Each untrimming result is converted to a list with the following structure: if the untrimming result has untrimmed Euclidean-space surfaces, they are gathered into a list. If the untrimming result has only UV-space data (surface patches or boundary curve pairs), then the first item in the list will be the containing surface, and the second will be a list of UV-space untrimmed objects. All the untrimmed surfaces are themselves gathered into a single list, which is returned by this function.

Parameters:

Untrimmed: A list of untrimming result structure.


Returned Value:

IPObjectStruct *: An object containing the untrimming results (as described above).


Keywords:

untrimming


TrimValidateNewTrimCntrs

(trimcntr.c:333)

Prototype:

  IPPolygonStruct *TrimValidateNewTrimCntrs(const CagdSrfStruct *Srf,
                                            const IPPolygonStruct *Cntrs)


Description:

Make sure the given trimming contours are in the surface domain. Points on contours that are found outside are coerced to be inside.

Parameters:

Srf: To trim into pieces.
Cntrs: Polylines to coerce to be inside.


Returned Value:

IPPolygonStruct *: New set of polylines that is guaranteed to be in Srf.


See Also:

TrimSrfsFromContours TrimSrfsFromContours2

Keywords: