A function class to calculate element matrices for the mass matrix. More...

#include <bf_identity.hh>

Inheritance diagram for hp2D::Identity< F >:
concepts::BilinearForm< Real, Real > hp2D::BilinearFormHelper_0_0< Real > concepts::Cloneable concepts::OutputOperator hp2D::LinearFormHelper_0< Real >

Public Member Functions

virtual Identity< F > * clone () const
 Virtual constructor. More...
 
virtual BilinearForm * clone () const=0
 Virtual constructor. More...
 
concepts::RCP< concepts::SharedJacobianDetdata () const
 Gets the pointer to the shared data. More...
 
void data (const concepts::RCP< concepts::SharedJacobianDet > d)
 Set the pointer to the shared data. More...
 
const concepts::ElementFormulaContainer< Real > formula () const
 
 Identity (const concepts::ElementFormulaContainer< F > frm=concepts::ElementFormulaContainer< F >(), bool all=false)
 Constructor. More...
 
void operator() (const concepts::Element< Real > &elmX, const concepts::Element< Real > &elmY, concepts::ElementMatrix< F > &em) const
 
virtual void operator() (const Element< Real > &elmX, const Element< Real > &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< Real > &elmX, const Element< Real > &elmY, ElementMatrix< Real > &em, const ElementPair< Real > &ep) const
 Evaluates the bilinear form for all shape functions on elmX and elmY and stores the result in the matrix em. More...
 
virtual ~Identity ()
 

Protected Member Functions

void computeIntermediate_ (const BaseQuad< concepts::Real > &elm) 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...
 
concepts::Array< Real > intermediateValue_
 Intermediate value (on each quadrature point) More...
 
concepts::RCP< concepts::SharedJacobianDetsharedData_
 Shared data for vectorial bilinear forms. More...
 

Private Member Functions

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

Private Attributes

bool all_
 For SumFactorization: More...
 
concepts::ElementMatrix< Real > mass1D_
 Local 1D mass matrix for Karniadakis basis. More...
 

Detailed Description

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

A function class to calculate element matrices for the mass matrix.

Test:
test::TestMatrices2D
Author
Kersten Schmidt, 2003
Examples
BGT_0.cc, exactDtN.cc, howToGetStarted.cc, hpFEM2d.cc, inhomDirichletBCs.cc, inhomDirichletBCsLagrange.cc, inhomNeumannBCs.cc, parallelizationTutorial.cc, and RobinBCs.cc.

Definition at line 57 of file bf_identity.hh.

Constructor & Destructor Documentation

◆ Identity()

template<class F = Real>
hp2D::Identity< F >::Identity ( const concepts::ElementFormulaContainer< F >  frm = concepts::ElementFormulaContainer< F >(),
bool  all = false 
)

Constructor.

◆ ~Identity()

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

Member Function Documentation

◆ clone() [1/2]

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

Virtual constructor.

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

Implements concepts::Cloneable.

◆ clone() [2/2]

virtual BilinearForm* concepts::BilinearForm< Real , Real >::clone ( ) const
pure virtualinherited

Virtual constructor.

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

◆ computeIntermediate_()

void hp2D::LinearFormHelper_0< Real >::computeIntermediate_ ( const BaseQuad< concepts::Real > &  elm) const
protectedinherited

Compute the intermediate data for element matrix computation.

This method is important for the derivated linear forms.

◆ data() [1/2]

Gets the pointer to the shared data.

◆ data() [2/2]

void hp2D::LinearFormHelper_0< Real >::data ( const concepts::RCP< concepts::SharedJacobianDet d)
inherited

Set the pointer to the shared data.

◆ formula()

const concepts::ElementFormulaContainer<Real > hp2D::BilinearFormHelper_0_0< Real >::formula
inlineinherited

Definition at line 57 of file bilinearFormHelper.hh.

◆ info()

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

Returns information in an output stream.

Reimplemented from concepts::BilinearForm< Real, Real >.

◆ operator()() [1/5]

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

◆ operator()() [2/5]

virtual void concepts::BilinearForm< Real , Real >::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

◆ operator()() [3/5]

virtual void concepts::BilinearForm< Real , Real >::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

Definition at line 57 of file bilinearForm.hh.

◆ operator()() [4/5]

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

Assembling for hp2D::InfiniteLaguerreQuad.

◆ operator()() [5/5]

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

Assembling for hp2D::Quad.

Member Data Documentation

◆ all_

template<class F = Real>
bool hp2D::Identity< F >::all_
private

For SumFactorization:

Definition at line 76 of file bf_identity.hh.

◆ frm_

concepts::ElementFormulaContainer<Real > hp2D::LinearFormHelper_0< Real >::frm_
protectedinherited

ElementFormula.

Definition at line 77 of file linearFormHelper.hh.

◆ intermediateValue_

concepts::Array<Real > hp2D::LinearFormHelper_0< Real >::intermediateValue_
mutableprotectedinherited

Intermediate value (on each quadrature point)

\[f(F_K(\xi))^\top \det J\]

Definition at line 74 of file linearFormHelper.hh.

◆ mass1D_

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

Local 1D mass matrix for Karniadakis basis.

Definition at line 74 of file bf_identity.hh.

◆ sharedData_

concepts::RCP<concepts::SharedJacobianDet> hp2D::LinearFormHelper_0< Real >::sharedData_
protectedinherited

Shared data for vectorial bilinear forms.

Definition at line 79 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