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

#include <bf_partialderiv.hh>

Inheritance diagram for hp2D::BilinearFormOnePartDeriv< F >:
concepts::BilinearForm< Real > hp2D::BilinearFormHelper_0_1_Part< Real > concepts::Cloneable concepts::OutputOperator

Public Member Functions

 BilinearFormOnePartDeriv (const enum partDerivType i, const concepts::ElementFormulaContainer< F > frm=concepts::ElementFormulaContainer< F >())
 Constructor. More...
 
virtual BilinearFormOnePartDeriv< F > * clone () const
 Virtual constructor. More...
 
virtual BilinearForm * clone () const=0
 Virtual constructor. 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 void operator() (const concepts::Element< Real > &elmX, const concepts::Element< Real > &elmY, concepts::ElementMatrix< F > &em) 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...
 
virtual ~BilinearFormOnePartDeriv ()
 

Protected Member Functions

void computeIntermediate_ (const BaseQuad< concepts::Real > &elm, const int i) const
 Compute the intermediate data for element matrix computation. More...
 
virtual std::ostream & info (std::ostream &os) const
 Returns information in an output stream. More...
 

Protected Attributes

concepts::ElementFormulaContainer< Real > frm_
 ElementFormula. More...
 
ArrayElementFormula< concepts::Point< Real, 2 > > intermediateVector_
 Intermediate vector (on each quadrature point) More...
 
concepts::RCP< concepts::SharedJacobianAdj< 2 > > sharedData_
 Shared data for vectorial bilinear forms. More...
 

Private Member Functions

bool operator() (const hp2D::Quad< Real > *elmX, const hp2D::Quad< Real > *elmY, concepts::ElementMatrix< F > &em) const
 Assembling for hp2D::Quad. More...
 
bool operator() (const InfiniteLaguerreQuad *elmX, const InfiniteLaguerreQuad *elmY, concepts::ElementMatrix< F > &em) const
 

Private Attributes

enum partDerivType i_
 Component of the test function. More...
 
concepts::ElementMatrix< Real > mass1D_
 
concepts::ElementMatrix< Real > stiff1D_
 Local 1D stiffness and mass matrices for Karniadakis basis. More...
 

Detailed Description

template<class F = Real>
class hp2D::BilinearFormOnePartDeriv< F >

A function class to calculate element matrices for the bilinear form related to a partial derivative of the 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) u \partial_i v \; d\xi = \int\limits_{\hat{K}} f(F_K\hat{\xi})\hat{u}(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 85 of file bf_partialderiv.hh.

Constructor & Destructor Documentation

◆ BilinearFormOnePartDeriv()

Constructor.

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

◆ ~BilinearFormOnePartDeriv()

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

Definition at line 98 of file bf_partialderiv.hh.

Member Function Documentation

◆ clone() [1/2]

template<class F = Real>
virtual BilinearFormOnePartDeriv<F>* hp2D::BilinearFormOnePartDeriv< F >::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 100 of file bf_partialderiv.hh.

◆ clone() [2/2]

virtual BilinearForm* concepts::BilinearForm< Real , typename Realtype<Real >::type >::clone
pure virtualinherited

Virtual constructor.

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

◆ computeIntermediate_()

void hp2D::BilinearFormHelper_0_1_Part< Real >::computeIntermediate_ ( const BaseQuad< concepts::Real > &  elm,
const int  i 
) const
protectedinherited

Compute the intermediate data for element matrix computation.

Parameters
itake only i-th column of Jacobian matrix (for test function)

This method is important for the derivated bilinear forms.

◆ data() [1/2]

Gets the pointer to the shared data.

◆ data() [2/2]

void hp2D::BilinearFormHelper_0_1_Part< Real >::data ( const concepts::RCP< concepts::SharedJacobianAdj< 2 > >  d)
inherited

Set the pointer to the shared data.

◆ info()

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

Returns information in an output stream.

Reimplemented from concepts::BilinearForm< Real >.

◆ operator()() [1/5]

template<class F = Real>
virtual void hp2D::BilinearFormOnePartDeriv< F >::operator() ( const concepts::Element< Real > &  elmX,
const concepts::Element< Real > &  elmY,
concepts::ElementMatrix< F > &  em 
) const
virtual

◆ operator()() [2/5]

virtual void concepts::BilinearForm< Real , typename Realtype<Real >::type >::operator() ( const Element< typename Realtype<Real >::type > &  elmX,
const Element< typename Realtype<Real >::type > &  elmY,
ElementMatrix< Real > &  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()() [3/5]

virtual void concepts::BilinearForm< Real , typename Realtype<Real >::type >::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
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.

◆ operator()() [4/5]

template<class F = Real>
bool hp2D::BilinearFormOnePartDeriv< F >::operator() ( const hp2D::Quad< Real > *  elmX,
const hp2D::Quad< Real > *  elmY,
concepts::ElementMatrix< F > &  em 
) const
private

Assembling for hp2D::Quad.

◆ operator()() [5/5]

template<class F = Real>
bool hp2D::BilinearFormOnePartDeriv< F >::operator() ( const InfiniteLaguerreQuad elmX,
const InfiniteLaguerreQuad elmY,
concepts::ElementMatrix< F > &  em 
) const
private

Member Data Documentation

◆ frm_

ElementFormula.

Definition at line 103 of file bilinearFormHelper.hh.

◆ i_

template<class F = Real>
enum partDerivType hp2D::BilinearFormOnePartDeriv< F >::i_
private

Component of the test function.

Definition at line 109 of file bf_partialderiv.hh.

◆ intermediateVector_

ArrayElementFormula<concepts::Point<Real ,2> > hp2D::BilinearFormHelper_0_1_Part< Real >::intermediateVector_
mutableprotectedinherited

Intermediate vector (on each quadrature point)

\[\underline{f}(F_K(\xi))(\mbox{adj}(J)_{\cdot,i})^\top\]

Definition at line 101 of file bilinearFormHelper.hh.

◆ mass1D_

template<class F = Real>
concepts::ElementMatrix<Real> hp2D::BilinearFormOnePartDeriv< F >::mass1D_
private

Definition at line 115 of file bf_partialderiv.hh.

◆ sharedData_

concepts::RCP<concepts::SharedJacobianAdj<2> > hp2D::BilinearFormHelper_0_1_Part< Real >::sharedData_
protectedinherited

Shared data for vectorial bilinear forms.

Definition at line 106 of file bilinearFormHelper.hh.

◆ stiff1D_

template<class F = Real>
concepts::ElementMatrix<Real> hp2D::BilinearFormOnePartDeriv< F >::stiff1D_
private

Local 1D stiffness and mass matrices for Karniadakis basis.

Definition at line 115 of file bf_partialderiv.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