hp3D::Laplace< F > Class Template Referenceabstract
A function class to calculate element matrices for the Laplacian. More...
#include <bilinearForm.hh>
Public Member Functions | |
virtual Laplace< F > * | clone () const |
Virtual constructor. More... | |
virtual BilinearForm * | clone () const=0 |
Virtual constructor. More... | |
Laplace () | |
Constructor. More... | |
virtual void | operator() (const concepts::Element< Real > &elmX, const concepts::Element< Real > &elmY, concepts::ElementMatrix< F > &em) const |
virtual void | operator() (const Element< typename Realtype< Real >::type > &elmX, const Element< typename Realtype< Real >::type > &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< typename Realtype< Real >::type > &elmX, const Element< typename Realtype< Real >::type > &elmY, ElementMatrix< Real > &em, const ElementPair< typename Realtype< Real >::type > &ep) const |
Evaluates the bilinear form for all shape functions on elmX and elmY and stores the result in the matrix em . More... | |
void | operator() (const Hexahedron &elmX, const Hexahedron &elmY, concepts::ElementMatrix< F > &em) const |
virtual | ~Laplace () |
Protected Member Functions | |
virtual std::ostream & | info (std::ostream &os) const |
Returns information in an output stream. More... | |
Private Attributes | |
concepts::Array< F > | jacobian_ |
Intermediate data for element matrix computation. More... | |
concepts::Array< concepts::MapReal3d > | jacobianInv_ |
const Hexahedron * | oldElm_ |
Detailed Description
template<class F = Real>
class hp3D::Laplace< F >
A function class to calculate element matrices for the Laplacian.
Definition at line 51 of file bilinearForm.hh.
Constructor & Destructor Documentation
◆ Laplace()
hp3D::Laplace< F >::Laplace | ( | ) |
Constructor.
◆ ~Laplace()
|
virtual |
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 62 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.
◆ info()
|
protectedvirtual |
Returns information in an output stream.
Reimplemented from concepts::BilinearForm< Real >.
◆ operator()() [1/4]
|
virtual |
◆ operator()() [2/4]
|
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/4]
|
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.
◆ operator()() [4/4]
void hp3D::Laplace< F >::operator() | ( | const Hexahedron & | elmX, |
const Hexahedron & | elmY, | ||
concepts::ElementMatrix< F > & | em | ||
) | const |
Member Data Documentation
◆ jacobian_
|
mutableprivate |
Intermediate data for element matrix computation.
Definition at line 67 of file bilinearForm.hh.
◆ jacobianInv_
|
mutableprivate |
Definition at line 68 of file bilinearForm.hh.
◆ oldElm_
|
mutableprivate |
Definition at line 70 of file bilinearForm.hh.
The documentation for this class was generated from the following file:
- hp3D/bilinearForm.hh