linDG1D::BoundaryIntStab Class Referenceabstract

Stabilizing boundary integral for the DG FEM in 1D. More...

#include <bilinearForm.hh>

Inheritance diagram for linDG1D::BoundaryIntStab:
concepts::BilinearForm< Real > concepts::Cloneable concepts::OutputOperator

Public Member Functions

virtual BoundaryIntStabclone () const
 Virtual constructor. More...
 
virtual BilinearForm * clone () const=0
 Virtual constructor. More...
 
virtual void operator() (const concepts::Element< Real > &elmX, const concepts::Element< Real > &elmY, concepts::ElementMatrix< Real > &em) const
 
virtual void operator() (const concepts::Element< Real > &elmX, const concepts::Element< Real > &elmY, concepts::ElementMatrix< Real > &em, const concepts::ElementPair< Real > &ep) 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

Stabilizing boundary integral for the DG FEM in 1D.

This bilinear form integrates

\[ \int_{K\cap K'} [\phi_i] \cdot [\phi_j] \, dx. \]

Author
Thomas Wihler, 2003
Examples
linearDG1d.cc.

Definition at line 49 of file bilinearForm.hh.

Member Function Documentation

◆ clone() [1/2]

virtual BoundaryIntStab* linDG1D::BoundaryIntStab::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 60 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()

◆ operator()() [1/4]

virtual void linDG1D::BoundaryIntStab::operator() ( const concepts::Element< Real > &  elmX,
const concepts::Element< Real > &  elmY,
concepts::ElementMatrix< Real > &  em 
) const
inlinevirtual

Definition at line 51 of file bilinearForm.hh.

◆ operator()() [2/4]

virtual void linDG1D::BoundaryIntStab::operator() ( const concepts::Element< Real > &  elmX,
const concepts::Element< Real > &  elmY,
concepts::ElementMatrix< Real > &  em,
const concepts::ElementPair< Real > &  ep 
) const
virtual

◆ 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