Medusa  1.1
Coordinate Free Mehless Method implementation
mm::NURBSPatch< vec_t, param_vec_t > Class Template Reference

#include <NURBSPatch_fwd.hpp>

Detailed Description

template<typename vec_t, typename param_vec_t>
class mm::NURBSPatch< vec_t, param_vec_t >

Class representing a single NURBS patch in an arbitrary dimensional space, defined on an arbitrary dimensional domain and generated by a tensor product of NURBS curves.

Warning
Currently only supports NURBS surfaces (2D domain) and NURBS curves (1D domain).

Usage example:

int p = 2;
Range<double> w({1, sqrt(3) / 2, 1, sqrt(3) / 2, 1, sqrt(3) / 2, 1, sqrt(3) / 2, 1});
Range<double> knots({0, 0, 0, PI / 2, PI / 2, PI, PI, 3 * PI / 2,
3 * PI / 2, 2 * PI, 2 * PI, 2 * PI});
Range<Vec3d> cp({Vec3d(1, 0, 0), Vec3d(1, 1, 0), Vec3d(0, 1, 0), Vec3d(-1, 1, 0),
Vec3d(-1, 0, 0), Vec3d(-1, -1, 0), Vec3d(0, -1, 0), Vec3d(1, -1, 0),
Vec3d(1, 0, 0)});
NURBSPatch<Vec3d, Vec1d> curve(cp, w, {knots}, {p});
Vec3d pt = curve.evaluate({PI});
Template Parameters
param_vec_tVector type of the parametric space.
vec_tVector type of ambient space.

Definition at line 72 of file NURBSPatch_fwd.hpp.

Public Member Functions

 NURBSPatch (const NdRange< proj_vec_t > &wcp, const std::array< Range< scalar_t >, param_dim > &ks, const std::array< int, param_dim > &in_p)
 Construct NURBS patch with weighted control points. More...
 
 NURBSPatch (const NdRange< vec_t > &cp, const NdRange< double > &w, const std::array< Range< scalar_t >, param_dim > &ks, const std::array< int, param_dim > &in_p)
 Construct NURBS patch with control points and weights. More...
 
 NURBSPatch (const NURBSPatch< vec_t, param_vec_t > &patch) noexcept
 Copy constructor. More...
 
NURBSPatch< vec_t, param_vec_t > & operator= (const NURBSPatch< vec_t, param_vec_t > &) noexcept
 Copy assignment. More...
 
 NURBSPatch (NURBSPatch< vec_t, param_vec_t > &&) noexcept=default
 Move constructor. More...
 
NURBSPatch< vec_t, param_vec_t > & operator= (NURBSPatch< vec_t, param_vec_t > &&) noexcept=default
 Move assignment. More...
 
void computeDerivativeStructure ()
 Calculate data needed for derivative evaluation. More...
 
vec_t evaluate (const param_vec_t &t, scalar_t *w) const
 Evaluate NURBS patch in one point along the second dimension. More...
 
vec_t evaluate (const param_vec_t &t) const
 Overload if weight is not required. More...
 
proj_vec_t evaluateWeighted (const param_vec_t &t) const
 Evaluate NURBS patch in one point. More...
 
Eigen::Matrix< scalar_t, dim, param_dimjacobian (const param_vec_t &t, const vec_t &pt, const scalar_t &w) const
 Evaluate NURBS jacobian matrix in one point. More...
 
Eigen::Matrix< scalar_t, dim, param_dimjacobian (const param_vec_t &t) const
 Overload that calculates the point automatically. More...
 
std::pair< vec_t, Eigen::Matrix< scalar_t, dim, param_dim > > evaluatePointAndJacobian (const param_vec_t &t) const
 Evaluate NURBS and its jacobian matrix in one parameter. More...
 
BoxShape< param_vec_t > getDomain () const
 Get domain of the NURBS patch. More...
 
Range< NURBSPatch< vec_t, Vec< scalar_t, param_dim - 1 > > > getBoundaries () const
 Get all boundaries of the NURBS patch as a Range of NURBS patches. More...
 
param_vec_t getPatchParameterFromBoundaryParameter (const Vec< scalar_t, param_dim - 1 > &t, int i, scalar_t epsilon=0) const
 Get vector from the parametric domain of the NURBS patch from the parameter of the NURBS patch representing one of the boundaries of the NURBS. More...
 
