G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsTensorBasis.h
Go to the documentation of this file.
1
14#pragma once
15
16#include <gsCore/gsBasis.h>
17#include <gsCore/gsBoundary.h>
18
19namespace gismo
20{
21
32template<short_t d, class T>
33class gsTensorBasis : public gsBasis<T>
34{
35public:
36 typedef memory::shared_ptr< gsTensorBasis > Ptr;
37 typedef memory::unique_ptr< gsTensorBasis > uPtr;
38
40
41 typedef gsBasis<T> Basis_t;
42
44
46 typedef T Scalar_t;
47
49 static const short_t Dim = d;
50
52 typedef Basis_t** iterator;
53 typedef Basis_t* const* const_iterator;
54
55 typedef typename Basis_t::domainIter domainIter;
56
57public:
58
59 gsTensorBasis() : m_bases() { }
60
61 virtual ~gsTensorBasis() { freeAll(m_bases, m_bases+d); }
62
63 gsTensorBasis( const gsTensorBasis & o);
64 gsTensorBasis& operator=( const gsTensorBasis & o);
65
66#if EIGEN_HAS_RVALUE_REFERENCES
67 gsTensorBasis(gsTensorBasis&& other)
68 {
69 util::copy(other.m_bases, other.m_bases+d, m_bases);
70 std::fill (other.m_bases, other.m_bases+d, nullptr);
71 }
72 gsTensorBasis & operator=(gsTensorBasis&&other)
73 {
74 freeAll(m_bases, m_bases+d);
75 util::copy(other.m_bases, other.m_bases+d, m_bases);
76 std::fill (other.m_bases, other.m_bases+d, nullptr);
77 return *this;
78 }
79#endif
80 bool isValid() const { return std::find(m_bases,m_bases+d,
81 static_cast<Basis_t*>(0)) == m_bases+d; }
82
84 gsTensorBasis( Basis_t* x, Basis_t* y);
85 // template<class U> gsTensorBasis(
86 // typename enable_if<d==2 && is_same<U,Basis_t*>::value,U>::type x, U y)
87 // { m_bases[0] = x; m_bases[1] = y; }
88
90 gsTensorBasis( Basis_t* x, Basis_t* y, Basis_t* z ) ;
91
93 gsTensorBasis( Basis_t* x, Basis_t* y, Basis_t* z, Basis_t* w ) ;
94
96 explicit gsTensorBasis(iterator it)
97 {
98 for (short_t i = 0; i < d; ++i)
99 m_bases[i] = *(it++);
100 }
101
102public:
103
104 // Returns the dimension of the basis
105 short_t domainDim() const { return Dim; }
106
108 index_t size() const
109 {
110 index_t r=1;
111 for (short_t i = 0; i < d; ++i)
112 r *= m_bases[i]->size();
113 return r;
114 }
115
116 // Look at gsBasis class for a description
117 size_t numTotalElements() const
118 {
119 size_t nElem = m_bases[0]->numElements();
120 for (short_t dim = 1; dim < d; ++dim)
121 nElem *= m_bases[dim]->numElements();
122 return nElem;
123 }
124
125 // Look at gsBasis class for a description
126 size_t numElements(boxSide const & s = boundary::none) const
127 {
128 if (0==s.index()) return this->numTotalElements();
129 const short_t dir = s.direction();
130 size_t nElem = 1;
131 for (short_t dim = 0; dim < d; ++dim)
132 {
133 if(dim == dir)
134 continue;
135 nElem *= m_bases[dim]->numElements();
136 }
137 return nElem;
138 }
139
140 // Look at gsBasis class for a description
141 size_t elementIndex(const gsVector<T> & u ) const
142 {
143 GISMO_ASSERT( u.rows() == d, "Wrong vector dimension");
144
145 size_t ElIndex = m_bases[d-1]->elementIndex( u.row(d-1) );
146 for ( short_t i=d-2; i>=0; --i )
147 ElIndex = ElIndex * m_bases[i]->numElements()
148 + m_bases[i]->elementIndex( u.row(i) );
149
150 return ElIndex;
151 }
152
153 // Look at gsBasis class for a description
155
158 {
159 result.resize(d);
160 for (short_t dim = 0; dim < d; ++dim)
161 result(dim) = static_cast<unsigned>(m_bases[dim]->numElements());
162 }
163
165 void anchors_into(gsMatrix<T>& result) const;
166
168 void anchor_into(index_t i, gsMatrix<T>& result) const;
169
170 // TODO: Why is this documentation not in gsBasis?
207 virtual void active_into(const gsMatrix<T> & u, gsMatrix<index_t>& result) const;
208
209 // Look at gsBasis class for documentation
210 bool isActive(const index_t i, const gsVector<T>& u) const;
211
217 gsVector<index_t,d>& upp ) const;
218
219 // Look at gsBasis class for documentation
220 virtual void connectivity(const gsMatrix<T> & nodes, gsMesh<T> & mesh) const;
221
225
228 gsMatrix<index_t> boundaryOffset(boxSide const & s, index_t offset) const;
229
230 index_t functionAtCorner(boxCorner const & c) const;
231
233 void getComponentsForSide(boxSide const & s, std::vector<Basis_t*> & rr) const;
234
235 // see gsBasis for doxygen documentation
236 // Returns a bounding box for the basis' domain
237 gsMatrix<T> support() const ;
238
239 // see gsBasis for doxygen documentation
240 // Returns a bounding box for the support of the ith basis function
241 gsMatrix<T> support(const index_t & i ) const ;
242
243 // see gsBasis for doxygen documentation
244 // Evaluates the non-zero basis functions (and optionally their
245 // first k derivatives) at value u into result
246 virtual void eval_into(const gsMatrix<T> & u, gsMatrix<T>& result) const ;
247
248 // see gsBasis for doxygen documentation
249 // Evaluate the i-th basis function at all columns of the matrix
250 // (or vector) u
251 virtual void evalSingle_into(index_t i, const gsMatrix<T> & u, gsMatrix<T>& result) const ;
252
254 virtual void eval_into(const gsMatrix<T> & u, const gsMatrix<T> & coefs, gsMatrix<T>& result ) const;
255
256 // see gsBasis for doxygen documentation
257 // Evaluate the nonzero basis functions and their derivatives up to
258 // order n at all columns of u
259 virtual void evalAllDers_into(const gsMatrix<T> & u, int n,
260 std::vector<gsMatrix<T> >& result,
261 bool sameElement = false) const;
262
263 // see gsBasis for doxygen documentation
264 // Evaluates the gradient the non-zero basis functions at value u.
265 virtual void deriv_into(const gsMatrix<T> & u, gsMatrix<T>& result ) const;
266
267 // Evaluates the second derivatives of the non-zero basis functions at value u.
268 virtual void deriv2_into(const gsMatrix<T> & u, gsMatrix<T>& result ) const;
269
270private:
271 // Internal function
272 //
273 // values: array of std::vectors of gsMatrix<T>
274 // values[i], i = 0,...,d-1, contains the result of
275 // the univariate evalAllDers_into, corresponding to coordinate direction i.
276 //
277 // values[i] is a std::vector< gsMatrix<T> >
278 // values[i][j] contains the j-th derivatives in coordinate direction i
279 //
280 // size: gsVector of length d, size[i] contains the number of basis functions
281 // in coordinate direction i.
282 static void deriv2_tp(const std::vector< gsMatrix<T> > values[],
284 gsMatrix<T>& result);
285
286public:
287 // see gsBasis for doxygen documentation
288 // Evaluate the i-th basis function derivative at all columns of
289 virtual void derivSingle_into(index_t i, const gsMatrix<T> & u, gsMatrix<T>& result) const ;
290
291 virtual void deriv2Single_into(index_t i, const gsMatrix<T> & u, gsMatrix<T>& result) const ;
292
293 // Evaluates the (partial) derivatives of an element given by coefs at (the columns of) u.
294 //void deriv_into(const gsMatrix<T> & u, const gsMatrix<T> & coefs, gsMatrix<T>& result ) const ;
295
296 // Look at gsBasis class for documentation
297 typename gsBasis<T>::domainIter makeDomainIterator() const;
298
299 // Look at gsBasis class for documentation
300 typename gsBasis<T>::domainIter makeDomainIterator(const boxSide & s) const;
301
302 // Look at gsBasis class for documentation
303 virtual typename gsGeometry<T>::uPtr interpolateAtAnchors(gsMatrix<T> const& vals) const;
304
309 std::vector<gsMatrix<T> >const& grid) const;
310
312 virtual std::ostream &print(std::ostream &os) const = 0;
313
314 // Look at gsBasis class for documentation
315 virtual void uniformRefine(int numKnots = 1, int mul = 1, short_t dir = -1)
316 {
317 if (-1==dir)
318 for (short_t j = 0; j < d; ++j)
319 m_bases[j]->uniformRefine(numKnots,mul);
320 else
321 m_bases[dir]->uniformRefine(numKnots,mul);
322 }
323
338 void refineElements(std::vector<index_t> const & elements);
339
342 void uniformRefine_withCoefs(gsMatrix<T>& coefs, int numKnots = 1, int mul = 1, short_t const dir = -1);
343
346 void uniformRefine_withTransfer(gsSparseMatrix<T,RowMajor> & transfer, int numKnots = 1, int mul = 1);
347
348 // Look at gsBasis class for documentation
349 virtual void uniformCoarsen(int numKnots = 1)
350 {
351 for (short_t j = 0; j < d; ++j)
352 m_bases[j]->uniformCoarsen(numKnots);
353 }
354
355 // Look at gsBasis class for documentation
356 void uniformCoarsen_withTransfer(gsSparseMatrix<T,RowMajor> & transfer, int numKnots = 1);
357
358 // Look at gsBasis class for documentation
359 virtual void degreeElevate(short_t const & i = 1, short_t const dir = -1)
360 {
361 if (dir == -1)
362 {
363 for (short_t j = 0; j < d; ++j)
364 m_bases[j]->degreeElevate(i);
365 }
366 else
367 {
368 GISMO_ASSERT( dir < this->dim(),
369 "Invalid basis component requested" );
370 m_bases[dir]->degreeElevate(i);
371 }
372 }
373
374 // Look at gsBasis class for documentation
375 virtual void degreeIncrease(short_t const & i = 1, short_t const dir = -1)
376 {
377 if (dir == -1)
378 {
379 for (short_t j = 0; j < d; ++j)
380 m_bases[j]->degreeIncrease(i);
381 }
382 else
383 {
384 GISMO_ASSERT( dir < this->dim(),
385 "Invalid basis component requested" );
386 m_bases[dir]->degreeIncrease(i);
387 }
388 }
389
390 // Look at gsBasis class for documentation
391 virtual void degreeDecrease(short_t const & i = 1, short_t const dir = -1)
392 {
393 if (dir == -1)
394 {
395 for (short_t j = 0; j < d; ++j)
396 m_bases[j]->degreeDecrease(i);
397 }
398 else
399 {
400 GISMO_ASSERT( dir < this->dim(),
401 "Invalid basis component requested" );
402 m_bases[dir]->degreeDecrease(i);
403 }
404 }
405
406 // Look at gsBasis class for documentation
407 virtual void degreeReduce(short_t const & i = 1, short_t const dir = -1)
408 {
409 GISMO_UNUSED(dir);
410 GISMO_ASSERT( dir < this->dim(),
411 "Invalid basis component requested" );
412 for (short_t j = 0; j < d; ++j)
413 m_bases[j]->degreeReduce(i);
414 }
415
418 const_iterator begin() const
419 { return &m_bases[0]; }
420
423 const_iterator end() const
424 { return &m_bases[d]; }
425
429 { return &m_bases[0]; }
430
434 { return &m_bases[d]; }
435
437 index_t size(short_t k) const { return m_bases[k]->size(); }
438
440 template<int s>
441 void size_cwise(gsVector<index_t,s> & result) const
442 {
443 result.resize(d);
444 for ( short_t k = 0; k!=d; ++k )
445 result[k] = m_bases[k]->size();
446 }
447
463
466 {
467 return m_bases[i]->degree(0);
468 }
469
471 {
472 short_t td = m_bases[0]->degree(0);
473 // take maximum of coordinate bases degrees
474 for (short_t k=1; k!=d; ++k)
475 td = math::max(td, m_bases[k]->degree(0));
476 return td;
477 }
478
480 {
481 short_t td = m_bases[0]->degree(0);
482 // take minimum of coordinate bases degrees
483 for (short_t k=1; k!=d; ++k)
484 td = math::min(td, m_bases[k]->degree(0));
485 return td;
486 }
487
489 {
490 short_t td = 0;
491 for (short_t k=0; k!=d; ++k)
492 td = + m_bases[k]->degree(0);
493 return td;
494 }
495
496 gsVector<int> cwiseDegree() const
497 {
498 gsVector<int> deg(d);
499 for (short_t k=0; k!=d; ++k)
500 deg[k] = m_bases[k]->degree(0);
501 return deg;
502 }
503
506 inline unsigned index(unsigned i, unsigned j, unsigned k=0) const;
507
509 inline unsigned stride(short_t dir) const;
510
512 void stride_cwise(gsVector<index_t,d> & result) const
513 {
514 //result.resize(d);
515 result[0] = 1;
516 for ( short_t i=1; i != d; ++i )
517 result[i] = result[i-1] * m_bases[i-1]->size();
518 }
519
522 inline index_t index(gsVector<index_t,d> const & v) const;
523 // inline unsigned index(gsVector<unsigned> & v) const;
524
528 {
530 int mm = m;
531 for (short_t i = 0; i<d; ++i )
532 {
533 ind(i)= mm % size(i);
534 mm -= ind(i);
535 mm /= size(i);
536 }
537 return ind;
538 }
539
540 void swapDirections(const unsigned i, const unsigned j)
541 {
542 GISMO_ASSERT( static_cast<int>(i) < Dim && static_cast<int>(j) < Dim,
543 "Invalid basis components "<<i<<" and "<<j<<" requested" );
544 std::swap(m_bases[i],m_bases[j]);
545 }
546
549 inline bool indexOnBoundary(const gsVector<index_t, d> & ind) const
550 {
551 for ( short_t i = 0; i < d; ++i )
552 if ( ind[i] == size(i)-1 )
553 return true;
554 return ((0 == ind.array()).any() );
555 }
556
559 inline bool indexOnBoundary(const index_t m) const
560 {
561 return ( indexOnBoundary( tensorIndex(m) ) );
562 }
563
564 // see gsBasis for documentation
565 void matchWith(const boundaryInterface & bi, const gsBasis<T> & other,
566 gsMatrix<index_t> & bndThis, gsMatrix<index_t> & bndOther) const;
567
568 // see gsBasis for documentation
569 void matchWith(const boundaryInterface & bi, const gsBasis<T> & other,
570 gsMatrix<index_t> & bndThis, gsMatrix<index_t> & bndOther, index_t offset) const;
571
573 virtual T getMinCellLength() const;
574
576 virtual T getMaxCellLength() const;
577
578 Basis_t& x() const
579 {
580 return *m_bases[0];
581 }
582
583 Basis_t& y() const {
584 if (d > 1) return *m_bases[1];
585 else
586 GISMO_ERROR("gsTensorBasis has no y component");
587 }
588
589 Basis_t& z() const {
590 if (d > 2)
591 return *m_bases[2];
592 else
593 GISMO_ERROR("gsTensorBasis has no z component");
594 }
595
597 {
598 GISMO_ASSERT( dir < Dim,
599 "Invalid basis component requested" );
600 return *m_bases[dir];
601 }
602
603 const Basis_t & component(short_t dir) const
604 {
605 GISMO_ASSERT( dir < Dim,
606 "Invalid basis component requested" );
607 return *m_bases[dir];
608 }
609
610 //inline int trueSize(int k) const { return m_bases[k]->trueSize(); }
611
612// Data members
613protected:
614
615 void swap(gsTensorBasis& o)
616 { std::swap_ranges(m_bases, m_bases+d, o.m_bases); }
617
618 Basis_t* m_bases[d];
619
620}; // class gsTensorBasis
621
622// Next line disallows instantization of gsTensorBasis<0,T>
623template<typename T> class gsTensorBasis<0,T>
624{using T::GISMO_ERROR_gsTensorBasis_cannot_have_dimension_zero;};
625
626/*
627 * @brief
628 * Class for a Tensor product spline space of dimension 1.
629 * This specialization is mainly for compatibility.
630 *
631 * \tparam T coefficient type
632 *
633 * \ingroup Tensor
634 */
635template<class T>
636class gsTensorBasis<1,T> : public gsBasis<T>
637{
638public:
639 typedef memory::shared_ptr< gsTensorBasis<1,T> > Ptr;
640 typedef memory::unique_ptr< gsTensorBasis<1,T> > uPtr;
641
642 static const short_t Dim = 1;
643
644 typedef gsBasis<T> Base;
645
646 typedef gsTensorBasis<1,T> Self_t;
647
648 typedef gsBasis<T> Basis_t;
649
651 typedef T Scalar_t;
652
653 typedef gsBasis<T> CoordinateBasis;
654
656 typedef Basis_t** iterator;
657 typedef Basis_t* const* const_iterator;
658
659 typedef gsBasis<T> ** base_iterator;
660public:
661
663 gsTensorBasis() : Base()
664 { m_address = this;}
665
667 gsTensorBasis( const gsTensorBasis & o)
668 : Basis_t(o)
669 {
670 m_address = this;
671 }
672
674 gsTensorBasis& operator=( const gsTensorBasis & o)
675 {
676 this->Base::operator=(o);
677 m_address = this;
678 return *this;
679 }
680
681 // Destructor
682 virtual ~gsTensorBasis()
683 {
684 m_address = NULL;
685 }
686
687 #if EIGEN_HAS_RVALUE_REFERENCES
688 gsTensorBasis(gsTensorBasis&& other)
689 { gsTensorBasis::operator=(std::forward<gsTensorBasis>(other)); }
690 gsTensorBasis & operator=(gsTensorBasis&&other)
691 {
692 other.m_address = m_address;
693 other.m_address = nullptr;
694 return *this;
695 }
696#endif
697 bool isValid() const { return static_cast<Basis_t*>(0) != m_address; }
698
701 explicit gsTensorBasis(Base * x)
702 : Base(*x)
703 {
704 m_address = this;
705 delete x;
706 x = NULL;
707 }
708
711 explicit gsTensorBasis(base_iterator it)
712 //: Basis_t(*static_cast<Basis_t*>(*it))
713 : Basis_t(**it)
714 {
715 m_address = this;
716 delete *it;
717 *it = NULL;
718 }
719
720public:
721
722 short_t dim() const { return 1;}
723
728 void active_cwise(const gsMatrix<T> & u,
729 gsVector<index_t,1>& low,
730 gsVector<index_t,1>& upp ) const
731 {
732 gsMatrix<index_t> act;
733 this->active_into(u, act);
734 low[0]= act(0,0);
735 upp[0]= act(act.size()-1, 0 );
736 }
737
739 void stride_cwise(gsVector<index_t,1> & result) const
740 {
741 result[0] = 1;
742 }
743
746 const_iterator begin() const
747 { return &m_address; }
748
751 const_iterator end() const
752 { return (&m_address)+1; }
753
757 { return &m_address; }
758
761 iterator end()
762 { return (&m_address)+1; }
763
764 // Unhide/forward gsBasis<T>::size(), since the following
765 // overload with size(k) automatically hides it in this class
766 // Note that MSVC 2010 produces compilation error if we just
767 // do a "using gsBasis<T>::size"
768 index_t size() const = 0;
769
772 index_t size(short_t k) const
773 {
774 GISMO_UNUSED(k);
775 GISMO_ASSERT(k==0, "Invalid direction");
776 return this->size();
777 }
778
781 void size_cwise(gsVector<index_t,1> & result) const
782 { result[0] = this->size(); }
783
784 gsVector<int> cwiseDegree() const
785 {
786 gsVector<int> deg(1);
787 deg[0] = this->degree(0);
788 return deg;
789 }
790
791 void swapDirections(const unsigned i, const unsigned j)
792 {
793 GISMO_UNUSED(i);
794 GISMO_UNUSED(j);
795 GISMO_ASSERT( static_cast<int>(i) == 0 && static_cast<int>(j) == 0,
796 "Invalid basis components "<<i<<" and "<<j<<" requested" );
797 }
798
800 gsMatrix<index_t> coefSlice(short_t dir, index_t k) const
801 {
802 GISMO_UNUSED(dir);
803 GISMO_UNUSED(k);
804 GISMO_ASSERT(dir == 0, "Invalid direction");
805 GISMO_ASSERT(k < this->size(), "Invalid index");
806 // return 0 or size()-1
808 }
809
810 inline unsigned index(unsigned i) const
811 { return i; }
812
814 inline unsigned index(unsigned , unsigned ) const
815 { GISMO_ERROR("The basis is 1D"); }
816
818 inline unsigned stride(short_t dir) const
819 {
820 GISMO_UNUSED(dir);
821 GISMO_ASSERT(dir==0,"Invalid direction");
822 return 1;
823 }
824
826 void getComponentsForSide(boxSide const &, std::vector<gsBasis<T>*> & rr) const
827 { rr.clear(); }
828
831 inline index_t index(gsVector<index_t,1> const & v) const
832 { return v[0]; }
833
836 inline gsVector<index_t,1> tensorIndex(const index_t& m) const
837 {
838 return gsVector<index_t,1>::Constant(1,m);
839 }
840
841 const Basis_t& x() const
842 {
843 return *this;
844 }
845
846 Basis_t & component(short_t i)
847 {
848 GISMO_UNUSED(i);
849 GISMO_ASSERT(i==0,"Invalid component requested");
850 return *this;
851 }
852
853 const Basis_t & component(short_t i) const
854 {
855 GISMO_UNUSED(i);
856 GISMO_ASSERT(i==0,"Invalid component requested");
857 return *this;
858 }
859
860private:
861
863 Basis_t * m_address;
864
865}; // class gsTensorBasis<1,T>
866
867/* ******************************************** */
868/* ******************************************** */
869
870template<short_t d, class Basis_t >
872{
873 index_t ind;
874
875 ind = v(d-1) ;//compute global index in the tensor product
876 for ( int i=d-2; i>=0; --i )
877 ind = ind * size(i) + v(i) ;
878 return ind;
879}
880
881template<short_t d, class Basis_t >
882inline unsigned gsTensorBasis<d,Basis_t>::index(unsigned i, unsigned j, unsigned k ) const
883{
884 return size(0) * (size(1) * k + j) + i;
885}
886
887
888template<short_t d, class Basis_t >
890{
891 GISMO_ASSERT( dir>=0 && dir< this->dim(),
892 "Something went wrong with requested direction." );
893 unsigned s(1);
894 for ( short_t i=0; i<dir; ++i )
895 s *= size(i);
896 return s;
897}
898
899
900} // namespace gismo
901
902#ifndef GISMO_BUILD_LIB
903#include GISMO_HPP_HEADER(gsTensorBasis.hpp)
904#endif
905
906
Struct which represents a certain side of a box.
Definition gsBoundary.h:85
A basis represents a family of scalar basis functions defined over a common parameter domain.
Definition gsBasis.h:79
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 gsBasis.hpp:603
virtual void degreeDecrease(short_t const &i=1, short_t const dir=-1)
Lower the degree of the basis by the given amount, preserving knots multiplicity.
Definition gsBasis.hpp:641
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
virtual void degreeIncrease(short_t const &i=1, short_t const dir=-1)
Elevate the degree of the basis by the given amount, preserve knots multiplicity.
Definition gsBasis.hpp:637
virtual void degreeElevate(short_t const &i=1, short_t const dir=-1)
Elevate the degree of the basis by the given amount, preserve smoothness.
Definition gsBasis.hpp:629
virtual size_t numElements(boxSide const &s=0) const
The number of elements on side s.
Definition gsBasis.hpp:551
virtual size_t elementIndex(const gsVector< T > &u) const
Returns an index for the element which contains point u.
Definition gsBasis.hpp:555
virtual index_t size() const
size
Definition gsFunctionSet.h:613
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
Abstract base class for tensor product bases.
Definition gsTensorBasis.h:34
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 gsTensorBasis.hpp:731
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
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 gsTensorBasis.hpp:876
unsigned stride(short_t dir) const
Returns the stride for dimension dir.
Definition gsTensorBasis.h:889
void refineElements(std::vector< index_t > const &elements)
Refine elements defined by their tensor-index.
Definition gsTensorBasis.hpp:795
void anchors_into(gsMatrix< T > &result) const
Returns the anchors (graville absissae) that represent the members of the basis.
Definition gsTensorBasis.hpp:124
virtual T getMaxCellLength() const
Get the maximum mesh size, as expected for approximation error estimates.
Definition gsTensorBasis.hpp:1054
short_t maxDegree() const
If the basis is of polynomial or piecewise polynomial type, then this function returns the maximum po...
Definition gsTensorBasis.h:470
gsBasis< T >::domainIter makeDomainIterator() const
Create a domain iterator for the computational mesh of this basis, that points to the first element o...
Definition gsTensorBasis.hpp:891
bool indexOnBoundary(const index_t m) const
Returns true iff the basis function indexed m is on the boundary.
Definition gsTensorBasis.h:559
gsMatrix< index_t > coefSlice(short_t dir, index_t k) const
Returns all the basis functions with tensor-numbering.
Definition gsTensorBasis.hpp:250
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
const_iterator begin() const
Definition gsTensorBasis.h:418
gsGeometry< T >::uPtr interpolateGrid(gsMatrix< T > const &vals, std::vector< gsMatrix< T > >const &grid) const
Definition gsTensorBasis.hpp:921
void uniformRefine_withTransfer(gsSparseMatrix< T, RowMajor > &transfer, int numKnots=1, int mul=1)
Definition gsTensorBasis.hpp:861
gsMatrix< index_t > allBoundary() const
Definition gsTensorBasis.hpp:283
gsTensorBasis(iterator it)
Constructor nD (takes ownership of the passed bases)
Definition gsTensorBasis.h:96
bool indexOnBoundary(const gsVector< index_t, d > &ind) const
Returns true iff the basis function with multi-index ind is on the boundary.
Definition gsTensorBasis.h:549
void numElements_cwise(gsVector< unsigned > &result) const
Returns the number of elements (component wise)
Definition gsTensorBasis.h:157
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 gsTensorBasis.hpp:410
void stride_cwise(gsVector< index_t, d > &result) const
Returns the strides for all dimensions.
Definition gsTensorBasis.h:512
virtual void connectivity(const gsMatrix< T > &nodes, gsMesh< T > &mesh) const
Definition gsTensorBasis.hpp:163
unsigned index(unsigned i, unsigned j, unsigned k=0) const
Definition gsTensorBasis.h:882
short_t totalDegree() const
If the basis is of polynomial or piecewise polynomial type, then this function returns the total poly...
Definition gsTensorBasis.h:488
void uniformRefine_withCoefs(gsMatrix< T > &coefs, int numKnots=1, int mul=1, short_t const dir=-1)
Definition gsTensorBasis.hpp:824
virtual std::ostream & print(std::ostream &os) const =0
Prints the object as a string, pure virtual function of gsTensorBasis.
Basis_t ** iterator
Iterators on coordinate bases.
Definition gsTensorBasis.h:52
short_t degree(short_t i) const
Returns the degree of the basis wrt variable i.
Definition gsTensorBasis.h:465
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 gsTensorBasis.hpp:634
T Scalar_t
Coefficient type.
Definition gsTensorBasis.h:46
virtual void eval_into(const gsMatrix< T > &u, gsMatrix< T > &result) const
Evaluates nonzero basis functions at point u into result.
Definition gsTensorBasis.hpp:502
index_t size() const
Returns the number of elements in the basis.
Definition gsTensorBasis.h:108
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 gsTensorBasis.hpp:204
Basis_t & component(short_t dir)
For a tensor product basis, return the 1-d basis for the i-th parameter component.
Definition gsTensorBasis.h:596
static const short_t Dim
Dimension of the parameter domain.
Definition gsTensorBasis.h:49
gsMatrix< index_t > boundaryOffset(boxSide const &s, index_t offset) const
Definition gsTensorBasis.hpp:304
short_t domainDim() const
Dimension of the (source) domain.
Definition gsTensorBasis.h:105
gsVector< index_t, d > tensorIndex(const index_t &m) const
Returns the tensor index of the basis function with global index m.
Definition gsTensorBasis.h:527
virtual gsGeometry< T >::uPtr interpolateAtAnchors(gsMatrix< T > const &vals) const
Applies interpolation of values pts using the anchors as parameter points. May be reimplemented in de...
Definition gsTensorBasis.hpp:908
void getComponentsForSide(boxSide const &s, std::vector< Basis_t * > &rr) const
Returns the components for a basis on the face s.
Definition gsTensorBasis.hpp:377
void anchor_into(index_t i, gsMatrix< T > &result) const
Returns the anchors (graville absissae) that represent the members of the basis.
Definition gsTensorBasis.hpp:148
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 gsTensorBasis.hpp:239
virtual void deriv_into(const gsMatrix< T > &u, gsMatrix< T > &result) const
Evaluates the first partial derivatives of the nonzero basis function.
Definition gsTensorBasis.hpp:594
virtual T getMinCellLength() const
Get the minimum mesh size, as expected for inverse inequalities.
Definition gsTensorBasis.hpp:1042
gsMatrix< T > support() const
Returns (a bounding box for) the domain of the whole basis.
Definition gsTensorBasis.hpp:390
iterator end()
Definition gsTensorBasis.h:433
virtual void degreeElevate(short_t const &i=1, short_t const dir=-1)
Elevate the degree of the basis by the given amount, preserve smoothness.
Definition gsTensorBasis.h:359
const_iterator end() const
Definition gsTensorBasis.h:423
virtual void degreeIncrease(short_t const &i=1, short_t const dir=-1)
Elevate the degree of the basis by the given amount, preserve knots multiplicity.
Definition gsTensorBasis.h:375
iterator begin()
Definition gsTensorBasis.h:428
index_t size(short_t k) const
The number of basis functions in the direction of the k-th parameter component.
Definition gsTensorBasis.h:437
short_t minDegree() const
If the basis is of polynomial or piecewise polynomial type, then this function returns the minimum po...
Definition gsTensorBasis.h:479
size_t numElements(boxSide const &s=boundary::none) const
The number of elements on side s.
Definition gsTensorBasis.h:126
gsMatrix< T > elementInSupportOf(index_t j) const
Returns (the coordinates of) an element in the support of basis function j.
Definition gsTensorBasis.hpp:1066
virtual void degreeDecrease(short_t const &i=1, short_t const dir=-1)
Lower the degree of the basis by the given amount, preserving knots multiplicity.
Definition gsTensorBasis.h:391
virtual void degreeReduce(short_t const &i=1, short_t const dir=-1)
Reduce the degree of the basis by the given amount, preserve smoothness.
Definition gsTensorBasis.h:407
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 gsTensorBasis.h:603
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 gsTensorBasis.hpp:454
size_t elementIndex(const gsVector< T > &u) const
Returns an index for the element which contains point u.
Definition gsTensorBasis.h:141
index_t index(gsVector< index_t, d > const &v) const
Definition gsTensorBasis.h:871
void size_cwise(gsVector< index_t, s > &result) const
The number of basis functions in the direction of the k-th parameter component.
Definition gsTensorBasis.h:441
void active_cwise(const gsMatrix< T > &u, gsVector< index_t, d > &low, gsVector< index_t, d > &upp) const
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 gsTensorBasis.hpp:429
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
Provides declaration of Basis abstract interface.
Provides structs and classes related to interfaces and boundaries.
#define short_t
Definition gsConfig.h:35
#define index_t
Definition gsConfig.h:32
#define GISMO_NO_IMPLEMENTATION
Definition gsDebug.h:129
#define GISMO_ERROR(message)
Definition gsDebug.h:118
#define GISMO_UNUSED(x)
Definition gsDebug.h:112
#define GISMO_ASSERT(cond, message)
Definition gsDebug.h:89
void copy(T begin, T end, U *result)
Small wrapper for std::copy mimicking std::copy for a raw pointer destination, copies n positions sta...
Definition gsMemory.h:391
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
Struct which represents a certain corner of a hyper-cube.
Definition gsBoundary.h:292