Class documentation of Concepts

Loading...
Searching...
No Matches
hp2D::GradLinearForm< F > Class Template Referenceabstract

#include <linearForm.hh>

Inheritance diagram for hp2D::GradLinearForm< F >:
concepts::LinearForm< F, G > hp2D::LinearFormHelper_1< F > concepts::LinearFormChoice concepts::OutputOperator

Public Member Functions

 GradLinearForm (const concepts::ElementFormulaContainer< F > frm1, const concepts::ElementFormulaContainer< F > frm2, bool ignoreMissingElem=false)
 
 GradLinearForm (const concepts::ElementFormulaContainer< concepts::Point< F, 2 > > frm, bool ignoreMissingElem=false)
 
virtual void operator() (const concepts::Element< Real > &elm, concepts::ElementMatrix< F > &em) const
 
virtual void operator() (const Element< G > &elm, ElementMatrix< F > &em) const =0
 
virtual void setBasis (Basis basis)
 

Public Attributes

Basis basis_
 

Protected Member Functions

virtual std::ostream & info (std::ostream &os) const
 Returns information in an output stream.
 
void computeIntermediate_ (const BaseQuad< concepts::Real > &elm) const
 

Protected Attributes

ArrayElementFormula< concepts::Point< F, 2 > > intermediateVector_
 
concepts::ElementFormulaContainer< concepts::Point< F, 2 > > frm_
 ElementFormula.
 

Detailed Description

template<class F = concepts::Real>
class hp2D::GradLinearForm< F >

Linear form in 2D.

This linear form computes

\[ \int_K \vec{f} \cdot \mbox{\bf grad}{v} \, dx \]

Currently only on quadrilaterals.

Author
Kersten Schmidt, 2005
Examples
inhomDirichletBCs.cc.

Definition at line 136 of file linearForm.hh.

Constructor & Destructor Documentation

◆ GradLinearForm() [1/2]

template<class F = concepts::Real>
hp2D::GradLinearForm< F >::GradLinearForm ( const concepts::ElementFormulaContainer< F >  frm1,
const concepts::ElementFormulaContainer< F >  frm2,
bool  ignoreMissingElem = false 
)

Constructor.

Parameters
frm1First component of the formula
frm2Second component of the formula

◆ GradLinearForm() [2/2]

template<class F = concepts::Real>
hp2D::GradLinearForm< F >::GradLinearForm ( const concepts::ElementFormulaContainer< concepts::Point< F, 2 > >  frm,
bool  ignoreMissingElem = false 
)

Constructor.

Parameters
frmVectorial formula
ignoreMissingElemdo not throw an exceptiont for incompatible elements, just ignore elements like GfemQuads

Member Function Documentation

◆ computeIntermediate_()

template<class F >
void hp2D::LinearFormHelper_1< F >::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.

◆ info()

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

Returns information in an output stream.

Reimplemented from concepts::LinearForm< F, G >.

Reimplemented in hp2D::PlCurlLinearForm< F >.

◆ operator()() [1/2]

template<class F = concepts::Real>
virtual void hp2D::GradLinearForm< F >::operator() ( const concepts::Element< Real > &  elm,
concepts::ElementMatrix< F > &  em 
) const
virtual

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/2]

template<class F , class G = typename Realtype<F>::type>
virtual void concepts::LinearForm< F, G >::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

Implemented in vectorial::LinearForm< F, G >.

◆ 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.

◆ frm_

template<class F >
concepts::ElementFormulaContainer<concepts::Point<F, 2> > hp2D::LinearFormHelper_1< F >::frm_
protectedinherited

ElementFormula.

Definition at line 139 of file linearFormHelper.hh.

◆ intermediateVector_

template<class F >
ArrayElementFormula<concepts::Point<F,2> > hp2D::LinearFormHelper_1< F >::intermediateVector_
mutableprotectedinherited

Intermediate vector (on each quadrature point)

\[\underline{f}(F_K(\xi))^\top \mbox{adj}(J)^\top\]

Definition at line 136 of file linearFormHelper.hh.


The documentation for this class was generated from the following file: