|
NGSolve
4.9
|
High Order Finite Element Space. More...
#include <h1hofespace.hpp>
Public Member Functions | |
| H1HighOrderFESpace (const MeshAccess &ama, const Flags &flags, bool checkflags=false) | |
| virtual string | GetClassName () const |
| virtual void | Update (LocalHeap &lh) |
| update dof-tables, old style | |
| virtual void | PrintReport (ostream &ost) |
| print report to stream | |
| virtual int | GetNDof () const |
| number of dofs of process | |
| virtual int | GetNDofLevel (int alevel) const |
| number of dofs on the level | |
| virtual const FiniteElement & | GetFE (int elnr, LocalHeap &lh) const |
| returns finite element. | |
| virtual const FiniteElement & | GetSFE (int elnr, LocalHeap &lh) const |
| returns surface element for boundary interals | |
| virtual void | GetDofNrs (int elnr, Array< int > &dnums) const |
| get dof-nrs of the element | |
| virtual void | GetVertexDofNrs (int vnr, Array< int > &dnums) const |
| get number of low-order dofs for node of type nt | |
| virtual void | GetEdgeDofNrs (int ednr, Array< int > &dnums) const |
| get dofs on edge enr | |
| virtual void | GetFaceDofNrs (int fanr, Array< int > &dnums) const |
| get dofs on face fnr | |
| virtual void | GetInnerDofNrs (int elnr, Array< int > &dnums) const |
| get dofs on element (=cell) elnr | |
| virtual void | GetSDofNrs (int selnr, Array< int > &dnums) const |
| returns dofs of sourface element | |
| virtual Table< int > * | CreateSmoothingBlocks (const Flags &precflags) const |
| virtual Array< int > * | CreateDirectSolverClusters (const Flags &flags) const |
| for anisotropic plane smoothing: | |
| int | GetFirstFaceDof (int i) const |
| int | GetFirstEdgeDof (int i) const |
| int | GetFirstElementDof (int i) const |
| void | UpdateDofTables () |
| virtual void | UpdateCouplingDofArray () |
| void | SetEdgeOrder (int enr, int eo) |
| void | SetFaceOrder (int fnr, int fo) |
| void | SetFaceOrder (int fnr, int ox, int oy) |
| void | SetElementOrder (int elnr, int elo) |
| void | SetElementOrder (int elnr, int ox, int oy, int oz) |
| virtual int | GetRelOrder () const |
| get relative (to mesh) order of finite elements | |
| virtual bool | VarOrder () const |
| IntRange | GetEdgeDofs (int nr) const |
| IntRange | GetFaceDofs (int nr) const |
| IntRange | GetElementDofs (int nr) const |
Protected Attributes | |
| int | level |
| Array< int > | first_edge_dof |
| Array< int > | first_face_dof |
| Array< int > | first_element_dof |
| int | rel_order |
| relative order to mesh-order | |
| bool | var_order |
| bool | fixed_order |
| bool | wb_loedge |
| Array< int > | order_edge |
| Array< INT< 2 > > | order_face |
| Array< INT< 3 > > | order_inner |
| Array< bool > | used_vertex |
| Array< bool > | used_edge |
| Array< bool > | used_face |
| int | ndof |
| int | uniform_order_inner |
| int | uniform_order_face |
| int | uniform_order_edge |
| int | uniform_order_quad |
| int | uniform_order_trig |
| Array< INT< 3 > > | dom_order_min |
| Array< INT< 3 > > | dom_order_max |
| int | smoother |
| Array< int > | ndlevel |
| bool | level_adapted_order |
| bool | nodalp2 |
High Order Finite Element Space.
| virtual void ngcomp::H1HighOrderFESpace::GetVertexDofNrs | ( | int | vnr, |
| Array< int > & | dnums | ||
| ) | const [virtual] |
get number of low-order dofs for node of type nt
get dofs on vertex vnr
Reimplemented from ngcomp::FESpace.
| virtual void ngcomp::H1HighOrderFESpace::Update | ( | LocalHeap & | lh | ) | [virtual] |
1.7.6.1