39template<
short_t d,
class T>
66 typedef memory::shared_ptr< gsTensorNurbs >
Ptr;
69 typedef memory::unique_ptr< gsTensorNurbs >
uPtr;
84 GISMO_ASSERT(d==2,
"Wrong dimension: tried to make a "<< d
85 <<
"D NURBS using 2 knot-vectors.");
92 this->
m_basis =
new Basis(tbasis) ;
95 "Coefficient matrix for the NURBS does not have "
96 "the expected number of control points (rows)." );
103 GISMO_ASSERT(d==2,
"Wrong dimension: tried to make a "<< d
104 <<
"D NURBS using 2 knot-vectors.");
112 "Coefficient matrix for the NURBS does not have "
113 "the expected number of control points (rows)." );
127 GISMO_ASSERT(d==3,
"Wrong dimension: tried to make a "<< d
128 <<
"D NURBS using 3 knot-vectors.");
136 "Coefficient matrix for the NURBS does not have "
137 "the expected number of control points (rows)." );
150 GISMO_ASSERT(d==3,
"Wrong dimension: tried to make a "<< d
151 <<
"D NURBS using 3 knot-vectors.");
159 "Coefficient matrix for the NURBS does not have "
160 "the expected number of control points (rows)." );
162 this->
m_basis =
new Basis(tbasis) ;
175 GISMO_ASSERT(d==4,
"Wrong dimension: tried to make a "<< d
176 <<
"D NURBS using 4 knot-vectors.");
178 std::vector< gsBSplineBasis<T>* > cbases;
187 "Coefficient matrix for the NURBS does not have "
188 "the expected number of control points (rows)." );
194 GISMO_BASIS_ACCESSORS
205 std::ostream &
print(std::ostream &os)
const
206 { os <<
"Tensor-NURBS geometry "<<
"R^"<< this->
parDim() <<
207 " --> R^"<< this->
geoDim()<<
", #control pnts= "<< this->
coefsSize() <<
": "
210 << this->
basis().weights().at(0) <<
" ... "
220 const KnotVectorType &
knots(
const int i)
const
221 {
return this->
basis().source().knots(i); }
223 KnotVectorType &
knots(
const int i)
224 {
return this->
basis().source().knots(i); }
230 GISMO_ASSERT( dir >= 0 &&
static_cast<unsigned>(dir) < d,
231 "Invalid basis component "<< dir <<
" requested for degree elevation" );
237 const index_t n = projectiveCoefs.cols();
245 const index_t nc = sz.template tail<static_cast<short_t>(d-1)>().prod();
246 projectiveCoefs.resize( sz[0], n * nc );
249 gsBoehm(tbs.knots(dir), projectiveCoefs , knot, i,
true );
250 sz[0] = projectiveCoefs.rows();
253 const index_t ncoef = sz.prod();
254 projectiveCoefs.resize(ncoef, n );
273 {
return this->
basis().source().component(i).degree(); }
287 str[0] = tbsbasis.
source().stride( k );
288 str[1] = tbsbasis.
source().stride(!k );
292 for (
int i=0; i< sz[0]; i++ )
293 for (
int j=0; j< sz[1]/2; j++ )
295 this->
m_coefs.row(i*str[0] + j*str[1] ).swap(
296 this->
m_coefs.row(i*str[0] + (sz[1]-j-1)*str[1] )
299 w.row(i*str[0] + j*str[1] ).swap(
300 w.row(i*str[0] + (sz[1]-j-1)*str[1] )
311 GISMO_ASSERT(d-1>=0,
"d must be greater or equal than 1");
312 GISMO_ASSERT(dir_fixed>=0 &&
static_cast<unsigned>(dir_fixed)<d,
"cannot fix a dir greater than dim or smaller than 0");
328 const int mult = this->
basis().knots(dir_fixed).multiplicity(par);
329 const int degree = this->
basis().degree(dir_fixed);
338 this->
basis().size_cwise(sizes);
339 constructCoefsForSlice<d, T>(dir_fixed, index, this->
coefs(), sizes,
coefs);
347 this->
basis().stride_cwise(intStrides);
349 clone->basis().knots(dir_fixed),
clone->coefs(),par,dir_fixed,
350 intStrides.template cast<unsigned>(),
degree-mult,
true);
356 clone->basis().size_cwise(sizes);
357 constructCoefsForSlice<d, T>(dir_fixed, index,
clone->coefs(), sizes,
coefs);
368 void swapDirections(
const unsigned i,
const unsigned j)
371 this->
basis().size_cwise(sz);
373 this->
basis().swapDirections(i,j);
386 std::vector<gsGeometry<T>*> result
389 for (
typename std::vector<
gsGeometry<T>*>::iterator it = result.begin();
390 it != result.end(); ++it )
422 tmp.
splitAt(dir, xi, leftBS, rightBS);
430 leftNurbs->
m_basis = leftBasis;
439 rightNurbs->
m_basis = rightBasis;
Struct which represents a certain side of a box.
Definition gsBoundary.h:85
A univariate B-spline basis.
Definition gsBSplineBasis.h:700
virtual const gsMatrix< T > & weights() const
Only for compatibility reasons, with gsRationalBasis. It returns an empty matrix.
Definition gsBasis.h:411
virtual const gsBasis & source() const
Definition gsBasis.h:710
uPtr clone()
Clone methode. Produceds a deep copy inside a uPtr.
Abstract base class representing a geometry map.
Definition gsGeometry.h:93
gsMatrix< T > & coefs()
Definition gsGeometry.h:340
gsMatrix< T >::RowXpr coef(index_t i)
Returns the i-th coefficient of the geometry as a row expression.
Definition gsGeometry.h:346
void eval_into(const gsMatrix< T > &u, gsMatrix< T > &result) const
Evaluate the function at points u into result.
Definition gsGeometry.hpp:166
short_t geoDim() const
Dimension n of the absent physical space.
Definition gsGeometry.h:292
virtual const gsBasis< T > & basis() const =0
Returns a const reference to the basis of the geometry.
short_t parDim() const
Dimension d of the parameter domain (same as domainDim()).
Definition gsGeometry.hpp:190
gsBasis< T > * m_basis
Pointer to the basis of this geometry.
Definition gsGeometry.h:632
gsMatrix< T > m_coefs
Coefficient matrix of size coefsSize() x geoDim()
Definition gsGeometry.h:629
index_t coefsSize() const
Return the number of coefficients (control points)
Definition gsGeometry.h:371
Class for representing a knot vector.
Definition gsKnotVector.h:80
int degree() const
Returns the degree of the knot vector.
Definition gsKnotVector.hpp:700
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
static void setFromProjectiveCoefs(const gsMatrix< T > &pr_coefs, gsMatrix< T > &coefs, gsMatrix< T > &weights)
Definition gsRationalBasis.h:356
const gsMatrix< T > & weights() const
Returns the weights of the rational basis.
Definition gsRationalBasis.h:300
A tensor product B-spline basis.
Definition gsTensorBSplineBasis.h:37
const Basis_t & component(short_t dir) const
For a tensor product basis, return the (const) 1-d basis for the i-th parameter component.
Definition gsTensorBSplineBasis.h:202
A tensor product of d B-spline functions, with arbitrary target dimension.
Definition gsTensorBSpline.h:45
void splitAt(index_t dir, T xi, gsTensorBSpline< d, T > &left, gsTensorBSpline< d, T > &right) const
Definition gsTensorBSpline.hpp:449
std::vector< gsGeometry< T > * > uniformSplit(index_t dir=-1) const
Definition gsTensorBSpline.hpp:386
index_t size() const
Returns the number of elements in the basis.
Definition gsTensorBasis.h:108
void size_cwise(gsVector< index_t, s > &result) const
The number of basis functions in the direction of the k-th parameter component.
Definition gsTensorBasis.h:441
A tensor product Non-Uniform Rational B-spline (NURBS) basis.
Definition gsTensorNurbsBasis.h:38
memory::unique_ptr< gsTensorNurbsBasis > uPtr
Unique pointer for gsTensorNurbsBasis.
Definition gsTensorNurbsBasis.h:70
A tensor product Non-Uniform Rational B-spline function (NURBS) of parametric dimension d,...
Definition gsTensorNurbs.h:41
void insertKnot(T knot, int dir, int i=1)
Inserts knot knot at direction dir, i times.
Definition gsTensorNurbs.h:227
gsTensorNurbs(gsKnotVector< T > const &KV1, gsKnotVector< T > const &KV2, gsKnotVector< T > const &KV3, gsMatrix< T > tcoefs)
Definition gsTensorNurbs.h:145
gsTensorNurbs(gsKnotVector< T > const &KV1, gsKnotVector< T > const &KV2, gsKnotVector< T > const &KV3, gsMatrix< T > tcoefs, gsMatrix< T > wgts)
Definition gsTensorNurbs.h:121
memory::unique_ptr< gsTensorNurbs > uPtr
Unique pointer for gsTensorNurbs.
Definition gsTensorNurbs.h:69
gsTensorNurbs(gsKnotVector< T > const &KV1, gsKnotVector< T > const &KV2, gsMatrix< T > tcoefs, gsMatrix< T > wgts)
Construct 2D tensor NURBS by knot vectors, degrees, weights and coefficient matrix.
Definition gsTensorNurbs.h:100
void reverse(unsigned k)
Toggle orientation wrt coordinate k.
Definition gsTensorNurbs.h:276
gsBSplineTraits< static_cast< short_t >(d-1), T >::RatGeometry BoundaryGeometryType
Associated boundary geometry type.
Definition gsTensorNurbs.h:59
void splitAt(index_t dir, T xi, gsTensorNurbs< d, T > &left, gsTensorNurbs< d, T > &right) const
Definition gsTensorNurbs.h:411
const gsMatrix< T > & weights() const
Returns the NURBS weights.
Definition gsTensorNurbs.h:266
std::vector< gsGeometry< T > * > uniformSplit(index_t dir=-1) const
Definition gsTensorNurbs.h:379
const KnotVectorType & knots(const int i) const
Returns a reference to the knot vector i.
Definition gsTensorNurbs.h:220
gsTensorNurbs()
Default empty constructor.
Definition gsTensorNurbs.h:74
gsTensorNurbs(gsKnotVector< T > const &KV1, gsKnotVector< T > const &KV2, gsMatrix< T > tcoefs)
Definition gsTensorNurbs.h:81
T & weight(int i) const
Access to i-th weight.
Definition gsTensorNurbs.h:263
void slice(index_t dir_fixed, T par, BoundaryGeometryType &result) const
Definition gsTensorNurbs.h:309
gsBSplineTraits< static_cast< short_t >(d-1), T >::RatBasis BoundaryBasisType
Associated boundary basis type.
Definition gsTensorNurbs.h:63
gsBSplineBasis< T > Family_t
Family type.
Definition gsTensorNurbs.h:53
gsMatrix< T > & weights()
Returns the NURBS weights as non-const reference.
Definition gsTensorNurbs.h:269
std::ostream & print(std::ostream &os) const
Prints the object as a string.
Definition gsTensorNurbs.h:205
gsTensorNurbs(gsKnotVector< T > const &KV1, gsKnotVector< T > const &KV2, gsKnotVector< T > const &KV3, gsKnotVector< T > const &KV4, gsMatrix< T > tcoefs, gsMatrix< T > wgts)
Definition gsTensorNurbs.h:168
memory::shared_ptr< gsTensorNurbs > Ptr
Shared pointer for gsTensorNurbs.
Definition gsTensorNurbs.h:66
short_t degree(unsigned i) const
Returns the degree of the basis wrt direction i.
Definition gsTensorNurbs.h:272
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
void gsTensorBoehm(KnotVectorType &knots, Mat &coefs, T val, int direction, gsVector< unsigned > str, int r=1, bool update_knots=true)
Definition gsBoehm.hpp:245
void swapTensorDirection(int k1, int k2, gsVector< index_t, d > &sz, gsMatrix< T > &coefs)
Definition gsTensorTools.h:129
#define short_t
Definition gsConfig.h:35
#define index_t
Definition gsConfig.h:32
#define GISMO_ASSERT(cond, message)
Definition gsDebug.h:89
Provides declaration of Geometry abstract interface.
Provides declaration of TensorNurbsBasis abstract interface.
The G+Smo namespace, containing all definitions for the library.
void gsBoehm(KnotVectorType &knots, Mat &coefs, T val, int r=1, bool update_knots=true)
Performs insertion of multiple knot on "knots" and coefficients "coefs".
Definition gsBoehm.hpp:29
S give(S &x)
Definition gsMemory.h:266
Traits for BSplineBasis in more dimensions.
Definition gsBSplineBasis.h:32