A 2D hp FEM space of edge elements with piecewise polynomial basis functions with two components, where the tangential component is continuous Currently, only quadrilaterials and conforming meshes are possible. More...

#include <space.hh>

Inheritance diagram for hp2Dedge::Space:
concepts::SpaceOnCells< Real > concepts::AdaptiveSpace< Real, concepts::AdaptiveAdjustP< 2 > > concepts::Space< Real >

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...
 
concepts::BoundaryConditionsboundary ()
 
uint dim ()
 Returns the dimension of the space. More...
 
uint dim () const
 
virtual uint getOutputDimension () const
 
virtual uint getOutputDimension () const
 Returns the default output dimension, when we consider plotting a real-valued operator on this space. More...
 
concepts::Scan< const concepts::Quad > * interScan ()
 scanner over the intermediate cells More...
 
concepts::Scan< const concepts::Quad > * interScan () const
 scanner over the intermediate cells More...
 
uint nelm ()
 Returns the number of elements in the space. More...
 
uint nelm () const
 
void rebuild ()
 Rebuilds the space. More...
 
void recomputeShapefunctions ()
 Recompute shape functions, e.g. More...
 
Scanscan ()
 scanner over the cells of the coarsest mesh More...
 
Scanscan () const
 scanner over the cells of the coarsest mesh More...
 
 Space (concepts::Mesh2 &msh, uint l, uint p, concepts::BoundaryConditions *bc=0, bool trunk=false)
 Constructor. More...
 
virtual ~Space ()
 

Protected Member Functions

virtual std::ostream & info (std::ostream &os) const
 

Private Member Functions

bool activeCell_ (concepts::Cell &cell)
 Checks if a cell is in active region. More...
 
concepts::TColumn< Real > * buildEdge_ (concepts::TColumn< Real > *T1, const concepts::Quad &cntr, const uint *nn)
 Builds the T columns for the edge degrees of freedom. More...
 
void buildElements_ (concepts::Quad2d &cell)
 Creates the elements and their T matrices. More...
 
concepts::TColumn< Real > * buildInterior_ (concepts::TColumn< Real > *T1, const concepts::Quad &cntr, const uint *nn)
 Builds the T columns for the interior degrees of freedom. More...
 
void buildTColumns_ (const concepts::Quad &cntr, ushort *Pmax, concepts::TColumn< Real > *T0=0)
 Builds the T columns, stores them together with the neccessary polynomial degrees in a list T_. More...
 
void computePmax_ (const concepts::Connector2 &cntr, ushort Pmax[2]) const
 Computes the maximal polynomial degree of cell in each direction. More...
 
void createCellList_ (const concepts::Connector &cntr, const concepts::CellData *father)
 Initially fills cellList_ (recursively) More...
 
void createEdgList_ (const concepts::Cell2 &cell)
 Initially fills edgeList_. More...
 
void deactivate_ (concepts::Edge &edg)
 
void deactivate_ (concepts::Vertex &vtx)
 
void edgeOrientation_ (const concepts::Connector2 &cntr, const uint i, uint &pIndex) const
 Helps distributing the polynomial degree to the edges. More...
 
void enforceBC_ (concepts::Quad2d &cell)
 Enforce Dirichlet boundary conditions. More...
 
void enrichElm_ (const concepts::Quad2d &cell, int *p)
 Enriches the polynomial degrees. More...
 
void matlabstatus_ ()
 
void meshAndPoly_ (concepts::Quad2d &cell, concepts::Level< 2 > &L, int *P)
 Adjusts the mesh and the polynomial degree. More...
 
void minimumRule_ (const concepts::Connector2 &cntr, const concepts::Connector1 &edge)
 Used to enforce the minimum rule. More...
 
void nonactiveRegion_ ()
 Set that cells in nonactive region to passive. More...
 
void recomputeSmatrices_ (const ushort *p, const uint *nn)
 Checks if the S matrices need to be recomputed and does so if necessary. More...
 
void status_ ()
 Gives the cells and edges with given polynomial degree (for debugging) More...
 

Private Attributes

concepts::DynArray< concepts::AdaptiveAdjustP< 2 > > adj_
 Hash table of the adjustment information of the elements. More...
 
