API Reference
Core.Union
— Method(basis::Bspline)(evalpoints)
Evaluate B-spline basis at all evalpoints.
Core.Union
— Method(basis::Bspline)(evalpoints)
Evaluate B-spline basis at all evalpoints.
NURBS.Bspline
— TypeBspline{F} <: Basis{F}
B-spline basis.
NURBS.Bspline
— Method(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).
NURBS.Bspline
— Method(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).
NURBS.BsplineCurve
— TypeBsplineCurve{F} <: Curve{F}
B-spline curve defined by the basis and the control points.
NURBS.BsplineCurve
— Method(curve::BsplineCurve)(uVector, k::Int)
Convenience function to compute points on all k derivatives of a B-spline curve.
NURBS.BsplineCurve
— Method(curve::BsplineCurve)(uVector)
Convenience function to compute points on a B-spline curve.
NURBS.BsplineSurface
— TypeBsplineSurface{F} <: Surface{F}
Surface defined by a B-spline basis and the control points.
NURBS.BsplineSurface
— Method(Patch::BsplineSurface)(uEvalpoints, vEvalpoints, k::Int)
Convenience function to compute points on all k derivatives of a B-spline surface.
NURBS.BsplineSurface
— Method(Patch::BsplineSurface)(uEvalpoints, vEvalpoints)
Convenience function to compute points on a B-spline surface.
NURBS.Connectivity
— MethodConnectivity(patches::Vector{<:Surface}, cU, cV; tol=1e-3)
Constructor for Connectivity structure.
NURBS.CurrySchoenberg
— TypeCurrySchoenberg{F} <: Basis{F}
Normalized B-spline basis.
NURBS.Interface
— TypeInterface
An interface is defined by two patches: their ID and wether the vectors normal to the interface point in the same direction.
NURBS.InterfacePatchwise
— TypeInterfacePatchwise
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.
NURBS.NURB
— TypeNURB{F} <: Basis{F}
NURBS basis.
NURBS.NURBScurve
— TypeNURBScurve{F} <: Curve{F}
B-spline curve defined by the basis and the control points.
NURBS.NURBScurve
— Method(curve::NURBScurve)(uVector, k::Int)
Convenience function to compute points on all k derivatives of a NURBS curve.
NURBS.NURBScurve
— Method(curve::NURBScurve)(uVector)
Convenience function to compute points on a NURBS curve.
NURBS.NURBSsurface
— TypeNURBSsurface{F} <: Surface{F}
Surface defined by a B-spline basis, the control points, and the weights.
NURBS.NURBSsurface
— Method(Patch::NURBSsurface)(uEvalpoints, vEvalpoints, k::Int, prealloc)
Convenience function to compute points on all k derivatives of a NURBSsurface, for preallocated memory.
NURBS.NURBSsurface
— Method(Patch::NURBSsurface)(uEvalpoints, vEvalpoints, k::Int)
Convenience function to compute points on all k derivatives of a NURBSsurface.
NURBS.NURBSsurface
— Method(Patch::NURBSsurface)(uEvalpoints, vEvalpoints)
Convenience function to compute points on a NURBSsurface.
NURBS.PatchInterface
— TypeLocal 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.
Base.split
— FunctionBase.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.
Base.split
— MethodBase.split(C::Curve, n::Int)
Split a curve in the parametric domain equally to obtain n curves.
Base.split
— MethodBase.split(C::Curve, n::Int)
Split a curve in the parametric domain equally to obtain n curves.
Base.split
— MethodBase.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.
Base.split
— MethodBase.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.
LinearAlgebra.rotate!
— Methodrotate!(shape::S, rotAxis::SVector{3,T}, angle::Real) where {S<:Shape,T}
Rotate a shape around the rotation axis by an angle (in rad).
LinearAlgebra.rotate!
— Methodrotate!(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).
NURBS.JacobiDet
— MethodJacobiDet(Ju, Jv)
Compute Jacobi determinant from the Jacobi matrix
NURBS.Jacobian
— MethodJacobian(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.
NURBS.anchors
— Methodanchors(kVec, degree::Int)
Return the anchors (as defined in [4]) corresponding to the given knotvector and the polynomial degree.
NURBS.bSplineNaive
— MethodbSplineNaive(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].
NURBS.bSplineNaive
— MethodbSplineNaive(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.
NURBS.bSplineNaive
— MethodbSplineNaive(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.
NURBS.bSplineNaiveDerivative
— MethodbSplineNaiveDer(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'.
NURBS.bSplineNaiveDerivative
— MethodbSplineNaiveDerivative(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.)
NURBS.basisFun!
— MethodbasisFun!(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.
NURBS.basisFun
— MethodbasisFun(knotSpan, uVector, basis::Basis)
Allocate memory and call basisFun!
NURBS.checkRange
— MethodcheckRange(point<:Real)
Check whether points is greater 0 and smaller 1.
NURBS.coarsen
— Functioncoarsen(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.
NURBS.curveDerivativesPoints
— MethodcurveDerivativesPoints(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).
NURBS.curveDerivativesPoints
— MethodcurvePoints(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]
NURBS.curvePoints
— MethodcurvePoints(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).
NURBS.curvePoints
— MethodcurvePoints(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]
NURBS.derBasisFun!
— MethodderBasisFun!(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.
NURBS.derBasisFun
— MethodderBasisFun(degree::Int, evalpoints, knotVector, numberDerivatives::Int)
Same as derBasisFun! but with memory allocation.
NURBS.evalNaive
— MethodevalNaive(basis::NURB, i::Int, evalpoints)
i-th NURB basis function evaluated at all 'evalpoints'.
NURBS.evalNaive
— MethodevalNaive(basis::Basis, i::Int, evalpoints)
i-th B-spline or Curry-Schoenberg basis function evaluated at all 'evalpoints'.
NURBS.evalNaiveDerivative
— MethodevalNaiveDerivative(basis::Bspline, i::Int, k::Int, evalpoints)
Compute the k-th derivative of i-th b-spline basis function evaluated at all 'evalpoints'.
NURBS.evalNaiveDerivative
— MethodevalNaiveDerivative(basis::NURB, i::Int, k::Int, evalpoints)
Compute the k-th derivative of i-th NURBS basis function evaluated at all 'evalpoints'.
NURBS.extendControlPoints!
— FunctionextendControlPoints!(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.
NURBS.extendControlPoints!
— MethodLowest order polynomial degree = 0.
Can multiplicity be maximum 1 and/or oldMult maximum 0? => simplify below algorithm further.
NURBS.extendKnotVector!
— MethodextendKnotVector!(knotVecOrig, newParametricPoint::Real, multiplicity::Int)
Insert the new value with a given mulitplicity (insert the point multiple times).
Modifies knotVecOrig.
NURBS.extractBracketList
— FunctionextractBracketList(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'.
NURBS.extractCtrlPoints
— FunctionextractCtrlPoints(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))
NURBS.extractKnotVecs
— FunctionextractKnotVecs(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'.
NURBS.findSpan!
— MethodfindSpan!(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)
NURBS.findSpan
— MethodfindSpan(b::T, uVec, kVec, p::T)
Allocate memory and call findSpan!
NURBS.generateKnotVec
— MethodgenerateKnotVec(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].
NURBS.getCartesianPoints
— FunctiongetCartesianPoints(stringVec)
Extract all Cartesian points and the offset between the first point and the #counting. (Assumption: all points appear in a consecutive list.)
NURBS.getCompleteString
— MethodgetCompleteString(stringVec::String, nextLineInd::Int)
Retreive all lines of the string vector 'stringVec' starting from 'nextLineInd' till a semicolon is encountered.
NURBS.getCompleteStringTill
— MethodgetCompleteStringTill(stringVec, nextLineInd::Int, stopWhenEncountered::String)
Retreive all lines of the string vector 'stringVec' starting from 'nextLineInd' till the next line starts with 'stopWhenEncountered'.
NURBS.getGlobalIndex
— MethodgetGlobalIndex(P, Plist; tol=1e-3)
Get the global index of the point P with respect to the ordering in the list of points Plist.
NURBS.getPatchInterfaces
— MethodgetPatchInterfaces(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
NURBS.greville
— Methodgreville(kVec, degree::Int)
Return the Greville sites (as defined in [3]) corresponding to the given knotvector and the polynomial degree.
NURBS.identifyInterfaces
— MethodidentifyInterfaces(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.
NURBS.insertKnot
— FunctioninsertKnot(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.
NURBS.insertKnot
— FunctioninsertKnot(C::Curve, newParametricPoint::Real, multiplicity::Int)
Convenience function to insert a knot into a curve.
NURBS.insertKnot!
— FunctioninsertKnot!(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.
NURBS.insertKnotU
— FunctioninsertKnotU(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.
NURBS.insertKnotV
— FunctioninsertKnotV(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.
NURBS.isValidKnotVector!
— MethodisValidKnotVector!(kVec)
Check whether the knot vector has only entries in [0, 1] and is in ascending order.
If not the knot vector is modified.
NURBS.isin
— Methodisin(P, Plist; tol=1e-3)
Check if 3D point P is in list of points Plist.
NURBS.limitMultiplicity
— MethodlimitMultiplicity(newMult::Int, degree::Int)
Limit the multiplicity to the maximum allowed value.
NURBS.load
— Functionload(f::File{format"DAT"}, T=Float64)
Load a .dat file via the FileIO package.
NURBS.load
— Functionload(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
NURBS.mirror!
— Methodtranslate!(shape::S, shift::SVector{3,T}) where {S<:Shape, T}
Mirror a shape through a plane defined by its normal and an anchor point.
NURBS.mirror!
— Methodmirror!(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.
NURBS.mirror
— Methodmirror(shape, rotAxis, angle)
Mirror through a plane defined by its normal and an anchor point.
NURBS.normalization
— Methodnormalization(basis::CurrySchoenberg, i::Int)
Normalization for the standard B-splines to obtain Curry-Schoenberg splines.
NURBS.normalization
— Methodnormalization(basis::Basis, i::Int)
For Bsplines there is no normalization.
NURBS.normalizationDerivative
— MethodnormalizationDerivative(knotVector, degree::Int, weights, u::Real, k::Int)
Return the k-th derivative of the weight function (21) in [2].
NURBS.normalize!
— Methodnormalize!(N, i::T, degree::T, knotVector, basis::CurrySchoenberg)
Normalize the standard B-splines to obtain Curry-Schoenberg splines.
NURBS.normalize!
— Methodnormalize!(N, i::T, degree::T, knotVector, basis::Bspline)
For Bsplines there is no normalization.
NURBS.numBasisFunctions
— MethodnumBasisFunctions(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.
NURBS.nurbsNaive
— MethodnurbsNaive(knotVector, i::Int, degree::Int, evalpoints, weights)
i-th NURB basis function of degree 'degree', the 'knotVector', and with 'weights' evaluated at all 'evalpoints'.
NURBS.nurbsNaive
— MethodnurbsNaive(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.
NURBS.nurbsNaiveDerivative
— MethodnurbsNaiveDerivative(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'.
NURBS.nurbsNaiveDerivative
— MethodnurbsNaiveDerivative(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].
NURBS.parseCtrlPoints
— FunctionparseCtrlPoints(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.
NURBS.parseLine
— FunctionparseLine(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'.
NURBS.parseSinglePatch
— MethodparseSinglePatch(stringVec, lineInd)
Take the string vector and the line index where "PATCH" stands and extract the information of a single patch
NURBS.parse_B_SPLINE_SURFACE_WITH_KNOTS
— Functionparse_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.);
NURBS.parse_NURBS_step
— Functionparse_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() );
NURBS.preAlloc
— MethodpreAlloc(degree::Int, uVector::Vector{T})
Allocate memory for basisFun.
NURBS.preAllocDer
— MethodpreAlloc(degree::Int, evalpoints, numberDerivatives::Int)
Allocate memory for derBasisFun.
NURBS.preAllocNURBSsurface
— MethodpreAllocNURBSsurface(uDegree::Int, vDegree::Int, uVector, vVector, k::Int)
Allocate all memory for surfaceDerivativesPoints!
NURBS.readMultipatch
— FunctionreadMultipatch(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 ...
NURBS.readStep
— FunctionreadStep(filename::String, T=Float64)
Read a step file.
So far only BSPLINESURFACEWITHKNOTS and BOUNDED_SURFACE() are supported.
NURBS.refine
— Functionrefine(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.
NURBS.refine
— Methodrefine(C::Curve, newParametricPoints::Vector)
Convenience function to refine a curve.
NURBS.refine
— Methodrefine(S::Surface; U=[], V=[])
NOTE: There are more efficient ways to do this.
NURBS.removeKnot
— FunctionremoveKnot(knotVec, controlPoints, degree::Int, pointToRemove::Real, multiplicity::Int, weights)
Remove a knot 'multiplicity' times.
NURBS.removeKnot!
— FunctionremoveKnot!(knotVec, controlPoints, degree::Int, pointToRemove::Real, multiplicity::Int, weights)
Modifies solely 'knotVec', 'controlPoints', and 'weights'.
NURBS.removeKnot
— MethodremoveKnot(C::Curve, pointToRemove::Real, multiplicity::Int)
Remove a knot from a curve ´multiplicity´ times.
NURBS.removeKnotU
— MethodremoveKnotU(S::Surface, pointToRemove::Real, multiplicity::Int)
Remove a knot in u-direction 'multiplicity' times.
TODO: replace by more efficient algorithm of Tiller.
NURBS.removeKnotV
— MethodremoveKnotV(S::Surface, pointToRemove::Real, multiplicity::Int)
Remove a knot in v-direction 'multiplicity' times.
TODO: replace by more efficient algorithm of Tiller.
NURBS.rotate
— Methodrotate(shape, rotAxis, angle)
Rotate a shape around the rotation axis by an angle (in rad).
NURBS.rotationMatrix
— MethodrotationMatrix(rotAxis, angle::Real)
Determine rotation matrix for a rotation axis and an angle (in rad).
NURBS.saveVtk
— MethodsaveVtk(filename::String, patches; resolution=0.01)
Save as an unstructured grid in a .vtu file to be displayed by ParaView.
NURBS.scale!
— Methodscale!(shapes::T, factor::Real) where {T<:Shape}
Scale a vector of shapes by a real factor.
NURBS.scale!
— Methodscale!(shape::T, factor::Real) where {T<:Shape}
Scale a shape by a real factor
NURBS.scale
— Methodscale(shape, factor)
Scale a shape by a real factor.
NURBS.similarCurve
— MethodsimilarCurve(curve::BsplineCurve, p::Int, kVec, cPts, w)
Construct B-spline curve from underlying data: ignore empty weights.
NURBS.similarCurve
— MethodsimilarCurve(curve::NURBScurve, p::Int, kVec, cPts, w)
Construct NURBS curve from underlying data.
NURBS.similarSurface
— MethodsimilarSurface(surface::BsplineSurface, pu::Int, pv::Int, uVec, vVec, cPts, w)
Construct B-spline surface from underlying data: ignore empty weights.
NURBS.similarSurface
— MethodsimilarSurface(surface::NURBSsurface, pu::Int, pv::Int, uVec, vVec, cPts, weights)
Construct NURBS surface from underlying data.
NURBS.spanRanges
— MethodspanRanges(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).
NURBS.splitControlPoints
— MethodsplitControlPoints(degree::Int, wOri, ctrlPtsOri, kVecOri, oldMultIndices, oldMult, limitedMult, splitPoint)
Split the control points (and the weights) by inserting a knot multiple times.
NURBS.splitData
— MethodsplitData(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.
NURBS.splitKnots
— MethodsplitKnots(degree::Int, kVecOri, splitPoint)
Split knot vector according to the split point.
NURBS.splitU
— MethodsplitU(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).
NURBS.splitV
— MethodsplitV(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).
NURBS.splitWeights
— MethodsplitWeights(weights::Vector{T}, Ind::Int) where {T}
Split the weights including the case when there are no weights.
NURBS.surfaceDerivativesPoints!
— MethodsurfaceDerivativesPoints!(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).
NURBS.surfaceDerivativesPoints
— MethodsurfaceDerivativesPoints(uDegree::Int, vDegree::Int, uKnotVector, vKnotVector, controlPoints, uVector, vVector, weights, k::Int)
Allocate memory and call surfaceDerivativesPoints!
NURBS.surfaceDerivativesPoints
— MethodsurfaceDerivativesPoints(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'
NURBS.surfaceDerivativesPointsUV
— MethodsurfaceDerivativesPointsUV(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'.
NURBS.surfacePoints
— MethodsurfacePoints(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
NURBS.surfacePoints
— MethodsurfacePoints(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).
NURBS.translate!
— Methodtranslate!(shape::S, shift::SVector{3,T}) where {S<:Shape, T}
Translate a shape into the direction given by the vector 'shift'.
NURBS.translate!
— Methodtranslate!(shapes::Vector{S}, shift::SVector{3, T}) where {S<:Shape, T}
Translate a vector of shapes into the direction given by the vector 'shift'.
NURBS.translate
— Methodtranslate(shape, shift)
Translate a shape into the direction given by the vector 'shift'.
NURBS.trimControlPoints!
— MethodtrimControlPoints!(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'.
NURBS.trimability
— Methodtrimability(knotVec, pointToRemove::Real, multiplicity::Int)
Analyze the knot vector: how many times is the knot actually contained, at which positions, and how many times.
NURBS.uniqueVertices
— MethoduniqueVertices(Patches::Vector{<:Surface{F}}; tol=1e-3) where {F}
Retreive the unique corner points of the patches in the provided list.
NURBS.vtk
— Functionvtk(patches, resolution=0.01)
Prepare the patches to be displayed by Paraview:
The output can be written to the VTK format using the WriteVTK.jl package.
using WriteVTK
vtk_grid("path/filename", x, y, z, cellV) do vtk end
NURBS.weights
— Functionweights(srfc::NURBSsurface)
Return the weights of a NURBS surface.
NURBS.weights
— Methodweights(basis::NURB)
Return the weights of a NURBS basis.
NURBS.weights
— Methodweights(basis::Basis)
Except for a NURBS basis all other bases have no weights.
NURBS.weights
— Methodweights(srfc::Surface{T}, i=0, j=0) where {T}
Except for a NURBS surfaces all other surfaces have no weights.
NURBS.weights
— Methodweights(w::Array{T}, i, j) where {T}
General case.