hp2D::LinearFormHelper_1< F > Class Template Reference

Helper class for linearforms l(v), where v is a one form. More...

#include <linearFormHelper.hh>

Inheritance diagram for hp2D::LinearFormHelper_1< F >:
hp2D::Advection< F >

Public Member Functions

 LinearFormHelper_1 (const concepts::ElementFormulaContainer< concepts::Point< F, 2 > > frm)
 
 LinearFormHelper_1 (const concepts::ElementFormulaContainer< F > frm1, const concepts::ElementFormulaContainer< F > frm2)
 
virtual ~LinearFormHelper_1 ()
 

Protected Member Functions

void computeIntermediate_ (const BaseQuad< concepts::Real > &elm) const
 Compute the intermediate data for element matrix computation. More...
 

Protected Attributes

concepts::ElementFormulaContainer< concepts::Point< F, 2 > > frm_
 ElementFormula. More...
 
ArrayElementFormula< concepts::Point< F, 2 > > intermediateVector_
 Intermediate vector (on each quadrature point) More...
 

Detailed Description

template<class F>
class hp2D::LinearFormHelper_1< F >

Helper class for linearforms l(v), where v is a one form.

\[\displaystyle l(v) = \int\limits_K \underline{f}(x)^\top\underline{v}\,dx = \int\limits_{\hat{K}} \underline{f}(F_K(\xi))^\top J^{-1} \hat{\underline{v}}\,\circ F_K^{-1} \det J\,d\xi = \int\limits_{\hat{K}} \underline{f}(F_K(\xi))^\top \mbox{adj}(J)^\top\hat{\underline{v}}\,\circ F_K^{-1}\,d\xi\]

Here, $F_K$ is the element mapping from reference element $\hat{K} = [0,1]^2$, $J$ is the Jacobian matrix.

Precomputes intermediate data for element matrix computation, this is

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

which is stored as a single vector intermediateVector_ for each quadrature point.

The class can be used as well for bilinear forms a(u,v) where u is a 0-form and v is a 1-form. One example is the bilinear form hp2D::Advection

\[ \int\limits_{K} \underline{k}^\top u \; \nabla{v} \; d\xi = \int\limits_{\hat{K}} \underline{k}^\top \hat{u} \; J^{-\top}\;\nabla{\hat{v}} \;|\det J| \; d\hat{\xi}. \]

Author
Kersten Schmidt, 2005

Definition at line 117 of file linearFormHelper.hh.

Constructor & Destructor Documentation

◆ LinearFormHelper_1() [1/2]

template<class F >
hp2D::LinearFormHelper_1< F >::LinearFormHelper_1 ( const concepts::ElementFormulaContainer< F >  frm1,
const concepts::ElementFormulaContainer< F >  frm2 
)

◆ LinearFormHelper_1() [2/2]

◆ ~LinearFormHelper_1()

template<class F >
virtual hp2D::LinearFormHelper_1< F >::~LinearFormHelper_1 ( )
virtual

Member Function Documentation

◆ computeIntermediate_()

template<class F >
void hp2D::LinearFormHelper_1< F >::computeIntermediate_ ( const BaseQuad< concepts::Real > &  elm) const
protected

Compute the intermediate data for element matrix computation.

This method is important for the derivated linear forms.

Member Data Documentation

◆ frm_

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

ElementFormula.

Definition at line 139 of file linearFormHelper.hh.

◆ intermediateVector_

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

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