30template<
short_t d, 
class T>
 
   42        const index_t shift = ( par ? -1 : 1);
 
   44        std::vector<index_t> temp;
 
   48        for(
size_t i = 0; i != m_xmatrix.size(); i++)
 
   50            GISMO_ASSERT(
static_cast<int>(offset)<this->m_bases[i]->size(k),
 
   51                         "Offset cannot be bigger than the amount of basis" 
   52                         "functions orthogonal to Boxside s!");
 
   54            index_t r = ( par ? this->m_bases[i]->size(k) - 1 : 0);
 
   55            for (
typename CMatrix::const_iterator it = m_xmatrix[i].begin();
 
   56                 it != m_xmatrix[i].end(); it++)
 
   58                ind = this->m_bases[i]->tensorIndex(*it);
 
   62                    hi = this->flatTensorIndexToHierachicalIndex(this->m_bases[i]->index(ind),i);
 
   64                    GISMO_ASSERT(hi!=-1,
"Neightbouring basis function with coordinates "<<ind.transpose()<<
" of level "<<i<<
" does not exist.");
 
   70        return makeMatrix<index_t>(temp.begin(),temp.size(),1 );
 
 
   75template<
short_t d, 
class T>
 
   78    GISMO_ASSERT(d-1>=0,
"d must be greater or equal than 1");
 
   79    GISMO_ASSERT(dir_fixed>=0 && 
static_cast<index_t>(dir_fixed)<d,
"cannot fix a dir greater than dim or smaller than 0");
 
   80    const boxSide side(dir_fixed,0);
 
   82        this->m_bases[0]->boundaryBasis(side);
 
   88        std::vector<index_t> boxes;
 
   89        this->getBoxesAlongSlice(dir_fixed,par,boxes);
 
 
   95template<
short_t d, 
class T>
 
   99    this->m_is_truncated.resize(this->size());
 
  100    m_presentation.clear();
 
  104    for (
index_t j = 0; j < this->size(); ++j)
 
  106        index_t level = this->levelOf(j);
 
  107        index_t tensor_index = this->flatTensorIndexOf(j, level);
 
  110        this->m_bases[level]->elementSupport_into(tensor_index, element_ind);
 
  113        low = element_ind.col(0); 
 
  114        high = element_ind.col(1); 
 
  117            this->_knotIndexToDiadicIndex(level,low);
 
  118            this->_knotIndexToDiadicIndex(level,high);
 
  125        index_t clevel = this->m_tree.query4(low, high, level);
 
  129            this->m_tree.computeFinestIndex(low, level, low);
 
  130            this->m_tree.computeFinestIndex(high, level, high);
 
  132            this->m_is_truncated[j] = clevel;
 
  133            _representBasisFunction(j, clevel, low, high);
 
  137            this->m_is_truncated[j] = -1;
 
 
  142template<
short_t d, 
class T>
 
  145    const unsigned pres_level,
 
  149    const unsigned cur_level = this->levelOf(j);
 
  153    act_size_of_coefs.fill(1);
 
  156    unsigned nmb_of_coefs = _updateSizeOfCoefs(cur_level, pres_level,
 
  157                                               finest_low, finest_high,
 
  161    coefs.row(0).setOnes();
 
  166    vec_nmb_of_coefs.fill(1);
 
  168    unsigned tensor_index = this->flatTensorIndexOf(j, cur_level);
 
  172        this->m_bases[cur_level]->tensorIndex(tensor_index);
 
  176    std::vector<gsKnotVector<T> > vector_of_kv(d);
 
  180    cur_size_of_coefs.fill(1);
 
  182    for (
unsigned level = cur_level; level < pres_level; ++level)
 
  184        _updateSizeOfCoefs(level, level + 1, finest_low,
 
  185                           finest_high, cur_size_of_coefs);
 
  192        this->m_tree.computeLevelIndex(finest_low, level, clow);
 
  193        this->m_tree.computeLevelIndex(finest_high, level, chigh);
 
  196            this->_diadicIndexToKnotIndex(level,clow);
 
  197            this->_diadicIndexToKnotIndex(level,chigh);
 
  199        this->m_tree.computeLevelIndex(finest_low, level + 1, flow);
 
  200        this->m_tree.computeLevelIndex(finest_high, level + 1, fhigh);
 
  203            this->_diadicIndexToKnotIndex(level + 1,flow);
 
  204            this->_diadicIndexToKnotIndex(level + 1,fhigh);
 
  207        std::vector<T> knots;
 
  209        for (
unsigned dim = 0; dim < d; ++dim)
 
  214            if (level == cur_level)
 
  215                vector_of_kv[dim] = ckv;
 
  218            std::set_symmetric_difference(ckv.
beginAt(clow[dim]), ckv.
endAt(chigh[dim]),
 
  220                                          std::back_inserter(knots));
 
  225                                     typename std::vector<T>::const_iterator>
 
  226                (vector_of_kv[dim], bspl_vec_ti[dim], coefs, vec_nmb_of_coefs,
 
  228                 cur_size_of_coefs, dim, knots.begin(), knots.end(),
 
  232        _truncate(coefs, act_size_of_coefs, cur_size_of_coefs,
 
  233                  level + 1, bspl_vec_ti, cur_level, finest_low);
 
  235    _saveNewBasisFunPresentation(coefs, act_size_of_coefs,
 
  236                                 j, pres_level, finest_low);
 
 
  240template<
short_t d, 
class T>
 
  245    const unsigned pres_level,
 
  248    const unsigned level = this->levelOf(j);
 
  249    const unsigned tensor_index = this->flatTensorIndexOf(j, level);
 
  252        this->m_bases[level]->tensorIndex(tensor_index);
 
  255    const unsigned f_ten_index = _basisFunIndexOnLevel(bspl_vec_ti, level,
 
  256                                                       finest_low, pres_level);
 
  259    bspline::buildCoeffsStrides<d>(act_size_of_coefs, act_coefs_strides);
 
  268    bspline::getLastIndexLocal<d>(act_size_of_coefs, last_point);
 
  270    this->m_presentation[j] =
 
  282        unsigned ten_index = f_ten_index;
 
  283        for (
unsigned dim = 0; dim < d; dim++)
 
  285            ten_index += position(dim) *
 
  286                this->m_bases[pres_level]->stride(
static_cast<int>(dim));
 
  289        unsigned coef_index = bspline::getIndex<d>(act_coefs_strides, position);
 
  291        if (coefs(coef_index) != 0)
 
  292            presentation(ten_index) = coefs(coef_index);
 
 
  299template<
short_t d, 
class T>
 
  302    const unsigned level,
 
  304    const unsigned new_level)
 
  307    this->m_tree.computeLevelIndex(fin_low, level, low);
 
  310    this->m_tree.computeLevelIndex(fin_low, new_level, flow);
 
  314        this->_diadicIndexToKnotIndex(level,low);
 
  315        this->_diadicIndexToKnotIndex(new_level,flow);
 
  320    for (
unsigned dim = 0; dim < d; dim++)
 
  323            this->m_bases[level]->knots(dim);
 
  326            this->m_bases[new_level]->knots(dim);
 
  334    return this->m_bases[new_level]->index(new_index);
 
 
  338template<
short_t d, 
class T>
 
  343    const unsigned level,
 
  345    const unsigned bspl_vec_ti_level,
 
  349    if (this->m_xmatrix[level].size() == 0)
 
  353    const unsigned const_ten_index = _basisFunIndexOnLevel(bspl_vec_ti,
 
  354                                                           bspl_vec_ti_level, finest_low, level);
 
  356    bspline::buildCoeffsStrides<d>(act_size_of_coefs, act_coefs_strides);
 
  360    bspline::getLastIndexLocal<d>(size_of_coefs, last_point);
 
  369    unsigned xmatrix_index = 0;
 
  370    unsigned tensor_active_index = this->m_xmatrix[level][0];
 
  372    unsigned numb_of_point = size_of_coefs[0];
 
  383        unsigned ten_index = const_ten_index;
 
  384        for (
unsigned dim = 1; dim < d; dim++)
 
  386            ten_index += position(dim) *
 
  387                this->m_bases[level]->stride(
static_cast<int>(dim));
 
  390        unsigned coef_index = bspline::getIndex<d>(act_coefs_strides, position);
 
  392        for (
unsigned index = 0; index < numb_of_point; ++index)
 
  394            if (tensor_active_index < ten_index)
 
  396                while(tensor_active_index < ten_index)
 
  399                    if (xmatrix_index == this->m_xmatrix[level].size())
 
  407                    tensor_active_index =
 
  408                        this->m_xmatrix[level][xmatrix_index];
 
  412            if (ten_index == tensor_active_index) 
 
  413                coefs(coef_index + index, 0) = 0;
 
 
  423template<
short_t d, 
class T>
 
  425    const unsigned clevel,
 
  426    const unsigned flevel,
 
  433    this->m_tree.computeLevelIndex(finest_low, clevel, clow);
 
  434    this->m_tree.computeLevelIndex(finest_high, clevel, chigh);
 
  437    this->m_tree.computeLevelIndex(finest_low, flevel, flow);
 
  438    this->m_tree.computeLevelIndex(finest_high, flevel, fhigh);
 
  442        this->_diadicIndexToKnotIndex(clevel,clow);
 
  443        this->_diadicIndexToKnotIndex(clevel,chigh);
 
  444        this->_diadicIndexToKnotIndex(flevel,flow);
 
  445        this->_diadicIndexToKnotIndex(flevel,fhigh);
 
  448    unsigned nmb_of_coefs = 1;
 
  450    for (
unsigned dim = 0; dim < d; ++dim)
 
  453            this->m_bases[clevel]->knots(dim);
 
  455            this->m_bases[flevel]->knots(dim);
 
  463        size_of_coefs(dim) += fnmb_knts - cnmb_knts;
 
  464        nmb_of_coefs *= size_of_coefs(dim);
 
 
  470template<
short_t d, 
class T>
 
  472typename util::enable_if<dd==2,void>::type
 
  482    const unsigned loc2glob = ( 1<< (this->maxLevel() - level) );
 
  483    if( b1[0]%loc2glob != 0 ) b1[0] -= b1[0]%loc2glob;
 
  484    if( b1[1]%loc2glob != 0 ) b1[1] -= b1[1]%loc2glob;
 
  485    if( b2[0]%loc2glob != 0 ) b2[0] += loc2glob -(b2[0]%loc2glob);
 
  486    if( b2[1]%loc2glob != 0 ) b2[1] += loc2glob -(b2[1]%loc2glob);
 
  490    this->m_tree.computeLevelIndex( b1, level, b1_outputs );
 
  491    this->m_tree.computeLevelIndex( b2, level, b2_outputs );
 
  492    int i0 = b1_outputs(0);
 
  493    int i1 = b2_outputs(0);
 
  494    int j0 = b1_outputs(1);
 
  495    int j1 = b2_outputs(1);
 
  496    i0 = m_bases[level]->knots(0).lastKnotIndex(i0) - m_deg[0];
 
  497    i1 = m_bases[level]->knots(0).firstKnotIndex(i1) - 1;
 
  498    j0 = m_bases[level]->knots(1).lastKnotIndex(j0) - m_deg[1];
 
  499    j1 = m_bases[level]->knots(1).firstKnotIndex(j1) - 1;
 
  501    const index_t sz0   = m_bases[level]->size(0);
 
  502    const index_t newSz = (i1 - i0 + 1)*(j1 - j0 + 1);
 
  503    cp.resize(newSz, geom_coef.cols());
 
  506    globalRefinement(geom_coef, level, temp);
 
  509    for(
int j = j0; j <= j1; j++)
 
  510        for(
int k = i0; k <= i1; k++)
 
  511            cp.row(cc++) = temp.row(j*sz0+k);
 
  515                         m_bases[level]->knots(0).begin() + i1 + m_deg[0] + 2);
 
  517                         m_bases[level]->knots(1).begin() + j1 + m_deg[1] + 2);
 
  521template<