Vec< scalar_t, param_dim - 1 > getBoundaryParameterFromPatchParameter (const param_vec_t &t, int i) const
 Get parameter from the parametric domain of the NURBS patch representing the boundary of the NURBS patch corresponding to a parameter from the parametric domain of the NURBS surface. More...
 

Public Types

enum  { dim = vec_t::dim, proj_dim = dim + 1, param_dim = param_vec_t::dim }
 Store dimensionality. More...
 
typedef vec_t::scalar_t scalar_t
 Scalar type. More...
 
typedef Vec< scalar_t, proj_dimproj_vec_t
 Vector type before projection. More...
 
template<typename T >
using NdRange = typename nurbs_patch_internal::NestedRange< param_dim, T >::type
 Range of param_dim dimensions. More...
 

Friends

struct nurbs_patch_internal::NURBSPatchHelper< vec_t, param_vec_t >
 

Private Attributes

std::array< int, param_dimp
 Underlying B-spline orders. More...
 
NdRange< proj_vec_tweighted_control_points
 Multidimensional Range of control points. More...
 
std::array< Range< scalar_t >, param_dimknots
 Array of Ranges of knots (knot vector). More...
 
std::unique_ptr< std::array< NURBSPatch< vec_t, param_vec_t >, param_dim > > der_structure
 NURBS surfaces needed for derivative calculation. More...
 

Member Enumeration Documentation

◆ anonymous enum

template<typename vec_t , typename param_vec_t >
anonymous enum

Store dimensionality.

Enumerator
dim 

Dimensionality of space after projection.

proj_dim 

Dimensionality of space before projection.

param_dim 

Dimensionality of parametric space.

Definition at line 78 of file NURBSPatch_fwd.hpp.

Constructor & Destructor Documentation

◆ NURBSPatch() [1/4]

template<typename vec_t , typename param_vec_t >
mm::NURBSPatch< vec_t, param_vec_t >::NURBSPatch ( const NdRange< proj_vec_t > &  wcp,
const std::array< Range< scalar_t >, param_dim > &  ks,
const std::array< int, param_dim > &  in_p 
)

Construct NURBS patch with weighted control points.

Parameters
wcpMultidimensional Range of weighted control points, e. g. Vec<double, 4>(w * x, w * y, w * z, w).
ksArray of Ranges of knots (pair of knot vectors).
in_pArray of underlying B-spline orders (starting from order 0).

Definition at line 17 of file NURBSPatch.hpp.

◆ NURBSPatch() [2/4]

template<typename vec_t , typename param_vec_t >
mm::NURBSPatch< vec_t, param_vec_t >::NURBSPatch ( const NdRange< vec_t > &  cp,
const NdRange< double > &  w,
const std::array< Range< scalar_t >, param_dim > &  ks,
const std::array< int, param_dim > &  in_p 
)

Construct NURBS patch with control points and weights.

Parameters
cpMultidimensional Range of control points.
wMultidimensional Range of weights.
ksArray of Ranges of knots (pair of knot vectors).
in_pArray of underlying B-spline orders (starting from order 0).

Definition at line 26 of file NURBSPatch.hpp.

◆ NURBSPatch() [3/4]

template<typename vec_t , typename param_vec_t >
mm::NURBSPatch< vec_t, param_vec_t >::NURBSPatch ( const NURBSPatch< vec_t, param_vec_t > &  patch)
noexcept

Copy constructor.

Definition at line 37 of file NURBSPatch.hpp.

◆ NURBSPatch() [4/4]

template<typename vec_t , typename param_vec_t >
mm::NURBSPatch< vec_t, param_vec_t >::NURBSPatch ( NURBSPatch< vec_t, param_vec_t > &&  )
defaultnoexcept

Move constructor.

Member Function Documentation

◆ computeDerivativeStructure()

template<typename vec_t , typename param_vec_t >
void mm::NURBSPatch< vec_t, param_vec_t >::computeDerivativeStructure

Calculate data needed for derivative evaluation.

Definition at line 64 of file NURBSPatch.hpp.

◆ evaluate() [1/2]

template<typename vec_t , typename param_vec_t >
vec_t mm::NURBSPatch< vec_t, param_vec_t >::evaluate ( const param_vec_t &  t) const

