hp1D::IdentityParallel< F > Class Template Referenceabstract
A function class to calculate element matrices for the mass matrix for elements that are shifted. More...
#include <bilinearForm.hh>
Public Member Functions | |
virtual IdentityParallel< F > * | clone () const |
Virtual constructor. More... | |
virtual BilinearForm * | clone () const=0 |
Virtual constructor. More... | |
IdentityParallel (const concepts::ElementFormulaContainer< F > frm=concepts::ElementFormulaContainer< F >(), const concepts::Real3d shiftX=concepts::Real3d(), const concepts::Real3d shiftY=concepts::Real3d(), const bool negShift=false) | |
Constructor. More... | |
void | operator() (const concepts::Element< Real > &elmX, const concepts::Element< Real > &elmY, concepts::ElementMatrix< F > &em) const |
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 | ~IdentityParallel () |
Protected Member Functions | |
virtual std::ostream & | info (std::ostream &os) const |
Returns information in an output stream. More... | |
Private Attributes | |
const bool | negShift_ |
const concepts::Point< Real, 3 > | shiftX_ |
const concepts::Point< Real, 3 > | shiftY_ |
Detailed Description
template<class F = Real>
class hp1D::IdentityParallel< F >
A function class to calculate element matrices for the mass matrix for elements that are shifted.
There are two independent coordinate shifts for the trial and test function. In contrast to hp1D::Identity the elements do not need to have overlapping support but their shifted start and end points have to be identical. Note that the orientation of the two elements may also be of reverse direction. If negShift
is set to true it is checked as well whether a coordinate shift with opposite sign matches the criterion.
Definition at line 93 of file bilinearForm.hh.
Constructor & Destructor Documentation
◆ IdentityParallel()
hp1D::IdentityParallel< F >::IdentityParallel | ( | const concepts::ElementFormulaContainer< F > | frm = concepts::ElementFormulaContainer< F >() , |
const concepts::Real3d | shiftX = concepts::Real3d() , |
||
const concepts::Real3d | shiftY = concepts::Real3d() , |
||
const bool | negShift = false |
||
) |
Constructor.
- Parameters
-
frm Formula of the bilinear form shiftX Coordinate shift of the trial function shiftY Coordinate shift of the test function negShift Boolean indicating whether shift can be understood also in negative direction
◆ ~IdentityParallel()
|
virtual |
Member Function Documentation
◆ clone() [1/2]
|
virtual |
Virtual constructor.
Returns a pointer to a copy of itself. The caller is responsible to destroy this copy.
Implements concepts::Cloneable.
◆ clone() [2/2]
|
pure virtualinherited |
Virtual constructor.
Returns a pointer to a copy of itself. The caller is responsible to destroy this copy.
◆ info()
|
protectedvirtual |
Returns information in an output stream.
Reimplemented from concepts::BilinearForm< Real, Real >.
◆ operator()() [1/3]
void hp1D::IdentityParallel< F >::operator() | ( | const concepts::Element< Real > & | elmX, |
const concepts::Element< Real > & | elmY, | ||
concepts::ElementMatrix< F > & | em | ||
) | const |
◆ 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
◆ negShift_
|
private |
Definition at line 122 of file bilinearForm.hh.
◆ shiftX_
|
private |
Definition at line 120 of file bilinearForm.hh.
◆ shiftY_
|
private |
Definition at line 121 of file bilinearForm.hh.
The documentation for this class was generated from the following file:
- hp1D/bilinearForm.hh