short_t d, 
class T>
 
  526    this->m_tree.getBoxes(b1,b2,level); 
 
  527    int nboxes = level.size();
 
  532    p1.resize(this->dim());
 
  533    p2.resize(this->dim());
 
  536    nvertices.resize(nboxes,this->dim());
 
  538    for (
int i = 0; i < nboxes; i++)
 
  540        p1 = b1.row(i).transpose();
 
  541        p2 = b2.row(i).transpose();
 
  543        this->getBsplinePatchGlobal(p1, p2, level[i], geom_coef, temp1, cku, ckv);
 
  551            int cprows = cp.rows();
 
  552            temp2.resize(cp.rows()+temp1.rows(),cp.cols());
 
  554            for (
int j = 0; j < cp.rows(); j++)
 
  556                for (
int k = 0; k < cp.cols(); k++)
 
  558                    temp2(j,k) = cp(j,k);
 
  562            for (
int j = 0; j < temp1.rows(); j++)
 
  564                for (
int k = 0; k < temp1.cols(); k++)
 
  566                    temp2(cprows+j,k) = temp1(j,k);
 
 
  578template<
short_t d, 
class T>
 
  587    this->m_tree.getBoxes(b1,b2,level); 
 
  589    const int nboxes = level.size();
 
  594    for (
int i = 0; i < nboxes; i++) 
 
  596        p1 = b1.row(i).transpose();
 
  597        p2 = b2.row(i).transpose();
 
  599        this->getBsplinePatchGlobal(p1, p2, level[i], geom_coef, temp1, cku, ckv);
 
 
  606template<
short_t d, 
class T>
 
  608    std::vector<std::vector<std::vector< std::vector<index_t> > > >& connectedComponents, 
gsVector<index_t>& level)
 const 
  612    for(
unsigned int i = 0; i < this->m_xmatrix.size(); i++)
 
  614        if(this->m_xmatrix[i].size()>0)
 
  620    gsDebug<<
"new min level"<<
"\n";
 
  621    std::vector< std::vector< std::vector< std::vector<index_t> > > > res; 
 
  622    std::vector< std::vector< std::vector<index_t > > > aabb;
 
  623    std::vector< std::vector<index_t > > boxes;
 
  624    aabb = this->domainBoundariesIndices(res);
 
  625    for(
unsigned int i = 0; i < aabb.size(); i++)
 
  628        for(
unsigned int j = 0; j < aabb[i].size(); j++)
 
  630            bool is_boundary = 
true;
 
  631            for(
unsigned int k = 0; k < aabb[i].size(); k++)
 
  635                    if(     (aabb[i][j][0] > aabb[i][k][0]) &&
 
  636                            (aabb[i][j][2] < aabb[i][k][2]) &&
 
  637                            (aabb[i][j][1] > aabb[i][k][1]) &&
 
  638                            (aabb[i][j][3] < aabb[i][k][3]) )
 
  640                        is_boundary= !is_boundary;
 
  646                boxes.push_back(aabb[i][j]);
 
  647                boxes[boxes.size()-1].push_back(i+first_level);
 
  648                connectedComponents.resize(connectedComponents.size()+1);
 
  649                connectedComponents[connectedComponents.size()-1].push_back(res[i][j]);
 
  655    level.resize(boxes.size());
 
  656    for(
size_t i = 0; i < boxes.size(); i++){
 
  657        level[i] = boxes[i][2*d];
 
  661    for(
size_t l = 0; l < aabb.size();l++) 
 
  663        for(
size_t i = 0; i < aabb[l].size(); i++) 
 
  665            int closest_box = -1;
 
  666            for(
size_t j = 0; j < boxes.size(); j++)
 
  668                if(l == 
static_cast<size_t>(boxes[j][4]) )
 
  670                    if(     (aabb[l][i][0] > boxes[j][0]) &&
 
  671                            (aabb[l][i][1] > boxes[j][1]) &&
 
  672                            (aabb[l][i][2] < boxes[j][2]) &&
 
  673                            (aabb[l][i][3] < boxes[j][3]))
 
  678                                   (boxes[closest_box][0] > boxes[j][0]) &&
 
  679                                   (boxes[closest_box][1] > boxes[j][1]) &&
 
  680                                   (boxes[closest_box][2] < boxes[j][2]) &&
 
  681                                   (boxes[closest_box][3] < boxes[j][3])))
 
  692                    else if(    (aabb[l][i][0] == boxes[j][0]) &&
 
  693                                (aabb[l][i][1] == boxes[j][1]) &&
 
  694                                (aabb[l][i][2] == boxes[j][2]) &&
 
  695                                (aabb[l][i][3] == boxes[j][3]))
 
  707                connectedComponents[closest_box].push_back(res[l][i]);
 
  716template<
short_t d, 
class T>
 
  724    std::vector<std::vector<std::vector< std::vector<T> > > >& trim_curves)
 const 
  736    std::vector< std::vector< std::vector< std::vector< T > > > > res; 
 
  737    std::vector< std::vector< std::vector<index_t> > > aabb;
 
  738    std::vector< std::vector<index_t> > boxes;
 
  739    aabb = this->domainBoundariesParams(res);
 
  740    for(
unsigned int i = 0; i < aabb.size(); i++)
 
  743        for(
unsigned int j = 0; j < aabb[i].size(); j++)
 
  745            bool is_boundary = 