concepts::BoundaryConditionsbc_
 Boundary conditions. More...
 
std::map< uint, concepts::CellDatacellList_
 List of cells. More...
 
concepts::DynArray< concepts::AdaptiveControl<> > ctrl0_
 Hash table of the control information for the vertices. More...
 
concepts::DynArray< 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...
 
uint dim_
 Dimension of the FE space. More...
 
std::map< uint, concepts::EdgeDataedgeList_
 List of edges. More...
 
concepts::Joiner< hp2D::Element< Real > *, 1 > * elm_
 Linked list of the elements. More...
 
concepts::Joiner< const concepts::Quad *, 1 > * interCells_
 
std::set< uint > interEdges_
 
concepts::Mesh2msh_
 Mesh. More...
 
uint nelm_
 Number of elements currently active in the mesh. More...
 
bool rebuild_
 If true: the elements have to be rebuilt. More...
 
std::unique_ptr< concepts::SMatrix1DS1left_n_
 
std::unique_ptr< concepts::SMatrix1DS1left_t_
 S matrices in 1D. More...
 
std::unique_ptr< concepts::SMatrix1DS1right_n_
 
std::unique_ptr< concepts::SMatrix1DS1right_t_
 
std::unique_ptr< concepts::SMatrixBase< Real > > Smatrices2H_ [2]
 S matrices for horizontal subdivision. More...
 
std::unique_ptr< concepts::SMatrixBase< Real > > Smatrices2H_X_ [2]
 S matrices for tangential and normal component for different subdivisions. More...
 
std::unique_ptr< concepts::SMatrixBase< Real > > Smatrices2H_Y_ [2]
 
std::unique_ptr< concepts::SMatrixBase< Real > > Smatrices2V_ [2]
 S matrices for vertical subdivision. More...
 
std::unique_ptr< concepts::SMatrixBase< Real > > Smatrices2V_X_ [2]
 
std::unique_ptr< concepts::SMatrixBase< Real > > Smatrices2V_Y_ [2]
 
std::unique_ptr< concepts::SMatrixBase< Real > > Smatrices4_ [4]
 S matrices for subdivision into 4 quads. More...
 
concepts::TColumnSet< Real, 2 > T_
 List of TColumns (for the member of the space cells) More...
 
bool trunk_
 Build trunk space. More...
 

Detailed Description

A 2D hp FEM space of edge elements with piecewise polynomial basis functions with two components, where the tangential component is continuous Currently, only quadrilaterials and conforming meshes are possible.

Author
Kersten Schmidt, 2004

Definition at line 50 of file space.hh.

Member Typedef Documentation

◆ Scan

Definition at line 54 of file space.hh.

◆ Scanner

typedef Scan<ElementWithCell<Real > > concepts::SpaceOnCells< Real >::Scanner
inherited

Definition at line 84 of file space.hh.

◆ type

typedef ElementWithCell<Real > concepts::SpaceOnCells< Real >::type
inherited

Definition at line 83 of file space.hh.

Constructor & Destructor Documentation

◆ Space()

hp2Dedge::Space::Space ( concepts::Mesh2 msh,
uint  l,
uint  p,
concepts::BoundaryConditions bc = 0,
bool  trunk = false 
)

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
mshThe domain of interest partitioned into a mesh.
lLevel of refinement
pDegree of the polynomials to be used.
bcBoundary conditions
trunkCreate trunk space ie., use the the rules p+q <= pIntMax + 1 in the interior instead of allowing all p, q combinations.

◆ ~Space()

virtual hp2Dedge::Space::~Space ( )
virtual

Member Function Documentation

◆ activeCell_()

bool hp2Dedge::Space::activeCell_ ( concepts::Cell cell)
private

Checks if a cell is in active region.

◆ adjust() [1/2]

virtual void hp2Dedge::Space::adjust ( const concepts::Element< Real > &  elm,
const concepts::AdaptiveAdjustP< 2 > &  a 
)
virtual

◆ adjust() [2/2]

virtual void concepts::AdaptiveSpace< Real , concepts::AdaptiveAdjustP< 2 > >::adjust ( const Element< F > &  elm,
const Tadj &  a 
)
pure virtualinherited

