G+Smo
24.08.0
Geometry + Simulation Modules
|
A B-spline function of one argument, with arbitrary target dimension.
This is the geometry type associated with gsBSplineBasis.
T | coefficient type |
Public Types | |
typedef memory::shared_ptr < gsBSpline > | Ptr |
Shared pointer for gsBSpline. | |
typedef memory::unique_ptr < gsBSpline > | uPtr |
Unique pointer for gsBSpline. | |
Public Member Functions | |
gsMatrix< index_t > | active (const gsMatrix< T > &u) const |
Returns the indices of active (nonzero) functions at points u, as a list of indices. More... | |
void | active_into (const gsMatrix< T > &u, gsMatrix< index_t > &result) const |
Indices of active (non-zero) function(s) for each point. More... | |
const gsBasis< T > & | basis (const index_t k) const |
Helper which casts and returns the k-th piece of this function set as a gsBasis. | |
uPtr | clone () |
Clone methode. Produceds a deep copy inside a uPtr. | |
T | closestPointTo (const gsVector< T > &pt, gsVector< T > &result, const T accuracy=1e-6, const bool useInitialPoint=false) const |
virtual void | computeMap (gsMapData< T > &InOut) const |
Computes map function data. More... | |
bool | contains (gsMatrix< T > const &p, T const &tol=1e-6) |
Returns true iff the point p is contained (approximately) on the curve, with the given tolerance. | |
short_t | degree (short_t i=0) const |
Returns the degree of the B-spline. | |
void | degreeElevate (short_t const i=1, short_t const dir=-1) |
Elevate the degree by the given amount i for the direction dir. If dir is -1 then degree elevation is done for all directions. Uses gsBasis<T>::degreeElevate. | |
gsMatrix< T > | deriv (const gsMatrix< T > &u) const |
Evaluate the derivatives,. More... | |
gsMatrix< T > | deriv2 (const gsMatrix< T > &u) const |
Evaluates the second derivatives of active (i.e., non-zero) functions at points u. More... | |
T | directedHausdorffDistance (const gsGeometry &other, const index_t nsamples=1000, const T accuracy=1e-6) const |
virtual T | distanceL2 (gsFunction< T > const &) const |
Computes the L2-distance between this function and the field and a function func. | |
T | domainEnd () const |
Returns the end value of the domain of the basis. | |
T | domainStart () const |
Returns the starting value of the domain of the basis. | |
gsMatrix< T > | eval (const gsMatrix< T > &u) const |
Evaluate the function,. More... | |
std::vector< gsMatrix< T > > | evalAllDers (const gsMatrix< T > &u, int n, bool sameElement=false) const |
Evaluate all derivatives upto order n,. More... | |
void | findCorner (const gsMatrix< T > &v, gsVector< index_t, 1 > &curr, T tol=1e-3) |
returns the tensor-index curr of the corner control point v, or an invalid index if the corner is not found within the tolerance tol | |
const gsFunction< T > & | function (const index_t k) const |
Helper which casts and returns the k-th piece of this function set as a gsFunction. | |
gsGeometrySlice< T > | getIsoParametricSlice (index_t dir_fixed, T par) const |
gsBSpline () | |
Default empty constructor. | |
gsBSpline (const Basis &basis, gsMatrix< T > coefs) | |
Construct B-Spline by basis and coefficient matrix. | |
gsBSpline (KnotVectorType KV, gsMatrix< T > coefs, bool periodic=false) | |
Construct B-Spline by a knot vector and coefficient matrix. | |
gsBSpline (T u0, T u1, unsigned interior, int degree, gsMatrix< T > coefs, unsigned mult_interior=1, bool periodic=false) | |
Construct a B-spline from an interval and knot vector specification. More... | |
T | HausdorffDistance (const gsGeometry &other, const index_t nsamples=1000, const T accuracy=1e-6, const bool directed=false) const |
Computes the Hausdorff distance between *this to other. | |
size_t | id () const |
Returns the patch index for this patch. | |
void | insertKnot (T knot, index_t i=1) |
template<class It > | |
void | insertKnots (It inBegin, It inEnd) |
std::vector < internal::gsCurveIntersectionResult < T > > | intersect (const gsBSpline< T > &other, T tolerance=1e-5, T curvatureTolerance=1+1e-6) const |
virtual void | invertPoints (const gsMatrix< T > &points, gsMatrix< T > &result, const T accuracy=1e-6, const bool useInitialPoint=false) const |
bool | isClosed (T tol=1e-10) const |
Returns true if this curve is closed. | |
bool | isOn (gsMatrix< T > const &u, T tol=1e-3) const |
Return true if point u is on the curve with tolerance tol. | |
bool | isPatchCorner (gsMatrix< T > const &v, T tol=1e-3) const |
Return true if point u is an endpoint (corner) of the curve with tolerance tol. | |
KnotVectorType & | knots (const int i=0) |
Returns a reference to the knot vector. | |
const KnotVectorType & | knots (const int i=0) const |
Returns a (const )reference to the knot vector. | |
void | merge (gsGeometry< T > *otherG) |
Merge other B-spline into this one. | |
int | newtonRaphson (const gsVector< T > &value, gsVector< T > &arg, bool withSupport=true, const T accuracy=1e-6, int max_loop=100, T damping_factor=1) const |
virtual index_t | nPieces () const |
Number of pieces in the domain of definition. | |
int | orientation () const |
Evaluates if the geometry orientation coincide with the ambient orientation. This is computed in the center of the parametrization and will fail to be meaningful if the geometry is singular. returns one if codimension is not zero. | |
virtual gsMatrix< T > | parameterCenter () const |
Returns a "central" point inside inside the parameter domain. | |
gsMatrix< T > | parameterCenter (const boxCorner &bc) const |
Get coordinates of the boxCorner bc in the parameter domain. | |
gsMatrix< T > | parameterCenter (const boxSide &bs) const |
Get coordinates of the midpoint of the boxSide bs in the parameter domain. | |
virtual const gsFunction & | piece (const index_t k) const |
Returns the piece(s) of the function(s) at subdomain k. | |
std::ostream & | print (std::ostream &os) const |
Prints the object as a string. | |
T | pseudoCurvature () const |
void | recoverPoints (gsMatrix< T > &xyz, gsMatrix< T > &uv, index_t k, const T accuracy=1e-6) const |
gsMatrix< T > | sample (int npoints=50) const |
Sample npoints uniformly distributed (in parameter domain) points on the curve. | |
gsBSpline< T > | segmentFromTo (T u0, T u1, T tolerance=1e-15) const |
void | setFurthestCorner (gsMatrix< T > const &v) |
Modifies the parameterization such that the point v is the ending point of the curve. Assumes that v is either the starting or the ending point of the curve. | |
void | setId (const size_t i) |
Sets the patch index for this patch. | |
void | setOriginCorner (gsMatrix< T > const &v) |
Modifies the parameterization such that the point v is the starting point of the curve. Assumes that v is either the starting or the ending point of the curve. | |
void | setPeriodic (bool flag=true) |
Tries to convert the curve into periodic. | |
index_t | size () const |
size More... | |
void | splitAt (T u0, gsBSpline< T > &left, gsBSpline< T > &right, T tolerance=1e-15) const |
gsMultiPatch< T > | toBezier (T tolerance=1e-15) const |
Evaluation functions | |
re-implemented from gsFunction, see also Evaluation members | |
void | eval_into (const gsMatrix< T > &u, gsMatrix< T > &result) const |
Evaluate the function at points u into result. More... | |
virtual void | deriv_into (const gsMatrix< T > &u, gsMatrix< T > &result) const |
Evaluate derivatives of the function \(f:\mathbb{R}^d\rightarrow\mathbb{R}^n\) at points u into result. More... | |
virtual void | deriv2_into (const gsMatrix< T > &u, gsMatrix< T > &result) const |
Evaluate second derivatives of the function at points u into result. More... | |
virtual void | evalAllDers_into (const gsMatrix< T > &u, int n, std::vector< gsMatrix< T > > &result, bool sameElement=false) const |
Evaluate the nonzero functions and their derivatives up to order n at points u into result. More... | |
virtual void | compute (const gsMatrix< T > &in, gsFuncData< T > &out) const |
Computes function data. More... | |
Accessors | |
virtual const gsBasis< T > & | basis () const =0 |
Returns a const reference to the basis of the geometry. More... | |
virtual gsBasis< T > & | basis ()=0 |
Returns a reference to the basis of the geometry. More... | |
short_t | targetDim () const |
Dimension of the ambient physical space (overriding gsFunction::targetDim()) | |
short_t | coefDim () const |
Dimension n of the coefficients (control points) | |
short_t | geoDim () const |
Dimension n of the absent physical space. | |
virtual short_t | domainDim () const |
Dimension d of the parameter domain (overriding gsFunction::domainDim()). | |
short_t | parDim () const |
Dimension d of the parameter domain (same as domainDim()). | |
short_t | coDim () const |
Co-dimension of the geometric object. | |
gsMatrix< T > | support () const |
Returns the range of parameters (same as parameterRange()) | |
gsMatrix< T > | parameterRange () const |
Returns the range of parameters as a matrix with two columns, [lower upper]. | |
std::vector< boxSide > | locateOn (const gsMatrix< T > &u, gsVector< bool > &onG2, gsMatrix< T > &preIm, bool lookForBoundary=false, real_t tol=1.e-6) const |
Get back the side of point u. More... | |
Coefficient access functions | |
These functions allow direct access to the coefficient matrix of the geometry. | |
gsMatrix< T > & | coefs () |
const gsMatrix< T > & | coefs () const |
Returns the coefficient matrix of the geometry. | |
gsMatrix< T >::RowXpr | coef (index_t i) |
Returns the i-th coefficient of the geometry as a row expression. | |
gsMatrix< T >::ConstRowXpr | coef (index_t i) const |
Returns the i-th coefficient of the geometry as a row expression. | |
T & | coef (index_t i, index_t j) |
Returns the j-th coordinate of the i-th coefficient of the geometry. | |
const T | coef (index_t i, index_t j) const |
Returns the j-th coordinate of the i-th coefficient of the geometry. | |
void | setCoefs (gsMatrix< T > cc) |
Set the coefficient matrix of the geometry, taking ownership of the matrix. | |
index_t | coefsSize () const |
Return the number of coefficients (control points) | |
Transformation functions | |
These functions apply various linear and affine transformations to the coefficients. | |
void | linearTransform (const gsMatrix< T > &mat) |
Apply the given square matrix to every control point. | |
void | rotate (T angle, const gsVector< T, 3 > &axis) |
Apply 3D Rotation by angle radians around axis axis. | |
void | rotate (T angle) |
Apply 2D Rotation by angle radians. | |
void | scale (T s, int coord=-1) |
Apply Scaling by factor s. | |
void | scale (gsVector< T > const &v) |
Apply Scaling coord-wise by a vector v. | |
void | translate (gsVector< T > const &v) |
Apply translation by vector v. | |
gsMatrix< T >::RowXpr | coefAtCorner (boxCorner const &c) |
Returns the control point at corner c. | |
gsMatrix< T >::ConstRowXpr | coefAtCorner (boxCorner const &c) const |
Returns the control point at corner c. | |
Other miscellaneous functions | |
virtual void | uniformRefine (int numKnots=1, int mul=1, int dir=-1) |
Refine the geometry uniformly, inserting numKnots new knots into each knot span. | |
virtual void | uniformCoarsen (int numKnots=1) |
Coarsen the geometry uniformly, removing numKnots new knots into each knot span. | |
void | refineElements (std::vector< index_t > const &boxes) |
Refines the basis and adjusts the coefficients to keep the geometry the same. More... | |
void | unrefineElements (std::vector< index_t > const &boxes) |
gsGeometry::uPtr | coord (const index_t c) const |
void | embed3d () |
Embeds coefficients in 3D. | |
void | embed (index_t N, bool pad_right=true) |
Embeds coefficients in N dimensions For the new dimensions zeros are added (or removed) on the right (if pad_right is true) or on the left (if pad_right is false) | |
short_t | degree (const short_t &i) const |
Returns the degree wrt direction i. | |
virtual void | degreeIncrease (short_t const i=1, short_t const dir=-1) |
Elevate the degree by the given amount i for the direction dir. If dir is -1 then degree elevation is done for all directions. Uses gsBasis<T>::degreeIncrease. | |
virtual void | degreeReduce (short_t const i=1, short_t const dir=-1) |
Reduces the degree by the given amount i for the direction dir. If dir is -1 then degree reduction is done for all directions. Uses gsBasis<T>::degreeReduce. | |
virtual void | degreeDecrease (short_t const i=1, short_t const dir=-1) |
Reduces the degree by the given amount i for the direction dir. If dir is -1 then degree reduction is done for all directions. Uses gsBasis<T>::degreeDecrease. | |
virtual void | hessian_into (const gsMatrix< T > &u, gsMatrix< T > &result, index_t coord) const |
void | controlNet (gsMesh< T > &mesh) const |
Return the control net of the geometry. | |
void | outerNormal_into (const gsMatrix< T > &u, gsMatrix< T > &result) const |
Computes the outer normals at parametric points u. More... | |
std::vector< gsGeometry * > | boundary () const |
Get boundary of this geometry as a vector of new gsGeometry instances. | |
virtual gsGeometry::uPtr | boundary (boxSide const &s) const |
Get parametrization of boundary side s as a new gsGeometry uPtr. | |
virtual gsGeometry::uPtr | iface (const boundaryInterface &bi, const gsGeometry &other) const |
Computes and returns the interface with other as a new geometry. | |
gsGeometry::uPtr | component (boxComponent const &bc) const |
Get parametrization of box component bc as a new gsGeometry uPtr. | |
virtual void | merge (gsGeometry *other) |
Merge the given other geometry into this one. | |
virtual void | toMesh (gsMesh< T > &msh, int npoints) const |
gsGeometry::uPtr | approximateLinearly () const |
void | evaluateMesh (gsMesh< T > &mesh) const |
virtual std::vector < gsGeometry< T > * > | uniformSplit (index_t dir=-1) const |
Evaluation functions | |
These functions allow one to evaluate the function as well as its derivatives at one or more points in the parameter space. See also Evaluation members. | |
virtual void | eval_component_into (const gsMatrix< T > &u, const index_t comp, gsMatrix< T > &result) const |
Evaluate the function for component comp in the target dimension at points u into result. | |
virtual void | jacobian_into (const gsMatrix< T > &u, gsMatrix< T > &result) const |
Computes for each point u a block of result containing the Jacobian matrix. | |
void | div_into (const gsMatrix< T > &u, gsMatrix< T > &result) const |
Computes for each point u a block of result containing the divergence matrix. | |
gsMatrix< T > | jacobian (const gsMatrix< T > &u) const |
virtual gsMatrix< T > | hessian (const gsMatrix< T > &u, index_t coord=0) const |
virtual gsMatrix< T > | laplacian (const gsMatrix< T > &u) const |
Evaluate the Laplacian at points u. More... | |
Protected Attributes | |
gsBasis< T > * | m_basis |
Pointer to the basis of this geometry. | |
gsMatrix< T > | m_coefs |
Coefficient matrix of size coefsSize() x geoDim() | |
size_t | m_id |
Private Member Functions | |
virtual void | insertKnot (T knot, index_t dir, index_t i=1) |
Inserts knot knot at direction dir, i times. | |
Returns the indices of active (nonzero) functions at points u, as a list of indices.
|
inlinevirtualinherited |
Indices of active (non-zero) function(s) for each point.
The columns are sorted in increasing order, if on a point there are less active then the number of rows in the result matrix (some other point has more actives) then the rest of the column is filled with 0s.
u | |
result |
Reimplemented from gsFunctionSet< T >.
|
pure virtualinherited |
Returns a const reference to the basis of the geometry.
Implemented in gsConstantFunction< T >.
|
pure virtualinherited |
Returns a reference to the basis of the geometry.
Implemented in gsConstantFunction< T >.
|
inherited |
Returns the parameters of closest point to pt as an argument, and the Euclidean distance as a return value
|
inlineinherited |
Returns the coefficient matrix of the geometry Coefficient matrix of size coefsSize() x geoDim()
|
virtualinherited |
Computes function data.
This function evaluates the functions and their derivatives at the points in and writes them in the corresponding fields of out. Which field to write (and what to compute) is controlled by the out.flags (see also gsFuncData).
The input points in are expected to be compatible with the implementation/representation of the function, i.e. they should be points inside the domain of definitition of the function
[in] | in | |
[out] | out |
Reimplemented from gsFunctionSet< T >.
Reimplemented in gsConstantFunction< T >.
|
virtualinherited |
Computes map function data.
This function evaluates the functions and their derivatives at the points InOut.points and writes them in the corresponding fields of InOut. Which field to write (and what to compute) is controlled by the InOut.flags (see also gsMapData). This is intended for parametrizations only and it works on functions sets of cardinality 1 only.
[in,out] | InOut |
Evaluate the derivatives,.
Evaluates the second derivatives of active (i.e., non-zero) functions at points u.
See documentation for deriv2_into() (the one without input parameter coefs) for details.
[in] | u | Evaluation points in columns. |
Evaluate second derivatives of the function at points u into result.
Let d be the dimension of the source space ( d = domainDim() ).
Let n be the dimension of the image/target space ( n = targetDim() ).
Let N denote the number of evaluation points.
[in] | u | gsMatrix of size d x N, where each column of u represents one evaluation point. |
[out] | result | gsMatrix of size (S*n) x N, where S=d*(d+1)/2. Each column in result corresponds to one point (i.e., one column in u) and contains the following values (for d=3, n=3): \( (\partial_{xx} f^{(1)}, \partial_{yy} f^{(1)}, \partial_{zz} f^{(1)}, \partial_{xy} f^{(1)}, \partial_{xz} f^{(1)}, \partial_{yz} f^{(1)}, \partial_{xx} f^{(2)},\ldots,\partial_{yz} f^{(3)} )^T\) |
Reimplemented from gsFunction< T >.
Reimplemented in gsConstantFunction< T >.
Evaluate derivatives of the function \(f:\mathbb{R}^d\rightarrow\mathbb{R}^n\) at points u into result.
Let d be the dimension of the source space ( d = domainDim() ).
Let n be the dimension of the image/target space ( n = targetDim() ).
Let N denote the number of evaluation points.
Let \( f:\mathbb R^2 \rightarrow \mathbb R^3 \), i.e., \( f(x,y) = ( f^{(1)}(x,y), f^{(2)}(x,y), f^{(3)}(x,y) )^T\),
and let \( u = ( u_1, \ldots, u_N) = ( (x_1,y_1)^T, \ldots, (x_N, y_N)^T )\).
Then, result is of the form
\[ \left[ \begin{array}{cccc} \partial_x f^{(1)}(u_1) & \partial_x f^{(1)}(u_2) & \ldots & \partial_x f^{(1)}(u_N) \\ \partial_y f^{(1)}(u_1) & \partial_y f^{(1)}(u_2) & \ldots & \partial_y f^{(1)}(u_N) \\ \partial_x f^{(2)}(u_1) & \partial_x f^{(2)}(u_2) & \ldots & \partial_x f^{(2)}(u_N) \\ \partial_y f^{(2)}(u_1) & \partial_y f^{(2)}(u_2) & \ldots & \partial_x f^{(2)}(u_N) \\ \partial_x f^{(3)}(u_1) & \partial_x f^{(3)}(u_2) & \ldots & \partial_x f^{(3)}(u_N)\\ \partial_y f^{(3)}(u_1) & \partial_y f^{(3)}(u_2) & \ldots & \partial_y f^{(3)}(u_N) \end{array} \right] \]
[in] | u | gsMatrix of size d x N, where each column of u represents one evaluation point. |
[out] | result | gsMatrix of size (d * n) x N. Each row of result corresponds to one component in the target space and contains the gradients for each evaluation point, as row vectors, one after the other (see above for details on the format). |
Reimplemented from gsFunction< T >.
Reimplemented in gsConstantFunction< T >.
|
inherited |
Computes the Hausdorff distance in a single direction from *this to other. The Hausdorff distance is computed by taking the maximum of the shortest distances between points of this and other.
Evaluate the function,.
Evaluate the function at points u into result.
Let d be the dimension of the source space ( d = domainDim() ).
Let n be the dimension of the image/target space ( n = targetDim() ).
Let N denote the number of evaluation points.
[in] | u | gsMatrix of size d x N, where each column of u represents one evaluation point. |
[out] | result | gsMatrix of size n x N, where each column of u represents the result of the function at the respective valuation point. |
Implements gsFunction< T >.
|
inherited |
Evaluate all derivatives upto order n,.
|
virtualinherited |
Evaluate the nonzero functions and their derivatives up to order n at points u into result.
The derivatives (the 0-th derivative is the function value) are stored in a result. result is a std::vector, where result[i] is a gsMatrix which contains the i-th derivatives.
The entries in result[0], result[1], and result[2] are ordered as in eval_into(), deriv_into(), and deriv2_into(), respectively. For i > 2, the derivatives are stored in lexicographical order, e.g. for order i = 3 and dimension 2 the derivatives are stored as follows: \( \partial_{xxx}, \, \partial_{xxy}, \, \partial_{xyy}, \, \partial_{yyy}.\, \)
[in] | u | Evaluation points, each column corresponds to one evaluation point. |
[in] | n | All derivatives up to order n are computed and stored in result. |
[in,out] | result | See above for format. |
Reimplemented from gsFunctionSet< T >.
Reimplemented in gsConstantFunction< T >.
|
inherited |
Updates the vertices of input mesh by evaluating the geometry at vertices. Vertices of the new mesh are
{ geom(v) | v vertex of input mesh }
|
inherited |
Gives back an isoParametric slice of the geometry with fixed par in direction dim_fixed as an gsGeometrySlice object.
|
inlinevirtualinherited |
Evaluates the Hessian (matrix of second partial derivatives) of coordinate coord at points u.
|
virtualinherited |
Compute the Hessian matrix of the coordinate coord evaluated at points u
Reimplemented from gsFunction< T >.
void insertKnot | ( | T | knot, |
index_t | i = 1 |
||
) |
Insert the given new knot (multiplicity i) without changing the curve.
|
inline |
Insert the given new knots in the range [inBegin .. inEend ) without changing the curve.
inBegin | iterator pointing to the first knot to be inserted |
inEnd | iterator pointing to one position after the last knot to be inserted |
std::vector< internal::gsCurveIntersectionResult< T > > intersect | ( | const gsBSpline< T > & | other, |
T | tolerance = 1e-5 , |
||
T | curvatureTolerance = 1+1e-6 |
||
) | const |
Compute all intersections with another B-Spline curve
other | the other B-Spline curve |
tolerance | the tolerance for computing the intersection |
curvatureTolerance | the curvature tolerance for the splitting stage |
|
virtualinherited |
Takes the physical points and computes the corresponding parameter values. If the point cannot be inverted (eg. is not part of the geometry) the corresponding parameter values will be undefined
Evaluate the Laplacian at points u.
By default uses central finite differences with h=0.00001
Reimplemented in gsFunctionExpr< T >.
|
inherited |
Get back the side of point u.
Check if points u also lie on the geometry and if required computes the if the points in u lie on one of the boundaries of the geometry
u | Matrix of points of the form geoDim() x #points |
onG2 | gsVector of booleans which indicate if the point is in the domain or outside |
preIm | Matrix of preimages of the points u |
lookForBoundary | if required the boundaries are computed the points in u lie on |
|
inherited |
Newton-Raphson method to find a solution of the equation f(arg) = value with starting vector arg. If the point cannot be inverted the corresponding parameter values will be undefined
Computes the outer normals at parametric points u.
Assumes that u is a list of points on the boundary of the geometry.
T pseudoCurvature | ( | ) | const |
calculates curvature approximation for a B-spline curve by dividing the total length of the polyline formed by its control points by the direct distance between the first and last control point
|
inherited |
Recovers a point on the (geometry) together with its parameters uv, assuming that the k-th coordinate of the point xyz is not known (and has a random value as input argument).
|
inherited |
Refines the basis and adjusts the coefficients to keep the geometry the same.
The syntax of boxes depends on the implementation in the underlying basis. See gsBasis::refineElements_withCoefs() for details.
gsBSpline< T > segmentFromTo | ( | T | u0, |
T | u1, | ||
T | tolerance = 1e-15 |
||
) | const |
Segment this BSpline curve between u0 and u1. Either of these values can be outside the parameter range of the curve, but u1 must be greater than u0. The degree of the curve is not modified.
u0 | starting parameter |
u1 | end parameter |
tolerance | proximity of the segment boundaries and B-spline knots to treat them as equal |
|
inlinevirtualinherited |
size
Reimplemented from gsFunctionSet< T >.
Reimplemented in gsPiecewiseFunction< T >.
Split this BSpline curve at u0. This value must be inside the parameter range of the curve.
u0 | split parameter |
left | curve segment from the left of parameter domain to u0 |
right | curve segment from u0 to the right of parameter domain |
tolerance | proximity of the segment boundaries and B-spline knots to treat them as equal |
gsMultiPatch< T > toBezier | ( | T | tolerance = 1e-15 | ) | const |
Convert BSpline curve into Bezier segments
tolerance | proximity of the segment boundaries and B-spline knots to treat them as equal |
|
virtualinherited |
Splits the geometry 2^d parts, where each direction is divided into two parts in in a uniform way, i.e., in the middle of the corresponding side. This method allocated new space for each new geometry, the original one stays unchanged.
Reimplemented in gsTensorNurbs< d, T >, gsTensorBSpline< d, T >, and gsTensorBSpline< domainDim+1, T >.
|
protectedinherited |
An auxiliary index for this geometry (eg. in case it is part of a multi-patch object)