Overload if weight is not required.

Parameters
tPoint of evaluation.
Returns
Vector to the point on the NURBS surface corresponding to parameter t.

Definition at line 81 of file NURBSPatch.hpp.

◆ evaluate() [2/2]

template<typename vec_t , typename param_vec_t >
vec_t mm::NURBSPatch< vec_t, param_vec_t >::evaluate ( const param_vec_t &  t,
scalar_t w 
) const

Evaluate NURBS patch in one point along the second dimension.

Parameters
tPoint of evaluation.
wVariable in which the weight is stored.
Returns
Vector to the point on the NURBS surface corresponding to parameter t.

Definition at line 69 of file NURBSPatch.hpp.

◆ evaluatePointAndJacobian()

template<typename vec_t , typename param_vec_t >
std::pair< vec_t, Eigen::Matrix< typename NURBSPatch< vec_t, param_vec_t >::scalar_t, NURBSPatch< vec_t, param_vec_t >::dim, NURBSPatch< vec_t, param_vec_t >::param_dim > > mm::NURBSPatch< vec_t, param_vec_t >::evaluatePointAndJacobian ( const param_vec_t &  t) const

Evaluate NURBS and its jacobian matrix in one parameter.

Parameters
tParameter of evaluation.
Returns
Pair of point on the NURBS and jacobian matrix of the NURBS at parameter t.

Definition at line 114 of file NURBSPatch.hpp.

◆ evaluateWeighted()

template<typename vec_t , typename param_vec_t >
NURBSPatch< vec_t, param_vec_t >::proj_vec_t mm::NURBSPatch< vec_t, param_vec_t >::evaluateWeighted ( const param_vec_t &  t) const

Evaluate NURBS patch in one point.

Parameters
tPoint of evaluation.
Returns
Vector to the weighted point on the NURBS patch corresponding to parameter t.

Definition at line 88 of file NURBSPatch.hpp.

◆ getBoundaries()

template<typename vec_t , typename param_vec_t >
Range< NURBSPatch< vec_t, Vec< typename NURBSPatch< vec_t, param_vec_t >::scalar_t, NURBSPatch< vec_t, param_vec_t >::param_dim - 1 > > > mm::NURBSPatch< vec_t, param_vec_t >::getBoundaries

Get all boundaries of the NURBS patch as a Range of NURBS patches.

Returns
An array of NURBS curves representing the boundaries.

Definition at line 135 of file NURBSPatch.hpp.

◆ getBoundaryParameterFromPatchParameter()

template<typename vec_t , typename param_vec_t >
Vec< typename NURBSPatch< vec_t, param_vec_t >::scalar_t, NURBSPatch< vec_t, param_vec_t >::param_dim - 1 > mm::NURBSPatch< vec_t, param_vec_t >::getBoundaryParameterFromPatchParameter ( const param_vec_t &  t,
int  i 
) const

Get parameter from the parametric domain of the NURBS patch representing the boundary of the NURBS patch corresponding to a parameter from the parametric domain of the NURBS surface.

Parameters
tParameter of the NURBS patch.
iIndex of the boundary.
Returns
Parameter from the parametric domain of the NURBS boundary curve corresponding to the input parameter.

Definition at line 149 of file NURBSPatch.hpp.

◆ getDomain()

template<typename vec_t , typename param_vec_t >
BoxShape< param_vec_t > mm::NURBSPatch< vec_t, param_vec_t >::getDomain

Get domain of the NURBS patch.

Returns
Domain of the NURBS in the form of a BoxShape.

Definition at line 122 of file NURBSPatch.hpp.

◆ getPatchParameterFromBoundaryParameter()

template<typename vec_t , typename param_vec_t >
param_vec_t mm::NURBSPatch< vec_t, param_vec_t >::getPatchParameterFromBoundaryParameter ( const Vec< scalar_t, param_dim - 1 > &  t,
int  i,
scalar_t  epsilon = 0 
) const

Get vector from the parametric domain of the NURBS patch from the parameter of the NURBS patch representing one of the boundaries of the NURBS.

Patch parameter can be generated epsilon times the length of the domain away from the boundary. This is useful when normals are not defined at the boundary of the patch.