true;
 
  746            for(
unsigned int k = 0; k < aabb[i].size(); k++)
 
  750                    if(     (aabb[i][j][0] > aabb[i][k][0]) && 
 
  751                            (aabb[i][j][2] < aabb[i][k][2]) &&
 
  752                            (aabb[i][j][1] > aabb[i][k][1]) &&
 
  753                            (aabb[i][j][3] < aabb[i][k][3]) )
 
  755                        is_boundary= !is_boundary;
 
  761                boxes.push_back(aabb[i][j]);
 
  762                boxes[boxes.size()-1].push_back(i);
 
  763                trim_curves.resize(trim_curves.size()+1);
 
  764                trim_curves[trim_curves.size()-1].push_back(res[i][j]);
 
  769    b1.resize(boxes.size(),d);
 
  770    b2.resize(boxes.size(),d);
 
  771    level.resize(boxes.size());
 
  772    for(
size_t i = 0; i < boxes.size(); i++){
 
  773        for(
unsigned j = 0; j < d; j++){
 
  776            b1(i,j) = boxes[i][j];
 
  777            b2(i,j) = boxes[i][j+d];
 
  779        level[i] = boxes[i][2*d];
 
  782    int nboxes = level.size();
 
  787    nvertices.resize(nboxes,this->dim());
 
  789    for (
int i = 0; i < nboxes; i++)
 
  795        p1 = b1.row(i).transpose();
 
  796        p2 = b2.row(i).transpose();
 
  798        this->getBsplinePatchGlobal(p1, p2, level[i], geom_coef, new_cp, cku, ckv);
 
  807            int cprows = cp.rows();
 
  811            bigger.resize(cp.rows()+new_cp.rows(), cp.cols());
 
  812            for(
int j=0; j< bigger.rows(); j++)
 
  816                    bigger.row(j) = cp.row(j);
 
  820                    bigger.row(j) = new_cp.row(j-cprows);
 
  831    for(
size_t l = 0; l < aabb.size();l++) 
 
  833        for(
size_t i = 0; i < aabb[l].size(); i++) 
 
  835            int closest_box = -1;
 
  836            for(
size_t j = 0; j < boxes.size(); j++)
 
  838                if(l == 
static_cast<size_t>(boxes[j][4]) )
 
  840                    if(     (aabb[l][i][0] > boxes[j][0]) &&
 
  841                            (aabb[l][i][1] > boxes[j][1]) &&
 
  842                            (aabb[l][i][2] < boxes[j][2]) &&
 
  843                            (aabb[l][i][3] < boxes[j][3]))
 
  849                                   (boxes[closest_box][0] > boxes[j][0]) &&
 
  850                                   (boxes[closest_box][1] > boxes[j][1]) &&
 
  851                                   (boxes[closest_box][2] < boxes[j][2]) &&
 
  852                                   (boxes[closest_box][3] < boxes[j][3])))
 
  863                    else if(    (aabb[l][i][0] == boxes[j][0]) &&
 
  864                                (aabb[l][i][1] == boxes[j][1]) &&
 
  865                                (aabb[l][i][2] == boxes[j][2]) &&
 
  866                                (aabb[l][i][3] == boxes[j][3]))
 
  878                trim_curves[closest_box].push_back(res[l][i]);
 
 
  886template<
short_t d, 
class T>
 
  889    std::vector<std::vector<std::vector< std::vector<T> > > >& trim_curves)
 const 
  897    for(
size_t i = 0; i < this->m_xmatrix.size(); i++)
 
  899        if(this->m_xmatrix[i].size()>0)
 
  905    gsDebug<<
"new min level"<<
"\n";
 
  906    std::vector< std::vector< std::vector< std::vector< T > > > > res; 
 
  907    std::vector< std::vector< std::vector<index_t > > > aabb;
 
  908    std::vector< std::vector<index_t > > boxes;
 
  909    aabb = this->domainBoundariesParams(res);
 
  910    for(
size_t i = 0; i < aabb.size(); i++)
 
  913        for(
size_t j = 0; j < aabb[i].size(); j++)
 
  915            bool is_boundary = 
true;
 
  916            for(
size_t k = 0; k < aabb[i].size(); k++)
 
  920                    if(     (aabb[i][j][0] > aabb[i][k][0]) &&
 
  921                            (aabb[i][j][2] < aabb[i][k][2]) &&
 
  922                            (aabb[i][j][1] > aabb[i][k][1]) &&
 
  923                            (aabb[i][j][3] < aabb[i][k][3]) )
 
  925                        is_boundary= !is_boundary;
 
  931                boxes.push_back(aabb[i][j]);
 
  932                boxes[boxes.size()-1].push_back(i+first_level);
 
  933                trim_curves.resize(trim_curves.size()+1);
 
  934                trim_curves[trim_curves.size()-1].push_back(res[i][j]);
 
  939    b1.resize(boxes.size(),d);
 
  940    b2.resize(boxes.size(),d);
 
  941    level.resize(boxes.size());
 
  942    for(
size_t i = 0; i < boxes.size(); i++)
 
  944        for(
unsigned j = 0; j < d; j++)
 
  948            b1(i,j) = boxes[i][j];
 
  949            b2(i,j) = boxes[i][j+d];
 
  951        level[i] = boxes[i][2*d];
 
  955    int nboxes = level.size();
 
  960    p1.resize(this->dim());
 
  961    p2.resize(this->dim());
 
  966    for (
int i = 0; i < nboxes; i++)
 
  968        p1 = b1.row(i).transpose();
 
  969        p2 = b2.row(i).transpose();
 
  971        this->getBsplinePatchGlobal(p1, p2, level[i], geom_coef, temp1, cku, ckv);
 
  977    for(
size_t l = 0; l < aabb.size();l++) 
 
  979        for(
size_t i = 0; i < aabb[l].size(); i++) 
 
  981            int closest_box = -1;
 
  982            for(
size_t j = 0; j < boxes.size(); j++)
 
  984                if(l == 
static_cast<size_t>(boxes[j][4]) )
 
  986                    if(     (aabb[l][i][0] > boxes[j][0]) &&
 
  987                            (aabb[l][i][1] > boxes[j][1]) &&
 
  988                            (aabb[l][i][2] < boxes[j][2]) &&
 
  989                            (aabb[l][i][3] < boxes[j][3]))
 
  994                                   (boxes[closest_box][0] > boxes[j][0]) &&
 
  995                                   (boxes[closest_box][1] > boxes[j][1]) &&
 
  996                                   (boxes[closest_box][2] < boxes[j][2]) &&
 
  997                                   (boxes[closest_box][3] < boxes[j][3])))
 
 1008                    else if(    (aabb[l][i][0] == boxes[j][0]) &&
 
 1009                                (aabb[l][i][1] == boxes[j][1]) &&
 
 1010                                (aabb[l][i][2] == boxes[j][2]) &&
 
 1011                                (aabb[l][i][3] == boxes[j][3]))
 
 1023                trim_curves[closest_box].push_back(res[l][i]);
 
 
 1036template<
short_t d, 
class T>
 
 1040    const index_t n = thbCoefs.cols();
 
 1043    lvlCoefs.setZero(m_bases[0]->size(), n);
 
 1044    for(cmatIterator it = m_xmatrix[0].begin(); it != m_xmatrix[0].end(); ++it)
 
 1046        const int hIndex = m_xmatrix_offset[0] + (it - m_xmatrix[0].begin());
 
 1047        lvlCoefs.row(*it) = thbCoefs.row(hIndex);
 
 1051    std::vector<T> knots_x, knots_y;
 
 1053    for(
int l = 1; l <=level; l++)
 
 1055        k1 = m_bases[l-1]->knots(0);
 
 1056        k2 = m_bases[l-1]->knots(1);
 
 1063        lvlCoefs.resize(m_bases[l-1]->size(0), n * m_bases[l-1]->size(1));
 
 1064        gsBoehmRefine(k1, lvlCoefs, m_deg[0], knots_x.begin(), knots_x.end(), 
false);
 
 1068        gsBoehmRefine(k2, lvlCoefs, m_deg[1], knots_y.begin(), knots_y.end(), 
false);
 
 1070        lvlCoefs.resize(m_bases[l]->size(), n); 
 
 1073        for(cmatIterator it = m_xmatrix[l].begin(); it != m_xmatrix[l].end(); ++it)
 
 1075            const int hIndex = m_xmatrix_offset[l] + (it - m_xmatrix[l].begin());
 
 1076            lvlCoefs.row(*it) = thbCoefs.row(hIndex);
 
 
 1086    for (
short_t i = d-1; i>=0; --i)
 
 1088        rvo[d-i-1] = m % str[i];
 
 1094template<
