#include <bilinearForm.hh>
Public Member Functions | |
| Advection (const concepts::ElementFormulaContainer< F > frm1, const concepts::ElementFormulaContainer< F > frm2, const concepts::ElementFormulaContainer< F > frm3, bool all=false) | |
| Advection (const concepts::ElementFormulaContainer< concepts::Point< F, 3 > > frm, bool all=false) | |
| virtual Advection< 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 |
Protected Member Functions | |
| virtual std::ostream & | info (std::ostream &os) const |
| Returns information in an output stream. | |
| void | computeIntermediate_ (const Hexahedron &elm) const |
Protected Attributes | |
| ArrayElementFormula< concepts::Point< F, 3 > > | intermediateVector_ |
| concepts::ElementFormulaContainer< concepts::Point< F, 3 > > | frm_ |
| ElementFormula. | |
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}.
\]](form_489.png)
Here k is an arbitrary vector-valued function with coefficients. For some k, the resulting matrix might be singular, e.g. k=[0,0].
Definition at line 186 of file bilinearForm.hh.
|
inline |
Definition at line 190 of file bilinearForm.hh.
|
inline |
Definition at line 197 of file bilinearForm.hh.
|
inlinevirtual |
Definition at line 203 of file bilinearForm.hh.
|
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 205 of file bilinearForm.hh.
|
protectedinherited |
Compute the intermediate data for element matrix computation
This method is important for the derivated linear forms.
|
protectedvirtual |
Returns information in an output stream.
Reimplemented from concepts::BilinearForm< F, G >.
|
pure virtualinherited |
Evaluates the bilinear form for all shape functions on elmX and elmY and stores the result in the matrix em.
em has the correct size. | elmX | Left element (test functions) |
| elmY | Right element (trial functions) |
| em | Return 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 >.
|
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.
em has the correct size. | elmX | Left element |
| elmY | Right element |
| em | Return element matrix |
| ep | Element pair holding more information on the pair elmX and elmY |
Reimplemented in vectorial::BilinearForm< F, G >.
Definition at line 57 of file bilinearForm.hh.
|
protectedinherited |
ElementFormula.
Definition at line 88 of file linearFormHelper.hh.
|
mutableprotectedinherited |