NURBS module
G+Smo (Geometry + Simulation Modules): NURBS module
- class pygismo.nurbs.gsBSpline
Bases:
gsGeometry
- coefDim(self: pygismo.nurbs.gsBSpline) int
Returns the number of coefficients defining this B-Spline
- degree(self: pygismo.nurbs.gsBSpline, direction: int = 0) int
Returns the degree of the B-Spline
- degreeElevate(self: pygismo.nurbs.gsBSpline, arg0: int, arg1: int) None
Elevate the degree
- domainEnd(self: pygismo.nurbs.gsBSpline) float
Returns the right end of the domain
- domainStart(self: pygismo.nurbs.gsBSpline) float
Returns the left end of the domain
- insertKnot(self: pygismo.nurbs.gsBSpline, knot: float, i: int = 1) None
Insert a knot with multiplicity i without changing the curve
- knots(*args, **kwargs)
Overloaded function.
knots(self: pygismo.nurbs.gsBSpline, direction: int = 0) -> pygismo.nurbs.gsKnotVector
Gets the knots
knots(self: pygismo.nurbs.gsBSpline, direction: int = 0) -> pygismo.nurbs.gsKnotVector
Gets the knots as a reference
- numCoefs(self: pygismo.nurbs.gsBSpline) int
Returns the number of coefficients
- sample(self: pygismo.nurbs.gsBSpline, arg0: int) numpy.ndarray[numpy.float64[m, n]]
Returns samples on the Bspline curve
- class pygismo.nurbs.gsBSplineBasis
Bases:
gsBasis
- component(*args, **kwargs)
Overloaded function.
component(self: pygismo.nurbs.gsBSplineBasis, arg0: int) -> pygismo.nurbs.gsBSplineBasis
Returns the basis component as a reference
component(self: pygismo.nurbs.gsBSplineBasis, arg0: int) -> pygismo.nurbs.gsBSplineBasis
Returns the basis component as a const reference
- dim(self: pygismo.nurbs.gsBSplineBasis) int
Returns the dimension
- eval(self: pygismo.nurbs.gsBSplineBasis, arg0: numpy.ndarray[numpy.float64[m, n]]) numpy.ndarray[numpy.float64[m, n]]
Evaluates points into a matrix
- evalSingle(self: pygismo.nurbs.gsBSplineBasis, arg0: int, arg1: numpy.ndarray[numpy.float64[m, n]]) numpy.ndarray[numpy.float64[m, n]]
Evaluates the basis function i
- evalSingle_into(self: pygismo.nurbs.gsBSplineBasis, arg0: int, arg1: numpy.ndarray[numpy.float64[m, n]], arg2: numpy.ndarray[numpy.float64[m, n]]) None
Evaluates the basis function i
- function(self: pygismo.nurbs.gsBSplineBasis, arg0: int) pygismo.core.gsBasisFun
Returns the basis function i
- knots(*args, **kwargs)
Overloaded function.
knots(self: pygismo.nurbs.gsBSplineBasis, arg0: int) -> pygismo.nurbs.gsKnotVector
Get the knot vector as a reference
knots(self: pygismo.nurbs.gsBSplineBasis, arg0: int) -> pygismo.nurbs.gsKnotVector
Get the knot vector as a const reference
- numElements(self: pygismo.nurbs.gsBSplineBasis, arg0: pygismo.core.boxSide) int
Returns the number of Elements
- size(self: pygismo.nurbs.gsBSplineBasis) int
Returns the size
- class pygismo.nurbs.gsKnotVector
Bases:
pybind11_object
- check(self: pygismo.nurbs.gsKnotVector) bool
Checks whether the knot vector is in a consistent state
- degree(self: pygismo.nurbs.gsKnotVector) int
Returns the degree of the knot vector
- first(self: pygismo.nurbs.gsKnotVector) float
Returns the first knot
- get(self: pygismo.nurbs.gsKnotVector) List[float]
Returns the knot vector data
- greville(self: pygismo.nurbs.gsKnotVector) numpy.ndarray[numpy.float64[m, n]]
Returns the Greville points
- iFind(self: pygismo.nurbs.gsKnotVector, arg0: float) __gnu_cxx::__normal_iterator<double const*, std::vector<double, std::allocator<double> > >
Returns pointer to the last occurrence of the knot at the beginning of the knot interval containing the knot
- inDomain(self: pygismo.nurbs.gsKnotVector, arg0: float) bool
Checks, whether the given value is inside the domain
- insert(self: pygismo.nurbs.gsKnotVector, arg0: float, arg1: int) None
- isConsistent(self: List[float], arg0: List[int]) bool
Sanity check
- knot(self: pygismo.nurbs.gsKnotVector, arg0: int) float
Returns the i-th knot
- last(self: pygismo.nurbs.gsKnotVector) float
Returns the last knot
- multiplicities(self: pygismo.nurbs.gsKnotVector) List[int]
Returns vector of multiplicities of the knots
- numElements(self: pygismo.nurbs.gsKnotVector) int
Returns the number of knot intervals inside the domain
- size(self: pygismo.nurbs.gsKnotVector) int
Returns the KnotVector number of knots including repetitions
- uFind(self: pygismo.nurbs.gsKnotVector, arg0: float) gismo::internal::gsUKnotIterator<double>
Returns poiter to the knot at the beginning of the knot interval containing the knot
- uSize(self: pygismo.nurbs.gsKnotVector) int
Returns the KnotVector number of knots without repetitions
- uValue(self: pygismo.nurbs.gsKnotVector, arg0: int) float
Returns the value of the i-th knot
- class pygismo.nurbs.gsNurbsCreator
Bases:
pybind11_object
- static BSplineAmoeba(r: float = 1, x: float = 0, y: float = 0) pygismo.nurbs.gsBSpline
- static BSplineAmoeba3degree(r: float = 1, x: float = 0, y: float = 0) pygismo.nurbs.gsBSpline
- static BSplineAmoebaBig(r: float = 1, x: float = 0, y: float = 0) pygismo.nurbs.gsBSpline
- static BSplineAustria(r: float = 1, x: float = 0, y: float = 0) pygismo.nurbs.gsBSpline
- static BSplineCube(*args, **kwargs)
Overloaded function.
BSplineCube(arg0: float, arg1: float, arg2: float, arg3: float) -> pygismo.nurbs.gsTensorBSpline3
Cube of side r, with lower left corner at (x,y,z)
BSplineCube(arg0: int) -> pygismo.nurbs.gsTensorBSpline3
The unit cube represented as a tensor B-spline of degree deg
- static BSplineCubeGrid(n: int, m: int, p: int, r: float = 1, lx: float = 0, ly: float = 0, lz: float = 0) pygismo.core.gsMultiPatch
Creates a m n X m m X p rectangle multipatch consisting of B-splines squares with lower left corner at at (lx,ly,lz)
- static BSplineE(r: float = 1, x: float = 0, y: float = 0) pygismo.nurbs.gsBSpline
- static BSplineFatCircle(r: float = 1, x: float = 0, y: float = 0) pygismo.nurbs.gsBSpline
Inexact circle using B-splines
- static BSplineFatDisk(r: float = 1, x: float = 0, y: float = 0) pygismo.nurbs.gsTensorBSpline2
Inexact disk using B-splines
- static BSplineFatQuarterAnnulus(r0: float = 0, r1: float = 1) pygismo.nurbs.gsTensorBSpline2
Fat annulus using B-splines, discarding the weights of the exact NURBS
- static BSplineFish(r: float = 1, x: float = 0, y: float = 0) pygismo.nurbs.gsBSpline
- static BSplineHalfCube(r: float = 1, x: float = 0, y: float = 0, z: float = 0) pygismo.nurbs.gsTensorBSpline3
- static BSplineLShapeMultiPatch_p2() pygismo.core.gsMultiPatch
L-Shaped domain represented as a multipatch (3 patches) tensor B-spline of degree 2. 1. Patch is the middel part, 2. Patch is the upper part, 3 Patch is the right part.
- static BSplineLShape_p1(r: float = 1) pygismo.nurbs.gsTensorBSpline2
L-Shaped domain represented as a tensor B-spline of degree 1
- static BSplineLShape_p2C0() pygismo.nurbs.gsTensorBSpline2
L-Shaped domain represented as a tensor B-spline of degree 2, with C0-continuity across the diagonal.
- static BSplineLShape_p2C1() pygismo.nurbs.gsTensorBSpline2
L-Shaped domain represented as a tensor B-spline of degree 2, with C1-continuity and double control points at the corners.
- static BSplineLineSegment(arg0: numpy.ndarray[numpy.float64[m, n]], arg1: numpy.ndarray[numpy.float64[m, n]]) pygismo.nurbs.gsBSpline
- static BSplineQuarterAnnulus(deg: int = 2) pygismo.core.gsGeometry
Inexact annulus using B-splines
- static BSplineRectangle(low_x: float = 0, low_y: float = 0, upp_x: float = 1, upp_y: float = 1, turn_def: float = 0) pygismo.nurbs.gsTensorBSpline2
2d-rectange [low_x,upp_x] x [low_y,upp_y], rotated by turndeg degrees.
- static BSplineRectangleWithPara(low_x: float = 0, low_y: float = 0, upp_x: float = 1, upp_y: float = 1) pygismo.nurbs.gsTensorBSpline2
Rectangle described by the identity mapping over the given parameter domain, using tensor product B-splines.
- static BSplineSaddle() gismo::gsTensorNurbs<3, double>
Saddle using B-splines
- static BSplineSegment(u0: float = 0, u1: float = 1) pygismo.nurbs.gsBSpline
- static BSplineSquare(*args, **kwargs)
Overloaded function.
BSplineSquare(arg0: numpy.ndarray[numpy.float64[m, n]]) -> pygismo.nurbs.gsTensorBSpline2
Square of side defined by two corners (columsn)
BSplineSquare(arg0: float, arg1: float, arg2: float) -> pygismo.nurbs.gsTensorBSpline2
Square of side r, with lower left corner at (x,y)
- static BSplineSquareDeg(deg: int, scale: float = 1) pygismo.nurbs.gsTensorBSpline2
The unit square represented as a tensor B-spline of degree deg
- static BSplineSquareGrid(n: int, m: int, r: float = 1, lx: float = 0, ly: float = 0) pygismo.core.gsMultiPatch
Creates a m n X m m rectangle multipatch consisting of B-splines squares with lower left corner at at (lx,ly)
- static BSplineTriangle(*args, **kwargs)
Overloaded function.
BSplineTriangle(H: float = 1, W: float = 0) -> pygismo.nurbs.gsTensorBSpline2
Makes a Isosceles triangle with height H and width W
BSplineTriangle(H: float = 1, W: float = 0) -> pygismo.nurbs.gsTensorBSpline2
Makes a star with N patches, outer radius R0 and inner radius R1
- static BSplineUnitInterval(arg0: int) pygismo.nurbs.gsBSpline
Creates a B-Spline unit interval
- static NurbsAmoebaFull(r: float = 1, x: float = 0, y: float = 0) gismo::gsNurbs<double>
- static NurbsAnnulus(r0: float = 0, r1: float = 1) gismo::gsTensorNurbs<2, double>
Exact full annulus using NURBS
- static NurbsBean(r: float = 1, x: float = 0, y: float = 0) gismo::gsNurbs<double>
- static NurbsCircle(r: float = 1, x: float = 0, y: float = 0) gismo::gsNurbs<double>
Circle using NURBS
- static NurbsCube(r: float = 1, x: float = 0, y: float = 0, z: float = 0) gismo::gsTensorNurbs<3, double>
Cube of side r, with lower left corner at (x,y,z) using NURBS
- static NurbsCurve1(r: float = 1, x: float = 0, y: float = 0) gismo::gsNurbs<double>
- static NurbsCurve2(r: float = 1, x: float = 0, y: float = 0) gismo::gsNurbs<double>
- static NurbsDisk(r: float = 1, x: float = 0, y: float = 0) gismo::gsTensorNurbs<2, double>
- static NurbsQrtPlateWHoleC0() pygismo.nurbs.gsTensorBSpline2
- static NurbsQuarterAnnulus(r0: float = 0, r1: float = 1) gismo::gsTensorNurbs<2, double>
Exact annulus using NURBS
- static NurbsSphere(r: float = 1, x: float = 0, y: float = 0, z: float = 0) gismo::gsTensorNurbs<2, double>
Sphere using NURBS
- static lift3D(*args, **kwargs)
Overloaded function.
lift3D(arg0: pygismo.nurbs.gsTensorBSpline2, arg1: float) -> pygismo.nurbs.gsTensorBSpline3
Lifts a 2D geometry into 3D space
lift3D(arg0: gismo::gsTensorNurbs<2, double>, arg1: float) -> gismo::gsTensorNurbs<3, double>
Lifts a 2D geometry into 3D space
- static lift4D(*args, **kwargs)
Overloaded function.
lift4D(arg0: pygismo.nurbs.gsTensorBSpline3, arg1: float) -> pygismo.nurbs.gsTensorBSpline4
Lifts a 3D geometry into 4D space
lift4D(arg0: gismo::gsTensorNurbs<3, double>, arg1: float) -> gismo::gsTensorNurbs<4, double>
Lifts a 3D geometry into 4D space
- class pygismo.nurbs.gsTensorBSpline2
Bases:
gsGeometry
- degree(self: pygismo.nurbs.gsTensorBSpline2, arg0: int) int
Returns the degree
- knots(*args, **kwargs)
Overloaded function.
knots(self: pygismo.nurbs.gsTensorBSpline2, arg0: int) -> pygismo.nurbs.gsKnotVector
Get the knot vector as a reference
knots(self: pygismo.nurbs.gsTensorBSpline2, arg0: int) -> pygismo.nurbs.gsKnotVector
Get the knot vector as a const reference
- class pygismo.nurbs.gsTensorBSpline3
Bases:
gsGeometry
- degree(self: pygismo.nurbs.gsTensorBSpline3, arg0: int) int
Returns the degree
- knots(*args, **kwargs)
Overloaded function.
knots(self: pygismo.nurbs.gsTensorBSpline3, arg0: int) -> pygismo.nurbs.gsKnotVector
Get the knot vector as a reference
knots(self: pygismo.nurbs.gsTensorBSpline3, arg0: int) -> pygismo.nurbs.gsKnotVector
Get the knot vector as a const reference
- class pygismo.nurbs.gsTensorBSpline4
Bases:
gsGeometry
- degree(self: pygismo.nurbs.gsTensorBSpline4, arg0: int) int
Returns the degree
- knots(*args, **kwargs)
Overloaded function.
knots(self: pygismo.nurbs.gsTensorBSpline4, arg0: int) -> pygismo.nurbs.gsKnotVector
Get the knot vector as a reference
knots(self: pygismo.nurbs.gsTensorBSpline4, arg0: int) -> pygismo.nurbs.gsKnotVector
Get the knot vector as a const reference
- class pygismo.nurbs.gsTensorBSplineBasis2
Bases:
gsBasis
- active(self: pygismo.nurbs.gsTensorBSplineBasis2, arg0: numpy.ndarray[numpy.float64[m, n]]) numpy.ndarray[numpy.int32[m, n]]
Gives actives at points into a matrix
- component(*args, **kwargs)
Overloaded function.
component(self: pygismo.nurbs.gsTensorBSplineBasis2, arg0: int) -> pygismo.nurbs.gsBSplineBasis
Returns the basis component as a reference
component(self: pygismo.nurbs.gsTensorBSplineBasis2, arg0: int) -> pygismo.nurbs.gsBSplineBasis
Returns the basis component as a const reference
- degree(self: pygismo.nurbs.gsTensorBSplineBasis2, arg0: int) int
Returns the degree
- deriv(self: pygismo.nurbs.gsTensorBSplineBasis2, arg0: numpy.ndarray[numpy.float64[m, n]]) numpy.ndarray[numpy.float64[m, n]]
Evaluates derivatives at points into a matrix
- deriv2(self: pygismo.nurbs.gsTensorBSplineBasis2, arg0: numpy.ndarray[numpy.float64[m, n]]) numpy.ndarray[numpy.float64[m, n]]
Evaluates second derivatives at points into a matrix
- dim(self: pygismo.nurbs.gsTensorBSplineBasis2) int
Returns the dimension
- eval(self: pygismo.nurbs.gsTensorBSplineBasis2, arg0: numpy.ndarray[numpy.float64[m, n]]) numpy.ndarray[numpy.float64[m, n]]
Evaluates points into a matrix
- evalSingle(self: pygismo.nurbs.gsTensorBSplineBasis2, arg0: int, arg1: numpy.ndarray[numpy.float64[m, n]]) numpy.ndarray[numpy.float64[m, n]]
Evaluates the basis function i
- evalSingle_into(self: pygismo.nurbs.gsTensorBSplineBasis2, arg0: int, arg1: numpy.ndarray[numpy.float64[m, n]], arg2: numpy.ndarray[numpy.float64[m, n]]) None
Evaluates the basis function i
- function(self: pygismo.nurbs.gsTensorBSplineBasis2, arg0: int) pygismo.core.gsBasisFun
Returns the basis function i
- knots(*args, **kwargs)
Overloaded function.
knots(self: pygismo.nurbs.gsTensorBSplineBasis2, arg0: int) -> pygismo.nurbs.gsKnotVector
Get the knot vector as a reference
knots(self: pygismo.nurbs.gsTensorBSplineBasis2, arg0: int) -> pygismo.nurbs.gsKnotVector
Get the knot vector as a const reference
- size(self: pygismo.nurbs.gsTensorBSplineBasis2) int
Returns the size
- class pygismo.nurbs.gsTensorBSplineBasis3
Bases:
gsBasis
- component(*args, **kwargs)
Overloaded function.
component(self: pygismo.nurbs.gsTensorBSplineBasis3, arg0: int) -> pygismo.nurbs.gsBSplineBasis
Returns the basis component as a reference
component(self: pygismo.nurbs.gsTensorBSplineBasis3, arg0: int) -> pygismo.nurbs.gsBSplineBasis
Returns the basis component as a const reference
- degree(self: pygismo.nurbs.gsTensorBSplineBasis3, arg0: int) int
Returns the degree
- dim(self: pygismo.nurbs.gsTensorBSplineBasis3) int
Returns the dimension
- eval(self: pygismo.nurbs.gsTensorBSplineBasis3, arg0: numpy.ndarray[numpy.float64[m, n]]) numpy.ndarray[numpy.float64[m, n]]
Evaluates points into a matrix
- evalSingle(self: pygismo.nurbs.gsTensorBSplineBasis3, arg0: int, arg1: numpy.ndarray[numpy.float64[m, n]]) numpy.ndarray[numpy.float64[m, n]]
Evaluates the basis function i
- evalSingle_into(self: pygismo.nurbs.gsTensorBSplineBasis3, arg0: int, arg1: numpy.ndarray[numpy.float64[m, n]], arg2: numpy.ndarray[numpy.float64[m, n]]) None
Evaluates the basis function i
- function(self: pygismo.nurbs.gsTensorBSplineBasis3, arg0: int) pygismo.core.gsBasisFun
Returns the basis function i
- knots(*args, **kwargs)
Overloaded function.
knots(self: pygismo.nurbs.gsTensorBSplineBasis3, arg0: int) -> pygismo.nurbs.gsKnotVector
Get the knot vector as a reference
knots(self: pygismo.nurbs.gsTensorBSplineBasis3, arg0: int) -> pygismo.nurbs.gsKnotVector
Get the knot vector as a const reference
- size(self: pygismo.nurbs.gsTensorBSplineBasis3) int
Returns the size
- class pygismo.nurbs.gsTensorBSplineBasis4
Bases:
gsBasis
- component(*args, **kwargs)
Overloaded function.
component(self: pygismo.nurbs.gsTensorBSplineBasis4, arg0: int) -> pygismo.nurbs.gsBSplineBasis
Returns the basis component as a reference
component(self: pygismo.nurbs.gsTensorBSplineBasis4, arg0: int) -> pygismo.nurbs.gsBSplineBasis
Returns the basis component as a const reference
- degree(self: pygismo.nurbs.gsTensorBSplineBasis4, arg0: int) int
Returns the degree
- dim(self: pygismo.nurbs.gsTensorBSplineBasis4) int
Returns the dimension
- eval(self: pygismo.nurbs.gsTensorBSplineBasis4, arg0: numpy.ndarray[numpy.float64[m, n]]) numpy.ndarray[numpy.float64[m, n]]
Evaluates points into a matrix
- evalSingle(self: pygismo.nurbs.gsTensorBSplineBasis4, arg0: int, arg1: numpy.ndarray[numpy.float64[m, n]]) numpy.ndarray[numpy.float64[m, n]]
Evaluates the basis function i
- evalSingle_into(self: pygismo.nurbs.gsTensorBSplineBasis4, arg0: int, arg1: numpy.ndarray[numpy.float64[m, n]], arg2: numpy.ndarray[numpy.float64[m, n]]) None
Evaluates the basis function i
- function(self: pygismo.nurbs.gsTensorBSplineBasis4, arg0: int) pygismo.core.gsBasisFun
Returns the basis function i
- knots(*args, **kwargs)
Overloaded function.
knots(self: pygismo.nurbs.gsTensorBSplineBasis4, arg0: int) -> pygismo.nurbs.gsKnotVector
Get the knot vector as a reference
knots(self: pygismo.nurbs.gsTensorBSplineBasis4, arg0: int) -> pygismo.nurbs.gsKnotVector
Get the knot vector as a const reference
- size(self: pygismo.nurbs.gsTensorBSplineBasis4) int
Returns the size