G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsRationalTHBSpline.h
Go to the documentation of this file.
1
14#pragma once
15
16#include <gsCore/gsGeometry.h>
19
20
21namespace gismo
22{
23
38template<short_t d, class T>
39class gsRationalTHBSpline : public gsGeoTraits<d,T>::GeometryBase
40{
41
42public:
43 typedef gsKnotVector<T> KnotVectorType;
44
45 typedef typename gsGeoTraits<d,T>::GeometryBase Base;
46
47 typedef T Scalar_t;
48
49 typedef gsTensorBSplineBasis<d,T> TBasis; // underlying tensor basis
50
53
54 // rational version of tensor basis (basis for this geometry)
56
58 // typedef typename gsBSplineTraits<static_cast<short_t>(d-1),T>::RatGeometry BoundaryGeometryType;
59 // typedef gsRationalTHBSpline<static_cast<short_t>(d-1),T> BoundaryGeometryType;
60 typedef typename
61 util::conditional<d==1, gsConstantFunction<T>, gsRationalTHBSpline<static_cast<short_t>(d-1),T>
63
66
68 typedef memory::shared_ptr< gsRationalTHBSpline > Ptr;
69
71 typedef memory::unique_ptr< gsRationalTHBSpline > uPtr;
72
73public:
74
77
78 gsRationalTHBSpline( const Basis & basis, gsMatrix<T> coefs ) :
79 Base( basis, give(coefs) ) { }
80
81 // gsRationalTHBSpline( gsTensorNurbs<d> nurbs): Base(new gsRationalTHBSplineBasis<d>(nurbs.basis().source(), nurbs.weights() ), nurbs.coefs() ){ }
83 explicit gsRationalTHBSpline( const gsTensorNurbs<d,T> & nurbs )
84 {
85 // gsTHBSplineBasis<d> thbBasis( nurbs.basis().source().clone() );
86 // this->m_basis = new Basis( *thbBasis.clone(), nurbs.weights() );
87 this->m_basis = new Basis( new gsTHBSplineBasis<d>( nurbs.basis().source() ), nurbs.weights() );
88 this->m_coefs = nurbs.coefs();
89 }
90
91
92
93 GISMO_BASIS_ACCESSORS
94
95 public:
96
97// ***********************************************
98// Virtual member functions required by the base class
99// ***********************************************
100
101 GISMO_CLONE_FUNCTION(gsRationalTHBSpline)
102
103
104 std::ostream &print(std::ostream &os) const
105 { os << "Tensor-NURBS geometry "<< "R^"<< this->parDim() <<
106 " --> R^"<< this->geoDim()<< ", #control pnts= "<< this->coefsSize() <<": "
107 << this->coef(0) <<" ... "<< this->coef(this->coefsSize()-1);
108 os << "\nweights: "
109 << this->basis().weights().at(0) <<" ... "
110 << this->basis().weights().at(this->coefsSize()-1)
111 <<"\n" ;
112 return os; }
113
114// ***********************************************
115// Additional members
116// ***********************************************
117
118 // /// Inserts knot \a knot at direction \a dir, \a i times
119 // void insertKnot( T knot, int dir, int i = 1)
120 // {
121 // GISMO_ASSERT( i>0, "multiplicity must be at least 1");
122 // GISMO_ASSERT( dir >= 0 && static_cast<unsigned>(dir) < d,
123 // "Invalid basis component "<< dir <<" requested for degree elevation" );
124
125 // // Combine control points and weights to n+1 dimensional "control points"
126 // gsMatrix<T> projectiveCoefs = basis().projectiveCoefs(m_coefs, weights());
127
128 // // Dimension of physical space + 1 for the weights
129 // const index_t n = projectiveCoefs.cols();
130
131 // Basis & tbs = this->basis().source();
132 // gsVector<index_t,d> sz;
133 // tbs.size_cwise(sz); // Size of basis per (parametric) direction
134
135 // swapTensorDirection(0, dir, sz, projectiveCoefs );
136
137 // const index_t nc = sz.template tail<static_cast<short_t>(d-1)>().prod();
138 // projectiveCoefs.resize( sz[0], n * nc );
139
140 // // Insert knot and update projective control points
141 // gsBoehm(tbs.knots(dir), projectiveCoefs , knot, i, true );
142 // sz[0] = projectiveCoefs.rows();
143
144 // // New number of coefficients
145 // const index_t ncoef = sz.prod();
146 // projectiveCoefs.resize(ncoef, n );
147
148 // swapTensorDirection(0, dir, sz, projectiveCoefs );
149
150 // // Unpack projective control points and update the coefs and weights of the original NURBS
151 // basis().setFromProjectiveCoefs(projectiveCoefs, m_coefs, weights() );
152 // }
153
155 T & weight(int i) const { return this->basis().weight(i); }
156
158 const gsMatrix<T> & weights() const { return this->basis().weights(); }
159
161 gsMatrix<T> & weights() { return this->basis().weights(); }
162
163 // /// Returns the degree of the basis wrt direction i
164 // short_t degree(unsigned i) const
165 // { return this->basis().source().component(i).degree(); }
166
167protected:
168 // todo: check function: check the coefficient number, degree, knot vector ...
169
170 using gsGeometry<T>::m_coefs;
171 using gsGeometry<T>::m_basis;
172
173}; // class gsRationalTHBSpline
174
175
176// ***********************************************
177// ***********************************************
178
179
180} // namespace gismo
A univariate B-spline basis.
Definition gsBSplineBasis.h:700
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
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
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
A rational Truncated Hierarchical B-Spline basis.
Definition gsRationalTHBSplineBasis.h:37
A rational truncated hierarchical B-Spline function of parametric dimension d, with arbitrary target ...
Definition gsRationalTHBSpline.h:40
memory::unique_ptr< gsRationalTHBSpline > uPtr
Unique pointer for gsRationalTHBSpline.
Definition gsRationalTHBSpline.h:71
util::conditional< d==1, gsConstantFunction< T >, gsRationalTHBSpline< static_cast< short_t >(d-1), T > >::type BoundaryGeometryType
Associated boundary geometry type.
Definition gsRationalTHBSpline.h:62
gsRationalTHBSpline()
Default empty constructor.
Definition gsRationalTHBSpline.h:76
const gsMatrix< T > & weights() const
Returns the NURBS weights.
Definition gsRationalTHBSpline.h:158
gsRationalTHBSplineBasis< d-1, T > BoundaryBasisType
Associated boundary basis type.
Definition gsRationalTHBSpline.h:65
T & weight(int i) const
Access to i-th weight.
Definition gsRationalTHBSpline.h:155
gsBSplineBasis< T > Family_t
Family type.
Definition gsRationalTHBSpline.h:52
gsMatrix< T > & weights()
Returns the NURBS weights as non-const reference.
Definition gsRationalTHBSpline.h:161
std::ostream & print(std::ostream &os) const
Prints the object as a string.
Definition gsRationalTHBSpline.h:104
gsRationalTHBSpline(const gsTensorNurbs< d, T > &nurbs)
Construct B-Spline from a Tensor B-Spline.
Definition gsRationalTHBSpline.h:83
memory::shared_ptr< gsRationalTHBSpline > Ptr
Shared pointer for gsRationalTHBSpline.
Definition gsRationalTHBSpline.h:68
Truncated hierarchical B-spline basis.
Definition gsTHBSplineBasis.h:36
A tensor product B-spline basis.
Definition gsTensorBSplineBasis.h:37
A tensor product Non-Uniform Rational B-spline function (NURBS) of parametric dimension d,...
Definition gsTensorNurbs.h:41
const gsMatrix< T > & weights() const
Returns the NURBS weights.
Definition gsTensorNurbs.h:266
Provides forward declarations of types and structs.
Provides declaration of Geometry abstract interface.
Provides declaration of RationalTHBSplineBasis abstract interface.
The G+Smo namespace, containing all definitions for the library.
S give(S &x)
Definition gsMemory.h:266