API Reference

Core.UnionMethod
(basis::Bspline)(evalpoints)

Evaluate B-spline basis at all evalpoints.

source
Core.UnionMethod
(basis::Bspline)(evalpoints)

Evaluate B-spline basis at all evalpoints.

source
NURBS.BsplineMethod
(basis::Bspline)(evalpoints, k::Int, prealloc)

Evaluate k-the derivative of B-spline basis at all evalpoints (all basis functions different from 0 at the evalpoints are evaluated).

source
NURBS.BsplineMethod
(basis::Bspline)(evalpoints, k::Int)

Evaluate k-the derivative of B-spline basis at all evalpoints (all basis functions different from 0 at the evalpoints are evaluated).

source
NURBS.BsplineCurveMethod
(curve::BsplineCurve)(uVector, k::Int)

Convenience function to compute points on all k derivatives of a B-spline curve.

source
NURBS.BsplineCurveMethod
(curve::BsplineCurve)(uVector)

Convenience function to compute points on a B-spline curve.

source
NURBS.BsplineSurfaceMethod
(Patch::BsplineSurface)(uEvalpoints, vEvalpoints, k::Int)

Convenience function to compute points on all k derivatives of a B-spline surface.

source
NURBS.BsplineSurfaceMethod
(Patch::BsplineSurface)(uEvalpoints, vEvalpoints)

Convenience function to compute points on a B-spline surface.

source
NURBS.ConnectivityMethod
Connectivity(patches::Vector{<:Surface}, cU, cV; tol=1e-3)

Constructor for Connectivity structure.

source
NURBS.InterfaceType
Interface

An interface is defined by two patches: their ID and wether the vectors normal to the interface point in the same direction.

source
NURBS.InterfacePatchwiseType
InterfacePatchwise

Save connectivity information for a single patch:

Each entry is a vector, e.g., localEdgesVec[i] contains the local edges being part of an interface of patch i.

source
NURBS.NURBScurveType
NURBScurve{F} <: Curve{F}

B-spline curve defined by the basis and the control points.

source
NURBS.NURBScurveMethod
(curve::NURBScurve)(uVector, k::Int)

Convenience function to compute points on all k derivatives of a NURBS curve.

source
NURBS.NURBScurveMethod
(curve::NURBScurve)(uVector)

Convenience function to compute points on a NURBS curve.

source
NURBS.NURBSsurfaceType
NURBSsurface{F} <: Surface{F}

Surface defined by a B-spline basis, the control points, and the weights.

source
NURBS.NURBSsurfaceMethod
(Patch::NURBSsurface)(uEvalpoints, vEvalpoints, k::Int, prealloc)

Convenience function to compute points on all k derivatives of a NURBSsurface, for preallocated memory.

source
NURBS.NURBSsurfaceMethod
(Patch::NURBSsurface)(uEvalpoints, vEvalpoints, k::Int)

Convenience function to compute points on all k derivatives of a NURBSsurface.

source
NURBS.NURBSsurfaceMethod
(Patch::NURBSsurface)(uEvalpoints, vEvalpoints)

Convenience function to compute points on a NURBSsurface.

source
NURBS.PatchInterfaceType

Local edge and vertex numbering

     2   2

1 o––––o––-> v-direction | | 1| |3 | | 4 o––––o 3 | 4 v u-direction

An interface of a patch.

source
Base.splitFunction
Base.split(C::Curve, splitPoint=0.5)

