A function class to calculate element matrices of the bilinear form. More...

#include <bilinearForm.hh>

Inheritance diagram for hp3D::Advection< F >:
concepts::BilinearForm< Real, Real > hp3D::LinearFormHelper_1< Real > concepts::Cloneable concepts::OutputOperator

Public Member Functions

 Advection (const concepts::ElementFormulaContainer< concepts::Point< F, 3 > > frm, bool all=false)
 
 Advection (const concepts::ElementFormulaContainer< F > frm1, const concepts::ElementFormulaContainer< F > frm2, const concepts::ElementFormulaContainer< F > frm3, bool all=false)
 
virtual Advection< F > * clone () const
 Virtual constructor. More...
 
virtual BilinearForm * clone () const=0
 Virtual constructor. More...
 
virtual 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 ~Advection ()
 

Protected Member Functions

void computeIntermediate_ (const Hexahedron &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< concepts::Point< Real, 3 > > frm_
 ElementFormula. More...
 
ArrayElementFormula< concepts::Point< Real, 3 > > intermediateVector_
 Intermediate vector (on each quadrature point) More...
 

Private Member Functions

bool operator() (const Hexahedron *elmX, const Hexahedron *elmY, concepts::ElementMatrix< F > &em) const
 Assembling for hp3D::Hexahedron. More...
 

Private Attributes

bool all_
 

Detailed Description

template<class F = Real>
class hp3D::Advection< F >

A function class to calculate element matrices of the bilinear form.

\[ \int\limits_{K} \underline{k}^\top u \; \nabla{v} \; d\xi = \int\limits_{\hat{K}} \underline{k}^\top \hat{u} \; J^{-\top}\;\nabla{\hat{v}} \;|\det J| \; d\hat{\xi}. \]

Here k is an arbitrary vector-valued function with coefficients. For some k, the resulting matrix might be singular, e.g. k=[0,0].

Author
Kersten Schmidt, 2014

Definition at line 186 of file bilinearForm.hh.

Constructor & Destructor Documentation

◆ Advection() [1/2]

template<class F = Real>
hp3D::Advection< F >::Advection ( const concepts::ElementFormulaContainer< F >  frm1,
const concepts::ElementFormulaContainer< F >  frm2,
const concepts::ElementFormulaContainer< F >  frm3,
bool  all = false 
)
inline

Definition at line 190 of file bilinearForm.hh.

◆ Advection() [2/2]

template<class F = Real>
hp3D::Advection< F >::Advection ( const concepts::ElementFormulaContainer< concepts::Point< F, 3 > >  frm,
bool  all = false 
)
inline

Definition at line 197 of file bilinearForm.hh.

◆ ~Advection()

template<class F = Real>
virtual hp3D::Advection< F >::~Advection ( )
inlinevirtual

Definition at line 203 of file bilinearForm.hh.

Member Function Documentation

◆ clone() [1/2]

template<class F = Real>
virtual Advection<F>* hp3D::Advection< 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 205 of file bilinearForm.hh.

◆ 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 hp3D::LinearFormHelper_1< Real >::computeIntermediate_ ( const Hexahedron elm) const
protectedinherited

Compute the intermediate data for element matrix computation.

This method is important for the derivated linear forms.

◆ info()

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

Returns information in an output stream.

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

◆ operator()() [1/4]

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

◆ operator()() [2/4]

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/4]

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/4]

template<class F = Real>
bool hp3D::Advection< F >::operator() ( const Hexahedron elmX,
const Hexahedron elmY,
concepts::ElementMatrix< F > &  em 
) const
private

Assembling for hp3D::Hexahedron.

Member Data Documentation

◆ all_

template<class F = Real>
bool hp3D::Advection< F >::all_
private

Definition at line 216 of file bilinearForm.hh.

◆ frm_

concepts::ElementFormulaContainer<concepts::Point<Real , 3> > hp3D::LinearFormHelper_1< Real >::frm_
protectedinherited

ElementFormula.

Definition at line 88 of file linearFormHelper.hh.

◆ intermediateVector_

ArrayElementFormula<concepts::Point<Real ,3> > hp3D::LinearFormHelper_1< Real >::intermediateVector_
mutableprotectedinherited

Intermediate vector (on each quadrature point)

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

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