Class documentation of Concepts

Loading...
Searching...
No Matches
hp2Dedge::GraduvBase< F, G > Class Template Reference

#include <bf_graduv.hh>

Inheritance diagram for hp2Dedge::GraduvBase< F, G >:
hp2D::BilinearFormHelper_1_1< F, G > hp2Dedge::Graduv< F > hp2Dedge::GraduvMatrix< F >

Public Types

typedef concepts::Combtype< F, G >::type value_type
 

Public Member Functions

 GraduvBase (const concepts::ElementFormulaContainer< F > frm, bool all=false)
 
 GraduvBase (const concepts::ElementFormulaContainer< concepts::Mapping< G, 2 > > frm, bool all=false)
 
void data (const concepts::RCP< concepts::SharedJacobianAdj< 2 > > d)
 Set the pointer to the shared data.
 
concepts::RCP< concepts::SharedJacobianAdj< 2 > > data () const
 Gets the pointer to the shared data.
 

Protected Member Functions

bool assemble_ (const hp2D::Quad< Real > *elmX, const Quad< Real > *elmY, concepts::ElementMatrix< value_type > &em) const
 
void computeIntermediate_ (const BaseQuad< Real > &elm, const int i=-1, const int j=-1) const
 

Protected Attributes

bool all_
 Parameter for the sum factorisation.
 
concepts::Array< F > intermediateValue_
 
concepts::Array< concepts::Mapping< G, 2 > > intermediateMatrix_
 
concepts::ElementFormulaContainer< F > frm_
 Element formula.
 
concepts::ElementFormulaContainer< concepts::Mapping< G, 2 > > frmM_
 Matrix element formula.
 

Detailed Description

template<class F = Real, class G = typename concepts::Realtype<F>::type>
class hp2Dedge::GraduvBase< F, G >

Base class to calculate element matrices for the (grad u, v)-bilinearform, for both scalar and matrix formulas.

This class is inspired by hp2D::LaplaceBase.

Author
Vinh Tho Ma, 2016

Definition at line 63 of file bf_graduv.hh.

Member Typedef Documentation

◆ value_type

template<class F = Real, class G = typename concepts::Realtype<F>::type>
typedef concepts::Combtype<F,G>::type hp2Dedge::GraduvBase< F, G >::value_type

Definition at line 66 of file bf_graduv.hh.

Constructor & Destructor Documentation

◆ GraduvBase()

template<class F = Real, class G = typename concepts::Realtype<F>::type>
hp2Dedge::GraduvBase< F, G >::GraduvBase ( const concepts::ElementFormulaContainer< F >  frm,
bool  all = false 
)

Constructor. The formula frm is evaluated in each quadrature point.

Member Function Documentation

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

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.

Member Data Documentation

◆ all_

template<class F = Real, class G = typename concepts::Realtype<F>::type>
bool hp2Dedge::GraduvBase< F, G >::all_
protected

Parameter for the sum factorisation.

Definition at line 79 of file bf_graduv.hh.

◆ frm_

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

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_
protectedinherited

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_
mutableprotectedinherited

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_
mutableprotectedinherited

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.


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