hp1D::Jump1Jump1< F > Class Template Referenceabstract
A function class to calculate element matrices for the Product of the jumps of the derivative. More...
#include <bilinearForm.hh>
Public Member Functions | |
virtual Jump1Jump1< F > * | clone () const |
Virtual constructor. More... | |
virtual BilinearForm * | clone () const=0 |
Virtual constructor. More... | |
const concepts::ElementFormula< Real > * | formula () const |
Getter for formula. More... | |
const concepts::Array< Real > | getFrmVal (const hp1D::Element< Real > &elm) const |
Returns the values of frm on {0,1}. More... | |
const concepts::Array< Real > | getInvJacobeanDet (const hp1D::Element< Real > &elm) const |
Returns the inverse of the jacobean determinant (i.J.d.) More... | |
const concepts::Array< Real > | getInvJacobeanDetSq (const hp1D::Element< Real > &elm) const |
Returns the inverse of the jacobean determinant squared (i.J.d.s.) More... | |
const concepts::Array< Real > | getTangentialSecondDerivative (const hp1D::Element< Real > &elm) const |
Returns the second tangential derivative (s.t.d.) of the inverse Parametrisation on {0,1}. More... | |
Jump1Jump1 (const concepts::ElementFormulaContainer< F > frm=concepts::ElementFormulaContainer< F >()) | |
Constructor. More... | |
void | operator() (const concepts::Element< Real > &elmX, const concepts::Element< Real > &elmY, concepts::ElementMatrix< F > &em) const |
Operator. More... | |
virtual void | operator() (const Element< Real > &elmX, const Element< Real > &elmY, ElementMatrix< Real > &em) const=0 |
Evaluates the bilinear form for all shape functions on elmX and elmY and stores the result in the matrix em . More... | |
virtual void | operator() (const Element< Real > &elmX, const Element< Real > &elmY, ElementMatrix< Real > &em, const ElementPair< Real > &ep) const |
Evaluates the bilinear form for all shape functions on elmX and elmY and stores the result in the matrix em . More... | |
virtual | ~Jump1Jump1 () |
Destructor. More... | |
Protected Member Functions | |
virtual std::ostream & | info (std::ostream &os) const |
Returns information in an output stream. More... | |
Protected Attributes | |
const concepts::ElementFormulaContainer< Real > | frm_ |
ElementFormula. More... | |
Detailed Description
template<class F = Real>
class hp1D::Jump1Jump1< F >
A function class to calculate element matrices for the Product of the jumps of the derivative.
Definition at line 165 of file bilinearForm.hh.
Constructor & Destructor Documentation
◆ Jump1Jump1()
|
inline |
Constructor.
The formula frm
is evaluated in each quadrature point.
Definition at line 172 of file bilinearForm.hh.
◆ ~Jump1Jump1()
|
virtual |
Destructor.
Member Function Documentation
◆ clone() [1/2]
|
inlinevirtual |
Virtual constructor.
Returns a pointer to a copy of itself. The caller is responsible to destroy this copy.
Implements concepts::Cloneable.
Definition at line 184 of file bilinearForm.hh.
◆ clone() [2/2]
|
pure virtualinherited |
Virtual constructor.
Returns a pointer to a copy of itself. The caller is responsible to destroy this copy.
◆ formula()
|
inlineinherited |
Getter for formula.
Definition at line 272 of file formula.hh.
◆ getFrmVal()
|
inherited |
Returns the values of frm on {0,1}.
- Returns
- Array of size 2, with the values frm on {0,1}
◆ getInvJacobeanDet()
|
inherited |
Returns the inverse of the jacobean determinant (i.J.d.)
- Returns
- Array of size 2, with the values i.J.d. on {0,1}
◆ getInvJacobeanDetSq()
|
inherited |
Returns the inverse of the jacobean determinant squared (i.J.d.s.)
- Returns
- Array of size 2, with the values of i.J.d.s. on {0,1}
◆ getTangentialSecondDerivative()
|
inherited |
Returns the second tangential derivative (s.t.d.) of the inverse Parametrisation on {0,1}.
- Returns
- Array of size 2, with the values of the s.t.d. on {0,1}
◆ info()
|
protectedvirtual |
Returns information in an output stream.
Reimplemented from concepts::BilinearForm< Real, Real >.
◆ operator()() [1/3]
void hp1D::Jump1Jump1< F >::operator() | ( | const concepts::Element< Real > & | elmX, |
const concepts::Element< Real > & | elmY, | ||
concepts::ElementMatrix< F > & | em | ||
) | const |
Operator.
◆ operator()() [2/3]
|
pure virtualinherited |
Evaluates the bilinear form for all shape functions on elmX
and elmY
and stores the result in the matrix em
.
- Postcondition
- The returned matrix
em
has the correct size.
- Parameters
-
elmX Left element (test functions) elmY Right element (trial functions) em Return element matrix
◆ operator()() [3/3]
|
inlinevirtualinherited |
Evaluates the bilinear form for all shape functions on elmX
and elmY
and stores the result in the matrix em
.
If this method is not reimplemented in a derived class, the default behaviour is to call the application operator without ep
.
- Postcondition
- The returned matrix
em
has the correct size.
- Parameters
-
elmX Left element elmY Right element em Return element matrix ep Element pair holding more information on the pair elmX
andelmY
Definition at line 57 of file bilinearForm.hh.
Member Data Documentation
◆ frm_
|
protectedinherited |
ElementFormula.
Definition at line 308 of file formula.hh.
The documentation for this class was generated from the following file:
- hp1D/bilinearForm.hh