34template<
short_t d,
class T>
37 GISMO_ENSURE(m_manualLevels,
"Add level only works when m_manualLevels==true");
39 std::vector<std::vector<index_t>> lvlIndices(d);
41 std::vector<T> difference, intersection;
42 std::vector<T> knots1, knots2;
43 std::vector<index_t> dirIndices;
45 for (
short_t dim=0; dim!=d; dim++)
47 dirIndices = m_uIndices.back()[dim];
49 tmpknots= tb2->knots(dim);
53 knots1 = tb2->knots(dim).unique();
54 knots2 = next_basis.knots(dim).unique();
60 std::set_intersection( knots2.begin(), knots2.end(),
61 knots1.begin(), knots1.end(),
62 std::back_inserter(intersection) );
64 std::set_difference( knots1.begin(), knots1.end(),
65 intersection.begin(), intersection.end(),
66 std::back_inserter(difference) );
67 GISMO_ASSERT(difference.size()==0,
"Knot vector is not nested!");
71 std::set_difference(knots2.begin(), knots2.end(),
72 knots1.begin(), knots1.end(),
73 std::back_inserter(difference) );
76 std::reverse(difference.begin(),difference.end());
78 std::transform(dirIndices.begin(), dirIndices.end(), dirIndices.begin(), [=](
index_t& i){return 2*i;});
81 while (difference.size()!=0)
84 typename std::vector<T>::iterator diff_ptr = difference.begin();
86 i = tmpknots.
uFind(*diff_ptr).uIndex();
89 for (n = 0, diff_ptr = difference.begin(); diff_ptr!=difference.end() && tmpknots.
uFind(*diff_ptr).uIndex()==i; n++, diff_ptr++);
91 GISMO_ENSURE(n < dirIndices[i+1] - dirIndices[i],
"Trying to insert more knots than available in the tensor structure");
97 diff_ptr = std::next(difference.begin(),(n-1)/2);
98 newIdx = (dirIndices[i] + dirIndices[i+1]) / 2;
103 diff_ptr = difference.begin();
104 newIdx = dirIndices[i] + 1;
108 dirIndices.insert(std::upper_bound( dirIndices.begin(), dirIndices.end(), newIdx ),newIdx);
110 tmpknots.
insert(*diff_ptr);
112 difference.erase(diff_ptr);
115 lvlIndices[dim] = dirIndices;
117 m_uIndices.push_back(lvlIndices);
120 m_bases.push_back( next_basis.
clone().release() );
123template<
short_t d,
class T>
126 return m_bases[0]->support();
130template<
short_t d,
class T>
134 int lvl = levelOf(i);
136 return m_bases[lvl]->support( this->flatTensorIndexOf(i) );
140template<
short_t d,
class T>
inline
143 GISMO_ASSERT(Pt.cols() == 1,
"Waiting for single point");
146 const int maxLevel = m_tree.getMaxInsLevel();
150 for(
int i =0; i < Dim; i++)
151 loIdx[i] = m_bases[maxLevel]->knots(i).uFind( Pt(i,0) ).uIndex();
154 this->_knotIndexToDiadicIndex(maxLevel,loIdx);
156 return m_tree.levelOf( loIdx, maxLevel);
159template<
short_t d,
class T>
inline
162 const int maxLevel = m_tree.getMaxInsLevel();
166 return m_tree.levelOf( Pt, maxLevel);
169template<
short_t d,
class T>
inline
174 lvl.resize( Pt.cols() );
175 loIdx.resize( Pt.rows(), Pt.cols() );
178 for(
index_t i = 0; i < Pt.cols(); i++)
180 lvl[i] = getLevelAtPoint( Pt.col(i) );
181 for(
index_t j = 0; j < Pt.rows(); j++)
182 loIdx(j,i) = m_bases[ lvl[i] ]->knots(j).uFind( Pt(j,i) ).uIndex() ;
186template<
short_t d,
class T>
inline
189 result.resize( u.cols() );
193 const int maxLevel = m_tree.getMaxInsLevel();
195 for(
index_t p = 0; p < u.cols(); p++ )
197 for(
short_t i = 0; i != d; ++i)
198 low[i] = m_bases[maxLevel]->knots(i).uFind(u(i,p)).uIndex();
201 const int lvl = m_tree.levelOf(low, maxLevel);
203 for(
int i = 0; i <= lvl; i++)
205 m_bases[i]->active_cwise(u.col(p), low, upp);
209 CMatrix::const_iterator it =
210 m_xmatrix[i].find_it_or_fail( m_bases[i]->index(cur) );
212 if( it != m_xmatrix[i].end() )
220template<
short_t d,
class T>
226 const CMatrix & cmat = m_xmatrix[lvl];
230 for (
index_t i = 0; i < d; ++i)
235 for (
index_t i = 0; i < d; ++i)
245 for (
index_t j = 0; j != end[i]; ++j)
247 if (cmat.bContains(k) && cmat.bContains(k + s))
250 const index_t kInd = m_xmatrix_offset[lvl] +
251 (std::lower_bound(cmat.begin(), cmat.end(), k)
255 const index_t kNextInd = m_xmatrix_offset[lvl] +
256 (std::lower_bound(cmat.begin(), cmat.end(), k + s)
259 mesh.addEdge(kInd, kNextInd);
267template<
short_t d,
class T>
274 for (
index_t i = 0; i< sz; ++i)
275 mesh.addVertex(nodes.row(i).transpose());
277 addConnectivity(level, mesh);
280template<
short_t d,
class T>
287 for(
index_t i = 0; i< sz; ++i )
288 mesh.addVertex( nodes.row(i).transpose() );
291 for(
size_t lvl = 0; lvl <= maxLevel(); lvl++)
293 addConnectivity(lvl, mesh);
297template<
short_t d,
class T>
300 return m_xmatrix_offset.back();
302template<
short_t d,
class T>
305 std::vector<gsSortedVector<index_t> > OX = m_xmatrix;
308 this->transfer(OX, transf);
309 gsDebug<<
"tranf orig:\n"<<transf<<std::endl;
310 coefs = transf*coefs;
313template<
short_t d,
class T>
317 this->refineElements_withTransfer(boxes,transf);
319 coefs = transf*coefs;
322template<
short_t d,
class T>
325 std::vector<gsSortedVector<index_t> > OX = m_xmatrix;
326 this->refineElements(boxes);
327 this->transfer(OX, tran);
330template<
short_t d,
class T>
331void gsHTensorBasis<d,T>::refineElements_withTransfer2(std::vector<index_t>
const & boxes, gsSparseMatrix<T> & tran)
333 std::vector<gsSortedVector<index_t> > OX = m_xmatrix;
334 this->refineElements(boxes);
335 this->transfer2(OX, tran);
338template<
short_t d,
class T>
339void gsHTensorBasis<d,T>::refineElements_withCoefs2(gsMatrix<T> & coefs,std::vector<index_t>
const & boxes)
341 std::vector<gsSortedVector<index_t> > OX = m_xmatrix;
342 refineElements(boxes);
343 gsSparseMatrix<> transf;
344 this->transfer2(OX, transf);
346 coefs = transf*coefs;
360template<
short_t d,
class T>
361void gsHTensorBasis<d,T>::unrefineElements_withCoefs(gsMatrix<T> & coefs,std::vector<index_t>
const & boxes)
363 gsSparseMatrix<> transf;
364 this->unrefineElements_withTransfer(boxes,transf);
367 typename gsSparseSolver<T>::QR solver(transf);
368 coefs=solver.solve(coefs);
371template<
short_t d,
class T>
372void gsHTensorBasis<d,T>::unrefineElements_withTransfer(std::vector<index_t>
const & boxes, gsSparseMatrix<T> & tran)
375 this->unrefineElements(boxes);
376 cp->transfer(this->m_xmatrix,tran);
379template<
short_t d,
class T>
385 std::vector<gsSortedVector<index_t> > OX = m_xmatrix;
393 std::vector<index_t> boxes;
395 for (
typename hdomain_type::literator it = m_tree.beginLeafIterator(); it.good(); it.next() )
401 lvl = it.level() + 1;
402 const point & l = it.lowerCorner();
403 const point & u = it.upperCorner();
405 boxes.push_back(lvl);
406 for(
short_t i = 0; i < d; i++)
407 boxes.push_back( l(i) * 2);
408 for(
short_t i = 0; i < d; i++)
409 boxes.push_back( u(i) * 2);
412 this->clone()->refineElements_withCoefs(coefs, boxes);
413 this->uniformRefine(numKnots, mul);
418template<
short_t d,
class T>
421 std::vector<gsSortedVector<index_t> > OX = m_xmatrix;
422 std::vector<index_t> boxes;
424 for (
typename hdomain_type::literator it = m_tree.beginLeafIterator(); it.good(); it.next() )
428 lvl = it.level() - 1;
429 const point & l = it.lowerCorner();
430 const point & u = it.upperCorner();
432 boxes.push_back(lvl);
433 for(
short_t i = 0; i < d; i++)
434 boxes.push_back( l(i) / 2);
435 for(
short_t i = 0; i < d; i++)
436 boxes.push_back( u(i)/2 + (
index_t)(u(i)%2!=0));
439 this->clone()->unrefineElements_withCoefs(coefs, boxes);
440 this->uniformCoarsen(numKnots);
444template<
short_t d,
class T>
447 GISMO_ASSERT(boxes.rows() == d,
"refine() needs d rows of boxes.");
448 GISMO_ASSERT(boxes.cols()%2 == 0,
"Each box needs two corners but you don't provide refine() with them.");
452 for(
int i = 0; i < boxes.cols()/2; i++)
454 for(
short_t j = 0; j < d; j++ )
457 "In refine() the first corner is outside the computational domain.");
459 "In refine() the second corner is outside the computational domain." );
468 this->refine( boxes );
471 needLevel( m_tree.getMaxInsLevel() );
476 std::vector<index_t> refVector = this->asElements(boxes, refExt);
479 this->refineElements( refVector );
486template<
short_t d,
class T>
489 GISMO_ASSERT(boxes.rows() == d,
"refine() needs d rows of boxes.");
490 GISMO_ASSERT(boxes.cols()%2 == 0,
"Each box needs two corners but you don't provide refine() with them.");
494 for(
int i = 0; i < boxes.cols()/2; i++)
496 for(
short_t j = 0; j < d; j++ )
499 "In refine() the first corner is outside the computational domain.");
501 "In refine() the second corner is outside the computational domain." );
507 std::vector<index_t> refVector = this->asElementsUnrefine(boxes, refExt);
510 this->unrefineElements( refVector );
540template<
short_t d,
class T>
541std::vector<index_t> gsHTensorBasis<d,T>::asElements(gsMatrix<T>
const & boxes,
int refExt)
const
548 const int offset = 2*d+1;
552 std::vector<index_t> refVector( offset * boxes.cols()/2 );
553 gsMatrix<T> ctr(d,1);
556 for(
index_t i = 0; i < boxes.cols()/2; i++)
558 ctr = ( boxes.col( 2*i ) + boxes.col( 2*i+1) )/2;
563 const int refLevel = getLevelAtPoint( ctr ) + 1;
566 needLevel( refLevel );
568 for(
index_t j = 0; j < boxes.rows();j++)
571 const gsKnotVector<T> & kv = m_bases[refLevel]->knots(j);
573 boxes(j,2*i ) ) - 1).uIndex();
575 boxes(j,2*i+1) ) - 1).uIndex();
587 _knotIndexToDiadicIndex(refLevel,j,k1);
588 _knotIndexToDiadicIndex(refLevel,j,k2);
590 _knotIndexToDiadicIndex(refLevel,j,maxKtIndex);
595 ( k1 < 2*refExt ? k1=0 : k1-=2*refExt );
596 ( k2 + 2*refExt >= maxKtIndex ? k2=maxKtIndex-1 : k2+=2*refExt);
599 refVector[i*offset] = refLevel;
600 refVector[i*offset+1+j] = k1;
601 refVector[i*offset+1+j+d] = k2;
612template<
short_t d,
class T>
613std::vector<index_t> gsHTensorBasis<d,T>::asElementsUnrefine(gsMatrix<T>
const & boxes,
int refExt)
const
620 const int offset = 2*d+1;
624 std::vector<index_t> refVector;
625 refVector.reserve( offset * boxes.cols()/2 );
626 gsMatrix<T> ctr(d,1);
632 for(
index_t i = 0; i < boxes.cols()/2; i++)
634 ctr = ( boxes.col( 2*i ) + boxes.col( 2*i+1) )/2;
639 const int refLevel = getLevelAtPoint( ctr ) - 1;
642 if (refLevel < 0)
continue;
644 refVector.resize(refVector.size() + offset);
646 for(
index_t j = 0; j < boxes.rows();j++)
649 const gsKnotVector<T> & kv = m_bases[refLevel+1]->knots(j);
651 boxes(j,2*i ) ) - 1).uIndex();
653 boxes(j,2*i+1) ) - 1).uIndex();
665 _knotIndexToDiadicIndex(refLevel,j,k1);
666 _knotIndexToDiadicIndex(refLevel,j,k2);
668 _knotIndexToDiadicIndex(refLevel,j,maxKtIndexd);
672 ( k1 < refExt ? k1=0 : k1-=refExt );
673 ( k2 + refExt >= maxKtIndexd ? k2=maxKtIndexd-1 : k2+=refExt);
677 k2 = k2/2 + (
index_t)(k2%2 != 0);
680 refVector[I*offset] = refLevel;
681 refVector[I*offset+1+j] = k1;
682 refVector[I*offset+1+j+d] = k2;
691template<
short_t d,
class T>
694 GISMO_ASSERT(boxes.rows() == d,
"refine() needs d rows of boxes.");
695 GISMO_ASSERT(boxes.cols()%2 == 0,
"Each box needs two corners but you don't provide refine() with them.");
699 for(
int i = 0; i < boxes.cols()/2; i++)
701 for(
short_t j = 0; j < d; j++ )
704 "In refine() the first corner is outside the computational domain.");
706 "In refine() the second corner is outside the computational domain." );
712 for(
index_t i = 0; i < boxes.cols()/2; i++)
715 const int fLevel = m_bases.size()-1;
717 for(
index_t j = 0; j < k1.size();j++)
721 boxes(j,2*i ) ) - 1).uIndex();
723 boxes(j,2*i+1) ) - 1).uIndex();
728 if (0!=k1[j]) {--k1[j];}
735 _knotIndexToDiadicIndex(fLevel,k1);
736 _knotIndexToDiadicIndex(fLevel,k2);
746 m_tree.sinkBox(k1, k2, fLevel);
748 needLevel( m_tree.getMaxInsLevel() );
811template<
short_t d,
class T>
815 const index_t lvl = this->levelOf(i);
818 m_bases[lvl]->elementSupport_into(m_xmatrix[lvl][ i - m_xmatrix_offset[lvl] ],
820 point low = elements.col(0);
821 point upp = elements.col(1);
823 for (
short_t k = 0; k!=d; ++k )
825 low[k] = low[k] << 1;
826 upp[k] = upp[k] << 1;
829 m_tree.insertBox(low,upp,lvl+1);
843template<
short_t d,
class T>
850 "The points did not define boxes properly. The boxes were not added to the basis.");
851 for(
size_t i = 0; i < (boxes.size())/(2*d+1); i++)
853 for(
short_t j = 0; j < d; j++ )
855 i1[j] = boxes[(i*(2*d+1))+j+1];
856 i2[j] = boxes[(i*(2*d+1))+d+j+1];
858 insert_box(i1,i2,boxes[i*(2*d+1)]);
864template<
short_t d,
class T>
871 "The points did not define boxes properly. The boxes were not added to the basis.");
873 for(
size_t i = 0; i < (boxes.size())/(2*d+1); i++)
875 for(
short_t j = 0; j < d; j++ )
877 i1[j] = boxes[(i*(2*d+1))+j+1];
878 i2[j] = boxes[(i*(2*d+1))+d+j+1];
881 m_tree.clearBox(i1,i2, boxes[i*(2*d+1)]);
887 auto leafIt = m_tree.beginLeafIterator();
888 for (; leafIt.good(); leafIt.next())
890 if ( leafIt.level()>0 )
892 leafIt.upperCorner(), leafIt.level() );
897 m_tree.computeMaxInsLevel();
901template<
short_t d,
class T>
904 const index_t dir = side.direction();
905 const index_t par = side.parameter();
907 rf(dir,!par) = rf(dir,par);
908 for (
index_t i = 0; i!=lvl; ++i)
912template<
short_t d,
class T>
919 if(
const Self_t * _other =
dynamic_cast<const Self_t*
>( &other) )
930 bndThis = this->boundaryOffset( bi.
first().side(), offset );
934 bndOther= _other->boundaryOffset( bi.
second().side(), offset );
936 "Input error, sizes do not match: "
937 <<bndThis.rows()<<
"!="<<bndOther.rows() );
941 for(
index_t i=0; i < bndThis.rows(); i++)
944 index_t L = this->levelOf( bndThis(i,0) );
947 index_t flat0 = this->flatTensorIndexOf( bndThis(i,0) );
949 tens0 = this->tensorLevel(L).tensorIndex( flat0 );
961 N[j] = _other->tensorLevel(L).component(j).size();
971 tens1[jj] = tens0[j];
979 tens1[jj] = N[jj]-(1+offset);
988 tens1[jj] = N[jj]-1 - tens1[jj];
993 flat1 = _other->tensorLevel(L).index( tens1 );
997 cont1 = _other->flatTensorIndexToHierachicalIndex( flat1, L );
999 bndOther( i, 0 ) = cont1;
1003 gsWarn<<
"Cannot match with "<< other <<
"\n";
1006template<
short_t d,
class T>
1012 matchWith(bi,other,bndThis,bndOther,0);
1021template<
short_t d,
class T>
1022void gsHTensorBasis<d,T>::set_activ1(
int level)
1029 CMatrix & cmat = m_xmatrix[level];
1036 if ( level >
static_cast<int>(m_tree.getMaxInsLevel() ) )
1040 gsVector<knotIter,d> starts, ends, curr;
1044 for(
short_t i = 0; i != d; ++i)
1047 starts[i] = m_bases[level]->knots(i).sbegin() ;
1049 ends [i] = m_bases[level]->knots(i).send()-m_deg[i]-1;
1055 for(
short_t i = 0; i != d; ++i)
1057 low[i] = curr[i].uIndex();
1058 upp[i] = (curr[i]+m_deg[i]+1).uIndex();
1059 ind[i] = curr[i].index();
1065 _knotIndexToDiadicIndex(level,low);
1066 _knotIndexToDiadicIndex(level,upp);
1069 if ( m_tree.query3(low, upp,level) == level)
1070 cmat.push_unsorted( m_bases[level]->index( ind ) );
1079template<
short_t d,
class T>
1081 const int level, point & actLow, point & actUpp)
1084 for(
short_t i = 0; i != d; ++i)
1086 actLow[i] = tb.knots(i).lastKnotIndex (boxLow[i]) - m_deg[i];
1087 actUpp[i] = tb.knots(i).firstKnotIndex(boxUpp[i]) - 1 ;
1095template<
short_t d,
class T>
1113 for (
typename hdomain_type::literator it = m_tree.beginLeafIterator();
1114 it.good(); it.next() )
1116 const int lvl = it.level();
1117 CMatrix & cmat = m_xmatrix[lvl];
1120 functionOverlap(it.lowerCorner(), it.upperCorner(), lvl, curr, actUpp);
1124 const index_t gi = m_bases[lvl]->index( curr );
1127 m_bases[lvl]->elementSupport_into(gi, elSupp);
1129 if ( (elSupp.col(0).array() >= it.lowerCorner().array()).all() &&
1130 (elSupp.col(1).array() <= it.upperCorner().array()).all() )
1133 cmat.push_unsorted( gi );
1138 if ( m_tree.query3(elSupp.col(0), elSupp.col(1), lvl) == lvl)
1139 cmat.push_unsorted( gi );
1145 for(
size_t lvl = 0; lvl != m_xmatrix.size(); ++lvl)
1147 m_xmatrix[lvl].sort();
1148 m_xmatrix[lvl].erase( std::unique( m_xmatrix[lvl].begin(), m_xmatrix[lvl].end() ),
1149 m_xmatrix[lvl].end() );
1153template<
short_t d,
class T>
1155 std::vector<CMatrix> & x_matrix_lvl)
const
1157 x_matrix_lvl.resize(level+1);
1164 for(
int j =0; j < level+1; j++)
1167 x_matrix_lvl[j].clear();
1169 for(
index_t i = 0; i != d; ++i)
1172 starts[i] = m_bases[j]->knots(i).sbegin() ;
1174 ends [i] = m_bases[j]->knots(i).send()-m_deg[i]-1;
1180 for(
index_t i = 0; i != d; ++i)
1182 low[i] = curr[i].uIndex();
1183 upp[i] = (curr[i]+m_deg[i]+1).uIndex();
1184 ind[i] = curr[i].index();
1188 if ( m_tree.query3(low, upp,j) == j)
1190 x_matrix_lvl[j].push_unsorted( m_bases[j]->index( ind ) );
1193 if ( m_tree.query3(low, upp,j) >= j)
1195 x_matrix_lvl[j].push_unsorted( m_bases[j]->index( ind ) );
1201 x_matrix_lvl[j].sort();
1207template<
short_t d,
class T>
inline
1209 typename gsHTensorBasis<d,T>::point
const & k2,
1215 m_tree.insertBox(k1,k2, lvl);
1216 needLevel( m_tree.getMaxInsLevel() );
1219template<
short_t d,
class T>
1225 while ( ! m_xmatrix_offset[1] )
1227 delete m_bases.front();
1228 m_bases.erase( m_bases.begin() );
1229 m_tree.decrementLevel();
1230 m_xmatrix.erase( m_xmatrix.begin() );
1231 m_xmatrix_offset.erase( m_xmatrix_offset.begin() );
1236template<
short_t d,
class T>
1239 GISMO_ASSERT(
static_cast<size_t>(level) < m_xmatrix.size(),
"Requested level does not exist.\n");
1241 CMatrix::const_iterator xmat_pointer = m_xmatrix[level].begin();
1242 CMatrix::const_iterator xmat_end = m_xmatrix[level].end();
1243 gsSortedVector< int >::iterator ind_pointer = indexes.begin();
1244 gsSortedVector< int >::iterator ind_end = indexes.end();
1246 while(ind_pointer!=ind_end&&xmat_pointer!=xmat_end)
1248 if(*ind_pointer<
static_cast<int>(*xmat_pointer))
1253 else if(*ind_pointer==
static_cast<int>(*xmat_pointer))
1255 (*ind_pointer)=m_xmatrix_offset[level]+index;
1266 while(ind_pointer!=ind_end)
1273template<
short_t d,
class T>
1276 if( m_xmatrix.size()<=
static_cast<size_t>(level) )
1278 CMatrix::const_iterator it = std::lower_bound(m_xmatrix[level].begin(), m_xmatrix[level].end(), index );
1279 if(it == m_xmatrix[level].end() || *it != index)
1282 return m_xmatrix_offset[level] + ( it - m_xmatrix[level].begin() );
1285template<
short_t d,
class T>
1291 const index_t sz = bound.rows();
1294 indexes.resize(sz,-1);
1295 if(level<=maxLevel())
1298 indexes[i]=(bound)(i,0);
1299 flatTensorIndexesToHierachicalIndexes(indexes,level);
1301 actives.resize(indexes.size(),
false);
1302 std::fill (actives.begin(),actives.end(),
false);
1303 for(
size_t i = 0;i<indexes.size();i++)
1308template<
short_t d,
class T>
1312 needLevel( m_tree.getMaxInsLevel() );
1316 m_xmatrix.resize( m_tree.getMaxInsLevel()+1 );
1318 m_tree.makeCompressed();
1320 for(
size_t i = 0; i != m_xmatrix.size(); i ++)
1327 m_xmatrix_offset.clear();
1328 m_xmatrix_offset.reserve(m_xmatrix.size()+1);
1329 m_xmatrix_offset.push_back(0);
1330 for (
size_t i = 0; i != m_xmatrix.size(); i++)
1332 m_xmatrix_offset.push_back(
1333 m_xmatrix_offset.back() + m_xmatrix[i].size() );
1337template<
short_t d,
class T>
1340 GISMO_ENSURE(!m_manualLevels || (
size_t)(maxLevel)<m_uIndices.size(),
"Maximum manual level reached, maxLevel = "<<maxLevel<<
", m_uIndices.size() = "<<m_uIndices.size());
1342 const int extraLevels = maxLevel + 1 - m_bases.size();
1343 for (
int i = 0; i < extraLevels; ++i )
1347 m_bases.push_back (next_basis);
1351template<
short_t d,
class T>
1357 for(
index_t i = 0; i < d; i++)
1358 m_deg[i] = tbasis.
degree(i);
1361 if (
const tensorBasis * tb2 =
1362 dynamic_cast<const tensorBasis*
>(&tbasis) )
1364 m_bases.push_back(tb2->clone().release());
1368 GISMO_ERROR(
"Cannot construct a Hierarchical basis from "<< tbasis );
1373 for (
index_t i = 0; i!=d; ++i )
1374 upp[i] = m_bases[0]->knots(i).numElements();
1383template<
short_t d,
class T>
1387 point low, upp, cur;
1388 const int maxLevel = m_tree.getMaxInsLevel();
1390 std::vector<std::vector<index_t> > temp_output;
1391 temp_output.resize( u.cols() );
1396 for(
index_t p = 0; p < u.cols(); p++)
1398 currPoint = u.col(p);
1399 for(
short_t i = 0; i != d; ++i)
1400 low[i] = m_bases[maxLevel]->knots(i).uFind( currPoint(i,0) ).uIndex();
1403 this->_knotIndexToDiadicIndex(maxLevel,low);
1406 const int lvl = m_tree.levelOf(low, maxLevel);
1407 for(
int i = 0; i <= lvl; i++)
1418 m_bases[i]->active_cwise(currPoint, low, upp);
1422 CMatrix::const_iterator it =
1423 m_xmatrix[i].find_it_or_fail( m_bases[i]->index(cur) );
1425 if( it != m_xmatrix[i].end() )
1427 temp_output[p].push_back(
1428 this->m_xmatrix_offset[i] + (it - m_xmatrix[i].begin() )
1436 if ( temp_output[p].size() > sz )
1437 sz = temp_output[p].size();
1440 result.resize(sz, u.cols() );
1441 for(
index_t i = 0; i < result.cols(); i++)
1443 result.col(i).topRows(temp_output[i].size())
1445 result.col(i).bottomRows(sz-temp_output[i].size()).setZero();
1449template<
short_t d,
class T>
1452 std::vector<index_t> temp;
1454 for(
size_t i = 0; i != m_xmatrix[i].size(); i++)
1455 for (CMatrix::const_iterator it = m_xmatrix[i].begin();
1456 it != m_xmatrix[i].end(); it++)
1458 ind = this->m_bases[i]->tensorIndex(*it);
1459 for (
unsigned j=0; j!=d; ++j )
1460 if ( (ind[j]==0) || (ind[j]==(this->m_bases[i]->size(j)-1)) )
1462 temp.push_back(m_xmatrix_offset[i] + (it-m_xmatrix[i].begin()) );
1466 return makeMatrix<index_t>(temp.begin(),temp.size(),1 );
1469template<
short_t d,
class T>
1477 std::vector<index_t> temp;
1480 GISMO_ASSERT(this->maxLevel() < this->m_bases.size(),
"Something went wrong: maxLevel() < m_bases.size(), "<<this->maxLevel()<<
" < "<<m_bases.size());
1481 needLevel(m_xmatrix.size()-1);
1483 for(
size_t i = 0; i != m_xmatrix.size(); i++)
1485 GISMO_ASSERT(
static_cast<int>(offset)<this->m_bases[i]->size(k),
1486 "Offset ("<<offset<<
") cannot be bigger than the amount of basis"
1487 "functions orthogonal to Boxside s! ("<<this->m_bases[i]->size(k)<<
")");
1489 index_t r = ( par ? this->m_bases[i]->size(k) - 1 -offset : offset);
1490 for (CMatrix::const_iterator it = m_xmatrix[i].begin();
1491 it != m_xmatrix[i].end(); it++)
1493 ind = this->m_bases[i]->tensorIndex(*it);
1496 m_xmatrix_offset[i] + (it - m_xmatrix[i].begin() )
1500 return makeMatrix<index_t>(temp.begin(),temp.size(),1 );
1503template<
short_t d,
class T>
1513template<
short_t d,
class T>
1515levelAtCorner(boxCorner
const & c)
const
1518 gsVector<bool> pars;
1519 c.parameters_into(d,pars);
1521 gsMatrix<T> supp = m_bases.front()->support();
1523 gsVector<T> vec(supp.rows());
1524 for (
index_t r = 0; r!=supp.rows(); r++)
1525 vec(r) = supp(r,pars(r));
1527 return getLevelAtPoint(vec);
1530template<
short_t d,
class T>
1532functionAtCorner(boxCorner
const & c)
const
1534 index_t lvl = this->levelAtCorner(c);
1537 index_t index = m_bases[lvl]->functionAtCorner(c);
1540 return this->flatTensorIndexToHierachicalIndex(index,lvl);
1543template<
short_t d,
class T>
1545functionAtCorner(boxCorner
const & c,
index_t level)
const
1548 index_t index = m_bases[level]->functionAtCorner(c);
1551 return this->flatTensorIndexToHierachicalIndex(index,level);
1564template<
short_t d,
class T>
1565void gsHTensorBasis<d,T>::uniformRefine(
int numKnots,
int mul,
int dir)
1568 GISMO_ASSERT(numKnots == 1,
"Only implemented for numKnots = 1");
1570 GISMO_ASSERT( m_tree.getMaxInsLevel() <
static_cast<unsigned>(m_bases.size()),
1571 "Problem with max inserted levels: "<< m_tree.getMaxInsLevel()
1572 <<
"<" << m_bases.size() <<
"\n");
1575 tensorBasis * last_basis = m_bases.back()->clone().release();
1576 last_basis->uniformRefine(1,mul,dir);
1577 m_bases.push_back( last_basis );
1580 delete m_bases.front();
1581 m_bases.erase( m_bases.begin() );
1584 m_tree.multiplyByTwo();
1589template<
short_t d,
class T>
1593 GISMO_ASSERT(numKnots == 1,
"Only implemented for numKnots = 1");
1597 m_bases.insert( m_bases.begin(), first_basis );
1600 delete m_bases.back();
1604 m_tree.divideByTwo();
1611template<
short_t d,
class T>
1614 std::vector< std::vector< std::vector< std::vector<index_t > > > > dummy;
1615 return domainBoundariesGeneric( dummy, result,
false );
1618template<
short_t d,
class T>
1621 std::vector< std::vector< std::vector< std::vector< T > > > > dummy;
1622 return domainBoundariesGeneric( result, dummy,
true );
1626template<
short_t d,
class T>
1628 std::vector< std::vector< std::vector< std::vector< T > > > >& params,
1629 bool indicesFlag )
const
1633 std::vector< std::vector< std::vector< int > > > res_aabb;
1634 std::vector< std::vector< std::vector<index_t > > > res_aabb_unsigned;
1635 std::vector< std::vector< std::vector< std::vector< index_t > > > > polylines;
1637 polylines = this->m_tree.getPolylines();
1638 res_aabb.resize( polylines.size() );
1645 indices.resize( polylines.size() );
1648 params.resize(polylines.size());
1650 int maxLevel =
static_cast<int>( this->maxLevel() );
1653 std::vector<T> x_dir(m_bases[maxLevel]->knots(0).unique());
1654 std::vector<T> y_dir(m_bases[maxLevel]->knots(1).unique());
1656 for(
unsigned int i0 = 0; i0 < polylines.size(); i0++)
1659 indices[i0].resize( polylines[i0].size() );
1661 params[i0].resize( polylines[i0].size() );
1663 res_aabb[i0].resize( polylines[i0].size() );
1664 for(
unsigned int i1 = 0; i1 < polylines[i0].size(); i1++)
1667 indices[i0][i1].resize( polylines[i0][i1].size() );
1669 params[i0][i1].resize( polylines[i0][i1].size() );
1671 res_aabb[i0][i1].resize( 4 );
1672 res_aabb[i0][i1][0] = 1000000;
1673 res_aabb[i0][i1][1] = 1000000;
1674 res_aabb[i0][i1][2] = -10000000;
1675 res_aabb[i0][i1][3] = -10000000;
1676 for(
unsigned int i2 = 0; i2 < polylines[i0][i1].size(); i2++)
1680 indices[i0][i1][i2].resize(4);
1681 indices[i0][i1][i2][0] = polylines[i0][i1][i2][0];
1682 indices[i0][i1][i2][1] = polylines[i0][i1][i2][1];
1683 indices[i0][i1][i2][2] = polylines[i0][i1][i2][2];
1684 indices[i0][i1][i2][3] = polylines[i0][i1][i2][3];
1689 params[i0][i1][i2].resize(4);
1691 params[i0][i1][i2][0] = (x_dir[polylines[i0][i1][i2][0]]);
1692 params[i0][i1][i2][1] = (y_dir[polylines[i0][i1][i2][1]]);
1693 params[i0][i1][i2][2] = (x_dir[polylines[i0][i1][i2][2]]);
1694 params[i0][i1][i2][3] = (y_dir[polylines[i0][i1][i2][3]]);
1696 if(res_aabb[i0][i1][0]>
signed(polylines[i0][i1][i2][0]))
1698 res_aabb[i0][i1][0] = polylines[i0][i1][i2][0];
1700 if(res_aabb[i0][i1][1]>
signed(polylines[i0][i1][i2][1]))
1702 res_aabb[i0][i1][1] = polylines[i0][i1][i2][1];
1704 if(res_aabb[i0][i1][2]<
signed(polylines[i0][i1][i2][2]))
1706 res_aabb[i0][i1][2] = polylines[i0][i1][i2][2];
1708 if(res_aabb[i0][i1][3]<
signed(polylines[i0][i1][i2][3]))
1710 res_aabb[i0][i1][3] = polylines[i0][i1][i2][3];
1716 res_aabb_unsigned.resize(res_aabb.size());
1717 for (
unsigned int i = 0; i < res_aabb.size(); i++)
1719 res_aabb_unsigned[i].resize(res_aabb[i].size());
1720 for (
unsigned int j = 0; j < res_aabb[i].size(); j++)
1722 res_aabb_unsigned[i][j].resize(res_aabb[i][j].size());
1723 for (
unsigned int k = 0; k < res_aabb[i][j].size(); k++)
1725 if(res_aabb[i][j][k]<0)
1726 gsWarn <<
"conversion form signed to unsigned\n";
1727 res_aabb_unsigned[i][j][k] = res_aabb[i][j][k];
1731 return res_aabb_unsigned;
1735template<
short_t d,
class T>
1739 needLevel( old.size() );
1743 std::vector< gsSparseMatrix<T,RowMajor> > transfer(m_bases.size()-1);
1744 std::vector<std::vector<T> > knots(d);
1746 for(
size_t i = 1; i < m_bases.size(); ++i)
1749 for(
unsigned dim = 0; dim != d; ++dim)
1766 while ( old.size() >= m_xmatrix.size() )
1769 result = this->coarsening_direct(old, m_xmatrix, transfer);
1776 while( m_xmatrix.back().size() == 0 )
1777 m_xmatrix.pop_back();
1780 const int sizeDiff =
static_cast<int>( m_bases.size() - m_xmatrix.size() );
1783 freeAll(m_bases.end() - sizeDiff, m_bases.end());
1784 m_bases.resize(m_xmatrix.size());
1787 result.makeCompressed();
1790template<
short_t d,
class T>
1794 needLevel( old.size() );
1796 tensorBasis T_0_copy = this->tensorLevel(0);
1797 std::vector< gsSparseMatrix<T,RowMajor> > transfer( m_bases.size()-1 );
1798 std::vector<std::vector<T> > knots(d);
1800 for(
size_t i = 1; i < m_bases.size(); ++i)
1803 for(
short_t dim = 0; dim != d; ++dim)
1815 T_0_copy.refine_withTransfer(transfer[i-1], knots);
1819 while ( old.size() >= m_xmatrix.size())
1820 m_xmatrix.push_back( gsSortedVector<index_t>() );
1822 result = this->coarsening_direct2(old, m_xmatrix, transfer);
1824 result.makeCompressed();
1828template<
short_t d,
class T>
1831 GISMO_ASSERT(
static_cast<size_t>(lvl) < m_xmatrix.size(),
1832 "Requested level does not exist.\n");
1834 if (m_bases[lvl]->knots(dir).has(knotValue))
1836 for(
unsigned int i =lvl;i < m_bases.size();i++)
1837 m_bases[i]->component(dir).insertKnot(knotValue,mult);
1841 gsWarn<<
"Knot value not in the given knot vector."<<std::endl;
1848template<
short_t d,
class T>
1851 for(
unsigned int k =0; k < knotValue.size(); ++k)
1853 if (m_bases[lvl]->knots(dir).has(knotValue[k]))
1855 for(
unsigned int i =lvl;i < m_bases.size(); ++i)
1856 m_bases[i]->component(dir).insertKnot(knotValue[k],mult);
1860 gsWarn<<
"Knot value not in the given knot vector."<<std::endl;
1866template<
short_t d,
class T>
1871 m_tree.getBoxesInLevelIndex(b1,b2,level);
1873 for(
index_t i = 0;i<level.rows();i++)
1878 const index_t par_index = m_bases[l]->knots(dir).uFind(par).uIndex();
1879 if(l>0 && (par_index>=min(dir)) && (par_index<=max(dir)))
1882 for(
index_t j=0;j<min.rows();++j)
1884 boxes.push_back(min(j));
1885 for(
index_t j=0;j<max.rows();++j)
1887 boxes.push_back(max(j));
1892template<
short_t d,
class T>
1895 GISMO_ASSERT(m_manualLevels,
"Only works for manual levels");
1896 knotIndex = m_uIndices[level][dir][knotIndex];
1899template<
short_t d,
class T>
1902 for (
index_t r = 0; r != d; r++)
1903 this->_knotIndexToDiadicIndex(level,r,knotIndex[r]);
1906template<
short_t d,
class T>
1909 GISMO_ASSERT(m_manualLevels,
"Only works for manual levels");
1910 typename std::vector<index_t>::const_iterator it = std::find_if(m_uIndices[level][dir].begin(),m_uIndices[level][dir].end(),[&diadicIndex](
const index_t i) {
return i >= diadicIndex; });
1911 GISMO_ASSERT(it!=m_uIndices[level][dir].end(),
"Index not found");
1912 diadicIndex = std::distance(m_uIndices[level][dir].begin(),it);
1915template<
short_t d,
class T>
1918 for (
index_t r = 0; r != d; r++)
1919 this->_diadicIndexToKnotIndex(level,r,diadicIndex[r]);
1923template<
short_t d,
class T>
1926 for (
size_t level=0;level<m_bases.size();++level)
1927 m_bases[level]->degreeElevate(i,dir);
1929 for(
unsigned c=0; c<d; ++c)
1930 m_deg[c]=m_bases[0]->degree(c);
1932 this->update_structure();
1935template<
short_t d,
class T>
1936void gsHTensorBasis<d,T>::degreeReduce(
int const & i,
int const dir)
1938 for (
size_t level=0;level<m_bases.size();++level)
1939 m_bases[level]->degreeReduce(i,dir);
1941 for(
unsigned c=0; c<d; ++c)
1942 m_deg[c]=m_bases[0]->degree(c);
1944 this->update_structure();
1947template<
short_t d,
class T>
1948void gsHTensorBasis<d,T>::degreeIncrease(
int const & i,
int const dir)
1950 for (
size_t level=0;level<m_bases.size();++level)
1951 m_bases[level]->degreeIncrease(i,dir);
1953 for(
unsigned c=0; c<d; ++c)
1954 m_deg[c]=m_bases[0]->degree(c);
1956 this->update_structure();
1959template<
short_t d,
class T>
1960void gsHTensorBasis<d,T>::degreeDecrease(
int const & i,
int const dir)
1962 for (
size_t level=0;level<m_bases.size();++level)
1963 m_bases[level]->degreeDecrease(i,dir);
1965 for(
unsigned c=0; c<d; ++c)
1966 m_deg[c]=m_bases[0]->degree(c);
1968 this->update_structure();
1971template<
short_t d,
class T>
1975 np.setConstant(npts);
1977 gsMatrix<T> grid = gsPointGrid<T>(supp.col(0),supp.col(1),np);
1979 this->eval_into(grid,res);
1981 return ((sums.array()>1-tol && sums.array()<1+tol).all());
1989template<
short_t d,
class T>
1995 GSXML_COMMON_FUNCTIONS(gsHTensorBasis<TMPLA2(d,T)>);
1996 static std::string tag () {
return "Basis"; }
1997 static std::string type () {
return ""; }
1999 static gsHTensorBasis<d,T> * get (gsXmlNode * node)
2001 gsXmlAttribute * btype = node->first_attribute(
"type");
2004 gsWarn<<
"Basis without a type in the xml file.\n";
2007 std::string s = btype->value() ;
2008 if ( s.compare(0, 9,
"HBSplineB" , 9 ) == 0 )
2009 return gsXml< gsHBSplineBasis<d,T> >::get(node);
2010 if ( s.compare(0, 10,
"THBSplineB", 10) == 0 )
2011 return gsXml< gsTHBSplineBasis<d,T> >::get(node);
2013 gsWarn<<
"gsXmlUtils: gsHTensorBasis: No known basis \""<<s<<
"\". Error.\n";
2017 static gsXmlNode * put (
const gsHTensorBasis<d,T> & obj,
2020 const gsBasis<T> * ptr = & obj;
2023 if (
const gsHBSplineBasis<d,T> * g =
2024 dynamic_cast<const gsHBSplineBasis<d,T> *
>( ptr ) )
2025 return gsXml< gsHBSplineBasis<d,T> >::put(*g,data);
2028 if (
const gsTHBSplineBasis<d,T> * g =
2029 dynamic_cast<const gsTHBSplineBasis<d,T> *
>( ptr ) )
2030 return gsXml< gsTHBSplineBasis<d,T> >::put(*g,data);
2032 gsWarn<<
"gsXmlUtils put: getBasis: No known basis \""<<obj<<
"\". Error.\n";
Struct which represents a certain side of a box.
Definition gsBoundary.h:85
bool parameter() const
Returns the parameter value (false=0=start, true=1=end) that corresponds to this side.
Definition gsBoundary.h:128
short_t direction() const
Returns the parametric direction orthogonal to this side.
Definition gsBoundary.h:113
Creates a mapped object or data pointer to a const vector without copying data.
Definition gsAsMatrix.h:285
A basis represents a family of scalar basis functions defined over a common parameter domain.
Definition gsBasis.h:79
virtual short_t degree(short_t i) const
Degree with respect to the i-th variable. If the basis is a tensor product of (piecewise) polynomial ...
Definition gsBasis.hpp:699
uPtr clone()
Clone methode. Produceds a deep copy inside a uPtr.
Class with a hierarchical domain structure represented by a box k-d-tree.
Definition gsHDomain.h:76
void insertBox(point const &lower, point const &upper, node *_node, int lvl)
The insert function which insert box defined by points lower and upper to level lvl.
Definition gsHDomain.hpp:143
Class representing a (scalar) hierarchical tensor basis of functions .
Definition gsHTensorBasis.h:75
memory::unique_ptr< gsHTensorBasis > uPtr
Unique pointer for gsHTensorBasis.
Definition gsHTensorBasis.h:81
index_t getLevelAtPoint(const gsMatrix< T > &Pt) const
Returns the level(s) at point(s) in the parameter domain.
Definition gsHTensorBasis.hpp:141
virtual void uniformCoarsen(int numKnots=1)
Coarsen the basis uniformly by removing groups of numKnots consecutive knots, each knot removed mul t...
Definition gsHTensorBasis.hpp:1590
void _knotIndexToDiadicIndex(const index_t level, const index_t dir, index_t &knotIndex) const
Transfers the knotIndex in the knot span in direction dir on level level to diadic indices.
Definition gsHTensorBasis.hpp:1893
virtual void unrefineElements(std::vector< index_t > const &boxes)
Clear the given boxes into the quadtree.
Definition gsHTensorBasis.hpp:865
void numActive_into(const gsMatrix< T > &u, gsVector< index_t > &result) const
The number of active basis functions at points u.
Definition gsHTensorBasis.hpp:187
void refineSide(const boxSide side, index_t lvl)
Refines all the cells on the side side up to level lvl.
Definition gsHTensorBasis.hpp:902
void flatTensorIndexesToHierachicalIndexes(gsSortedVector< int > &indexes, const int level) const
transformes a sortedVector indexes of flat tensor index of the bspline basis of level to hierachical ...
Definition gsHTensorBasis.hpp:1237
gsMatrix< index_t > allBoundary() const
Returns the indices of the basis functions that are nonzero at the domain boundary.
Definition gsHTensorBasis.hpp:1450
void getBoxesAlongSlice(int dir, T par, std::vector< index_t > &boxes) const
Definition gsHTensorBasis.hpp:1867
virtual void refine(gsMatrix< T > const &boxes, int refExt)
Refine the basis to levels and in the areas defined by boxes with an extension.
Definition gsHTensorBasis.hpp:445
std::vector< std::vector< std::vector< index_t > > > domainBoundariesIndices(std::vector< std::vector< std::vector< std::vector< index_t > > > > &result) const
Gives polylines on the boundaries between different levels of the mesh.
Definition gsHTensorBasis.hpp:1619
virtual void connectivity(const gsMatrix< T > &nodes, gsMesh< T > &mesh) const
Definition gsHTensorBasis.hpp:281
void activeBoundaryFunctionsOfLevel(const unsigned level, const boxSide &s, std::vector< bool > &actives) const
Definition gsHTensorBasis.hpp:1286
void uniformRefine_withCoefs(gsMatrix< T > &coefs, int numKnots=1, int mul=1, short_t const dir=-1)
Refine the basis uniformly.
Definition gsHTensorBasis.hpp:380
void makeCompressed()
Cleans the basis, removing any inactive levels.
Definition gsHTensorBasis.hpp:1220
void transfer(const std::vector< gsSortedVector< index_t > > &old, gsSparseMatrix< T > &result)
Returns transfer matrix between the hirarchical spline given by the characteristic matrix "old" and t...
Definition gsHTensorBasis.hpp:1736
void _diadicIndexToKnotIndex(const index_t level, gsVector< index_t, d > &diadicIndex) const
Transfers the diadicIndex in the knot span in direction on level level to knot indices.
Definition gsHTensorBasis.hpp:1916
void addLevel(const gsTensorBSplineBasis< d, T > &next_basis)
Adds a level, only if manual levels are activated.
Definition gsHTensorBasis.hpp:35
index_t size() const
The number of basis functions in this basis.
Definition gsHTensorBasis.hpp:298
bool testPartitionOfUnity(const index_t npts=100, const T tol=1e-12) const
Test the partition of unity.
Definition gsHTensorBasis.hpp:1972
virtual void increaseMultiplicity(index_t lvl, int dir, T knotValue, int mult=1)
Increases the multiplicity of a knot with the value knotValue in level lvl in direction dir by mult....
Definition gsHTensorBasis.hpp:1829
virtual void active_into(const gsMatrix< T > &u, gsMatrix< index_t > &result) const
Returns the indices of active basis functions at points u, as a list of indices, in result....
Definition gsHTensorBasis.hpp:1384
void insert_box(point const &k1, point const &k2, int lvl)
Inserts a domain into the basis.
Definition gsHTensorBasis.hpp:1208
void functionOverlap(const point &boxLow, const point &boxUpp, const int level, point &actLow, point &actUpp)
Returns the basis functions of level which have support on box, represented as an index box.
Definition gsHTensorBasis.hpp:1080
void refineBasisFunction(const index_t i)
Refines the basis function with (hierarchical) index i.
Definition gsHTensorBasis.hpp:812
int flatTensorIndexToHierachicalIndex(index_t index, const int level) const
takes a flat tensor index of the bspline basis of level and gives back the hierachical index....
Definition gsHTensorBasis.hpp:1274
virtual gsMatrix< index_t > boundaryOffset(boxSide const &s, index_t offset) const
Definition gsHTensorBasis.hpp:1471
std::vector< std::vector< std::vector< index_t > > > domainBoundariesGeneric(std::vector< std::vector< std::vector< std::vector< index_t > > > > &indices, std::vector< std::vector< std::vector< std::vector< T > > > > ¶ms, bool indicesFlag) const
Implementation of the features common to domainBoundariesParams and domainBoundariesIndices....
Definition gsHTensorBasis.hpp:1627
index_t getLevelAtIndex(const point &Pt) const
Returns the level(s) at indexes in the parameter domain.
Definition gsHTensorBasis.hpp:160
void uniformCoarsen_withCoefs(gsMatrix< T > &coefs, int numKnots=1)
Coarsen the basis uniformly.
Definition gsHTensorBasis.hpp:419
gsMatrix< T > support() const
Returns the boundary basis for side s.
Definition gsHTensorBasis.hpp:124
void needLevel(int maxLevel) const
Makes sure that there are numLevels grids computed in the hierarachy.
Definition gsHTensorBasis.hpp:1338
void refineElements_withCoefs(gsMatrix< T > &coefs, std::vector< index_t > const &boxes)
Definition gsHTensorBasis.hpp:314
void setActiveToLvl(int level, std::vector< CMatrix > &x_matrix_lvl) const
Creates characteristic matrices for basis where "level" is the maximum level i.e. ignoring higher lev...
Definition gsHTensorBasis.hpp:1154
virtual void update_structure()
Updates the basis structure (eg. charact. matrices, etc), to be called after any modifications.
Definition gsHTensorBasis.hpp:1309
void getLevelUniqueSpanAtPoints(const gsMatrix< T > &Pt, gsVector< index_t > &lvl, gsMatrix< index_t > &loIdx) const
Returns the level(s) and knot span(s) at point(s) in the parameter domain.
Definition gsHTensorBasis.hpp:170
virtual void refineElements(std::vector< index_t > const &boxes)
Insert the given boxes into the quadtree.
Definition gsHTensorBasis.hpp:844
std::vector< std::vector< std::vector< index_t > > > domainBoundariesParams(std::vector< std::vector< std::vector< std::vector< T > > > > &result) const
Gives polylines on the boundaries between different levels of the mesh.
Definition gsHTensorBasis.hpp:1612
Class for representing a knot vector.
Definition gsKnotVector.h:80
size_t size() const
Number of knots (including repetitions).
Definition gsKnotVector.h:242
mult_t maxInteriorMultiplicity() const
Returns the maximum multiplicity in the interior.
Definition gsKnotVector.hpp:724
internal::gsKnotIterator< T > smart_iterator
Definition gsKnotVector.h:107
int degree() const
Returns the degree of the knot vector.
Definition gsKnotVector.hpp:700
uiterator domainUBegin() const
Returns a unique iterator pointing to the starting knot of the domain.
Definition gsKnotVector.h:359
void insert(T knot, mult_t mult=1)
Inserts knot knot into the knot vector with multiplicity mult.
Definition gsKnotVector.hpp:325
uiterator domainUEnd() const
Returns a unique iterator pointing to the ending knot of the domain.
Definition gsKnotVector.h:367
uiterator uFind(const T u) const
Returns the uiterator pointing to the knot at the beginning of the knot interval containing u....
Definition gsKnotVector.hpp:747
void symDifference(const gsKnotVector< T > &other, std::vector< T > &result) const
Computes the symmetric difference between this knot-vector and other.
Definition gsKnotVector.hpp:172
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
Class Representing a triangle mesh with 3D vertices.
Definition gsMesh.h:32
This class is derived from std::vector, and adds sort tracking.
Definition gsSortedVector.h:110
Sparse matrix class, based on gsEigen::SparseMatrix.
Definition gsSparseMatrix.h:139
A tensor product B-spline basis.
Definition gsTensorBSplineBasis.h:37
void refine_withTransfer(gsSparseMatrix< T, RowMajor > &transfer, const std::vector< std::vector< T > > &refineKnots)
Takes a vector of coordinate wise knot values and inserts these values to the basis.
Definition gsTensorBSplineBasis.hpp:66
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
virtual void uniformCoarsen(int numKnots=1)
Coarsen the basis uniformly by removing groups of numKnots consecutive knots, each knot removed mul t...
Definition gsTensorBasis.h:349
unsigned stride(short_t dir) const
Returns the stride for dimension dir.
Definition gsTensorBasis.h:889
virtual void uniformRefine(int numKnots=1, int mul=1, short_t dir=-1)
Refine the basis uniformly by inserting numKnots new knots with multiplicity mul on each knot span.
Definition gsTensorBasis.h:315
unsigned index(unsigned i, unsigned j, unsigned k=0) const
Definition gsTensorBasis.h:882
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
bool nextLexicographic(Vec &cur, const Vec &size)
Iterates through a tensor lattice with the given size. Updates cur and returns true if another entry ...
Definition gsCombinatorics.h:196
bool nextCubePoint(Vec &cur, const Vec &end)
Iterate in lexigographic order through the points of the integer lattice contained in the cube [0,...
Definition gsCombinatorics.h:327
#define short_t
Definition gsConfig.h:35
#define index_t
Definition gsConfig.h:32
#define gsDebug
Definition gsDebug.h:61
#define GISMO_ERROR(message)
Definition gsDebug.h:118
#define gsWarn
Definition gsDebug.h:50
#define GISMO_UNUSED(x)
Definition gsDebug.h:112
#define GISMO_ENSURE(cond, message)
Definition gsDebug.h:102
#define GISMO_ASSERT(cond, message)
Definition gsDebug.h:89
Provides declaration of the Mesh class.
Provides functions to generate structured point data.
Provides implementation of generic XML functions.
Provides declaration of input/output XML utilities struct.
The G+Smo namespace, containing all definitions for the library.
void freeAll(It begin, It end)
Frees all pointers in the range [begin end)
Definition gsMemory.h:312
Struct which represents an interface between two patches.
Definition gsBoundary.h:650
patchSide & first()
first, returns the first patchSide of this interface
Definition gsBoundary.h:776
patchSide & second()
second, returns the second patchSide of this interface
Definition gsBoundary.h:782