class Iter, 
class Vec> 
 
 1096bool findNextMatch(Iter & it, Iter stop, Vec & cur,
 
 1097                     const Vec & low, 
const Vec & upp,
 
 1106        else if ( ci < *it )
 
 1133            it = std::lower_bound(it, stop, ci);
 
 1141template<
class T, 
class Vec>
 
 1142bool findNextMatch(
const gsSparseVector<T> & sv, 
index_t & ii,
 
 1143                     Vec & cur, 
const Vec & low, 
const Vec & upp,
 
 1147    while ( ii<sv.data().size() )
 
 1150        if ( ci == sv.data().index(ii) )
 
 1152        else if ( ci < sv.data().index(ii) )
 
 1157            ii = sv.data().searchLowerIndex(ii, sv.data().size(), ci);
 
 1248template<
short_t d, 
class T>
 
 1249index_t gsTHBSplineBasis<d,T>::numActiveMax(
const gsMatrix<T> & u,
 
 1250                                            gsMatrix<index_t> & offset)
 const 
 1252    point low, upp, cur, mstr, str, ll, uu, cc;
 
 1255    const int maxLevel = this->m_tree.getMaxInsLevel();
 
 1256    gsMatrix<T> currPoint;
 
 1257    offset.setZero(maxLevel+1, u.cols()); 
 
 1259    std::vector<std::vector<std::pair<index_t,index_t> > > tfunction(maxLevel+1);
 
 1261    for (
index_t i = 0; i != u.cols(); i++) 
 
 1263        currPoint = u.col(i);
 
 1264        for(
short_t k = 0; k != d; ++k)
 
 1265            low[k] = m_bases[maxLevel]->knots(k).uFind( currPoint.at(k) ).uIndex();
 
 1267            this->_knotIndexToDiadicIndex(maxLevel,low);
 
 1269        index_t lvl = std::min(this->m_tree.levelOf(low, maxLevel),(
int) m_xmatrix.size()-1);
 
 1272        typename CMatrix::const_iterator it;
 
 1273        for(
int level = 0; level <= maxLevel; level++) 
 
 1275            if (level>lvl && tfunction[level].empty()) 
 
 1278            m_bases[level]->active_cwise(currPoint, low, upp);
 
 1279            m_bases[level]->stride_cwise(mstr);
 
 1283                it = m_xmatrix[level].begin();
 
 1285                while ( findNextMatch(it, m_xmatrix[level].end(), cur, low, upp, mstr) )
 
 1287                    const index_t act = this->m_xmatrix_offset[level] + (it - m_xmatrix[level].begin());
 
 1288                    if ( isTruncated(act) )
 
 1291                        const int rlvl = m_is_truncated.at(act);
 
 1292                        tfunction[rlvl].push_back( std::make_pair(act,level) );
 
 1302            if ( !tfunction[level].empty() )
 
 1304                for ( std::pair<index_t,index_t> & q : tfunction[level] )
 
 1306                    const gsSparseVector<T>& coefs = getCoefs(q.first);
 
 1309                    if ( findNextMatch(coefs, ii, cur, low, upp, mstr) )
 
 1310                        ++offset(q.second, i);
 
 1312                tfunction[level].clear();
 
 1321    return offset.sum();
 
 1326template<
short_t d, 
class T>
 
 1327void gsTHBSplineBasis<d,T>::activeAtLevel_into(
index_t lvl, 
const gsMatrix<T>& u,
 
 1328                                               std::vector<index_t> & result)
 const 
 1330    gsMatrix<index_t> ind;
 
 1332    point low, upp, cur, mstr, str, ll, uu, cc;
 
 1333    m_bases[lvl]->active_cwise(u, low, upp);
 
 1335    this->m_bases[lvl]->stride_cwise(mstr);
 
 1337    typename CMatrix::const_iterator it = m_xmatrix[lvl].begin();
 
 1338    typename CMatrix::const_iterator end = m_xmatrix[lvl].end();
 
 1339    while ( findNextMatch(it, end, cur, low, upp, mstr) )
 
 1341        const index_t act = this->m_xmatrix_offset[lvl] + (it - m_xmatrix[lvl].begin());
 
 1342        if ( isTruncated(act) )
 
 1344            const gsSparseVector<T>& coefs = getCoefs(act);
 
 1345            const gsTensorBSplineBasis<d, T>& base =
 
 1346                *this->m_bases[this->m_is_truncated[act]];
 
 1348            base.active_cwise(u, ll, uu);
 
 1349            base.stride_cwise(str);
 
 1352            if ( findNextMatch(coefs, ii, cc, ll, uu, str) )
 
 1353                result.push_back(act);
 
 1410            result.push_back(act);
 
 1416template<
short_t d, 
class T>
 
 1421    const int maxLevel = this->m_tree.getMaxInsLevel();
 
 1426    std::vector<std::vector<index_t> > temp_output;
 
 1427    temp_output.resize( u.cols() );
 
 1428    size_t ov = 1, sz = 0;
 
 1429    for(
short_t i = 0; i != d; ++i)
 
 1430        ov *= this->m_bases.front()->component(i).degree() + 1;
 
 1432    for(
index_t p = 0; p < u.cols(); p++) 
 
 1434        temp_output[p].reserve(ov+2);
 
 1435        currPoint = u.col(p);
 
 1436        for(
short_t i = 0; i != d; ++i)
 
 1437            low[i] = m_bases[maxLevel]->knots(i).uFind( currPoint(i,0) ).uIndex();
 
 1440            this->_knotIndexToDiadicIndex(maxLevel,low);
 
 1443        const int lvl = std::min(this->m_tree.levelOf(low, maxLevel),(
int) m_xmatrix.size()-1);
 
 1445        for(
int i = 0; i <= lvl; i++)
 
 1446            activeAtLevel_into(i,currPoint,temp_output[p]);
 
 1449        if ( temp_output[p].size() > sz )
 
 1450            sz = temp_output[p].size();
 
 1453    result.resize(sz, u.cols() );
 
 1454    for(
index_t i = 0; i < u.cols(); i++)
 
 1456        std::copy(temp_output[i].begin(), temp_output[i].end(), result.col(i).data() );
 
 1457        result.col(i).bottomRows(sz-temp_output[i].size()).setZero();
 
 
 1461template<
short_t d, 
class T>
 
 1466    if (this->m_is_truncated[i] == -1)  
 
 1468        unsigned level = this->levelOf(i);
 
 1469        unsigned tensor_index = flatTensorIndexOf(i, level);
 
 1470        this->m_bases[level]->evalSingle_into(tensor_index, u, result);
 
 1475        unsigned level = this->m_is_truncated[i];
 
 1478            *this->m_bases[level];
 
 1481            (u, base, coefs, result);
 
 
 1485template<
short_t d, 
class T>
 
 1491    if (this->m_is_truncated[i] == -1) 
 
 1493        const unsigned level = this->levelOf(i);
 
 1494        const unsigned fl_tensor_index = flatTensorIndexOf(i, level);
 
 1495        this->m_bases[level]->deriv2Single_into(fl_tensor_index, u, result);
 
 1499        const unsigned level = this->m_is_truncated[i];
 
 1502            *this->m_bases[level];
 
 1504        gsTensorDeriv2_into<d, T, gsKnotVector<T>,
 
 
 1509template<
short_t d, 
class T>
 
 1521    this->active_into(u, indices);
 
 1523    result.setZero(indices.rows(), u.cols());
 
 1525    for (
index_t i = 0; i < indices.cols(); i++)
 
 1527        for (
index_t j = 0; j < indices.rows(); j++)
 
 1529            const index_t index = indices(j, i);
 
 1530            if (j != 0 && index == 0)
 
 1533            evalSingle_into(index, u.col(i), res);
 
 1534            result(j, i) = res.value();
 
 
 1540template<
short_t d, 
class T>
 
 1544    this->active_into(u, indices);
 
 1546    static const short_t numDers = (d * (d + 1)) / 2;
 
 1549    result.setZero(indices.rows() * numDers, u.cols());
 
 1551    for (
int i = 0; i < indices.cols(); i++)
 
 1553        for (
int j = 0; j < indices.rows(); j++)
 
 1555            const index_t index = indices(j, i);
 
 1556            if (j != 0 && index == 0)
 
 1559            this->deriv2Single_into(index, u.col(i), res);
 
 1561            result.template block<numDers, 1>(j * numDers, i) = res;
 
 
 1567template<
short_t d, 
class T>
 
 1571    this->active_into(u, indices);
 
 1574    result.setZero(indices.rows() * d, u.cols());
 
 1576    for (
int i = 0; i < indices.cols(); i++)
 
 1578        for (
int j = 0; j < indices.rows(); j++)
 
 1581            const unsigned index = indices(j, i);
 
 1582            if (j != 0 && index == 0)
 
 1585            this->derivSingle_into(index, u.col(i), res);
 
 1586            result.template block<d,1>(j * d, i) = res;
 
 
 1592template<
short_t d, 
class T>
 
 1598    if (this->m_is_truncated[i] == -1) 
 
 1600        unsigned level = this->levelOf(i);
 
 1601        unsigned fl_tensor_index = flatTensorIndexOf(i, level);
 
 1602        this->m_bases[level]->derivSingle_into(fl_tensor_index, u, result);
 
 1606        unsigned level = this->m_is_truncated[i];
 
 1609            *this->m_bases[level];
 
 1610        gsTensorDeriv_into<d, T, gsKnotVector<T>,
 
 
 1616template<
short_t d, 
class T>
 
 1617inline void eval_tp(
const std::vector<std::vector<
gsMatrix<T> > > & cw,
 
 1625        result[0](m,i) = cw[0][0].at(ti.
at(0));
 
 1627            result[0](m,i) *= cw[k][0].at(ti.
at(k));
 
 1632        T * acc = result[1].col(i).data() + str[1]*m;
 
 1636            *acc = cw[k][1].
at(ti.
at(k));   
 
 1638                *acc *= cw[l][0].at(ti.
at(l));
 
 1639            for ( 
short_t l=k+1; l<d; ++l)
 
 1640                *acc *= cw[l][0].at(ti.
at(l));
 
 1648        T * acc = result[2].col(i).data() + str[2]*m;
 
 1651            *(acc+k) = cw[k][2].at(ti.
at(k)); 
 
 1653                *(acc+k) *= cw[l][0].
at(ti.
at(l));
 
 1654            for ( 
short_t l=k+1; l<d; ++l)
 
 1656                *(acc+k) *= cw[l][0].at(ti.
at(l));
 
 1658                *(acc+mc) = cw[k][1].at(ti.
at(k)) * cw[l][1].at(ti.
at(l));
 
 1660                    *(acc+mc) *= cw[q][0].at(ti.
at(q));
 
 1661                for ( 
short_t q=k+1; q<l; ++q)
 
 1662                    *(acc+mc) *= cw[q][0].at(ti.
at(q));
 
 1663                for ( 
short_t q=l+1; q<d; ++q)
 
 1664                    *(acc+mc) *= cw[q][0].at(ti.
at(q));
 
 1672template<
short_t d, 
class T>
 
 1675                                             bool sameElement)
 const 
 1679    if (0==u.cols()) 
return;
 
 1681    const int maxLevel = this->m_tree.getMaxInsLevel();
 
 1684    std::vector<std::vector<std::pair<index_t,index_t> > > tfunction(maxLevel+1);
 
 1687    std::vector<std::vector<gsMatrix<T> > > cw(d);
 
 1688    std::vector<gsMatrix<T> > temp(n+1), cwa(d);
 
 1689    point low(d), ki(d), kil, til, tlow, tupp, tstr, tcur;
 
 1693        active_into(u.col(0), act);
 
 1695        active_into(u, act);
 
 1699    index_t result_size = act.rows();
 
 1705    for (
int l = 0; l <= n; l++)
 
 1708        result[l].setZero(result_size * str[l], u.cols());
 
 1709        temp[l].resize(str[l], 1);
 
 1716    std::vector<std::vector<index_t> >thbact;
 
 1719    for (
index_t i = 0; i != u.cols(); i++) 
 
 1721        if (!sameElement || 0==i)
 
 1724            for(
short_t k = 0; k != d; ++k)
 
 1726                auto kit = m_bases[maxLevel]->knots(k).uFind( u(k,i) );
 
 1727                ki [k] = kit.lastAppearance() - m_bases[maxLevel]->degree(k);
 
 1728                low[k] = kit.uIndex();
 
 1731                this->_knotIndexToDiadicIndex(maxLevel,low);
 
 1732            lvl = std::min(this->m_tree.levelOf(low, maxLevel),(
int) m_xmatrix.size()-1);
 
 1736        thbact.resize(lvl+1);
 
 1738        for(
int level = 0; level <= maxLevel; level++) 
 
 1740            if ( level<=lvl && (!sameElement || 0==i) )
 
 1742                thbact[level].clear();
 
 1743                activeAtLevel_into(level,u.col(i),thbact[level]);
 
 1746            if ( (level>lvl || thbact[level].empty()) && tfunction[level].empty()) 
 
 1751                *tmp.data() = u(k,i);
 
 1752                this->m_bases[level]->component(k).evalAllDers_into(tmp, n, cw[k], sameElement);
 
 1753                kil[k] = this->m_bases[level]->component(k).firstActive( u(k,i) ); 
 
 1757                for (
size_t j = 0; j!=thbact[level].size(); ++j, ++m) 
 
 1759                    index = thbact[level][j];
 
 1760                    if ( !isTruncated(index) ) 
 
 1762                        til = this->m_bases[level]->tensorIndex(flatTensorIndexOf(index, level));
 
 1764                        eval_tp(cw,til,n,str,m,i,result);
 
 1766                    else if (!sameElement || 0==i)
 
 1772                        const int rlvl = m_is_truncated.at(index);
 
 1773                        tfunction[rlvl].push_back( std::make_pair(index,m) );
 
 1778            if ( !tfunction[level].empty() )
 
 1780                for(
short_t k = 0; k!=d; ++k) til[k] = kil[k] + this->m_bases[level]->degree(k);
 
 1781                this->m_bases[level]->stride_cwise(tstr);
 
 1784                for ( std::pair<index_t,index_t> & q : tfunction[level] )
 
 1791                    while ( findNextMatch(coefs, ii, tcur, kil, til, tstr) )
 
 1794                        eval_tp(cw,low,n,str,0,0,temp);
 
 1795                        for ( 
index_t l = 0; l <= n; ++l)
 
 1797                            auto acc = result[l].col(i).segment(q.second*str[l],str[l]);
 
 1798                            acc += coefs.data().value(ii) * temp[l];
 
 1825                    tfunction[level].clear();
 
 
 1950template<
short_t d, 
class T>
 
 1952    typename gsTHBSplineBasis<d, T>::AxisAlignedBoundingBox& boundaryAABB,
 
 1953    typename gsTHBSplineBasis<d, T>::TrimmingCurves& trimCurves)
 const 
 1955    Polylines polylines;
 
 1956    AxisAlignedBoundingBox aabb;
 
 1958    aabb = this->domainBoundariesParams(polylines);
 
 1959    breakCycles(aabb, polylines);
 
 1961    int numBoundaryBoxes = 0;
 
 1962    for (
unsigned level = 0; level != aabb.size(); level++)
 
 1964        boundaryAABB.push_back(std::vector< std::vector<index_t> >());
 
 1965        trimCurves.push_back(std::vector< std::vector< std::vector< std::vector<T> > > >());
 
 1968        for (
unsigned boxI = 0; boxI != aabb[level].size(); boxI++)
 
 1970            bool isBoundaryBox = 
true;
 
 1971            for (
unsigned boxJ = 0; boxJ != aabb[level].size(); boxJ++)
 
 1975                    if (isFirstBoxCompletelyInsideSecond(aabb[level][boxI], aabb[level][boxJ]))
 
 1977                        isBoundaryBox = !isBoundaryBox;
 
 1985                boundaryAABB[level].push_back(aabb[level][boxI]);
 
 1988                trimCurves[level].push_back(std::vector< std::vector< std::vector<T> > >());
 
 1989                trimCurves[level][trimCurves[level].
size() - 1].push_back(polylines[level][boxI]);
 
 1994    for (
unsigned level = 0; level != aabb.size(); level++)
 
 1996        for (
unsigned box = 0; box != aabb[level].size(); box++)
 
 1998            int closestBox = -1;
 
 1999            for (
unsigned boundBox = 0;
 
 2000                 boundBox != boundaryAABB[level].
size();
 
 2003                if (isFirstBoxCompletelyInsideSecond(aabb[level][box], boundaryAABB[level][boundBox]))
 
 2005                    if (closestBox == -1 ||
 
 2006                        !isFirstBoxCompletelyInsideSecond(boundaryAABB[level][closestBox],
 
 2007                                                          boundaryAABB[level][boundBox]))
 
 2009                        closestBox = boundBox;
 
 2012                else if (areBoxesTheSame(aabb[level][box], boundaryAABB[level][boundBox]))
 
 2019            if (-1 < closestBox)
 
 2021                trimCurves[level][closestBox].push_back(polylines[level][box]);
 
 
 2028template<
short_t d, 
class T>
 
 2030typename util::enable_if<dd==2,gsTensorBSpline<d,T> >::type
 
 2032                                            const unsigned level,
 
 2036    for (
unsigned dim = 0; dim != d; dim++)
 
 2038        low(dim) = boundingBox[dim];
 
 2039        upp(dim) = boundingBox[d + dim];
 
 2041    this->m_tree.computeLevelIndex(low, level, low);
 
 2042    this->m_tree.computeLevelIndex(upp, level, upp);
 
 2044    const gsKnotVector<T> & knots0 = m_bases[level]->knots(0);
 
 2045    const gsKnotVector<T> & knots1 = m_bases[level]->knots(1);
 
 2049        this->_diadicIndexToKnotIndex(level,low);
 
 2050        this->_diadicIndexToKnotIndex(level,low);
 
 2053    const int lowIndex0 = knots0.lastKnotIndex (low(0)) - m_deg[0];
 
 2054    const int uppIndex0 = knots0.firstKnotIndex(upp(0)) - 1;
 
 2055    const int lowIndex1 = knots1.lastKnotIndex (low(1)) - m_deg[1];
 
 2056    const int uppIndex1 = knots1.firstKnotIndex(upp(1)) - 1;
 
 2058    const int numDirection0 = uppIndex0 - lowIndex0 + 1;
 
 2059    const int numDirection1 = uppIndex1 - lowIndex1 + 1;
 
 2060    const int numNewCoefs = numDirection0 * numDirection1;
 
 2061    gsMatrix<T> newCoefs(numNewCoefs, geomCoefs.cols());
 
 2062    const index_t sz0 = m_bases[level]->size(0);
 
 2065    globalRefinement(geomCoefs, level, coefs);
 
 2067    for (
int j = lowIndex1; j <= uppIndex1; j++)
 
 2068        for (
int i = lowIndex0; i <= uppIndex0; i++)
 
 2069            newCoefs.row(cc++) =  coefs.row(j*sz0+i);
 
 2071    std::vector<gsKnotVector<T> > kv(2);
 
 2073    kv[0] = gsKnotVector<T>(m_deg[0], knots0.begin() + lowIndex0,
 
 2074                            knots0.begin() + uppIndex0 + m_deg[0] + 2);
 
 2076    kv[1] = gsKnotVector<T>(m_deg[1], knots1.begin() + lowIndex1,
 
 2077                            knots1.begin() + uppIndex1 + m_deg[1] + 2);
 
 2079    tensorBasis basis(kv);
 
 2081    return gsTensorBSpline<d, T> (basis, newCoefs);
 
 2084template<
short_t d, 
class T>
 
 2091{ getBsplinePatchGlobal_impl<d>(b1,b2,level,geom_coef,cp,k1,k2); };
 
 
 2093template<
short_t d, 
class T>
 
 2095                                         const unsigned level,
 
 2097{ 
return getBSplinePatch_impl<d>(boundingBox, level, geomCoefs); }
 
 
 2100template<
short_t d, 
class T>
 
 2102    typename gsTHBSplineBasis<d, T>::AxisAlignedBoundingBox& aabb,
 
 2103    typename gsTHBSplineBasis<d, T>::Polylines& polylines)
 const 
 2105    for (
size_t level = 0; level != polylines.
size(); level++)
 
 2107        for (
size_t line = 0; line != polylines[level].
size(); line++)
 
 2110            index_t segment = identifyCycle(polylines[level][line], pt);
 
 2114                std::vector< std::vector<T> > part1, part2;
 
 2115                breakPolylineIntoTwoParts(polylines[level][line], segment, pt,
 
 2118                polylines[level][line] = part1;
 
 2119                polylines[level].push_back(part2);
 
 2121                std::vector<index_t> aabb1, aabb2;
 
 2122                findNewAABB(part1, aabb1);
 
 2123                findNewAABB(part2, aabb2);
 
 2125                aabb[level][line] = aabb1;
 
 2126                aabb[level].push_back(aabb2);
 
 
 2140template<
short_t d, 
class T>
 
 2142                                              std::pair<T, T>& pt)
 const 
 2144    std::map< std::pair<T, T>, 
index_t > times;
 
 2145    std::map< std::pair<T, T>, 
index_t > index;
 
 2147    for (
size_t seg = 0; seg != line.size(); seg++)
 
 2149        const size_t seg1 = (seg + 1) % line.
size();
 
 2151        std::pair<T, T> currentPt( line[seg][0], line[seg][1] );
 
 2152        if (!((currentPt.first == line[seg1][0] && currentPt.second == line[seg1][1]) ||
 
 2153              (currentPt.first == line[seg1][2] && currentPt.second == line[seg1][3])))
 
 2155            currentPt.first = line[seg][2];
 
 2156            currentPt.second = line[seg][3];
 
 2159        size_t count = times.count(currentPt);
 
 2162            times[currentPt] += 1;
 
 2166            times[currentPt] = 1;
 
 2167            index[currentPt] = seg1;
 
 2171    typedef typename std::map< std::pair<T, T>, 
index_t>::iterator iterator;
 
 2172    for (iterator it = times.begin(); it != times.end(); it++)
 
 2174        if (it->second == 2)
 
 2177            return index[it->first];
 
 2181            GISMO_ENSURE(2 > it->second, 
"Internal error. Check the polylines from the domainBoundariesParam." );
 
 
 2188template<
short_t d, 
class T>
 
 2190    const std::vector< std::vector< T> >& line,
 
 2192    const std::pair<T, T>& meetingPt,
 
 2193    std::vector< std::vector<T> >& part1,
 
 2194    std::vector< std::vector<T> >& part2)
 const 
 2200    for (
index_t i = 0; i != length; i++)
 
 2202        const index_t seg = (i + segment) % length;
 
 2207            part1.push_back(line[seg]);
 
 2212            if ((meetingPt.first == line[seg][0] && meetingPt.second == line[seg][1]) ||
 
 2213                (meetingPt.first == line[seg][2] && meetingPt.second == line[seg][3]))
 
 2217                    part1.push_back(line[seg]);
 
 2223                    part2.push_back(line[seg]);
 
 2230                    part1.push_back(line[seg]);
 
 2234                    part2.push_back(line[seg]);
 
 
 2242template<