Adjusts the space in the next rebuild step for this element.

◆ boundary()

concepts::BoundaryConditions* hp2Dedge::Space::boundary ( )
inline

Definition at line 95 of file space.hh.

◆ buildEdge_()

concepts::TColumn<Real>* hp2Dedge::Space::buildEdge_ ( concepts::TColumn< Real > *  T1,
const concepts::Quad cntr,
const uint *  nn 
)
private

Builds the T columns for the edge degrees of freedom.

◆ buildElements_()

void hp2Dedge::Space::buildElements_ ( concepts::Quad2d cell)
private

Creates the elements and their T matrices.

This routine is called recursively until the finest level is reached. With the TColumns and polynomial degrees, which are taken from T_, the elements are created and the T matrix is created on the finest level.

◆ buildInterior_()

concepts::TColumn<Real>* hp2Dedge::Space::buildInterior_ ( concepts::TColumn< Real > *  T1,
const concepts::Quad cntr,
const uint *  nn 
)
private

Builds the T columns for the interior degrees of freedom.

◆ buildTColumns_()

void hp2Dedge::Space::buildTColumns_ ( const concepts::Quad cntr,
ushort *  Pmax,
concepts::TColumn< Real > *  T0 = 0 
)
private

Builds the T columns, stores them together with the neccessary polynomial degrees in a list T_.

  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 T columns and polynomial degrees are stored.

Parameters
cntrThe current cell (geometric)
PmaxRequest for maximal polynomial degree, array with 2 entries
T0Precomputed T columns which belong to a parent of the cell
Precondition
Pmax is an array of 2 integers

◆ computePmax_()

void hp2Dedge::Space::computePmax_ ( const concepts::Connector2 cntr,
ushort  Pmax[2] 
) const
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.

◆ createCellList_()

void hp2Dedge::Space::createCellList_ ( const concepts::Connector cntr,
const concepts::CellData father 
)
private

Initially fills cellList_ (recursively)

◆ createEdgList_()

void hp2Dedge::Space::createEdgList_ ( const concepts::Cell2 cell)
private

Initially fills edgeList_.

◆ deactivate_() [1/2]

void hp2Dedge::Space::deactivate_ ( concepts::Edge edg)
private

◆ deactivate_() [2/2]

void hp2Dedge::Space::deactivate_ ( concepts::Vertex vtx)
private

◆ dim() [1/2]

uint hp2Dedge::Space::dim ( )
inlinevirtual

Returns the dimension of the space.

Implements concepts::SpaceOnCells< Real >.

Definition at line 366 of file space.hh.

◆ dim() [2/2]

uint hp2Dedge::Space::dim ( ) const
inline

Definition at line 346 of file space.hh.

◆ edgeOrientation_()

void hp2Dedge::Space::edgeOrientation_ ( const concepts::Connector2 cntr,
const uint  i,
uint &  pIndex 
) const
private

Helps distributing the polynomial degree to the edges.

Parameters
cntrConnector
iIndex of the face [0..5]
pIndexReturn value, the index of the polynomial degree of this edge in the array of the element

◆ enforceBC_()

void hp2Dedge::Space::enforceBC_ ( concepts::Quad2d cell)
private

Enforce Dirichlet boundary conditions.

This routine is called by meshAndPoly_ on every cell.

Parameters
cellCell

◆ enrichElm_()

void hp2Dedge::Space::enrichElm_ ( const concepts::Quad2d cell,
int *  p 
)
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() [1/2]

virtual uint hp2Dedge::Space::getOutputDimension ( ) const
inlinevirtual

Definition at line 74 of file space.hh.

◆ getOutputDimension() [2/2]

virtual uint concepts::Space< Real >::getOutputDimension ( ) const
inlinevirtualinherited

Returns the default output dimension, when we consider plotting a real-valued operator on this space.

Definition at line 50 of file space.hh.

◆ info()

virtual std::ostream& hp2Dedge::Space::info ( std::ostream &  os) const
protectedvirtual

Reimplemented from concepts::SpaceOnCells< Real >.

◆ interScan() [1/2]

concepts::Scan< const concepts::Quad > * hp2Dedge::Space::interScan ( )
inline

