Dune Core Modules (2.11.0)
#include <dune/vtk/gridcreators/lagrangegridcreator.hh>
Public Member Functions | |
| void | insertVerticesImpl (std::vector< GlobalCoordinate > const &points, std::vector< std::uint64_t > const &) |
Implementation of the interface function insertVertices() | |
| void | insertElementsImpl (std::vector< std::uint8_t > const &types, std::vector< std::int64_t > const &offsets, std::vector< std::int64_t > const &connectivity) |
Implementation of the interface function insertElements() | |
| LocalParametrization | localParametrization (unsigned int insertionIndex) const |
| Construct an element parametrization. More... | |
| LocalParametrization | localParametrization (Element const &element) const |
| Construct an element parametrization. More... | |
| LocalGeometry | localGeometry (Element const &element) const |
| Construct a transformation of local element coordinates. More... | |
| int | order (GeometryType type, std::size_t nNodes) const |
| Determine lagrange order from number of points. | |
| int | order () const |
| Determine lagrange order from number of points from the first element parametrization. | |
| EntitySet | entitySet () const |
| Dummy function returning a placeholder entityset. | |
| GlobalCoordinate | operator() (GlobalCoordinate const &) const |
| Dummy function returning a placeholder entityset. | |
| GridFactory< Grid > & | factory () |
| Return the associated GridFactory. | |
| GridFactory< Grid > const & | factory () const |
| Return the associated (const) GridFactory. | |
| void | insertVertices (std::vector< GlobalCoordinate > const &points, std::vector< std::uint64_t > const &point_ids) |
| Insert all points as vertices into the factory. | |
| void | insertElements (std::vector< std::uint8_t > const &types, std::vector< std::int64_t > const &offsets, std::vector< std::int64_t > const &connectivity) |
| Create elements based on type and connectivity description. | |
| void | insertPieces (std::vector< std::string > const &pieces) |
| Insert part of a grid stored in file into factory. | |
| std::unique_ptr< Grid > | createGrid () const |
| Construct the actual grid using the GridFactory. | |
| auto | comm () const |
| Return the mpi collective communicator. | |
Friends | |
| LocalFunction | localFunction (LagrangeGridCreator &gridCreator) |
| Local function representing the parametrization of the grid. More... | |
Detailed Description
class Dune::Vtk::LagrangeGridCreator< GridType >
The grid is created from the first nodes of a cell parametrization, representing the corner vertices. Thus a piecewise "flat" grid is constructed. The parametrization is 1. passed as a local element parametrization to the insertElement() function of a gridFactory to allow the grid itself to handle the parametrization and 2. is stored internally that can be accessed by using this GridCreator object as a grid function, or by extracting locally the parametrization on each existing grid element after creation of the grid.
So, the LagrangeGridCreator models both, a GridCreator and a GridFunction.
Member Function Documentation
◆ localGeometry()
|
inline |
Construct a transformation of local element coordinates.
An element might have a different local coordinate system than the coordinate system used to defined the element parametrization. Thus coordinate transform of the local parametrization is needed for element-local evaluations. This local geometry transform is obtained by figuring out the permutation of corners in the element corresponding to the inserted corner vertices.
References Dune::Vtk::LagrangeGridCreator< GridType >::factory(), Dune::GridFactoryInterface< GridType >::insertionIndex(), Dune::Vtk::LagrangeGridCreator< GridType >::localGeometry(), and VTK_ASSERT.
Referenced by Dune::Vtk::LagrangeGridCreator< GridType >::localGeometry(), and Dune::Vtk::LagrangeGridCreator< GridType >::localParametrization().
◆ localParametrization() [1/2]
|
inline |
Construct an element parametrization.
The returned LocalParametrization is a mapping GlobalCoordinate(LocalCoordinate) where LocalCoordinate is w.r.t. the local coordinate system in the passed element andGlobalCoordinate` a world coordinate in the parametrized grid.
References Dune::Vtk::LagrangeGridCreator< GridType >::factory(), Dune::GridFactoryInterface< GridType >::insertionIndex(), Dune::Vtk::LagrangeGridCreator< GridType >::localGeometry(), Dune::Vtk::LagrangeGridCreator< GridType >::order(), and VTK_ASSERT.
◆ localParametrization() [2/2]
|
inline |
Construct an element parametrization.
The returned LocalParametrization is a mapping GlobalCoordinate(LocalCoordinate) where LocalCoordinate is w.r.t. the local coordinate system in an element with giveninsertionIndex<tt>(defined by the inserted corner vertices) and GlobalCoordinate` a world coordinate in the parametrized grid.
References Dune::Vtk::LagrangeGridCreator< GridType >::order(), and VTK_ASSERT.
Referenced by Dune::Vtk::LagrangeGridCreator< GridType >::insertElementsImpl().
Friends And Related Function Documentation
◆ localFunction
|
friend |
Local function representing the parametrization of the grid.
The returned object models Functions::Concept::LocalFunction and can thus be bound to an element of the created grid and evaluated in the local coordinates of the bound element.
It is implemented in terms of the LocalParametrization function returned by the method localParametrization(element). See comments there for further details.
Note, this methods requires the GridCreator to be based by lvalue-reference. This is necessary, since we want to guarantee that all internal storage is preserved while evaluating the local function.
The documentation for this class was generated from the following file:
- dune/vtk/gridcreators/lagrangegridcreator.hh
|
Legal Statements / Impressum |
Hosted by TU Dresden & Uni Heidelberg |
generated with Hugo v0.111.3
(Feb 14, 23:39, 2026)