short_t d, 
class T>
 
 2244                                         std::vector<index_t>& aabb)
 const 
 2246    T minX = polyline[0][0];
 
 2247    T minY = polyline[0][1];
 
 2248    T maxX = polyline[0][2];
 
 2249    T maxY = polyline[0][3];
 
 2252    for (
size_t seg = 0; seg != polyline.size(); seg++)
 
 2254        if (polyline[seg][0] < minX)
 
 2256            minX = polyline[seg][0];
 
 2258        if (polyline[seg][1] < minY)
 
 2260            minY = polyline[seg][1];
 
 2262        if (maxX < polyline[seg][2])
 
 2264            maxX = polyline[seg][2];
 
 2266        if (maxY < polyline[seg][3])
 
 2268            maxY = polyline[seg][3];
 
 2272    unsigned maxLevel = this->maxLevel();
 
 2277    for (
unsigned i = 0; i != kv0.
uSize(); i++)
 
 2290    for (
unsigned i = 0; i != kv1.
uSize(); i++)
 
 
 2310template<
short_t d, 
class T>
 
 2316    this->m_tree.getBoxesInLevelIndex(b1,b2,level);
 
 2318    std::vector< gsSparseMatrix<T,RowMajor> > transfer;
 
 2319    transfer.resize(this->maxLevel() );
 
 2320    std::vector<std::vector<T> > knots(d);
 
 2322    for(
unsigned i = 0; i < this->maxLevel(); ++i)
 
 2325        for(
short_t dim = 0; dim < d; dim++)
 
 2341    std::vector< gsSparseMatrix<T,RowMajor> > temp_transf;
 
 2342    for(
unsigned j = 0; j < this->maxLevel(); ++j)
 
 2344        std::vector<CMatrix> x_mat_old_0, x_matrix_lvl;
 
 2345        this->setActiveToLvl(j,x_mat_old_0);
 
 2346        this->setActiveToLvl(j+1,x_matrix_lvl);
 
 2347        temp_transf.push_back(transfer[j]);
 
 2349        gsSparseMatrix<T> crs = this->coarsening_direct(x_mat_old_0, x_matrix_lvl, temp_transf);
 
 2350        result.push_back(crs);
 
 
 2355template<