scanner over the intermediate cells

Definition at line 381 of file space.hh.

◆ interScan() [2/2]

concepts::Scan< const concepts::Quad > * hp2Dedge::Space::interScan ( ) const
inline

scanner over the intermediate cells

Definition at line 361 of file space.hh.

◆ matlabstatus_()

void hp2Dedge::Space::matlabstatus_ ( )
private

◆ meshAndPoly_()

void hp2Dedge::Space::meshAndPoly_ ( concepts::Quad2d cell,
concepts::Level< 2 > &  L,
int *  P 
)
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_

  @param cell Cell
Parameters
LThe desired level
PThe desired polynomial degree
Precondition
P is an array of 3 integers
See also
enforceBC_

◆ minimumRule_()

void hp2Dedge::Space::minimumRule_ ( const concepts::Connector2 cntr,
const concepts::Connector1 edge 
)
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]

uint hp2Dedge::Space::nelm ( )
inlinevirtual

Returns the number of elements in the space.

Implements concepts::SpaceOnCells< Real >.

Definition at line 371 of file space.hh.

◆ nelm() [2/2]

uint hp2Dedge::Space::nelm ( ) const
inline

Definition at line 351 of file space.hh.

◆ nonactiveRegion_()

void hp2Dedge::Space::nonactiveRegion_ ( )
private

Set that cells in nonactive region to passive.

A whole region can be deactivated by setting the attribute of the cells in the region to be of Dirichlet boundary condition. That not includes the edges of the region.

If the edges of the region can be deactivated by usual setting of Dirichlet boundary condition.

◆ rebuild()

void hp2Dedge::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.

◆ recomputeShapefunctions()

void hp2Dedge::Space::recomputeShapefunctions ( )

Recompute shape functions, e.g.

for other abscissas redefined through setIntegrationRule

◆ recomputeSmatrices_()

void hp2Dedge::Space::recomputeSmatrices_ ( const ushort *  p,
const uint *  nn 
)
private

Checks if the S matrices need to be recomputed and does so if necessary.

Parameters
ppolynomial degrees of the element
nnmaximal polynomial degree (in x- and y-direction) for the two vector components, i.e. array with 4 entries, 1st and 2nd are for first vector component, 3rd and 4th are for the second vector component.

◆ scan() [1/2]

Space::Scan * hp2Dedge::Space::scan ( )
inlinevirtual

scanner over the cells of the coarsest mesh

Implements concepts::SpaceOnCells< Real >.

Definition at line 376 of file space.hh.

◆ scan() [2/2]

Space::Scan * hp2Dedge::Space::scan ( ) const
inline

scanner over the cells of the coarsest mesh

Definition at line 356 of file space.hh.

◆ status_()

void hp2Dedge::Space::status_ ( )
private

Gives the cells and edges with given polynomial degree (for debugging)

Member Data Documentation

◆ adj_

concepts::DynArray<concepts::AdaptiveAdjustP<2> > hp2Dedge::Space::adj_
private

Hash table of the adjustment information of the elements.

Definition at line 169 of file space.hh.

◆ bc_

concepts::BoundaryConditions* hp2Dedge::Space::bc_
private

Boundary conditions.

Definition at line 124 of file space.hh.

◆ cellList_

std::map<uint, concepts::CellData> hp2Dedge::Space::cellList_
private

List of cells.

Definition at line 151 of file space.hh.

◆ ctrl0_

concepts::DynArray<concepts::AdaptiveControl<> > hp2Dedge::Space::ctrl0_
private

Hash table of the control information for the vertices.

Definition at line 144 of file space.hh.

◆ ctrl1_

concepts::DynArray<concepts::AdaptiveControlP<1> > hp2Dedge::Space::ctrl1_
private

Hash table of the control information for the edges.

Definition at line 146 of file space.hh.

◆ ctrl2_

std::unordered_map<uint, concepts::AdaptiveControlP<2> > hp2Dedge::Space::ctrl2_
private

Hash table of the control information for the 2D elements.

Definition at line 148 of file space.hh.

◆ dim_

uint hp2Dedge::Space::dim_
private

Dimension of the FE space.

