3 #ifndef DUNE_FUNCTIONS_FUNCTIONSPACEBASES_BSPLINEBASIS_HH
4 #define DUNE_FUNCTIONS_FUNCTIONSPACEBASES_BSPLINEBASIS_HH
14 #include <dune/common/dynmatrix.hh>
16 #include <dune/localfunctions/common/localbasis.hh>
17 #include <dune/common/diagonalmatrix.hh>
18 #include <dune/localfunctions/common/localkey.hh>
19 #include <dune/localfunctions/common/localfiniteelementtraits.hh>
20 #include <dune/geometry/type.hh>
21 #include <dune/functions/functionspacebases/nodes.hh>
22 #include <dune/functions/functionspacebases/defaultglobalbasis.hh>
23 #include <dune/functions/functionspacebases/leafprebasismixin.hh>
30 template<
typename GV,
typename R>
31 class BSplineLocalFiniteElement;
34 class BSplinePreBasis;
45 template<
class GV,
class R>
50 typedef typename GV::ctype D;
51 enum {dim = GV::dimension};
55 typedef LocalBasisTraits<D,dim,FieldVector<D,dim>,R,1,FieldVector<R,1>,
64 : preBasis_(preBasis),
72 std::vector<FieldVector<R,1> >& out)
const
74 FieldVector<D,dim> globalIn = offset_;
75 scaling_.umv(in,globalIn);
77 preBasis_.evaluateFunction(globalIn, out, lFE_.currentKnotSpan_);
84 std::vector<FieldMatrix<D,1,dim> >& out)
const
86 FieldVector<D,dim> globalIn = offset_;
87 scaling_.umv(in,globalIn);
89 preBasis_.evaluateJacobian(globalIn, out, lFE_.currentKnotSpan_);
91 for (
size_t i=0; i<out.size(); i++)
92 for (
int j=0; j<dim; j++)
93 out[i][0][j] *= scaling_[j][j];
98 inline void evaluate (
const typename std::array<int,k>& directions,
99 const typename Traits::DomainType& in,
100 std::vector<typename Traits::RangeType>& out)
const
109 FieldVector<D,dim> globalIn = offset_;
110 scaling_.umv(in,globalIn);
112 preBasis_.evaluate(directions, globalIn, out, lFE_.currentKnotSpan_);
114 for (
size_t i=0; i<out.size(); i++)
115 out[i][0] *= scaling_[directions[0]][directions[0]];
120 FieldVector<D,dim> globalIn = offset_;
121 scaling_.umv(in,globalIn);
123 preBasis_.evaluate(directions, globalIn, out, lFE_.currentKnotSpan_);
125 for (
size_t i=0; i<out.size(); i++)
126 out[i][0] *= scaling_[directions[0]][directions[0]]*scaling_[directions[1]][directions[1]];
130 DUNE_THROW(NotImplemented,
"B-Spline derivatives of order " << k <<
" not implemented yet!");
143 return *std::max_element(preBasis_.order_.begin(), preBasis_.order_.end());
160 FieldVector<D,dim> offset_;
161 DiagonalMatrix<D,dim> scaling_;
181 std::array<unsigned int,dim> multiindex (
unsigned int i)
const
183 std::array<unsigned int,dim> alpha;
184 for (
int j=0; j<dim; j++)
186 alpha[j] = i % sizes_[j];
193 void setup1d(std::vector<unsigned int>& subEntity)
204 unsigned lastIndex=0;
205 subEntity[lastIndex++] = 0;
206 for (
unsigned i = 0; i < sizes_[0] - 2; ++i)
207 subEntity[lastIndex++] = 0;
209 subEntity[lastIndex++] = 1;
211 assert(
size()==lastIndex);
214 void setup2d(std::vector<unsigned int>& subEntity)
216 unsigned lastIndex=0;
230 subEntity[lastIndex++] = 0;
231 for (
unsigned i = 0; i < sizes_[0]-2; ++i)
232 subEntity[lastIndex++] = 2;
234 subEntity[lastIndex++] = 1;
237 for (
unsigned e = 0; e < sizes_[1]-2; ++e)
239 subEntity[lastIndex++] = 0;
240 for (
unsigned i = 0; i < sizes_[0]-2; ++i)
241 subEntity[lastIndex++] = 0;
242 subEntity[lastIndex++] = 1;
246 subEntity[lastIndex++] = 2;
247 for (
unsigned i = 0; i < sizes_[0]-2; ++i)
248 subEntity[lastIndex++] = 3;
250 subEntity[lastIndex++] = 3;
252 assert(
size()==lastIndex);
257 void init(
const std::array<unsigned,dim>& sizes)
264 std::vector<unsigned int> codim(li_.size());
266 for (std::size_t i=0; i<codim.size(); i++)
271 std::array<unsigned int,dim> mIdx = multiindex(i);
272 for (
int j=0; j<dim; j++)
273 if (mIdx[j]==0 or mIdx[j]==sizes[j]-1)
282 std::vector<unsigned int> index(
size());
284 for (std::size_t i=0; i<index.size(); i++)
288 std::array<unsigned int,dim> mIdx = multiindex(i);
290 for (
int j=dim-1; j>=0; j--)
291 if (mIdx[j]>0 and mIdx[j]<sizes[j]-1)
292 index[i] = (sizes[j]-1)*index[i] + (mIdx[j]-1);
296 std::vector<unsigned int> subEntity(li_.size());
298 if (subEntity.size() > 0)
304 }
else if (dim==2 and sizes_[0]>1 and sizes_[1]>1) {
311 for (
size_t i=0; i<li_.size(); i++)
312 li_[i] = LocalKey(subEntity[i], codim[i], index[i]);
318 return std::accumulate(sizes_.begin(), sizes_.end(), 1, std::multiplies<unsigned int>());
330 std::array<unsigned, dim> sizes_;
332 std::vector<LocalKey> li_;
339 template<
int dim,
class LB>
344 template<
typename F,
typename C>
347 DUNE_THROW(NotImplemented,
"BSplineLocalInterpolation::interpolate");
361 template<
class GV,
class R>
364 typedef typename GV::ctype D;
365 enum {dim = GV::dimension};
371 typedef LocalFiniteElementTraits<BSplineLocalBasis<GV,R>,
378 : preBasis_(preBasis),
379 localBasis_(preBasis,*this)
385 : preBasis_(other.preBasis_),
386 localBasis_(preBasis_,*this)
395 void bind(
const std::array<unsigned,dim>& elementIdx)
398 for (
size_t i=0; i<elementIdx.size(); i++)
400 currentKnotSpan_[i] = 0;
403 while (preBasis_.knotVectors_[i][currentKnotSpan_[i]+1] < preBasis_.knotVectors_[i][currentKnotSpan_[i]]+1e-8)
404 currentKnotSpan_[i]++;
406 for (
size_t j=0; j<elementIdx[i]; j++)
408 currentKnotSpan_[i]++;
411 while (preBasis_.knotVectors_[i][currentKnotSpan_[i]+1] < preBasis_.knotVectors_[i][currentKnotSpan_[i]]+1e-8)
412 currentKnotSpan_[i]++;
416 localBasis_.offset_[i] = preBasis_.knotVectors_[i][currentKnotSpan_[i]];
417 localBasis_.scaling_[i][i] = preBasis_.knotVectors_[i][currentKnotSpan_[i]+1] - preBasis_.knotVectors_[i][currentKnotSpan_[i]];
421 std::array<unsigned int, dim> sizes;
422 for (
size_t i=0; i<dim; i++)
424 localCoefficients_.init(sizes);
436 return localCoefficients_;
442 return localInterpolation_;
449 for (
int i=0; i<dim; i++)
458 return GeometryTypes::cube(dim);
466 const auto& order = preBasis_.order_;
467 unsigned int r = order[i]+1;
468 if (currentKnotSpan_[i]<order[i])
469 r -= (order[i] - currentKnotSpan_[i]);
470 if ( order[i] > (preBasis_.knotVectors_[i].size() - currentKnotSpan_[i] - 2) )
471 r -= order[i] - (preBasis_.knotVectors_[i].size() - currentKnotSpan_[i] - 2);
482 std::array<unsigned,dim> currentKnotSpan_;
486 template<
typename GV>
498 template<
typename GV>
504 static const int dim = GV::dimension;
507 class MultiDigitCounter
514 MultiDigitCounter(
const std::array<unsigned int,dim>& limits)
517 std::fill(counter_.begin(), counter_.end(), 0);
521 MultiDigitCounter& operator++()
523 for (
int i=0; i<dim; i++)
528 if (counter_[i] < limits_[i])
537 const unsigned int& operator[](
int i)
const
543 unsigned int cycle()
const
546 for (
int i=0; i<dim; i++)
554 const std::array<unsigned int,dim> limits_;
557 std::array<unsigned int,dim> counter_;
565 using size_type = std::size_t;
567 using Node = BSplineNode<GV>;
591 const std::vector<double>& knotVector,
593 bool makeOpen =
true)
601 assert( std::accumulate(
elements_.begin(),
elements_.end(), 1, std::multiplies<unsigned>()) == gridView_.size(0) );
603 for (
int i=0; i<dim; i++)
608 for (
unsigned int j=0; j<order; j++)
614 for (
unsigned int j=0; j<order; j++)
643 const FieldVector<double,dim>& lowerLeft,
644 const FieldVector<double,dim>& upperRight,
645 const std::array<unsigned int,dim>& elements,
647 bool makeOpen =
true)
653 assert( std::accumulate(
elements_.begin(),
elements_.end(), 1, std::multiplies<unsigned>()) == gridView_.size(0) );
655 for (
int i=0; i<dim; i++)
660 for (
unsigned int j=0; j<order; j++)
664 for (
size_t j=0; j<elements[i]+1; j++)
665 knotVectors_[i].push_back(lowerLeft[i] + j*(upperRight[i]-lowerLeft[i]) / elements[i]);
668 for (
unsigned int j=0; j<order; j++)
702 size_type result = 1;
703 for (
int i=0; i<dim; i++)
709 template<
typename It>
714 std::array<unsigned int, dim> localSizes;
715 for (
int i=0; i<dim; i++)
716 localSizes[i] = node.finiteElement().size(i);
717 for (size_type i = 0, end = node.size() ; i < end ; ++i, ++it)
719 std::array<unsigned int,dim> localIJK =
getIJK(i, localSizes);
721 const auto currentKnotSpan = node.finiteElement().currentKnotSpan_;
722 const auto order =
order_;
724 std::array<unsigned int,dim> globalIJK;
725 for (
int i=0; i<dim; i++)
726 globalIJK[i] = std::max((
int)currentKnotSpan[i] - (
int)order[i], 0) + localIJK[i];
729 size_type globalIdx = globalIJK[dim-1];
731 for (
int i=dim-2; i>=0; i--)
732 globalIdx = globalIdx *
size(i) + globalIJK[i];
742 unsigned int result = 1;
743 for (
size_t i=0; i<dim; i++)
749 unsigned int size (
size_t d)
const
759 std::vector<FieldVector<R,1> >& out,
760 const std::array<unsigned,dim>& currentKnotSpan)
const
763 std::array<std::vector<R>, dim> oneDValues;
765 for (
size_t i=0; i<dim; i++)
768 std::array<unsigned int, dim> limits;
769 for (
int i=0; i<dim; i++)
770 limits[i] = oneDValues[i].
size();
772 MultiDigitCounter ijkCounter(limits);
774 out.resize(ijkCounter.cycle());
776 for (
size_t i=0; i<out.size(); i++, ++ijkCounter)
779 for (
size_t j=0; j<dim; j++)
780 out[i] *= oneDValues[j][ijkCounter[j]];
790 std::vector<FieldMatrix<R,1,dim> >& out,
791 const std::array<unsigned,dim>& currentKnotSpan)
const
794 std::array<unsigned int, dim> limits;
795 for (
int i=0; i<dim; i++)
798 if (currentKnotSpan[i]<
order_[i])
799 limits[i] -= (
order_[i] - currentKnotSpan[i]);
805 std::array<unsigned int, dim> offset;
806 for (
int i=0; i<dim; i++)
807 offset[i] = std::max((
int)(currentKnotSpan[i] -
order_[i]),0);
810 std::array<std::vector<R>, dim> oneDValues;
813 std::array<std::vector<R>, dim> lowOrderOneDValues;
815 std::array<DynamicMatrix<R>, dim> values;
817 for (
size_t i=0; i<dim; i++)
821 for (
size_t j=0; j<oneDValues[i].size(); j++)
822 oneDValues[i][j] = values[i][
order_[i]][j];
827 for (
size_t j=0; j<lowOrderOneDValues[i].size(); j++)
828 lowOrderOneDValues[i][j] = values[i][
order_[i]-1][j];
834 std::array<std::vector<R>, dim> oneDDerivatives;
835 for (
size_t i=0; i<dim; i++)
837 oneDDerivatives[i].resize(limits[i]);
840 std::fill(oneDDerivatives[i].begin(), oneDDerivatives[i].end(), R(0.0));
843 for (
size_t j=offset[i]; j<offset[i]+limits[i]; j++)
848 if (std::isnan(derivativeAddend1))
849 derivativeAddend1 = 0;
850 if (std::isnan(derivativeAddend2))
851 derivativeAddend2 = 0;
852 oneDDerivatives[i][j-offset[i]] =
order_[i] * ( derivativeAddend1 - derivativeAddend2 );
859 std::array<std::vector<R>, dim> oneDValuesShort;
861 for (
int i=0; i<dim; i++)
863 oneDValuesShort[i].resize(limits[i]);
865 for (
size_t j=0; j<limits[i]; j++)
866 oneDValuesShort[i][j] = oneDValues[i][offset[i] + j];
872 MultiDigitCounter ijkCounter(limits);
874 out.resize(ijkCounter.cycle());
877 for (
size_t i=0; i<out.size(); i++, ++ijkCounter)
878 for (
int j=0; j<dim; j++)
881 for (
int k=0; k<dim; k++)
882 out[i][0][j] *= (j==k) ? oneDDerivatives[k][ijkCounter[k]]
883 : oneDValuesShort[k][ijkCounter[k]];
889 template <
size_type k>
890 void evaluate(
const typename std::array<int,k>& directions,
891 const FieldVector<typename GV::ctype,dim>& in,
892 std::vector<FieldVector<R,1> >& out,
893 const std::array<unsigned,dim>& currentKnotSpan)
const
895 if (k != 1 && k != 2)
896 DUNE_THROW(RangeError,
"Differentiation order greater than 2 is not supported!");
899 std::array<std::vector<R>, dim> oneDValues;
900 std::array<std::vector<R>, dim> oneDDerivatives;
901 std::array<std::vector<R>, dim> oneDSecondDerivatives;
905 for (
size_t i=0; i<dim; i++)
906 evaluateAll(in[i], oneDValues[i],
true, oneDDerivatives[i],
false, oneDSecondDerivatives[i],
knotVectors_[i],
order_[i], currentKnotSpan[i]);
908 for (
size_t i=0; i<dim; i++)
912 std::array<unsigned int, dim> offset;
913 for (
int i=0; i<dim; i++)
914 offset[i] = std::max((
int)(currentKnotSpan[i] -
order_[i]),0);
917 std::array<unsigned int, dim> limits;
918 for (
int i=0; i<dim; i++)
923 if (currentKnotSpan[i]<
order_[i])
924 limits[i] -= (
order_[i] - currentKnotSpan[i]);
931 std::array<std::vector<R>, dim> oneDValuesShort;
933 for (
int i=0; i<dim; i++)
935 oneDValuesShort[i].resize(limits[i]);
937 for (
size_t j=0; j<limits[i]; j++)
938 oneDValuesShort[i][j] = oneDValues[i][offset[i] + j];
942 MultiDigitCounter ijkCounter(limits);
944 out.resize(ijkCounter.cycle());
949 for (
size_t i=0; i<out.size(); i++, ++ijkCounter)
952 for (
int l=0; l<dim; l++)
953 out[i][0] *= (directions[0]==l) ? oneDDerivatives[l][ijkCounter[l]]
954 : oneDValuesShort[l][ijkCounter[l]];
961 for (
size_t i=0; i<out.size(); i++, ++ijkCounter)
964 for (
int j=0; j<dim; j++)
966 if (directions[0] != directions[1])
967 if (directions[0] == j || directions[1] == j)
968 out[i][0] *= oneDDerivatives[j][ijkCounter[j]];
970 out[i][0] *= oneDValuesShort[j][ijkCounter[j]];
972 if (directions[0] == j)
973 out[i][0] *= oneDSecondDerivatives[j][ijkCounter[j]];
975 out[i][0] *= oneDValuesShort[j][ijkCounter[j]];
986 static std::array<unsigned int,dim>
getIJK(
typename GridView::IndexSet::IndexType idx, std::array<unsigned int,dim> elements)
988 std::array<unsigned,dim> result;
989 for (
int i=0; i<dim; i++)
991 result[i] = idx%elements[i];
1006 const std::vector<R>& knotVector,
1008 unsigned int currentKnotSpan)
1010 std::size_t outSize = order+1;
1011 if (currentKnotSpan<order)
1012 outSize -= (order - currentKnotSpan);
1013 if ( order > (knotVector.size() - currentKnotSpan - 2) )
1014 outSize -= order - (knotVector.size() - currentKnotSpan - 2);
1015 out.resize(outSize);
1018 DynamicMatrix<R> N(order+1, knotVector.size()-1);
1026 for (
size_t i=0; i<knotVector.size()-1; i++)
1027 N[0][i] = (i == currentKnotSpan);
1029 for (
size_t r=1; r<=order; r++)
1030 for (
size_t i=0; i<knotVector.size()-r-1; i++)
1032 R factor1 = ((knotVector[i+r] - knotVector[i]) > 1e-10)
1033 ? (in - knotVector[i]) / (knotVector[i+r] - knotVector[i])
1035 R factor2 = ((knotVector[i+r+1] - knotVector[i+1]) > 1e-10)
1036 ? (knotVector[i+r+1] - in) / (knotVector[i+r+1] - knotVector[i+1])
1038 N[r][i] = factor1 * N[r-1][i] + factor2 * N[r-1][i+1];
1045 for (
size_t i=0; i<out.size(); i++) {
1046 out[i] = N[order][std::max((
int)(currentKnotSpan - order),0) + i];
1063 DynamicMatrix<R>& out,
1064 const std::vector<R>& knotVector,
1066 unsigned int currentKnotSpan)
1069 DynamicMatrix<R>& N = out;
1071 N.resize(order+1, knotVector.size()-1);
1079 for (
size_t i=0; i<knotVector.size()-1; i++)
1080 N[0][i] = (i == currentKnotSpan);
1082 for (
size_t r=1; r<=order; r++)
1083 for (
size_t i=0; i<knotVector.size()-r-1; i++)
1085 R factor1 = ((knotVector[i+r] - knotVector[i]) > 1e-10)
1086 ? (in - knotVector[i]) / (knotVector[i+r] - knotVector[i])
1088 R factor2 = ((knotVector[i+r+1] - knotVector[i+1]) > 1e-10)
1089 ? (knotVector[i+r+1] - in) / (knotVector[i+r+1] - knotVector[i+1])
1091 N[r][i] = factor1 * N[r-1][i] + factor2 * N[r-1][i+1];
1106 std::vector<R>& out,
1108 bool evaluateHessian, std::vector<R>& outHess,
1109 const std::vector<R>& knotVector,
1111 unsigned int currentKnotSpan)
1116 if (currentKnotSpan<order)
1117 limit -= (order - currentKnotSpan);
1118 if ( order > (knotVector.size() - currentKnotSpan - 2) )
1119 limit -= order - (knotVector.size() - currentKnotSpan - 2);
1122 unsigned int offset;
1123 offset = std::max((
int)(currentKnotSpan - order),0);
1126 DynamicMatrix<R> values;
1130 out.resize(knotVector.size()-order-1);
1131 for (
size_t j=0; j<out.size(); j++)
1132 out[j] = values[order][j];
1135 std::vector<R> lowOrderOneDValues;
1139 lowOrderOneDValues.resize(knotVector.size()-(order-1)-1);
1140 for (
size_t j=0; j<lowOrderOneDValues.size(); j++)
1141 lowOrderOneDValues[j] = values[order-1][j];
1145 std::vector<R> lowOrderTwoDValues;
1147 if (order>1 && evaluateHessian)
1149 lowOrderTwoDValues.resize(knotVector.size()-(order-2)-1);
1150 for (
size_t j=0; j<lowOrderTwoDValues.size(); j++)
1151 lowOrderTwoDValues[j] = values[order-2][j];
1157 outJac.resize(limit);
1160 std::fill(outJac.begin(), outJac.end(), R(0.0));
1163 for (
size_t j=offset; j<offset+limit; j++)
1165 R derivativeAddend1 = lowOrderOneDValues[j] / (knotVector[j+order]-knotVector[j]);
1166 R derivativeAddend2 = lowOrderOneDValues[j+1] / (knotVector[j+order+1]-knotVector[j+1]);
1168 if (std::isnan(derivativeAddend1))
1169 derivativeAddend1 = 0;
1170 if (std::isnan(derivativeAddend2))
1171 derivativeAddend2 = 0;
1172 outJac[j-offset] = order * ( derivativeAddend1 - derivativeAddend2 );
1178 if (evaluateHessian)
1180 outHess.resize(limit);
1183 std::fill(outHess.begin(), outHess.end(), R(0.0));
1186 for (
size_t j=offset; j<offset+limit; j++)
1188 assert(j+2 < lowOrderTwoDValues.size());
1189 R derivativeAddend1 = lowOrderTwoDValues[j] / (knotVector[j+order]-knotVector[j]) / (knotVector[j+order-1]-knotVector[j]);
1190 R derivativeAddend2 = lowOrderTwoDValues[j+1] / (knotVector[j+order]-knotVector[j]) / (knotVector[j+order]-knotVector[j+1]);
1191 R derivativeAddend3 = lowOrderTwoDValues[j+1] / (knotVector[j+order+1]-knotVector[j+1]) / (knotVector[j+order]-knotVector[j+1]);
1192 R derivativeAddend4 = lowOrderTwoDValues[j+2] / (knotVector[j+order+1]-knotVector[j+1]) / (knotVector[j+1+order]-knotVector[j+2]);
1195 if (std::isnan(derivativeAddend1))
1196 derivativeAddend1 = 0;
1197 if (std::isnan(derivativeAddend2))
1198 derivativeAddend2 = 0;
1199 if (std::isnan(derivativeAddend3))
1200 derivativeAddend3 = 0;
1201 if (std::isnan(derivativeAddend4))
1202 derivativeAddend4 = 0;
1203 outHess[j-offset] = order * (order-1) * ( derivativeAddend1 - derivativeAddend2 -derivativeAddend3 + derivativeAddend4 );
1224 template<
typename GV>
1226 public LeafBasisNode
1228 static const int dim = GV::dimension;
1232 using size_type = std::size_t;
1233 using Element =
typename GV::template Codim<0>::Entity;
1237 preBasis_(preBasis),
1238 finiteElement_(*preBasis)
1242 const Element& element()
const
1251 const FiniteElement& finiteElement()
const
1253 return finiteElement_;
1257 void bind(
const Element& e)
1260 auto elementIndex = preBasis_->gridView().indexSet().index(e);
1261 finiteElement_.bind(preBasis_->getIJK(elementIndex,preBasis_->elements_));
1262 this->setSize(finiteElement_.size());
1267 const BSplinePreBasis<GV>* preBasis_;
1269 FiniteElement finiteElement_;
1275 namespace BasisFactory {
1283 inline auto bSpline(
const std::vector<double>& knotVector,
1285 bool makeOpen =
true)
1287 return [&knotVector, order, makeOpen](
const auto& gridView) {
1288 return BSplinePreBasis<std::decay_t<decltype(gridView)>>(gridView, knotVector, order, makeOpen);
1304 template<
typename GV>
LocalBasis class in the sense of dune-localfunctions, presenting the restriction of a B-spline patch ...
Definition: bsplinebasis.hh:47
LocalBasisTraits< D, dim, FieldVector< D, dim >, R, 1, FieldVector< R, 1 >, FieldMatrix< R, 1, dim > > Traits
export type traits for function signature
Definition: bsplinebasis.hh:56
unsigned int order() const
Polynomial order of the shape functions.
Definition: bsplinebasis.hh:141
std::size_t size() const
Return the number of basis functions on the current knot span.
Definition: bsplinebasis.hh:148
void evaluate(const typename std::array< int, k > &directions, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions and derivatives of any order.
Definition: bsplinebasis.hh:98
void evaluateFunction(const FieldVector< D, dim > &in, std::vector< FieldVector< R, 1 > > &out) const
Evaluate all shape functions.
Definition: bsplinebasis.hh:71
void evaluateJacobian(const FieldVector< D, dim > &in, std::vector< FieldMatrix< D, 1, dim > > &out) const
Evaluate Jacobian of all shape functions.
Definition: bsplinebasis.hh:83
BSplineLocalBasis(const BSplinePreBasis< GV > &preBasis, const BSplineLocalFiniteElement< GV, R > &lFE)
Constructor with a given B-spline patch.
Definition: bsplinebasis.hh:62
Attaches a shape function to an entity.
Definition: bsplinebasis.hh:179
const LocalKey & localKey(std::size_t i) const
get i'th index
Definition: bsplinebasis.hh:322
std::size_t size() const
number of coefficients
Definition: bsplinebasis.hh:316
LocalFiniteElement in the sense of dune-localfunctions, for the B-spline basis on tensor-product grid...
Definition: bsplinebasis.hh:363
BSplineLocalFiniteElement(const BSplineLocalFiniteElement &other)
Copy constructor.
Definition: bsplinebasis.hh:384
const BSplineLocalBasis< GV, R > & localBasis() const
Hand out a LocalBasis object.
Definition: bsplinebasis.hh:428
LocalFiniteElementTraits< BSplineLocalBasis< GV, R >, BSplineLocalCoefficients< dim >, BSplineLocalInterpolation< dim, BSplineLocalBasis< GV, R > > > Traits
Export various types related to this LocalFiniteElement.
Definition: bsplinebasis.hh:373
BSplineLocalFiniteElement(const BSplinePreBasis< GV > &preBasis)
Constructor with a given B-spline basis.
Definition: bsplinebasis.hh:377
void bind(const std::array< unsigned, dim > &elementIdx)
Bind LocalFiniteElement to a specific knot span of the spline patch.
Definition: bsplinebasis.hh:395
const BSplineLocalCoefficients< dim > & localCoefficients() const
Hand out a LocalCoefficients object.
Definition: bsplinebasis.hh:434
GeometryType type() const
Return the reference element that the local finite element is defined on (here, a hypercube)
Definition: bsplinebasis.hh:456
unsigned size() const
Number of shape functions in this finite element.
Definition: bsplinebasis.hh:446
unsigned int size(int i) const
Number of degrees of freedom for one coordinate direction.
Definition: bsplinebasis.hh:464
const BSplineLocalInterpolation< dim, BSplineLocalBasis< GV, R > > & localInterpolation() const
Hand out a LocalInterpolation object.
Definition: bsplinebasis.hh:440
Local interpolation in the sense of dune-localfunctions, for the B-spline basis on tensor-product gri...
Definition: bsplinebasis.hh:341
void interpolate(const F &f, std::vector< C > &out) const
Local interpolation of a function.
Definition: bsplinebasis.hh:345
Pre-basis for B-spline basis.
Definition: bsplinebasis.hh:501
std::array< unsigned, dim > elements_
Number of grid elements in the different coordinate directions.
Definition: bsplinebasis.hh:1217
static std::array< unsigned int, dim > getIJK(typename GridView::IndexSet::IndexType idx, std::array< unsigned int, dim > elements)
Compute integer element coordinates from the element index.
Definition: bsplinebasis.hh:986
const GridView & gridView() const
Obtain the grid view that the basis is defined on.
Definition: bsplinebasis.hh:680
static void evaluateFunctionFull(const typename GV::ctype &in, DynamicMatrix< R > &out, const std::vector< R > &knotVector, unsigned int order, unsigned int currentKnotSpan)
Evaluate all one-dimensional B-spline functions for a given coordinate direction.
Definition: bsplinebasis.hh:1062
void evaluateFunction(const FieldVector< typename GV::ctype, dim > &in, std::vector< FieldVector< R, 1 > > &out, const std::array< unsigned, dim > ¤tKnotSpan) const
Evaluate all B-spline basis functions at a given point.
Definition: bsplinebasis.hh:758
std::array< unsigned int, dim > order_
Order of the B-spline for each space dimension.
Definition: bsplinebasis.hh:1211
static void evaluateAll(const typename GV::ctype &in, std::vector< R > &out, bool evaluateJacobian, std::vector< R > &outJac, bool evaluateHessian, std::vector< R > &outHess, const std::vector< R > &knotVector, unsigned int order, unsigned int currentKnotSpan)
Evaluate the second derivatives of all one-dimensional B-spline functions for a given coordinate dire...
Definition: bsplinebasis.hh:1105
static void evaluateFunction(const typename GV::ctype &in, std::vector< R > &out, const std::vector< R > &knotVector, unsigned int order, unsigned int currentKnotSpan)
Evaluate all one-dimensional B-spline functions for a given coordinate direction.
Definition: bsplinebasis.hh:1005
unsigned int size(size_t d) const
Number of shape functions in one direction.
Definition: bsplinebasis.hh:749
GV GridView
The grid view that the FE space is defined on.
Definition: bsplinebasis.hh:564
It indices(const Node &node, It it) const
Maps from subtree index set [0..size-1] to a globally unique multi index in global basis.
Definition: bsplinebasis.hh:710
void evaluate(const typename std::array< int, k > &directions, const FieldVector< typename GV::ctype, dim > &in, std::vector< FieldVector< R, 1 > > &out, const std::array< unsigned, dim > ¤tKnotSpan) const
Evaluate Derivatives of all B-spline basis functions.
Definition: bsplinebasis.hh:890
void update(const GridView &gv)
Update the stored grid view, to be called if the grid has changed.
Definition: bsplinebasis.hh:686
unsigned int dimension() const
Total number of B-spline basis functions.
Definition: bsplinebasis.hh:740
void initializeIndices()
Initialize the global indices.
Definition: bsplinebasis.hh:676
Node makeNode() const
Create tree node.
Definition: bsplinebasis.hh:694
BSplinePreBasis(const GridView &gridView, const std::vector< double > &knotVector, unsigned int order, bool makeOpen=true)
Construct a B-spline basis for a given grid view and set of knot vectors.
Definition: bsplinebasis.hh:590
void evaluateJacobian(const FieldVector< typename GV::ctype, dim > &in, std::vector< FieldMatrix< R, 1, dim > > &out, const std::array< unsigned, dim > ¤tKnotSpan) const
Evaluate Jacobian of all B-spline basis functions.
Definition: bsplinebasis.hh:789
std::array< std::vector< double >, dim > knotVectors_
The knot vectors, one for each space dimension.
Definition: bsplinebasis.hh:1214
size_type maxNodeSize() const
Get the maximal number of DOFs associated to node for any element.
Definition: bsplinebasis.hh:700
BSplinePreBasis(const GridView &gridView, const FieldVector< double, dim > &lowerLeft, const FieldVector< double, dim > &upperRight, const std::array< unsigned int, dim > &elements, unsigned int order, bool makeOpen=true)
Construct a B-spline basis for a given grid view with uniform knot vectors.
Definition: bsplinebasis.hh:642
Global basis for given pre-basis.
Definition: defaultglobalbasis.hh:46
A generic MixIn class for PreBasis.
Definition: leafprebasismixin.hh:32
size_type size() const
Get the total dimension of the space spanned by this basis.
Definition: leafprebasismixin.hh:56
auto bSpline(const std::vector< double > &knotVector, unsigned int order, bool makeOpen=true)
Create a pre-basis factory that can create a B-spline pre-basis.
Definition: bsplinebasis.hh:1283
Definition: polynomial.hh:13