Linear form in 2D. More...

#include <linearForm.hh>

Inheritance diagram for hp2D::Riesz< F >:
concepts::LinearForm< Real > hp2D::LinearFormHelper_0< Real > concepts::LinearFormChoice

Public Member Functions

concepts::RCP< concepts::SharedJacobianDetdata () 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::SharedJacobianDetsharedData_
 Shared data for vectorial bilinear forms. More...
 

Private Member Functions

void operator() (const Quad< Real > &elm, concepts::ElementMatrix< F > &em) const
 

Private Attributes

concepts::BoundaryConditionsbc_
 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

\[ \int_K f v \, dx + \int_{\partial K \cap \Gamma_N} g v \, \mathrm{d}s. \]

Currently only on quadrilaterals. The boundary integral is only correct for straigt boundaries.

Author
Kersten Schmidt, 2003
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]

template<class F = Real>
hp2D::Riesz< F >::Riesz ( const concepts::ElementFormulaContainer< F >  frm,
concepts::BoundaryConditions bc = 0,
bool  ignoreMissingElem = false 
)

Constructor.

Parameters
frmThe formula, given elementwise
bcBoundary conditions, defaults to homogeneous

◆ Riesz() [2/2]

template<class F = Real>
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
frmThe formula, given as string

◆ ~Riesz()

template<class F = Real>
virtual hp2D::Riesz< F >::~Riesz ( )
virtual

Member Function Documentation

◆ computeIntermediate_()

void hp2D::LinearFormHelper_0< Real >::computeIntermediate_ ( const BaseQuad< concepts::Real > &  elm) const
protectedinherited

Compute the intermediate data for element matrix computation.

This method is important for the derivated linear forms.

◆ data() [1/2]

Gets the pointer to the shared data.

◆ data() [2/2]

void hp2D::LinearFormHelper_0< Real >::data ( const concepts::RCP< concepts::SharedJacobianDet d)
inherited

Set the pointer to the shared data.

◆ info()

template<class F = Real>
virtual std::ostream& hp2D::Riesz< F >::info ( std::ostream &  os) const
protectedvirtual

Reimplemented from concepts::LinearForm< Real >.

◆ operator()() [1/3]

template<class F = Real>
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
elmThe element for which the load vector should be computed.
emThe load vector

◆ operator()() [2/3]

virtual void concepts::LinearForm< Real , typename Realtype<Real >::type >::operator() ( const Element< G > &  elm,
ElementMatrix< F > &  em 
) const
pure virtualinherited

Computes the element contribution to the function.

Parameters
elmElement on which the computations should be performed
emThe local matrix

◆ operator()() [3/3]

template<class F = Real>
void hp2D::Riesz< F >::operator() ( const Quad< Real > &  elm,
concepts::ElementMatrix< F > &  em 
) const
private

◆ setBasis()

virtual void concepts::LinearFormChoice::setBasis ( Basis  basis)
inlinevirtualinherited

Definition at line 68 of file linearForm.hh.

Member Data Documentation

◆ basis_

Basis concepts::LinearFormChoice::basis_
mutableinherited

Definition at line 71 of file linearForm.hh.

◆ bc_

template<class F = Real>
concepts::BoundaryConditions* hp2D::Riesz< F >::bc_
private

Reference to the boundary conditions.

Definition at line 106 of file linearForm.hh.

◆ frm_

concepts::ElementFormulaContainer<Real > hp2D::LinearFormHelper_0< Real >::frm_
protectedinherited

ElementFormula.

Definition at line 77 of file linearFormHelper.hh.

◆ ignoreMissingElem

template<class F = Real>
bool hp2D::Riesz< F >::ignoreMissingElem
private

Definition at line 107 of file linearForm.hh.

◆ intermediateValue_

concepts::Array<Real > hp2D::LinearFormHelper_0< Real >::intermediateValue_
mutableprotectedinherited

Intermediate value (on each quadrature point)

\[f(F_K(\xi))^\top \det J\]

Definition at line 74 of file linearFormHelper.hh.

◆ jacobian_

template<class F = Real>
concepts::Array<Real> hp2D::Riesz< F >::jacobian_
mutableprivate

Intermediate data for element matrix computation.

Definition at line 100 of file linearForm.hh.

◆ sharedData_

concepts::RCP<concepts::SharedJacobianDet> hp2D::LinearFormHelper_0< Real >::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:
Page URL: http://wiki.math.ethz.ch/bin/view/Concepts/WebHome
21 August 2020
© 2020 Eidgenössische Technische Hochschule Zürich