G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsTHBSplineBasis.hpp
Go to the documentation of this file.
1
14#pragma once
15
16#include <gsCore/gsMultiPatch.h>
17
18#include <gsNurbs/gsBoehm.h>
19#include <gsNurbs/gsDeboor.hpp>
20
21#include <gsIO/gsXml.h>
23
25
26namespace gismo
27{
28
29
30template<short_t d, class T>
32boundaryOffset(boxSide const & s,index_t offset) const
33{
34 if (1!=offset)
36 else
37 {
38 //get information on the side
39 const index_t k = s.direction();
40 const bool par = s.parameter();
41
42 const index_t shift = ( par ? -1 : 1);
43
44 std::vector<index_t> temp;
46 index_t hi;
47 // i goes through all levels of the hierarchical basis
48 for(size_t i = 0; i != m_xmatrix.size(); i++)
49 {
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!");
53
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++)
57 {
58 ind = this->m_bases[i]->tensorIndex(*it);
59 if ( ind[k]==r )
60 {
61 ind[k]+=shift;
62 hi = this->flatTensorIndexToHierachicalIndex(this->m_bases[i]->index(ind),i);
63
64 GISMO_ASSERT(hi!=-1,"Neightbouring basis function with coordinates "<<ind.transpose()<<" of level "<<i<<" does not exist.");
65
66 temp.push_back(hi);
67 }
68 }
69 }
70 return makeMatrix<index_t>(temp.begin(),temp.size(),1 );
71 }
72
73}
74
75template<short_t d, class T>
77{
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);
81 const typename gsTensorBSplineBasis<d,T>::BoundaryBasisType::uPtr bBSplineBasis =
82 this->m_bases[0]->boundaryBasis(side);
84 new typename gsTHBSplineBasis<d,T>::BoundaryBasisType(*bBSplineBasis);//,this->m_tree.getMaxInsLevel()+1);
85
86 if(d!=1)
87 {
88 std::vector<index_t> boxes;
89 this->getBoxesAlongSlice(dir_fixed,par,boxes);
90 bBasis->refineElements(boxes);
91 }
92 return bBasis;
93}
94
95template<short_t d, class T>
97{
98 // Cleanup previous basis
99 this->m_is_truncated.resize(this->size());
100 m_presentation.clear();
101
102 gsMatrix<index_t, d, 2> element_ind(d, 2);
103 gsVector<index_t, d > low, high;
104 for (index_t j = 0; j < this->size(); ++j)
105 {
106 index_t level = this->levelOf(j);
107 index_t tensor_index = this->flatTensorIndexOf(j, level);
108
109 // element indices
110 this->m_bases[level]->elementSupport_into(tensor_index, element_ind);
111
112 // I tried with block, I can not trick the compiler to use references
113 low = element_ind.col(0); //block<d, 1>(0, 0);
114 high = element_ind.col(1); //block<d, 1>(0, 1);
115 if (m_manualLevels)
116 {
117 this->_knotIndexToDiadicIndex(level,low);
118 this->_knotIndexToDiadicIndex(level,high);
119 }
120
121 // Finds coarsest level that function, with supports given with
122 // support indices of the coarsest level (low & high), has presentation
123 // based only on B-Splines (and not THB-Splines).
124 // this is not the same as query 3
125 index_t clevel = this->m_tree.query4(low, high, level);
126
127 if (level != clevel) // we must compute its presentation
128 {
129 this->m_tree.computeFinestIndex(low, level, low);
130 this->m_tree.computeFinestIndex(high, level, high);
131
132 this->m_is_truncated[j] = clevel;
133 _representBasisFunction(j, clevel, low, high);
134 }
135 else
136 {
137 this->m_is_truncated[j] = -1;
138 }
139 }
140}
141
142template<short_t d, class T>
144 const unsigned j,
145 const unsigned pres_level,
146 const gsVector<index_t, d>& finest_low,
147 const gsVector<index_t, d>& finest_high)
148{
149 const unsigned cur_level = this->levelOf(j);
150
151 // actual size of the coefficients
152 gsVector<index_t, d> act_size_of_coefs(d);
153 act_size_of_coefs.fill(1);
154
155 // number of new coefficients
156 unsigned nmb_of_coefs = _updateSizeOfCoefs(cur_level, pres_level,
157 finest_low, finest_high,
158 act_size_of_coefs);
159 gsMatrix<T> coefs(nmb_of_coefs, 1);
160 coefs.fill(0);
161 coefs.row(0).setOnes();
162
163 // vector of the numbers of the coefficients (in each dimension)
164 // stored in coefs
165 gsVector<index_t, d> vec_nmb_of_coefs(d);
166 vec_nmb_of_coefs.fill(1);
167
168 unsigned tensor_index = this->flatTensorIndexOf(j, cur_level);
169
170 // B-Spline vector tensor index
171 gsVector<index_t, d> bspl_vec_ti =
172 this->m_bases[cur_level]->tensorIndex(tensor_index);
173
174 // we need to separately save knot vectors because we will modify
175 // them, when we proceed from one level on another
176 std::vector<gsKnotVector<T> > vector_of_kv(d);
177
178 // size of the coefficients that are affected in individual iteration
179 gsVector<index_t, d> cur_size_of_coefs(d);
180 cur_size_of_coefs.fill(1);
181
182 for (unsigned level = cur_level; level < pres_level; ++level)
183 {
184 _updateSizeOfCoefs(level, level + 1, finest_low,
185 finest_high, cur_size_of_coefs);
186
187
188 // index of a support of the j-th basis function (l_low, l_high
189 // on level, and l1_high, l1_low on level + 1)
190 gsVector<index_t, d> clow, chigh, fhigh, flow;
191
192 this->m_tree.computeLevelIndex(finest_low, level, clow);
193 this->m_tree.computeLevelIndex(finest_high, level, chigh);
194 if (m_manualLevels)
195 {
196 this->_diadicIndexToKnotIndex(level,clow);
197 this->_diadicIndexToKnotIndex(level,chigh);
198 }
199 this->m_tree.computeLevelIndex(finest_low, level + 1, flow);
200 this->m_tree.computeLevelIndex(finest_high, level + 1, fhigh);
201 if (m_manualLevels)
202 {
203 this->_diadicIndexToKnotIndex(level + 1,flow);
204 this->_diadicIndexToKnotIndex(level + 1,fhigh);
205 }
206
207 std::vector<T> knots;
208
209 for (unsigned dim = 0; dim < d; ++dim)
210 {
211 const gsKnotVector<T>& ckv = m_bases[level ]->knots(dim);
212 const gsKnotVector<T>& fkv = m_bases[level+1]->knots(dim);
213
214 if (level == cur_level)
215 vector_of_kv[dim] = ckv;
216
217 knots.clear();
218 std::set_symmetric_difference(ckv.beginAt(clow[dim]), ckv.endAt(chigh[dim]),
219 fkv.beginAt(flow[dim]), fkv.endAt(fhigh[dim]),
220 std::back_inserter(knots));
221
225 typename std::vector<T>::const_iterator>
226 (vector_of_kv[dim], bspl_vec_ti[dim], coefs, vec_nmb_of_coefs,
227 act_size_of_coefs,
228 cur_size_of_coefs, dim, knots.begin(), knots.end(),
229 true);
230 }
231
232 _truncate(coefs, act_size_of_coefs, cur_size_of_coefs,
233 level + 1, bspl_vec_ti, cur_level, finest_low);
234 }
235 _saveNewBasisFunPresentation(coefs, act_size_of_coefs,
236 j, pres_level, finest_low);
237}
238
239
240template<short_t d, class T>
242 const gsMatrix<T>& coefs,
243 const gsVector<index_t, d>& act_size_of_coefs,
244 const unsigned j,
245 const unsigned pres_level,
246 const gsVector<index_t, d>& finest_low)
247{
248 const unsigned level = this->levelOf(j);
249 const unsigned tensor_index = this->flatTensorIndexOf(j, level);
250
251 gsVector<index_t, d> bspl_vec_ti =
252 this->m_bases[level]->tensorIndex(tensor_index);
253
254 // finer tensor index
255 const unsigned f_ten_index = _basisFunIndexOnLevel(bspl_vec_ti, level,
256 finest_low, pres_level);
257
258 gsVector<index_t, d> act_coefs_strides(d);
259 bspline::buildCoeffsStrides<d>(act_size_of_coefs, act_coefs_strides);
260
261
262 gsVector<index_t, d> position(d);
263 position.fill(0);
264
265
266 gsVector<index_t, d> first_point(position);
267 gsVector<index_t, d> last_point(d);
268 bspline::getLastIndexLocal<d>(act_size_of_coefs, last_point);
269
270 this->m_presentation[j] =
271 gsSparseVector<T>(this->m_bases[pres_level]->size());
272
273 gsSparseVector<T>& presentation = this->m_presentation[j];
274 //presentation.reserve(coefs.rows()); // reserve memory ?
275
276 do
277 {
278 // ten_index - (tensor) index of a bspline function with respect to
279 // the coef at position "position"
280 // coef_index - (local) index of a ceofficient at position
281
282 unsigned ten_index = f_ten_index;
283 for (unsigned dim = 0; dim < d; dim++)
284 {
285 ten_index += position(dim) *
286 this->m_bases[pres_level]->stride(static_cast<int>(dim));
287 }
288
289 unsigned coef_index = bspline::getIndex<d>(act_coefs_strides, position);
290
291 if (coefs(coef_index) != 0)
292 presentation(ten_index) = coefs(coef_index);
293
294 } while(nextCubePoint<gsVector<index_t, d> > (position, first_point,
295 last_point));
296}
297
298
299template<short_t d, class T>
301 const gsVector<index_t, d>& index,
302 const unsigned level,
303 const gsVector<index_t, d>& fin_low,
304 const unsigned new_level)
305{
307 this->m_tree.computeLevelIndex(fin_low, level, low);
308
309 gsVector<index_t, d> flow(d);
310 this->m_tree.computeLevelIndex(fin_low, new_level, flow);
311
312 if (m_manualLevels)
313 {
314 this->_diadicIndexToKnotIndex(level,low);
315 this->_diadicIndexToKnotIndex(new_level,flow);
316 }
317
318 gsVector<index_t, d> new_index(d);
319
320 for (unsigned dim = 0; dim < d; dim++)
321 {
322 const gsKnotVector<T>& ckv =
323 this->m_bases[level]->knots(dim);
324
325 const gsKnotVector<T>& fkv =
326 this->m_bases[new_level]->knots(dim);
327
328
329 unsigned mult = index[dim] - ckv.firstKnotIndex(low[dim]);
330
331 new_index(dim) = fkv.firstKnotIndex(flow[dim]) + mult;
332 }
333
334 return this->m_bases[new_level]->index(new_index);
335}
336
337
338template<short_t d, class T>
340 gsMatrix<T>& coefs,
341 const gsVector<index_t, d>& act_size_of_coefs,
342 const gsVector<index_t, d>& size_of_coefs,
343 const unsigned level,
344 const gsVector<index_t, d>& bspl_vec_ti,
345 const unsigned bspl_vec_ti_level,
346 const gsVector<index_t, d>& finest_low)
347{
348 // if we dont have any active function in this level, we do not truncate
349 if (this->m_xmatrix[level].size() == 0)
350 return;
351
352 // global tensor index
353 const unsigned const_ten_index = _basisFunIndexOnLevel(bspl_vec_ti,
354 bspl_vec_ti_level, finest_low, level);
355 gsVector<index_t, d> act_coefs_strides(d);
356 bspline::buildCoeffsStrides<d>(act_size_of_coefs, act_coefs_strides);
357
358
359 gsVector<index_t, d> last_point(d);
360 bspline::getLastIndexLocal<d>(size_of_coefs, last_point);
361 last_point(0) = 0;
362
363
364 gsVector<index_t, d> position(d);
365 position.fill(0);
366
367 gsVector<index_t, d> first_point(position);
368
369 unsigned xmatrix_index = 0;
370 unsigned tensor_active_index = this->m_xmatrix[level][0];
371
372 unsigned numb_of_point = size_of_coefs[0];
373
374 do
375 {
376 // indices
377 // ten_index - (tensor) index of a bspline function with respect to
378 // the coef at position "position"
379 // coef_index - (local) index of a ceofficient at position
380 // tensor_active_index - (tensor) index of a active functions
381
382
383 unsigned ten_index = const_ten_index;
384 for (unsigned dim = 1; dim < d; dim++)
385 {
386 ten_index += position(dim) *
387 this->m_bases[level]->stride(static_cast<int>(dim));
388 }
389
390 unsigned coef_index = bspline::getIndex<d>(act_coefs_strides, position);
391
392 for (unsigned index = 0; index < numb_of_point; ++index)
393 {
394 if (tensor_active_index < ten_index)
395 {
396 while(tensor_active_index < ten_index)
397 {
398 xmatrix_index++;
399 if (xmatrix_index == this->m_xmatrix[level].size())
400 {
401 // we don't have any active basis function,
402 // so all the rest basis function in our
403 // representation should not be truncated
404 return;
405 }
406
407 tensor_active_index =
408 this->m_xmatrix[level][xmatrix_index];
409 }
410 // ten_index <= tensor_active_index holds
411 }
412 if (ten_index == tensor_active_index) // truncate
413 coefs(coef_index + index, 0) = 0;
414
415 ten_index++;
416 }
417
418 } while(nextCubePoint<gsVector<index_t, d> >(position, first_point,
419 last_point));
420}
421
422
423template<short_t d, class T>
425 const unsigned clevel,
426 const unsigned flevel,
427 const gsVector<index_t, d>& finest_low,
428 const gsVector<index_t, d>& finest_high,
429 gsVector<index_t, d>& size_of_coefs)
430{
431 gsVector<index_t, d> clow, chigh;
432
433 this->m_tree.computeLevelIndex(finest_low, clevel, clow);
434 this->m_tree.computeLevelIndex(finest_high, clevel, chigh);
435
436 gsVector<index_t, d> flow, fhigh;
437 this->m_tree.computeLevelIndex(finest_low, flevel, flow);
438 this->m_tree.computeLevelIndex(finest_high, flevel, fhigh);
439
440 if (m_manualLevels)
441 {
442 this->_diadicIndexToKnotIndex(clevel,clow);
443 this->_diadicIndexToKnotIndex(clevel,chigh);
444 this->_diadicIndexToKnotIndex(flevel,flow);
445 this->_diadicIndexToKnotIndex(flevel,fhigh);
446 }
447 // number of new coefficients
448 unsigned nmb_of_coefs = 1;
449
450 for (unsigned dim = 0; dim < d; ++dim)
451 {
452 const gsKnotVector<T>& ckv =
453 this->m_bases[clevel]->knots(dim);
454 const gsKnotVector<T>& fkv =
455 this->m_bases[flevel]->knots(dim);
456
457 // Number of knots in the coarse knot vector
458 unsigned cnmb_knts = ckv.lastKnotIndex(chigh[dim]) -ckv.firstKnotIndex(clow[dim]);
459
460 // Number of knots in the fine knot vector
461 unsigned fnmb_knts = fkv.lastKnotIndex(fhigh[dim]) -fkv.firstKnotIndex(flow[dim]);
462
463 size_of_coefs(dim) += fnmb_knts - cnmb_knts;
464 nmb_of_coefs *= size_of_coefs(dim);
465 }
466 return nmb_of_coefs;
467}
468
469// return the B-spline representation of a THB-spline subpatch
470template<short_t d, class T>
471template<short_t dd>
472typename util::enable_if<dd==2,void>::type
475 unsigned level,
476 const gsMatrix<T>& geom_coef,
477 gsMatrix<T>& cp,
478 gsKnotVector<T>& k1,
479 gsKnotVector<T>& k2) const
480{
481 // check if the indices in b1, and b2 are correct with respect to the given level
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);
487
488 // select the indices of all B-splines of the given level acting on the given box
489 gsVector<index_t,d> b1_outputs, b2_outputs;
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;
500
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());
504
505 gsMatrix<T> temp;
506 globalRefinement(geom_coef, level, temp);
507
508 index_t cc = 0;
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);
512
513 // compute the new vectors for the B-spline patch
514 k1 = gsKnotVector<T>(m_deg[0], m_bases[level]->knots(0).begin() + i0 ,
515 m_bases[level]->knots(0).begin() + i1 + m_deg[0] + 2);
516 k2 = gsKnotVector<T>(m_deg[1], m_bases[level]->knots(1).begin() + j0 ,
517 m_bases[level]->knots(1).begin() + j1 + m_deg[1] + 2);
518}
519
520// returns the list of B-spline patches to represent a THB-spline geometry
521template<short_t d, class T>
524 gsVector<index_t>& level, gsMatrix<index_t>& nvertices) const
525{
526 this->m_tree.getBoxes(b1,b2,level); // splitting based on the quadtree
527 int nboxes = level.size();
528 //------------------------------------------------------------------------------------------------------------------------------
529 // iteration on the boxes to call getBsplinePatchGlobal()
530 //------------------------------------------------------------------------------------------------------------------------------
531 gsVector<index_t> p1, p2;
532 p1.resize(this->dim());
533 p2.resize(this->dim());
534 gsMatrix<T> temp1, temp2;
535 gsKnotVector<T> cku, ckv;
536 nvertices.resize(nboxes,this->dim());
537
538 for (int i = 0; i < nboxes; i++)
539 {
540 p1 = b1.row(i).transpose();
541 p2 = b2.row(i).transpose();
542
543 this->getBsplinePatchGlobal(p1, p2, level[i], geom_coef, temp1, cku, ckv);
544
545 if (i == 0)
546 {
547 cp = temp1;
548 }
549 else
550 {
551 int cprows = cp.rows();
552 temp2.resize(cp.rows()+temp1.rows(),cp.cols());
553
554 for (int j = 0; j < cp.rows(); j++)
555 {
556 for (int k = 0; k < cp.cols(); k++)
557 {
558 temp2(j,k) = cp(j,k);
559 }
560 }
561
562 for (int j = 0; j < temp1.rows(); j++)
563 {
564 for (int k = 0; k < temp1.cols(); k++)
565 {
566 temp2(cprows+j,k) = temp1(j,k);
567 }
568 }
569 cp = temp2;
570 }
571
572 nvertices(i,0) = cku.size()-cku.degree()-1;
573 nvertices(i,1) = ckv.size()-ckv.degree()-1;
574 }
575}
576
577// returns the list of B-spline patches to represent a THB-spline geometry
578template<short_t d, class T>
580{
581 GISMO_ASSERT(d==2,"Dim must be 2 for now");
582
583 gsMultiPatch<T> result;
584
585 gsMatrix<index_t> b1, b2;
586 gsVector<index_t> level;
587 this->m_tree.getBoxes(b1,b2,level); // splitting based on the quadtree
588
589 const int nboxes = level.size();
590 gsVector<index_t> p1, p2;
591 gsMatrix<T> temp1;
592 gsKnotVector<T> cku, ckv;
593
594 for (int i = 0; i < nboxes; i++) // for all boxes
595 {
596 p1 = b1.row(i).transpose();
597 p2 = b2.row(i).transpose();
598
599 this->getBsplinePatchGlobal(p1, p2, level[i], geom_coef, temp1, cku, ckv);
600 result.addPatch(typename gsTensorBSpline<2, T>::uPtr(new gsTensorBSpline<2, T>(cku, ckv, give(temp1))));
601 }
602
603 return result;
604}
605
606template<short_t d, class T>
608 std::vector<std::vector<std::vector< std::vector<index_t> > > >& connectedComponents, gsVector<index_t>& level) const
609{
610 //identify the outer polylines- conected components
611 int first_level = 0;
612 for(unsigned int i = 0; i < this->m_xmatrix.size(); i++)
613 {
614 if(this->m_xmatrix[i].size()>0)
615 {
616 first_level = i;
617 break;
618 }
619 }
620 gsDebug<<"new min level"<<"\n";
621 std::vector< std::vector< std::vector< std::vector<index_t> > > > res; //things to assign to trim_curves
622 std::vector< std::vector< std::vector<index_t > > > aabb;//axis aligned bounding box
623 std::vector< std::vector<index_t > > boxes;
624 aabb = this->domainBoundariesIndices(res);
625 for(unsigned int i = 0; i < aabb.size(); i++)//level
626 {
627 //compare every aabb with the others
628 for(unsigned int j = 0; j < aabb[i].size(); j++)//aabb
629 {
630 bool is_boundary = true;
631 for(unsigned int k = 0; k < aabb[i].size(); k++)//aabb
632 {
633 if(j!=k)
634 {
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]) )
639 {
640 is_boundary= !is_boundary;
641 }
642 }
643 }
644 if(is_boundary)
645 {
646 boxes.push_back(aabb[i][j]);
647 boxes[boxes.size()-1].push_back(i+first_level);//adding the level information
648 connectedComponents.resize(connectedComponents.size()+1);
649 connectedComponents[connectedComponents.size()-1].push_back(res[i][j]);//trim_curves
650 }
651 }
652 }
653 //create the matrices b1,b2 and vector level
654
655 level.resize(boxes.size());
656 for(size_t i = 0; i < boxes.size(); i++){
657 level[i] = boxes[i][2*d];
658 }
659
660 // identify holes
661 for(size_t l = 0; l < aabb.size();l++) //level
662 {
663 for(size_t i = 0; i < aabb[l].size(); i++) //bb
664 {
665 int closest_box = -1;
666 for(size_t j = 0; j < boxes.size(); j++)
667 {
668 if(l == static_cast<size_t>(boxes[j][4]) )
669 {
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]))
674 {
675 if(closest_box!=-1){
676 //test with previous closest_box
677 if(!(
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])))
682 {
683 closest_box = j;
684 }
685 }
686 else
687 {
688 closest_box = j;
689 }
690
691 }
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]))
696 {
697 //boxes[j] == aabb[l][i]
698 closest_box = -1;
699 break;
700 }
701 }
702
703 }
704 //prirad s res do trim_curves
705 if(closest_box>-1)
706 {
707 connectedComponents[closest_box].push_back(res[l][i]);
708 }
709 }
710 }
711
712
713}
714
715//return data for trimming in parasolid
716template<short_t d, class T>
718 const gsMatrix<T>& geom_coef,
719 gsMatrix<T>& cp,
722 gsVector<index_t>& level,
723 gsMatrix<index_t>& nvertices,
724 std::vector<std::vector<std::vector< std::vector<T> > > >& trim_curves) const
725{
726 //identify the outer polylines- conected components
727// int first_level = 0;
728// for(unsigned int i = 0; i < this->m_xmatrix.size(); i++){
729// if(this->m_xmatrix[i].size()>0){
730// first_level = i - 1;
731// break;
732// }
733// }
734
735 // gsDebug<<"new min level"<< first_level << "\n";
736 std::vector< std::vector< std::vector< std::vector< T > > > > res; //things to assign to trim_curves
737 std::vector< std::vector< std::vector<index_t> > > aabb;//axis aligned bounding box
738 std::vector< std::vector<index_t> > boxes;
739 aabb = this->domainBoundariesParams(res);
740 for(unsigned int i = 0; i < aabb.size(); i++)//level
741 {
742 //compare every aabb with the others
743 for(unsigned int j = 0; j < aabb[i].size(); j++)//aabb
744 {
745 bool is_boundary = true;
746 for(unsigned int k = 0; k < aabb[i].size(); k++)//aabb
747 {
748 if(j!=k)
749 {
750 if( (aabb[i][j][0] > aabb[i][k][0]) && // aabb[i][j] is completely inside aabb[i][k]
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]) )
754 {
755 is_boundary= !is_boundary;
756 }
757 }
758 }
759 if(is_boundary)
760 {
761 boxes.push_back(aabb[i][j]);
762 boxes[boxes.size()-1].push_back(i);//adding the level information
763 trim_curves.resize(trim_curves.size()+1);
764 trim_curves[trim_curves.size()-1].push_back(res[i][j]);//trim_curves
765 }
766 }
767 }
768 //create the matrices b1,b2 and vector level
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++){
774 //convert from param spece to index space in highest level
775 // (!) next lines: Conversion form double to int, possible loss of data
776 b1(i,j) = boxes[i][j];
777 b2(i,j) = boxes[i][j+d];
778 }
779 level[i] = boxes[i][2*d];
780 }
781 //use getBsplinePatches
782 int nboxes = level.size();
783 //------------------------------------------------------------------------------------------------------------------------------
784 // iteration on the boxes to call getBsplinePatchGlobal()
785 //------------------------------------------------------------------------------------------------------------------------------
786
787 nvertices.resize(nboxes,this->dim());
788
789 for (int i = 0; i < nboxes; i++)
790 {
791 gsMatrix<T> new_cp;
792 gsVector<index_t> p1, p2;
793 gsKnotVector<T> cku, ckv;
794
795 p1 = b1.row(i).transpose();
796 p2 = b2.row(i).transpose();
797
798 this->getBsplinePatchGlobal(p1, p2, level[i], geom_coef, new_cp, cku, ckv);
799
800 if (i == 0)
801 {
802 cp = new_cp;
803 }
804 else
805 {
806 gsMatrix<T> bigger;
807 int cprows = cp.rows();
808 /*cp.conservativeResize( cp.rows() + temp_cp.rows(), gsEigen::NoChange );
809 for( int j=cprows; j < cp.rows(); j++ )
810 cp.row(j) = temp2.row(j-cprows);*/
811 bigger.resize(cp.rows()+new_cp.rows(), cp.cols());
812 for(int j=0; j< bigger.rows(); j++)
813 {
814 if(j<cprows)
815 {
816 bigger.row(j) = cp.row(j);
817 }
818 else
819 {
820 bigger.row(j) = new_cp.row(j-cprows);
821 }
822 }
823 cp=bigger;
824 }
825
826 nvertices(i,0) = cku.size()-cku.degree()-1;
827 nvertices(i,1) = ckv.size()-ckv.degree()-1;
828
829 }
830 // identify holes
831 for(size_t l = 0; l < aabb.size();l++) //level
832 {
833 for(size_t i = 0; i < aabb[l].size(); i++) //bb
834 {
835 int closest_box = -1;
836 for(size_t j = 0; j < boxes.size(); j++)
837 {
838 if(l == static_cast<size_t>(boxes[j][4]) )
839 {
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]))
844 {
845 if(closest_box!=-1)
846 {
847 //test with previous closest_box
848 if(!(
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])))
853 {
854 closest_box = j;
855 }
856 }
857 else
858 {
859 closest_box = j;
860 }
861
862 }
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]))
867 {
868 //boxes[j] == aabb[l][i]
869 closest_box = -1;
870 break;
871 }
872 }
873
874 }
875 //prirad s res do trim_curves
876 if(closest_box>-1)
877 {
878 trim_curves[closest_box].push_back(res[l][i]);
879 }
880 }
881 }
882}
883
884
885//return data for trimming in parasolid
886template<short_t d, class T>
888 const gsMatrix<T>& geom_coef,
889 std::vector<std::vector<std::vector< std::vector<T> > > >& trim_curves) const
890{
893 gsVector<index_t> level;
894 gsMultiPatch<T> result;
895 //identify the outer polylines- conected components
896 int first_level = 0;
897 for(size_t i = 0; i < this->m_xmatrix.size(); i++)
898 {
899 if(this->m_xmatrix[i].size()>0)
900 {
901 first_level = i;
902 break;
903 }
904 }
905 gsDebug<<"new min level"<<"\n";
906 std::vector< std::vector< std::vector< std::vector< T > > > > res; //things to assign to trim_curves
907 std::vector< std::vector< std::vector<index_t > > > aabb;//axis aligned bounding box
908 std::vector< std::vector<index_t > > boxes;
909 aabb = this->domainBoundariesParams(res);
910 for(size_t i = 0; i < aabb.size(); i++)//level
911 {
912 //compare every aabb with the others
913 for(size_t j = 0; j < aabb[i].size(); j++)//aabb
914 {
915 bool is_boundary = true;
916 for(size_t k = 0; k < aabb[i].size(); k++)//aabb
917 {
918 if(j!=k)
919 {
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]) )
924 {
925 is_boundary= !is_boundary;
926 }
927 }
928 }
929 if(is_boundary)
930 {
931 boxes.push_back(aabb[i][j]);
932 boxes[boxes.size()-1].push_back(i+first_level);//adding the level information
933 trim_curves.resize(trim_curves.size()+1);
934 trim_curves[trim_curves.size()-1].push_back(res[i][j]);//trim_curves
935 }
936 }
937 }
938 //create the matrices b1,b2 and vector level
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++)
943 {
944 for(unsigned j = 0; j < d; j++)
945 {
946 //convert from param spece to index space in highest level
947 // (!) next lines: Conversion form double to int, possible loss of data
948 b1(i,j) = boxes[i][j];
949 b2(i,j) = boxes[i][j+d];
950 }
951 level[i] = boxes[i][2*d];
952 }
953 //use getBsplinePatches
954
955 int nboxes = level.size();
956 //------------------------------------------------------------------------------------------------------------------------------
957 // iteration on the boxes to call getBsplinePatchGlobal()
958 //------------------------------------------------------------------------------------------------------------------------------
959 gsVector<index_t> p1, p2;
960 p1.resize(this->dim());
961 p2.resize(this->dim());
962 gsMatrix<T> temp1;
963 gsKnotVector<T> cku, ckv;
964
965
966 for (int i = 0; i < nboxes; i++)
967 {
968 p1 = b1.row(i).transpose();
969 p2 = b2.row(i).transpose();
970
971 this->getBsplinePatchGlobal(p1, p2, level[i], geom_coef, temp1, cku, ckv);
972 gsTensorBSplineBasis<2, T> tbasis(cku, ckv);
973 result.addPatch(typename gsTensorBSpline<2, T>::uPtr(new gsTensorBSpline<2, T>(tbasis, give(temp1))));
974
975 }
976 // identify holes
977 for(size_t l = 0; l < aabb.size();l++) //level
978 {
979 for(size_t i = 0; i < aabb[l].size(); i++) //bb
980 {
981 int closest_box = -1;
982 for(size_t j = 0; j < boxes.size(); j++)
983 {
984 if(l == static_cast<size_t>(boxes[j][4]) )
985 {
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]))
990 {
991 if(closest_box!=-1){
992 //test with previous closest_box
993 if(!(
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])))
998 {
999 closest_box = j;
1000 }
1001 }
1002 else
1003 {
1004 closest_box = j;
1005 }
1006
1007 }
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]))
1012 {
1013 //boxes[j] == aabb[l][i]
1014 closest_box = -1;
1015 break;
1016 }
1017 }
1018
1019 }
1020 //prirad s res do trim_curves
1021 if(closest_box>-1)
1022 {
1023 trim_curves[closest_box].push_back(res[l][i]);
1024 }
1025 }
1026 }
1027 return result;
1028}
1029
1030
1031/*
1032 Private functions
1033*/
1034
1035
1036template<short_t d, class T>
1038 int level, gsMatrix<T> & lvlCoefs) const
1039{
1040 const index_t n = thbCoefs.cols();
1041
1042 // Initialize level 0 coefficients
1043 lvlCoefs.setZero(m_bases[0]->size(), n);
1044 for(cmatIterator it = m_xmatrix[0].begin(); it != m_xmatrix[0].end(); ++it)
1045 {
1046 const int hIndex = m_xmatrix_offset[0] + (it - m_xmatrix[0].begin());
1047 lvlCoefs.row(*it) = thbCoefs.row(hIndex);
1048 }
1049
1050 gsKnotVector<T> k1, k2;// fixme: boehm refine needs non-const kv
1051 std::vector<T> knots_x, knots_y;
1052
1053 for(int l = 1; l <=level; l++)
1054 {
1055 k1 = m_bases[l-1]->knots(0);
1056 k2 = m_bases[l-1]->knots(1);
1057
1058 // global dyadic refinement with respect to previous level
1059 k1.getUniformRefinementKnots(1,knots_x);
1060 k2.getUniformRefinementKnots(1,knots_y);
1061
1062 // refine direction 0
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);
1065
1066 // refine direction 1
1067 lvlCoefs.blockTransposeInPlace(m_bases[l-1]->size(1));
1068 gsBoehmRefine(k2, lvlCoefs, m_deg[1], knots_y.begin(), knots_y.end(), false);
1069 lvlCoefs.blockTransposeInPlace(m_bases[l]->size(0));
1070 lvlCoefs.resize(m_bases[l]->size(), n); //lvlCoefs: control points at level \a l
1071
1072 // overwrite with the THB coefficients of level \a l
1073 for(cmatIterator it = m_xmatrix[l].begin(); it != m_xmatrix[l].end(); ++it)
1074 {
1075 const int hIndex = m_xmatrix_offset[l] + (it - m_xmatrix[l].begin());
1076 lvlCoefs.row(*it) = thbCoefs.row(hIndex);
1077 }
1078 }
1079}
1080
1081template<index_t d>
1082inline gsVector<index_t, d> tensorIndex(index_t m,
1083 const gsVector<index_t, d> & str)
1084{
1086 for (short_t i = d-1; i>=0; --i)
1087 {
1088 rvo[d-i-1] = m % str[i];
1089 m -= rvo[i];
1090 }
1091 return rvo;
1092}
1093
1094template<class Iter, class Vec> //raname as: next intersection
1095// rename: findEqualIndex
1096bool findNextMatch(Iter & it, Iter stop, Vec & cur,
1097 const Vec & low, const Vec & upp,
1098 const Vec & str)
1099{
1100 index_t ci;
1101 while (it!=stop)
1102 {
1103 ci = cur.dot(str);
1104 if (ci == *it)
1105 return true;
1106 else if ( ci < *it )
1107 {
1108 // problem: we should be inside [low,upp].
1109 // but what happens if we are not ?
1110 /*
1111 cur = tensorIndex(*it, str);
1112 for (index_t i = cur.size()-1; i>=0; ++i)
1113 {
1114 if (cur[i]>upp[i]) //impossible for i=d-1
1115 cur[i] = upp[i];
1116 else if (cur[i]<low[i])
1117 {
1118 cur[i] = low[i];
1119 if (i+1<cur.size())
1120 cur[i+1] += 1;
1121 else
1122 return false;
1123 }
1124 }
1125 // HOW TO BREAK ?
1126 //if ( (cur.array() > upp.array()).any() ) break;
1127 if (cur.dot(str) > upp.dot(str) ) break;
1128 //*/
1129
1130 if (!nextCubePoint(cur,low,upp)) break; //move by one..
1131 }
1132 else // *it < ci
1133 it = std::lower_bound(it, stop, ci);
1134 }
1135
1136 return false;
1137}
1138
1139// finds nonzero index in the iterator it.
1140// assumes [low,upp] and iterator have stride str
1141template<class T, class Vec>
1142bool findNextMatch(const gsSparseVector<T> & sv, index_t & ii,
1143 Vec & cur, const Vec & low, const Vec & upp,
1144 const Vec & str)
1145{
1146 index_t ci;
1147 while ( ii<sv.data().size() )
1148 {
1149 ci = cur.dot(str);
1150 if ( ci == sv.data().index(ii) )
1151 return true;
1152 else if ( ci < sv.data().index(ii) )
1153 {
1154 if (!nextCubePoint(cur,low,upp)) break;
1155 }
1156 else // ii < ci
1157 ii = sv.data().searchLowerIndex(ii, sv.data().size(), ci);
1158 }
1159 return false;
1160}
1161
1162 /*
1163//template<short_t d, class T>
1164//index_t gsTHBSplineBasis<d,T>::
1165template<class T>
1166active_detail(const gsMatrix<T> & u,
1167 gsMatrix<index_t> &actives,//sorted active functions
1168 gsMatrix<index_t> & offset,//level offset
1169 // ???
1170 std::vector<std:vector<index_t> > & repLvl, // [rlvl]->list(act)
1171 gsMatrix<index_t & repLevel> // for each function: [rlvl] -> []
1172 ) const
1173{
1174 point low, upp, cur, mstr, str, ll, uu, cc;
1175 index_t ii;
1176 // todo: store active_cwise/stride_cwise
1177 const int maxLevel = this->m_tree.getMaxInsLevel();
1178 gsMatrix<T> currPoint;
1179 offset.setZero(maxLevel+1, u.cols()); //trim at the end
1180 // stores [rlvl]->(act,level)
1181 std::vector<std::vector<std::pair<index_t,index_t> > > tfunction(maxLevel+1);
1182
1183 for (index_t i = 0; i != u.cols(); i++) // for all points
1184 {
1185 currPoint = u.col(i);
1186 for(short_t k = 0; k != d; ++k)
1187 low[k] = m_bases[maxLevel]->knots(k).uFind( currPoint.at(k) ).uIndex();
1188 if (m_manualLevels)
1189 this->_knotIndexToDiadicIndex(maxLevel,low);
1190 // Identify the level of the point
1191 index_t lvl = std::min(this->m_tree.levelOf(low, maxLevel),(int) m_xmatrix.size()-1);
1192 //offset.col(i).setZero();
1193
1194 typename CMatrix::const_iterator it;
1195 for(int level = 0; level <= maxLevel; level++) // for all relevant levels
1196 {
1197 if (level>lvl && tfunction[level].empty()) //nothing to do here
1198 continue; // 0 offset
1199
1200 m_bases[level]->active_cwise(currPoint, low, upp);//my be improved: start from finest lvl
1201 m_bases[level]->stride_cwise(mstr);
1202
1203 if ( level<= lvl)
1204 {
1205 it = m_xmatrix[level].begin();
1206 cur = low;
1207 while ( findNextMatch(it, m_xmatrix[level].end(), cur, low, upp, mstr) )
1208 {
1209 const index_t act = this->m_xmatrix_offset[level] + (it - m_xmatrix[level].begin());//flattensortohier ?
1210 if ( isTruncated(act) )
1211 {
1212 //Record truncated functions with their representation level
1213 const int rlvl = m_is_truncated.at(act);
1214 tfunction[rlvl].push_back( std::make_pair(act,level) );
1215 }
1216 else
1217 ++offset(level, i);
1218 ++it;//advance both
1219 nextCubePoint(cur,low,upp);
1220 }
1221 }
1222
1223 // count truncated active functions
1224 if ( !tfunction[level].empty() )
1225 {
1226 for ( std::pair<index_t,index_t> & q : tfunction[level] )
1227 {
1228 const gsSparseVector<T>& coefs = getCoefs(q.first);
1229 typename gsSparseVector<T>::InnerIterator it(coefs);
1230 cur = low;
1231 ii = 0;
1232 if ( findNextMatch(coefs, ii, cur, low, upp, mstr) )
1233 ++offset(q.second, i);
1234 }
1235 tfunction[level].clear();
1236 }
1237 }// end level
1238 }//end points
1239
1240 // ii =maxLevel;
1241 // while (offset[ii]==) --ii;
1242 // offset.conservativeresize(ii,gsEigen::NoChange);
1243
1244 return offset.sum();
1245}
1246// */
1247
1248template<short_t d, class T>
1249index_t gsTHBSplineBasis<d,T>::numActiveMax(const gsMatrix<T> & u,
1250 gsMatrix<index_t> & offset) const
1251{
1252 point low, upp, cur, mstr, str, ll, uu, cc;
1253 index_t ii;
1254 // todo: store active_cwise/stride_cwise
1255 const int maxLevel = this->m_tree.getMaxInsLevel();
1256 gsMatrix<T> currPoint;
1257 offset.setZero(maxLevel+1, u.cols()); //trim at the end
1258 // stores [rlvl]->(act,level)
1259 std::vector<std::vector<std::pair<index_t,index_t> > > tfunction(maxLevel+1);
1260
1261 for (index_t i = 0; i != u.cols(); i++) // for all points
1262 {
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();
1266 if (m_manualLevels)
1267 this->_knotIndexToDiadicIndex(maxLevel,low);
1268 // Identify the level of the point
1269 index_t lvl = std::min(this->m_tree.levelOf(low, maxLevel),(int) m_xmatrix.size()-1);
1270 //offset.col(i).setZero();
1271
1272 typename CMatrix::const_iterator it;
1273 for(int level = 0; level <= maxLevel; level++) // for all relevant levels
1274 {
1275 if (level>lvl && tfunction[level].empty()) //nothing to do here
1276 continue; // 0 offset
1277
1278 m_bases[level]->active_cwise(currPoint, low, upp);//my be improved: start from finest lvl
1279 m_bases[level]->stride_cwise(mstr);
1280
1281 if ( level<= lvl)
1282 {
1283 it = m_xmatrix[level].begin();
1284 cur = low;
1285 while ( findNextMatch(it, m_xmatrix[level].end(), cur, low, upp, mstr) )
1286 {
1287 const index_t act = this->m_xmatrix_offset[level] + (it - m_xmatrix[level].begin());//flattensortohier ?
1288 if ( isTruncated(act) )
1289 {
1290 //Record truncated functions with their representation level
1291 const int rlvl = m_is_truncated.at(act);
1292 tfunction[rlvl].push_back( std::make_pair(act,level) );
1293 }
1294 else
1295 ++offset(level, i);
1296 ++it;//advance both
1297 nextCubePoint(cur,low,upp);
1298 }
1299 }
1300
1301 // count truncated active functions
1302 if ( !tfunction[level].empty() )
1303 {
1304 for ( std::pair<index_t,index_t> & q : tfunction[level] )
1305 {
1306 const gsSparseVector<T>& coefs = getCoefs(q.first);
1307 cur = low;
1308 ii = 0;
1309 if ( findNextMatch(coefs, ii, cur, low, upp, mstr) )
1310 ++offset(q.second, i);
1311 }
1312 tfunction[level].clear();
1313 }
1314 }// end level
1315 }//end points
1316
1317 // ii =maxLevel;
1318 // while (offset[ii]==) --ii;
1319 // offset.conservativeresize(ii,gsEigen::NoChange);
1320
1321 return offset.sum();
1322}
1323
1324// returns all actives at \a u from level \a lvl only
1325// Does not clear result, only appends data
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
1329{
1330 gsMatrix<index_t> ind;
1331 GISMO_ASSERT(1==u.cols(), "Expecting single point");
1332 point low, upp, cur, mstr, str, ll, uu, cc;
1333 m_bases[lvl]->active_cwise(u, low, upp);//my be improved: start from finest lvl
1334 index_t ii;
1335 this->m_bases[lvl]->stride_cwise(mstr);
1336 cur = low;
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) )
1340 {
1341 const index_t act = this->m_xmatrix_offset[lvl] + (it - m_xmatrix[lvl].begin());
1342 if ( isTruncated(act) )
1343 {
1344 const gsSparseVector<T>& coefs = getCoefs(act);
1345 const gsTensorBSplineBasis<d, T>& base =
1346 *this->m_bases[this->m_is_truncated[act]];
1347
1348 base.active_cwise(u, ll, uu);
1349 base.stride_cwise(str);
1350 cc = ll;
1351 ii = 0;
1352 if ( findNextMatch(coefs, ii, cc, ll, uu, str) )
1353 result.push_back(act);
1354
1355 //*/
1356
1357 /* // known to work
1358 base.active_into(u, ind);
1359 for (index_t k = 0; k < ind.rows(); ++k)
1360 if (coefs(ind.at(k)) != 0)
1361 {
1362 gsInfo<< "1.GOT value indeed "<< coefs.data().value(ii)<< " at " << ind.at(k) <<".\n";
1363 if (!rr) gsInfo << "Missed "<< ind.at(k) <<"\n";
1364 result.push_back(act);
1365 break;
1366 }
1367 if (rr)
1368 {
1369 gsInfo << "False pos "<< act <<"\n";
1370 //*/
1371
1372 /* // TESTING
1373 gsInfo<< "---- Comparing \n" << base.active(u).transpose()
1374 <<"\n";
1375 typename gsSparseVector<T>::InnerIterator it(coefs);
1376 while (it) { gsInfo<< " " << it.index(); ++it;}
1377 gsInfo<<"\n";
1378 typename gsSparseVector<T>::InnerIterator it2(coefs);
1379 while (it2) { gsInfo<< " " << it2.value(); ++it2;}
1380 gsInfo<<"\n";
1381 base.active_cwise(u, ll, uu);
1382 base.stride_cwise(str);
1383 cc = ll;
1384 ii = 0;
1385
1386 while ( findNextMatch(coefs, ii, cc, ll, uu, str) )
1387 {
1388 result.push_back(act);
1389
1390 //coefs.data().index(ii);
1391 //coefs.data().value(ii);
1392 gsInfo<< "got value "<< coefs.data().value(ii)<< " at " << cc.dot(str) <<" = "<< coefs.data().index(ii)<<".\n";
1393 //result.push_back(act);
1394
1395 break;
1396 ii++;
1397 //nextCubePoint(cc,ll,uu);
1398 }
1399
1400 for (index_t k = 0; k < ind.rows(); ++k)
1401 if (coefs(ind.at(k)) != 0)
1402 {
1403 gsInfo<< "2.GOT value indeed "<< coefs.data().value(ii)<< " at " << ind.at(k) <<".\n";
1404 break;
1405 }
1406 }
1407 //*/
1408 }
1409 else
1410 result.push_back(act);
1411 ++it;//advance both
1412 nextCubePoint(cur,low,upp);
1413 }
1414}
1415
1416template<short_t d, class T>
1418{
1419 gsMatrix<T> currPoint; // HIGHLY INEFFICIENT
1420 point low;
1421 const int maxLevel = this->m_tree.getMaxInsLevel();
1422
1423// gsMatrix<index_t> offset(lvl,u.cols());
1424// std::vector<std::vector<index_t> > offset_tmp(u.cols());
1425
1426 std::vector<std::vector<index_t> > temp_output;//collects the outputs
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;
1431
1432 for(index_t p = 0; p < u.cols(); p++) //for all input points
1433 {
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();
1438
1439 if (m_manualLevels)
1440 this->_knotIndexToDiadicIndex(maxLevel,low);
1441
1442 // Identify the level of the point
1443 const int lvl = std::min(this->m_tree.levelOf(low, maxLevel),(int) m_xmatrix.size()-1);
1444
1445 for(int i = 0; i <= lvl; i++)
1446 activeAtLevel_into(i,currPoint,temp_output[p]);
1447
1448 // find maximum overload sz
1449 if ( temp_output[p].size() > sz )
1450 sz = temp_output[p].size();
1451 }
1452
1453 result.resize(sz, u.cols() );
1454 for(index_t i = 0; i < u.cols(); i++)
1455 {
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();
1458 }
1459}
1460
1461template<short_t d, class T>
1463 const gsMatrix<T>& u,
1464 gsMatrix<T>& result) const
1465{
1466 if (this->m_is_truncated[i] == -1) // basis function not truncated
1467 {
1468 unsigned level = this->levelOf(i);
1469 unsigned tensor_index = flatTensorIndexOf(i, level);
1470 this->m_bases[level]->evalSingle_into(tensor_index, u, result);
1471 }
1472 else
1473 {
1474
1475 unsigned level = this->m_is_truncated[i];
1476 const gsSparseVector<T>& coefs = getCoefs(i);
1477 const gsTensorBSplineBasis<d, T>& base =
1478 *this->m_bases[level];
1479
1480 gsTensorDeboor<d, T, gsKnotVector<T>, gsSparseVector<T> >
1481 (u, base, coefs, result);
1482 }
1483}
1484
1485template<short_t d, class T>
1487 const gsMatrix<T>& u,
1488 gsMatrix<T>& result) const
1489{
1490
1491 if (this->m_is_truncated[i] == -1) // basis function not truncated
1492 {
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);
1496 }
1497 else
1498 {
1499 const unsigned level = this->m_is_truncated[i];
1500 const gsSparseVector<T>& coefs = this->getCoefs(i);
1501 const gsTensorBSplineBasis<d, T> & base =
1502 *this->m_bases[level];
1503
1504 gsTensorDeriv2_into<d, T, gsKnotVector<T>,
1505 gsSparseVector<T> >(u, base, coefs, result);
1506 }
1507}
1508
1509template<short_t d, class T>
1511{
1512 /*
1513 // slightly slower currently (!sameElement)
1514 std::vector<gsMatrix<T> > tmp(1);
1515 evalAllDers_into(u, 0, tmp);
1516 tmp[0].swap(result);
1517 return;
1518 //*/
1519 gsMatrix<index_t> indices;
1520 gsMatrix<T> res(1, 1);
1521 this->active_into(u, indices);//note: evaluating too many separate points together is inefficient
1522
1523 result.setZero(indices.rows(), u.cols());
1524
1525 for (index_t i = 0; i < indices.cols(); i++)
1526 {
1527 for (index_t j = 0; j < indices.rows(); j++)
1528 {
1529 const index_t index = indices(j, i);
1530 if (j != 0 && index == 0)
1531 break;
1532
1533 evalSingle_into(index, u.col(i), res);
1534 result(j, i) = res.value();
1535 }
1536 }
1537}
1538
1539
1540template<short_t d, class T>
1542{
1543 gsMatrix<index_t> indices;
1544 this->active_into(u, indices);
1545
1546 static const short_t numDers = (d * (d + 1)) / 2;
1547 gsMatrix<T> res(numDers, 1); // result of deriv2Single_into
1548
1549 result.setZero(indices.rows() * numDers, u.cols());
1550
1551 for (int i = 0; i < indices.cols(); i++)
1552 {
1553 for (int j = 0; j < indices.rows(); j++)
1554 {
1555 const index_t index = indices(j, i);
1556 if (j != 0 && index == 0)
1557 break;
1558
1559 this->deriv2Single_into(index, u.col(i), res);
1560
1561 result.template block<numDers, 1>(j * numDers, i) = res;
1562 }
1563 }
1564}
1565
1566
1567template<short_t d, class T>
1569{
1570 gsMatrix<index_t> indices;
1571 this->active_into(u, indices);
1572 gsMatrix<T> res(d, 1);
1573
1574 result.setZero(indices.rows() * d, u.cols());
1575
1576 for (int i = 0; i < indices.cols(); i++)
1577 {
1578 for (int j = 0; j < indices.rows(); j++)
1579 {
1580
1581 const unsigned index = indices(j, i);
1582 if (j != 0 && index == 0)
1583 break;
1584
1585 this->derivSingle_into(index, u.col(i), res);
1586 result.template block<d,1>(j * d, i) = res;
1587 }
1588 }
1589}
1590
1591
1592template<short_t d, class T>
1594 const gsMatrix<T> & u,
1595 gsMatrix<T>& result) const
1596{
1597
1598 if (this->m_is_truncated[i] == -1) // basis function not truncated
1599 {
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);
1603 }
1604 else
1605 {
1606 unsigned level = this->m_is_truncated[i];
1607 const gsSparseVector<T>& coefs = this->getCoefs(i);
1608 const gsTensorBSplineBasis<d,T>& base =
1609 *this->m_bases[level];
1610 gsTensorDeriv_into<d, T, gsKnotVector<T>,
1611 gsSparseVector<T> >(u, base, coefs, result);
1612 }
1613}
1614
1615
1616template<short_t d, class T>
1617inline void eval_tp(const std::vector<std::vector<gsMatrix<T> > > & cw,
1618 const gsVector<index_t,d> & ti,
1619 int n, const gsVector<index_t> & str,
1620 index_t m, index_t i,
1621 std::vector<gsMatrix<T> > & result)
1622{
1623 if (n>-1)
1624 {
1625 result[0](m,i) = cw[0][0].at(ti.at(0));//[dir][value].at(func)
1626 for ( short_t k=1; k<d; ++k)
1627 result[0](m,i) *= cw[k][0].at(ti.at(k));
1628 }
1629
1630 if (n>0)
1631 {
1632 T * acc = result[1].col(i).data() + str[1]*m;
1633 for ( short_t k=0; k<d; ++k)
1634 {
1635 // derivative w.r.t. k-th variable
1636 *acc = cw[k][1].at(ti.at(k)); //cw[k][1].at
1637 for ( short_t l=0; l<k; ++l)
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));
1641 ++acc;
1642 }
1643 }
1644
1645 if (n>1)
1646 {
1647 short_t mc = ti.size();
1648 T * acc = result[2].col(i).data() + str[2]*m;
1649 for ( short_t k=0; k<d; ++k)
1650 {
1651 *(acc+k) = cw[k][2].at(ti.at(k)); // pure 2nd derivative w.r.t. k-th variable
1652 for ( short_t l=0; l<k; ++l)
1653 *(acc+k) *= cw[l][0].at(ti.at(l));
1654 for ( short_t l=k+1; l<d; ++l)
1655 {
1656 *(acc+k) *= cw[l][0].at(ti.at(l));
1657 // Then all mixed derivatives follow in lex order
1658 *(acc+mc) = cw[k][1].at(ti.at(k)) * cw[l][1].at(ti.at(l));
1659 for ( short_t q=0; q<k; ++q)
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));
1665 ++mc;
1666 }
1667 }
1668 //++
1669 }
1670}
1671
1672template<short_t d, class T>
1674 std::vector<gsMatrix<T> >& result,
1675 bool sameElement) const
1676{
1677 //gsBasis<T>::evalAllDers_into(u,n,result); return;
1678 result.resize(n+1);
1679 if (0==u.cols()) return;
1680
1681 const int maxLevel = this->m_tree.getMaxInsLevel();
1682 // stores [rlvl]->[act,m]
1683 //std::vector<index_t,std::list<std::pair<index_t,index_t> > > tfunction;
1684 std::vector<std::vector<std::pair<index_t,index_t> > > tfunction(maxLevel+1);
1685
1686 gsMatrix<T> tmp(1,1);
1687 std::vector<std::vector<gsMatrix<T> > > cw(d);//, cwt(d);
1688 std::vector<gsMatrix<T> > temp(n+1), cwa(d);
1689 point low(d), ki(d), kil, til, tlow, tupp, tstr, tcur;
1690
1691 gsMatrix<index_t> tact, act;
1692 if (sameElement)
1693 active_into(u.col(0), act);
1694 else
1695 active_into(u, act);
1696
1697// gsMatrix<index_t> astr;
1698// index_t result_size(sameElement ? numActiveMax(u.col(0),astr) : numActiveMax(u,astr));
1699 index_t result_size = act.rows();//to try: conservativeResize on the go?
1700
1701 // yes: BETTER COMPUTE ACTIVES (maybe also m_bases->actCWise) HERE ONCE AND FOR ALL
1702 // together with active offsets and maybe also tfunction
1703
1704 gsVector<index_t> str(n+1);
1705 for (int l = 0; l <= n; l++)
1706 {
1707 str[l] = numCompositions(l,d);
1708 result[l].setZero(result_size * str[l], u.cols());//zeros in case of sameElement=false
1709 temp[l].resize(str[l], 1);
1710 }
1711
1712 index_t m, index(0);
1713
1714// /* // NEW IMPL
1715
1716 std::vector<std::vector<index_t> >thbact;
1717
1718 int lvl(0);
1719 for (index_t i = 0; i != u.cols(); i++) // for all points
1720 {
1721 if (!sameElement || 0==i)
1722 {
1723 // Identify the level of the point
1724 for(short_t k = 0; k != d; ++k)
1725 {
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();//unique index
1729 }
1730 if (m_manualLevels)
1731 this->_knotIndexToDiadicIndex(maxLevel,low);
1732 lvl = std::min(this->m_tree.levelOf(low, maxLevel),(int) m_xmatrix.size()-1);
1733 }
1734
1735 m = 0; // (!) problem with !sameElement -> numbering in result: counting relies on thbact.size()
1736 thbact.resize(lvl+1);
1737
1738 for(int level = 0; level <= maxLevel; level++) // for all relevant levels(..how about act_offset)
1739 {
1740 if ( level<=lvl && (!sameElement || 0==i) )
1741 {
1742 thbact[level].clear();
1743 activeAtLevel_into(level,u.col(i),thbact[level]);
1744 }
1745
1746 if ( (level>lvl || thbact[level].empty()) && tfunction[level].empty()) //nothing to do here
1747 continue;
1748
1749 for(short_t k = 0; k!=d; ++k)
1750 {
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) ); // get it differently?
1754 }
1755
1756 if (level<=lvl)
1757 for (size_t j = 0; j!=thbact[level].size(); ++j, ++m) // for all actives
1758 {
1759 index = thbact[level][j];
1760 if ( !isTruncated(index) ) // basis function not truncated
1761 {
1762 til = this->m_bases[level]->tensorIndex(flatTensorIndexOf(index, level));
1763 til -= kil; //numbering in element-local position
1764 eval_tp(cw,til,n,str,m,i,result);
1765 }
1766 else if (!sameElement || 0==i)// function is truncated
1767 {
1768 //during active search we also know tfunction.
1769 //If remembered we can evan evaluate here instead.
1770
1771 //Record truncated functions with their representation level
1772 const int rlvl = m_is_truncated.at(index);
1773 tfunction[rlvl].push_back( std::make_pair(index,m) );
1774 }
1775 }
1776
1777 // Evaluate any truncated functions
1778 if ( !tfunction[level].empty() )
1779 {
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);//todo: active cwise - we already jave kil (=lower)
1782 index_t ii;
1783 //-------
1784 for ( std::pair<index_t,index_t> & q : tfunction[level] )
1785 {
1786 const gsSparseVector<T>& coefs = getCoefs(q.first);
1787
1788// /*
1789 ii = 0;
1790 tcur = kil;
1791 while ( findNextMatch(coefs, ii, tcur, kil, til, tstr) )
1792 {
1793 low = tcur-kil;
1794 eval_tp(cw,low,n,str,0,0,temp);
1795 for ( index_t l = 0; l <= n; ++l)
1796 {
1797 auto acc = result[l].col(i).segment(q.second*str[l],str[l]);
1798 acc += coefs.data().value(ii) * temp[l];
1799 }
1800
1801 ++ii;//advance both
1802 nextCubePoint(tcur,kil,til);
1803 }
1804
1805 //*/
1806
1807 /* // known to work
1808 this->m_bases[level]->active_into(u.col(i), tact);// expensive
1809 for ( index_t k = 0; k < tact.rows(); ++k)
1810 {
1811 if ( 0!=coefs[tact.at(k)] )
1812 {
1813 til = this->m_bases[level]->tensorIndex(tact.at(k)) - kil;//reused numbering in element-local position
1814 eval_tp(cw,til,n,str,0,0,temp);
1815 for ( index_t l = 0; l <= n; ++l)
1816 {
1817 auto acc = result[l].col(i).segment(q.second*str[l],str[l]);
1818 acc += coefs[tact.at(k)] * temp[l];
1819 }
1820 }
1821 }
1822 //*/
1823 }
1824 if(!sameElement)
1825 tfunction[level].clear();
1826 }
1827
1828 }//end level
1829 }//end point i
1830
1831 //------------------------------------------
1832 //*/
1833 // known to work
1834 /*
1835 std::vector<gsMatrix<T> > result0;
1836 result0 = result;
1837 for (int l = 0; l <= n; l++)
1838 result[l].setZero();
1839
1840 index_t level;
1841 for (index_t j = 0; j < act.rows(); j++) // for all active functions
1842 {
1843 if (sameElement)
1844 {
1845 index = act.at(j);
1846 level = this->levelOf(index);
1847 til = this->m_bases[level]->tensorIndex(flatTensorIndexOf(index, level));
1848 }
1849 for (index_t i = 0; i != u.cols(); i++) // for all points
1850 {
1851 if (!sameElement)
1852 {
1853 index = act(j, i);
1854 if (j != 0 && index == 0) break;
1855 level = this->levelOf(index);//
1856 til = this->m_bases[level]->tensorIndex(flatTensorIndexOf(index, level));//
1857 }
1858
1859 if ( !isTruncated(index) ) // basis function not truncated
1860 {
1861 for(short_t k = 0; k!=d; ++k)
1862 {
1863 *tmp.data() = u(k,i);
1864 this->m_bases[level]->component(k).evalAllDersSingle_into(til[k], tmp, n, cwa[k]);
1865 }
1866 // -- Start eval
1867 if (n>-1)
1868 {
1869 result[0](j,i) = cwa[0].at(0);
1870 for ( short_t k=1; k<d; ++k)
1871 result[0](j,i) *= cwa[k].at(0);
1872 }
1873
1874 if (n>0)
1875 {
1876 T * acc = result[1].col(i).data() + str[1]*j;
1877 for ( short_t k=0; k<d; ++k)
1878 {
1879 // derivative w.r.t. k-th variable
1880 *acc = cwa[k].at(1);
1881 for ( short_t l=0; l<k; ++l)
1882 *acc *= cwa[l].at(0);
1883 for ( short_t l=k+1; l<d; ++l)
1884 *acc *= cwa[l].at(0);
1885 ++acc;
1886 }
1887 }
1888
1889 if (n>1)
1890 {
1891 short_t m = d;
1892 T * acc = result[2].col(i).data() + str[2]*j;
1893 for ( short_t k=0; k<d; ++k)
1894 {
1895 *(acc+k) = cwa[k].at(2); // pure 2nd derivative w.r.t. k-th variable
1896 for ( short_t l=0; l<k; ++l)
1897 *(acc+k) *= cwa[l].at(0);
1898 for ( short_t l=k+1; l<d; ++l)
1899 {
1900 *(acc+k) *= cwa[l].at(0);
1901 // Then all mixed derivatives follow in lex order
1902 *(acc+m) = cwa[k].at(1) * cwa[l].at(1);
1903 for ( short_t q=0; q<k; ++q)
1904 *(acc+m) *= cwa[q].at(0);
1905 for ( short_t q=k+1; q<l; ++q)
1906 *(acc+m) *= cwa[q].at(0);
1907 for ( short_t q=l+1; q<d; ++q)
1908 *(acc+m) *= cwa[q].at(0);
1909 ++m;
1910 }
1911 }
1912 //TODO higher n here
1913 }
1914 // -- End eval
1915 }
1916 else// function is truncated
1917 {
1918 const int rlvl = m_is_truncated.at(index);
1919 const gsSparseVector<T>& coefs = getCoefs(index);
1920 this->m_bases[rlvl]->active_into(u.col(i), tact);
1921 this->m_bases[rlvl]->evalAllDers_into(u.col(i), n, temp, sameElement);
1922 for (int l = 0; l <= n; l++)
1923 {
1924 auto acc = result[l].col(i).segment(j*str[l],str[l]);
1925 for ( index_t q = 0; q < tact.rows(); ++q) // for all non-zero basis functions
1926 acc += coefs[tact.at(q)] * temp[l].col(0).segment(str[l]*q, str[l]);
1927 }
1928 }
1929 }
1930 }
1931
1932 for (int l = 0; l <= n; l++)
1933 {
1934 if ( (result[l] - result0[l]).norm() > 1e-5)
1935 {
1936 gsInfo <<"l ======== "<<l<<"\n";
1937 gsInfo <<"Old:\n"<< result[l] <<"\n";
1938 gsInfo <<"New:\n"<< result0[l] <<"\n";
1939 gsInfo <<"diff:\n"<< (result[l] - result0[l]) <<"\n";
1940 for (m = 0; m!= act.rows(); ++m)
1941 gsInfo<< (isTruncated(act.at(m)) ? "yes : " : "" )
1942 << m_is_truncated.at(index) <<"\n";
1943 GISMO_ERROR("not Ok");
1944 }
1945 }
1946 //*/
1947
1948}
1949
1950template<short_t d, class T>
1952 typename gsTHBSplineBasis<d, T>::AxisAlignedBoundingBox& boundaryAABB,
1953 typename gsTHBSplineBasis<d, T>::TrimmingCurves& trimCurves) const
1954{
1955 Polylines polylines;
1956 AxisAlignedBoundingBox aabb;
1957
1958 aabb = this->domainBoundariesParams(polylines);
1959 breakCycles(aabb, polylines);
1960
1961 int numBoundaryBoxes = 0;
1962 for (unsigned level = 0; level != aabb.size(); level++)
1963 {
1964 boundaryAABB.push_back(std::vector< std::vector<index_t> >());
1965 trimCurves.push_back(std::vector< std::vector< std::vector< std::vector<T> > > >());
1966
1967 //compare every aabb with the others
1968 for (unsigned boxI = 0; boxI != aabb[level].size(); boxI++)
1969 {
1970 bool isBoundaryBox = true;
1971 for (unsigned boxJ = 0; boxJ != aabb[level].size(); boxJ++)
1972 {
1973 if (boxI != boxJ)
1974 {
1975 if (isFirstBoxCompletelyInsideSecond(aabb[level][boxI], aabb[level][boxJ]))
1976 {
1977 isBoundaryBox = !isBoundaryBox;
1978 }
1979 }
1980 }
1981
1982 if (isBoundaryBox)
1983 {
1984 numBoundaryBoxes++;
1985 boundaryAABB[level].push_back(aabb[level][boxI]);
1986
1987 // make new componenet
1988 trimCurves[level].push_back(std::vector< std::vector< std::vector<T> > >());
1989 trimCurves[level][trimCurves[level].size() - 1].push_back(polylines[level][boxI]);
1990 }
1991 }
1992 }
1993
1994 for (unsigned level = 0; level != aabb.size(); level++)
1995 {
1996 for (unsigned box = 0; box != aabb[level].size(); box++)
1997 {
1998 int closestBox = -1;
1999 for (unsigned boundBox = 0;
2000 boundBox != boundaryAABB[level].size();
2001 boundBox++)
2002 {
2003 if (isFirstBoxCompletelyInsideSecond(aabb[level][box], boundaryAABB[level][boundBox]))
2004 {
2005 if (closestBox == -1 ||
2006 !isFirstBoxCompletelyInsideSecond(boundaryAABB[level][closestBox],
2007 boundaryAABB[level][boundBox]))
2008 {
2009 closestBox = boundBox;
2010 }
2011 }
2012 else if (areBoxesTheSame(aabb[level][box], boundaryAABB[level][boundBox]))
2013 {
2014 closestBox = -1;
2015 break;
2016 }
2017 }
2018
2019 if (-1 < closestBox)
2020 {
2021 trimCurves[level][closestBox].push_back(polylines[level][box]);
2022 }
2023 }
2024 }
2025}
2026
2027
2028template<short_t d, class T>
2029template<short_t dd>
2030typename util::enable_if<dd==2,gsTensorBSpline<d,T> >::type
2031gsTHBSplineBasis<d,T>::getBSplinePatch_impl(const std::vector<index_t>& boundingBox,
2032 const unsigned level,
2033 const gsMatrix<T>& geomCoefs) const
2034{
2035 gsVector<index_t, d> low, upp;
2036 for (unsigned dim = 0; dim != d; dim++)
2037 {
2038 low(dim) = boundingBox[dim];
2039 upp(dim) = boundingBox[d + dim];
2040 }
2041 this->m_tree.computeLevelIndex(low, level, low);
2042 this->m_tree.computeLevelIndex(upp, level, upp);
2043
2044 const gsKnotVector<T> & knots0 = m_bases[level]->knots(0);
2045 const gsKnotVector<T> & knots1 = m_bases[level]->knots(1);
2046
2047 if (m_manualLevels)
2048 {
2049 this->_diadicIndexToKnotIndex(level,low);
2050 this->_diadicIndexToKnotIndex(level,low);
2051 }
2052
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;
2057
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);
2063
2064 gsMatrix<T> coefs;
2065 globalRefinement(geomCoefs, level, coefs);
2066 index_t cc = 0;
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);
2070
2071 std::vector<gsKnotVector<T> > kv(2);
2072
2073 kv[0] = gsKnotVector<T>(m_deg[0], knots0.begin() + lowIndex0,
2074 knots0.begin() + uppIndex0 + m_deg[0] + 2);
2075
2076 kv[1] = gsKnotVector<T>(m_deg[1], knots1.begin() + lowIndex1,
2077 knots1.begin() + uppIndex1 + m_deg[1] + 2);
2078
2079 tensorBasis basis(kv);
2080
2081 return gsTensorBSpline<d, T> (basis, newCoefs);
2082}
2083
2084template<short_t d, class T>
2087 unsigned level,
2088 const gsMatrix<T>& geom_coef,
2090 gsKnotVector<T>& k2) const
2091{ getBsplinePatchGlobal_impl<d>(b1,b2,level,geom_coef,cp,k1,k2); };
2092
2093template<short_t d, class T>
2095 const unsigned level,
2096 const gsMatrix<T>& geomCoefs) const
2097{ return getBSplinePatch_impl<d>(boundingBox, level, geomCoefs); }
2098
2099
2100template<short_t d, class T>
2102 typename gsTHBSplineBasis<d, T>::AxisAlignedBoundingBox& aabb,
2103 typename gsTHBSplineBasis<d, T>::Polylines& polylines) const
2104{
2105 for (size_t level = 0; level != polylines.size(); level++)
2106 {
2107 for (size_t line = 0; line != polylines[level].size(); line++)
2108 {
2109 std::pair<T, T> pt; // point
2110 index_t segment = identifyCycle(polylines[level][line], pt);
2111
2112 if (-1 < segment)
2113 {
2114 std::vector< std::vector<T> > part1, part2;
2115 breakPolylineIntoTwoParts(polylines[level][line], segment, pt,
2116 part1, part2);
2117
2118 polylines[level][line] = part1;
2119 polylines[level].push_back(part2);
2120
2121 std::vector<index_t> aabb1, aabb2;
2122 findNewAABB(part1, aabb1);
2123 findNewAABB(part2, aabb2);
2124
2125 aabb[level][line] = aabb1;
2126 aabb[level].push_back(aabb2);
2127
2128 // very important, this will check current line again if it has more cycles
2129 line--;
2130 }
2131 }
2132 }
2133}
2134
2135// ''''''''''''''''''''''''''''''''''''''''''''''''''''''''''''''''''''''
2136// utility funcitions for breakCycles
2137// ......................................................................
2138
2139// utility funcition for breakCycles
2140template<short_t d, class T>
2141index_t gsTHBSplineBasis<d, T>::identifyCycle(const std::vector< std::vector< T> >& line,
2142 std::pair<T, T>& pt) const
2143{
2144 std::map< std::pair<T, T>, index_t > times;
2145 std::map< std::pair<T, T>, index_t > index;
2146
2147 for (size_t seg = 0; seg != line.size(); seg++)
2148 {
2149 const size_t seg1 = (seg + 1) % line.size();
2150
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])))
2154 {
2155 currentPt.first = line[seg][2];
2156 currentPt.second = line[seg][3];
2157 }
2158
2159 size_t count = times.count(currentPt);
2160 if (0 < count)
2161 {
2162 times[currentPt] += 1;
2163 }
2164 else
2165 {
2166 times[currentPt] = 1;
2167 index[currentPt] = seg1;
2168 }
2169 }
2170
2171 typedef typename std::map< std::pair<T, T>, index_t>::iterator iterator;
2172 for (iterator it = times.begin(); it != times.end(); it++)
2173 {
2174 if (it->second == 2)
2175 {
2176 pt = it->first;
2177 return index[it->first];
2178 }
2179 else
2180 {
2181 GISMO_ENSURE(2 > it->second, "Internal error. Check the polylines from the domainBoundariesParam." );
2182 }
2183 }
2184 return -1;
2185}
2186
2187// utility funcition for breakCycles
2188template<short_t d, class T>
2190 const std::vector< std::vector< T> >& line,
2191 const index_t segment,
2192 const std::pair<T, T>& meetingPt,
2193 std::vector< std::vector<T> >& part1,
2194 std::vector< std::vector<T> >& part2) const
2195{
2196 bool p1 = false; // inside part 1
2197 bool p2 = false; // inside part 2
2198
2199 index_t length = static_cast<index_t> (line.size());
2200 for (index_t i = 0; i != length; i++)
2201 {
2202 const index_t seg = (i + segment) % length;
2203
2204 if (!p1 && !p2) // start
2205 {
2206 p1 = true;
2207 part1.push_back(line[seg]);
2208 }
2209 else // not start
2210 {
2211 // we hit the meeting point again
2212 if ((meetingPt.first == line[seg][0] && meetingPt.second == line[seg][1]) ||
2213 (meetingPt.first == line[seg][2] && meetingPt.second == line[seg][3]))
2214 {
2215 if (p1) // end of part 1
2216 {
2217 part1.push_back(line[seg]);
2218 p1 = false;
2219 p2 = true;
2220 }
2221 else if (p2) // start or finish
2222 {
2223 part2.push_back(line[seg]);
2224 }
2225 }
2226 else
2227 {
2228 if (p1)
2229 {
2230 part1.push_back(line[seg]);
2231 }
2232 else if (p2)
2233 {
2234 part2.push_back(line[seg]);
2235 }
2236 }
2237 }
2238 }
2239}
2240
2241// utility funcition for breakCycles
2242template<short_t d, class T>
2243void gsTHBSplineBasis<d, T>::findNewAABB(const std::vector< std::vector<T> >& polyline,
2244 std::vector<index_t>& aabb) const
2245{
2246 T minX = polyline[0][0];
2247 T minY = polyline[0][1];
2248 T maxX = polyline[0][2];
2249 T maxY = polyline[0][3];
2250
2251
2252 for (size_t seg = 0; seg != polyline.size(); seg++)
2253 {
2254 if (polyline[seg][0] < minX)
2255 {
2256 minX = polyline[seg][0];
2257 }
2258 if (polyline[seg][1] < minY)
2259 {
2260 minY = polyline[seg][1];
2261 }
2262 if (maxX < polyline[seg][2])
2263 {
2264 maxX = polyline[seg][2];
2265 }
2266 if (maxY < polyline[seg][3])
2267 {
2268 maxY = polyline[seg][3];
2269 }
2270 }
2271
2272 unsigned maxLevel = this->maxLevel();
2273 const gsKnotVector<T>& kv0 = this->m_bases[maxLevel]->knots(0);
2274 const gsKnotVector<T>& kv1 = this->m_bases[maxLevel]->knots(1);
2275
2276 aabb.resize(4);
2277 for (unsigned i = 0; i != kv0.uSize(); i++)
2278 {
2279 if (kv0(i) <= minX)
2280 {
2281 aabb[0] = i;
2282 }
2283 if (maxX <= kv0(i))
2284 {
2285 aabb[2] = i;
2286 break;
2287 }
2288 }
2289
2290 for (unsigned i = 0; i != kv1.uSize(); i++)
2291 {
2292 if (kv1(i) <= minY)
2293 {
2294 aabb[1] = i;
2295 }
2296 if (maxY <= kv1(i))
2297 {
2298 aabb[3] = i;
2299 break;
2300 }
2301 }
2302}
2303
2304
2305
2306// --------------------------------------------------------------------------------
2307// Code for hierarchical coarsening
2308// --------------------------------------------------------------------------------
2309
2310template<short_t d, class T>
2312{
2313 result.clear();
2314 gsVector<index_t> level;
2315 gsMatrix<index_t> b1, b2;//boxes in highes level numbering
2316 this->m_tree.getBoxesInLevelIndex(b1,b2,level);//return boxes in level indices
2317 tensorBasis T_0_copy = this->tensorLevel(0);
2318 std::vector< gsSparseMatrix<T,RowMajor> > transfer;
2319 transfer.resize(this->maxLevel() );
2320 std::vector<std::vector<T> > knots(d);
2321
2322 for(unsigned i = 0; i < this->maxLevel(); ++i)
2323 {
2324 //T_0_copy.uniformRefine_withTransfer(transfer[i], 1);
2325 for(short_t dim = 0; dim < d; dim++)
2326 {
2327 const gsKnotVector<T> & ckv = m_bases[i]->knots(dim);
2328 const gsKnotVector<T> & fkv = m_bases[i + 1]->knots(dim);
2329 ckv.symDifference(fkv, knots[dim]);
2330
2331 //gsDebug << "level: " << i << "\n"
2332 // << "direction: " << dim << "\n"
2333 // << "dirknots:\n"<<gsAsMatrix<T>(dirKnots)<<std::endl;
2334 }
2335 T_0_copy.refine_withTransfer(transfer[i], knots);
2336
2337 // Must use refine_withTransfer(gsSparseMatrix<T,RowMajor> & transfer, const std::vector<T>& knots)
2338 // must correctly find the knots to insert
2339 // T_0_copy.uniformRefine_withTransfer(transfer[i], 1);
2340 }
2341 std::vector< gsSparseMatrix<T,RowMajor> > temp_transf;
2342 for(unsigned j = 0; j < this->maxLevel(); ++j)
2343 {
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]);
2348
2349 gsSparseMatrix<T> crs = this->coarsening_direct(x_mat_old_0, x_matrix_lvl, temp_transf);
2350 result.push_back(crs);
2351 }
2352}
2353
2354//todo remove
2355template<short_t d, class T>
2357{
2358 int size1= 0, size2 = 0;
2359 int glob_numb = 0;//continous numbering of hierarchical basis
2360 for(unsigned int i =0; i< old.size();i++)
2361 {//count the number of basis functions in old basis
2362 size1 += old[i].size();
2363 }
2364 for(unsigned int i =0; i< n.size();i++)
2365 {//count the number of basis functions in new basis
2366 size2 += n[i].size();
2367 }
2368 gsSparseMatrix<T> result(size2,size1);
2369
2370 gsMatrix<T> transferDense = transfer;
2371
2372 for (unsigned int i = 0; i < old.size(); i++)//iteration through the levels of the old basis
2373 {
2374 // find starting index of level i in new basis
2375 int start_lv_i = 0;
2376 for(unsigned int l =0; l < i; l++)
2377 {
2378 start_lv_i += n[l].size();
2379 }
2380
2381 for (unsigned int j = 0; j < old[i].size();j++)//iteration through the basis functions in the given level
2382 {
2383 start_lv_i = 0;
2384 for(unsigned int l =0; l < i; l++)
2385 {
2386 start_lv_i += n[l].size();
2387 }
2388 const unsigned old_ij = old[i][j]; // tensor product index
2389
2390 if( n[i].bContains(old_ij) )//it he basis function was not refined
2391 {
2392 result(start_lv_i + std::distance(n[i].begin(), n[i].find_it_or_fail(old_ij) ),glob_numb ) = 1;//settign the coefficient of the not refined basis function to 1
2393 for(int k = 0; k < transferDense.rows(); k++)//basis function was refined->looking for the coeficinets from the global transfer matrix
2394 {
2395 //TODO test if the basis function is in the basis
2396 if(transferDense(k, old_ij) != 0)//if the coefficient is non zero we find the coresponding function in n
2397 {
2398 if(n[i+1].bContains(k))
2399 {
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);
2402 }
2403 }
2404 }
2405 }
2406 else
2407 {
2408 for(int k = 0; k < transferDense.rows(); k++)//basis function was refined->looking for the coeficinets from the global transfer matrix
2409 {
2410 if(transferDense(k, old_ij) != 0)//if the coefficient is non zero we find the coresponding function in n
2411 {
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);
2414 }
2415 }
2416 }
2417 glob_numb++;
2418 }
2419 }
2420 return result;
2421}
2422template<short_t d, class T>
2424 const std::vector<gsSortedVector<index_t> >& n,
2425 const std::vector<gsSparseMatrix<T,RowMajor> >& transfer) const
2426{
2427 int size1= 0;int size2 = 0;
2428 int glob_numb = 0;//continous numbering of hierarchical basis
2429 for(unsigned int i =0; i< old.size();i++){//count the number of basis functions in old basis
2430 size1 += old[i].size();
2431 }
2432 for(unsigned int i =0; i< n.size();i++){//count the number of basis functions in new basis
2433 size2 += n[i].size();
2434 }
2435 gsSparseMatrix<T> result(size2,size1);
2436
2437
2438 //std::vector<gsMatrix<T> > transferDense;// = transfer;
2439 //transferDense.resize(transfer.size());
2440 //for (unsigned int i = 0; i < transfer.size();i++){
2441 // transferDense[i] = transfer[i];
2442 //}
2443 std::vector<gsSparseMatrix<T,ColMajor> > temptransfer;// = transfer;
2444 temptransfer.resize(transfer.size());
2445 for (unsigned int i = 0; i < transfer.size();i++){
2446 temptransfer[i] = transfer[i];
2447 //gsDebug<<"transfer"<<i<<"\n"<<transfer[i]<<std::endl;
2448 }
2449
2450 for (unsigned int i = 0; i < old.size(); i++)//iteration through the levels of the old basis
2451 {
2452 // find starting index of level i in new basis
2453 int start_lv_i = 0;
2454 for(unsigned int l =0; l < i; l++)
2455 {
2456 start_lv_i += n[l].size();
2457 }
2458
2459 //gsDebug<<old[i].size()<<std::endl;
2460 for (unsigned int j = 0; j < old[i].size();j++)//iteration through the basis functions in the given level
2461 {
2462 //gsDebug<<"j = "<<j<<std::endl;
2463 start_lv_i = 0;
2464 for(unsigned int l =0; l < i; l++)
2465 {
2466 start_lv_i += n[l].size();
2467 }
2468 const unsigned old_ij = old[i][j]; // tensor product index
2469 gsMatrix<index_t, d, 2> supp(d, 2);
2470 this->m_bases[i]->elementSupport_into(old_ij, supp);//this->support(start_lv_i+old_ij);
2471 //gsDebug<<"supp "<< supp<<std::endl;
2472 //unsigned max_lvl =
2473 // math::min<index_t>( this->m_tree.query4(supp.col(0),supp.col(1), i), transfer.size() ) ;
2474 gsSparseVector<T,RowMajor> t(this->m_bases[i]->size());
2475 t.setZero();
2476 t[old_ij] = 1;
2477 for(unsigned int k = i; k < n.size();k++){
2478 //gsDebug<<"i: "<<i<<" j: "<<j<<" k: "<<k<<std::endl;
2479 if(k > i)
2480 {
2481 //compare with old matrix
2482 for(int l = 0 ; l < t.size();l++){
2483 if(t[l]!=0)
2484 if( (k) < (old.size()))
2485 if(old[k].bContains(l)){
2486 t[l]=0;
2487 }
2488
2489 }
2490 }
2491 //gsDebug<<"ksize:"<<n[k].size()<<std::endl;
2492 if(k!=0)
2493 {
2494 start_lv_i = 0;
2495 for(unsigned int l =0; l < k-1; l++)
2496 {
2497 start_lv_i += n[l].size();
2498 }
2499 }
2500 //for all non zero in t comapre with new
2501 for(int l = 0 ; l < t.size();l++)
2502 {
2503 //gsDebug<<"i: "<<i<<" j: "<<j<<" l: "<<l<<"nsize"<<n.size()<<std::endl;
2504 if(t[l]!=0)
2505 if(n[k].bContains(l))
2506 {
2507 //gsDebug<<"j: "<<j<<" "<<"oldij"<<old_ij<<" "<<"l:"<<l<<" ";
2508 //gsDebug<<"k "<<k<<" j: "<<j<<" "<<"globnumb "<<glob_numb<<" oldij "<<old_ij<<" "<<"l:"<<l<<" "<<t[l]<<" ";
2509 int p = 0;
2510 if(k!=0){
2511 p = n[k-1].size();
2512 }
2513 const int pos = start_lv_i + p + std::distance(n[k].begin(), n[k].find_it_or_fail(l));
2514 //const int pos = start_lv_i + n[k].size() + std::distance(n[k+1].begin(), n[k+1].find_it_or_fail(l));
2515 //gsDebug<<"pos:"<<pos<<std::endl;
2516 //double ppp = transferDense[coeffs[0].lvl](k, coeffs[0].pos);
2517 result(pos,glob_numb) = t[l];//transferDense[coeffs[0].lvl](k, coeffs[0].pos);
2518 //t[l] = 0;
2519 //const int pos = start_lv_i + n[coeffs[0].lvl].size() + std::distance(n[coeffs[0].lvl+1].begin(), n[coeffs[0].lvl+1].find_it_or_fail(k.row()));
2520 //double ppp = transferDense[coeffs[0].lvl](k, coeffs[0].pos);
2521 //result(pos,glob_numb) += coeffs[0].coef * temptransfer[coeffs[0].lvl](k.row(), coeffs[0].pos);//transferDense[coeffs[0].lvl](k, coeffs[0].pos);
2522 }
2523 }
2524 //gsDebug<<"b"<<std::endl;
2525 gsSparseVector<T,RowMajor> temp(this->m_bases[i]->size());
2526 //gsDebug<<"c"<<std::endl;
2527 temp.setZero();
2528 if(k<temptransfer.size())
2529 temp = temptransfer[k] * t.transpose();
2530 t = temp;
2531 //refine
2532 }
2533
2534
2535// if( n[i].bContains(old_ij) )//it he basis function was not refined
2536// {
2537// result(start_lv_i + std::distance(n[i].begin(), n[i].find_it_or_fail(old_ij) ),glob_numb ) = 1;//settign the coefficient of the not refined basis function to 1
2538// }
2539// else
2540// {
2541
2542// gsMatrix<index_t, d, 2> supp(d, 2);
2543// this->m_bases[i]->elementSupport_into(old_ij, supp);//this->support(start_lv_i+old_ij);
2544// //gsDebug<<<"supp "<< supp<<std::endl;
2545// unsigned max_lvl =
2546// math::min<index_t>( this->m_tree.query4(supp.col(0),supp.col(1), i), transfer.size() ) ;
2547// gsSparseVector<T,RowMajor> t(this->m_bases[i]->size());
2548// t.setZero();
2549// t[old_ij] = 1;
2550// //gsDebug<<"j = "<<j<<std::endl;
2551// //gsDebug<<"nsize"<<n.size()<<std::endl;
2552// for(unsigned int k = i+1; k < n.size();k++){
2553// start_lv_i = 0;
2554// for(unsigned int l =0; l < k-1; l++)
2555// {
2556// start_lv_i += n[l].size();
2557// }
2558// // gsDebug<<"k "<<k<<std::endl;
2559// // gsDebug<<"nk"<<std::endl;
2560// // for(int a = 0; a < n[k].size();a++){
2561// // gsDebug<<n[k][a]<<" ";
2562// // }
2563// //gsDebug<<std::endl;
2564// gsSparseVector<T,RowMajor> M;
2565// M.setZero();
2566// M = temptransfer[k-1] * t.transpose();
2567
2568// //gsDebug<<"M\n"<<M<<std::endl;
2569// for(int l = 0 ; l < M.size();l++)
2570// //for(typename gsSparseVector<T,RowMajor>::InnerIterator l(M); l; ++l)//basis function was refined->looking for the coeficinets from the global transfer matrix
2571// {
2572// if(M[l]!=0)
2573// if(n[k].bContains(l))
2574// {
2575// //gsDebug<<"j: "<<j<<" "<<"oldij"<<old_ij<<" "<<"l:"<<l<<" ";
2576// //gsDebug<<"k "<<k<<" j: "<<j<<" "<<"globnumb "<<glob_numb<<" "<<"l:"<<l<<" "<<M[l]<<" ";
2577// const int pos = start_lv_i + n[k-1].size() + std::distance(n[k].begin(), n[k].find_it_or_fail(l));
2578// //gsDebug<<"pos:"<<pos<<std::endl;
2579// //double ppp = transferDense[coeffs[0].lvl](k, coeffs[0].pos);
2580// result(pos,glob_numb) = M[l];//transferDense[coeffs[0].lvl](k, coeffs[0].pos);
2581// M[l] = 0;
2582// //const int pos = start_lv_i + n[coeffs[0].lvl].size() + std::distance(n[coeffs[0].lvl+1].begin(), n[coeffs[0].lvl+1].find_it_or_fail(k.row()));
2583// //double ppp = transferDense[coeffs[0].lvl](k, coeffs[0].pos);
2584// //result(pos,glob_numb) += coeffs[0].coef * temptransfer[coeffs[0].lvl](k.row(), coeffs[0].pos);//transferDense[coeffs[0].lvl](k, coeffs[0].pos);
2585// }
2586// }
2587// t = M;
2588// //gsDebug<<"M after \n"<<M<<std::endl;
2589// }
2590// }
2591 glob_numb++;
2592 }
2593
2594 }
2595
2596 return result;
2597}
2598
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
2603{
2604 GISMO_ASSERT(old.size() < n.size(), "old,n problem in coarsening.");
2605
2606 int size1= 0;int size2 = 0;
2607 int glob_numb = 0;//continous numbering of hierarchical basis
2608 for(unsigned int i =0; i< old.size();i++)
2609 {//count the number of basis functions in old basis
2610 size1 += old[i].size();
2611 }
2612 for(unsigned int i =0; i< n.size();i++)
2613 {//count the number of basis functions in new basis
2614 size2 += n[i].size();
2615 }
2616 gsSparseMatrix<T> result(size2,size1);
2617
2618
2619// std::vector<gsMatrix<T> > transferDense;// = transfer;
2620// transferDense.resize(transfer.size());
2621// for (unsigned int i = 0; i < transfer.size();i++)
2622// {
2623// transferDense[i] = transfer[i];
2624// }
2625 std::vector<gsSparseMatrix<T,ColMajor> > temptransfer;// = transfer;
2626 temptransfer.resize(transfer.size());
2627 for (unsigned int i = 0; i < transfer.size();i++)
2628 {
2629 temptransfer[i] = transfer[i];
2630 }
2631 //gsDebug<<"temp:\n"<< temptransfer[0]<<std::endl;
2632
2633 for (unsigned int i = 0; i < old.size(); i++)//iteration through the levels of the old basis
2634 {
2635 // find starting index of level i in new basis
2636 int start_lv_i = 0;
2637 for(unsigned int l =0; l < i; l++)
2638 {
2639 start_lv_i += n[l].size();
2640 }
2641
2642 for (unsigned int j = 0; j < old[i].size();j++)//iteration through the basis functions in the given level
2643 {
2644 //gsDebug<<"j...."<<j<<endl;
2645 //gsDebug<<"i = "<< i<< " j= "<<j<<std::endl;
2646 start_lv_i = 0;
2647 for(unsigned int l =0; l < i; l++)
2648 {
2649 start_lv_i += n[l].size();
2650 }
2651 const unsigned old_ij = old[i][j]; // tensor product index
2652 //gsDebug<<"oldij = "<< old_ij<<std::endl;
2653 if( n[i].bContains(old_ij) )//it he basis function was not refined
2654 {
2655 result(start_lv_i + std::distance(n[i].begin(), n[i].find_it_or_fail(old_ij) ),glob_numb ) = 1;//settign the coefficient of the not refined basis function to 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);//this->support(start_lv_i+old_ij);
2659 gsVector<index_t, d > low = supp.col(0);
2660 gsVector<index_t, d > upp = supp.col(1);
2661 if (m_manualLevels)
2662 {
2663 this->_knotIndexToDiadicIndex(i,low);
2664 this->_knotIndexToDiadicIndex(i,upp);
2665 }
2666 unsigned max_lvl = math::min<index_t>( this->m_tree.query4(low,upp, i), transfer.size()) ;//transfer.size();//
2667 //gsDebug<<"transfer size "<< transfer.size()<<" max lvl"<< this->m_tree.query4(supp.col(0),supp.col(1), this->levelOf(start_lv_i+j))<<" lvl of"<< this->levelOf(start_lv_i+j)<<" support\n"<<supp<<std::endl;
2668 lvl_coef temp;
2669 temp.pos = old_ij;
2670 temp.coef = 1;
2671 temp.lvl = i;
2672 if(temp.lvl<n.size()-1)//in case of no refinement
2673 {
2674 coeffs.push_back(temp);
2675 }
2676 for (size_t coeff_index = 0; coeff_index < coeffs.size(); ++coeff_index)
2677 {
2678 const lvl_coef coeff = coeffs[coeff_index];
2679
2680 start_lv_i = 0;
2681 for(unsigned int l =0; l < coeff.lvl; l++)
2682 {
2683 start_lv_i += n[l].size();
2684 }
2685
2686 for(typename gsSparseMatrix<T,ColMajor>::InnerIterator k(temptransfer[coeff.lvl],coeff.pos); k; ++k)//basis function was refined->looking for the coeficinets from the global transfer matrix
2687 {
2688// if(transferDense[coeff.lvl](k, coeff.pos) != 0)
2689// {
2690 //gsDebug<<"first while "<< coeffs.size()<<std::endl;
2691 bool p = true;
2692 if( (coeff.lvl+1) < old.size())
2693 {
2694 if(old[coeff.lvl+1].bContains(k.row()))
2695 {
2696 p = false;
2697 }
2698 }
2699 if(n[coeff.lvl+1].bContains(k.row()))
2700 {
2701 if(p)
2702 {
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()));
2704
2705 result(pos,glob_numb) += coeff.coef * temptransfer[coeff.lvl](k.row(), coeff.pos);//transferDense[coeff.lvl](k, coeff.pos);
2706 if(coeff.lvl + 1 < max_lvl)//transfer.size()
2707 {
2708 temp.pos = k.row();
2709 temp.coef = temptransfer[coeff.lvl](k.row(), coeff.pos) * coeff.coef;
2710 temp.lvl = coeff.lvl+1;
2711 coeffs.push_back(temp);
2712 }
2713 }
2714 }else
2715 {
2716 //TODO test if level would not be higher than the q4 for this function- or max inserted level
2717 if(p)
2718 {
2719 if( coeff.lvl + 1< max_lvl)
2720 {
2721 temp.pos = k.row();
2722 temp.coef = temptransfer[coeff.lvl](k.row(), coeff.pos) * coeff.coef;
2723 temp.lvl = coeff.lvl+1;
2724 coeffs.push_back(temp);
2725 }
2726 }
2727 }
2728 //}
2729 }
2730
2731// for(unsigned int ii = 0; ii < coeffs.size();ii++){
2732// gsDebug<<"( "<< coeffs[ii].pos<<" , "<<coeffs[ii].coef<<" , "<<coeffs[ii].lvl<<" )"<<std::endl;
2733// }
2734// gsDebug<<"j ="<<j<<std::endl;
2735 }
2736 }
2737 else
2738 {
2739 gsMatrix<index_t, d, 2> supp(d, 2);
2740 this->m_bases[i]->elementSupport_into(old_ij, supp);//this->support(start_lv_i+old_ij);
2741 //gsDebug<<"supp "<< supp<<std::endl;
2742 gsVector<index_t, d > low = supp.col(0);
2743 gsVector<index_t, d > upp = supp.col(1);
2744 if (m_manualLevels)
2745 {
2746 this->_knotIndexToDiadicIndex(i,low);
2747 this->_knotIndexToDiadicIndex(i,upp);
2748 }
2749
2750 unsigned max_lvl =
2751 math::min<index_t>( this->m_tree.query4(low,upp, i), transfer.size() ) ;
2752 std::vector<lvl_coef> coeffs;
2753 lvl_coef temp;
2754 temp.pos = old_ij;
2755 temp.coef = 1;
2756 temp.lvl = i;
2757 //gsDebug<<"temp.pos: "<<temp.pos<<" temp.coef: "<<temp.coef<<" temp.lvl "<< temp.lvl<<std::endl;
2758 coeffs.push_back(temp);
2759 for (size_t coeff_index = 0; coeff_index < coeffs.size(); ++coeff_index)
2760 {
2761 const lvl_coef coeff = coeffs[coeff_index];
2762
2763 start_lv_i = 0;
2764 for(unsigned int l =0; l < coeff.lvl; l++)
2765 {
2766 start_lv_i += n[l].size();
2767 }
2768
2769 for(typename gsSparseMatrix<T,ColMajor>::InnerIterator k(temptransfer[coeff.lvl],coeff.pos); k; ++k)//basis function was refined->looking for the coeficinets from the global transfer matrix
2770 {
2771 //if(transferDense[coeff.lvl](k, coeff.pos) != 0)
2772 //{
2773 //gsDebug<<"second while "<< coeffs.size()<<std::endl;
2774 bool p = true;
2775 if( (coeff.lvl+1) < (old.size()))
2776 {
2777 if(old[coeff.lvl+1].bContains(k.row()))
2778 {
2779 p = false;
2780 }
2781 }
2782 //gsDebug<<"k.row "<<k.row()<<std::endl;
2783 if(n[coeff.lvl+1].bContains(k.row()))
2784 {
2785 if(p)
2786 {
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()));
2788 //T ppp = transferDense[coeff.lvl](k, coeff.pos);
2789 //gsDebug<<"pos "<<pos<<" oldij "<<old_ij<<" "<< "coeflvl "<<coeff.lvl<<" ";
2790 //gsDebug<<"inserted coef "<< coeff.coef<<"*"<< temptransfer[coeff.lvl](k.row(), coeff.pos)<<std::endl;
2791 result(pos,glob_numb) += coeff.coef * temptransfer[coeff.lvl](k.row(), coeff.pos);//transferDense[coeff.lvl](k, coeff.pos);
2792 if(coeff.lvl < max_lvl-1)
2793 {
2794 temp.pos = k.row();
2795 temp.coef = temptransfer[coeff.lvl](k.row(), coeff.pos) * coeff.coef;
2796 temp.lvl = coeff.lvl+1;
2797 //gsDebug<<"temp.pos: "<<temp.pos<<" temp.coef: "<<temp.coef<<" temp.lvl "<< temp.lvl<<std::endl;
2798 coeffs.push_back(temp);
2799 }
2800 }
2801 }
2802 else
2803 {
2804 if(p)
2805 {
2806 if(coeff.lvl < max_lvl-1)
2807 {
2808 temp.pos = k.row();
2809 temp.coef = temptransfer[coeff.lvl](k.row(), coeff.pos) * coeff.coef;
2810 temp.lvl = coeff.lvl+1;
2811 //gsDebug<<"temp.pos: "<<temp.pos<<" temp.coef: "<<temp.coef<<" temp.lvl "<< temp.lvl<<std::endl;
2812 coeffs.push_back(temp);
2813 }
2814 }
2815 }
2816 //}
2817 }
2818
2819// for(unsigned int ii = 0; ii < coeffs.size();ii++){
2820// gsDebug<<"( "<< coeffs[ii].pos<<" , "<<coeffs[ii].coef<<" , "<<coeffs[ii].lvl<<" )"<<std::endl;
2821// }
2822 }
2823 }
2824 glob_numb++;
2825
2826 }
2827 }
2828 return result;
2829}
2830
2831
2832namespace internal
2833{
2834
2836template<short_t d, class T>
2837class gsXml< gsTHBSplineBasis<d,T> >
2838{
2839private:
2840 gsXml() { }
2841public:
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):""); }
2845
2846 static gsTHBSplineBasis<d,T> * get (gsXmlNode * node)
2847 {
2848 return getHTensorBasisFromXml< gsTHBSplineBasis<d,T> > (node);
2849 }
2850
2851 static gsXmlNode * put (const gsTHBSplineBasis<d,T> & obj,
2852 gsXmlTree & data )
2853 {
2854 return putHTensorBasisToXml< gsTHBSplineBasis<d,T> > (obj, data);
2855 }
2856};
2857
2858
2859}
2860
2861} // namespace gismo
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.
Utility functions related to tensor-structured objects.
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