Helper class for bilinearforms a(u,v), where u and v are 1-forms, which computes intermediate data for element matrix computation. More...

#include <bilinearFormHelper.hh>

Public Member Functions

 BilinearFormHelper_1_1 (const concepts::ElementFormulaContainer< concepts::Mapping< G, 2 > > frmM)
 Constructor for matrix valued function. More...
 
 BilinearFormHelper_1_1 (const concepts::ElementFormulaContainer< F > frm)
 Constructor for scalar function. More...
 
concepts::RCP< concepts::SharedJacobianAdj< 2 > > data () const
 Gets the pointer to the shared data. More...
 
void data (const concepts::RCP< concepts::SharedJacobianAdj< 2 > > d)
 Set the pointer to the shared data. More...
 
virtual ~BilinearFormHelper_1_1 ()
 Destructor. More...
 

Protected Member Functions

void computeIntermediate_ (const BaseQuad< Real > &elm, const int i=-1, const int j=-1) const
 Compute the intermediate data for element matrix computation. More...
 

Protected Attributes

concepts::ElementFormulaContainer< F > frm_
 Element formula. More...
 
concepts::ElementFormulaContainer< concepts::Mapping< G, 2 > > frmM_
 Matrix element formula. More...
 
concepts::Array< concepts::Mapping< G, 2 > > intermediateMatrix_
 Intermediate matrix. More...
 
concepts::Array< F > intermediateValue_
 Intermediate value. More...
 

Private Member Functions

void computeadjJ_adjJT_rank1_ (concepts::Array< concepts::Mapping< Real, 2 > > &intermediateMatrix, const int i, const int j) const
 Computes the either adj(J)*adj(J)^T or in the case of partial derivatives (i > NONE, j > NONE) the rank-1 product of j-th column of adj(J) with i-th row of adj(J)^T. More...
 
void computeJacobianMatrix_ (const BaseQuad< Real > &elm, concepts::Array< concepts::Mapping< Real, 2 > > &J, concepts::Array< Real > &detJ_inv) const
 Compute the Jacobian matrix and the inverse of its determinant on each quadrature point. More...
 

Private Attributes

concepts::RCP< concepts::SharedJacobianAdj< 2 > > sharedData_
 Shared data for vectorial bilinear forms. More...
 

Detailed Description

template<class F, class G = typename concepts::Realtype<F>::type>
class hp2D::BilinearFormHelper_1_1< F, G >

Helper class for bilinearforms a(u,v), where u and v are 1-forms, which computes intermediate data for element matrix computation.

It works with scalar functions and matrix valued functions.

For scalar function $f(x)$ we have

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

For matrix valued function $A(x)$ we have

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

For scalar functions also the intermediate data for a bilinear form with partial derivatives can be stored, where the intermediate matrix is the rank-1 product of j-th column of $adj(J)$ with i-th row of $adj(J)^T$.

Here, J is the Jacobian matrix.

The adjunct of the Jacobian matrix and the inverse of the determinant of the Jacobian matrix can be shared between different instances of the class. This is useful for bilinear forms on vectorial spaces.

Author
Kersten Schmidt, 2007

Definition at line 140 of file bilinearFormHelper.hh.

Constructor & Destructor Documentation

◆ BilinearFormHelper_1_1() [1/2]

template<class F , class G = typename concepts::Realtype<F>::type>
hp2D::BilinearFormHelper_1_1< F, G >::BilinearFormHelper_1_1 ( const concepts::ElementFormulaContainer< F >  frm)

Constructor for scalar function.

◆ BilinearFormHelper_1_1() [2/2]

template<class F , class G = typename concepts::Realtype<F>::type>
hp2D::BilinearFormHelper_1_1< F, G >::BilinearFormHelper_1_1 ( const concepts::ElementFormulaContainer< concepts::Mapping< G, 2 > >  frmM)

Constructor for matrix valued function.

◆ ~BilinearFormHelper_1_1()

template<class F , class G = typename concepts::Realtype<F>::type>
virtual hp2D::BilinearFormHelper_1_1< F, G >::~BilinearFormHelper_1_1 ( )
virtual

Destructor.

Member Function Documentation

◆ computeadjJ_adjJT_rank1_()

