hp2D::Space Class Referenceabstract
A 2D hp FEM space with continuous, piecewise polynomial basis functions. More...
#include <space.hh>
Public Types | |
typedef concepts::Scan< hp2D::Element< Real > > | Scan |
typedef Scan< ElementWithCell< Real > > | Scanner |
typedef ElementWithCell< Real > | type |
Public Member Functions | |
virtual void | adjust (const concepts::Element< Real > &elm, const concepts::AdaptiveAdjustP< 2 > &a) |
virtual void | adjust (const Element< Real > &elm, const concepts::AdaptiveAdjustP< 2 > &a)=0 |
Adjusts the space in the next rebuild step for this element. More... | |
Scan * | constScan () const |
Returns constant scanner over elements of last build, In contrast to scan(), this method also workes, if the space would need a rebuild, * i.e. More... | |
uint | dim () |
Returns the dimension of the space. More... | |
virtual uint | dim () const |
virtual uint | getOutputDimension () const |
Returns the default output dimension, when we consider plotting a real-valued operator on this space. More... | |
uint | nelm () |
Returns the number of elements in the space. More... | |
virtual uint | nelm () const |
Space | rebuild () |
Rebuilds the space. More... | |
Scan * | scan () |
Returns a scanner to iterate over the elements of the space. More... | |
virtual Scan * | scan () const |
Space (concepts::Mesh2 &msh, uint l, uint p, concepts::BoundaryConditions *bc=0) | |
Constructor. More... | |
Space (const Space &spc) | |
Copy constructor. More... | |
virtual | ~Space () |
Friends | |
class | BuildDofsBase |
Strategies to Build the Degrees of Freedom | |
concepts::Mesh2 & | msh_ |
Mesh. More... | |
concepts::BoundaryConditions * | bc_ |
Boundary conditions. More... | |
bool | rebuild_ |
If true: the elements have to be rebuilt. More... | |
uint | dim_ |
Dimension of the FE space. More... | |
uint | nelm_ |
Number of elements currently active in the mesh. More... | |
concepts::Joiner< Element< Real > *, 1 > * | elm_ |
Linked list of the elements. More... | |
std::unordered_map< uint, concepts::AdaptiveControl<> > | ctrl0_ |
Hash table of the control information for the vertices. More... | |
std::unique_ptr< std::unordered_map< uint, concepts::AdaptiveControlP< 1 > > > | ctrl1_ |
Hash table of the control information for the edges. More... | |
std::unordered_map< uint, concepts::AdaptiveControlP< 2 > > | ctrl2_ |
Hash table of the control information for the 2D elements. More... | |
std::unordered_map< uint, concepts::AdaptiveControl<> > | ctrlc0_ |
Hash table for connected vertex dofs (mapping from attribute) More... | |
std::unordered_map< uint, const concepts::Attribute * > | vtxConnect_ |
Hash table for connected vertex dofs (mapping from cell key) More... | |
std::map< uint, concepts::CellData > | cellList_ |
List of cells. More... | |
std::map< uint, concepts::EdgeData > | edgeList_ |
List of edges. More... | |
std::map< uint, concepts::VertexData > | vertexList_ |
List of vertices. More... | |
std::unordered_map< uint, concepts::AdaptiveAdjustP< 2 > > | adj_ |
Hash table of the adjustment information of the elements. More... | |
std::unique_ptr< concepts::SMatrix1D > | S1left_ |
S matrices in 1D. More... | |
std::unique_ptr< concepts::SMatrix1D > | S1right_ |
Mesh. More... | |
std::unique_ptr< concepts::SMatrixBase< Real > > | Smatrices2H_ [2] |
S matrices for horizontal subdivision. More... | |
std::unique_ptr< concepts::SMatrixBase< Real > > | Smatrices2V_ [2] |
S matrices for vertical subdivision. More... | |
std::unique_ptr< concepts::SMatrixBase< Real > > | Smatrices4_ [4] |
S matrices for subdivision into 4 quads. More... | |
std::unique_ptr< BuildDofsBase > | buildVertexDofs_ |
Strategy to build the vertex degrees of freedom. More... | |
std::unique_ptr< BuildDofsBase > | buildEdgeDofs_ |
Strategy to build the edge degrees of freedom. More... | |
std::unique_ptr< BuildDofsBase > | buildInnerDofs_ |
Strategy to build the inner degrees of freedom. More... | |
void | buildVertexDofs (const BuildDofsBase &b) |
Change the strategy how the degrees of freedom for the vertices are built. More... | |
void | buildEdgeDofs (const BuildDofsBase &b) |
Change the strategy how the degrees of freedom for the edge are built. More... | |
void | buildInnerDofs (const BuildDofsBase &b) |
Change the strategy how the degrees of freedom for the interior are built. More... | |
uint | connectedIdx (const concepts::Attribute &a) const |
Gives the number of connected dofs for one attribute. More... | |
void | recomputeShapefunctions () |
Recompute shape functions, e.g. More... | |
virtual std::ostream & | info (std::ostream &os) const |
Mesh. More... | |
void | meshAndPoly_ (concepts::Quad2d &cell, concepts::Level< 2 > &L, int *P) |
Adjusts the mesh and the polynomial degree. More... | |
void | enforceBC_ (const concepts::Connector2 &cell) |
Enforce Dirichlet boundary conditions. More... | |
void | createCellList_ (const concepts::Connector &cntr, const concepts::CellData *father) |
Initially fills cellList_ (recusively) More... | |
void | createVtxEdgList_ (const concepts::Connector2 &cell) |
Initially fills edgeListVtx_ and vertexList_ . More... | |
void | enrichElm_ (const concepts::Quad2d &cell, int *p) |
Enriches the polynomial degrees. More... | |
void | minimumRule_ (const concepts::Connector2 &cntr, const concepts::Connector1 &edge) |
Used to enforce the minimum rule. More... | |
void | buildElements_ (concepts::Quad2d &cell, ushort *Pmax, concepts::TColumn< Real > *T0=0) |
Creates the elements and their T matrices. More... | |
void | buildEmptyElements_ (concepts::Quad2d &cell) |
Mesh. More... | |
void | computePmax_ (const concepts::Connector2 &cntr, ushort Pmax[2]) const |
Computes the maximal polynomial degree of cell in each direction. More... | |
void | recomputeSmatrices_ (const Quad< Real > &elm) |
Checks if the S matrices need to be recomputed and does so if necessary. More... | |
void | deactivate_ (const concepts::Edge &edg) |
Mesh. More... | |
void | deactivate_ (const concepts::Vertex &vtx) |
Mesh. More... | |
bool | activeCell_ (concepts::Cell &cell) |
Checks if a cell is in active region. More... | |
void | status_ () |
Gives the cells and edges with given polynomial degree (for debugging) More... | |
void | matlabstatus_ () |
Mesh. More... | |
Space (const Space &spc, uint i) | |
Copy constructor used to return a snapshot of the old space from rebuild() More... | |
Detailed Description
A 2D hp FEM space with continuous, piecewise polynomial basis functions.
Currently, only hexahedrons and conforming meshes are possible.
- Adaptivity
- You can adaptively refine this space by using adjust() on some elements to refine them or increase their polynomial degree. Please note that these extensions can be locally varying (ie. non-uniform) and anisotropic.
- Trunk Space
- Classically, the whole tensor product polynomial space is used on each element to build the degress of freedom. However, a considerable amount of degrees of freedom can be saved, when only a trunk of this space is used. In fact, only very little is lost in accuracy while the gain in terms of degrees of freedom is noteworthy. The methods buildVertexDofs(), buildEdgeDofs() and buildInnerDofs() can be used to change the way the degrees of freedom on each part of each element are built.
- Examples
- BGT_0.cc, hpFEM2d-simple.cc, and hpFEM2d.cc.
Member Typedef Documentation
◆ Scan
typedef concepts::Scan<hp2D::Element<Real> > hp2D::Space::Scan |
◆ Scanner
|
inherited |
◆ type
|
inherited |
Constructor & Destructor Documentation
◆ Space() [1/3]
hp2D::Space::Space | ( | concepts::Mesh2 & | msh, |
uint | l, | ||
uint | p, | ||
concepts::BoundaryConditions * | bc = 0 |
||
) |
Constructor.
Scans the mesh and sets the cells in the mesh active and the level of refinement and the polynomial degree in all cells to the given values. rebuild_ is set to true, ie. if the mesh is used it will firstly be rebuilt.
- Parameters
-
msh The domain of interest partitioned into a mesh. l Level of refinement p Degree of the polynomials to be used. bc Boundary conditions
◆ Space() [2/3]
hp2D::Space::Space | ( | const Space & | spc | ) |
Copy constructor.
Copies mesh, boundary conditions and ctrl2_
(information about refinements and polynomial degrees of the elements).
◆ ~Space()
|
virtual |
◆ Space() [3/3]
|
private |
Copy constructor used to return a snapshot of the old space from rebuild()
- Parameters
-
spc To be copied i Dummy argument to select this copy constructor
Member Function Documentation
◆ activeCell_()
|
private |
Checks if a cell is in active region.
A whole region can be deactivated by setting the attribute of the cells in the region to be of Dirichlet boundary condition. In this region empty elements. Thats necessary for vectorial spaces, which live on the same mesh, so that the spaces in all dimension have the same number of elements.
The edges of the nonactive regions could be active, that means that the shape functions lie only on one side of the edge. The edges of the region can be deactivated by usual setting of Dirichlet boundary condition.
◆ adjust() [1/2]
|
virtual |
◆ adjust() [2/2]
|
pure virtualinherited |
Adjusts the space in the next rebuild step for this element.
◆ buildEdgeDofs()
void hp2D::Space::buildEdgeDofs | ( | const BuildDofsBase & | b | ) |
Change the strategy how the degrees of freedom for the edge are built.
The default strategy is BuildEdgeDofs. You can change this strategy any time you chose.
buildEdgeDofs_ is reset with a clone of b
.
- Parameters
-
b New strategy
◆ buildElements_()
|
private |
Creates the elements and their T matrices.
This routine is called recursively until the finest level is reached. If there have degrees of freedom to be assigned on higher level, the according T columns are adapted using the S matrices. On the finest level, the elements are created and the T matrix is created from the T columns.
- Parameters
-
cell The current cell Pmax Request for maximal polynomial degree, array with 2 entries T0 Precomputed T columns which belong to a parent of the cell
- Precondition
Pmax
is an array of 2 integers
◆ buildEmptyElements_()
|
private |
Mesh.
◆ buildInnerDofs()
void hp2D::Space::buildInnerDofs | ( | const BuildDofsBase & | b | ) |
Change the strategy how the degrees of freedom for the interior are built.
The default strategy is BuildInnerDofsLinTrunk. You can change this strategy any time you chose.
buildInnerDofs_ is reset with a clone of b
.
- Parameters
-
b New strategy
◆ buildVertexDofs()
void hp2D::Space::buildVertexDofs | ( | const BuildDofsBase & | b | ) |
Change the strategy how the degrees of freedom for the vertices are built.
The default strategy is BuildVertexDofs. You can change this strategy any time you chose.
buildVertexDofs_ is reset with a clone of b
.
- Parameters
-
b New strategy
◆ computePmax_()
|
private |
Computes the maximal polynomial degree of cell
in each direction.
Takes into account higher polynomial degrees on faces and edges compared to the interior.
◆ connectedIdx()
uint hp2D::Space::connectedIdx | ( | const concepts::Attribute & | a | ) | const |
Gives the number of connected dofs for one attribute.
◆ constScan()
|
inline |
◆ createCellList_()
|
private |
Initially fills cellList_
(recusively)
◆ createVtxEdgList_()
|
private |
Initially fills edgeListVtx_
and vertexList_
.
◆ deactivate_() [1/2]
|
private |
Mesh.
◆ deactivate_() [2/2]
|
private |
Mesh.
◆ dim() [1/2]
|
inlinevirtual |
Returns the dimension of the space.
Implements concepts::SpaceOnCells< Real >.
◆ dim() [2/2]
|
inlinevirtual |
- Examples
- BGT_0.cc, and hpFEM2d.cc.
◆ enforceBC_()
|
private |
Enforce Dirichlet boundary conditions.
This routine is called by meshAndPoly_
on every cell.
- Parameters
-
cell Cell
◆ enrichElm_()
|
private |
Enriches the polynomial degrees.
This enrichement is necessary for the edges to be able to create conforming global basis functions.
In refined cells, some of the small edges might have a polynomial degree which is too low to represent a part of a basis function which is defined on a larger element.
To achieve this, the cells are traversed top down until the level of the elements is reached. During the traversal, the needed polynomial degree for the small elements is computed and enforced on the element level (only on the edges).
- Precondition
p
is an array of 6 integers
◆ getOutputDimension()
|
inlinevirtualinherited |
◆ info()
|
protectedvirtual |
Mesh.
Reimplemented from concepts::SpaceOnCells< Real >.
◆ matlabstatus_()
|
private |
Mesh.
◆ meshAndPoly_()
|
private |
Adjusts the mesh and the polynomial degree.
On the old mesh, the adjustment information is taken. This leads to mesh refinements or coarsening and a certain polynomial degree distribution.
The polynomial degree is fixed on the finest level (where the elements are eventually built) and given to the coarser levels (for edges and the interior of the cell). This is done bottom up.
Finally, the boundary conditions are enforced on every cell using enforceBC_
- Parameters
-
cell Cell L The desired level P The desired polynomial degree
- Precondition
P
is an array of 3 integers
- See also
- enforceBC_
◆ minimumRule_()
|
private |
Used to enforce the minimum rule.
The minimum rule states that the polynomial degree on a edge has to be the minimum of the adjacent elements.
This routine recursively calls itself and asserts that all the elements adjacent to edge
enforce their polynomial degree on edge
. This is done if the cell
is member of the space. Otherwise, minimumRule_
is called with the children of cell
.
◆ nelm() [1/2]
|
inlinevirtual |
Returns the number of elements in the space.
Implements concepts::SpaceOnCells< Real >.
◆ nelm() [2/2]
◆ rebuild()
Space hp2D::Space::rebuild | ( | ) |
Rebuilds the space.
First the old list of Finite Elements in the space is removed (not the elements in the topology!). Then the control information for the vertices, edges and faces is cleared (the information on the 2D elements is retained). For every cell in the original mesh, rebuild0_ and rebuild1_ are called. From the mesh, they build the space and store the elements in elm_. An essential part of the generated information are the T matrices. They describe how the local shape functions are glued together to form the global shape functions.
rebuild0_ marks the elements with the right values and rebuild1_ really does the work.
- See also
- TMatrix for more information on T matrices and how they are stored.
- Returns
- Copy of the space on the old level
- Examples
- BGT_0.cc, hpFEM2d-simple.cc, and hpFEM2d.cc.
◆ recomputeShapefunctions()
void hp2D::Space::recomputeShapefunctions | ( | ) |
Recompute shape functions, e.g.
for other abscissas redefined through setIntegrationRule
- Examples
- BGT_0.cc.
◆ recomputeSmatrices_()
|
private |
Checks if the S matrices need to be recomputed and does so if necessary.
- Parameters
-
elm Element which the S matrices are needed for (the polynomial degrees are important)
◆ scan() [1/2]
|
inlinevirtual |
Returns a scanner to iterate over the elements of the space.
Implements concepts::SpaceOnCells< Real >.
◆ scan() [2/2]
|
inlinevirtual |
◆ status_()
|
private |
Gives the cells and edges with given polynomial degree (for debugging)
Friends And Related Function Documentation
◆ BuildDofsBase
|
friend |
Member Data Documentation
◆ adj_
|
private |
◆ bc_
|
private |
◆ buildEdgeDofs_
|
private |
◆ buildInnerDofs_
|
private |
◆ buildVertexDofs_
|
private |
◆ cellList_
|
private |
◆ ctrl0_
|
private |
◆ ctrl1_
|
private |
◆ ctrl2_
|
private |
◆ ctrlc0_
|
private |
◆ dim_
◆ edgeList_
|
private |
◆ elm_
|
private |
◆ msh_
|
private |
◆ nelm_
|
private |
◆ rebuild_
|
private |
◆ S1left_
|
private |
◆ S1right_
|
private |
◆ Smatrices2H_
|
private |
S matrices for horizontal subdivision.
- See also
- concepts::Quad2dSubdiv2H
◆ Smatrices2V_
|
private |
S matrices for vertical subdivision.
- See also
- concepts::Quad2dSubdiv2V
◆ Smatrices4_
|
private |
S matrices for subdivision into 4 quads.
- See also
- concepts::Quad2dSubdiv4
◆ vertexList_
|
private |
◆ vtxConnect_
|
private |
The documentation for this class was generated from the following file:
- hp2D/space.hh