lagrangebasis.hh 19.5 KB
Newer Older
Carsten Gräser's avatar
Carsten Gräser committed
1 2 3 4 5
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
#ifndef DUNE_FUNCTIONS_FUNCTIONSPACEBASES_LAGRANGEBASIS_HH
#define DUNE_FUNCTIONS_FUNCTIONSPACEBASES_LAGRANGEBASIS_HH

6 7
#include <dune/common/exceptions.hh>

8
#include <dune/functions/common/lagrangelfecache.hh>
9
#include <dune/functions/functionspacebases/nodes.hh>
Carsten Gräser's avatar
Carsten Gräser committed
10 11
#include <dune/functions/functionspacebases/flatmultiindex.hh>
#include <dune/functions/functionspacebases/defaultglobalbasis.hh>
12 13 14
#include <dune/functions/functionspacebases/geometrytypeprovider.hh>


Carsten Gräser's avatar
Carsten Gräser committed
15 16 17 18 19


namespace Dune {
namespace Functions {

20 21 22 23 24 25 26 27 28 29 30 31
// *****************************************************************************
// This is the reusable part of the LagrangeBasis. It contains
//
//   LagrangePreBasis
//   LagrangeNodeIndexSet
//   LagrangeNode
//
// The pre-basis allows to create the others and is the owner of possible shared
// state. These three components do _not_ depend on the global basis or index
// set and can be used without a global basis.
// *****************************************************************************

32
template<typename GV, int k, typename TP, typename GeometryTypeProvider=AutoGeometryTypeProvider>
33 34
class LagrangeNode;

35
template<typename GV, int k, class MI, class TP, typename GeometryTypeProvider=AutoGeometryTypeProvider>
36 37
class LagrangeNodeIndexSet;

38
template<typename GV, int k, class MI, typename GeometryTypeProvider=AutoGeometryTypeProvider>
39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56
class LagrangePreBasis;



/**
 * \brief A pre-basis for a PQ-lagrange bases with given order
 *
 * \ingroup FunctionSpaceBasesImplementations
 *
 * \tparam GV  The grid view that the FE basis is defined on
 * \tparam k   The polynomial order of ansatz functions
 * \tparam MI  Type to be used for multi-indices
 *
 * \note This only works for certain grids.  The following restrictions hold
 * - If k is no larger than 2, then the grids can have any dimension
 * - If k is larger than 3 then the grid must be two-dimensional
 * - If k is 3, then the grid can be 3d *if* it is a simplex grid
 */
57
template<typename GV, int k, class MI, typename GeometryTypeProvider>
58 59 60 61 62 63 64 65 66 67 68 69 70 71
class LagrangePreBasis
{
  static const int dim = GV::dimension;

public:

  //! The grid view that the FE basis is defined on
  using GridView = GV;

  //! Type used for indices and size information
  using size_type = std::size_t;

private:

72
  template<typename, int, class, class, class>
73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96
  friend class LagrangeNodeIndexSet;

  // Precompute the number of dofs per entity type
  const static size_type dofsPerVertex =
      k == 0 ? (dim == 0 ? 1 : 0) : 1;
  const static size_type dofsPerEdge =
      k == 0 ? (dim == 1 ? 1 : 0) : k-1;
  const static size_type dofsPerTriangle =
      k == 0 ? (dim == 2 ? 1 : 0) : (k-1)*(k-2)/2;
  const static size_type dofsPerQuad =
      k == 0 ? (dim == 2 ? 1 : 0) : (k-1)*(k-1);
  const static size_type dofsPerTetrahedron =
      k == 0 ? (dim == 3 ? 1 : 0) : (k-3)*(k-2)*(k-1)/6;
  const static size_type dofsPerPrism =
      k == 0 ? (dim == 3 ? 1 : 0) : (k-1)*(k-1)*(k-2)/2;
  const static size_type dofsPerHexahedron =
      k == 0 ? (dim == 3 ? 1 : 0) : (k-1)*(k-1)*(k-1);
  const static size_type dofsPerPyramid =
      k == 0 ? (dim == 3 ? 1 : 0) : (k-2)*(k-1)*(2*k-3)/6;

public:

  //! Template mapping root tree path to type of created tree node
  template<class TP>
97
  using Node = LagrangeNode<GV, k, TP, GeometryTypeProvider>;
98 99 100

  //! Template mapping root tree path to type of created tree node index set
  template<class TP>
101
  using IndexSet = LagrangeNodeIndexSet<GV, k, MI, TP, GeometryTypeProvider>;
102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119

  //! Type used for global numbering of the basis vectors
  using MultiIndex = MI;