Parameters
tParameter of the NURBS patch representing the boundary.
iIndex of the boundary.
epsilonEpsilon.
Returns
Vector from the parametric domain of the NURBS patch corresponding to the input parameter.

Definition at line 140 of file NURBSPatch.hpp.

◆ jacobian() [1/2]

template<typename vec_t , typename param_vec_t >
Eigen::Matrix< typename NURBSPatch< vec_t, param_vec_t >::scalar_t, NURBSPatch< vec_t, param_vec_t >::dim, NURBSPatch< vec_t, param_vec_t >::param_dim > mm::NURBSPatch< vec_t, param_vec_t >::jacobian ( const param_vec_t &  t) const

Overload that calculates the point automatically.

Warning
Less efficient.
Parameters
tPoint of evaluation.
Returns
Jacobian matrix of the NURBS patch at point t.
Warning
computeDerivativeStructure() must be called on the object before calling this function.

Definition at line 103 of file NURBSPatch.hpp.

◆ jacobian() [2/2]

template<typename vec_t , typename param_vec_t >
Eigen::Matrix< typename NURBSPatch< vec_t, param_vec_t >::scalar_t, NURBSPatch< vec_t, param_vec_t >::dim, NURBSPatch< vec_t, param_vec_t >::param_dim > mm::NURBSPatch< vec_t, param_vec_t >::jacobian ( const param_vec_t &  t,
const vec_t &  pt,
const scalar_t w 
) const

Evaluate NURBS jacobian matrix in one point.

Parameters
tPoint of evaluation.
ptPoint on the NURBS patch corresponding to t.
wWeight of pt.
Returns
Jacobian matrix of the NURBS patch at point t.
Warning
computeDerivativeStructure() must be called on the object before calling this function.

Definition at line 95 of file NURBSPatch.hpp.

◆ operator=() [1/2]

template<typename vec_t , typename param_vec_t >
NURBSPatch< vec_t, param_vec_t > & mm::NURBSPatch< vec_t, param_vec_t >::operator= ( const NURBSPatch< vec_t, param_vec_t > &  other)
noexcept

Copy assignment.

Definition at line 47 of file NURBSPatch.hpp.

◆ operator=() [2/2]

template<typename vec_t , typename param_vec_t >
NURBSPatch<vec_t, param_vec_t>& mm::NURBSPatch< vec_t, param_vec_t >::operator= ( NURBSPatch< vec_t, param_vec_t > &&  )
defaultnoexcept

Move assignment.

Member Data Documentation

◆ der_structure

template<typename vec_t , typename param_vec_t >
std::unique_ptr<std::array<NURBSPatch<vec_t, param_vec_t>, param_dim> > mm::NURBSPatch< vec_t, param_vec_t >::der_structure
private

NURBS surfaces needed for derivative calculation.

Definition at line 93 of file NURBSPatch_fwd.hpp.

◆ knots

template<typename vec_t , typename param_vec_t >
std::array<Range<scalar_t>, param_dim> mm::NURBSPatch< vec_t, param_vec_t >::knots
private

Array of Ranges of knots (knot vector).

Definition at line 90 of file NURBSPatch_fwd.hpp.

◆ p

template<typename vec_t , typename param_vec_t >
std::array<int, param_dim> mm::NURBSPatch< vec_t, param_vec_t >::p
private

Underlying B-spline orders.

Definition at line 88 of file NURBSPatch_fwd.hpp.

◆ weighted_control_points

template<typename vec_t , typename param_vec_t >
NdRange<proj_vec_t> mm::NURBSPatch< vec_t, param_vec_t >::weighted_control_points
private

Multidimensional Range of control points.

Definition at line 89 of file NURBSPatch_fwd.hpp.


The documentation for this class was generated from the following files:
mm::NURBSPatch::p
std::array< int, param_dim > p
Underlying B-spline orders.
Definition: NURBSPatch_fwd.hpp:88
mm::NURBSPatch::knots
std::array< Range< scalar_t >, param_dim > knots
Array of Ranges of knots (knot vector).
Definition: NURBSPatch_fwd.hpp:90
mm::Vec3d
Vec< double, 3 > Vec3d
Convenience typedef for 3d vector of doubles.
Definition: Vec_fwd.hpp:35
mm::PI
static const double PI
Mathematical constant pi in double precision.
Definition: Config.hpp:44