hp2D::Riesz< F > Class Template Referenceabstract
Linear form in 2D. More...
#include <linearForm.hh>
Public Member Functions | |
concepts::RCP< concepts::SharedJacobianDet > | data () const |
Gets the pointer to the shared data. More... | |
void | data (const concepts::RCP< concepts::SharedJacobianDet > d) |
Set the pointer to the shared data. More... | |
void | operator() (const concepts::Element< Real > &elm, concepts::ElementMatrix< F > &em) const |
Computes the element load vector. More... | |
virtual void | operator() (const Element< typename Realtype< Real >::type > &elm, ElementMatrix< Real > &em) const=0 |
Computes the element contribution to the function. More... | |
Riesz (const concepts::ElementFormulaContainer< F > frm, concepts::BoundaryConditions *bc=0, bool ignoreMissingElem=false) | |
Constructor. More... | |
Riesz (const std::string &str, concepts::BoundaryConditions *bc=0, bool ignoreMissingElem=false) | |
Contructor via a string of real valued Riesz. More... | |
virtual void | setBasis (Basis basis) |
virtual | ~Riesz () |
Public Attributes | |
Basis | basis_ |
Protected Member Functions | |
void | computeIntermediate_ (const BaseQuad< concepts::Real > &elm) const |
Compute the intermediate data for element matrix computation. More... | |
virtual std::ostream & | info (std::ostream &os) const |
Protected Attributes | |
concepts::ElementFormulaContainer< Real > | frm_ |
ElementFormula. More... | |
concepts::Array< Real > | intermediateValue_ |
Intermediate value (on each quadrature point) More... | |
concepts::RCP< concepts::SharedJacobianDet > | sharedData_ |
Shared data for vectorial bilinear forms. More... | |
Private Member Functions | |
void | operator() (const Quad< Real > &elm, concepts::ElementMatrix< F > &em) const |
Private Attributes | |
concepts::BoundaryConditions * | bc_ |
Reference to the boundary conditions. More... | |
bool | ignoreMissingElem |
concepts::Array< Real > | jacobian_ |
Intermediate data for element matrix computation. More... | |
Detailed Description
template<class F = Real>
class hp2D::Riesz< F >
Linear form in 2D.
This linear form computes
Currently only on quadrilaterals. The boundary integral is only correct for straigt boundaries.
- Examples
- elasticity2D_tutorial.cc, howToGetStarted.cc, hpFEM2d-simple.cc, hpFEM2d.cc, inhomDirichletBCs.cc, inhomDirichletBCsLagrange.cc, and parallelizationTutorial.cc.
Definition at line 60 of file linearForm.hh.
Constructor & Destructor Documentation
◆ Riesz() [1/2]
hp2D::Riesz< F >::Riesz | ( | const concepts::ElementFormulaContainer< F > | frm, |
concepts::BoundaryConditions * | bc = 0 , |
||
bool | ignoreMissingElem = false |
||
) |
Constructor.
- Parameters
-
frm The formula, given elementwise bc Boundary conditions, defaults to homogeneous
◆ Riesz() [2/2]
hp2D::Riesz< F >::Riesz | ( | const std::string & | str, |
concepts::BoundaryConditions * | bc = 0 , |
||
bool | ignoreMissingElem = false |
||
) |
Contructor via a string of real valued Riesz.
This will be interpreted as a concepts::ParsedFormula
- Parameters
-
frm The formula, given as string
◆ ~Riesz()
|
virtual |
Member Function Documentation
◆ computeIntermediate_()
|
protectedinherited |
Compute the intermediate data for element matrix computation.
This method is important for the derivated linear forms.
◆ data() [1/2]
|
inherited |
Gets the pointer to the shared data.
◆ data() [2/2]
|
inherited |
Set the pointer to the shared data.
◆ info()
|
protectedvirtual |
Reimplemented from concepts::LinearForm< Real >.
◆ operator()() [1/3]
void hp2D::Riesz< F >::operator() | ( | const concepts::Element< Real > & | elm, |
concepts::ElementMatrix< F > & | em | ||
) | const |
Computes the element load vector.
As for the computation of an element stiffness matrix, there are the loops over all quadrature points and the loops over all shape functions.
- Parameters
-
elm The element for which the load vector should be computed. em The load vector
◆ operator()() [2/3]
|
pure virtualinherited |
Computes the element contribution to the function.
- Parameters
-
elm Element on which the computations should be performed em The local matrix
◆ operator()() [3/3]
|
private |
◆ setBasis()
|
inlinevirtualinherited |
Definition at line 68 of file linearForm.hh.
Member Data Documentation
◆ basis_
|
mutableinherited |
Definition at line 71 of file linearForm.hh.
◆ bc_
|
private |
Reference to the boundary conditions.
Definition at line 106 of file linearForm.hh.
◆ frm_
|
protectedinherited |
ElementFormula.
Definition at line 77 of file linearFormHelper.hh.
◆ ignoreMissingElem
|
private |
Definition at line 107 of file linearForm.hh.
◆ intermediateValue_
|
mutableprotectedinherited |
◆ jacobian_
|
mutableprivate |
Intermediate data for element matrix computation.
Definition at line 100 of file linearForm.hh.
◆ sharedData_
|
protectedinherited |
Shared data for vectorial bilinear forms.
Definition at line 79 of file linearFormHelper.hh.
The documentation for this class was generated from the following file:
- hp2D/linearForm.hh