  //! Type used for prefixes handed to the size() method
  using SizePrefix = Dune::ReservedVector<size_type, 1>;

  //! Constructor for a given grid view object
  LagrangePreBasis(const GridView& gv) :
    gridView_(gv)
  {}

  //! Initialize the global indices
  void initializeIndices()
  {
    vertexOffset_        = 0;
    edgeOffset_          = vertexOffset_          + dofsPerVertex * ((size_type)gridView_.size(dim));

120 121 122 123 124 125
    if (dim>=2)
    {
      triangleOffset_      = edgeOffset_            + dofsPerEdge * ((size_type) gridView_.size(dim-1));

      quadrilateralOffset_ = triangleOffset_        + dofsPerTriangle * ((size_type)gridView_.size(Dune::GeometryTypes::triangle));
    }
126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245

    if (dim==3) {
      tetrahedronOffset_   = quadrilateralOffset_ + dofsPerQuad * ((size_type)gridView_.size(Dune::GeometryTypes::quadrilateral));

      prismOffset_         = tetrahedronOffset_   +   dofsPerTetrahedron * ((size_type)gridView_.size(Dune::GeometryTypes::tetrahedron));

      hexahedronOffset_    = prismOffset_         +   dofsPerPrism * ((size_type)gridView_.size(Dune::GeometryTypes::prism));

      pyramidOffset_       = hexahedronOffset_    +   dofsPerHexahedron * ((size_type)gridView_.size(Dune::GeometryTypes::hexahedron));
    }
  }

  //! Obtain the grid view that the basis is defined on
  const GridView& gridView() const
  {
    return gridView_;
  }

  //! Update the stored grid view, to be called if the grid has changed
  void update (const GridView& gv)
  {
    gridView_ = gv;
  }

  /**
   * \brief Create tree node with given root tree path
   *
   * \tparam TP Type of root tree path
   * \param tp Root tree path
   *
   * By passing a non-trivial root tree path this can be used
   * to create a node suitable for being placed in a tree at
   * the position specified by the root tree path.
   */
  template<class TP>
  Node<TP> node(const TP& tp) const
  {
    return Node<TP>{tp};
  }

  /**
   * \brief Create tree node index set with given root tree path
   *
   * \tparam TP Type of root tree path
   * \param tp Root tree path
   *
   * Create an index set suitable for the tree node obtained
   * by node(tp).
   */
  template<class TP>
  IndexSet<TP> indexSet() const
  {
    return IndexSet<TP>{*this};
  }

  //! Same as size(prefix) with empty prefix
  size_type size() const
  {
    switch (dim)
    {
      case 1:
        return dofsPerVertex * ((size_type)gridView_.size(dim))
          + dofsPerEdge*((size_type)gridView_.size(dim-1));
      case 2:
      {
        return dofsPerVertex * ((size_type)gridView_.size(dim))
          + dofsPerEdge * ((size_type)gridView_.size(dim-1))
          + dofsPerTriangle * ((size_type)gridView_.size(Dune::GeometryTypes::triangle))
          + dofsPerQuad * ((size_type)gridView_.size(Dune::GeometryTypes::quadrilateral));
      }
      case 3:
      {
        return dofsPerVertex * ((size_type)gridView_.size(dim))
          + dofsPerEdge * ((size_type)gridView_.size(dim-1))
          + dofsPerTriangle * ((size_type)gridView_.size(Dune::GeometryTypes::triangle))
          + dofsPerQuad * ((size_type)gridView_.size(Dune::GeometryTypes::quadrilateral))
          + dofsPerTetrahedron * ((size_type)gridView_.size(Dune::GeometryTypes::tetrahedron))
          + dofsPerPyramid * ((size_type)gridView_.size(Dune::GeometryTypes::pyramid))
          + dofsPerPrism * ((size_type)gridView_.size(Dune::GeometryTypes::prism))
          + dofsPerHexahedron * ((size_type)gridView_.size(Dune::GeometryTypes::hexahedron));
      }
    }
    DUNE_THROW(Dune::NotImplemented, "No size method for " << dim << "d grids available yet!");
  }

  //! Return number of possible values for next position in multi index
  size_type size(const SizePrefix prefix) const
  {
    assert(prefix.size() == 0 || prefix.size() == 1);
    return (prefix.size() == 0) ? size() : 0;
  }

  //! Get the total dimension of the space spanned by this basis
  size_type dimension() const
  {
    return size();
  }

  //! Get the maximal number of DOFs associated to node for any element
  size_type maxNodeSize() const
  {
    return StaticPower<(k+1),GV::dimension>::power;
  }

protected:
  GridView gridView_;

