Enumerations | |
| enum | partDerivType { NO_DERIV = 0 , X_DERIV = 1 , Y_DERIV = 2 , Z_DERIV = 3 } | 
Functions | |
| template<class F > | |
| void | setupIdentity (vectorial::BilinearForm< F, typename concepts::Realtype< F >::type > &bf) | 
| template<class F > | |
| void | setupIdentity (vectorial::BilinearForm< F, typename concepts::Realtype< F >::type > &bf, const concepts::ElementFormulaContainer< concepts::Mapping< F, 3 > > frm, bool transposed=false) | 
| template<class F > | |
| void | setupAdvection (vectorial::BilinearForm< F, typename concepts::Realtype< F >::type > &bf, const concepts::ElementFormulaContainer< concepts::Point< F, 3 > > frm) | 
| template<class F > | |
| ShapeFunction3D< F > | makeShapeFunction3D (const Hexahedron &elm) | 
| std::ostream & | operator<< (std::ostream &os, const TrivialWeight &p) | 
| std::ostream & | operator<< (std::ostream &os, const ShortestDist &p) | 
| std::ostream & | operator<< (std::ostream &os, const PostprocessSquare &p) | 
| std::ostream & | operator<< (std::ostream &os, const PostprocessRoot &p) | 
| std::ostream & | operator<< (std::ostream &os, const PostprocessRoot4 &p) | 
| std::ostream & | operator<< (std::ostream &os, const Postprocess34 &p) | 
| std::ostream & | operator<< (std::ostream &os, const Postprocess14 &p) | 
| std::ostream & | operator<< (std::ostream &os, const Postprocess4 &p) | 
| std::ostream & | operator<< (std::ostream &os, const Postprocess3 &p) | 
| template<typename DistClass , typename Function > | |
| std::ostream & | operator<< (std::ostream &os, const DistancePost< DistClass, Function > &p) | 
| std::ostream & | operator<< (std::ostream &os, const ProductOfAll &p) | 
| std::ostream & | operator<< (std::ostream &os, const DaugeWeight &p) | 
| std::ostream & | operator<< (std::ostream &os, const ShortestDistLimited &p) | 
3D hp-FEM for H1-conforming elements.
The Space can be built using full tensor product polynomial spaces in the elements or trunk spaces (and much more) by changing the way the degrees of freedom are built in the BuildDofsBase class (and its specializations).
| enum hp3D::partDerivType | 
Definition at line 31 of file bf3D_partialDeriv.hh.
      
  | 
  inline | 
Definition at line 175 of file shortestDist.hh.
| void hp3D::setupAdvection | ( | vectorial::BilinearForm< F, typename concepts::Realtype< F >::type > & | bf, | 
| const concepts::ElementFormulaContainer< concepts::Point< F, 3 > > | frm | ||
| ) | 
Function to setup a bilinear form related to the vector Advection, namely
![\[\sum_{i,j=1}^3 \int_K w_j u_i\frac{\partial v_i}{\partial_{x_j}}\mathrm{d}x
  = \sum_{i=1}^3 \int_K \underline{w}\cdot\nabla v_i u_i \mathrm{d}x
  = \int_K \underline{w}^\top\nabla\underline{v}\underline{u}\mathrm{d}x\]](form_573.png)
on vectorial spaces, where $\underline{u}$ and $\underline{v}$ are trial and test functions respectively, $\underline{w}$ is a given vectorial formula, and
![\[
\nabla \underline{v} = \left(\begin{array}{ccc}
\frac{\partial u_1}{\partial x_1} & 
\frac{\partial u_1}{\partial x_2} &
\frac{\partial u_1}{\partial x_3} \\
\frac{\partial u_2}{\partial x_1} & 
\frac{\partial u_2}{\partial x_2} &
\frac{\partial u_2}{\partial x_3} \\
\frac{\partial u_3}{\partial x_1} & 
\frac{\partial u_3}{\partial x_2} &
\frac{\partial u_3}{\partial x_3}.
\end{array}\right)
\]](form_574.png)
| void hp3D::setupIdentity | ( | vectorial::BilinearForm< F, typename concepts::Realtype< F >::type > & | bf | ) | 
| void hp3D::setupIdentity | ( | vectorial::BilinearForm< F, typename concepts::Realtype< F >::type > & | bf, | 
| const concepts::ElementFormulaContainer< concepts::Mapping< F, 3 > > | frm, | ||
| bool | transposed = false  | 
        ||
| ) | 
Function to setup a bilinear form related to the vector Identity, namely
![\[\int_K \underline{u}^\top A\underline{v}\mathrm{d}x\]](form_572.png)
wih a matrix function A on vectorial spaces.
If transposed is true then the transposed of A is taken instead.