Dune Core Modules (2.9.0)

bsplinebasis.hh
Go to the documentation of this file.
1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 #ifndef DUNE_FUNCTIONS_FUNCTIONSPACEBASES_BSPLINEBASIS_HH
4 #define DUNE_FUNCTIONS_FUNCTIONSPACEBASES_BSPLINEBASIS_HH
5 
10 #include <array>
11 #include <numeric>
12 
14 #include <dune/common/dynmatrix.hh>
15 
16 #include <dune/localfunctions/common/localbasis.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 
24 namespace Dune
25 {
26 namespace Functions {
27 
28 // A maze of dependencies between the different parts of this. We need a few forward declarations
29 template<typename GV, typename R>
30 class BSplineLocalFiniteElement;
31 
32 template<typename GV>
33 class BSplinePreBasis;
34 
35 
44 template<class GV, class R>
46 {
47  friend class BSplineLocalFiniteElement<GV,R>;
48 
49  typedef typename GV::ctype D;
50  enum {dim = GV::dimension};
51 public:
52 
56 
63  : preBasis_(preBasis),
64  lFE_(lFE)
65  {}
66 
71  std::vector<FieldVector<R,1> >& out) const
72  {
73  FieldVector<D,dim> globalIn = offset_;
74  scaling_.umv(in,globalIn);
75 
76  preBasis_.evaluateFunction(globalIn, out, lFE_.currentKnotSpan_);
77  }
78 
83  std::vector<FieldMatrix<D,1,dim> >& out) const
84  {
85  FieldVector<D,dim> globalIn = offset_;
86  scaling_.umv(in,globalIn);
87 
88  preBasis_.evaluateJacobian(globalIn, out, lFE_.currentKnotSpan_);
89 
90  for (size_t i=0; i<out.size(); i++)
91  for (int j=0; j<dim; j++)
92  out[i][0][j] *= scaling_[j][j];
93  }
94 
96  template<size_t k>
97  inline void evaluate (const typename std::array<int,k>& directions,
98  const typename Traits::DomainType& in,
99  std::vector<typename Traits::RangeType>& out) const
100  {
101  switch(k)
102  {
103  case 0:
104  evaluateFunction(in, out);
105  break;
106  case 1:
107  {
108  FieldVector<D,dim> globalIn = offset_;
109  scaling_.umv(in,globalIn);
110 
111  preBasis_.evaluate(directions, globalIn, out, lFE_.currentKnotSpan_);
112 
113  for (size_t i=0; i<out.size(); i++)
114  out[i][0] *= scaling_[directions[0]][directions[0]];
115  break;
116  }
117  case 2:
118  {
119  FieldVector<D,dim> globalIn = offset_;
120  scaling_.umv(in,globalIn);
121 
122  preBasis_.evaluate(directions, globalIn, out, lFE_.currentKnotSpan_);
123 
124  for (size_t i=0; i<out.size(); i++)
125  out[i][0] *= scaling_[directions[0]][directions[0]]*scaling_[directions[1]][directions[1]];
126  break;
127  }
128  default:
129  DUNE_THROW(NotImplemented, "B-Spline derivatives of order " << k << " not implemented yet!");
130  }
131  }
132 
140  unsigned int order () const
141  {
142  return *std::max_element(preBasis_.order_.begin(), preBasis_.order_.end());
143  }
144 
147  std::size_t size() const
148  {
149  return lFE_.size();
150  }
151 
152 private:
153  const BSplinePreBasis<GV>& preBasis_;
154 
156 
157  // Coordinates in a single knot span differ from coordinates on the B-spline patch
158  // by an affine transformation. This transformation is stored in offset_ and scaling_.
159  FieldVector<D,dim> offset_;
160  DiagonalMatrix<D,dim> scaling_;
161 };
162 
176 template<int dim>
178 {
179  // Return i as a d-digit number in the (k+1)-nary system
180  std::array<unsigned int,dim> multiindex (unsigned int i) const
181  {
182  std::array<unsigned int,dim> alpha;
183  for (int j=0; j<dim; j++)
184  {
185  alpha[j] = i % sizes_[j];
186  i = i/sizes_[j];
187  }
188  return alpha;
189  }
190 
192  void setup1d(std::vector<unsigned int>& subEntity)
193  {
194  if (sizes_[0]==1)
195  {
196  subEntity[0] = 0;
197  return;
198  }
199 
200  /* edge and vertex numbering
201  0----0----1
202  */
203  unsigned lastIndex=0;
204  subEntity[lastIndex++] = 0; // corner 0
205  for (unsigned i = 0; i < sizes_[0] - 2; ++i)
206  subEntity[lastIndex++] = 0; // inner dofs of element (0)
207 
208  subEntity[lastIndex++] = 1; // corner 1
209 
210  assert(size()==lastIndex);
211  }
212 
213  void setup2d(std::vector<unsigned int>& subEntity)
214  {
215  unsigned lastIndex=0;
216 
217  // LocalKey: entity number , entity codim, dof indices within each entity
218  /* edge and vertex numbering
219  2----3----3
220  | |
221  | |
222  0 1
223  | |
224  | |
225  0----2----1
226  */
227 
228  // lower edge (2)
229  subEntity[lastIndex++] = 0; // corner 0
230  for (unsigned i = 0; i < sizes_[0]-2; ++i)
231  subEntity[lastIndex++] = 2; // inner dofs of lower edge (2)
232 
233  subEntity[lastIndex++] = 1; // corner 1
234 
235  // iterate from bottom to top over inner edge dofs
236  for (unsigned e = 0; e < sizes_[1]-2; ++e)
237  {
238  subEntity[lastIndex++] = 0; // left edge (0)
239  for (unsigned i = 0; i < sizes_[0]-2; ++i)
240  subEntity[lastIndex++] = 0; // face dofs
241  subEntity[lastIndex++] = 1; // right edge (1)
242  }
243 
244  // upper edge (3)
245  subEntity[lastIndex++] = 2; // corner 2
246  for (unsigned i = 0; i < sizes_[0]-2; ++i)
247  subEntity[lastIndex++] = 3; // inner dofs of upper edge (3)
248 
249  subEntity[lastIndex++] = 3; // corner 3
250 
251  assert(size()==lastIndex);
252  }
253 
254 
255 public:
256  void init(const std::array<unsigned,dim>& sizes)
257  {
258  sizes_ = sizes;
259 
260  li_.resize(size());
261 
262  // Set up array of codimension-per-dof-number
263  std::vector<unsigned int> codim(li_.size());
264 
265  for (std::size_t i=0; i<codim.size(); i++)
266  {
267  codim[i] = 0;
268  // Codimension gets increased by 1 for each coordinate direction
269  // where dof is on boundary
270  std::array<unsigned int,dim> mIdx = multiindex(i);
271  for (int j=0; j<dim; j++)
272  if (mIdx[j]==0 or mIdx[j]==sizes[j]-1)
273  codim[i]++;
274  }
275 
276  // Set up index vector (the index of the dof in the set of dofs of a given subentity)
277  // Algorithm: the 'index' has the same ordering as the dof number 'i'.
278  // To make it consecutive we interpret 'i' in the (k+1)-adic system, omit all digits
279  // that correspond to axes where the dof is on the element boundary, and transform the
280  // rest to the (k-1)-adic system.
281  std::vector<unsigned int> index(size());
282 
283  for (std::size_t i=0; i<index.size(); i++)
284  {
285  index[i] = 0;
286 
287  std::array<unsigned int,dim> mIdx = multiindex(i);
288 
289  for (int j=dim-1; j>=0; j--)
290  if (mIdx[j]>0 and mIdx[j]<sizes[j]-1)
291  index[i] = (sizes[j]-1)*index[i] + (mIdx[j]-1);
292  }
293 
294  // Set up entity and dof numbers for each (supported) dimension separately
295  std::vector<unsigned int> subEntity(li_.size());
296 
297  if (subEntity.size() > 0)
298  {
299  if (dim==1) {
300 
301  setup1d(subEntity);
302 
303  } else if (dim==2 and sizes_[0]>1 and sizes_[1]>1) {
304 
305  setup2d(subEntity);
306 
307  }
308  }
309 
310  for (size_t i=0; i<li_.size(); i++)
311  li_[i] = LocalKey(subEntity[i], codim[i], index[i]);
312  }
313 
315  std::size_t size () const
316  {
317  return std::accumulate(sizes_.begin(), sizes_.end(), 1, std::multiplies<unsigned int>());
318  }
319 
321  const LocalKey& localKey (std::size_t i) const
322  {
323  return li_[i];
324  }
325 
326 private:
327 
328  // Number of shape functions on this element per coordinate direction
329  std::array<unsigned, dim> sizes_;
330 
331  std::vector<LocalKey> li_;
332 };
333 
338 template<int dim, class LB>
340 {
341 public:
343  template<typename F, typename C>
344  void interpolate (const F& f, std::vector<C>& out) const
345  {
346  DUNE_THROW(NotImplemented, "BSplineLocalInterpolation::interpolate");
347  }
348 
349 };
350 
360 template<class GV, class R>
362 {
363  typedef typename GV::ctype D;
364  enum {dim = GV::dimension};
365  friend class BSplineLocalBasis<GV,R>;
366 public:
367 
373 
377  : preBasis_(preBasis),
378  localBasis_(preBasis,*this)
379  {}
380 
384  : preBasis_(other.preBasis_),
385  localBasis_(preBasis_,*this)
386  {}
387 
394  void bind(const std::array<unsigned,dim>& elementIdx)
395  {
396  /* \todo In the long run we need to precompute a table for this */
397  for (size_t i=0; i<elementIdx.size(); i++)
398  {
399  currentKnotSpan_[i] = 0;
400 
401  // Skip over degenerate knot spans
402  while (preBasis_.knotVectors_[i][currentKnotSpan_[i]+1] < preBasis_.knotVectors_[i][currentKnotSpan_[i]]+1e-8)
403  currentKnotSpan_[i]++;
404 
405  for (size_t j=0; j<elementIdx[i]; j++)
406  {
407  currentKnotSpan_[i]++;
408 
409  // Skip over degenerate knot spans
410  while (preBasis_.knotVectors_[i][currentKnotSpan_[i]+1] < preBasis_.knotVectors_[i][currentKnotSpan_[i]]+1e-8)
411  currentKnotSpan_[i]++;
412  }
413 
414  // Compute the geometric transformation from knotspan-local to global coordinates
415  localBasis_.offset_[i] = preBasis_.knotVectors_[i][currentKnotSpan_[i]];
416  localBasis_.scaling_[i][i] = preBasis_.knotVectors_[i][currentKnotSpan_[i]+1] - preBasis_.knotVectors_[i][currentKnotSpan_[i]];
417  }
418 
419  // Set up the LocalCoefficients object
420  std::array<unsigned int, dim> sizes;
421  for (size_t i=0; i<dim; i++)
422  sizes[i] = size(i);
423  localCoefficients_.init(sizes);
424  }
425 
428  {
429  return localBasis_;
430  }
431 
434  {
435  return localCoefficients_;
436  }
437 
440  {
441  return localInterpolation_;
442  }
443 
445  unsigned size () const
446  {
447  std::size_t r = 1;
448  for (int i=0; i<dim; i++)
449  r *= size(i);
450  return r;
451  }
452 
456  {
457  return GeometryTypes::cube(dim);
458  }
459 
460 //private:
461 
463  unsigned int size(int i) const
464  {
465  const auto& order = preBasis_.order_;
466  unsigned int r = order[i]+1; // The 'normal' value
467  if (currentKnotSpan_[i]<order[i]) // Less near the left end of the knot vector
468  r -= (order[i] - currentKnotSpan_[i]);
469  if ( order[i] > (preBasis_.knotVectors_[i].size() - currentKnotSpan_[i] - 2) )
470  r -= order[i] - (preBasis_.knotVectors_[i].size() - currentKnotSpan_[i] - 2);
471  return r;
472  }
473 
474  const BSplinePreBasis<GV>& preBasis_;
475 
476  BSplineLocalBasis<GV,R> localBasis_;
477  BSplineLocalCoefficients<dim> localCoefficients_;
479 
480  // The knot span we are bound to
481  std::array<unsigned,dim> currentKnotSpan_;
482 };
483 
484 
485 template<typename GV>
486 class BSplineNode;
487 
497 template<typename GV>
499 {
500  static const int dim = GV::dimension;
501 
503  class MultiDigitCounter
504  {
505  public:
506 
510  MultiDigitCounter(const std::array<unsigned int,dim>& limits)
511  : limits_(limits)
512  {
513  std::fill(counter_.begin(), counter_.end(), 0);
514  }
515 
517  MultiDigitCounter& operator++()
518  {
519  for (int i=0; i<dim; i++)
520  {
521  ++counter_[i];
522 
523  // no overflow?
524  if (counter_[i] < limits_[i])
525  break;
526 
527  counter_[i] = 0;
528  }
529  return *this;
530  }
531 
533  const unsigned int& operator[](int i) const
534  {
535  return counter_[i];
536  }
537 
539  unsigned int cycle() const
540  {
541  unsigned int r = 1;
542  for (int i=0; i<dim; i++)
543  r *= limits_[i];
544  return r;
545  }
546 
547  private:
548 
550  const std::array<unsigned int,dim> limits_;
551 
553  std::array<unsigned int,dim> counter_;
554 
555  };
556 
557 public:
558 
560  using GridView = GV;
561  using size_type = std::size_t;
562 
563  using Node = BSplineNode<GV>;
564 
565  static constexpr size_type maxMultiIndexSize = 1;
566  static constexpr size_type minMultiIndexSize = 1;
567  static constexpr size_type multiIndexBufferSize = 1;
568 
569  // Type used for function values
570  using R = double;
571 
591  const std::vector<double>& knotVector,
592  unsigned int order,
593  bool makeOpen = true)
594  : gridView_(gridView)
595  {
596  // \todo Detection of duplicate knots
597  std::fill(elements_.begin(), elements_.end(), knotVector.size()-1);
598 
599  // Mediocre sanity check: we don't know the number of grid elements in each direction.
600  // but at least we know the total number of elements.
601  assert( std::accumulate(elements_.begin(), elements_.end(), 1, std::multiplies<unsigned>()) == gridView_.size(0) );
602 
603  for (int i=0; i<dim; i++)
604  {
605  // Prepend the correct number of additional knots to open the knot vector
607  if (makeOpen)
608  for (unsigned int j=0; j<order; j++)
609  knotVectors_[i].push_back(knotVector[0]);
610 
611  knotVectors_[i].insert(knotVectors_[i].end(), knotVector.begin(), knotVector.end());
612 
613  if (makeOpen)
614  for (unsigned int j=0; j<order; j++)
615  knotVectors_[i].push_back(knotVector.back());
616  }
617 
618  std::fill(order_.begin(), order_.end(), order);
619  }
620 
643  const FieldVector<double,dim>& lowerLeft,
644  const FieldVector<double,dim>& upperRight,
645  const std::array<unsigned int,dim>& elements,
646  unsigned int order,
647  bool makeOpen = true)
648  : elements_(elements),
649  gridView_(gridView)
650  {
651  // Mediocre sanity check: we don't know the number of grid elements in each direction.
652  // but at least we know the total number of elements.
653  assert( std::accumulate(elements_.begin(), elements_.end(), 1, std::multiplies<unsigned>()) == gridView_.size(0) );
654 
655  for (int i=0; i<dim; i++)
656  {
657  // Prepend the correct number of additional knots to open the knot vector
659  if (makeOpen)
660  for (unsigned int j=0; j<order; j++)
661  knotVectors_[i].push_back(lowerLeft[i]);
662 
663  // Construct the actual knot vector
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]);
666 
667  if (makeOpen)
668  for (unsigned int j=0; j<order; j++)
669  knotVectors_[i].push_back(upperRight[i]);
670  }
671 
672  std::fill(order_.begin(), order_.end(), order);
673  }
674 
677  {}
678 
680  const GridView& gridView() const
681  {
682  return gridView_;
683  }
684 
686  void update(const GridView& gv)
687  {
688  gridView_ = gv;
689  }
690 
694  Node makeNode() const
695  {
696  return Node{this};
697  }
698 
699  // Ideally this method should be implemented as
700  //
701  // template<class SizePrefix>
702  // size_type size(const SizePrefix& prefix) const
703  //
704  // But leads to ambiguity with the other size method:
705  //
706  // unsigned int size (size_t d) const
707  //
708  // Once the latter is removed, this implementation should be changed.
709 
711  template<class ST, int i>
712  size_type size(const Dune::ReservedVector<ST, i>& prefix) const
713  {
714  assert(prefix.size() == 0 || prefix.size() == 1);
715  return (prefix.size() == 0) ? size() : 0;
716  }
717 
719  size_type dimension() const
720  {
721  return size();
722  }
723 
725  size_type maxNodeSize() const
726  {
727  size_type result = 1;
728  for (int i=0; i<dim; i++)
729  result *= order_[i]+1;
730  return result;
731  }
732 
734  template<typename It>
735  It indices(const Node& node, It it) const
736  {
737  // Local degrees of freedom are arranged in a lattice.
738  // We need the lattice dimensions to be able to compute lattice coordinates from a local index
739  std::array<unsigned int, dim> localSizes;
740  for (int i=0; i<dim; i++)
741  localSizes[i] = node.finiteElement().size(i);
742  for (size_type i = 0, end = node.size() ; i < end ; ++i, ++it)
743  {
744  std::array<unsigned int,dim> localIJK = getIJK(i, localSizes);
745 
746  const auto currentKnotSpan = node.finiteElement().currentKnotSpan_;
747  const auto order = order_;
748 
749  std::array<unsigned int,dim> globalIJK;
750  for (int i=0; i<dim; i++)
751  globalIJK[i] = std::max((int)currentKnotSpan[i] - (int)order[i], 0) + localIJK[i]; // needs to be a signed type!
752 
753  // Make one global flat index from the globalIJK tuple
754  size_type globalIdx = globalIJK[dim-1];
755 
756  for (int i=dim-2; i>=0; i--)
757  globalIdx = globalIdx * size(i) + globalIJK[i];
758 
759  *it = {{globalIdx}};
760  }
761  return it;
762  }
763 
765  unsigned int size () const
766  {
767  unsigned int result = 1;
768  for (size_t i=0; i<dim; i++)
769  result *= size(i);
770  return result;
771  }
772 
774  unsigned int size (size_t d) const
775  {
776  return knotVectors_[d].size() - order_[d] - 1;
777  }
778 
782  std::vector<FieldVector<R,1> >& out,
783  const std::array<unsigned,dim>& currentKnotSpan) const
784  {
785  // Evaluate
786  std::array<std::vector<R>, dim> oneDValues;
787 
788  for (size_t i=0; i<dim; i++)
789  evaluateFunction(in[i], oneDValues[i], knotVectors_[i], order_[i], currentKnotSpan[i]);
790 
791  std::array<unsigned int, dim> limits;
792  for (int i=0; i<dim; i++)
793  limits[i] = oneDValues[i].size();
794 
795  MultiDigitCounter ijkCounter(limits);
796 
797  out.resize(ijkCounter.cycle());
798 
799  for (size_t i=0; i<out.size(); i++, ++ijkCounter)
800  {
801  out[i] = R(1.0);
802  for (size_t j=0; j<dim; j++)
803  out[i] *= oneDValues[j][ijkCounter[j]];
804  }
805  }
806 
813  std::vector<FieldMatrix<R,1,dim> >& out,
814  const std::array<unsigned,dim>& currentKnotSpan) const
815  {
816  // How many shape functions to we have in each coordinate direction?
817  std::array<unsigned int, dim> limits;
818  for (int i=0; i<dim; i++)
819  {
820  limits[i] = order_[i]+1; // The 'standard' value away from the boundaries of the knot vector
821  if (currentKnotSpan[i]<order_[i])
822  limits[i] -= (order_[i] - currentKnotSpan[i]);
823  if ( order_[i] > (knotVectors_[i].size() - currentKnotSpan[i] - 2) )
824  limits[i] -= order_[i] - (knotVectors_[i].size() - currentKnotSpan[i] - 2);
825  }
826 
827  // The lowest knot spans that we need values from
828  std::array<unsigned int, dim> offset;
829  for (int i=0; i<dim; i++)
830  offset[i] = std::max((int)(currentKnotSpan[i] - order_[i]),0);
831 
832  // Evaluate 1d function values (needed for the product rule)
833  std::array<std::vector<R>, dim> oneDValues;
834 
835  // Evaluate 1d function values of one order lower (needed for the derivative formula)
836  std::array<std::vector<R>, dim> lowOrderOneDValues;
837 
838  std::array<DynamicMatrix<R>, dim> values;
839 
840  for (size_t i=0; i<dim; i++)
841  {
842  evaluateFunctionFull(in[i], values[i], knotVectors_[i], order_[i], currentKnotSpan[i]);
843  oneDValues[i].resize(knotVectors_[i].size()-order_[i]-1);
844  for (size_t j=0; j<oneDValues[i].size(); j++)
845  oneDValues[i][j] = values[i][order_[i]][j];
846 
847  if (order_[i]!=0)
848  {
849  lowOrderOneDValues[i].resize(knotVectors_[i].size()-(order_[i]-1)-1);
850  for (size_t j=0; j<lowOrderOneDValues[i].size(); j++)
851  lowOrderOneDValues[i][j] = values[i][order_[i]-1][j];
852  }
853  }
854 
855 
856  // Evaluate 1d function derivatives
857  std::array<std::vector<R>, dim> oneDDerivatives;
858  for (size_t i=0; i<dim; i++)
859  {
860  oneDDerivatives[i].resize(limits[i]);
861 
862  if (order_[i]==0) // order-zero functions are piecewise constant, hence all derivatives are zero
863  std::fill(oneDDerivatives[i].begin(), oneDDerivatives[i].end(), R(0.0));
864  else
865  {
866  for (size_t j=offset[i]; j<offset[i]+limits[i]; j++)
867  {
868  R derivativeAddend1 = lowOrderOneDValues[i][j] / (knotVectors_[i][j+order_[i]]-knotVectors_[i][j]);
869  R derivativeAddend2 = lowOrderOneDValues[i][j+1] / (knotVectors_[i][j+order_[i]+1]-knotVectors_[i][j+1]);
870  // The two previous terms may evaluate as 0/0. This is to be interpreted as 0.
871  if (std::isnan(derivativeAddend1))
872  derivativeAddend1 = 0;
873  if (std::isnan(derivativeAddend2))
874  derivativeAddend2 = 0;
875  oneDDerivatives[i][j-offset[i]] = order_[i] * ( derivativeAddend1 - derivativeAddend2 );
876  }
877  }
878  }
879 
880  // Working towards computing only the parts that we really need:
881  // Let's copy them out into a separate array
882  std::array<std::vector<R>, dim> oneDValuesShort;
883 
884  for (int i=0; i<dim; i++)
885  {
886  oneDValuesShort[i].resize(limits[i]);
887 
888  for (size_t j=0; j<limits[i]; j++)
889  oneDValuesShort[i][j] = oneDValues[i][offset[i] + j];
890  }
891 
892 
893 
894  // Set up a multi-index to go from consecutive indices to integer coordinates
895  MultiDigitCounter ijkCounter(limits);
896 
897  out.resize(ijkCounter.cycle());
898 
899  // Complete Jacobian is given by the product rule
900  for (size_t i=0; i<out.size(); i++, ++ijkCounter)
901  for (int j=0; j<dim; j++)
902  {
903  out[i][0][j] = 1.0;
904  for (int k=0; k<dim; k++)
905  out[i][0][j] *= (j==k) ? oneDDerivatives[k][ijkCounter[k]]
906  : oneDValuesShort[k][ijkCounter[k]];
907  }
908 
909  }
910 
912  template <size_type k>
913  void evaluate(const typename std::array<int,k>& directions,
915  std::vector<FieldVector<R,1> >& out,
916  const std::array<unsigned,dim>& currentKnotSpan) const
917  {
918  if (k != 1 && k != 2)
919  DUNE_THROW(RangeError, "Differentiation order greater than 2 is not supported!");
920 
921  // Evaluate 1d function values (needed for the product rule)
922  std::array<std::vector<R>, dim> oneDValues;
923  std::array<std::vector<R>, dim> oneDDerivatives;
924  std::array<std::vector<R>, dim> oneDSecondDerivatives;
925 
926  // Evaluate 1d function derivatives
927  if (k==1)
928  for (size_t i=0; i<dim; i++)
929  evaluateAll(in[i], oneDValues[i], true, oneDDerivatives[i], false, oneDSecondDerivatives[i], knotVectors_[i], order_[i], currentKnotSpan[i]);
930  else
931  for (size_t i=0; i<dim; i++)
932  evaluateAll(in[i], oneDValues[i], true, oneDDerivatives[i], true, oneDSecondDerivatives[i], knotVectors_[i], order_[i], currentKnotSpan[i]);
933 
934  // The lowest knot spans that we need values from
935  std::array<unsigned int, dim> offset;
936  for (int i=0; i<dim; i++)
937  offset[i] = std::max((int)(currentKnotSpan[i] - order_[i]),0);
938 
939  // Set up a multi-index to go from consecutive indices to integer coordinates
940  std::array<unsigned int, dim> limits;
941  for (int i=0; i<dim; i++)
942  {
943  // In a proper implementation, the following line would do
944  //limits[i] = oneDValues[i].size();
945  limits[i] = order_[i]+1; // The 'standard' value away from the boundaries of the knot vector
946  if (currentKnotSpan[i]<order_[i])
947  limits[i] -= (order_[i] - currentKnotSpan[i]);
948  if ( order_[i] > (knotVectors_[i].size() - currentKnotSpan[i] - 2) )
949  limits[i] -= order_[i] - (knotVectors_[i].size() - currentKnotSpan[i] - 2);
950  }
951 
952  // Working towards computing only the parts that we really need:
953  // Let's copy them out into a separate array
954  std::array<std::vector<R>, dim> oneDValuesShort;
955 
956  for (int i=0; i<dim; i++)
957  {
958  oneDValuesShort[i].resize(limits[i]);
959 
960  for (size_t j=0; j<limits[i]; j++)
961  oneDValuesShort[i][j] = oneDValues[i][offset[i] + j];
962  }
963 
964 
965  MultiDigitCounter ijkCounter(limits);
966 
967  out.resize(ijkCounter.cycle());
968 
969  if (k == 1)
970  {
971  // Complete Jacobian is given by the product rule
972  for (size_t i=0; i<out.size(); i++, ++ijkCounter)
973  {
974  out[i][0] = 1.0;
975  for (int l=0; l<dim; l++)
976  out[i][0] *= (directions[0]==l) ? oneDDerivatives[l][ijkCounter[l]]
977  : oneDValuesShort[l][ijkCounter[l]];
978  }
979  }
980 
981  if (k == 2)
982  {
983  // Complete derivation by deriving the tensor product
984  for (size_t i=0; i<out.size(); i++, ++ijkCounter)
985  {
986  out[i][0] = 1.0;
987  for (int j=0; j<dim; j++)
988  {
989  if (directions[0] != directions[1]) //derivation in two different variables
990  if (directions[0] == j || directions[1] == j) //the spline has to be derived (once) in this direction
991  out[i][0] *= oneDDerivatives[j][ijkCounter[j]];
992  else //no derivation in this direction
993  out[i][0] *= oneDValuesShort[j][ijkCounter[j]];
994  else //spline is derived two times in the same direction
995  if (directions[0] == j) //the spline is derived two times in this direction
996  out[i][0] *= oneDSecondDerivatives[j][ijkCounter[j]];
997  else //no derivation in this direction
998  out[i][0] *= oneDValuesShort[j][ijkCounter[j]];
999  }
1000  }
1001  }
1002  }
1003 
1004 
1009  static std::array<unsigned int,dim> getIJK(typename GridView::IndexSet::IndexType idx, std::array<unsigned int,dim> elements)
1010  {
1011  std::array<unsigned,dim> result;
1012  for (int i=0; i<dim; i++)
1013  {
1014  result[i] = idx%elements[i];
1015  idx /= elements[i];
1016  }
1017  return result;
1018  }
1019 
1028  static void evaluateFunction (const typename GV::ctype& in, std::vector<R>& out,
1029  const std::vector<R>& knotVector,
1030  unsigned int order,
1031  unsigned int currentKnotSpan)
1032  {
1033  std::size_t outSize = order+1; // The 'standard' value away from the boundaries of the knot vector
1034  if (currentKnotSpan<order) // Less near the left end of the knot vector
1035  outSize -= (order - currentKnotSpan);
1036  if ( order > (knotVector.size() - currentKnotSpan - 2) )
1037  outSize -= order - (knotVector.size() - currentKnotSpan - 2);
1038  out.resize(outSize);
1039 
1040  // It's not really a matrix that is needed here, a plain 2d array would do
1041  DynamicMatrix<R> N(order+1, knotVector.size()-1);
1042 
1043  // The text books on splines use the following geometric condition here to fill the array N
1044  // (see for example Cottrell, Hughes, Bazilevs, Formula (2.1). However, this condition
1045  // only works if splines are never evaluated exactly on the knots.
1046  //
1047  // for (size_t i=0; i<knotVector.size()-1; i++)
1048  // N[0][i] = (knotVector[i] <= in) and (in < knotVector[i+1]);
1049  for (size_t i=0; i<knotVector.size()-1; i++)
1050  N[0][i] = (i == currentKnotSpan);
1051 
1052  for (size_t r=1; r<=order; r++)
1053  for (size_t i=0; i<knotVector.size()-r-1; i++)
1054  {
1055  R factor1 = ((knotVector[i+r] - knotVector[i]) > 1e-10)
1056  ? (in - knotVector[i]) / (knotVector[i+r] - knotVector[i])
1057  : 0;
1058  R factor2 = ((knotVector[i+r+1] - knotVector[i+1]) > 1e-10)
1059  ? (knotVector[i+r+1] - in) / (knotVector[i+r+1] - knotVector[i+1])
1060  : 0;
1061  N[r][i] = factor1 * N[r-1][i] + factor2 * N[r-1][i+1];
1062  }
1063 
1068  for (size_t i=0; i<out.size(); i++) {
1069  out[i] = N[order][std::max((int)(currentKnotSpan - order),0) + i];
1070  }
1071  }
1072 
1085  static void evaluateFunctionFull(const typename GV::ctype& in,
1086  DynamicMatrix<R>& out,
1087  const std::vector<R>& knotVector,
1088  unsigned int order,
1089  unsigned int currentKnotSpan)
1090  {
1091  // It's not really a matrix that is needed here, a plain 2d array would do
1092  DynamicMatrix<R>& N = out;
1093 
1094  N.resize(order+1, knotVector.size()-1);
1095 
1096  // The text books on splines use the following geometric condition here to fill the array N
1097  // (see for example Cottrell, Hughes, Bazilevs, Formula (2.1). However, this condition
1098  // only works if splines are never evaluated exactly on the knots.
1099  //
1100  // for (size_t i=0; i<knotVector.size()-1; i++)
1101  // N[0][i] = (knotVector[i] <= in) and (in < knotVector[i+1]);
1102  for (size_t i=0; i<knotVector.size()-1; i++)
1103  N[0][i] = (i == currentKnotSpan);
1104 
1105  for (size_t r=1; r<=order; r++)
1106  for (size_t i=0; i<knotVector.size()-r-1; i++)
1107  {
1108  R factor1 = ((knotVector[i+r] - knotVector[i]) > 1e-10)
1109  ? (in - knotVector[i]) / (knotVector[i+r] - knotVector[i])
1110  : 0;
1111  R factor2 = ((knotVector[i+r+1] - knotVector[i+1]) > 1e-10)
1112  ? (knotVector[i+r+1] - in) / (knotVector[i+r+1] - knotVector[i+1])
1113  : 0;
1114  N[r][i] = factor1 * N[r-1][i] + factor2 * N[r-1][i+1];
1115  }
1116  }
1117 
1118 
1128  static void evaluateAll(const typename GV::ctype& in,
1129  std::vector<R>& out,
1130  bool evaluateJacobian, std::vector<R>& outJac,
1131  bool evaluateHessian, std::vector<R>& outHess,
1132  const std::vector<R>& knotVector,
1133  unsigned int order,
1134  unsigned int currentKnotSpan)
1135  {
1136  // How many shape functions to we have in each coordinate direction?
1137  unsigned int limit;
1138  limit = order+1; // The 'standard' value away from the boundaries of the knot vector
1139  if (currentKnotSpan<order)
1140  limit -= (order - currentKnotSpan);
1141  if ( order > (knotVector.size() - currentKnotSpan - 2) )
1142  limit -= order - (knotVector.size() - currentKnotSpan - 2);
1143 
1144  // The lowest knot spans that we need values from
1145  unsigned int offset;
1146  offset = std::max((int)(currentKnotSpan - order),0);
1147 
1148  // Evaluate 1d function values (needed for the product rule)
1149  DynamicMatrix<R> values;
1150 
1151  evaluateFunctionFull(in, values, knotVector, order, currentKnotSpan);
1152 
1153  out.resize(knotVector.size()-order-1);
1154  for (size_t j=0; j<out.size(); j++)
1155  out[j] = values[order][j];
1156 
1157  // Evaluate 1d function values of one order lower (needed for the derivative formula)
1158  std::vector<R> lowOrderOneDValues;
1159 
1160  if (order!=0)
1161  {
1162  lowOrderOneDValues.resize(knotVector.size()-(order-1)-1);
1163  for (size_t j=0; j<lowOrderOneDValues.size(); j++)
1164  lowOrderOneDValues[j] = values[order-1][j];
1165  }
1166 
1167  // Evaluate 1d function values of two order lower (needed for the (second) derivative formula)
1168  std::vector<R> lowOrderTwoDValues;
1169 
1170  if (order>1 && evaluateHessian)
1171  {
1172  lowOrderTwoDValues.resize(knotVector.size()-(order-2)-1);
1173  for (size_t j=0; j<lowOrderTwoDValues.size(); j++)
1174  lowOrderTwoDValues[j] = values[order-2][j];
1175  }
1176 
1177  // Evaluate 1d function derivatives
1178  if (evaluateJacobian)
1179  {
1180  outJac.resize(limit);
1181 
1182  if (order==0) // order-zero functions are piecewise constant, hence all derivatives are zero
1183  std::fill(outJac.begin(), outJac.end(), R(0.0));
1184  else
1185  {
1186  for (size_t j=offset; j<offset+limit; j++)
1187  {
1188  R derivativeAddend1 = lowOrderOneDValues[j] / (knotVector[j+order]-knotVector[j]);
1189  R derivativeAddend2 = lowOrderOneDValues[j+1] / (knotVector[j+order+1]-knotVector[j+1]);
1190  // The two previous terms may evaluate as 0/0. This is to be interpreted as 0.
1191  if (std::isnan(derivativeAddend1))
1192  derivativeAddend1 = 0;
1193  if (std::isnan(derivativeAddend2))
1194  derivativeAddend2 = 0;
1195  outJac[j-offset] = order * ( derivativeAddend1 - derivativeAddend2 );
1196  }
1197  }
1198  }
1199 
1200  // Evaluate 1d function second derivatives
1201  if (evaluateHessian)
1202  {
1203  outHess.resize(limit);
1204 
1205  if (order<2) // order-zero functions are piecewise constant, hence all derivatives are zero
1206  std::fill(outHess.begin(), outHess.end(), R(0.0));
1207  else
1208  {
1209  for (size_t j=offset; j<offset+limit; j++)
1210  {
1211  assert(j+2 < lowOrderTwoDValues.size());
1212  R derivativeAddend1 = lowOrderTwoDValues[j] / (knotVector[j+order]-knotVector[j]) / (knotVector[j+order-1]-knotVector[j]);
1213  R derivativeAddend2 = lowOrderTwoDValues[j+1] / (knotVector[j+order]-knotVector[j]) / (knotVector[j+order]-knotVector[j+1]);
1214  R derivativeAddend3 = lowOrderTwoDValues[j+1] / (knotVector[j+order+1]-knotVector[j+1]) / (knotVector[j+order]-knotVector[j+1]);
1215  R derivativeAddend4 = lowOrderTwoDValues[j+2] / (knotVector[j+order+1]-knotVector[j+1]) / (knotVector[j+1+order]-knotVector[j+2]);
1216  // The two previous terms may evaluate as 0/0. This is to be interpreted as 0.
1217 
1218  if (std::isnan(derivativeAddend1))
1219  derivativeAddend1 = 0;
1220  if (std::isnan(derivativeAddend2))
1221  derivativeAddend2 = 0;
1222  if (std::isnan(derivativeAddend3))
1223  derivativeAddend3 = 0;
1224  if (std::isnan(derivativeAddend4))
1225  derivativeAddend4 = 0;
1226  outHess[j-offset] = order * (order-1) * ( derivativeAddend1 - derivativeAddend2 -derivativeAddend3 + derivativeAddend4 );
1227  }
1228  }
1229  }
1230  }
1231 
1232 
1234  std::array<unsigned int, dim> order_;
1235 
1237  std::array<std::vector<double>, dim> knotVectors_;
1238 
1240  std::array<unsigned,dim> elements_;
1241 
1242  GridView gridView_;
1243 };
1244 
1245 
1246 
1247 template<typename GV>
1248 class BSplineNode :
1249  public LeafBasisNode
1250 {
1251  static const int dim = GV::dimension;
1252 
1253 public:
1254 
1255  using size_type = std::size_t;
1256  using Element = typename GV::template Codim<0>::Entity;
1257  using FiniteElement = BSplineLocalFiniteElement<GV,double>;
1258 
1259  BSplineNode(const BSplinePreBasis<GV>* preBasis) :
1260  preBasis_(preBasis),
1261  finiteElement_(*preBasis)
1262  {}
1263 
1265  const Element& element() const
1266  {
1267  return element_;
1268  }
1269 
1274  const FiniteElement& finiteElement() const
1275  {
1276  return finiteElement_;
1277  }
1278 
1280  void bind(const Element& e)
1281  {
1282  element_ = e;
1283  auto elementIndex = preBasis_->gridView().indexSet().index(e);
1284  finiteElement_.bind(preBasis_->getIJK(elementIndex,preBasis_->elements_));
1285  this->setSize(finiteElement_.size());
1286  }
1287 
1288 protected:
1289 
1290  const BSplinePreBasis<GV>* preBasis_;
1291 
1292  FiniteElement finiteElement_;
1293  Element element_;
1294 };
1295 
1296 
1297 
1298 namespace BasisFactory {
1299 
1306 inline auto bSpline(const std::vector<double>& knotVector,
1307  unsigned int order,
1308  bool makeOpen = true)
1309 {
1310  return [&knotVector, order, makeOpen](const auto& gridView) {
1311  return BSplinePreBasis<std::decay_t<decltype(gridView)>>(gridView, knotVector, order, makeOpen);
1312  };
1313 }
1314 
1315 } // end namespace BasisFactory
1316 
1317 // *****************************************************************************
1318 // This is the actual global basis implementation based on the reusable parts.
1319 // *****************************************************************************
1320 
1327 template<typename GV>
1329 
1330 
1331 } // namespace Functions
1332 
1333 } // namespace Dune
1334 
1335 #endif // DUNE_FUNCTIONS_FUNCTIONSPACEBASES_BSPLINEBASIS_HH
Construct a matrix with a dynamic size.
Definition: dynmatrix.hh:61
void resize(size_type r, size_type c, value_type v=value_type())
resize matrix to r × c
Definition: dynmatrix.hh:106
A dense n x m matrix.
Definition: fmatrix.hh:117
vector space out of a tensor product of fields.
Definition: fvector.hh:95
LocalBasis class in the sense of dune-localfunctions, presenting the restriction of a B-spline patch ...
Definition: bsplinebasis.hh:46
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:55
unsigned int order() const
Polynomial order of the shape functions.
Definition: bsplinebasis.hh:140
std::size_t size() const
Return the number of basis functions on the current knot span.
Definition: bsplinebasis.hh:147
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:97
void evaluateFunction(const FieldVector< D, dim > &in, std::vector< FieldVector< R, 1 > > &out) const
Evaluate all shape functions.
Definition: bsplinebasis.hh:70
void evaluateJacobian(const FieldVector< D, dim > &in, std::vector< FieldMatrix< D, 1, dim > > &out) const
Evaluate Jacobian of all shape functions.
Definition: bsplinebasis.hh:82
BSplineLocalBasis(const BSplinePreBasis< GV > &preBasis, const BSplineLocalFiniteElement< GV, R > &lFE)
Constructor with a given B-spline patch.
Definition: bsplinebasis.hh:61
Attaches a shape function to an entity.
Definition: bsplinebasis.hh:178
const LocalKey & localKey(std::size_t i) const
get i'th index
Definition: bsplinebasis.hh:321
std::size_t size() const
number of coefficients
Definition: bsplinebasis.hh:315
LocalFiniteElement in the sense of dune-localfunctions, for the B-spline basis on tensor-product grid...
Definition: bsplinebasis.hh:362
BSplineLocalFiniteElement(const BSplineLocalFiniteElement &other)
Copy constructor.
Definition: bsplinebasis.hh:383
const BSplineLocalBasis< GV, R > & localBasis() const
Hand out a LocalBasis object.
Definition: bsplinebasis.hh:427
LocalFiniteElementTraits< BSplineLocalBasis< GV, R >, BSplineLocalCoefficients< dim >, BSplineLocalInterpolation< dim, BSplineLocalBasis< GV, R > > > Traits
Export various types related to this LocalFiniteElement.
Definition: bsplinebasis.hh:372
BSplineLocalFiniteElement(const BSplinePreBasis< GV > &preBasis)
Constructor with a given B-spline basis.
Definition: bsplinebasis.hh:376
void bind(const std::array< unsigned, dim > &elementIdx)
Bind LocalFiniteElement to a specific knot span of the spline patch.
Definition: bsplinebasis.hh:394
const BSplineLocalCoefficients< dim > & localCoefficients() const
Hand out a LocalCoefficients object.
Definition: bsplinebasis.hh:433
GeometryType type() const
Return the reference element that the local finite element is defined on (here, a hypercube)
Definition: bsplinebasis.hh:455
unsigned size() const
Number of shape functions in this finite element.
Definition: bsplinebasis.hh:445
unsigned int size(int i) const
Number of degrees of freedom for one coordinate direction.
Definition: bsplinebasis.hh:463
const BSplineLocalInterpolation< dim, BSplineLocalBasis< GV, R > > & localInterpolation() const
Hand out a LocalInterpolation object.
Definition: bsplinebasis.hh:439
Local interpolation in the sense of dune-localfunctions, for the B-spline basis on tensor-product gri...
Definition: bsplinebasis.hh:340
void interpolate(const F &f, std::vector< C > &out) const
Local interpolation of a function.
Definition: bsplinebasis.hh:344
Pre-basis for B-spline basis.
Definition: bsplinebasis.hh:499
std::array< unsigned, dim > elements_
Number of grid elements in the different coordinate directions.
Definition: bsplinebasis.hh:1240
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:1009
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:1085
void evaluateFunction(const FieldVector< typename GV::ctype, dim > &in, std::vector< FieldVector< R, 1 > > &out, const std::array< unsigned, dim > &currentKnotSpan) const
Evaluate all B-spline basis functions at a given point.
Definition: bsplinebasis.hh:781
std::array< unsigned int, dim > order_
Order of the B-spline for each space dimension.
Definition: bsplinebasis.hh:1234
size_type dimension() const
Get the total dimension of the space spanned by this basis.
Definition: bsplinebasis.hh:719
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:1128
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:1028
unsigned int size(size_t d) const
Number of shape functions in one direction.
Definition: bsplinebasis.hh:774
GV GridView
The grid view that the FE space is defined on.
Definition: bsplinebasis.hh:560
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:735
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 > &currentKnotSpan) const
Evaluate Derivatives of all B-spline basis functions.
Definition: bsplinebasis.hh:913
void update(const GridView &gv)
Update the stored grid view, to be called if the grid has changed.
Definition: bsplinebasis.hh:686
void initializeIndices()
Initialize the global indices.
Definition: bsplinebasis.hh:676
size_type size(const Dune::ReservedVector< ST, i > &prefix) const
Return number of possible values for next position in multi index.
Definition: bsplinebasis.hh:712
unsigned int size() const
Total number of B-spline basis functions.
Definition: bsplinebasis.hh:765
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 > &currentKnotSpan) const
Evaluate Jacobian of all B-spline basis functions.
Definition: bsplinebasis.hh:812
std::array< std::vector< double >, dim > knotVectors_
The knot vectors, one for each space dimension.
Definition: bsplinebasis.hh:1237
size_type maxNodeSize() const
Get the maximal number of DOFs associated to node for any element.
Definition: bsplinebasis.hh:725
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
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:125
Describe position of one degree of freedom.
Definition: localkey.hh:23
Default exception for dummy implementations.
Definition: exceptions.hh:263
Default exception class for range errors.
Definition: exceptions.hh:254
A Vector class with statically reserved memory.
Definition: reservedvector.hh:47
constexpr size_type size() const noexcept
Returns number of elements in the vector.
Definition: reservedvector.hh:387
This file implements a quadratic diagonal matrix of fixed size.
This file implements a dense matrix with dynamic numbers of rows and columns.
void umv(const X &x, Y &y) const
y += A x
Definition: diagonalmatrix.hh:297
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
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:1306
constexpr GeometryType cube(unsigned int dim)
Returns a GeometryType representing a hypercube of dimension dim.
Definition: type.hh:472
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:291
auto max(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::max()
Definition: defaults.hh:81
constexpr HybridTreePath< T..., std::size_t > push_back(const HybridTreePath< T... > &tp, std::size_t i)
Appends a run time index to a HybridTreePath.
Definition: treepath.hh:281
Dune namespace.
Definition: alignedallocator.hh:13
Static tag representing a codimension.
Definition: dimension.hh:24
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:34
D DomainType
domain type
Definition: localbasis.hh:42
traits helper struct
Definition: localfiniteelementtraits.hh:13
A unique label for each type of element that can occur in a grid.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (Apr 25, 22:37, 2024)