  size_type vertexOffset_;
  size_type edgeOffset_;
  size_type triangleOffset_;
  size_type quadrilateralOffset_;
  size_type tetrahedronOffset_;
  size_type pyramidOffset_;
  size_type prismOffset_;
  size_type hexahedronOffset_;

};



246
template<typename GV, int k, typename TP, typename GeometryTypeProvider>
247 248 249 250 251 252 253
class LagrangeNode :
  public LeafBasisNode<std::size_t, TP>
{
  static const int dim = GV::dimension;
  static const int maxSize = StaticPower<(k+1),GV::dimension>::power;

  using Base = LeafBasisNode<std::size_t,TP>;
254 255

  using FiniteElementCache = typename Dune::Functions::LagrangeFiniteElementCache<typename GV::ctype, double, dim, k>;
256 257 258 259 260 261

public:

  using size_type = std::size_t;
  using TreePath = TP;
  using Element = typename GV::template Codim<0>::Entity;
262
  using FiniteElement = typename FiniteElementCache::template FiniteElement<typename GeometryTypeProvider::template Type<GV>>;
263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288

  LagrangeNode(const TreePath& treePath) :
    Base(treePath),
    finiteElement_(nullptr),
    element_(nullptr)
  {}

  //! Return current element, throw if unbound
  const Element& element() const
  {
    return *element_;
  }

  /** \brief Return the LocalFiniteElement for the element we are bound to
   *
   * The LocalFiniteElement implements the corresponding interfaces of the dune-localfunctions module
   */
  const FiniteElement& finiteElement() const
  {
    return *finiteElement_;
  }

  //! Bind to element.
  void bind(const Element& e)
  {
    element_ = &e;
289 290
    const auto& geometryType = GeometryTypeProvider::template type<GV>(*element_);
    finiteElement_ = &(cache_.get(geometryType));
291 292 293 294 295 296 297 298 299 300 301 302
    this->setSize(finiteElement_->size());
  }

protected:

  FiniteElementCache cache_;
  const FiniteElement* finiteElement_;
  const Element* element_;
};



303
template<typename GV, int k, class MI, class TP, typename GeometryTypeProvider>
304 305 306 307 308 309 310 311 312 313 314
class LagrangeNodeIndexSet
{
  enum {dim = GV::dimension};

public:

  using size_type = std::size_t;

  /** \brief Type used for global numbering of the basis vectors */
  using MultiIndex = MI;

315
  using PreBasis = LagrangePreBasis<GV, k, MI, GeometryTypeProvider>;
316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358

  using Node = typename PreBasis::template Node<TP>;

  LagrangeNodeIndexSet(const PreBasis& preBasis) :
    preBasis_(&preBasis),
    node_(nullptr)
  {}

  /** \brief Bind the view to a grid element
   *
   * Having to bind the view to an element before being able to actually access any of its data members
   * offers to centralize some expensive setup code in the 'bind' method, which can save a lot of run-time.
   */
  void bind(const Node& node)
  {
    node_ = &node;
  }

  /** \brief Unbind the view
   */
  void unbind()
  {
    node_ = nullptr;
  }

  /** \brief Size of subtree rooted in this node (element-local)
   */
  size_type size() const
  {
    assert(node_ != nullptr);
    return node_->finiteElement().size();
  }