template<class F , class G = typename concepts::Realtype<F>::type>
void hp2D::BilinearFormHelper_1_1< F, G >::computeadjJ_adjJT_rank1_ ( concepts::Array< concepts::Mapping< Real, 2 > > &  intermediateMatrix,
const int  i,
const int  j 
) const
private

Computes the either adj(J)*adj(J)^T or in the case of partial derivatives (i > NONE, j > NONE) the rank-1 product of j-th column of adj(J) with i-th row of adj(J)^T.

◆ computeIntermediate_()

template<class F , class G = typename concepts::Realtype<F>::type>
void hp2D::BilinearFormHelper_1_1< F, G >::computeIntermediate_ ( const BaseQuad< Real > &  elm,
const int  i = -1,
const int  j = -1 
) const
protected

Compute the intermediate data for element matrix computation.

Parameters
iif i=0 or 1, then take only i-th column of Jacobian matrix (for test function)
jif j=0 or 1, then take only j-th column of Jacobian matrix (for trial function)

The Jacobian matrices have to been taken both full (i,j = -1) or both partial (i,j = 0 or 1).

Matrix formulas and complex valued scalar formulas are only implemented for full Jacobians.

◆ computeJacobianMatrix_()

template<class F , class G = typename concepts::Realtype<F>::type>
void hp2D::BilinearFormHelper_1_1< F, G >::computeJacobianMatrix_ ( const BaseQuad< Real > &  elm,
concepts::Array< concepts::Mapping< Real, 2 > > &  J,
concepts::Array< Real > &  detJ_inv 
) const
private

Compute the Jacobian matrix and the inverse of its determinant on each quadrature point.

◆ data() [1/2]

template<class F , class G = typename concepts::Realtype<F>::type>
concepts::RCP<concepts::SharedJacobianAdj<2> > hp2D::BilinearFormHelper_1_1< F, G >::data ( ) const

Gets the pointer to the shared data.

◆ data() [2/2]

template<class F , class G = typename concepts::Realtype<F>::type>
void hp2D::BilinearFormHelper_1_1< F, G >::data ( const concepts::RCP< concepts::SharedJacobianAdj< 2 > >  d)

Set the pointer to the shared data.

Member Data Documentation

◆ frm_

template<class F , class G = typename concepts::Realtype<F>::type>
concepts::ElementFormulaContainer<F> hp2D::BilinearFormHelper_1_1< F, G >::frm_
protected

Element formula.

Definition at line 193 of file bilinearFormHelper.hh.

◆ frmM_

template<class F , class G = typename concepts::Realtype<F>::type>
concepts::ElementFormulaContainer<concepts::Mapping<G,2> > hp2D::BilinearFormHelper_1_1< F, G >::frmM_
protected

Matrix element formula.

Definition at line 195 of file bilinearFormHelper.hh.

◆ intermediateMatrix_

template<class F , class G = typename concepts::Realtype<F>::type>
concepts::Array<concepts::Mapping<G,2> > hp2D::BilinearFormHelper_1_1< F, G >::intermediateMatrix_
mutableprotected

Intermediate matrix.

In case of a scalar formula:

\[\mbox{adj}(J) \mbox{adj}(J)^\top\]

In case of a matrix formula $M$:

\[\mbox{adj}(J) M \mbox{adj}(J)^\top\]

In case of partial Jacobian:

\[\mbox{adj}(J)_{\cdot,j} (\mbox{adj}(J)_{\cdot,i})^\top\]

Definition at line 191 of file bilinearFormHelper.hh.

◆ intermediateValue_

template<class F , class G = typename concepts::Realtype<F>::type>
concepts::Array<F> hp2D::BilinearFormHelper_1_1< F, G >::intermediateValue_
mutableprotected

Intermediate value.

In case of a scalar formula:

\[\frac{f(F_K(\xi))}{\det J}\]

In case of a matrix formula:

\[\frac{1}{\det J}\]

Definition at line 179 of file bilinearFormHelper.hh.

◆ sharedData_

template<class F , class G = typename concepts::Realtype<F>::type>
concepts::RCP<concepts::SharedJacobianAdj<2> > hp2D::BilinearFormHelper_1_1< F, G >::sharedData_
private

Shared data for vectorial bilinear forms.

Definition at line 212 of file bilinearFormHelper.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