Class documentation of Concepts

Loading...
Searching...
No Matches
hp2D::BilinearFormTwoPartDeriv< F > Class Template Referenceabstract

#include <bf_partialderiv.hh>

Inheritance diagram for hp2D::BilinearFormTwoPartDeriv< F >:
concepts::BilinearForm< F, G > hp2D::BilinearFormHelper_1_1< F, G > concepts::Cloneable concepts::OutputOperator

Public Member Functions

 BilinearFormTwoPartDeriv (const enum partDerivType i, const enum partDerivType j, const concepts::ElementFormulaContainer< F > frm=concepts::ElementFormulaContainer< F >())
 
virtual BilinearFormTwoPartDeriv< F > * clone () const
 
virtual void operator() (const concepts::Element< Real > &elmX, const concepts::Element< Real > &elmY, concepts::ElementMatrix< F > &em) const
 
virtual void operator() (const Element< G > &elmX, const Element< G > &elmY, ElementMatrix< F > &em) const =0
 
virtual void operator() (const Element< G > &elmX, const Element< G > &elmY, ElementMatrix< F > &em, const ElementPair< G > &ep) const
 
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

virtual std::ostream & info (std::ostream &os) const
 Returns information in an output stream.
 
void computeIntermediate_ (const BaseQuad< Real > &elm, const int i=-1, const int j=-1) const
 

Protected Attributes

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 hp2D::BilinearFormTwoPartDeriv< F >

A function class to calculate element matrices for the bilinear form related to a partial derivative of the trial and (!) test functions (scalar).

The shape functions on the physical cell $K$ are defined by standard mapping local shape functions defined on the reference square $\hat{K}$.

\[
\int\limits_{K} f(\xi)\partial_j u \partial_i v \; d\xi
= 
\int\limits_{\hat{K}} f(F_K\hat{\xi})\nabla{\hat{u}} (J^{-1})_{\cdot,j}(J^{-\top})_{i,\cdot}\nabla{\hat{v}} \;|\det J| \; d\hat{\xi}  
\]

The class can be used to construct other bilinear forms, also for vectorial functions.

See also
vectorial::BilinearForm
Author
Kersten Schmidt, 2010

Definition at line 151 of file bf_partialderiv.hh.

Constructor & Destructor Documentation

◆ BilinearFormTwoPartDeriv()

template<class F = Real>
hp2D::BilinearFormTwoPartDeriv< F >::BilinearFormTwoPartDeriv ( const enum partDerivType  i,
const enum partDerivType  j,
const concepts::ElementFormulaContainer< F >  frm = concepts::ElementFormulaContainer< F >() 
)

Constructor

Parameters
iDirection of the partial derivative of the test function, X for $partial_x$ and Y for $partial_y$.
jDirection of the partial derivative of the trial function.

◆ ~BilinearFormTwoPartDeriv()

template<class F = Real>
virtual hp2D::BilinearFormTwoPartDeriv< F >::~BilinearFormTwoPartDeriv ( )
inlinevirtual

Definition at line 166 of file bf_partialderiv.hh.

Member Function Documentation

◆ clone()

template<class F = Real>
virtual BilinearFormTwoPartDeriv< F > * hp2D::BilinearFormTwoPartDeriv< F >::clone ( ) const
inlinevirtual

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

Implements concepts::BilinearForm< F, G >.

Definition at line 168 of file bf_partialderiv.hh.

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

◆ info()

template<class F = Real>
virtual std::ostream & hp2D::BilinearFormTwoPartDeriv< F >::info ( std::ostream &  os) const
protectedvirtual

Returns information in an output stream.

Reimplemented from concepts::BilinearForm< F, G >.

◆ operator()() [1/2]

template<class F , class G = typename Realtype<F>::type>
virtual void concepts::BilinearForm< F, G >::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

Implemented in vectorial::BilinearForm< F, G >, concepts::BilinearFormLiCo< F, G >, concepts::BilinearFormContainer< F, G >, concepts::BilinearF_Sum< F, H, J, G >, and concepts::BilinearF_W< F, H, J, G >.

◆ operator()() [2/2]

template<class F , class G = typename Realtype<F>::type>
virtual void concepts::BilinearForm< F, G >::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

Reimplemented in vectorial::BilinearForm< F, G >.

Definition at line 57 of file bilinearForm.hh.

Member Data Documentation

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