  //! Maps from subtree index set [0..size-1] to a globally unique multi index in global basis
  template<typename It>
  It indices(It it) const
  {
    assert(node_ != nullptr);
    for (size_type i = 0, end = node_->finiteElement().size() ; i < end ; ++it, ++i)
      {
        Dune::LocalKey localKey = node_->finiteElement().localCoefficients().localKey(i);
        const auto& gridIndexSet = preBasis_->gridView().indexSet();
        const auto& element = node_->element();
359
        const auto& geometryType = GeometryTypeProvider::template type<GV>(element);
360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385

        // The dimension of the entity that the current dof is related to
        auto dofDim = dim - localKey.codim();

        // Test for a vertex dof
        // The test for k==1 is redundant, but having it here allows the compiler to conclude
        // at compile-time that the dofDim==0 case is the only one that will ever happen.
        // This leads to measurable speed-up: see
        //   https://gitlab.dune-project.org/staging/dune-functions/issues/30
        if (k==1 || dofDim==0) {
          *it = {{ (size_type)(gridIndexSet.subIndex(element,localKey.subEntity(),dim)) }};
          continue;
        }

        if (dofDim==1)
          {  // edge dof
            if (dim==1)  // element dof -- any local numbering is fine
              {
                *it = {{ preBasis_->edgeOffset_
                         + preBasis_->dofsPerEdge * ((size_type)gridIndexSet.subIndex(element,0,0))
                         + localKey.index() }};
                continue;
              }
            else
              {
                const auto refElement
386
                  = Dune::referenceElement<double,dim>(geometryType);
387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407

                // we have to reverse the numbering if the local triangle edge is
                // not aligned with the global edge
                auto v0 = (size_type)gridIndexSet.subIndex(element,refElement.subEntity(localKey.subEntity(),localKey.codim(),0,dim),dim);
                auto v1 = (size_type)gridIndexSet.subIndex(element,refElement.subEntity(localKey.subEntity(),localKey.codim(),1,dim),dim);
                bool flip = (v0 > v1);
                *it = {{ (flip)
                         ? preBasis_->edgeOffset_
                         + preBasis_->dofsPerEdge*((size_type)gridIndexSet.subIndex(element,localKey.subEntity(),localKey.codim()))
                         + (preBasis_->dofsPerEdge-1)-localKey.index()
                         : preBasis_->edgeOffset_
                         + preBasis_->dofsPerEdge*((size_type)gridIndexSet.subIndex(element,localKey.subEntity(),localKey.codim()))
                         + localKey.index() }};
                continue;
              }
          }

        if (dofDim==2)
          {
            if (dim==2)   // element dof -- any local numbering is fine
              {
408
                if (geometryType.isTriangle())
409 410 411 412 413
                  {
                    const int interiorLagrangeNodesPerTriangle = (k-1)*(k-2)/2;
                    *it = {{ preBasis_->triangleOffset_ + interiorLagrangeNodesPerTriangle*((size_type)gridIndexSet.subIndex(element,0,0)) + localKey.index() }};
                    continue;
                  }
414
                else if (geometryType.isQuadrilateral())
415 416 417 418 419 420 421 422 423 424
                  {
                    const int interiorLagrangeNodesPerQuadrilateral = (k-1)*(k-1);
                    *it = {{ preBasis_->quadrilateralOffset_ + interiorLagrangeNodesPerQuadrilateral*((size_type)gridIndexSet.subIndex(element,0,0)) + localKey.index() }};
                    continue;
                  }
                else
                  DUNE_THROW(Dune::NotImplemented, "2d elements have to be triangles or quadrilaterals");
              } else
              {
                const auto refElement
425
                  = Dune::referenceElement<double,dim>(geometryType);
426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441

                if (k>3)
                  DUNE_THROW(Dune::NotImplemented, "LagrangeNodalBasis for 3D grids is only implemented if k<=3");

                if (k==3 and !refElement.type(localKey.subEntity(), localKey.codim()).isTriangle())
                  DUNE_THROW(Dune::NotImplemented, "LagrangeNodalBasis for 3D grids with k==3 is only implemented if the grid is a simplex grid");

                *it = {{ preBasis_->triangleOffset_ + ((size_type)gridIndexSet.subIndex(element,localKey.subEntity(),localKey.codim())) }};
                continue;
              }
          }

        if (dofDim==3)
          {
            if (dim==3)   // element dof -- any local numbering is fine
              {
442
                if (geometryType.isTetrahedron())
443 444 445 446
                  {
                    *it = {{ preBasis_->tetrahedronOffset_ + PreBasis::dofsPerTetrahedron*((size_type)gridIndexSet.subIndex(element,0,0)) + localKey.index() }};
                    continue;
                  }
447
                else if (geometryType.isHexahedron())
448 449 450 451
                  {
                    *it = {{ preBasis_->hexahedronOffset_ + PreBasis::dofsPerHexahedron*((size_type)gridIndexSet.subIndex(element,0,0)) + localKey.index() }};
                    continue;
                  }
452
                else if (geometryType.isPrism())
453 454 455 456
                  {
                    *it = {{ preBasis_->prismOffset_ + PreBasis::dofsPerPrism*((size_type)gridIndexSet.subIndex(element,0,0)) + localKey.index() }};
                    continue;
                  }
457
                else if (geometryType.isPyramid())
458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477
                  {
                    *it = {{ preBasis_->pyramidOffset_ + PreBasis::dofsPerPyramid*((size_type)gridIndexSet.subIndex(element,0,0)) + localKey.index() }};
                    continue;
                  }
                else
                  DUNE_THROW(Dune::NotImplemented, "3d elements have to be tetrahedra, hexahedra, prisms, or pyramids");
              } else
              DUNE_THROW(Dune::NotImplemented, "Grids of dimension larger than 3 are no supported");
          }
        DUNE_THROW(Dune::NotImplemented, "Grid contains elements not supported for the LagrangeNodalBasis");
      }
    return it;
  }

protected:
  const PreBasis* preBasis_;