Definition at line 135 of file space.hh.

◆ edgeList_

std::map<uint, concepts::EdgeData> hp2Dedge::Space::edgeList_
private

List of edges.

Definition at line 153 of file space.hh.

◆ elm_

concepts::Joiner<hp2D::Element<Real>*, 1>* hp2Dedge::Space::elm_
private

Linked list of the elements.

Definition at line 141 of file space.hh.

◆ interCells_

concepts::Joiner<const concepts::Quad*, 1>* hp2Dedge::Space::interCells_
private

Definition at line 161 of file space.hh.

◆ interEdges_

std::set<uint> hp2Dedge::Space::interEdges_
private

Definition at line 166 of file space.hh.

◆ msh_

concepts::Mesh2& hp2Dedge::Space::msh_
private

Mesh.

Definition at line 122 of file space.hh.

◆ nelm_

uint hp2Dedge::Space::nelm_
private

Number of elements currently active in the mesh.

Definition at line 138 of file space.hh.

◆ rebuild_

bool hp2Dedge::Space::rebuild_
private

If true: the elements have to be rebuilt.

Some elements have been marked for refinement or for coarsening and the elements are outdated.

Definition at line 132 of file space.hh.

◆ S1left_n_

std::unique_ptr<concepts::SMatrix1D> hp2Dedge::Space::S1left_n_
private

Definition at line 173 of file space.hh.

◆ S1left_t_

std::unique_ptr<concepts::SMatrix1D> hp2Dedge::Space::S1left_t_
private

S matrices in 1D.

Definition at line 173 of file space.hh.

◆ S1right_n_

std::unique_ptr<concepts::SMatrix1D> hp2Dedge::Space::S1right_n_
private

Definition at line 173 of file space.hh.

◆ S1right_t_

std::unique_ptr<concepts::SMatrix1D> hp2Dedge::Space::S1right_t_
private

Definition at line 173 of file space.hh.

◆ Smatrices2H_

std::unique_ptr<concepts::SMatrixBase<Real> > hp2Dedge::Space::Smatrices2H_[2]
private

S matrices for horizontal subdivision.

See also
concepts::Quad2dSubdiv2H

Definition at line 187 of file space.hh.

◆ Smatrices2H_X_

std::unique_ptr<concepts::SMatrixBase<Real> > hp2Dedge::Space::Smatrices2H_X_[2]
private

S matrices for tangential and normal component for different subdivisions.

Definition at line 179 of file space.hh.

◆ Smatrices2H_Y_

std::unique_ptr<concepts::SMatrixBase<Real> > hp2Dedge::Space::Smatrices2H_Y_[2]
private

Definition at line 180 of file space.hh.

◆ Smatrices2V_

std::unique_ptr<concepts::SMatrixBase<Real> > hp2Dedge::Space::Smatrices2V_[2]
private

S matrices for vertical subdivision.

See also
concepts::Quad2dSubdiv2V

Definition at line 191 of file space.hh.

◆ Smatrices2V_X_

std::unique_ptr<concepts::SMatrixBase<Real> > hp2Dedge::Space::Smatrices2V_X_[2]
private

Definition at line 181 of file space.hh.

◆ Smatrices2V_Y_

std::unique_ptr<concepts::SMatrixBase<Real> > hp2Dedge::Space::Smatrices2V_Y_[2]
private

Definition at line 182 of file space.hh.

◆ Smatrices4_

std::unique_ptr<concepts::SMatrixBase<Real> > hp2Dedge::Space::Smatrices4_[4]
private

S matrices for subdivision into 4 quads.

See also
concepts::Quad2dSubdiv4

Definition at line 195 of file space.hh.

◆ T_

concepts::TColumnSet<Real,2> hp2Dedge::Space::T_
private

List of TColumns (for the member of the space cells)

Definition at line 156 of file space.hh.

◆ trunk_

bool hp2Dedge::Space::trunk_
private

Build trunk space.

Definition at line 126 of file space.hh.


The documentation for this class was generated from the following file:
Page URL: http://wiki.math.ethz.ch/bin/view/Concepts/WebHome
21 August 2020
© 2020 Eidgenössische Technische Hochschule Zürich