Split a curve at a single parametric point in the range ]0,1[ by inserting a single knot multiple times.

source
Base.splitMethod
Base.split(C::Curve, n::Int)

Split a curve in the parametric domain equally to obtain n curves.

source
Base.splitMethod
Base.split(C::Curve, n::Int)

Split a curve in the parametric domain equally to obtain n curves.

source
Base.splitMethod
Base.split(S::SurfaceT; U=[], V=[]) where {SurfaceT<:Surface}

Split a surface at multiple parametric points in U and V, each in the range ]0,1[, by inserting knots.

source
Base.splitMethod
Base.split(C::CurveT, splits::Vector) where {CurveT<:Curve}

Split a curve at multiple parametric points, each in the range ]0,1[, by inserting knots.

The splits vector needs to be sorted and each value should only occur once.

source
LinearAlgebra.rotate!Method
rotate!(shape::S, rotAxis::SVector{3,T}, angle::Real) where {S<:Shape,T}

Rotate a shape around the rotation axis by an angle (in rad).

source
LinearAlgebra.rotate!Method
rotate!(shapes::Vector{S}, rotAxis::SVector{3,T}, angle::Real) where {S<:Shape,T}

Rotate a vector of shapes around the rotation axis by an angle (in rad).

source
NURBS.JacobianMethod
Jacobian(Patch::Surface, uEvalpoints, vEvalpoints)

Compute the Jacobian matrix and its (generalized) determinant at the parametric points 'uEvalpoints' and 'vEvalpoints'.

Return - J 2-dimensional vector: first for the derivative w.r.t 'u', second w.r.t 'v' each vector entry contains a matrix of size (uEvalpoints, vEvalpoints) each entry of the matrix is an SVector with the derivatives: SVector(∂x/∂u, ∂y/∂u, ∂y/∂u)

- dJ    matrix of size (uEvalpoints, vEvalpoints) where each entry is the Jacobi determinant evaluated at the points 'u' and 'v'.

Note: surface points are evaluated but thrown away: maybe change this/make use of it.

source
NURBS.anchorsMethod
anchors(kVec, degree::Int)

Return the anchors (as defined in [4]) corresponding to the given knotvector and the polynomial degree.

source
NURBS.bSplineNaiveMethod
bSplineNaive(knotVector, i::Int, degree::Int, evalpoints)

i-th b-spline basis function of degree 'degree' evaluated at all 'evalpoints'.

The knotvector is assumed to be normalized to [1, 0].

source
NURBS.bSplineNaiveMethod
bSplineNaive(knotVector, i::Int, degree::Int, u::Real, N)

i-th b-spline basis function of degree 'degree' evaluated at u.

Formula (2.5) of 'The NURBS Book' p. 50.

source
NURBS.bSplineNaiveMethod
bSplineNaive(knotVector, i::Int, u::Real)

i-th b-spline basis function of degree 0 evaluated at u.

Formula (2.5) of 'The NURBS Book' p. 50.

source
NURBS.bSplineNaiveDerivativeMethod
bSplineNaiveDer(knotVector, i::Int, degree::Int, evalpoints, k::Int; normalize=true)

Compute the k-th derivative of i-th b-spline basis function of degree 'degree' evaluated at all 'evalpoints'.

source
NURBS.bSplineNaiveDerivativeMethod
bSplineNaiveDerivative(knotVector, i::Int, degree::Int, u::Real, k::Int)

k-th derivative of the i-th b-spline basis function of degree 'degree' evaluated at (single) point 'u'.

Formula (2.9) of 'The NURBS Book' p. 61. (Recursive implementation to avoid the faculties in the (2.10) formula.)

source
NURBS.basisFun!Method
basisFun!(prealloc::pAlloc, knotSpan, uVector, basis::Basis)

Compute the nonvanishing B-spline basis functions of degree 'degree' at the parametric points defined by 'uVector'

Return the B-spline basis functions vector of size length(uVector) * (degree + 1).

Adapted from Algorithm A2.2 from 'The NURBS Book' p. 70.

source
NURBS.basisFunMethod
basisFun(knotSpan, uVector, basis::Basis)

Allocate memory and call basisFun!

source
NURBS.coarsenFunction
coarsen(knotVec, controlPoints, degree::Int, pointsToRemove::Vector, weights=[])

Coarsen a curve by removing parametric points from a curve's knot vector.

NOTE: There are more efficient ways to do this.

source
NURBS.curveDerivativesPointsMethod
curveDerivativesPoints(nbasisFun::Int, degree::Int, knotVector, controlPoints, uVector, weights, k::Int)

Compute points on the k-th derivatives of a NURBS curve: given the 'knotVector', the 'controlPoints', the 'degree', and the 'weights', the curve is evaluated at the points given in 'uVector'.

Using (4.8) on page 125 of 'The NURBS Book'.

Example for the controlPoints:

P1 = SVector(0.0, 0.0, 0.0) P2 = SVector(0.1, 0.25, 0.0) P3 = SVector(0.25, 0.3, 0.0)

controlPoints = [P1, P2, P3]

Note: the efficient evaluation via the B-spline basis is employed (no use of the naive evaluation of the NURBS basis).

source
NURBS.curveDerivativesPointsMethod
curvePoints(nbasisFun::Int, degree::Int, knotVector, controlPoints, uVector)

Compute a points on the k-th derivatives of a B-spline curve: given the 'knotVector', the 'controlPoints', and the 'degree', the curve is evaluated at the points given in 'uVector'.

Example for the controlPoints:

P1 = SVector(0.0, 0.0, 0.0) P2 = SVector(0.1, 0.25, 0.0) P3 = SVector(0.25, 0.3, 0.0)

controlPoints = [P1, P2, P3]

source
NURBS.curvePointsMethod
curvePoints(nbasisFun::Int, degree::Int, knotVector, controlPoints, uVector, weights)

Compute points on a NURBS curve: given the 'knotVector', the 'controlPoints', the 'degree', and the 'weights', the curve is evaluated at the points given in 'uVector'.

Example for the controlPoints:

P1 = SVector(0.0, 0.0, 0.0) P2 = SVector(0.1, 0.25, 0.0) P3 = SVector(0.25, 0.3, 0.0)

controlPoints = [P1, P2, P3]

Note: the efficient evaluation via the B-spline basis is employed (no use of the naive evaluation of the NURBS basis).

source
NURBS.curvePointsMethod
curvePoints(basis::Basis, controlPoints, uVector)

Compute a 1D B-spline curve: given the 'knotVector', the 'controlPoints', and the 'degree', the curve is evaluated at the points given in 'uVector'.

Example for the controlPoints:

P1 = SVector(0.0, 0.0, 0.0) P2 = SVector(0.1, 0.25, 0.0) P3 = SVector(0.25, 0.3, 0.0)

controlPoints = [P1, P2, P3]

source
NURBS.derBasisFun!Method
derBasisFun!(knotSpan, degree::Int, evalpoints, knotVector, numberDerivatives::Int)

Compute the nonvanishing B-spline basis functions and its derivatives of degree 'degree' at the parametric points defined by 'uVector'.

Organization of output: dersv[n, k, :] contains (k-1)-th derivative at n-th point.

Adapted from Algorithm A2.3 from 'The NURBS BOOK' p. 72.

source
NURBS.derBasisFunMethod
derBasisFun(degree::Int, evalpoints, knotVector, numberDerivatives::Int)

Same as derBasisFun! but with memory allocation.

source
NURBS.evalNaiveMethod
evalNaive(basis::NURB, i::Int, evalpoints)

i-th NURB basis function evaluated at all 'evalpoints'.

source
NURBS.evalNaiveMethod
evalNaive(basis::Basis, i::Int, evalpoints)

i-th B-spline or Curry-Schoenberg basis function evaluated at all 'evalpoints'.

source
NURBS.evalNaiveDerivativeMethod
evalNaiveDerivative(basis::Bspline, i::Int, k::Int, evalpoints)

Compute the k-th derivative of i-th b-spline basis function evaluated at all 'evalpoints'.

source
NURBS.evalNaiveDerivativeMethod
evalNaiveDerivative(basis::NURB, i::Int, k::Int, evalpoints)

Compute the k-th derivative of i-th NURBS basis function evaluated at all 'evalpoints'.

source
NURBS.extendControlPoints!Function
extendControlPoints!(controlPoints, knotVecOrig, degree::Int, pos::Int, uNew::Real, multiplicity::Int, oldMult::Int, weights)

Insert the new control points (and optionally the weights) corresponding to the new values in the knot vector 'multiplicity' times.

Adaption of Algorithm A5.1 from 'The NURBS Book' p. 151.

Modifies controlPoints and weights.

source
NURBS.extendControlPoints!Method

Lowest order polynomial degree = 0.

Can multiplicity be maximum 1 and/or oldMult maximum 0? => simplify below algorithm further.

source
NURBS.extendKnotVector!Method
extendKnotVector!(knotVecOrig, newParametricPoint::Real, multiplicity::Int)

Insert the new value with a given mulitplicity (insert the point multiple times).

Modifies knotVecOrig.

source
NURBS.extractBracketListFunction
extractBracketList(line::String, nextInd::Int, T=Float64)

Given a string with a list of numbers between two brackets (e.g., (2,3,1,5,...,6)), extract the numbers.

As second output the index after ')' is returned.

'nextInd' is the index of the '(' in 'line'.

source
NURBS.extractCtrlPointsFunction
extractCtrlPoints(line::String, firstCommaInd::Int, points, offset::Int, T=Float64)

Extract a list of points between '((' and '))' in the line starting from 'firstCommaInd' of the format

    ((#349,#350,#351,#352,#353),(#354,#355,#356,#357,#358),(#359,#360,#361,#362,#363),(#364,#365,#366,#367,#368),(#369,#370,#371,#372,#373))
source
NURBS.extractKnotVecsFunction
extractKnotVecs(line::String, startInd::Int, T=Float64)

Extract the knot vectors from the 'line' of the format 'SOMESTRING(5,5),(5,5),(0.,1.),(0.,1.)SOMESTRING' where 'startInd' points to the 'G'.

source
NURBS.findSpan!Method
findSpan!(spanVec, b::T, uVec, kVec, p::T)

Find the spans of a B-spline knot vector at the parametric points 'u', where 'b' is the number of basis functions (control points).

Span: the intervall index in which a point lies. E.g., knotVector = [0, 1, 2, 3, 4]. Hence, there are 4 intervalls. u=1.2 lies in the second intervall.

Modification of Algorithm A2.1 from 'The NURBS Book' p. 68.

Assumption that the knotVector is open! (the first and last knot are repeated degree + 1 times)

source
NURBS.generateKnotVecMethod
generateKnotVec(b::Int, degree::Int)

Convenience function to generate a knot vector for 'b' basis functions and a certain 'degree':

The first and last entry are repeated 'degree'+1 times. Normalized to [0, 1].

source
NURBS.getCartesianPointsFunction
getCartesianPoints(stringVec)

Extract all Cartesian points and the offset between the first point and the #counting. (Assumption: all points appear in a consecutive list.)

source
NURBS.getCompleteStringMethod
getCompleteString(stringVec::String, nextLineInd::Int)

Retreive all lines of the string vector 'stringVec' starting from 'nextLineInd' till a semicolon is encountered.

source
NURBS.getCompleteStringTillMethod
getCompleteStringTill(stringVec, nextLineInd::Int, stopWhenEncountered::String)

Retreive all lines of the string vector 'stringVec' starting from 'nextLineInd' till the next line starts with 'stopWhenEncountered'.

source
NURBS.getGlobalIndexMethod
getGlobalIndex(P, Plist; tol=1e-3)

Get the global index of the point P with respect to the ordering in the list of points Plist.

source
NURBS.getPatchInterfacesMethod
getPatchInterfaces(patches::Vector{<:Surface{F}}, interfaces, commonVtxs) where {F}

Determine for each patch the local edges where there is an interface and whether the orientation has to be flipped.

localEdgesVec[p] contains the local edges of the p-th patch that are attached to an interface
source
NURBS.grevilleMethod
greville(kVec, degree::Int)

Return the Greville sites (as defined in [3]) corresponding to the given knotvector and the polynomial degree.

source
NURBS.identifyInterfacesMethod
identifyInterfaces(patches::Vector{<:Surface{F}}; tol=1e-3) where {F}

Determine the interfaces between patches including the property if the basis functions normal to the edge point in the same direction.

source
NURBS.insertKnotFunction
insertKnot(C::Curve, newParametricPoint::Real, multiplicity::Int)

Convenience function to insert a knot into a curve.

source
NURBS.insertKnotFunction
insertKnot(knotVecOrig, controlPointsOrig, p::Int, newParametricPoint::Real, multiplicity::Int)

Insert a new value with a given mulitplicity (insert the value multiple times) into a knot vector for the polynomial degree 'degree'. Return the resutling knot vector and the new control points.

Adaption of Algorithm A5.1 from 'The NURBS Book' p. 151.

The provided knot vector and control-point vector are NOT modified.

source
NURBS.insertKnot!Function
insertKnot!(knotVecOrig, controlPointsOrig, p::Int, newParametricPoint::Real, multiplicity::Int)

Insert a new value with a given mulitplicity (insert the value multiple times) into a knot vector for the polynomial degree 'degree'. Return the resutling knot vector and the new control points.

Adaption of Algorithm A5.1 from 'The NURBS Book' p. 151.

The provided knot vector and control-point vector are modified.

source
NURBS.insertKnotUFunction
insertKnotU(S::Surface, newParametricPoint::Real, multiplicity::Int=1)

Insert a knot in u-direction with the given multiplicity.

TODO: the computation of the alphas is redundant.

source
NURBS.insertKnotVFunction
insertKnotV(S::Surface, newParametricPoint::Real, multiplicity::Int=1)

Insert a knot in v-direction with the given multiplicity.

TODO: the computation of the alphas is redundant.

source
NURBS.isValidKnotVector!Method
isValidKnotVector!(kVec)

Check whether the knot vector has only entries in [0, 1] and is in ascending order.

If not the knot vector is modified.

source
NURBS.isinMethod
isin(P, Plist; tol=1e-3)

Check if 3D point P is in list of points Plist.

source
NURBS.loadFunction
load(f::File{format"DAT"}, T=Float64)

Load a step file via the FileIO package.

If the file has no ending the magic byte is used to identify the file format.

Supported endings: .stp, .step, .stpnc, .p21, .210

source
NURBS.loadFunction
load(f::File{format"DAT"}, T=Float64)

Load a .dat file via the FileIO package.

source
NURBS.mirror!Method
translate!(shape::S, shift::SVector{3,T}) where {S<:Shape, T}

Mirror a shape through a plane defined by its normal and an anchor point.

source
NURBS.mirror!Method
mirror!(shapes::Vector{S}, normal::SVector{3,T}, anchor::SVector{3,T}) where {S<:Shape,T}

Mirror a vector of shapes through a plane defined by its normal and an anchor point.

source
NURBS.mirrorMethod
mirror(shape, rotAxis, angle)

Mirror through a plane defined by its normal and an anchor point.

source
NURBS.normalizationMethod
normalization(basis::CurrySchoenberg, i::Int)

Normalization for the standard B-splines to obtain Curry-Schoenberg splines.

source
NURBS.normalizationDerivativeMethod
normalizationDerivative(knotVector, degree::Int, weights, u::Real, k::Int)

Return the k-th derivative of the weight function (21) in [2].

source
NURBS.normalize!Method
normalize!(N, i::T, degree::T, knotVector, basis::CurrySchoenberg)

Normalize the standard B-splines to obtain Curry-Schoenberg splines.

source
NURBS.normalize!Method
normalize!(N, i::T, degree::T, knotVector, basis::Bspline)

For Bsplines there is no normalization.

source
NURBS.numBasisFunctionsMethod
numBasisFunctions(basis::Basis)

The number of basis functions is fixed by the knot vector and the degree.

Assumption: the first and last knot vector entry has mulitplicity degree + 1.

source
NURBS.nurbsNaiveMethod
nurbsNaive(knotVector, i::Int, degree::Int, evalpoints, weights)

i-th NURB basis function of degree 'degree', the 'knotVector', and with 'weights' evaluated at all 'evalpoints'.

source
NURBS.nurbsNaiveMethod
nurbsNaive(knotVector, i::Int, degree::Int, u::Real, weights)

i-th NURB basis function of degree 'degree', the 'knotVector', and with 'weights' evaluated at single point 'u'.

Formula (4.2) of 'The NURBS Book' p. 118.

source
NURBS.nurbsNaiveDerivativeMethod
nurbsNaiveDerivative(knotVector, i::Int, degree::Int, weights, evalpoints, k::Int; normalize=true)

Compute the k-th derivative of i-th NURBS basis function of degree 'degree' evaluated at all 'evalpoints'.

source
NURBS.nurbsNaiveDerivativeMethod
nurbsNaiveDerivative(knotVector, i::Int, degree::Int, weights, u::Real, k::Int)

k-th derivative of the i-th NURBS basis function of degree 'degree' evaluated at (single) point 'u'.

Formula (52) of [2].

source
NURBS.parseCtrlPointsFunction
parseCtrlPoints(pointDims, stringVec, T=Float64)

Take the 3 strings in 'stringVec' and parse it into a pointDims[1] x pointDims[2] matrix where each entry is a SVector for a controlpoint.

source
NURBS.parseLineFunction
parseLine(line::String, T=Float64)

Take the string 'line' of the form '0.1 0.2 0.6 ...' and parse it to a vector with eltype 'T'.

source
NURBS.parseSinglePatchMethod
parseSinglePatch(stringVec, lineInd)

Take the string vector and the line index where "PATCH" stands and extract the information of a single patch

source
NURBS.parse_B_SPLINE_SURFACE_WITH_KNOTSFunction
parse_B_SPLINE_SURFACE_WITH_KNOTS(line::String, points, offset::Int, T=Float64)

Parse a BSPLINESURFACEWITHKNOTS

Example:

#13764=BSPLINESURFACEWITHKNOTS('',3,3,((#16247,#16248,#16249,#16250, #16251),(#16252,#16253,#16254,#16255,#16256),(#16257,#16258,#16259,#16260, #16261),(#16262,#16263,#16264,#16265,#16266),(#16267,#16268,#16269,#16270, #16271)),.UNSPECIFIED.,.F.,.F.,.F.,(4,1,4),(4,1,4),(0.,1.,2.),(0.,1.,2.), .UNSPECIFIED.);

source
NURBS.parse_NURBS_stepFunction

parse_NURBS(lineCtlPts::String, lineKnVecs::String, lineWeight::String, points, offset::Int, T=Float64)

Example:

#163=( BOUNDEDSURFACE() BSPLINESURFACE(4,4,((#349,#350,#351,#352,#353),(#354,#355,#356,#357,#358), (#359,#360,#361,#362,#363),(#364,#365,#366,#367,#368),(#369,#370,#371,#372, #373)),.UNSPECIFIED.,.F.,.F.,.F.) BSPLINESURFACEWITHKNOTS((5,5),(5,5),(0.,1.),(0.,1.),.UNSPECIFIED.) GEOMETRICREPRESENTATIONITEM() RATIONALBSPLINESURFACE(((5.07179676972449,4.52004210360334,4.35726558990816, 4.52004210360334,5.07179676972449),(4.52004210360334,3.86602540378444,3.64492370567392, 3.86602540378444,4.52004210360334),(4.35726558990816,3.64492370567392,3.40455735015306, 3.64492370567392,4.35726558990816),(4.52004210360334,3.86602540378444,3.64492370567392, 3.86602540378444,4.52004210360334),(5.07179676972449,4.52004210360334,4.35726558990816, 4.52004210360334,5.07179676972449))) REPRESENTATION_ITEM('') SURFACE() );

source
NURBS.preAllocMethod
preAlloc(degree::Int, uVector::Vector{T})

Allocate memory for basisFun.

source
NURBS.preAllocDerMethod
preAlloc(degree::Int, evalpoints, numberDerivatives::Int)

Allocate memory for derBasisFun.

source
NURBS.preAllocNURBSsurfaceMethod
preAllocNURBSsurface(uDegree::Int, vDegree::Int, uVector, vVector, k::Int)

Allocate all memory for surfaceDerivativesPoints!

source
NURBS.readMultipatchFunction
readMultipatch(filename::String)

Read multipatch file organized as:

PATCH a 4 4 -> degree in u and v 5 5 -> number of control points in u and v 0.0 0.0 ... -> knot vector in u 0.0 0.0 ... -> knot vector in v 0.1 0.2 ... | 0.1 0.2 ... > xyz components of the control points (normalized with the weigths -> we remove this weighting when reading in the data) 0.1 0.2 ... | 1.0 1.0 ... -> weights

PATCH b ...

source
NURBS.readStepFunction
readStep(filename::String, T=Float64)

Read a step file.

So far only BSPLINESURFACEWITHKNOTS and BOUNDED_SURFACE() are supported.

source
NURBS.refineFunction
refine(kVec::Vector, controlPts::Vector, degree::Int, newParametricPoints::Vector, weights=[])

Refine a curve by inserting new parametric points into the curve's knot vector.

NOTE: There are more efficient ways to do this. See, e.g., 'The NURBS Book' p. 162.

source
NURBS.refineMethod
refine(C::Curve, newParametricPoints::Vector)

Convenience function to refine a curve.

source
NURBS.refineMethod
refine(S::Surface; U=[], V=[])

NOTE: There are more efficient ways to do this.

source
NURBS.removeKnotFunction
removeKnot(knotVec, controlPoints, degree::Int, pointToRemove::Real, multiplicity::Int, weights)

Remove a knot 'multiplicity' times.

source
NURBS.removeKnot!Function
removeKnot!(knotVec, controlPoints, degree::Int, pointToRemove::Real, multiplicity::Int, weights)

Modifies solely 'knotVec', 'controlPoints', and 'weights'.

source
NURBS.removeKnotMethod
removeKnot(C::Curve, pointToRemove::Real, multiplicity::Int)

Remove a knot from a curve ´multiplicity´ times.

source
NURBS.removeKnotUMethod
removeKnotU(S::Surface, pointToRemove::Real, multiplicity::Int)

Remove a knot in u-direction 'multiplicity' times.

TODO: replace by more efficient algorithm of Tiller.

source
NURBS.removeKnotVMethod
removeKnotV(S::Surface, pointToRemove::Real, multiplicity::Int)

Remove a knot in v-direction 'multiplicity' times.

TODO: replace by more efficient algorithm of Tiller.

source
NURBS.rotateMethod
rotate(shape, rotAxis, angle)

Rotate a shape around the rotation axis by an angle (in rad).

source
NURBS.rotationMatrixMethod
rotationMatrix(rotAxis, angle::Real)

Determine rotation matrix for a rotation axis and an angle (in rad).

source
NURBS.scale!Method
scale!(shapes::T, factor::Real) where {T<:Shape}

Scale a vector of shapes by a real factor.

source
NURBS.scale!Method
scale!(shape::T, factor::Real) where {T<:Shape}

Scale a shape by a real factor

source
NURBS.similarCurveMethod
similarCurve(curve::BsplineCurve, p::Int, kVec, cPts, w)

Construct B-spline curve from underlying data: ignore empty weights.

source
NURBS.similarCurveMethod
similarCurve(curve::NURBScurve, p::Int, kVec, cPts, w)

Construct NURBS curve from underlying data.

source
NURBS.similarSurfaceMethod
similarSurface(surface::BsplineSurface, pu::Int, pv::Int, uVec, vVec, cPts, w)

Construct B-spline surface from underlying data: ignore empty weights.

source
NURBS.similarSurfaceMethod
similarSurface(surface::NURBSsurface, pu::Int, pv::Int, uVec, vVec, cPts, weights)

Construct NURBS surface from underlying data.

source
NURBS.spanRangesMethod
spanRanges(Bspl::Bspline, points)

Determine the ranges of the points which lie in each span of the B-spline (assuming normalized open knot vectors).

Return a vector of ranges (one entry per span).

source
NURBS.splitControlPointsMethod
splitControlPoints(degree::Int, wOri, ctrlPtsOri, kVecOri, oldMultIndices, oldMult, limitedMult, splitPoint)

Split the control points (and the weights) by inserting a knot multiple times.

source
NURBS.splitDataMethod
splitData(degree::Int, wOri, kVecOri, ctrlPtsOri, splitPoint=0.5)

Split the underlying data of a curve at a single parametric point in the range ]0,1[ by inserting a single knot multiple times.

source
NURBS.splitKnotsMethod
splitKnots(degree::Int, kVecOri, splitPoint)

Split knot vector according to the split point.

source
NURBS.splitUMethod
splitU(S::SurfaceT, splits::Vector) where {SurfaceT<:Surface}

Split a surface along the u-direction at multiple parametric points, each in the range ]0,1[, by inserting knots.

The splits vector needs to be sorted and each value should only occur once (both are not checked since this function is meant for internal use).

source
NURBS.splitVMethod
splitV(S::SurfaceT, splits::Vector) where {SurfaceT<:Surface}

Split a surface along the v-direction at multiple parametric points, each in the range ]0,1[, by inserting knots.

The splits vector needs to be sorted and each value should only occur once (both are not checked since this function is meant for internal use).

source
NURBS.splitWeightsMethod
splitWeights(weights::Vector{T}, Ind::Int) where {T}

Split the weights including the case when there are no weights.

source
NURBS.surfaceDerivativesPoints!Method
surfaceDerivativesPoints!(prealloc, uDegree::Int, vDegree::Int, uKnotVector, vKnotVector, controlPoints, uVector, vVector, weights, k::Int)

Compute NURBS surface: given the knotvectors and the degrees in 'u' and 'v' direction, the surface is evaluated at the evaluation points (uVector, vVector).

Using (4.20) on page 136 of 'The NURBS Book'.

Control points ordering P_(xi,yj):

P11 ––- P12 ––- P13 –-> y / v direction | | | | | | P21 ––- P22 ––- P23 | | | | | | P31 ––- P32 ––- P_33 | x / u direction

Returns a (k x k) matrix where each entry is a matrix of size (uKnotVector x vKnotVector): surfaces[q, p] is the matrix for the (q-1)-th derivative in u-direction and the (p-1)-th derivative in v-direction.

Note: the efficient evaluation via the B-spline basis is employed (no use of the naive evaluation of the NURBS basis).

source
NURBS.surfaceDerivativesPointsMethod
surfaceDerivativesPoints(uDegree::Int, vDegree::Int, uKnotVector, vKnotVector, controlPoints, uVector, vVector, weights, k::Int)

Allocate memory and call surfaceDerivativesPoints!

source
NURBS.surfaceDerivativesPointsMethod
surfaceDerivativesPoints(uDegree::Int, vDegree::Int, uKnotVector, vKnotVector, controlPoints, uVector, vVector, k::Int)

Compute B-spline surface and its derivatives: given the knotvectors and the degrees in 'u' and 'v' direction, the surface and its derivatives are evaluated at the evaluation points (uVector, vVector).

Control points ordering P_(xi,yj):

P11 ––- P12 ––- P13 –-> y / v direction | | | | | | P21 ––- P22 ––- P23 | | | | | | P31 ––- P32 ––- P_33 | x / u direction

Returns a (k x k) matrix where each entry is a matrix of size (uKnotVector x vKnotVector): surfaces[q, p] is the matrix for the (q-1)-th derivative in u-direction and the (p-1)-th derivative in v-direction.

Note: not the most efficient implementation. TODO: implement algorithm A.37 and A.38 of 'The Nurbs book'

source
NURBS.surfaceDerivativesPointsUVMethod
surfaceDerivativesPointsUV(uDegree::Int, vDegree::Int, controlPoints, uVector, vVector, Nu, Nv, q::Int, p::Int, uSpan, vSpan)

Compute the q-th derivative along 'u' and the p-th derivative along 'v'.

source
NURBS.surfacePointsMethod
surfacePoints(uBasis::Basis, vBasis::Basis, controlPoints, uVector, vVector)

Compute B-spline surface: given the knotvectors and the degrees in 'u' and 'v' direction, the surface is evaluated at the evaluation points (uVector, vVector).

Control points ordering P_(xi,yj):

P11 ––- P12 ––- P13 –-> y / v direction | | | | | | P21 ––- P22 ––- P23 | | | | | | P31 ––- P32 ––- P_33 | x / u direction

source
NURBS.surfacePointsMethod
surfacePoints(uBasis::Basis, vBasis::Basis, controlPoints, uVector, vVector, weights)

Compute NURBS surface: given the knotvectors and the degrees in 'u' and 'v' direction, the surface is evaluated at the evaluation points (uVector, vVector).

Control points ordering P_(xi,yj):

P11 ––- P12 ––- P13 –-> y / v direction | | | | | | P21 ––- P22 ––- P23 | | | | | | P31 ––- P32 ––- P_33 | x / u direction

Note: the efficient evaluation via the B-spline basis is employed (no use of the naive evaluation of the NURBS basis).

source
NURBS.translate!Method
translate!(shape::S, shift::SVector{3,T}) where {S<:Shape, T}

Translate a shape into the direction given by the vector 'shift'.

source
NURBS.translate!Method
translate!(shapes::Vector{S}, shift::SVector{3, T}) where {S<:Shape, T}

Translate a vector of shapes into the direction given by the vector 'shift'.

source
NURBS.translateMethod
translate(shape, shift)

Translate a shape into the direction given by the vector 'shift'.

source
NURBS.trimControlPoints!Method
trimControlPoints!(cPts, kVecOri, degree::Int, pos::Int, pointToRemove::Real, limitedMult::Int, oldMult::Int, weights)

Adaption of Algorithm A5.8 from 'The NURBS Book' p. 185.

Modifies solely 'cPts' and 'weights'.

source
NURBS.trimabilityMethod
trimability(knotVec, pointToRemove::Real, multiplicity::Int)

Analyze the knot vector: how many times is the knot actually contained, at which positions, and how many times.

source
NURBS.uniqueVerticesMethod
uniqueVertices(Patches::Vector{<:Surface{F}}; tol=1e-3) where {F}

Retreive the unique corner points of the patches in the provided list.

source
NURBS.weightsFunction
weights(srfc::NURBSsurface)

Return the weights of a NURBS surface.

source
NURBS.weightsMethod
weights(basis::Basis)

Except for a NURBS basis all other bases have no weights.

source
NURBS.weightsMethod
weights(srfc::Surface{T}, i=0, j=0) where {T}

Except for a NURBS surfaces all other surfaces have no weights.

source