  const Node* node_;
};

Carsten Gräser's avatar
Carsten Gräser committed
478 479


480
namespace BasisFactory {
Carsten Gräser's avatar
Carsten Gräser committed
481

482 483
namespace Imp {

484
template<std::size_t k, class GeometryTypeProvider>
485 486 487 488 489 490 491 492
class LagrangePreBasisFactory
{
public:
  static const std::size_t requiredMultiIndexSize = 1;

  template<class MultiIndex, class GridView>
  auto makePreBasis(const GridView& gridView) const
  {
493
    return LagrangePreBasis<GridView, k, MultiIndex, GeometryTypeProvider>(gridView);
494 495 496 497 498 499 500 501
  }

};

} // end namespace BasisFactory::Imp



502
/**
503
 * \brief Create a pre-basis factory that can create a  Lagrange pre-basis
504
 *
505 506 507 508 509 510 511 512 513 514 515
 * The GeometryTypeProvider allows to explicitly fix how geometry types
 * are derived. The default version is AutoGeometryTypeProvider which
 * uses a fixed finite element implementation if the mesh is statically
 * known to only contain a single geometry type. Otherwise is uses
 * the polymorphic VirtualFiniteElement interface.
 *
 * If your grid type allows for mixed geometry types but you know,
 * that your grid only contains a single one, then you can explicitly
 * provide a SimplexGeometryTypeProvider or CubeGeometryTypeProvider
 * to avoid the polymorphic interface.
 *
516 517 518
 * \ingroup FunctionSpaceBasesImplementations
 *
 * \tparam k   The polynomial order of ansatz functions
519
 * \tparam GeometryTypeProvider Helper struct to determine geometry types
520
 */
521
template<std::size_t k, class GeometryTypeProvider=AutoGeometryTypeProvider>
522
auto lagrange()
Carsten Gräser's avatar
Carsten Gräser committed
523
{
524
  return Imp::LagrangePreBasisFactory<k, GeometryTypeProvider>();
Carsten Gräser's avatar
Carsten Gräser committed
525 526
}

527
} // end namespace BasisFactory
Carsten Gräser's avatar
Carsten Gräser committed
528 529 530 531



/** \brief Nodal basis of a scalar k-th-order Lagrangean finite element space
532 533
 *
 * \ingroup FunctionSpaceBasesImplementations
Carsten Gräser's avatar
Carsten Gräser committed
534 535 536 537 538 539
 *
 * \note This only works for certain grids.  The following restrictions hold
 * - If k is no larger than 2, then the grids can have any dimension
 * - If k is larger than 3 then the grid must be two-dimensional
 * - If k is 3, then the grid can be 3d *if* it is a simplex grid
 *
540
 * All arguments passed to the constructor will be forwarded to the constructor
541
 * of LagrangePreBasis.
542
 *
543 544 545 546 547 548 549 550 551 552 553
 * The GeometryTypeProvider allows to explicitly fix how geometry types
 * are derived. The default version is AutoGeometryTypeProvider which
 * uses a fixed finite element implementation if the mesh is statically
 * known to only contain a single geometry type. Otherwise is uses
 * the polymorphic VirtualFiniteElement interface.
 *
 * If your grid type allows for mixed geometry types but you know,
 * that your grid only contains a single one, then you can explicitly
 * provide a SimplexGeometryTypeProvider or CubeGeometryTypeProvider
 * to avoid the polymorphic interface.
 *
Carsten Gräser's avatar
Carsten Gräser committed
554 555
 * \tparam GV The GridView that the space is defined on
 * \tparam k The order of the basis
556
 * \tparam GeometryTypeProvider Helper struct to determine geometry types
Carsten Gräser's avatar
Carsten Gräser committed
557
 */
558 559
template<typename GV, int k, class GeometryTypeProvider=AutoGeometryTypeProvider>
using LagrangeBasis = DefaultGlobalBasis<LagrangePreBasis<GV, k, FlatMultiIndex<std::size_t>, GeometryTypeProvider> >;
Carsten Gräser's avatar
Carsten Gräser committed
560 561 562 563 564 565 566 567 568





} // end namespace Functions
} // end namespace Dune


Carsten Gräser's avatar
Carsten Gräser committed
569
#endif // DUNE_FUNCTIONS_FUNCTIONSPACEBASES_LAGRANGEBASIS_HH