hp2Dedge::EdgeIdentity Class Referenceabstract

Bilinear form in 2D, integral over an edge. More...

#include <bilinearForm.hh>

Inheritance diagram for hp2Dedge::EdgeIdentity:
concepts::BilinearForm< Real > concepts::Cloneable concepts::OutputOperator

Public Member Functions

virtual EdgeIdentityclone () const
 Virtual constructor. More...
 
virtual BilinearForm * clone () const=0
 Virtual constructor. More...
 
void operator() (const concepts::Element< Real > &elmX, const concepts::Element< Real > &elmY, concepts::ElementMatrix< Real > &em) const
 Computes the element mass matrix. More...
 
void operator() (const Edge< Real > &elmX, const Edge< Real > &elmY, concepts::ElementMatrix< Real > &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...
 

Protected Member Functions

virtual std::ostream & info (std::ostream &os) const
 Returns information in an output stream. More...
 

Detailed Description

Bilinear form in 2D, integral over an edge.

This linear form computes

\[\int_{\partial K \cap \Gamma_N} g \vec{t} \vec v \, ds. \]

where $\vec{t} = \vec{n}^\bot = {-n_y \choose n_x} $ is the tangential vector, $\vec{n}$ the outer normal vector and $g$ the tangential component of the Neumann boundary.

Author
Kersten Schmidt, 2004

Definition at line 118 of file bilinearForm.hh.

Member Function Documentation

◆ clone() [1/2]

virtual EdgeIdentity* hp2Dedge::EdgeIdentity::clone ( ) const
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 133 of file bilinearForm.hh.

◆ clone() [2/2]

virtual BilinearForm* concepts::BilinearForm< Real , typename Realtype<Real >::type >::clone ( ) const
pure virtualinherited

Virtual constructor.

Returns a pointer to a copy of itself. The caller is responsible to destroy this copy.

◆ info()

virtual std::ostream& hp2Dedge::EdgeIdentity::info ( std::ostream &  os) const
protectedvirtual

Returns information in an output stream.

Reimplemented from concepts::BilinearForm< Real >.

◆ operator()() [1/4]

void hp2Dedge::EdgeIdentity::operator() ( const concepts::Element< Real > &  elmX,
const concepts::Element< Real > &  elmY,
concepts::ElementMatrix< Real > &  em 
) const

Computes the element mass matrix.

There are the loops over all quadrature points and the loops over all shape functions.

Parameters
elmXThe edge with the test functions
elmYThe edge with the trial functions
emThe element matrix

◆ operator()() [2/4]

void hp2Dedge::EdgeIdentity::operator() ( const Edge< Real > &  elmX,
const Edge< Real > &  elmY,
concepts::ElementMatrix< Real > &  em 
) const

◆ operator()() [3/4]

virtual void concepts::BilinearForm< Real , typename Realtype<Real >::type >::operator() ( const Element< G > &  elmX,
const Element< G > &  elmY,
ElementMatrix< F > &  em 
) const
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
elmXLeft element (test functions)
elmYRight element (trial functions)
emReturn element matrix

◆ operator()() [4/4]

virtual void concepts::BilinearForm< Real , typename Realtype<Real >::type >::operator() ( const Element< G > &  elmX,
const Element< G > &  elmY,
ElementMatrix< F > &  em,
const ElementPair< G > &  ep 
) const
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
elmXLeft element
elmYRight element
emReturn element matrix
epElement pair holding more information on the pair elmX and elmY

Definition at line 57 of file bilinearForm.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