short_t d, 
class T>
 
 2358    int size1= 0, size2 = 0;
 
 2360    for(
unsigned int i =0; i< old.size();i++)
 
 2362        size1 += old[i].size();
 
 2364    for(
unsigned int i =0; i< n.size();i++)
 
 2366        size2 += n[i].size();
 
 2372    for (
unsigned int i = 0; i < old.size(); i++)
 
 2376        for(
unsigned int l =0; l < i; l++)
 
 2378            start_lv_i += n[l].size();
 
 2381        for (
unsigned int j = 0; j < old[i].size();j++)
 
 2384            for(
unsigned int l =0; l < i; l++)
 
 2386                start_lv_i += n[l].size();
 
 2388            const unsigned old_ij = old[i][j];  
 
 2390            if( n[i].bContains(old_ij) )
 
 2392                result(start_lv_i + std::distance(n[i].begin(), n[i].find_it_or_fail(old_ij) ),glob_numb ) = 1;
 
 2393                for(
int k = 0; k < transferDense.rows(); k++)
 
 2396                    if(transferDense(k, old_ij) != 0)
 
 2398                        if(n[i+1].bContains(k))
 
 2400                            const int pos = start_lv_i + n[i].size() + std::distance(n[i+1].begin(), n[i+1].find_it_or_fail(k));
 
 2401                            result(pos,glob_numb) = transferDense(k, old_ij);
 
 2408                for(
int k = 0; k < transferDense.rows(); k++)
 
 2410                    if(transferDense(k, old_ij) != 0)
 
 2412                        const int pos = start_lv_i + n[i].size() + std::distance(n[i+1].begin(), n[i+1].find_it_or_fail(k));
 
 2413                        result(pos,glob_numb) = transferDense(k, old_ij);
 
 
 2422template<
short_t d, 
class T>
 
 2427    int size1= 0;
int size2 = 0;
 
 2429    for(
unsigned int i =0; i< old.size();i++){
 
 2430        size1 += old[i].
size();
 
 2432    for(
unsigned int i =0; i< n.size();i++){
 
 2433        size2 += n[i].size();
 
 2435    gsSparseMatrix<T> result(size2,size1);
 
 2443    std::vector<gsSparseMatrix<T,ColMajor> > temptransfer;
 
 2444    temptransfer.resize(transfer.size());
 
 2445    for (
unsigned int i = 0; i < transfer.size();i++){
 
 2446        temptransfer[i] = transfer[i];
 
 2450    for (
unsigned int i = 0; i < old.size(); i++)
 
 2454        for(
unsigned int l =0; l < i; l++)
 
 2456            start_lv_i += n[l].size();
 
 2460        for (
unsigned int j = 0; j < old[i].size();j++)
 
 2464            for(
unsigned int l =0; l < i; l++)
 
 2466                start_lv_i += n[l].size();
 
 2468            const unsigned old_ij = old[i][j];  
 
 2469            gsMatrix<index_t, d, 2> supp(d, 2);
 
 2470            this->m_bases[i]->elementSupport_into(old_ij, supp);
 
 2474            gsSparseVector<T,RowMajor> t(this->m_bases[i]->size());
 
 2477            for(
unsigned int k = i; k < n.size();k++){
 
 2482                    for(
int l = 0 ; l < t.size();l++){
 
 2484                            if( (k) < (old.size()))
 
 2485                                if(old[k].bContains(l)){
 
 2495                    for(
unsigned int l =0; l < k-1; l++)
 
 2497                        start_lv_i += n[l].size();
 
 2501                for(
int l = 0 ; l < t.size();l++)
 
 2505                        if(n[k].bContains(l))
 
 2513                            const int pos = start_lv_i + p + std::distance(n[k].begin(), n[k].find_it_or_fail(l));
 
 2517                            result(pos,glob_numb) = t[l];
 
 2525                gsSparseVector<T,RowMajor> temp(this->m_bases[i]->size());
 
 2528                if(k<temptransfer.size())
 
 2529                    temp = temptransfer[k] * t.transpose();
 
 2599template<
short_t d, 
class T>
 
 2600gsSparseMatrix<T> gsTHBSplineBasis<d,T>::coarsening_direct( 
const std::vector<gsSortedVector<index_t> >& old,
 
 2601                                                      const std::vector<gsSortedVector<index_t> >& n,
 
 2602                                                      const std::vector<gsSparseMatrix<T,RowMajor> >& transfer)
 const 
 2604    GISMO_ASSERT(old.size() < n.size(), 
"old,n problem in coarsening.");
 
 2606    int size1= 0;
int size2 = 0;
 
 2608    for(
unsigned int i =0; i< old.size();i++)
 
 2610        size1 += old[i].size();
 
 2612    for(
unsigned int i =0; i< n.size();i++)
 
 2614        size2 += n[i].size();
 
 2616    gsSparseMatrix<T> result(size2,size1);
 
 2625    std::vector<gsSparseMatrix<T,ColMajor> > temptransfer;
 
 2626    temptransfer.resize(transfer.size());
 
 2627    for (
unsigned int i = 0; i < transfer.size();i++)
 
 2629        temptransfer[i] = transfer[i];
 
 2633    for (
unsigned int i = 0; i < old.size(); i++)
 
 2637        for(
unsigned int l =0; l < i; l++)
 
 2639            start_lv_i += n[l].size();
 
 2642        for (
unsigned int j = 0; j < old[i].size();j++)
 
 2647            for(
unsigned int l =0; l < i; l++)
 
 2649                start_lv_i += n[l].size();
 
 2651            const unsigned old_ij = old[i][j];  
 
 2653            if( n[i].bContains(old_ij) )
 
 2655                result(start_lv_i + std::distance(n[i].begin(), n[i].find_it_or_fail(old_ij) ),glob_numb ) = 1;
 
 2656                std::vector<lvl_coef> coeffs;
 
 2657                gsMatrix<index_t, d, 2> supp(d, 2);
 
 2658                this->m_bases[i]->elementSupport_into(old_ij, supp);
 
 2659                gsVector<index_t, d   > low = supp.col(0);
 
 2660                gsVector<index_t, d   > upp = supp.col(1);
 
 2663                    this->_knotIndexToDiadicIndex(i,low);
 
 2664                    this->_knotIndexToDiadicIndex(i,upp);
 
 2666                unsigned max_lvl = math::min<index_t>( this->m_tree.query4(low,upp, i), transfer.size()) ;
 
 2672                if(temp.lvl<n.size()-1)
 
 2674                    coeffs.push_back(temp);
 
 2676                for (
size_t coeff_index = 0; coeff_index < coeffs.size(); ++coeff_index)
 
 2678                    const lvl_coef coeff = coeffs[coeff_index];
 
 2681                    for(
unsigned int l =0; l < coeff.lvl; l++)
 
 2683                        start_lv_i += n[l].size();
 
 2686                    for(
typename gsSparseMatrix<T,ColMajor>::InnerIterator k(temptransfer[coeff.lvl],coeff.pos); k; ++k)
 
 2692                        if( (coeff.lvl+1) < old.size())
 
 2694                            if(old[coeff.lvl+1].bContains(k.row()))
 
 2699                        if(n[coeff.lvl+1].bContains(k.row()))
 
 2703                                const int pos = start_lv_i + n[coeff.lvl].size() + std::distance(n[coeff.lvl+1].begin(), n[coeff.lvl+1].find_it_or_fail(k.row()));
 
 2705                                result(pos,glob_numb) += coeff.coef * temptransfer[coeff.lvl](k.row(), coeff.pos);
 
 2706                                if(coeff.lvl + 1 < max_lvl)
 
 2709                                    temp.coef = temptransfer[coeff.lvl](k.row(), coeff.pos) * coeff.coef;
 
 2710                                    temp.lvl = coeff.lvl+1;
 
 2711                                    coeffs.push_back(temp);
 
 2719                                if( coeff.lvl + 1< max_lvl)
 
 2722                                    temp.coef = temptransfer[coeff.lvl](k.row(), coeff.pos) * coeff.coef;
 
 2723                                    temp.lvl = coeff.lvl+1;
 
 2724                                    coeffs.push_back(temp);
 
 2739                gsMatrix<index_t, d, 2> supp(d, 2);
 
 2740                this->m_bases[i]->elementSupport_into(old_ij, supp);
 
 2742                gsVector<index_t, d   > low = supp.col(0);
 
 2743                gsVector<index_t, d   > upp = supp.col(1);
 
 2746                    this->_knotIndexToDiadicIndex(i,low);
 
 2747                    this->_knotIndexToDiadicIndex(i,upp);
 
 2751                    math::min<index_t>( this->m_tree.query4(low,upp, i), transfer.size() ) ;
 
 2752                std::vector<lvl_coef> coeffs;
 
 2758                coeffs.push_back(temp);
 
 2759                for (
size_t coeff_index = 0; coeff_index < coeffs.size(); ++coeff_index)
 
 2761                    const lvl_coef coeff = coeffs[coeff_index];
 
 2764                    for(
unsigned int l =0; l < coeff.lvl; l++)
 
 2766                        start_lv_i += n[l].size();
 
 2769                    for(
typename gsSparseMatrix<T,ColMajor>::InnerIterator k(temptransfer[coeff.lvl],coeff.pos); k; ++k)
 
 2775                        if( (coeff.lvl+1) < (old.size()))
 
 2777                            if(old[coeff.lvl+1].bContains(k.row()))
 
 2783                        if(n[coeff.lvl+1].bContains(k.row()))
 
 2787                                const int pos = start_lv_i + n[coeff.lvl].size() + std::distance(n[coeff.lvl+1].begin(), n[coeff.lvl+1].find_it_or_fail(k.row()));
 
 2791                                result(pos,glob_numb) += coeff.coef * temptransfer[coeff.lvl](k.row(), coeff.pos);
 
 2792                                if(coeff.lvl < max_lvl-1)
 
 2795                                    temp.coef = temptransfer[coeff.lvl](k.row(), coeff.pos) * coeff.coef;
 
 2796                                    temp.lvl = coeff.lvl+1;
 
 2798                                    coeffs.push_back(temp);
 
 2806                                if(coeff.lvl < max_lvl-1)
 
 2809                                    temp.coef = temptransfer[coeff.lvl](k.row(), coeff.pos) * coeff.coef;
 
 2810                                    temp.lvl = coeff.lvl+1;
 
 2812                                    coeffs.push_back(temp);
 
 2836template<
short_t d, 
class T>
 
 2837class gsXml< gsTHBSplineBasis<d,T> >
 
 2842    GSXML_COMMON_FUNCTIONS(gsTHBSplineBasis<TMPLA2(d,T)>);
 
 2843    static std::string tag () { 
return "Basis"; }
 
 2844    static std::string type () { 
return "THBSplineBasis"+ (d>1 ? 
to_string(d):
""); }
 
 2846    static gsTHBSplineBasis<d,T> * get (gsXmlNode * node)
 
 2848        return getHTensorBasisFromXml< gsTHBSplineBasis<d,T> > (node);
 
 2851    static gsXmlNode * put (
const gsTHBSplineBasis<d,T> & obj,
 
 2854        return putHTensorBasisToXml< gsTHBSplineBasis<d,T> > (obj, data);
 
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
index_t size() const
The number of basis functions in this basis.
Definition gsHTensorBasis.hpp:298
virtual gsMatrix< index_t > boundaryOffset(boxSide const &s, index_t offset) const
Definition gsHTensorBasis.hpp:1471
virtual void refineElements(std::vector< index_t > const &boxes)
Insert the given boxes into the quadtree.
Definition gsHTensorBasis.hpp:844
Class for representing a knot vector.
Definition gsKnotVector.h:80
unsigned lastKnotIndex(const size_t &i) const
Definition gsKnotVector.h:844
void getUniformRefinementKnots(mult_t knotsPerSpan, knotContainer &result, mult_t mult=1) const
Definition gsKnotVector.hpp:1048
size_t size() const
Number of knots (including repetitions).
Definition gsKnotVector.h:242
iterator beginAt(mult_t upos) const
Definition gsKnotVector.hpp:131
int degree() const
Returns the degree of the knot vector.
Definition gsKnotVector.hpp:700
size_t uSize() const
Number of unique knots (i.e., without repetitions).
Definition gsKnotVector.h:245
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
iterator endAt(mult_t upos) const
Definition gsKnotVector.hpp:139
unsigned firstKnotIndex(const size_t &i) const
Definition gsKnotVector.h:835
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
void blockTransposeInPlace(const index_t colBlock)
Transposes in place the matrix block-wise. The matrix is.
Definition gsMatrix.h:468
Container class for a set of geometry patches and their topology, that is, the interface connections ...
Definition gsMultiPatch.h:100
index_t addPatch(typename gsGeometry< T >::uPtr g)
Add a patch from a gsGeometry<T>::uPtr.
Definition gsMultiPatch.hpp:211
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
Sparse vector class, based on gsEigen::SparseVector.
Definition gsSparseVector.h:35
Truncated hierarchical B-spline basis.
Definition gsTHBSplineBasis.h:36
void deriv2_into(const gsMatrix< T > &u, gsMatrix< T > &result) const
Evaluate the second derivatives of all active basis function at points u.
Definition gsTHBSplineBasis.hpp:1541
BoundaryBasisType * basisSlice(index_t dir_fixed, T par) const
Gives back the basis at a slice in dir_fixed at par.
Definition gsTHBSplineBasis.hpp:76
unsigned _updateSizeOfCoefs(const unsigned clevel, const unsigned flevel, const gsVector< index_t, d > &finest_low, const gsVector< index_t, d > &finest_high, gsVector< index_t, d > &size_of_coefs)
We get current size of the coefficients. Function updates this sizes accordingly to the refinement fr...
Definition gsTHBSplineBasis.hpp:424
void decomposeDomain(typename gsTHBSplineBasis::AxisAlignedBoundingBox &boundaryAABB, typename gsTHBSplineBasis::TrimmingCurves &trimCurves) const
Decomposes domain of the THB-Spline-Basis into partitions.
Definition gsTHBSplineBasis.hpp:1951
void _truncate(gsMatrix< T > &coefs, const gsVector< index_t, d > &act_size_of_coefs, const gsVector< index_t, d > &size_of_coefs, const unsigned level, const gsVector< index_t, d > &bspl_vec_ti, const unsigned bspl_vec_ti_level, const gsVector< index_t, d > &finest_low)
Performs truncation.
Definition gsTHBSplineBasis.hpp:339
void _saveNewBasisFunPresentation(const gsMatrix< T > &coefs, const gsVector< index_t, d > &act_size_of_coefs, const unsigned j, const unsigned pres_level, const gsVector< index_t, d > &finest_low)
Saves a presentation of the j-th basis function. Presentation is given by the coefficients coefs....
Definition gsTHBSplineBasis.hpp:241
void getBsplinePatches_trimming(const gsMatrix< T > &geom_coef, gsMatrix< T > &cp, gsMatrix< index_t > &b1, gsMatrix< index_t > &b2, gsVector< index_t > &level, gsMatrix< index_t > &nvertices, std::vector< std::vector< std::vector< std::vector< T > > > > &trim_curves) const
Return the list of B-spline patches to represent a THB-spline geometry.
Definition gsTHBSplineBasis.hpp:717
gsSparseMatrix< T > coarsening(const std::vector< gsSortedVector< index_t > > &old, const std::vector< gsSortedVector< index_t > > &n, const gsSparseMatrix< T, RowMajor > &transfer) const
returns a transfer matrix using the characteristic matrix of the old and new basis
Definition gsTHBSplineBasis.hpp:2356
void evalSingle_into(index_t i, const gsMatrix< T > &u, gsMatrix< T > &result) const
Evaluate the i-th basis function at points u into result.
Definition gsTHBSplineBasis.hpp:1462
void _representBasisFunction(const unsigned j, const unsigned pres_level, const gsVector< index_t, d > &finest_low, const gsVector< index_t, d > &finest_high)
Computes representation of j-th basis function on pres_level and saves it.
Definition gsTHBSplineBasis.hpp:143
gsMultiPatch< T > getBsplinePatchesToMultiPatch_trimming(const gsMatrix< T > &geom_coef, std::vector< std::vector< std::vector< std::vector< T > > > > &trim_curves) const
Return a multipatch structure of B-splines.
Definition gsTHBSplineBasis.hpp:887
util::conditional< d==1, gsConstantBasis< T >, gsTHBSplineBasis< static_cast< short_t >(d-1), T > >::type BoundaryBasisType
Associated Boundary basis type.
Definition gsTHBSplineBasis.h:57
void globalRefinement(const gsMatrix< T > &thbCoefs, int level, gsMatrix< T > &lvlCoefs) const
Returns a representation of thbCoefs as tensor-product B-spline coefficientes lvlCoefs at level level...
Definition gsTHBSplineBasis.hpp:1037
void representBasis()
Computes and saves representation of all basis functions.
Definition gsTHBSplineBasis.hpp:96
unsigned _basisFunIndexOnLevel(const gsVector< index_t, d > &index, const unsigned level, const gsVector< index_t, d > &fin_low, const unsigned new_level)
Computes tensor index of a basis function on a finer level (new_level) which is presented with tensor...
Definition gsTHBSplineBasis.hpp:300
void evalAllDers_into(const gsMatrix< T > &u, int n, std::vector< gsMatrix< T > > &result, bool sameElement=false) const
Evaluate the nonzero functions and their derivatives up to order n at points u into result.
Definition gsTHBSplineBasis.hpp:1673
void eval_into(const gsMatrix< T > &u, gsMatrix< T > &result) const
Evaluates nonzero basis functions at point u into result.
Definition gsTHBSplineBasis.hpp:1510
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 gsTHBSplineBasis.hpp:1417
void transferbyLvl(std::vector< gsSparseMatrix< T > > &result)
returns transfer matrices betweend the levels of the given hierarchical spline
Definition gsTHBSplineBasis.hpp:2311
void findNewAABB(const std::vector< std::vector< T > > &polyline, std::vector< index_t > &aabb) const
Finds new axis aligned bounding box for given polyline.
Definition gsTHBSplineBasis.hpp:2243
gsTensorBSpline< d, T > getBSplinePatch(const std::vector< index_t > &boundingBox, const unsigned level, const gsMatrix< T > &geomCoefs) const
Returns a tensor B-Spline patch defined by boundingBox.
Definition gsTHBSplineBasis.hpp:2094
void breakCycles(typename gsTHBSplineBasis::AxisAlignedBoundingBox &aabb, typename gsTHBSplineBasis::Polylines &polylines) const
Breaks the cycles of polylines and returns updated polylines.
Definition gsTHBSplineBasis.hpp:2101
gsMultiPatch< T > getBsplinePatchesToMultiPatch(const gsMatrix< T > &geom_coef) const
Return a multipatch structure of B-splines.
Definition gsTHBSplineBasis.hpp:579
void breakPolylineIntoTwoParts(const std::vector< std::vector< T > > &polyline, const index_t segment, const std::pair< T, T > &pt, std::vector< std::vector< T > > &part1, std::vector< std::vector< T > > &part2) const
Breaks polyline into two parts.
Definition gsTHBSplineBasis.hpp:2189
gsMatrix< index_t > boundaryOffset(boxSide const &s, index_t offset) const
Definition gsTHBSplineBasis.hpp:32
void getBsplinePatches(const gsMatrix< T > &geom_coef, gsMatrix< T > &cp, gsMatrix< index_t > &b1, gsMatrix< index_t > &b2, gsVector< index_t > &level, gsMatrix< index_t > &nvertices) const
Return the list of B-spline patches to represent a THB-spline geometry.
Definition gsTHBSplineBasis.hpp:522
index_t identifyCycle(const std::vector< std::vector< T > > &polyline, std::pair< T, T > &pt) const
Identify if the polyline can be split into two cycles.
Definition gsTHBSplineBasis.hpp:2141
void deriv_into(const gsMatrix< T > &u, gsMatrix< T > &result) const
Evaluates the first partial derivatives of the nonzero basis function.
Definition gsTHBSplineBasis.hpp:1568
GISMO_MAKE_GEOMETRY_NEW void getBsplinePatchGlobal(gsVector< index_t > b1, gsVector< index_t > b2, unsigned level, const gsMatrix< T > &geom_coef, gsMatrix< T > &cp, gsKnotVector< T > &k1, gsKnotVector< T > &k2) const
Returns the B-spline representation of a THB-spline subpatch.
Definition gsTHBSplineBasis.hpp:2085
void deriv2Single_into(index_t i, const gsMatrix< T > &u, gsMatrix< T > &result) const
Evaluate the (partial) derivatives of the i-th basis function at points u into result.
Definition gsTHBSplineBasis.hpp:1486
void derivSingle_into(index_t i, const gsMatrix< T > &u, gsMatrix< T > &result) const
Evaluates the (partial) derivatives of the i-th basis function at points u into result.
Definition gsTHBSplineBasis.hpp:1593
A tensor product B-spline basis.
Definition gsTensorBSplineBasis.h:37
memory::unique_ptr< Self_t > uPtr
Smart pointer for gsTensorBSplineBasis.
Definition gsTensorBSplineBasis.h:69
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
A tensor product of d B-spline functions, with arbitrary target dimension.
Definition gsTensorBSpline.h:45
memory::unique_ptr< gsTensorBSpline > uPtr
Unique pointer for gsTensorBSpline.
Definition gsTensorBSpline.h:62
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
T at(index_t i) const
Returns the i-th element of the vector.
Definition gsVector.h:177
void gsTensorBoehmRefineLocal(KnotVectorType &knots, const unsigned index, Mat &coefs, gsVector< index_t, d > &nmb_of_coefs, const gsVector< index_t, d > &act_size_of_coefs, const gsVector< index_t, d > &size_of_coefs, const unsigned direction, ValIt valBegin, ValIt valEnd, const bool update_knots)
Local refinement algorithm.
Definition gsBoehm.hpp:501
unsigned numCompositions(int sum, short_t dim)
Number of compositions of sum into dim integers.
Definition gsCombinatorics.h:691
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
Boehm's algorithm for knot insertion.
#define short_t
Definition gsConfig.h:35
#define index_t
Definition gsConfig.h:32
Implementation of deBoor and tensor deBoor algorithm.
#define gsDebug
Definition gsDebug.h:61
#define GISMO_ENSURE(cond, message)
Definition gsDebug.h:102
#define GISMO_ASSERT(cond, message)
Definition gsDebug.h:89
Provides declaration of the MultiPatch class.
Provides implementation of generic XML functions.
Provides declaration of input/output XML utilities struct.
std::string to_string(const unsigned &i)
Helper to convert small unsigned to string.
Definition gsXml.cpp:74
The G+Smo namespace, containing all definitions for the library.
S give(S &x)
Definition gsMemory.h:266
void gsBoehmRefine(KnotVectorType &knots, Mat &coefs, int p, ValIt valBegin, ValIt valEnd, bool update_knots=true)
Definition gsBoehm.hpp:163