G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsBSplineBasis.hpp
Go to the documentation of this file.
1
14#pragma once
15
16#include <gsNurbs/gsBSpline.h>
18
19#include <gsNurbs/gsDeboor.hpp>
20#include <gsNurbs/gsBoehm.h>
22
24
26
27#include <gsIO/gsXml.h>
28
29namespace gismo
30{
31
32
33template <class T>
34typename gsTensorBSplineBasis<1,T>::Self_t *
35gsTensorBSplineBasis<1,T>::New(std::vector<gsBasis<T>*> & bb )
36{
37 GISMO_ASSERT( bb.size() == 1, "Expecting one component");
38 Self_t * c = dynamic_cast<Self_t*>(bb.front());
39 if ( NULL != c )
40 {
41 bb.clear();
42 return c;
43 }
44 else
45 {
46 gsWarn<<"Something went wrong during BSpline construction.\n";
47 return NULL;
48 }
49}
50
51template <class T>
53{
54 return m_knots.uFind(u(0,0)).uIndex();
55}
56
57template <class T>
59{
60 return m_knots.uFind( u ).uIndex();
61}
62
63template <class T>
65{
66 typename KnotVectorType::smart_iterator it = m_knots.sbegin() + j;
67 index_t v = ( (it+m_knots.degree()+1).uIndex() - it.uIndex() ) / 2 ;
68 it.uAdd(v);
69 gsMatrix<T> rvo(1,2);
70 rvo(0,0) = it.value();
71 it.uNext();
72 rvo(0,1) = it.value();
73 // // take indexed support
74 // gsMatrix<index_t,1,2> isup(1,2);
75 // this->elementSupport_into(i, isup);
76 // // take the middle element
77 // isup.at(0) = (isup.at(0)+isup.at(1))/2;
78 // isup.at(1) = isup.at(0)+1;
79 // // take coordinates
80 // gsMatrix<T> rvo(1,2);
81 // rvo(0,0) = m_knots.uValue(isup.at(0));
82 // rvo(0,1) = m_knots.uValue(isup.at(1));
83 return rvo;
84}
85
86template <class T>
88 gsMesh<T> & mesh) const
89{
90 const index_t sz = size();
91 GISMO_ASSERT( nodes.rows() == sz, "Invalid input.");
92
93 // Add first vertex
94 if ( sz != 0 )
95 mesh.addVertex( nodes.row(0).transpose() );
96
97 // Add adges and vertices
98 for(index_t i = 1; i< sz; ++i )
99 {
100 mesh.addVertex( nodes.row(i).transpose() );
101 mesh.addEdge( i-1, i );
102 }
103
104 if ( m_periodic )
105 mesh.addEdge( sz-1, 0);
106}
107
108template <class T>
110 const gsBasis<T> & other,
111 gsMatrix<index_t> & bndThis,
112 gsMatrix<index_t> & bndOther) const
113{
114 if ( const TensorSelf_t * _other = dynamic_cast<const TensorSelf_t*>(&other) )
115 {
116 bndThis .resize(1,1);
117 bndOther.resize(1,1);
118 bndThis (0,0) = bi.first() .side() == 1 ? 0 : m_knots.size() -m_p-2;
119 bndOther(0,0) = bi.second().side() == 1 ? 0 : _other->m_knots.size()-_other->m_p-2;
120 return;
121 }
122
123 gsWarn<< "Cannot match with "<< other <<"\n";
124}
125
126template <class T>
128 gsMatrix<index_t>& result ) const
129{
130 result.resize(m_p+1, u.cols());
131
132 if ( m_periodic )
133 {
134 // We want to keep the non-periodic case unaffected wrt
135 // complexity, therefore we keep the modulo operation of the
136 // periodic case separate
137 const index_t s = size();
138 for (index_t j = 0; j < u.cols(); ++j)
139 {
140 unsigned first = firstActive(u(0,j));
141 for (int i = 0; i != m_p+1; ++i)
142 result(i,j) = (first++) % s;
143 }
144 }
145 else
146 {
147 for (index_t j = 0; j < u.cols(); ++j)
148 {
149 unsigned first = firstActive(u(0,j));
150 for (int i = 0; i != m_p+1; ++i)
151 result(i,j) = first++;
152 }
153 }
154}
155
156template <class T>
158{
159 GISMO_ASSERT( u.rows() == 1, "Invalid input.");
160 // Note: right end of the support will be considered active
161 return( (u.value() >= m_knots[i]) && (u.value() <= m_knots[i+m_p+1]) );
162}
163
164template <class T>
166{
167 if( m_periodic ) // Periodic basis does not have such things as boundaries.
168 {
169 gsWarn << "Periodic basis does not have such things as boundaries.\n";
170 // return NULL;
171 gsMatrix<index_t> matrix;
172 return matrix;
173 }
174 else
175 {
176 gsMatrix<index_t> res(2,1);
177 res(0,0) = 0;
178 res(1,0) = m_knots.size()-m_p-2;
179 return res;
180 }
181}
182
183
184template <class T>
186 index_t offset ) const
187{
188 if( m_periodic )
189 gsWarn << "Periodic basis does not have such things as boundaries.\n";
190
191 gsMatrix<index_t> res(1,1);
192 GISMO_ASSERT(offset+m_p+1 < static_cast<index_t>(m_knots.size()),
193 "Offset cannot be bigger than the amount of basis functions orthogonal to Boxside s!");
194 switch (s) {
195 case boundary::left : // left
196 res(0,0)= offset;
197 break;
198 case boundary::right : // right
199 res(0,0)= m_knots.size()-m_p-2-offset;
200 break;
201 default:
202 GISMO_ERROR("gsBSplineBasis: valid sides is left(west) and right(east).");
203 };
204 return res;
205}
206
208template <class T>
210{
211 return new gsConstantBasis<T>(1.0);
212}
214
215template <class T>
217{
218 // The support of a the whole B-spline basis is the interval
219 // between m_knots[m_p] and m_knots[i+m_p+1].
220 // First m_p and last m_p knots are the "ghost" knots and are not
221 // part of the domain
222 gsMatrix<T> res(1,2);
223 res << domainStart() , domainEnd() ;
224 return res ;
225}
226
227template <class T>
229{
230 // Note: in the periodic case last index is
231 // m_knots.size() - m_p - 1 - m_periodic
232 // but we may still accept the original index of the basis function.
233 // One issue is that in the periodic case the basis function
234 // support has two connected components, so probably one should
235 // call this function with twice for the two twins
236
237 GISMO_ASSERT( static_cast<size_t>(i) < m_knots.size()-m_p-1,
238 "Invalid index of basis function." );
239 gsMatrix<T> res(1,2);
240 res << ( i > m_p ? m_knots[i] : m_knots[m_p] ),
241 ( static_cast<size_t>(i) < (m_knots.size()-2*m_p-2) ? m_knots[i+m_p+1] :
242 m_knots[m_knots.size()-m_p-1] );
243 return res ;
244}
245
246template <class T>
248{
249 if( m_periodic == 0 )
250 return i;
251 const index_t s = size();
252 if ( i < static_cast<index_t>(m_periodic) )
253 i += s;
254 else if ( i > s )
255 i -= s;
256 return i;
257}
258
259template <class T>
261{
262 result.resize(m_p+1, u.cols() );
263
264#if (FALSE)
265
266 typename KnotVectorType::const_iterator kspan;
267
268 // Low degree specializations
269 switch (m_p)
270 {
271 case 0: // constant B-splines
272 result.setOnes();
273 return;
274
275 case 1:
276 for (index_t v = 0; v < u.cols(); ++v) // for all columns of u
277 if ( ! inDomain( u(0,v) ) )
278 result.col(v).setZero();
279 else
280 {
281 kspan = m_knots.findspanIter ( u(0,v) );
282 bspline::evalDeg1Basis( u(0,v), kspan, m_p, result.col(v) );
283 }
284 return;
285 case 2:
286 for (index_t v = 0; v < u.cols(); ++v) // for all columns of u
287 if ( ! inDomain( u(0,v) ) )
288 result.col(v).setZero();
289 else
290 {
291 kspan = m_knots.findspanIter ( u(0,v) );
292 bspline::evalDeg2Basis( u(0,v), kspan, m_p, result.col(v) );
293 }
294 return;
295 case 3:
296 for (index_t v = 0; v < u.cols(); ++v) // for all columns of u
297 if ( ! inDomain( u(0,v) ) )
298 result.col(v).setZero();
299 else
300 {
301 kspan = m_knots.findspanIter ( u(0,v) );
302 bspline::evalDeg3Basis( u(0,v), kspan, m_p, result.col(v) );
303 }
304 return;
305 };
306
307 STACK_ARRAY(T, left, m_p + 1);
308 STACK_ARRAY(T, right, m_p + 1);
309
310 for (index_t v = 0; v < u.cols(); ++v) // for all columns of u
311 {
312 // Check if the point is in the domain
313 if ( ! inDomain( u(0,v) ) )
314 {
315 // gsWarn<< "Point "<< u(0,v) <<" not in the BSpline domain.\n";
316 result.col(v).setZero();
317 }
318 else
319 { // Locate the point in the knot-vector
320 kspan = m_knots.findspanIter ( u(0,v) );
321 // Run evaluation algorithm
322 bspline::evalBasis( u(0,v), kspan, m_p, left, right, result.col(v) );
323 }
324 }
325
326#endif
327//#if (FALSE)
328 STACK_ARRAY(T, left, m_p + 1);
329 STACK_ARRAY(T, right, m_p + 1);
330
331 for (index_t v = 0; v < u.cols(); ++v) // for all columns of u
332 {
333 // Check if the point is in the domain
334 if ( ! inDomain( u(0,v) ) )
335 {
336 // gsWarn<< "Point "<< u(0,v) <<" not in the BSpline domain.\n";
337 result.col(v).setZero();
338 continue;
339 }
340
341 // Run evaluation algorithm
342
343 // Get span of absissae
344 unsigned span = m_knots.iFind( u(0,v) ) - m_knots.begin() ;
345
346 //ndu[0] = (T)(1); // 0-th degree function value
347 result(0,v)= (T)(1); // 0-th degree function value
348
349 for(int j=1; j<= m_p; ++j) // For all degrees ( ndu column)
350 {
351 left[j] = u(0,v) - m_knots[span+1-j];
352 right[j] = m_knots[span+j] - u(0,v);
353 T saved = (T)(0) ;
354
355 for(int r=0; r!=j ; ++r) // For all (except the last) basis functions of degree j ( ndu row)
356 {
357 // Strictly lower triangular part: Knot differences of distance j
358 //ndu[j*p1 + r] = right[r+1]+left[j-r] ;
359
360 //const T temp = ndu[r*p1 + j-1] / ndu[j*p1 + r] ;
361 const T temp = result(r,v) / ( right[r+1] + left[j-r] );
362 // Upper triangular part: Basis functions of degree j
363 //ndu[r*p1 + j] = saved + right[r+1] * temp ;// r-th function value of degree j
364 result(r,v) = saved + right[r+1] * temp ;// r-th function value of degree j
365 saved = left[j-r] * temp ;
366 }
367 //ndu[j*p1 + j] = saved ;// Diagonal: j-th (last) function value of degree j
368 result(j,v) = saved;
369 }
370
371 }// end for all columns v
372//#endif
373}
374
375
376template <class T>
378 const gsMatrix<T> & u,
379 gsMatrix<T>& result) const
380{
381 GISMO_ASSERT( static_cast<size_t>(i) < m_knots.size()-m_p-1,"Invalid index of basis function." );
382
383 result.resize(1, u.cols() );
384 STACK_ARRAY(T, N, m_p + 1);
385
386 for (index_t s=0;s<u.cols(); ++s)
387 {
388 //Special cases
389
390 // Periodic basis ?
391 // Note: for periodic, i is understood modulo the basis size
392 if( m_periodic )
393 {
394 // Basis is periodic and i could be among the crossing functions.
395 if( static_cast<int>(i) <= m_periodic )
396 {
397 gsMatrix<T> supp = support(i);
398 if( u(0,s) < supp(0) || u(0,s) > supp(1)) // u is outside the support of B_i
399 i += size(); // we use the basis function ignored ny the periodic construction
400 }
401 }
402
403 // Special case of C^{-1} on right end of support
404 if ( (static_cast<size_t>(i) == m_knots.size()-m_p-2) &&
405 (u(0,s) == m_knots.last()) && (u(0,s)== m_knots[m_knots.size()-m_p-1]) )
406 {
407 result(0,s)= (T)(1.0);
408 continue;
409 }
410
411 // Locality property
412 if ( (u(0,s) < m_knots[i]) || (u(0,s) >= m_knots[i+m_p+1]) )
413 {
414 result(0,s)= (T)(0.0);
415 continue;
416 }
417
418 //bspline::deBoorTriangle( u(0,s), m_knots.begin() + i, m_p, N );
419 //result(0,s) = N[ m_p ];
420
421 // Initialize zeroth degree functions
422 for (int j=0;j<=m_p; ++j)
423 if ( u(0,s) >= m_knots[i+j] && u(0,s) < m_knots[i+j+1] )
424 N[j] = (T)(1.0);
425 else
426 N[j] = (T)(0.0);
427 // Compute according to the trangular table
428 for (int k=1;k<=m_p; ++k)
429 {
430 T saved;
431 if (N[0] == (T)(0.0))
432 saved = (T)(0.0);
433 else
434 saved= ((u(0,s) - m_knots[i] )* N[0]) / (m_knots[k+i] - m_knots[i]);
435 for (int j=0;j<m_p-k+1; ++j)
436 {
437 const T kleft = m_knots[i+j+1];
438 const T kright = m_knots[i+j+k+1];
439 if ( N[j+1] == (T)(0.0) )
440 {
441 N[j] = saved;
442 saved = (T)(0.0);
443 }
444 else
445 {
446 const T temp = N[j+1] / ( kright - kleft );
447 N[j]= saved + ( kright-u(0,s) )*temp;
448 saved= (u(0,s) -kleft ) * temp;
449 }
450 }
451 }
452 result(0,s) = N[0];
453 }
454}
455
456template <class T>
458 const gsMatrix<T> & u,
459 int n,
460 gsMatrix<T>& result) const
461{
462 GISMO_ASSERT( u.rows() == 1 , "gsBSplineBasis accepts points with one coordinate.");
463 // This is a copy of the interior of evalAllDerSingle_into() except for the last for cycle with k.
464
465 // Notation from the NURBS book:
466 // p - degree
467 // U = {u[0],...,u[m]} the knot vector
468 // i - id of the basis function (hopefully the same ;)
469 // u - parameter (where to evaluate)
470
471 // Note that N[i][k] = N[ i*table_size + k ].
472
473 T saved;
474 T Uleft;
475 T Uright;
476 T temp;
477 const unsigned int table_size = m_p + 1;
478 STACK_ARRAY(T, N, table_size*table_size);
479 STACK_ARRAY(T, ND, table_size);
480
481 result.resize(1, u.cols()); // (1 derivative) times number of points of evaluation.
482
483 // From here on, we follow the book.
484 for( index_t s = 0; s < u.cols(); s++ )
485 {
486 // Special cases
487 // Periodic basis
488 if( (int)(i) <= m_periodic ) // Basis is periodic and i could be among the crossing functions.
489 {
490 gsMatrix<T> supp = support(i);
491 if( u(0,s) < supp(0) || u(0,s) > supp(1)) // u is outside the support of B_i
492 i += size(); // we use the basis function ignored ny the periodic construction
493 }
494
495 if( u(0,s) == m_knots.last() )
496 gsWarn << "evalDerSingle_into sometimes gives strange results for the last knot.\n";
497
498 if( u(0,s) < m_knots[i] || u(0,s) >= m_knots[i+m_p+1]) // Local property
499 {
500 // If we are outside the support, the value and all the
501 // derivatives there are equal to zero.
502 result(0,s ) = 0;
503 continue;// to the next point
504 }
505
506 bspline::deBoorTriangle( u(0,s), m_knots.begin() + i, m_p, N );
507
508 /*
509 for( int j = 0; j <= m_p; j++ ) // Initialize zeroth-degree functions
510 if( u(0,s) >= m_knots[i+j] && u(0,s) < m_knots[i+j+1] )
511 N[ j*table_size ] = 1;
512 else
513 N[ j*table_size ] = 0;
514 for( int k = 1; k <= m_p; k++ ) // Compute full triangular table
515 {
516 if( N[(k-1)] == 0 )
517 saved = 0;
518 else
519 saved = ((u(0,s)-m_knots[i])*N[ k-1 ])/(m_knots[i+k]-m_knots[i]);
520 for( int j = 0; j < m_p-k+1; j++)
521 {
522 Uleft = m_knots[i+j+1];
523 Uright = m_knots[i+j+k+1];
524 if( N[ (j+1)*table_size + (k-1) ] == 0 )
525 {
526 N[ j*table_size + k ] = saved;
527 saved = 0;
528 }
529 else
530 {
531 temp = N[ (j+1)*table_size + (k-1) ]/(Uright-Uleft);
532 N[ j*table_size + k ] = saved + (Uright-u(0,s))*temp;
533 saved = (u(0,s)-Uleft)*temp;
534 }
535 }
536 }*/
537
538 // Not tested few lines:
539 if( n == 0 )
540 {
541 result(0,s) = N[ m_p ]; // The function value
542 continue;
543 }
544
545 for( int j = 0; j <= n; j++ ) // Load appropriate column
546 ND[j] = N[ j*table_size + (m_p-n) ];
547 for( int jj = 1; jj <= n; jj++ ) // Compute table of width n
548 {
549 if( ND[0] == 0 )
550 saved = 0;
551 else
552 saved = ND[0]/(m_knots[ i+m_p-n+jj ] - m_knots[i]);
553 for( int j = 0; j < n-jj+1; j++ )
554 {
555 Uleft = m_knots[i+j+1];
556 Uright = m_knots[i+j+m_p-n+jj+1];
557 if( ND[j+1] == 0 )
558 {
559 ND[j] = static_cast<T>(m_p-n+jj)*saved;
560 saved = 0;
561 }
562 else
563 {
564 temp = ND[j+1]/(Uright - Uleft);
565 ND[j] = static_cast<T>(m_p-n+jj)*(saved - temp);
566 saved = temp;
567 }
568 }
569 }
570 result(0,s) = ND[0]; // n-th derivative
571 }
572 // Changes from the book: notation change, stack arrays instead of normal C arrays, changed all 0.0s and 1.0s into 0 and 1. No return in the first check
573 // Warning: Segmentation fault if n > m_p.
574
575}
576
577template <class T> inline
579 const gsMatrix<T> & coefs,
580 gsMatrix<T>& result) const
581{
582 GISMO_ASSERT( u.rows() == 1 , "gsBSplineBasis accepts points with one coordinate (got "
583 <<u.rows()<<").");
584 if( m_periodic == 0 )
585 gsDeboor(u, m_knots, m_p, coefs, result); // De Boor's algorithm
586 else
587 {
588 gsDeboor(u, m_knots, m_p, perCoefs(coefs), result); // Make sure that the ghost coefficients agree with the beginning of the curve.
589 }
590 // return gsBasis<T>::eval(u,coefs); // defalult gsBasis implementation
591}
592
593
594template <class T> inline
596{
597 GISMO_ASSERT( u.rows() == 1 , "gsBSplineBasis accepts points with one coordinate.");
598
599 const int pk = m_p-1 ;
600 const int p1 = m_p + 1; // degree plus one
601 STACK_ARRAY(T, ndu , m_p);
602 STACK_ARRAY(T, left , p1 );
603 STACK_ARRAY(T, right, p1 );
604
605 result.resize( m_p + 1, u.cols() ) ;
606
607 for (index_t v = 0; v < u.cols(); ++v) // for all columns of u
608 {
609 // Check if the point is in the domain
610 if ( ! inDomain( u(0,v) ) )
611 {
612 // gsWarn<< "Point "<< u(0,v) <<" not in the BSpline domain.\n";
613 result.col(v).setZero();
614 continue;
615 }
616
617 // Run evaluation algorithm and keep first derivative
618
619 // Get span of absissae
620 typename KnotVectorType::iterator span = m_knots.iFind( u(0,v) );
621
622 ndu[0] = (T)(1); // 0-th degree function value
623 left[0] = 0;
624
625 for(int j=1; j<m_p ; j++) // For all degrees
626 {
627 // Compute knot splits
628 left[j] = u(0,v) - *(span+1-j);
629 right[j] = *(span+j) - u(0,v);
630
631 // Compute Basis functions of degree m_p-1 ( ndu[] )
632 T saved = (T)(0) ;
633 for(int r=0; r<j ; r++) // For all (except the last) basis functions of degree 1..j
634 {
635 const T temp = ndu[r] / ( right[r+1] + left[j-r] ) ;
636 ndu[r] = saved + right[r+1] * temp ;// r-th function value of degree 1..j
637 saved = left[j-r] * temp ;
638 }
639 ndu[j] = saved ;// last function value of degree 1..j
640 }
641
642 // Compute last knot split
643 left[m_p] = u(0,v) - *(span+1-m_p);
644 right[m_p] = *(span+m_p) - u(0,v);
645
646 // Compute the first derivatives (using ndu[] and left+right)
647 right[0] = right[1]+left[m_p] ;
648 result(0 , v) = - static_cast<T>(m_p) * ndu[0] / right[0] ;
649 for(int r = 1; r < m_p; r++)
650 {
651 // Compute knot difference r of distance m_p (overwrite right[])
652 right[r] = right[r+1]+left[m_p-r] ;
653 result(r, v) = static_cast<T>(m_p) * ( ndu[r-1] / right[r-1] - ndu[r] / right[r] );
654 }
655 result(m_p, v) = static_cast<T>(m_p) * ndu[pk] / right[pk];
656
657 }// end for all columns v
658}
659
660template <class T> inline
662 gsMatrix<T>& result ) const
663{
664 std::vector<gsMatrix<T> > ev;
665 this->evalAllDers_into(u, 2, ev);
666 result.swap(ev[2]);
667}
668
669template <class T> inline
671 const gsMatrix<T> & u,
672 gsMatrix<T>& result ) const
673{
674 // \todo Redo an efficient implementation p. 76, Alg. A2.5 Nurbs book
675 result.resize(1, u.cols() );
676 gsMatrix<T> tmp;
678
679 for (index_t j = 0; j < u.cols(); ++j)
680 {
681 const index_t first = firstActive(u(0,j));
682 if ( (i>= first) && (i<= first + m_p) )
683 result(0,j) = tmp(i-first,j);
684 else
685 result(0,j) = (T)(0.0);
686 }
687}
688
689template <class T> inline
690void
692 const gsMatrix<T> & u,
693 int n,
694 gsMatrix<T>& result) const
695{
696 GISMO_ASSERT( u.rows() == 1 , "gsBSplineBasis accepts points with one coordinate.");
697 // Notation from the NURBS book:
698 // p - degree
699 // U = {u[0],...,u[m]} the knot vector
700 // i - id of the basis function (hopefully the same ;)
701 // u - parameter (where to evaluate)
702
703 T saved, Uleft, Uright, temp;
704 const unsigned int p1 = m_p + 1;
705 // Note that N[i][k] = N[ i*p1 + k ].
706 STACK_ARRAY(T, N, p1*p1);
707 STACK_ARRAY(T, ND, p1);
708
709 result.setZero(n + 1, u.cols());
710 n = ( n>m_p ? m_p : n );
711
712 for( index_t s = 0; s < u.cols(); s++ )
713 {
714 // Periodic basis
715 if( (int)(i) <= m_periodic ) // Basis is periodic and i could be among the crossing functions.
716 {
717 gsMatrix<T> supp = support(i);
718 if( u(0,s) < supp(0) || u(0,s) > supp(1)) // u is outside the support of B_i
719 i += size(); // we use the basis function ignored ny the periodic construction
720 }
721
722 if( u(0,s) < m_knots[i] || u(0,s) >= m_knots[i+m_p+1]) // Local property
723 if( u( 0,s ) != m_knots[m_knots.size()-m_p-1] )
724 continue; //zeros
725
726 //bspline::deBoorTriangle( u(0,s), m_knots.begin() + i, m_p, N );
727 if ( u(0,s) == m_knots[m_knots.size()-m_p-1] )
728 { // Initialize zeroth-degree functions
729 for( int j = 0; j <= m_p; j++ )
730 N[ j*p1 ] = (T)( u(0,s) > m_knots[i+j] && u(0,s) <= m_knots[i+j+1] );
731 }
732 else // not at right boundary
733 {
734 for( int j = 0; j <= m_p; j++ )
735 N[ j*p1 ] = (T)( u(0,s) >= m_knots[i+j] && u(0,s) < m_knots[i+j+1] );
736 }
737
738 for( int k = 1; k <= m_p; k++ ) // Compute full triangular table
739 {
740 saved = (N[(k-1)] == 0 ? (T)(0) :
741 ((u(0,s)-m_knots[i])*N[ k-1 ])/(m_knots[i+k]-m_knots[i]) );
742
743 for( int j = 0; j < m_p-k+1; j++)
744 {
745 Uleft = m_knots[i+j+1];
746 Uright = m_knots[i+j+k+1];
747 if( N[ (j+1)*p1 + (k-1) ] == 0 )
748 {
749 N[ j*p1 + k ] = saved;
750 saved = 0;
751 }
752 else
753 {
754 temp = N[ (j+1)*p1 + (k-1) ]/(Uright-Uleft);
755 N[ j*p1 + k ] = saved + (Uright-u(0,s))*temp;
756 saved = (u(0,s)-Uleft)*temp;
757 }
758 }
759 }
760
761 result(0,s) = N[ m_p ]; // The function value
762
763 for( int k = 1; k <= n; k++ ) // Compute the derivatives
764 {
765 for( int j = 0; j <= k; j++ ) // Load appropriate column
766 ND[j] = N[ j*p1 + (m_p-k) ];
767 for( int jj = 1; jj <= k; jj++ ) // Compute table of width k
768 {
769 if( ND[0] == 0 )
770 saved = 0;
771 else
772 saved = ND[0]/(m_knots[ i+m_p-k+jj ] - m_knots[i]);
773 for( int j = 0; j < k-jj+1; j++ )
774 {
775 Uleft = m_knots[i+j+1];
776 Uright = m_knots[i + m_p-k+jj + j+1];
777 if( ND[j+1] == 0 )
778 {
779 ND[j] = static_cast<T>(m_p-k+jj)*saved;
780 saved = 0;
781 }
782 else
783 {
784 temp = ND[j+1]/(Uright - Uleft);
785 ND[j] = static_cast<T>(m_p-k+jj)*(saved - temp);
786 saved = temp;
787 }
788 }
789 }
790 result(k,s) = ND[0]; // k-th derivative
791 }
792 }
793}
794
795template <class T> inline
796void gsTensorBSplineBasis<1,T>::deriv_into(const gsMatrix<T> & u, const gsMatrix<T> & coefs, gsMatrix<T>& result ) const
797{
798 // TO DO specialized computation for gsBSplineBasis
799 if( m_periodic == 0 )
800 gsBasis<T>::derivFunc_into(u, coefs, result);
801 else
802 gsBasis<T>::derivFunc_into(u, perCoefs(coefs), result);
803}
804
805template <class T> inline
806void gsTensorBSplineBasis<1,T>::deriv2_into(const gsMatrix<T> & u, const gsMatrix<T> & coefs, gsMatrix<T>& result ) const
807{
808 // TO DO specialized computation for gsBSplineBasis
809 if( m_periodic == 0 )
810 gsBasis<T>::deriv2Func_into(u, coefs, result);
811 else
812 gsBasis<T>::deriv2Func_into(u, perCoefs(coefs), result);
813}
814
815template <class T> inline
817{
818 // \todo Redo an efficient implementation p. 76, Alg. A2.5 Nurbs book
819 result.resize(1, u.cols() );
820 gsMatrix<T> tmp;
822
823 for (index_t j = 0; j < u.cols(); ++j)
824 {
825 const index_t first = firstActive(u(0,j));
826 if ( (i>= first) && (i<= first + m_p) )
827 result(0,j) = tmp(i-first,j);
828 else
829 result(0,j) = (T)(0.0);
830 }
831}
832
833template <class T> inline
835{
836 std::vector<gsMatrix<T> > ev;
837 this->evalAllDers_into(u, 2, ev);
838 return ev[2].colwise().sum();
839}
840
841template <class T>
842typename gsBasis<T>::uPtr
844{
845 typename Self_t::uPtr ptr1 = memory::convert_ptr<Self_t>(other.clone());
846
847 if ( ptr1 )
848 return typename gsBasis<T>::uPtr(new gsTensorBSplineBasis<2,T>( ptr1.release(), this->clone().release() ));
849 else
850 {
851 gsWarn<<"gsTensorBSplineBasis::tensorize: Invalid basis "<< other <<"\n";
852 return typename gsBasis<T>::uPtr();
853 }
854}
855
856template <class T>
857memory::unique_ptr<gsGeometry<T> > gsBSplineBasis<T>::makeGeometry( gsMatrix<T> coefs ) const
858{
859 return typename gsGeometry<T>::uPtr(new GeometryType(*this, give(coefs)));
860}
861
862template <class T>
864evalAllDers_into(const gsMatrix<T> & u, int n,
865 std::vector<gsMatrix<T> >& result,
866 bool sameElement) const
867{
868 // TO DO : Use less memory proportionally to n
869 // Only last n+1 columns and last n rows of ndu are needed
870 // Also a's size is proportional to n
871 GISMO_ASSERT( u.rows() == 1 , "gsBSplineBasis accepts points with one coordinate.");
872
873 const int p1 = m_p + 1; // degree plus one
874
875 STACK_ARRAY(T, ndu, p1 * p1 );
876 STACK_ARRAY(T, left, p1);
877 STACK_ARRAY(T, right, p1);
878 STACK_ARRAY(T, a, 2 * p1);
879
880 result.resize(n+1);
881 for(int k=0; k<=n; k++)
882 result[k].resize(m_p + 1, u.cols());
883
884#if FALSE
885
886 const int pn = m_p - n;
887 unsigned span;
888
889 for (index_t v = 0; v < u.cols(); ++v) // for all columns of u
890 {
891 if (!sameElement || 0==v)
892 {
893
894 // Check if the point is in the domain
895 if ( ! inDomain( u(0,v) ) )
896 {
897 // gsWarn<< "Point "<< u(0,v) <<" not in the BSpline domain.\n";
898 for(int k=0; k<=n; k++)
899 result[k].col(v).setZero();
900 continue;
901 }
902
903 // Run evaluation algorithm and keep the function values triangle & the knot differences
904 span = m_knots.findspan( u(0,v) ) ; // Get span of absissae.
905 }
906
907 for(int j=1; j<= m_p; j++) // For all degrees ( ndu column)
908 {
909 // Compute knot splits
910 left[j] = u(0,v) - m_knots[span+1-j];
911 right[j] = m_knots[span+j] - u(0,v);
912 }
913
914 ndu[0] = (T)(1) ; // 0-th degree function value
915 T saved = (T)(0) ;
916 // Compute Basis functions of degree m_p-n ( ndu[] )
917 for(int r=0; r<pn ; r++) // For all (except the last) basis functions of degree m_p-n
918 {
919 const T temp = ndu[r*p1 + pn -1] / ( right[r+1] + left[pn-r] ) ;
920 ndu[r*p1 + pn] = saved + right[r+1] * temp ;// r-th function value of degree j
921 saved = left[pn-r] * temp ;
922 }
923 ndu[pn*p1 + pn] = saved ;
924
925 // Compute n-th derivative and continue to n-1 ... 0
926
927 for(int r=0; r<pn ; r++) // For all (except the last) basis functions of degree j ( ndu row)
928 {
929 // Strictly lower triangular part: Knot differences of distance j
930 ndu[j*p1 + r] = right[r+1]+left[j-r] ;
931 const T temp = ndu[r*p1 + j-1] / ndu[j*p1 + r] ;
932 // Upper triangular part: Basis functions of degree j
933 ndu[r*p1 + j] = saved + right[r+1] * temp ;// r-th function value of degree j
934 saved = left[j-r] * temp ;
935 }
936 // Diagonal: j-th (last) function value of degree j
937 ndu[j*p1 + j] = saved ;
938 }
939
940#endif
941
942 typename KnotVectorType::iterator span;
943
944 for (index_t v = 0; v < u.cols(); ++v) // for all columns of u
945 {
946
947 if (!sameElement || 0==v)
948 {
949 // Check if the point is in the domain
950 if ( ! inDomain( u(0,v) ) )
951 {
952 //gsDebug<< "Point "<< u(0,v) <<" not in the BSpline domain ["
953 // << *(m_knots.begin()+m_p)<< ", "<<*(m_knots.end()-m_p-1)<<"].\n";
954 for(int k=0; k<=n; k++)
955 result[k].col(v).setZero();
956 continue;
957 }
958
959 // Run evaluation algorithm and keep the function values triangle & the knot differences
960 span = m_knots.iFind( u(0,v) );
961 }
962
963 ndu[0] = (T)(1) ; // 0-th degree function value
964 for(int j=1; j<= m_p; j++) // For all degrees ( ndu column)
965 {
966 // Compute knot splits
967 left[j] = u(0,v) - *(span+1-j);
968 right[j] = *(span+j) - u(0,v);
969
970 T saved = (T)(0) ;
971
972 for(int r=0; r<j ; r++) // For all (except the last) basis functions of degree j ( ndu row)
973 {
974 // Strictly lower triangular part: Knot differences of distance j
975 ndu[j*p1 + r] = right[r+1]+left[j-r] ;
976 const T temp = ndu[r*p1 + j-1] / ndu[j*p1 + r] ;
977 // Upper triangular part: Basis functions of degree j
978 ndu[r*p1 + j] = saved + right[r+1] * temp ;// r-th function value of degree j
979 saved = left[j-r] * temp ;
980 }
981 // Diagonal: j-th (last) function value of degree j
982 ndu[j*p1 + j] = saved ;
983 }
984
985 // Assign 0-derivative equal to function values
986 //result.front().block(0,v, p1,1) = ndu.col(m_p);
987 for (int j=0; j <= m_p ; ++j )
988 result.front()(j,v) = ndu[j*p1 + m_p];
989
990 // Compute the derivatives
991 for(int r = 0; r <= m_p; r++)
992 {
993 // alternate rows in array a
994 T* a1 = &a[0];
995 T* a2 = &a[p1];
996
997 a1[0] = (T)(1) ;
998
999 // Compute the k-th derivative of the r-th basis function
1000 for(int k=1; k<=n; k++)
1001 {
1002 int rk,pk,j1,j2 ;
1003 T d(0) ;
1004 rk = r-k ; pk = m_p-k ;
1005
1006 if(r >= k)
1007 {
1008 a2[0] = a1[0] / ndu[ (pk+1)*p1 + rk] ;
1009 d = a2[0] * ndu[rk*p1 + pk] ;
1010 }
1011
1012 j1 = ( rk >= -1 ? 1 : -rk );
1013 j2 = ( r-1 <= pk ? k-1 : m_p - r );
1014
1015 for(int j = j1; j <= j2; j++)
1016 {
1017 a2[j] = (a1[j] - a1[j-1]) / ndu[(pk+1)*p1 + rk+j] ;
1018 d += a2[j] * ndu[(rk+j)*p1 + pk] ;
1019 }
1020
1021 if(r <= pk)
1022 {
1023 a2[k] = -a1[k-1] / ndu[(pk+1)*p1 + r] ;
1024 d += a2[k] * ndu[r*p1 + pk] ;
1025 }
1026
1027 result[k](r, v) = d;
1028
1029 std::swap(a1, a2); // Switch rows
1030 }
1031 }
1032 }// end for all columns v
1033
1034 // Multiply through by the factor factorial(m_p)/factorial(m_p-k)
1035 int r = m_p ;
1036 for(int k=1; k<=n; k++)
1037 {
1038 result[k].array() *= (T)(r) ;
1039 r *= m_p - k ;
1040 }
1041}
1042
1043
1044template <class T>
1045void gsTensorBSplineBasis<1,T>::refine_withCoefs(gsMatrix<T>& coefs, const std::vector<T>& knots)
1046{
1047 // Probably does not work with periodic basis. We would need to shift the knots in the ghost areas by the active length,
1048 // sort the knot vector, swap the coefs accordingly, call this with perCoefs (or coefs?) and stretch outer knots.
1049 gsBoehmRefine(this->knots(), coefs, m_p, knots.begin(), knots.end());
1050}
1051
1052
1053template <class T>
1055{
1056 // See remark about periodic basis in refine_withCoefs, please.
1057 gsFiberMatrix<T> trans;
1058 trans.setIdentity( this->size() );
1059 gsBoehmRefine(this->knots(), trans, m_p, knots.begin(), knots.end());
1060 trans.toSparseMatrix( transfer );
1061}
1062
1063
1064template <class T>
1066{
1067 // See remark about periodic basis in refine_withCoefs, please.
1068 std::vector<T> newKnots;
1069 this->knots().getUniformRefinementKnots(numKnots, newKnots,mul);
1070 this->refine_withCoefs(coefs, newKnots);
1071}
1072
1073
1074template <class T>
1076{
1077 // See remark about periodic basis in refine_withCoefs, please.
1078 std::vector<T> newKnots;
1079 this->knots().getUniformRefinementKnots(numKnots, newKnots,mul);
1080 this->refine_withTransfer(transfer, newKnots);
1081}
1082
1083template <class T>
1085{
1086 // Simple implementation: coarsen and refine again.
1087 // Could be done more efficiently if needed.
1088 std::vector<T> removedKnots = m_knots.coarsen(numKnots);
1089 this->clone()->refine_withTransfer( transfer, removedKnots );
1090}
1091
1092template <class T>
1094{
1095 GISMO_ASSERT(c<3,"Invalid corner for 1D basis.");
1096 return ( c == 1 ? 0 : this->size()-1);
1097}
1098
1099template <class T>
1100std::ostream & gsTensorBSplineBasis<1,T>::print(std::ostream &os) const
1101{
1102 os << "BSplineBasis: deg=" <<degree()
1103 << ", size="<< size() << ", knot vector:\n";
1104 os << this->knots();
1105 if( m_periodic > 0 )
1106 os << ",\n m_periodic= " << m_periodic;
1107 return os;
1108}
1109template <class T>
1111{
1112 std::stringstream os;
1113 os << "BSplineBasis: deg=" <<degree()
1114 << ", size="<< size() << ", knot vector:\n";
1115 os << this->knots().detail();
1116 if( m_periodic > 0 )
1117 os << ",\n m_periodic= " << m_periodic;
1118 return os.str();
1119}
1120
1121
1122template<class T>
1124{
1125 int borderKnotMult = this->borderKnotMult(); // Otherwise complains that error: 'borderKnotMult' cannot be used as a function
1126 std::vector<T> newKnots(m_knots.begin(), m_knots.end());
1127 for( int i = 0; i <= m_p - borderKnotMult; i++ )
1128 {
1129 newKnots[i] = newKnots[ newKnots.size() - 2 * m_p - 2 + i + borderKnotMult ] - _activeLength();
1130 newKnots[ newKnots.size() - i - 1 ] = newKnots[ 2 * m_p - borderKnotMult + 1 - i ] + _activeLength();
1131 }
1132 m_knots=KnotVectorType(m_p, newKnots.begin(), newKnots.end());
1133}
1134
1135template <class T>
1137{
1138 gsWarn << "gsBSplineBasis: Converting basis to periodic"<< *this<< "\n";
1139
1140 int borderKnotMult = this->borderKnotMult();
1141 if( m_knots.size() < static_cast<size_t>(2 * m_p + 2) ) // We need at least one internal knot span.
1142 {
1143 gsWarn << "Your basis cannot be changed into periodic:\n Not enough internal control points for our construction of periodic basis.\n";
1144 m_periodic = 0;
1145 }
1146
1147 else if( isClamped() )
1148 {
1149 // Clamped knots are OK, if BOTH of them are clamped.
1150 _stretchEndKnots();
1151 m_periodic = m_p + 2 - borderKnotMult; // +2 and not +1 because borderKnotMult has been computed before stretching
1152 }
1153
1154 else
1155 {
1156 m_periodic = m_p + 1 - borderKnotMult; // We suppose it will be possible to convert the basis to periodic.
1157
1158 T i1, i2;
1159 for( int i = 2; i <= 2 * m_p + 1 - borderKnotMult; i ++ )
1160 {
1161 // Compare the knot intervals in the beginning with the corresponding knot intervals at the end.
1162 i1 = m_knots[i] - m_knots[i-1];
1163 i2 = m_knots[m_knots.size() - (2*m_p) + i - 2 + borderKnotMult] - m_knots[m_knots.size() - (2*m_p) + i - 3 + borderKnotMult];
1164 if( math::abs( i1 - i2 ) > (T)(1e-8) )
1165 {
1166 gsWarn << "Your basis cannot be changed into periodic:\n Trouble stretching interior knots.\n";
1167 //std::cerr << "i: " << i << ", i1: " << i1 << ", i2: " << i2 << std::endl;
1168 m_periodic = 0;
1169 return;
1170 }
1171 }
1172 }
1173}
1174
1175template <class T>
1177{
1178 // Terminology: Blue knots are the m_p + 1st from the beginning and from the end of knot vector.
1179 if( isClamped() )
1180 return m_p + 1;
1181
1182 int multiFirst = 0;
1183 int multiLast = 0;
1184 int lastBlueId = m_knots.size() - m_p - 1; // Index of the m_p + 1st knot from the end.
1185
1186 for( int i = 0; i < m_p; i++ )
1187 {
1188 if( m_knots[m_p - i] == m_knots[m_p])
1189 multiFirst++;
1190 else
1191 break;
1192 }
1193
1194 for( int i = 0; i < m_p; i++ )
1195 {
1196 if( m_knots[lastBlueId + i] == m_knots[ lastBlueId ])
1197 multiLast++;
1198 else
1199 break;
1200 }
1201
1202 if( multiFirst == multiLast )
1203 return multiLast;
1204 else
1205 {
1206 gsWarn << "Different multiplicity of the blue knots.\n";
1207 return 0;
1208 }
1209
1210}
1211
1212template <class T>
1214{
1215 GISMO_ASSERT( isClamped(), "_stretchEndKnots() is intended for use only to knot vectors with clamped end knots.");
1216 T curFirst=m_knots[0];
1217 T curLast=m_knots[m_knots.size() - 1];
1218 m_knots.remove(curFirst); m_knots.insert(m_knots[m_knots.size() - m_p - 2] - _activeLength());
1219 m_knots.remove(curLast); m_knots.insert(m_knots[m_p + 1] + _activeLength());
1220}
1221
1222/* ********************************************** */
1223
1224template <class T>
1226{
1227 GISMO_UNUSED(i);
1228 GISMO_ASSERT(i==0,"gsBSplineBasis has only one component");
1229 return const_cast<gsBSplineBasis&>(*this);
1230}
1231
1232template <class T>
1234{
1235 GISMO_UNUSED(i);
1236 GISMO_ASSERT(i==0,"gsBSplineBasis has only one component");
1237 return const_cast<gsBSplineBasis&>(*this);
1238}
1239
1240template <class T>
1241typename gsBasis<T>::uPtr
1242gsBSplineBasis<T>::create(std::vector<KnotVectorType> cKV)
1243{
1244 typedef typename gsBasis<T>::uPtr basisPtr;
1245
1246 const index_t dd = cKV.size();
1247 switch (dd)
1248 {
1249 case 1:
1250 return basisPtr(new gsBSplineBasis<T>(give(cKV)));
1251 break;
1252 case 2:
1253 return basisPtr(new gsTensorBSplineBasis<2,T>(give(cKV)));
1254 break;
1255 case 3:
1256 return basisPtr(new gsTensorBSplineBasis<3,T>(give(cKV)));
1257 break;
1258 case 4:
1259 return basisPtr(new gsTensorBSplineBasis<4,T>(give(cKV)));
1260 break;
1261 }
1262 GISMO_ERROR("Dimension should be between 1 and 4.");
1263}
1264
1265namespace internal
1266{
1267
1269template<class T>
1270class gsXml< gsBSplineBasis<T> >
1271{
1272private:
1273 gsXml() { }
1274public:
1275 GSXML_COMMON_FUNCTIONS(gsBSplineBasis<T>);
1276 static std::string tag () { return "Basis"; }
1277 static std::string type () { return "BSplineBasis"; }
1278
1279 static gsBSplineBasis<T> * get (gsXmlNode * node)
1280 {
1281 GISMO_ASSERT( !strcmp( node->name(),"Basis"), "Wrong tag." );
1282
1283 if (!strcmp(node->first_attribute("type")->value(),"TensorBSplineBasis1"))
1284 node = node->first_node("BSplineBasis");
1285
1286 GISMO_ASSERT( !strcmp(node->first_attribute("type")->value(),"BSplineBasis"),
1287 "Wrong XML type, expected BSplineBasis." );
1288
1289 gsXmlNode * tmp = node->first_node("KnotVector");
1290 // if type: == Plain, == Compact ..
1291 GISMO_ASSERT(tmp, "Did not find a KnotVector tag in the Xml file.");
1292 gsKnotVector<T> kv;
1293 gsXml<gsKnotVector<T> >::get_into(tmp, kv);
1294
1295 return new gsBSplineBasis<T>( kv );
1296 }
1297
1298 static gsXmlNode * put (const gsBSplineBasis<T> & obj,
1299 gsXmlTree & data)
1300 {
1301 // Add a new node for obj (without data)
1302 gsXmlNode* bs_node = internal::makeNode("Basis" , data);
1303 bs_node->append_attribute( makeAttribute("type", "BSplineBasis", data) );
1304
1305 // Write the knot vector
1306 gsXmlNode* tmp =
1307 internal::gsXml<gsKnotVector<T> >::put(obj.knots(), data );
1308 bs_node->append_node(tmp);
1309
1310 // All set, return the BSPlineBasis node
1311 return bs_node;
1312 }
1313
1314};
1315
1316
1317}// namespace internal
1318
1319} // namespace gismo
Struct which represents a certain side of a box.
Definition gsBoundary.h:85
A univariate B-spline basis.
Definition gsBSplineBasis.h:700
memory::unique_ptr< gsGeometry< T > > makeGeometry(gsMatrix< T > coefs) const
Create a gsGeometry of proper type for this basis with the given coefficient matrix.
Definition gsBSplineBasis.hpp:857
Self_t & component(short_t i)
For a tensor product basis, return the 1-d basis for the i-th parameter component.
Definition gsBSplineBasis.hpp:1225
A basis represents a family of scalar basis functions defined over a common parameter domain.
Definition gsBasis.h:79
virtual void deriv2_into(const gsMatrix< T > &u, gsMatrix< T > &result) const
Evaluate the second derivatives of all active basis function at points u.
Definition gsBasis.hpp:484
virtual gsMatrix< T > laplacian(const gsMatrix< T > &u) const
Compute the Laplacian of all nonzero basis functions at points u.
Definition gsBasis.hpp:183
virtual gsBasis::uPtr tensorize(const gsBasis &other) const
Return a tensor basis of this and other.
Definition gsBasis.hpp:537
virtual void uniformCoarsen_withTransfer(gsSparseMatrix< T, RowMajor > &transfer, int numKnots=1)
Coarsen the basis uniformly and produce a sparse matrix which maps coarse coefficient vectors to refi...
Definition gsBasis.hpp:624
virtual void evalFunc_into(const gsMatrix< T > &u, const gsMatrix< T > &coefs, gsMatrix< T > &result) const
Evaluate the function described by coefs at points u, i.e., evaluates a linear combination of coefs x...
Definition gsBasis.hpp:37
virtual void uniformRefine_withTransfer(gsSparseMatrix< T, RowMajor > &transfer, int numKnots=1, int mul=1)
Refine the basis uniformly.
Definition gsBasis.hpp:611
virtual gsMatrix< index_t > allBoundary() const
Returns the indices of the basis functions that are nonzero at the domain boundary.
Definition gsBasis.hpp:334
virtual void evalAllDersSingle_into(index_t i, const gsMatrix< T > &u, int n, std::vector< gsMatrix< T > > &result) const
Evaluate the basis function i and its derivatives up to order n at points u into result.
Definition gsBasis.hpp:494
virtual 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 gsBasis.hpp:470
virtual void matchWith(const boundaryInterface &bi, const gsBasis< T > &other, gsMatrix< index_t > &bndThis, gsMatrix< index_t > &bndOther, index_t offset=0) const
Computes the indices of DoFs that match on the interface bi. The interface is assumed to be a common ...
Definition gsBasis.hpp:707
virtual void connectivity(const gsMatrix< T > &nodes, gsMesh< T > &mesh) const
Definition gsBasis.hpp:312
virtual gsBasis::uPtr create() const
Create an empty basis of the derived type and return a pointer to it.
Definition gsBasis.hpp:532
virtual void uniformRefine_withCoefs(gsMatrix< T > &coefs, int numKnots=1, int mul=1, short_t const dir=-1)
Refine the basis uniformly.
Definition gsBasis.hpp:607
virtual void eval_into(const gsMatrix< T > &u, gsMatrix< T > &result) const
Evaluates nonzero basis functions at point u into result.
Definition gsBasis.hpp:466
virtual void derivFunc_into(const gsMatrix< T > &u, const gsMatrix< T > &coefs, gsMatrix< T > &result) const
Evaluate the derivatives of the function described by coefs at points u.
Definition gsBasis.hpp:84
virtual void deriv2Func_into(const gsMatrix< T > &u, const gsMatrix< T > &coefs, gsMatrix< T > &result) const
Evaluates the second derivatives of the function described by coefs at points u.
Definition gsBasis.hpp:101
memory::unique_ptr< gsBasis > uPtr
Unique pointer for gsBasis.
Definition gsBasis.h:89
virtual gsMatrix< index_t > boundaryOffset(boxSide const &s, index_t offset) const
Definition gsBasis.hpp:339
virtual bool isActive(const index_t i, const gsVector< T > &u) const
Returns true if there the point u with non-zero value or derivatives when evaluated at the basis func...
Definition gsBasis.hpp:320
virtual void deriv_into(const gsMatrix< T > &u, gsMatrix< T > &result) const
Evaluates the first partial derivatives of the nonzero basis function.
Definition gsBasis.hpp:474
virtual gsMatrix< T > support() const
Returns (a bounding box for) the domain of the whole basis.
Definition gsBasis.hpp:454
virtual gsMatrix< T > elementInSupportOf(index_t j) const
Returns (the coordinates of) an element in the support of basis function j.
Definition gsBasis.hpp:559
virtual 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 gsBasis.hpp:488
virtual size_t elementIndex(const gsVector< T > &u) const
Returns an index for the element which contains point u.
Definition gsBasis.hpp:555
virtual void evalDerSingle_into(index_t i, const gsMatrix< T > &u, int n, gsMatrix< T > &result) const
Evaluate the (partial) derivative(s) of order n the i-th basis function at points u into result.
Definition gsBasis.hpp:525
virtual std::string detail() const
Prints the object as a string with extended details.
Definition gsBasis.h:733
virtual 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 gsBasis.hpp:478
Class defining a dummy basis of constant functions. This is used for compatibility reasons.
Definition gsConstantBasis.h:35
A specialized sparse matrix class which stores each row as a separate sparse vector.
Definition gsFiberMatrix.h:31
uPtr clone()
Clone methode. Produceds a deep copy inside a uPtr.
virtual 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 gsFunctionSet.hpp:81
memory::unique_ptr< gsGeometry > uPtr
Unique pointer for gsGeometry.
Definition gsGeometry.h:100
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
Sparse matrix class, based on gsEigen::SparseMatrix.
Definition gsSparseMatrix.h:139
index_t size() const
Returns the number of elements in the basis.
Definition gsBSplineBasis.h:165
A tensor product B-spline basis.
Definition gsTensorBSplineBasis.h:37
void refine_withCoefs(gsMatrix< T > &coefs, 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:83
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 gsTensorBSplineBasis.hpp:176
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
std::ostream & print(std::ostream &os) const
Prints the object as a string.
Definition gsTensorBSplineBasis.h:223
A tensor product of d B-spline functions, with arbitrary target dimension.
Definition gsTensorBSpline.h:45
index_t size() const
Returns the number of elements in the basis.
Definition gsTensorBasis.h:108
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
Implementation of common algorithms for B-splines.
Represents a B-spline curve/function with one parameter.
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 GISMO_ERROR(message)
Definition gsDebug.h:118
#define gsWarn
Definition gsDebug.h:50
#define GISMO_UNUSED(x)
Definition gsDebug.h:112
#define GISMO_ASSERT(cond, message)
Definition gsDebug.h:89
A specialized sparse matrix class which stores separately each fiber.
Provides declaration of the Mesh class.
Provides declaration of TensorBSplineBasis abstract interface.
Provides declaration of input/output XML utilities struct.
void evalDeg3Basis(const T &u, const KnotIterator knot, MatrixType &result)
Evaluation for degree 3 B-spline basis.
Definition gsBSplineAlgorithms.h:94
void evalDeg1Basis(const T &u, const KnotIterator knot, gsEigen::MatrixBase< Derived > const &result)
Evaluation for degree 1 B-spline basis.
Definition gsBSplineAlgorithms.h:72
void evalBasis(T u, KnotIterator knot, int deg, gsEigen::MatrixBase< Derived > const &result)
Definition gsBSplineAlgorithms.h:34
void deBoorTriangle(T u, KnotIterator knot, int deg, T N[])
Definition gsBSplineAlgorithms.h:107
void evalDeg2Basis(const T &u, const KnotIterator knot, gsEigen::MatrixBase< Derived > const &result)
Evaluation for degree 2 B-spline basis.
Definition gsBSplineAlgorithms.h:85
gsXmlNode * makeNode(const std::string &name, gsXmlTree &data)
Helper to allocate XML node.
Definition gsXml.cpp:54
gsXmlAttribute * makeAttribute(const std::string &name, const std::string &value, gsXmlTree &data)
Helper to allocate XML attribute.
Definition gsXml.cpp:37
The G+Smo namespace, containing all definitions for the library.
void gsDeboor(const gsMatrix< T > &u, const KnotVectorType &knots, int deg, const gsMatrix< T > &coefs, gsMatrix< T > &result)
Definition gsDeboor.hpp:30
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
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
Struct which represents a certain corner of a hyper-cube.
Definition gsBoundary.h:292