G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsHTensorBasis.hpp
Go to the documentation of this file.
1
14#pragma once
15
17
18#include <gsUtils/gsPointGrid.h>
19
20#include <gsIO/gsXml.h>
22
23namespace gismo
24{
25 // Central element implementation
26// template<short_t d, class T>
27// gsMatrix<T> gsHTensorBasis<d,T>::elementInSupportOf(index_t i) const
28// {
29// index_t lvl = levelOf(i);
30// index_t j = flatTensorIndexOf(i);
31// return m_bases[lvl]->elementInSupportOf(j);
32// }
33
34template<short_t d, class T>
36{
37 GISMO_ENSURE(m_manualLevels,"Add level only works when m_manualLevels==true");
38 // Fill the unique indices
39 std::vector<std::vector<index_t>> lvlIndices(d);
40 const tensorBasis * tb2 = dynamic_cast<const tensorBasis*>(m_bases.back());
41 std::vector<T> difference, intersection;
42 std::vector<T> knots1, knots2;
43 std::vector<index_t> dirIndices;
44 gsKnotVector<T> tmpknots;
45 for (short_t dim=0; dim!=d; dim++)
46 {
47 dirIndices = m_uIndices.back()[dim];
48 // Take the difference of the knot vectors
49 tmpknots= tb2->knots(dim);
50
51 GISMO_ASSERT(tmpknots.maxInteriorMultiplicity()<tmpknots.degree()+1,"Knot vector cannot have knot multiplicities higher than p="<<tmpknots.degree());
52
53 knots1 = tb2->knots(dim).unique();
54 knots2 = next_basis.knots(dim).unique();
55
56 // Check nestedness.
57 intersection.clear();
58 difference.clear();
59 // The unique knots of the new basis must contain the ones of the previous level
60 std::set_intersection( knots2.begin(), knots2.end(),
61 knots1.begin(), knots1.end(),
62 std::back_inserter(intersection) );
63 // The difference of the two must not contain any knot in knots1!
64 std::set_difference( knots1.begin(), knots1.end(),
65 intersection.begin(), intersection.end(),
66 std::back_inserter(difference) );
67 GISMO_ASSERT(difference.size()==0,"Knot vector is not nested!");
68
69 difference.clear();
70 // The difference of the two must not contain any knot in knots1!
71 std::set_difference(knots2.begin(), knots2.end(),
72 knots1.begin(), knots1.end(),
73 std::back_inserter(difference) );
74
75 // We reverse since later on we loop and pop the elements on the back
76 std::reverse(difference.begin(),difference.end());
77 // Double all indices in dirIndices. These correspond to the indices that come from the nested bases
78 std::transform(dirIndices.begin(), dirIndices.end(), dirIndices.begin(), [=](index_t& i){return 2*i;});
79
80 // Find the extra indices for the nodes in difference
81 while (difference.size()!=0)
82 {
83 index_t i, n;
84 typename std::vector<T>::iterator diff_ptr = difference.begin();
85 // Find the index of the highest knot below the different knot
86 i = tmpknots.uFind(*diff_ptr).uIndex(); // index of the knot span in which *it is located
87
88 // Count how many difference knots lay in the same span
89 for (n = 0, diff_ptr = difference.begin(); diff_ptr!=difference.end() && tmpknots.uFind(*diff_ptr).uIndex()==i; n++, diff_ptr++);
90 // Count the distance between two element indices. It should be larger than n
91 GISMO_ENSURE(n < dirIndices[i+1] - dirIndices[i],"Trying to insert more knots than available in the tensor structure");
92
93 index_t newIdx;
94 if (n % 2 != 0)
95 {
96 // We add the middle knot
97 diff_ptr = std::next(difference.begin(),(n-1)/2);
98 newIdx = (dirIndices[i] + dirIndices[i+1]) / 2;
99 }
100 else
101 {
102 // We add the first knot
103 diff_ptr = difference.begin();
104 newIdx = dirIndices[i] + 1;
105 }
106
107 // Add the correct index
108 dirIndices.insert(std::upper_bound( dirIndices.begin(), dirIndices.end(), newIdx ),newIdx);
109 // Insert the removed difference knot in the temporary knot vector
110 tmpknots.insert(*diff_ptr);
111 // Remove the difference knot, since it is treated
112 difference.erase(diff_ptr);
113 }
114
115 lvlIndices[dim] = dirIndices;
116 }
117 m_uIndices.push_back(lvlIndices);
118
119 // m_bases.push_back( new gsTensorBSplineBasis<d, T>( give( next_basis ) ));
120 m_bases.push_back( next_basis.clone().release() );
121}
122
123template<short_t d, class T>
125{
126 return m_bases[0]->support();
127}
128
129//i in continuous numbering
130template<short_t d, class T>
132{
133 // Get the level
134 int lvl = levelOf(i);
135 // Return the the support
136 return m_bases[lvl]->support( this->flatTensorIndexOf(i) );
137}
138
139// S.K.
140template<short_t d, class T> inline
142{
143 GISMO_ASSERT(Pt.cols() == 1, "Waiting for single point");
144 point loIdx;
145
146 const int maxLevel = m_tree.getMaxInsLevel();
147
148 needLevel(maxLevel);
149
150 for( int i =0; i < Dim; i++)
151 loIdx[i] = m_bases[maxLevel]->knots(i).uFind( Pt(i,0) ).uIndex();
152
153 if (m_manualLevels)
154 this->_knotIndexToDiadicIndex(maxLevel,loIdx);
155
156 return m_tree.levelOf( loIdx, maxLevel);
157}
158
159template<short_t d, class T> inline
161{
162 const int maxLevel = m_tree.getMaxInsLevel();
163
164 needLevel(maxLevel);
165
166 return m_tree.levelOf( Pt, maxLevel);
167}
168
169template<short_t d, class T> inline
171 gsVector<index_t> & lvl,
172 gsMatrix<index_t> & loIdx ) const
173{
174 lvl.resize( Pt.cols() );
175 loIdx.resize( Pt.rows(), Pt.cols() );
176 lvl.setZero();
177 loIdx.setZero();
178 for( index_t i = 0; i < Pt.cols(); i++)
179 {
180 lvl[i] = getLevelAtPoint( Pt.col(i) );
181 for( index_t j = 0; j < Pt.rows(); j++)
182 loIdx(j,i) = m_bases[ lvl[i] ]->knots(j).uFind( Pt(j,i) ).uIndex() ;
183 }
184}
185
186template<short_t d, class T> inline
188{
189 result.resize( u.cols() );
190 result.setZero();
191
192 point low, upp, cur;
193 const int maxLevel = m_tree.getMaxInsLevel();
194
195 for(index_t p = 0; p < u.cols(); p++ ) //for all input points
196 {
197 for(short_t i = 0; i != d; ++i)
198 low[i] = m_bases[maxLevel]->knots(i).uFind(u(i,p)).uIndex();
199
200 // Identify the level of the point
201 const int lvl = m_tree.levelOf(low, maxLevel);
202
203 for(int i = 0; i <= lvl; i++)
204 {
205 m_bases[i]->active_cwise(u.col(p), low, upp);
206 cur = low;
207 do //iterate over all points in [low,upp]
208 {
209 CMatrix::const_iterator it =
210 m_xmatrix[i].find_it_or_fail( m_bases[i]->index(cur) );
211
212 if( it != m_xmatrix[i].end() )// if index is found
213 result[p]++;
214 }
215 while( nextCubePoint(cur,low,upp) );
216 }
217 }
218}
219
220template<short_t d, class T>
221void gsHTensorBasis<d, T>::addConnectivity(int lvl, gsMesh<T> & mesh) const
222{
224
225 const gsTensorBSplineBasis<d, T> & bb = *m_bases[lvl];
226 const CMatrix & cmat = m_xmatrix[lvl];
227
228 // Last tensor-index in level lvl
230 for (index_t i = 0; i < d; ++i)
231 end(i) = bb.component(i).size() - 1;
232
233 index_t k, s;
235 for (index_t i = 0; i < d; ++i) // For all axes
236 {
237 s = bb.stride(i);
238 v = low;
239 upp = end;
240 upp[i] = 0; // suppress to face v[i]==0
241
242 do // Insert all edges normal to axis i
243 {
244 k = bb.index(v);
245 for (index_t j = 0; j != end[i]; ++j)
246 {
247 if (cmat.bContains(k) && cmat.bContains(k + s))
248 {
249 // inefficient for now
250 const index_t kInd = m_xmatrix_offset[lvl] +
251 (std::lower_bound(cmat.begin(), cmat.end(), k)
252 - cmat.begin());
253
254 // inefficient for now
255 const index_t kNextInd = m_xmatrix_offset[lvl] +
256 (std::lower_bound(cmat.begin(), cmat.end(), k + s)
257 - cmat.begin());
258
259 mesh.addEdge(kInd, kNextInd);
260 }
261 k += s;
262 }
263 } while (nextCubePoint(v, low, upp));
264 }
265}
266
267template<short_t d, class T>
268void gsHTensorBasis<d, T>::connectivity(const gsMatrix<T> & nodes, int level, gsMesh<T> & mesh) const
269{
270 const index_t sz = size();
271 GISMO_ASSERT(nodes.rows() == sz, "Invalid input.");
272
273 // Add all vertices
274 for (index_t i = 0; i< sz; ++i)
275 mesh.addVertex(nodes.row(i).transpose());
276
277 addConnectivity(level, mesh);
278}
279
280template<short_t d, class T>
282{
283 const index_t sz = size();
284 GISMO_ASSERT( nodes.rows() == sz, "Invalid input.");
285
286 // Add vertices
287 for(index_t i = 0; i< sz; ++i )
288 mesh.addVertex( nodes.row(i).transpose() );
289
290 // For all levels
291 for(size_t lvl = 0; lvl <= maxLevel(); lvl++)
292 {
293 addConnectivity(lvl, mesh);
294 }
295}
296
297template<short_t d, class T>
299{
300 return m_xmatrix_offset.back();
301}
302template<short_t d, class T>
304{
305 std::vector<gsSortedVector<index_t> > OX = m_xmatrix;
306 refine(boxes);
307 gsSparseMatrix<> transf;
308 this->transfer(OX, transf);
309 gsDebug<<"tranf orig:\n"<<transf<<std::endl;
310 coefs = transf*coefs;
311}
312
313template<short_t d, class T>
314void gsHTensorBasis<d,T>::refineElements_withCoefs(gsMatrix<T> & coefs,std::vector<index_t> const & boxes)
315{
316 gsSparseMatrix<> transf;
317 this->refineElements_withTransfer(boxes,transf);
318 //gsDebug<<"tranf orig:\n"<<transf<<std::endl;
319 coefs = transf*coefs;
320}
321
322template<short_t d, class T>
323void gsHTensorBasis<d,T>::refineElements_withTransfer(std::vector<index_t> const & boxes, gsSparseMatrix<T> & tran)
324{
325 std::vector<gsSortedVector<index_t> > OX = m_xmatrix;
326 this->refineElements(boxes);
327 this->transfer(OX, tran);
328}
329
330template<short_t d, class T>
331void gsHTensorBasis<d,T>::refineElements_withTransfer2(std::vector<index_t> const & boxes, gsSparseMatrix<T> & tran)
332{
333 std::vector<gsSortedVector<index_t> > OX = m_xmatrix;
334 this->refineElements(boxes);
335 this->transfer2(OX, tran);
336}
337
338template<short_t d, class T>
339void gsHTensorBasis<d,T>::refineElements_withCoefs2(gsMatrix<T> & coefs,std::vector<index_t> const & boxes)
340{
341 std::vector<gsSortedVector<index_t> > OX = m_xmatrix;
342 refineElements(boxes);
343 gsSparseMatrix<> transf;
344 this->transfer2(OX, transf);
345 //gsDebug<<"tranf 2:\n"<<transf<<std::endl;
346 coefs = transf*coefs;
347}
348
349// template<short_t d, class T>
350// void gsHTensorBasis<d,T>::unrefine_withCoefs(gsMatrix<T> & coefs, gsMatrix<T> const & boxes)
351// {
352// std::vector<gsSortedVector<index_t> > OX = m_xmatrix;
353// unrefine(boxes);
354// gsSparseMatrix<> transf;
355// this->transfer(OX, transf);
356// gsDebug<<"tranf orig:\n"<<transf<<std::endl;
357// coefs = transf*coefs;
358// }
359
360template<short_t d, class T>
361void gsHTensorBasis<d,T>::unrefineElements_withCoefs(gsMatrix<T> & coefs,std::vector<index_t> const & boxes)
362{
363 gsSparseMatrix<> transf;
364 this->unrefineElements_withTransfer(boxes,transf);
365 //gsDebug<<"tranf orig:\n"<<transf<<std::endl;
366
367 typename gsSparseSolver<T>::QR solver(transf);
368 coefs=solver.solve(coefs);
369}
370
371template<short_t d, class T>
372void gsHTensorBasis<d,T>::unrefineElements_withTransfer(std::vector<index_t> const & boxes, gsSparseMatrix<T> & tran)
373{
374 typename gsHTensorBasis<d,T>::uPtr cp = this->clone();
375 this->unrefineElements(boxes);
376 cp->transfer(this->m_xmatrix,tran);
377}
378
379template<short_t d, class T>
380void gsHTensorBasis<d,T>::uniformRefine_withCoefs(gsMatrix<T>& coefs, int numKnots, int mul, int dir)
381{
382 GISMO_UNUSED(dir);
383 GISMO_ASSERT(dir==-1,"Direction is not implemented");
384
385 std::vector<gsSortedVector<index_t> > OX = m_xmatrix;
386 //uniformRefine(numKnots);
387 //gsMatrix<> transf;
388 //CMatrix temp;
389 //this->m_xmatrix.insert(this->m_xmatrix.begin(), temp);
390 //gsVector<index_t> level;
391 //gsMatrix<index_t> p1, p2;
392 //this->m_tree.getBoxes(p1,p2,level);
393 std::vector<index_t> boxes;
394 index_t lvl;
395 for ( typename hdomain_type::literator it = m_tree.beginLeafIterator(); it.good(); it.next() )
396 {
397 // gsDebug <<" level : "<< it.level() <<"\n";
398 // gsDebug <<" lower : "<< it.lowerCorner() <<"\n";
399 // gsDebug <<" upper : "<< it.upperCorner() <<"\n";
400
401 lvl = it.level() + 1;
402 const point & l = it.lowerCorner();
403 const point & u = it.upperCorner();
404
405 boxes.push_back(lvl);
406 for( short_t i = 0; i < d; i++)
407 boxes.push_back( l(i) * 2);
408 for( short_t i = 0; i < d; i++)
409 boxes.push_back( u(i) * 2);
410 }
411
412 this->clone()->refineElements_withCoefs(coefs, boxes);
413 this->uniformRefine(numKnots, mul);
414 //this->m_xmatrix.erase(this->m_xmatrix.begin(),this->m_xmatrix.begin()+1);
415 //coefs = transf*coefs;
416}
417
418template<short_t d, class T>
420{
421 std::vector<gsSortedVector<index_t> > OX = m_xmatrix;
422 std::vector<index_t> boxes;
423 index_t lvl;
424 for ( typename hdomain_type::literator it = m_tree.beginLeafIterator(); it.good(); it.next() )
425 {
426 if (it.level() == 0)
427 continue;
428 lvl = it.level() - 1;
429 const point & l = it.lowerCorner();
430 const point & u = it.upperCorner();
431
432 boxes.push_back(lvl);
433 for( short_t i = 0; i < d; i++)
434 boxes.push_back( l(i) / 2);
435 for( short_t i = 0; i < d; i++)
436 boxes.push_back( u(i)/2 + (index_t)(u(i)%2!=0));
437 }
438
439 this->clone()->unrefineElements_withCoefs(coefs, boxes);
440 this->uniformCoarsen(numKnots);
441}
442
443
444template<short_t d, class T>
445void gsHTensorBasis<d,T>::refine(gsMatrix<T> const & boxes, int refExt)
446{
447 GISMO_ASSERT(boxes.rows() == d, "refine() needs d rows of boxes.");
448 GISMO_ASSERT(boxes.cols()%2 == 0, "Each box needs two corners but you don't provide refine() with them.");
449
450#ifndef NDEBUG
451 gsMatrix<T> para = support();
452 for(int i = 0; i < boxes.cols()/2; i++)
453 {
454 for( short_t j = 0; j < d; j++ )
455 {
456 GISMO_ASSERT( para(j,0) <= boxes(j, 2*i) ,
457 "In refine() the first corner is outside the computational domain.");
458 GISMO_ASSERT( para(j,1) >= boxes(j, 2*i+1),
459 "In refine() the second corner is outside the computational domain." );
460 }
461 }
462#endif
463
464 if( refExt == 0 )
465 {
466 // If there is no refinement-extension, just use the
467 // "regular" refinement function refine( gsMatrix )
468 this->refine( boxes );
469
470 // Make sure there are enough levels
471 needLevel( m_tree.getMaxInsLevel() );
472 }
473 else
474 {
475 // Make an element vector
476 std::vector<index_t> refVector = this->asElements(boxes, refExt);//std::vector<unsigned> refVector = this->asElements(boxes, refExt);
477
478 // ...and refine
479 this->refineElements( refVector );
480 }
481
482 // Update the basis (already done by now)
483 //update_structure();
484}
485
486template<short_t d, class T>
487void gsHTensorBasis<d,T>::unrefine(gsMatrix<T> const & boxes, int refExt)
488{
489 GISMO_ASSERT(boxes.rows() == d, "refine() needs d rows of boxes.");
490 GISMO_ASSERT(boxes.cols()%2 == 0, "Each box needs two corners but you don't provide refine() with them.");
491
492#ifndef NDEBUG
493 gsMatrix<T> para = support();
494 for(int i = 0; i < boxes.cols()/2; i++)
495 {
496 for( short_t j = 0; j < d; j++ )
497 {
498 GISMO_ASSERT( para(j,0) <= boxes(j, 2*i) ,
499 "In refine() the first corner is outside the computational domain.");
500 GISMO_ASSERT( para(j,1) >= boxes(j, 2*i+1),
501 "In refine() the second corner is outside the computational domain." );
502 }
503 }
504#endif
505
506 // Make an element vector
507 std::vector<index_t> refVector = this->asElementsUnrefine(boxes, refExt);//std::vector<unsigned> refVector = this->asElements(boxes, refExt);
508
509 // ...and refine
510 this->unrefineElements( refVector );
511
512 // if( refExt == 0 )
513 // {
514 // // If there is no refinement-extension, just use the
515 // // "regular" unrefinement function unrefine( gsMatrix )
516 // // this->unrefine( boxes );
517
518 // // Make an element vector
519 // std::vector<index_t> refVector = this->asElements(boxes, 0);//std::vector<unsigned> refVector = this->asElements(boxes, refExt);
520
521 // // ...and refine
522 // this->unrefineElements( refVector );
523
524 // // Make an element vector
525 // // this->unrefine( boxes );
526 // }
527 // else
528 // {
529 // // Make an element vector
530 // std::vector<index_t> refVector = this->asElements(boxes, refExt);//std::vector<unsigned> refVector = this->asElements(boxes, refExt);
531
532 // // ...and refine
533 // this->unrefineElements( refVector );
534 // }
535
536 // Update the basis (already done by now)
537 //update_structure();
538}
539
540template<short_t d, class T>
541std::vector<index_t> gsHTensorBasis<d,T>::asElements(gsMatrix<T> const & boxes, int refExt) const
542{
543 // If there is a refinement-extension, we will have to use
544 // refineElements( std::vector )
545 //
546 // Each box will be represented by 2*d+1 entries specifying
547 // <level to be refined to>,<lower corner>,<upper corner>
548 const int offset = 2*d+1;
549
550 // Initialize vector of size
551 // "entries per box" times "number of boxes":
552 std::vector<index_t> refVector( offset * boxes.cols()/2 );
553 gsMatrix<T> ctr(d,1);
554
555 // Loop over all boxes:
556 for(index_t i = 0; i < boxes.cols()/2; i++)
557 {
558 ctr = ( boxes.col( 2*i ) + boxes.col( 2*i+1) )/2;
559
560 // Compute the level we want to refine to.
561 // Note that, if the box extends over several elements,
562 // the level at the centerpoint will be taken for reference
563 const int refLevel = getLevelAtPoint( ctr ) + 1;
564
565 // Make sure there are enough levels
566 needLevel( refLevel );
567
568 for(index_t j = 0; j < boxes.rows();j++)
569 {
570 // Convert the parameter coordinates to (unique) knot indices
571 const gsKnotVector<T> & kv = m_bases[refLevel]->knots(j);
572 index_t k1 = (std::upper_bound(kv.domainUBegin(), kv.domainUEnd(),
573 boxes(j,2*i ) ) - 1).uIndex();
574 index_t k2 = (std::upper_bound(kv.domainUBegin(), kv.domainUEnd()+1,
575 boxes(j,2*i+1) ) - 1).uIndex();
576
577 // Trivial boxes trigger some refinement
578 if ( k1 == k2)
579 {
580 if (0!=k1) {--k1;}
581 ++k2;
582 }
583
584 index_t maxKtIndex = kv.size();
585 if (m_manualLevels)
586 {
587 _knotIndexToDiadicIndex(refLevel,j,k1);
588 _knotIndexToDiadicIndex(refLevel,j,k2);
589
590 _knotIndexToDiadicIndex(refLevel,j,maxKtIndex);
591 }
592 // If applicable, add the refinement extension.
593 // Note that extending by one cell on level L means
594 // extending by two cells in level L+1
595 ( k1 < 2*refExt ? k1=0 : k1-=2*refExt );
596 ( k2 + 2*refExt >= maxKtIndex ? k2=maxKtIndex-1 : k2+=2*refExt);
597
598 // Store the data...
599 refVector[i*offset] = refLevel;
600 refVector[i*offset+1+j] = k1;
601 refVector[i*offset+1+j+d] = k2;
602 }
603 }
604 // gsDebug<<"begin\n";
605 // for (std::vector<unsigned>::const_iterator i = refVector.begin(); i != refVector.end(); ++i)
606 // std::cout << *i << ' ';
607 // gsDebug<<"end\n";
608
609 return refVector;
610}
611
612template<short_t d, class T>
613std::vector<index_t> gsHTensorBasis<d,T>::asElementsUnrefine(gsMatrix<T> const & boxes, int refExt) const
614{
615 // If there is a refinement-extension, we will have to use
616 // refineElements( std::vector )
617 //
618 // Each box will be represented by 2*d+1 entries specifying
619 // <level to be refined to>,<lower corner>,<upper corner>
620 const int offset = 2*d+1;
621
622 // Initialize vector of size
623 // "entries per box" times "number of boxes":
624 std::vector<index_t> refVector;
625 refVector.reserve( offset * boxes.cols()/2 );
626 gsMatrix<T> ctr(d,1);
627
628 // Index for the actual boxes that can be unrefined
629 // (so not the ones with level < 0 )
630 index_t I = 0;
631 // Loop over all boxes:
632 for(index_t i = 0; i < boxes.cols()/2; i++)
633 {
634 ctr = ( boxes.col( 2*i ) + boxes.col( 2*i+1) )/2;
635
636 // Compute the level we want to refine to.
637 // Note that, if the box extends over several elements,
638 // the level at the centerpoint will be taken for reference
639 const int refLevel = getLevelAtPoint( ctr ) - 1;
640
641 // If the level is 0, we cannot coarsen
642 if (refLevel < 0) continue;
643
644 refVector.resize(refVector.size() + offset);
645
646 for(index_t j = 0; j < boxes.rows();j++)
647 {
648 // Convert the parameter coordinates to (unique) knot indices
649 const gsKnotVector<T> & kv = m_bases[refLevel+1]->knots(j);
650 index_t k1 = (std::upper_bound(kv.domainUBegin(), kv.domainUEnd(),
651 boxes(j,2*i ) ) - 1).uIndex();
652 index_t k2 = (std::upper_bound(kv.domainUBegin(), kv.domainUEnd()+1,
653 boxes(j,2*i+1) ) - 1).uIndex();
654
655 // Trivial boxes trigger some refinement
656 if ( k1 == k2)
657 {
658 if (0!=k1) {--k1;}
659 ++k2;
660 }
661
662 index_t maxKtIndexd = kv.size();
663 if (m_manualLevels)
664 {
665 _knotIndexToDiadicIndex(refLevel,j,k1);
666 _knotIndexToDiadicIndex(refLevel,j,k2);
667
668 _knotIndexToDiadicIndex(refLevel,j,maxKtIndexd);
669 }
670
671 // If applicable, add the refinement extension.
672 ( k1 < refExt ? k1=0 : k1-=refExt );
673 ( k2 + refExt >= maxKtIndexd ? k2=maxKtIndexd-1 : k2+=refExt);
674
675 // go one level up;
676 k1 /= 2;
677 k2 = k2/2 + (index_t)(k2%2 != 0);
678
679 // Store the data...
680 refVector[I*offset] = refLevel;
681 refVector[I*offset+1+j] = k1;
682 refVector[I*offset+1+j+d] = k2;
683 }
684
685 I++;
686 }
687
688 return refVector;
689}
690
691template<short_t d, class T>
693{
694 GISMO_ASSERT(boxes.rows() == d, "refine() needs d rows of boxes.");
695 GISMO_ASSERT(boxes.cols()%2 == 0, "Each box needs two corners but you don't provide refine() with them.");
696
697#ifndef NDEBUG
698 gsMatrix<T> para = support();
699 for(int i = 0; i < boxes.cols()/2; i++)
700 {
701 for( short_t j = 0; j < d; j++ )
702 {
703 GISMO_ASSERT( para(j,0) <= boxes(j, 2*i) ,
704 "In refine() the first corner is outside the computational domain.");
705 GISMO_ASSERT( para(j,1) >= boxes(j, 2*i+1),
706 "In refine() the second corner is outside the computational domain." );
707 }
708 }
709#endif
710
711 gsVector<index_t,d> k1, k2;
712 for(index_t i = 0; i < boxes.cols()/2; i++)
713 {
714 // 1. Get a small cell containing the box
715 const int fLevel = m_bases.size()-1;
716
717 for(index_t j = 0; j < k1.size();j++)
718 {
719 const gsKnotVector<T> & kv = m_bases.back()->knots(j);
720 k1[j] = (std::upper_bound(kv.domainUBegin(), kv.domainUEnd(),
721 boxes(j,2*i ) ) - 1).uIndex();
722 k2[j] = (std::upper_bound(kv.domainUBegin(), kv.domainUEnd()+1,
723 boxes(j,2*i+1) ) - 1).uIndex();
724
725 // Trivial boxes trigger some refinement
726 if ( k1[j] == k2[j])
727 {
728 if (0!=k1[j]) {--k1[j];}
729 ++k2[j];
730 }
731 }
732
733 if (m_manualLevels)
734 {
735 _knotIndexToDiadicIndex(fLevel,k1);
736 _knotIndexToDiadicIndex(fLevel,k2);
737 }
738
739 // 2. Find the smallest level in which the box is completely contained
740 //const int level = m_tree.query3(k1,k2,fLevel) + 1;
741 // make sure that the grid is computed ( needLevel(level) )
742 //const tensorBasis & tb = tensorLevel(level);
743 //GISMO_UNUSED(tb);
744
745 // Sink box
746 m_tree.sinkBox(k1, k2, fLevel);
747 // Make sure we have enough levels
748 needLevel( m_tree.getMaxInsLevel() );
749 }
750
751 // Update the basis
752 update_structure();
753}
754
755// template<short_t d, class T>
756// void gsHTensorBasis<d,T>::unrefine(gsMatrix<T> const & boxes)
757// {
758// GISMO_ASSERT(boxes.rows() == d, "unrefine() needs d rows of boxes.");
759// GISMO_ASSERT(boxes.cols()%2 == 0, "Each box needs two corners but you don't provide unrefine() with them.");
760
761// #ifndef NDEBUG
762// gsMatrix<T> para = support();
763// for(int i = 0; i < boxes.cols()/2; i++)
764// {
765// for( short_t j = 0; j < d; j++ )
766// {
767// GISMO_ASSERT( para(j,0) <= boxes(j, 2*i) ,
768// "In unrefine() the first corner is outside the computational domain.");
769// GISMO_ASSERT( para(j,1) >= boxes(j, 2*i+1),
770// "In unrefine() the second corner is outside the computational domain." );
771// }
772// }
773// #endif
774
775// gsVector<index_t,d> k1, k2;
776// for(index_t i = 0; i < boxes.cols()/2; i++)
777// {
778// // 1. Get a small cell containing the box
779// const int fLevel = m_bases.size()-1;
780
781// for(index_t j = 0; j < k1.size();j++)
782// {
783// const gsKnotVector<T> & kv = m_bases.back()->knots(j);
784// k1[j] = (std::upper_bound(kv.domainUBegin(), kv.domainUEnd(),
785// boxes(j,2*i ) ) - 1).uIndex();
786// k2[j] = (std::upper_bound(kv.domainUBegin(), kv.domainUEnd()+1,
787// boxes(j,2*i+1) ) - 1).uIndex();
788
789// // Trivial boxes trigger some refinement
790// if ( k1[j] == k2[j])
791// {
792// if (0!=k1[j]) {--k1[j];}
793// ++k2[j];
794// }
795// }
796
797// // 2. Find the smallest level in which the box is completely contained
798// //const int level = m_tree.query3(k1,k2,fLevel) + 1;
799// // make sure that the grid is computed ( needLevel(level) )
800// //const tensorBasis & tb = tensorLevel(level);
801// //GISMO_UNUSED(tb);
802
803// // Sink box
804// m_tree.raiseBox(k1, k2, fLevel);
805// }
806
807// // Update the basis
808// // update_structure();
809// }
810
811template<short_t d, class T>
813{
814 // Get current level
815 const index_t lvl = this->levelOf(i);
816 // Get the support endpoints
818 m_bases[lvl]->elementSupport_into(m_xmatrix[lvl][ i - m_xmatrix_offset[lvl] ],
819 elements);
820 point low = elements.col(0);
821 point upp = elements.col(1);
822 // Advance the indices to one level deeper
823 for ( short_t k = 0; k!=d; ++k )
824 {
825 low[k] = low[k] << 1;
826 upp[k] = upp[k] << 1;
827 }
828 // Insert the domain to the lvl+1 nested domain
829 m_tree.insertBox(low,upp,lvl+1);
830 // Update the basis
831 update_structure();
832}
833
834
835/*
836 template<short_t d, class T>
837 void gsHTensorBasis<d,T>::refine(gsDomainIterator<T> const & boxes)
838 {
839
840 }
841*/
842
843template<short_t d, class T>
844void gsHTensorBasis<d,T>::refineElements(std::vector<index_t> const & boxes)
845{
846 point i1;
847 point i2;
848
849 GISMO_ASSERT( (boxes.size()%(2*d + 1))==0,
850 "The points did not define boxes properly. The boxes were not added to the basis.");
851 for(size_t i = 0; i < (boxes.size())/(2*d+1); i++)
852 {
853 for( short_t j = 0; j < d; j++ )
854 {
855 i1[j] = boxes[(i*(2*d+1))+j+1];
856 i2[j] = boxes[(i*(2*d+1))+d+j+1];
857 }
858 insert_box(i1,i2,boxes[i*(2*d+1)]);
859 }
860
861 update_structure();
862}
863
864template<short_t d, class T>
865void gsHTensorBasis<d,T>::unrefineElements(std::vector<index_t> const & boxes)
866{
867 point i1;
868 point i2;
869
870 GISMO_ASSERT( (boxes.size()%(2*d + 1))==0,
871 "The points did not define boxes properly. The boxes were not added to the basis.");
872
873 for(size_t i = 0; i < (boxes.size())/(2*d+1); i++)
874 {
875 for( short_t j = 0; j < d; j++ )
876 {
877 i1[j] = boxes[(i*(2*d+1))+j+1];
878 i2[j] = boxes[(i*(2*d+1))+d+j+1];
879 }
880
881 m_tree.clearBox(i1,i2, boxes[i*(2*d+1)]);
882 // needLevel( m_tree.getMaxInsLevel() );
883 }
884
885 // reconstruct the whole tree to fix alignment
886 gsHDomain<d> newtree( m_tree.upperCornerIndex() );
887 auto leafIt = m_tree.beginLeafIterator();
888 for (; leafIt.good(); leafIt.next())
889 {
890 if ( leafIt.level()>0 )
891 newtree.insertBox(leafIt.lowerCorner(),
892 leafIt.upperCorner(), leafIt.level() );
893 }
894 m_tree = newtree;
895
896 //recompute max-ins-level
897 m_tree.computeMaxInsLevel();
898 update_structure();
899}
900
901template<short_t d, class T>
903{
904 const index_t dir = side.direction();
905 const index_t par = side.parameter();
906 gsMatrix<T> rf = this->support();
907 rf(dir,!par) = rf(dir,par);
908 for (index_t i = 0; i!=lvl; ++i) // lazy impl., this can be more efficient
909 this->refine(rf);
910}
911
912template<short_t d, class T>
914 const gsBasis<T> & other,
915 gsMatrix<index_t> & bndThis,
916 gsMatrix<index_t> & bndOther,
917 index_t offset) const
918{
919 if( const Self_t * _other = dynamic_cast<const Self_t*>( &other) )
920 {
921 // tens1 will store the tensor-index on side second(),...
922 gsVector<index_t,d> N, tens0, tens1;
923
924 // see if the orientation is preserved on side second()
925 const gsVector<bool> dirOrient = bi.dirOrientation();
926 const gsVector<index_t> dirMap = bi.dirMap();
927
928 // get the global indices of the basis functions which are
929 // active on the interface
930 bndThis = this->boundaryOffset( bi.first().side(), offset );
931
932 // this is only for checking whether, at least, both involved
933 // bases have the same number of DOF on the interface.
934 bndOther= _other->boundaryOffset( bi.second().side(), offset );
935 GISMO_ASSERT( bndThis.rows() == bndOther.rows(),
936 "Input error, sizes do not match: "
937 <<bndThis.rows()<<"!="<<bndOther.rows() );
938 // bndOther gets overwritten completely, so here is the setZero():
939 bndOther.setZero();
940
941 for( index_t i=0; i < bndThis.rows(); i++)
942 {
943 // get the level of the basis function on side first()
944 index_t L = this->levelOf( bndThis(i,0) );
945 // get the flat tensor index
946 // (i.e., the single-number-index on level L)...
947 index_t flat0 = this->flatTensorIndexOf( bndThis(i,0) );
948 // ... and change it to the tensor-index.
949 tens0 = this->tensorLevel(L).tensorIndex( flat0 );
950
951 // ...flat1 the corresponding flat index
952 // (single-number on level)...
953 index_t flat1 = 0;
954 // ...and cont1 the corresponding continued (global) index.
955 index_t cont1 = 0;
956
957 // get the sizes of the components of the tensor-basis on this level,
958 // i.e., the sizes of the univariate bases corresponding
959 // to the respective coordinate directions
960 for( short_t j=0; j < d; j++)
961 N[j] = _other->tensorLevel(L).component(j).size();
962
963 // get the tensor-index of the basis function on level L on
964 // second() that should be matched with flatp/tens0
965 for( short_t j=0; j<d; j++)
966 {
967 // coordinate direction j on first() gets
968 // mapped to direction jj on second()
969 index_t jj = dirMap[j];
970 // store the respective component of the tensor-index
971 tens1[jj] = tens0[j];
972
973 if( jj == bi.second().direction() )
974 {
975 // if jj is the direction() of the interface,
976 // however, we need either the first
977 // or last basis function
978 if( bi.second().parameter() ) // true = 1 = end
979 tens1[jj] = N[jj]-(1+offset);
980 else
981 tens1[jj] = offset;
982 }
983 else
984 {
985 // otherwise, check if the orientation is
986 // preserved. If necessary, flip it.
987 if( !dirOrient[j] )
988 tens1[jj] = N[jj]-1 - tens1[jj];
989 }
990
991 }
992
993 flat1 = _other->tensorLevel(L).index( tens1 );
994
995 // compute the "continuous" index on second(), i.e., the index
996 // in the numbering which is global over all levels.
997 cont1 = _other->flatTensorIndexToHierachicalIndex( flat1, L );
998 // this is the index that has to be matched with bndThis(i,0)
999 bndOther( i, 0 ) = cont1;
1000 }
1001 return;
1002 }
1003 gsWarn<<"Cannot match with "<< other <<"\n";
1004}
1005
1006template<short_t d, class T>
1008 const gsBasis<T> & other,
1009 gsMatrix<index_t> & bndThis,
1010 gsMatrix<index_t> & bndOther) const
1011{
1012 matchWith(bi,other,bndThis,bndOther,0);
1013}
1014
1015
1016//protected functions
1017
1018// Construct the characteristic matrix of level \a level ; i.e., set
1019// all the matrix entries corresponding to active functions to one and
1020// the rest to zero.
1021template<short_t d, class T>
1022void gsHTensorBasis<d,T>::set_activ1(int level)
1023{
1024 typedef typename gsKnotVector<T>::smart_iterator knotIter;
1025
1026 //gsDebug<<" Setting level "<< level <<"\n";
1027 point low, upp;
1028
1029 CMatrix & cmat = m_xmatrix[level];
1030
1031 // Clear previous entries
1032 cmat.clear();
1033
1034 // If a level is to be checked which is larger than
1035 // the maximum inserted level, nothing need so to be done
1036 if ( level > static_cast<int>(m_tree.getMaxInsLevel() ) )
1037 return;
1038
1039
1040 gsVector<knotIter,d> starts, ends, curr;
1041 point ind;
1042 ind[0] = 0; // for d==1: warning: may be used uninitialized in this function (snap-ci)
1043
1044 for(short_t i = 0; i != d; ++i)
1045 {
1046 // beginning of the iteration in i-th direction
1047 starts[i] = m_bases[level]->knots(i).sbegin() ;
1048 // end of the iteration in i-th direction
1049 ends [i] = m_bases[level]->knots(i).send()-m_deg[i]-1;
1050 }
1051
1052 curr = starts;// start iteration
1053 do
1054 {
1055 for(short_t i = 0; i != d; ++i)
1056 {
1057 low[i] = curr[i].uIndex(); // lower left corner of the support of the function
1058 upp[i] = (curr[i]+m_deg[i]+1).uIndex(); // upper right corner of the support
1059 ind[i] = curr[i].index(); // index of the function in the matrix
1060 }
1061
1062 //Here: get knot indices in some standard indexing (eg. dyadic)
1063 if (m_manualLevels)
1064 {
1065 _knotIndexToDiadicIndex(level,low);
1066 _knotIndexToDiadicIndex(level,upp);
1067 }
1068
1069 if ( m_tree.query3(low, upp,level) == level) //if active ????
1070 cmat.push_unsorted( m_bases[level]->index( ind ) );
1071
1072 }
1073 while ( nextLexicographic(curr,starts,ends) ); // while there are some functions (i.e., some combinations of iterators) left
1074
1075 cmat.sort();
1076
1077}
1078
1079template<short_t d, class T>
1080void gsHTensorBasis<d,T>::functionOverlap(const point & boxLow, const point & boxUpp,
1081 const int level, point & actLow, point & actUpp)
1082{
1083 const tensorBasis & tb = *m_bases[level];
1084 for(short_t i = 0; i != d; ++i)
1085 {
1086 actLow[i] = tb.knots(i).lastKnotIndex (boxLow[i]) - m_deg[i];
1087 actUpp[i] = tb.knots(i).firstKnotIndex(boxUpp[i]) - 1 ;
1088
1089 // Note aao:
1090 //actLow[i] = firstKnotIndex(boxLow[i]);
1091 //actUpp[i] = tb.knots(i).lastKnotIndex(boxUpp[i])-m_deg[i]-1;
1092 }
1093}
1094
1095template<short_t d, class T>
1097{
1098 // for(size_t lvl = 0; lvl != m_xmatrix.size(); lvl++)
1099 // set_activ1(lvl);
1100 // return;
1101
1102 // iterate over leaf-boxes
1103 // for all overlapping supports with the box
1104 // set obvious to active
1105 // for the rest candidates (supp. not fully contained in box ~ !query2)
1106 // (equiv: actives on the boundary cells of the box)
1107 // query3(supp,box.level) == level (min. is level: no coarser)
1108 // take care: duplicates from different leaves or adj. cells
1109 point curr, actUpp;
1110 gsMatrix<index_t,d,2> elSupp;
1111
1112 // try: iteration per level
1113 for ( typename hdomain_type::literator it = m_tree.beginLeafIterator();
1114 it.good(); it.next() )
1115 {
1116 const int lvl = it.level();
1117 CMatrix & cmat = m_xmatrix[lvl];
1118
1119 // Get candidate functions
1120 functionOverlap(it.lowerCorner(), it.upperCorner(), lvl, curr, actUpp);
1121
1122 do
1123 {
1124 const index_t gi = m_bases[lvl]->index( curr );
1125
1126 // Get element support
1127 m_bases[lvl]->elementSupport_into(gi, elSupp);
1128
1129 if ( (elSupp.col(0).array() >= it.lowerCorner().array()).all() &&
1130 (elSupp.col(1).array() <= it.upperCorner().array()).all() )
1131 {
1132 // to do: all-at-once
1133 cmat.push_unsorted( gi );
1134 }
1135 else
1136 {
1137 // Check if active (iff no overlap with level less than lvl)
1138 if ( m_tree.query3(elSupp.col(0), elSupp.col(1), lvl) == lvl)
1139 cmat.push_unsorted( gi );
1140 }
1141 }
1142 while( nextCubePoint(curr, actUpp) );
1143 }
1144
1145 for(size_t lvl = 0; lvl != m_xmatrix.size(); ++lvl)
1146 {
1147 m_xmatrix[lvl].sort();
1148 m_xmatrix[lvl].erase( std::unique( m_xmatrix[lvl].begin(), m_xmatrix[lvl].end() ),
1149 m_xmatrix[lvl].end() );
1150 }
1151}
1152
1153template<short_t d, class T>
1155 std::vector<CMatrix> & x_matrix_lvl) const
1156{
1157 x_matrix_lvl.resize(level+1);
1158
1159 // vectors of iterators. one iterator for each dimension
1162 ind[0] = 0; // for d==1: warning: may be used uninitialized in this function (snap-ci)
1163 gsVector<index_t,d> low, upp;
1164 for(int j =0; j < level+1; j++)
1165 {
1166 // Clear previous entries
1167 x_matrix_lvl[j].clear();
1168
1169 for(index_t i = 0; i != d; ++i)
1170 {
1171 // beginning of the iteration in i-th direction
1172 starts[i] = m_bases[j]->knots(i).sbegin() ;
1173 // end of the iteration in i-th direction
1174 ends [i] = m_bases[j]->knots(i).send()-m_deg[i]-1;
1175 }
1176
1177 curr = starts; // set start of iteration
1178 do
1179 {
1180 for(index_t i = 0; i != d; ++i)
1181 {
1182 low[i] = curr[i].uIndex(); // lower left corner of the support of the function
1183 upp[i] = (curr[i]+m_deg[i]+1).uIndex(); // upper right corner of the support
1184 ind[i] = curr[i].index(); // index of the function in the matrix
1185 }
1186 if(j < level)
1187 {
1188 if ( m_tree.query3(low, upp,j) == j)
1189 { //if active
1190 x_matrix_lvl[j].push_unsorted( m_bases[j]->index( ind ) );
1191 }
1192 }else{
1193 if ( m_tree.query3(low, upp,j) >= j)
1194 { //if active
1195 x_matrix_lvl[j].push_unsorted( m_bases[j]->index( ind ) );
1196 }
1197 }
1198 }
1199 while ( nextLexicographic(curr,starts,ends) ); // while there are some functions (i.e., some combinations of iterators) left
1200
1201 x_matrix_lvl[j].sort();
1202 }
1203
1204}
1205
1207template<short_t d, class T> inline
1208void gsHTensorBasis<d,T>::insert_box(typename gsHTensorBasis<d,T>::point const & k1,
1209 typename gsHTensorBasis<d,T>::point const & k2,
1210 int lvl )
1211{
1212 // Remember box in History (for debugging)
1213 // m_boxHistory.push_back( box(k1,k2,lvl) );
1214
1215 m_tree.insertBox(k1,k2, lvl);
1216 needLevel( m_tree.getMaxInsLevel() );
1217}
1218
1219template<short_t d, class T>
1221{
1222 // Compress the tree
1223 // m_tree.makeCompressed();
1224
1225 while ( ! m_xmatrix_offset[1] )
1226 {
1227 delete m_bases.front();
1228 m_bases.erase( m_bases.begin() );
1229 m_tree.decrementLevel();
1230 m_xmatrix.erase( m_xmatrix.begin() );
1231 m_xmatrix_offset.erase( m_xmatrix_offset.begin() );
1232 }
1233 // Note/to do: cleaning up empty levels at the end as well.
1234}
1235
1236template<short_t d, class T>
1238{
1239 GISMO_ASSERT( static_cast<size_t>(level) < m_xmatrix.size(), "Requested level does not exist.\n");
1240
1241 CMatrix::const_iterator xmat_pointer = m_xmatrix[level].begin();
1242 CMatrix::const_iterator xmat_end = m_xmatrix[level].end();
1243 gsSortedVector< int >::iterator ind_pointer = indexes.begin();
1244 gsSortedVector< int >::iterator ind_end = indexes.end();
1245 index_t index = 0;
1246 while(ind_pointer!=ind_end&&xmat_pointer!=xmat_end)
1247 {
1248 if(*ind_pointer<static_cast<int>(*xmat_pointer))
1249 {
1250 (*ind_pointer)=-1;
1251 ++ind_pointer;
1252 }
1253 else if(*ind_pointer==static_cast<int>(*xmat_pointer))
1254 {
1255 (*ind_pointer)=m_xmatrix_offset[level]+index;
1256 ++xmat_pointer;
1257 ++index;
1258 ++ind_pointer;
1259 }
1260 else
1261 {
1262 ++xmat_pointer;
1263 ++index;
1264 }
1265 }
1266 while(ind_pointer!=ind_end)
1267 {
1268 (*ind_pointer)=-1;
1269 ++ind_pointer;
1270 }
1271}
1272
1273template<short_t d, class T>
1275{
1276 if( m_xmatrix.size()<=static_cast<size_t>(level) )
1277 return -1;
1278 CMatrix::const_iterator it = std::lower_bound(m_xmatrix[level].begin(), m_xmatrix[level].end(), index );
1279 if(it == m_xmatrix[level].end() || *it != index)
1280 return -1;
1281 else
1282 return m_xmatrix_offset[level] + ( it - m_xmatrix[level].begin() );
1283}
1284
1285template<short_t d, class T>
1286void gsHTensorBasis<d,T>::activeBoundaryFunctionsOfLevel(const unsigned level,const boxSide & s,std::vector<bool>& actives) const
1287{
1288 needLevel( level );
1289
1290 const gsMatrix<index_t> bound = m_bases[level]->boundary(s);
1291 const index_t sz = bound.rows();
1292 //gsSortedVector< int > indexes(bound->data(),bound->data()+sz);
1293 gsSortedVector< int > indexes;
1294 indexes.resize(sz,-1);
1295 if(level<=maxLevel())
1296 {
1297 for(index_t i = 0;i<sz;++i)
1298 indexes[i]=(bound)(i,0);
1299 flatTensorIndexesToHierachicalIndexes(indexes,level);
1300 }
1301 actives.resize(indexes.size(),false);
1302 std::fill (actives.begin(),actives.end(),false);
1303 for(size_t i = 0;i<indexes.size();i++)
1304 if(indexes[i]!=-1)
1305 actives[i]=true;
1306}
1307
1308template<short_t d, class T>
1309void gsHTensorBasis<d,T>::update_structure() // to do: rename as updateHook
1310{
1311 // Make sure we have computed enough levels
1312 needLevel( m_tree.getMaxInsLevel() );
1313
1314 // Setup the characteristic matrices
1315 m_xmatrix.clear();
1316 m_xmatrix.resize( m_tree.getMaxInsLevel()+1 );
1317 // Compress the tree
1318 m_tree.makeCompressed();
1319
1320 for(size_t i = 0; i != m_xmatrix.size(); i ++)
1321 set_activ1(i);
1322
1323 // Store all indices of active basis functions to m_matrix
1324 //setActive();
1325
1326 // Compute offsets
1327 m_xmatrix_offset.clear();
1328 m_xmatrix_offset.reserve(m_xmatrix.size()+1);
1329 m_xmatrix_offset.push_back(0);
1330 for (size_t i = 0; i != m_xmatrix.size(); i++)
1331 {
1332 m_xmatrix_offset.push_back(
1333 m_xmatrix_offset.back() + m_xmatrix[i].size() );
1334 }
1335}
1336
1337template<short_t d, class T>
1338void gsHTensorBasis<d,T>::needLevel(int maxLevel) const
1339{
1340 GISMO_ENSURE(!m_manualLevels || (size_t)(maxLevel)<m_uIndices.size(),"Maximum manual level reached, maxLevel = "<<maxLevel<<", m_uIndices.size() = "<<m_uIndices.size());
1341 // +1 for the initial basis in m_bases
1342 const int extraLevels = maxLevel + 1 - m_bases.size();
1343 for ( int i = 0; i < extraLevels; ++i )
1344 {
1345 tensorBasis * next_basis = m_bases.back()->clone().release();
1346 next_basis->uniformRefine(1);
1347 m_bases.push_back (next_basis); //note: m_bases is mutable
1348 }
1349}
1350
1351template<short_t d, class T>
1353{
1354 // Degrees
1355 //m_deg = tbasis.cwiseDegree();
1356 m_deg.resize(d);
1357 for( index_t i = 0; i < d; i++)
1358 m_deg[i] = tbasis.degree(i);
1359
1360 // Construct the initial basis
1361 if ( const tensorBasis * tb2 =
1362 dynamic_cast<const tensorBasis*>(&tbasis) )
1363 {
1364 m_bases.push_back(tb2->clone().release());
1365 }
1366 else
1367 {
1368 GISMO_ERROR("Cannot construct a Hierarchical basis from "<< tbasis );
1369 }
1370
1371 // Initialize the binary tree
1372 point upp;
1373 for ( index_t i = 0; i!=d; ++i )
1374 upp[i] = m_bases[0]->knots(i).numElements();
1375
1376 m_tree.init(upp);
1377
1378 // Need one level at least, in case refine(gsMatrix<T> boxes) is called
1379 this->needLevel(1);
1380}
1381
1382
1383template<short_t d, class T>
1385{
1386 gsMatrix<T> currPoint;
1387 point low, upp, cur;
1388 const int maxLevel = m_tree.getMaxInsLevel();
1389
1390 std::vector<std::vector<index_t> > temp_output;//collects the outputs
1391 temp_output.resize( u.cols() );
1392 size_t sz = 0;
1393
1394 //gsMatrix<index_t> activesLvl;
1395
1396 for(index_t p = 0; p < u.cols(); p++) //for all input points
1397 {
1398 currPoint = u.col(p);
1399 for(short_t i = 0; i != d; ++i)
1400 low[i] = m_bases[maxLevel]->knots(i).uFind( currPoint(i,0) ).uIndex();
1401
1402 if (m_manualLevels)
1403 this->_knotIndexToDiadicIndex(maxLevel,low);
1404
1405 // Identify the level of the point
1406 const int lvl = m_tree.levelOf(low, maxLevel);
1407 for(int i = 0; i <= lvl; i++)
1408 {
1409 /*
1410 m_bases[i]->active_into(currPoint, activesLvl);
1411
1412 std::set_intersection(m_xmatrix[i].begin(), m_xmatrix[i].end(),
1413 activesLvl.data(), activesLvl.data() + activesLvl.size(),
1414 std::back_inserter( temp_output[p] ) );
1415 +++ Renumbering to H-basis indexing
1416 */
1417
1418 m_bases[i]->active_cwise(currPoint, low, upp);
1419 cur = low;
1420 do
1421 {
1422 CMatrix::const_iterator it =
1423 m_xmatrix[i].find_it_or_fail( m_bases[i]->index(cur) );
1424
1425 if( it != m_xmatrix[i].end() )// if index is found
1426 {
1427 temp_output[p].push_back(
1428 this->m_xmatrix_offset[i] + (it - m_xmatrix[i].begin() )
1429 );
1430 }
1431 }
1432 while( nextCubePoint(cur,low,upp) );
1433 }
1434
1435 // update result size
1436 if ( temp_output[p].size() > sz )
1437 sz = temp_output[p].size();
1438 }
1439
1440 result.resize(sz, u.cols() );
1441 for(index_t i = 0; i < result.cols(); i++)
1442 {
1443 result.col(i).topRows(temp_output[i].size())
1444 = gsAsConstVector<index_t>(temp_output[i]);
1445 result.col(i).bottomRows(sz-temp_output[i].size()).setZero();
1446 }
1447}
1448
1449template<short_t d, class T>
1451{
1452 std::vector<index_t> temp;
1454 for(size_t i = 0; i != m_xmatrix[i].size(); i++)
1455 for (CMatrix::const_iterator it = m_xmatrix[i].begin();
1456 it != m_xmatrix[i].end(); it++)
1457 {
1458 ind = this->m_bases[i]->tensorIndex(*it);
1459 for (unsigned j=0; j!=d; ++j )
1460 if ( (ind[j]==0) || (ind[j]==(this->m_bases[i]->size(j)-1)) )
1461 {
1462 temp.push_back(m_xmatrix_offset[i] + (it-m_xmatrix[i].begin()) );
1463 break;
1464 }
1465 }
1466 return makeMatrix<index_t>(temp.begin(),temp.size(),1 );
1467}
1468
1469template<short_t d, class T>
1471boundaryOffset(boxSide const & s,index_t offset) const
1472{
1473 //get information on the side
1474 const index_t k = s.direction();
1475 const bool par = s.parameter();
1476
1477 std::vector<index_t> temp;
1479 // i goes through all levels of the hierarchical basis
1480 GISMO_ASSERT(this->maxLevel() < this->m_bases.size(),"Something went wrong: maxLevel() < m_bases.size(), "<<this->maxLevel()<<" < "<<m_bases.size());
1481 needLevel(m_xmatrix.size()-1);
1482
1483 for(size_t i = 0; i != m_xmatrix.size(); i++)
1484 {
1485 GISMO_ASSERT(static_cast<int>(offset)<this->m_bases[i]->size(k),
1486 "Offset ("<<offset<<") cannot be bigger than the amount of basis"
1487 "functions orthogonal to Boxside s! ("<<this->m_bases[i]->size(k)<<")");
1488
1489 index_t r = ( par ? this->m_bases[i]->size(k) - 1 -offset : offset);
1490 for (CMatrix::const_iterator it = m_xmatrix[i].begin();
1491 it != m_xmatrix[i].end(); it++)
1492 {
1493 ind = this->m_bases[i]->tensorIndex(*it);
1494 if ( ind[k]==r )
1495 temp.push_back(
1496 m_xmatrix_offset[i] + (it - m_xmatrix[i].begin() )
1497 );
1498 }
1499 }
1500 return makeMatrix<index_t>(temp.begin(),temp.size(),1 );
1501}
1502
1503template<short_t d, class T>
1505boundaryOffsetLevel(boxSide const & s,index_t offset, index_t level) const
1506{
1507 //get information on the side
1508 //index_t k = s.direction();
1509 //bool par = s.parameter();
1510 return this->m_bases[level]->boundaryOffset(s, offset);
1511}
1512
1513template<short_t d, class T>
1514index_t gsHTensorBasis<d,T>::
1515levelAtCorner(boxCorner const & c) const
1516{
1517 // Get parametric points of corner
1518 gsVector<bool> pars;
1519 c.parameters_into(d,pars);
1520 // Get the support of the underlying tensor basis
1521 gsMatrix<T> supp = m_bases.front()->support();
1522 // Assign the extreme of the support depending on the parametric coordinate
1523 gsVector<T> vec(supp.rows());
1524 for (index_t r = 0; r!=supp.rows(); r++)
1525 vec(r) = supp(r,pars(r));
1526
1527 return getLevelAtPoint(vec);
1528}
1529
1530template<short_t d, class T>
1531index_t gsHTensorBasis<d,T>::
1532functionAtCorner(boxCorner const & c) const
1533{
1534 index_t lvl = this->levelAtCorner(c);
1535
1536 // Get the index of the corner on the level
1537 index_t index = m_bases[lvl]->functionAtCorner(c);
1538
1539 // Transform to (continuous) HBspline index.
1540 return this->flatTensorIndexToHierachicalIndex(index,lvl);
1541}
1542
1543template<short_t d, class T>
1544index_t gsHTensorBasis<d,T>::
1545functionAtCorner(boxCorner const & c, index_t level) const
1546{
1547 // Get the index of the corner on the level
1548 index_t index = m_bases[level]->functionAtCorner(c);
1549
1550 // Transform to (continuous) HBspline index.
1551 return this->flatTensorIndexToHierachicalIndex(index,level);
1552}
1553
1554/*
1555template<short_t d, class T>
1556void gsHTensorBasis<d,T>::evalAllDers_into(const gsMatrix<T> & u, int n,
1557 std::vector<gsMatrix<T> >& result, bool sameElement) const;
1558{
1559 result.resize(n+1);
1560
1561}
1562*/
1563
1564template<short_t d, class T>
1565void gsHTensorBasis<d,T>::uniformRefine(int numKnots, int mul, int dir)
1566{
1567 GISMO_UNUSED(numKnots);
1568 GISMO_ASSERT(numKnots == 1, "Only implemented for numKnots = 1");
1569
1570 GISMO_ASSERT( m_tree.getMaxInsLevel() < static_cast<unsigned>(m_bases.size()),
1571 "Problem with max inserted levels: "<< m_tree.getMaxInsLevel()
1572 <<"<" << m_bases.size() <<"\n");
1573
1574 // Keep consistency of finest level
1575 tensorBasis * last_basis = m_bases.back()->clone().release();
1576 last_basis->uniformRefine(1,mul,dir);
1577 m_bases.push_back( last_basis );
1578
1579 // Delete the first level
1580 delete m_bases.front();
1581 m_bases.erase( m_bases.begin() );
1582
1583 // Lift all indices in the tree by one level
1584 m_tree.multiplyByTwo();
1585
1586 update_structure();
1587}
1588
1589template<short_t d, class T>
1591{
1592 GISMO_UNUSED(numKnots);
1593 GISMO_ASSERT(numKnots == 1, "Only implemented for numKnots = 1");
1594
1595 tensorBasis * first_basis = m_bases.front()->clone().release();
1596 first_basis->uniformCoarsen(1);
1597 m_bases.insert( m_bases.begin(), first_basis );
1598
1599 // Delete the last level
1600 delete m_bases.back();
1601 m_bases.pop_back();
1602
1603 // Lift all indices in the tree by one level
1604 m_tree.divideByTwo();
1605 // What happens when zero interior knots???????
1606
1607 update_structure();
1608}
1609
1610
1611template<short_t d, class T>
1612std::vector< std::vector< std::vector<index_t > > > gsHTensorBasis<d,T>::domainBoundariesParams( std::vector< std::vector< std::vector< std::vector< T > > > >& result) const
1613{
1614 std::vector< std::vector< std::vector< std::vector<index_t > > > > dummy;
1615 return domainBoundariesGeneric( dummy, result, false );
1616}
1617
1618template<short_t d, class T>
1619std::vector< std::vector< std::vector<index_t > > > gsHTensorBasis<d,T>::domainBoundariesIndices( std::vector< std::vector< std::vector< std::vector<index_t > > > >& result ) const
1620{
1621 std::vector< std::vector< std::vector< std::vector< T > > > > dummy;
1622 return domainBoundariesGeneric( result, dummy, true );
1623}
1624
1625
1626template<short_t d, class T>
1627std::vector< std::vector< std::vector<index_t > > > gsHTensorBasis<d,T>::domainBoundariesGeneric(std::vector< std::vector< std::vector< std::vector<index_t > > > >& indices,
1628 std::vector< std::vector< std::vector< std::vector< T > > > >& params,
1629 bool indicesFlag ) const
1630{
1631 indices.clear();
1632 params.clear();
1633 std::vector< std::vector< std::vector< int > > > res_aabb;
1634 std::vector< std::vector< std::vector<index_t > > > res_aabb_unsigned;
1635 std::vector< std::vector< std::vector< std::vector< index_t > > > > polylines;
1636
1637 polylines = this->m_tree.getPolylines();
1638 res_aabb.resize( polylines.size() );
1639 // We cannot simply say
1640 // result = this->m_tree.getPolylines();
1641 // because the return value of getPolylines() are vectors of ints and we have no implicit conversion of int to T.
1642
1643 // We want indices/params to be of the same size as polylines. We achieve this here and in the for cycles.
1644 if( indicesFlag )
1645 indices.resize( polylines.size() );
1646
1647 else
1648 params.resize(polylines.size());
1649
1650 int maxLevel = static_cast<int>( this->maxLevel() );
1651 // We precompute the parameter values corresponding to indices of m_maxInsLevel
1652 // although we don't need them if indicesFlag == true.
1653 std::vector<T> x_dir(m_bases[maxLevel]->knots(0).unique());
1654 std::vector<T> y_dir(m_bases[maxLevel]->knots(1).unique());
1655
1656 for(unsigned int i0 = 0; i0 < polylines.size(); i0++)
1657 {
1658 if( indicesFlag )
1659 indices[i0].resize( polylines[i0].size() );
1660 else
1661 params[i0].resize( polylines[i0].size() );
1662
1663 res_aabb[i0].resize( polylines[i0].size() );
1664 for(unsigned int i1 = 0; i1 < polylines[i0].size(); i1++)
1665 {
1666 if( indicesFlag )
1667 indices[i0][i1].resize( polylines[i0][i1].size() );
1668 else
1669 params[i0][i1].resize( polylines[i0][i1].size() );
1670
1671 res_aabb[i0][i1].resize( 4 );
1672 res_aabb[i0][i1][0] = 1000000;
1673 res_aabb[i0][i1][1] = 1000000;
1674 res_aabb[i0][i1][2] = -10000000;
1675 res_aabb[i0][i1][3] = -10000000;
1676 for(unsigned int i2 = 0; i2 < polylines[i0][i1].size(); i2++)
1677 {
1678 if( indicesFlag )
1679 {
1680 indices[i0][i1][i2].resize(4);
1681 indices[i0][i1][i2][0] = polylines[i0][i1][i2][0];
1682 indices[i0][i1][i2][1] = polylines[i0][i1][i2][1];
1683 indices[i0][i1][i2][2] = polylines[i0][i1][i2][2];
1684 indices[i0][i1][i2][3] = polylines[i0][i1][i2][3];
1685
1686 }
1687 else
1688 {
1689 params[i0][i1][i2].resize(4);
1690 // We could as well successively push_back() them but this should be slightly more efficient.
1691 params[i0][i1][i2][0] = (x_dir[polylines[i0][i1][i2][0]]);
1692 params[i0][i1][i2][1] = (y_dir[polylines[i0][i1][i2][1]]);
1693 params[i0][i1][i2][2] = (x_dir[polylines[i0][i1][i2][2]]);
1694 params[i0][i1][i2][3] = (y_dir[polylines[i0][i1][i2][3]]);
1695 }
1696 if(res_aabb[i0][i1][0]>signed(polylines[i0][i1][i2][0]))
1697 {
1698 res_aabb[i0][i1][0] = polylines[i0][i1][i2][0];
1699 }
1700 if(res_aabb[i0][i1][1]>signed(polylines[i0][i1][i2][1]))
1701 {
1702 res_aabb[i0][i1][1] = polylines[i0][i1][i2][1];
1703 }
1704 if(res_aabb[i0][i1][2]<signed(polylines[i0][i1][i2][2]))
1705 {
1706 res_aabb[i0][i1][2] = polylines[i0][i1][i2][2];
1707 }
1708 if(res_aabb[i0][i1][3]<signed(polylines[i0][i1][i2][3]))
1709 {
1710 res_aabb[i0][i1][3] = polylines[i0][i1][i2][3];
1711 }
1712 }
1713 }
1714 }
1715
1716 res_aabb_unsigned.resize(res_aabb.size());
1717 for (unsigned int i = 0; i < res_aabb.size(); i++)
1718 {
1719 res_aabb_unsigned[i].resize(res_aabb[i].size());
1720 for (unsigned int j = 0; j < res_aabb[i].size(); j++)
1721 {
1722 res_aabb_unsigned[i][j].resize(res_aabb[i][j].size());
1723 for (unsigned int k = 0; k < res_aabb[i][j].size(); k++)
1724 {
1725 if(res_aabb[i][j][k]<0)
1726 gsWarn << "conversion form signed to unsigned\n";
1727 res_aabb_unsigned[i][j][k] = res_aabb[i][j][k];
1728 }
1729 }
1730 }
1731 return res_aabb_unsigned;
1732}
1733
1734
1735template<short_t d, class T>
1737{
1738 // Note: implementation assumes number of old + 1 m_bases exists in this basis
1739 needLevel( old.size() );
1740
1741 tensorBasis T_0_copy = this->tensorLevel(0);
1742
1743 std::vector< gsSparseMatrix<T,RowMajor> > transfer(m_bases.size()-1);
1744 std::vector<std::vector<T> > knots(d);
1745
1746 for(size_t i = 1; i < m_bases.size(); ++i)
1747 {
1748 //T_0_copy.uniformRefine_withTransfer(transfer[i-1], 1);
1749 for(unsigned dim = 0; dim != d; ++dim)
1750 {
1751 const gsKnotVector<T> & ckv = m_bases[i-1]->knots(dim);
1752 const gsKnotVector<T> & fkv = m_bases[i ]->knots(dim);
1753 ckv.symDifference(fkv, knots[dim]);
1754 // equivalent (dyadic ref.):
1755 // ckv.getUniformRefinementKnots(1, knots[dim]);
1756
1757 //gsDebug << "level: " << i << "\n"
1758 // << "direction: " << dim << "\n";
1759 //gsDebugVar(gsAsMatrix<T>(dirKnots));
1760 }
1761
1762 T_0_copy.refine_withTransfer(transfer[i-1], knots);
1763 }
1764
1765 // Add missing empty char. matrices
1766 while ( old.size() >= m_xmatrix.size() )
1767 m_xmatrix.push_back( gsSortedVector<index_t>() );
1768
1769 result = this->coarsening_direct(old, m_xmatrix, transfer);
1770
1771 // This function automatically adds additional characteristic matrices,
1772 // even if they are not needed.
1773 // Check whether the characteristic matrices corresponding to the finest
1774 // levels are actually used. If they are empty, i.e., if there are no
1775 // active functions on that level, drop them...
1776 while( m_xmatrix.back().size() == 0 )
1777 m_xmatrix.pop_back();
1778
1779 // ...similarly, erase all those fine bases which are actually not used.
1780 const int sizeDiff = static_cast<int>( m_bases.size() - m_xmatrix.size() );
1781 if( sizeDiff > 0 )
1782 {
1783 freeAll(m_bases.end() - sizeDiff, m_bases.end());
1784 m_bases.resize(m_xmatrix.size());
1785 }
1786
1787 result.makeCompressed();
1788}
1789
1790template<short_t d, class T>
1791void gsHTensorBasis<d,T>::transfer2(const std::vector<gsSortedVector<index_t> >& old, gsSparseMatrix<T>& result)
1792{
1793 // Note: implementation assumes number of old + 1 m_bases exists in this basis
1794 needLevel( old.size() );
1795
1796 tensorBasis T_0_copy = this->tensorLevel(0);
1797 std::vector< gsSparseMatrix<T,RowMajor> > transfer( m_bases.size()-1 );
1798 std::vector<std::vector<T> > knots(d);
1799
1800 for(size_t i = 1; i < m_bases.size(); ++i)
1801 {
1802 //T_0_copy.uniformRefine_withTransfer(transfer[i], 1);
1803 for(short_t dim = 0; dim != d; ++dim)
1804 {
1805 const gsKnotVector<T> & ckv = m_bases[i-1]->knots(dim);
1806 const gsKnotVector<T> & fkv = m_bases[i ]->knots(dim);
1807 ckv.symDifference(fkv, knots[dim]);
1808 // equivalent (dyadic ref.):
1809 // ckv.getUniformRefinementKnots(1, knots[dim]);
1810
1811 //gsDebug << "level: " << i << "\n"
1812 // << "direction: " << dim << "\n";
1813 //gsDebugVar(gsAsMatrix<T>(dirKnots));
1814 }
1815 T_0_copy.refine_withTransfer(transfer[i-1], knots);
1816 }
1817
1818 // Add missing empty char. matrices
1819 while ( old.size() >= m_xmatrix.size())
1820 m_xmatrix.push_back( gsSortedVector<index_t>() );
1821
1822 result = this->coarsening_direct2(old, m_xmatrix, transfer);
1823
1824 result.makeCompressed();
1825}
1826
1827
1828template<short_t d, class T>
1829void gsHTensorBasis<d,T>::increaseMultiplicity(index_t lvl, int dir, T knotValue, int mult)
1830{
1831 GISMO_ASSERT( static_cast<size_t>(lvl) < m_xmatrix.size(),
1832 "Requested level does not exist.\n");
1833
1834 if (m_bases[lvl]->knots(dir).has(knotValue))
1835 {
1836 for(unsigned int i =lvl;i < m_bases.size();i++)
1837 m_bases[i]->component(dir).insertKnot(knotValue,mult);
1838 }
1839 else
1840 {
1841 gsWarn<<"Knot value not in the given knot vector."<<std::endl;
1842 }
1843
1844 update_structure();
1845}
1846
1847
1848template<short_t d, class T>
1849void gsHTensorBasis<d,T>::increaseMultiplicity(index_t lvl, int dir, const std::vector<T> & knotValue, int mult)
1850{
1851 for(unsigned int k =0; k < knotValue.size(); ++k)
1852 {
1853 if (m_bases[lvl]->knots(dir).has(knotValue[k]))
1854 {
1855 for(unsigned int i =lvl;i < m_bases.size(); ++i)
1856 m_bases[i]->component(dir).insertKnot(knotValue[k],mult);
1857 }
1858 else
1859 {
1860 gsWarn<<"Knot value not in the given knot vector."<<std::endl;
1861 }
1862 }
1863 update_structure();
1864}
1865
1866template<short_t d, class T>
1867void gsHTensorBasis<d,T>::getBoxesAlongSlice( int dir, T par,std::vector<index_t>& boxes ) const
1868{
1869 gsMatrix<index_t> b1,b2;
1870 gsVector<index_t> level;
1871 m_tree.getBoxesInLevelIndex(b1,b2,level);
1872 gsVector<index_t> min,max;
1873 for(index_t i = 0;i<level.rows();i++)
1874 {
1875 min = b1.row(i);
1876 max = b2.row(i);
1877 const index_t l = level(i);
1878 const index_t par_index = m_bases[l]->knots(dir).uFind(par).uIndex();
1879 if(l>0 && (par_index>=min(dir)) && (par_index<=max(dir)))
1880 {
1881 boxes.push_back(l);
1882 for(index_t j=0;j<min.rows();++j)
1883 if(j!=dir)
1884 boxes.push_back(min(j));
1885 for(index_t j=0;j<max.rows();++j)
1886 if(j!=dir)
1887 boxes.push_back(max(j));
1888 }
1889 }
1890}
1891
1892template<short_t d, class T>
1893void gsHTensorBasis<d,T>::_knotIndexToDiadicIndex(const index_t level, const index_t dir, index_t & knotIndex) const
1894{
1895 GISMO_ASSERT(m_manualLevels,"Only works for manual levels");
1896 knotIndex = m_uIndices[level][dir][knotIndex];
1897}
1898
1899template<short_t d, class T>
1901{
1902 for (index_t r = 0; r != d; r++)
1903 this->_knotIndexToDiadicIndex(level,r,knotIndex[r]);
1904}
1905
1906template<short_t d, class T>
1907void gsHTensorBasis<d,T>::_diadicIndexToKnotIndex(const index_t level, const index_t dir, index_t & diadicIndex) const
1908{
1909 GISMO_ASSERT(m_manualLevels,"Only works for manual levels");
1910 typename std::vector<index_t>::const_iterator it = std::find_if(m_uIndices[level][dir].begin(),m_uIndices[level][dir].end(),[&diadicIndex](const index_t i) { return i >= diadicIndex; });
1911 GISMO_ASSERT(it!=m_uIndices[level][dir].end(),"Index not found");
1912 diadicIndex = std::distance(m_uIndices[level][dir].begin(),it);
1913}
1914
1915template<short_t d, class T>
1917{
1918 for (index_t r = 0; r != d; r++)
1919 this->_diadicIndexToKnotIndex(level,r,diadicIndex[r]);
1920}
1921
1922
1923template<short_t d, class T>
1924void gsHTensorBasis<d,T>::degreeElevate(int const & i, int const dir)
1925{
1926 for (size_t level=0;level<m_bases.size();++level)
1927 m_bases[level]->degreeElevate(i,dir);
1928
1929 for(unsigned c=0; c<d; ++c)
1930 m_deg[c]=m_bases[0]->degree(c);
1931
1932 this->update_structure();
1933}
1934
1935template<short_t d, class T>
1936void gsHTensorBasis<d,T>::degreeReduce(int const & i, int const dir)
1937{
1938 for (size_t level=0;level<m_bases.size();++level)
1939 m_bases[level]->degreeReduce(i,dir);
1940
1941 for(unsigned c=0; c<d; ++c)
1942 m_deg[c]=m_bases[0]->degree(c);
1943
1944 this->update_structure();
1945}
1946
1947template<short_t d, class T>
1948void gsHTensorBasis<d,T>::degreeIncrease(int const & i, int const dir)
1949{
1950 for (size_t level=0;level<m_bases.size();++level)
1951 m_bases[level]->degreeIncrease(i,dir);
1952
1953 for(unsigned c=0; c<d; ++c)
1954 m_deg[c]=m_bases[0]->degree(c);
1955
1956 this->update_structure();
1957}
1958
1959template<short_t d, class T>
1960void gsHTensorBasis<d,T>::degreeDecrease(int const & i, int const dir)
1961{
1962 for (size_t level=0;level<m_bases.size();++level)
1963 m_bases[level]->degreeDecrease(i,dir);
1964
1965 for(unsigned c=0; c<d; ++c)
1966 m_deg[c]=m_bases[0]->degree(c);
1967
1968 this->update_structure();
1969}
1970
1971template<short_t d, class T>
1972bool gsHTensorBasis<d,T>::testPartitionOfUnity(const index_t npts, const T tol) const
1973{
1974 gsVector<unsigned> np(d);
1975 np.setConstant(npts);
1976 gsMatrix<T> supp = this->support();
1977 gsMatrix<T> grid = gsPointGrid<T>(supp.col(0),supp.col(1),np);
1978 gsMatrix<T> res;
1979 this->eval_into(grid,res);
1980 gsVector<T> sums = res.colwise().sum();
1981 return ((sums.array()>1-tol && sums.array()<1+tol).all());
1982}
1983
1984
1985namespace internal
1986{
1987
1989template<short_t d, class T>
1990class gsXml< gsHTensorBasis<d,T> >
1991{
1992private:
1993 gsXml() { }
1994public:
1995 GSXML_COMMON_FUNCTIONS(gsHTensorBasis<TMPLA2(d,T)>);
1996 static std::string tag () { return "Basis"; }
1997 static std::string type () { return ""; } // tag ?
1998
1999 static gsHTensorBasis<d,T> * get (gsXmlNode * node)
2000 {
2001 gsXmlAttribute * btype = node->first_attribute("type");
2002 if ( ! btype )
2003 {
2004 gsWarn<< "Basis without a type in the xml file.\n";
2005 return NULL;
2006 }
2007 std::string s = btype->value() ;
2008 if ( s.compare(0, 9, "HBSplineB" , 9 ) == 0 ) // needs correct d as well
2009 return gsXml< gsHBSplineBasis<d,T> >::get(node);
2010 if ( s.compare(0, 10,"THBSplineB", 10) == 0 )
2011 return gsXml< gsTHBSplineBasis<d,T> >::get(node);
2012
2013 gsWarn<<"gsXmlUtils: gsHTensorBasis: No known basis \""<<s<<"\". Error.\n";
2014 return NULL;
2015 }
2016
2017 static gsXmlNode * put (const gsHTensorBasis<d,T> & obj,
2018 gsXmlTree & data )
2019 {
2020 const gsBasis<T> * ptr = & obj;
2021
2022 // Hier. B-splines
2023 if ( const gsHBSplineBasis<d,T> * g =
2024 dynamic_cast<const gsHBSplineBasis<d,T> *>( ptr ) )
2025 return gsXml< gsHBSplineBasis<d,T> >::put(*g,data);
2026
2027 // Truncated hier. B-splines
2028 if ( const gsTHBSplineBasis<d,T> * g =
2029 dynamic_cast<const gsTHBSplineBasis<d,T> *>( ptr ) )
2030 return gsXml< gsTHBSplineBasis<d,T> >::put(*g,data);
2031
2032 gsWarn<<"gsXmlUtils put: getBasis: No known basis \""<<obj<<"\". Error.\n";
2033 return NULL;
2034 }
2035};
2036
2037} // namespace internal
2038
2039
2040} // 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
Creates a mapped object or data pointer to a const vector without copying data.
Definition gsAsMatrix.h:285
A basis represents a family of scalar basis functions defined over a common parameter domain.
Definition gsBasis.h:79
virtual short_t degree(short_t i) const
Degree with respect to the i-th variable. If the basis is a tensor product of (piecewise) polynomial ...
Definition gsBasis.hpp:699
uPtr clone()
Clone methode. Produceds a deep copy inside a uPtr.
Class with a hierarchical domain structure represented by a box k-d-tree.
Definition gsHDomain.h:76
void insertBox(point const &lower, point const &upper, node *_node, int lvl)
The insert function which insert box defined by points lower and upper to level lvl.
Definition gsHDomain.hpp:143
Class representing a (scalar) hierarchical tensor basis of functions .
Definition gsHTensorBasis.h:75
memory::unique_ptr< gsHTensorBasis > uPtr
Unique pointer for gsHTensorBasis.
Definition gsHTensorBasis.h:81
index_t getLevelAtPoint(const gsMatrix< T > &Pt) const
Returns the level(s) at point(s) in the parameter domain.
Definition gsHTensorBasis.hpp:141
virtual void uniformCoarsen(int numKnots=1)
Coarsen the basis uniformly by removing groups of numKnots consecutive knots, each knot removed mul t...
Definition gsHTensorBasis.hpp:1590
void _knotIndexToDiadicIndex(const index_t level, const index_t dir, index_t &knotIndex) const
Transfers the knotIndex in the knot span in direction dir on level level to diadic indices.
Definition gsHTensorBasis.hpp:1893
virtual void unrefineElements(std::vector< index_t > const &boxes)
Clear the given boxes into the quadtree.
Definition gsHTensorBasis.hpp:865
void numActive_into(const gsMatrix< T > &u, gsVector< index_t > &result) const
The number of active basis functions at points u.
Definition gsHTensorBasis.hpp:187
void refineSide(const boxSide side, index_t lvl)
Refines all the cells on the side side up to level lvl.
Definition gsHTensorBasis.hpp:902
void flatTensorIndexesToHierachicalIndexes(gsSortedVector< int > &indexes, const int level) const
transformes a sortedVector indexes of flat tensor index of the bspline basis of level to hierachical ...
Definition gsHTensorBasis.hpp:1237
gsMatrix< index_t > allBoundary() const
Returns the indices of the basis functions that are nonzero at the domain boundary.
Definition gsHTensorBasis.hpp:1450
void getBoxesAlongSlice(int dir, T par, std::vector< index_t > &boxes) const
Definition gsHTensorBasis.hpp:1867
virtual void refine(gsMatrix< T > const &boxes, int refExt)
Refine the basis to levels and in the areas defined by boxes with an extension.
Definition gsHTensorBasis.hpp:445
std::vector< std::vector< std::vector< index_t > > > domainBoundariesIndices(std::vector< std::vector< std::vector< std::vector< index_t > > > > &result) const
Gives polylines on the boundaries between different levels of the mesh.
Definition gsHTensorBasis.hpp:1619
virtual void connectivity(const gsMatrix< T > &nodes, gsMesh< T > &mesh) const
Definition gsHTensorBasis.hpp:281
void activeBoundaryFunctionsOfLevel(const unsigned level, const boxSide &s, std::vector< bool > &actives) const
Definition gsHTensorBasis.hpp:1286
void uniformRefine_withCoefs(gsMatrix< T > &coefs, int numKnots=1, int mul=1, short_t const dir=-1)
Refine the basis uniformly.
Definition gsHTensorBasis.hpp:380
void makeCompressed()
Cleans the basis, removing any inactive levels.
Definition gsHTensorBasis.hpp:1220
void transfer(const std::vector< gsSortedVector< index_t > > &old, gsSparseMatrix< T > &result)
Returns transfer matrix between the hirarchical spline given by the characteristic matrix "old" and t...
Definition gsHTensorBasis.hpp:1736
void _diadicIndexToKnotIndex(const index_t level, gsVector< index_t, d > &diadicIndex) const
Transfers the diadicIndex in the knot span in direction on level level to knot indices.
Definition gsHTensorBasis.hpp:1916
void addLevel(const gsTensorBSplineBasis< d, T > &next_basis)
Adds a level, only if manual levels are activated.
Definition gsHTensorBasis.hpp:35
index_t size() const
The number of basis functions in this basis.
Definition gsHTensorBasis.hpp:298
bool testPartitionOfUnity(const index_t npts=100, const T tol=1e-12) const
Test the partition of unity.
Definition gsHTensorBasis.hpp:1972
virtual void increaseMultiplicity(index_t lvl, int dir, T knotValue, int mult=1)
Increases the multiplicity of a knot with the value knotValue in level lvl in direction dir by mult....
Definition gsHTensorBasis.hpp:1829
virtual void active_into(const gsMatrix< T > &u, gsMatrix< index_t > &result) const
Returns the indices of active basis functions at points u, as a list of indices, in result....
Definition gsHTensorBasis.hpp:1384
void insert_box(point const &k1, point const &k2, int lvl)
Inserts a domain into the basis.
Definition gsHTensorBasis.hpp:1208
void functionOverlap(const point &boxLow, const point &boxUpp, const int level, point &actLow, point &actUpp)
Returns the basis functions of level which have support on box, represented as an index box.
Definition gsHTensorBasis.hpp:1080
void refineBasisFunction(const index_t i)
Refines the basis function with (hierarchical) index i.
Definition gsHTensorBasis.hpp:812
int flatTensorIndexToHierachicalIndex(index_t index, const int level) const
takes a flat tensor index of the bspline basis of level and gives back the hierachical index....
Definition gsHTensorBasis.hpp:1274
virtual gsMatrix< index_t > boundaryOffset(boxSide const &s, index_t offset) const
Definition gsHTensorBasis.hpp:1471
std::vector< std::vector< std::vector< index_t > > > domainBoundariesGeneric(std::vector< std::vector< std::vector< std::vector< index_t > > > > &indices, std::vector< std::vector< std::vector< std::vector< T > > > > &params, bool indicesFlag) const
Implementation of the features common to domainBoundariesParams and domainBoundariesIndices....
Definition gsHTensorBasis.hpp:1627
index_t getLevelAtIndex(const point &Pt) const
Returns the level(s) at indexes in the parameter domain.
Definition gsHTensorBasis.hpp:160
void uniformCoarsen_withCoefs(gsMatrix< T > &coefs, int numKnots=1)
Coarsen the basis uniformly.
Definition gsHTensorBasis.hpp:419
gsMatrix< T > support() const
Returns the boundary basis for side s.
Definition gsHTensorBasis.hpp:124
void needLevel(int maxLevel) const
Makes sure that there are numLevels grids computed in the hierarachy.
Definition gsHTensorBasis.hpp:1338
void refineElements_withCoefs(gsMatrix< T > &coefs, std::vector< index_t > const &boxes)
Definition gsHTensorBasis.hpp:314
void setActiveToLvl(int level, std::vector< CMatrix > &x_matrix_lvl) const
Creates characteristic matrices for basis where "level" is the maximum level i.e. ignoring higher lev...
Definition gsHTensorBasis.hpp:1154
virtual void update_structure()
Updates the basis structure (eg. charact. matrices, etc), to be called after any modifications.
Definition gsHTensorBasis.hpp:1309
void getLevelUniqueSpanAtPoints(const gsMatrix< T > &Pt, gsVector< index_t > &lvl, gsMatrix< index_t > &loIdx) const
Returns the level(s) and knot span(s) at point(s) in the parameter domain.
Definition gsHTensorBasis.hpp:170
virtual void refineElements(std::vector< index_t > const &boxes)
Insert the given boxes into the quadtree.
Definition gsHTensorBasis.hpp:844
std::vector< std::vector< std::vector< index_t > > > domainBoundariesParams(std::vector< std::vector< std::vector< std::vector< T > > > > &result) const
Gives polylines on the boundaries between different levels of the mesh.
Definition gsHTensorBasis.hpp:1612
Class for representing a knot vector.
Definition gsKnotVector.h:80
size_t size() const
Number of knots (including repetitions).
Definition gsKnotVector.h:242
mult_t maxInteriorMultiplicity() const
Returns the maximum multiplicity in the interior.
Definition gsKnotVector.hpp:724
internal::gsKnotIterator< T > smart_iterator
Definition gsKnotVector.h:107
int degree() const
Returns the degree of the knot vector.
Definition gsKnotVector.hpp:700
uiterator domainUBegin() const
Returns a unique iterator pointing to the starting knot of the domain.
Definition gsKnotVector.h:359
void insert(T knot, mult_t mult=1)
Inserts knot knot into the knot vector with multiplicity mult.
Definition gsKnotVector.hpp:325
uiterator domainUEnd() const
Returns a unique iterator pointing to the ending knot of the domain.
Definition gsKnotVector.h:367
uiterator uFind(const T u) const
Returns the uiterator pointing to the knot at the beginning of the knot interval containing u....
Definition gsKnotVector.hpp:747
void symDifference(const gsKnotVector< T > &other, std::vector< T > &result) const
Computes the symmetric difference between this knot-vector and other.
Definition gsKnotVector.hpp:172
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
Class Representing a triangle mesh with 3D vertices.
Definition gsMesh.h:32
This class is derived from std::vector, and adds sort tracking.
Definition gsSortedVector.h:110
Sparse matrix class, based on gsEigen::SparseMatrix.
Definition gsSparseMatrix.h:139
A tensor product B-spline basis.
Definition gsTensorBSplineBasis.h:37
void refine_withTransfer(gsSparseMatrix< T, RowMajor > &transfer, const std::vector< std::vector< T > > &refineKnots)
Takes a vector of coordinate wise knot values and inserts these values to the basis.
Definition gsTensorBSplineBasis.hpp:66
const Basis_t & component(short_t dir) const
For a tensor product basis, return the (const) 1-d basis for the i-th parameter component.
Definition gsTensorBSplineBasis.h:202
virtual void uniformCoarsen(int numKnots=1)
Coarsen the basis uniformly by removing groups of numKnots consecutive knots, each knot removed mul t...
Definition gsTensorBasis.h:349
unsigned stride(short_t dir) const
Returns the stride for dimension dir.
Definition gsTensorBasis.h:889
virtual void uniformRefine(int numKnots=1, int mul=1, short_t dir=-1)
Refine the basis uniformly by inserting numKnots new knots with multiplicity mul on each knot span.
Definition gsTensorBasis.h:315
unsigned index(unsigned i, unsigned j, unsigned k=0) const
Definition gsTensorBasis.h:882
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
bool nextLexicographic(Vec &cur, const Vec &size)
Iterates through a tensor lattice with the given size. Updates cur and returns true if another entry ...
Definition gsCombinatorics.h:196
bool nextCubePoint(Vec &cur, const Vec &end)
Iterate in lexigographic order through the points of the integer lattice contained in the cube [0,...
Definition gsCombinatorics.h:327
#define short_t
Definition gsConfig.h:35
#define index_t
Definition gsConfig.h:32
#define gsDebug
Definition gsDebug.h:61
#define GISMO_ERROR(message)
Definition gsDebug.h:118
#define gsWarn
Definition gsDebug.h:50
#define GISMO_UNUSED(x)
Definition gsDebug.h:112
#define GISMO_ENSURE(cond, message)
Definition gsDebug.h:102
#define GISMO_ASSERT(cond, message)
Definition gsDebug.h:89
Provides declaration of the Mesh class.
Provides functions to generate structured point data.
Provides implementation of generic XML functions.
Provides declaration of input/output XML utilities struct.
The G+Smo namespace, containing all definitions for the library.
void freeAll(It begin, It end)
Frees all pointers in the range [begin end)
Definition gsMemory.h:312
Struct which represents an interface between two patches.
Definition gsBoundary.h:650
patchSide & first()
first, returns the first patchSide of this interface
Definition gsBoundary.h:776
patchSide & second()
second, returns the second patchSide of this interface
Definition gsBoundary.h:782