#include <quad.hh>
Public Types | |
typedef F | type |
Public Member Functions | |
InfiniteQuad (concepts::InfiniteQuad2d &cell, const ushort *p, concepts::TColumn< Real > *T0, concepts::TColumn< Real > *T1) | |
const concepts::Connector2 & | support () const |
concepts::Real3d | vertex (uint i) const |
const concepts::InfiniteQuad2d & | cell () const |
Returns the cell on which the element is built. | |
const ushort * | p () const |
const concepts::QuadratureRule1d * | integrationX () const |
Returns the integration rule in x direction. | |
const concepts::QuadratureRule1d * | integrationY () const |
Returns the integration rule in y direction. | |
const concepts::Karniadakis< 1, 0 > * | shpfctX () const |
Returns the shape functions in x direction. | |
const concepts::Karniadakis< 1, 1 > * | shpfctDX () const |
Returns the derivatives of the shape functions in x direction. | |
virtual const concepts::LaguerreBasis< 0 > * | shpfctY () const =0 |
Returns the shape functions in x direction. | |
virtual const concepts::LaguerreBasis< 1 > * | shpfctDY () const =0 |
Returns the derivatives of the shape functions in x direction. | |
virtual const concepts::ElementGraphics< Real > * | graphics () const =0 |
virtual const concepts::TMatrix< Real > & | T () const |
Returns the T matrix of the element. | |
void | appendT (concepts::TColumn< Real > *T) |
Appends the T columns to the T matrix. | |
Real3d | elemMap (const Real coord_local) const |
Real3d | elemMap (const Real2d &coord_local) const |
Real3d | elemMap (const Real3d &coord_local) const |
uint & | tag () |
Returns the tag. | |
Static Public Member Functions | |
static void | setIntegrationRuleX (enum concepts::intRule rule, bool constant, uint points) |
static void | setIntegrationRuleY (bool constant, uint points, const Real border) |
static void | resetIntegrationRule () |
Set the standard type of integration. | |
static std::string | integrationRule () |
Returns information on the settings of the quadrature rule. | |
static Real & | borderY () |
Largest value in local y coordinates. Default is 1.0. | |
Protected Member Functions | |
void | computeIntegrationRuleFromP_ (const ushort *p) |
void | computeIntegrationRule_ (const uint nq[2]) |
void | computeShapefunctionsX_ () |
gets the shapefunctions in x direction | |
virtual std::ostream & | info (std::ostream &os) const |
Returns information in an output stream. | |
Protected Attributes | |
ushort | p_ [2] |
Polynomial degree. | |
std::unique_ptr< concepts::QuadratureRule1d > | intX_ |
The integration rules. | |
std::unique_ptr< concepts::QuadratureRule1d > | intY_ |
concepts::TMatrix< Real > | T_ |
T matrix of the element. | |
|
inherited |
Definition at line 81 of file element.hh.
hp2D::InfiniteQuad::InfiniteQuad | ( | concepts::InfiniteQuad2d & | cell, |
const ushort * | p, | ||
concepts::TColumn< Real > * | T0, | ||
concepts::TColumn< Real > * | T1 | ||
) |
Constructor
cell | Cell on which the element is defined |
p | Polynomial degree (might be anisotropic) |
T0 | Part of the T matrix |
T1 | Part of the T matrix |
|
inlineinherited |
Appends the T columns to the T matrix.
Definition at line 43 of file element.hh.
|
inlinestatic |
|
inlinevirtual |
Returns the cell on which the element is built.
Implements hp2D::Element< Real >.
|
inlineinherited |
Definition at line 86 of file element.hh.
|
inlineinherited |
Definition at line 90 of file element.hh.
|
inlineinherited |
Definition at line 94 of file element.hh.
|
pure virtual |
Reimplemented from concepts::Element< F >.
|
protectedvirtual |
Returns information in an output stream.
Reimplemented from hp2D::Element< Real >.
Reimplemented in hp2D::InfiniteLaguerreQuad.
|
inline |
|
inline |
|
inline |
|
static |
Sets the type of integration rule in x direction to use.
rule | Determines the integration rule |
constant | Use constant number of points |
points | Number of points to use (for constant = true ) or points to add to approximation order (for constant = false ) |
|
static |
Sets the type of integration rule in y direction to use.
constant | Use constant number of points |
points | Number of points to use (for constant = true ) or points to add to approximation order (for constant = false ) |
border | Largest value in local y coordinates. |
|
inline |
|
pure virtual |
Returns the derivatives of the shape functions in x direction.
Implemented in hp2D::InfiniteLaguerreQuad.
|
inline |
|
pure virtual |
Returns the shape functions in x direction.
Implemented in hp2D::InfiniteLaguerreQuad.
|
inlinevirtual |
Returns the topological support of the element. Possible supports for an element are quadrilaterals and triangles.
Implements hp2D::Element< Real >.
|
inlinevirtualinherited |
Returns the T matrix of the element.
Implements concepts::ElementWithCell< F >.
Definition at line 40 of file element.hh.
|
inlineinherited |
Returns the tag.
Definition at line 66 of file element.hh.
|
inlinevirtual |
Returns the coordinates of the ith vertex of this element.
i | Index of the vertex |
Implements hp2D::Element< Real >.
|
protected |
|
protected |
|
protected |
|
protectedinherited |
T matrix of the element.
Definition at line 57 of file element.hh.