In file fespace.hpp: Base class for finite element space.
Documentation
Base class for finite element space.
Provides finite elements, global degrees of freedom,
and transformations of element-matrices and element-vectors
Inheritance:
Public Fields
-
mutable ARRAY<SpecialElement*> specialelements
- special elements for hacks (used for contact, periodic-boundary-penalty-constraints,
Public Methods
-
FESpace(const MeshAccess & ama, int aorder, int adim, bool acomplex, bool parseflags=false)
- old constructor type
-
FESpace(const MeshAccess & ama, const Flags & flags, bool parseflags=false)
- Constructor.
-
virtual ~FESpace()
-
virtual void Update()
- update dof-tables
-
virtual void Update(LocalHeap & lh)
- update dof-table, preferred form
-
virtual void PrintReport(ostream & ost)
- print report to stream
-
int GetOrder() const
- order of finite elements
-
int GetDimension() const
- how many components
-
bool IsComplex() const
- complex space ?
-
virtual const char* GetType()
- will be replaced by getclassname !
-
virtual int GetNDof() const = 0
- number of global dofs
-
virtual int GetNDofLevel(int level) const
- number of global dofs on the level
-
virtual const FiniteElement& GetFE(int elnr, LocalHeap & lh) const
- returs finite element.
-
virtual void GetDofNrs(int elnr, ARRAY<int> & dnums) const = 0
- get dof-nrs of the element
-
virtual void GetExternalDofNrs(int elnr, ARRAY<int> & dnums) const
- get remaining dofs after static condensation
-
virtual void GetWireBasketDofNrs(int vnr, ARRAY<int> & dnums) const
- experiments with new preconditioners
-
virtual void GetVertexDofNrs(int vnr, ARRAY<int> & dnums) const
- experiments with new preconditioners
-
virtual void GetEdgeDofNrs(int ednr, ARRAY<int> & dnums) const
- experiments with new preconditioners
-
virtual void GetFaceDofNrs(int fanr, ARRAY<int> & dnums) const
- experiments with new preconditioners
-
virtual void GetInnerDofNrs(int elnr, ARRAY<int> & dnums) const
- experiments with new preconditioners
-
virtual const FiniteElement& GetSFE(int selnr, LocalHeap & lh) const
- returns surface element for boundary interals
-
virtual void GetSDofNrs(int selnr, ARRAY<int> & dnums) const = 0
- returns dofs of sourface element
-
bool DefinedOn(int elnr) const
- is the FESpace defined for this sub-domain nr ?
-
bool DefinedOnBoundary(int belnr) const
-
void SetDefinedOn(const BitArray & defon)
-
void SetDefinedOnBoundary(const BitArray & defon)
-
void SetDirichletBoundaries(const BitArray & dirbnds)
-
const FiniteElement& GetFE(ELEMENT_TYPE type) const
- Get reference element for tet, prism, trig, etc
-
FESpace& LowOrderFESpace()
- according low-order FESpace (if available)
-
const FESpace& LowOrderFESpace() const
- according low-order FESpace (if available)
-
virtual void LockSomeDofs(BaseMatrix & mat) const
-
virtual Table<int> * CreateSmoothingBlocks(int type = 0) const
-
virtual Table<int> * CreateSmoothingBlocks(const Flags & flags) const
-
virtual BitArray* CreateIntermediatePlanes(int type = 0) const
- for anisotropic plane smoothing
-
virtual const ngmg::Prolongation* GetProlongation() const
- Returns multigrid-prolongation
-
void SetProlongation(ngmg::Prolongation* aprol)
- Set multigrid prolongation
-
const BilinearFormIntegrator* GetEvaluator() const
- returns function-evaluator
-
const BilinearFormIntegrator* GetBoundaryEvaluator() const
- returns function-evaluator for boundary values
-
virtual MatrixGraph* GetGraph(int level, bool symmetric)
- generates matrix graph
Protected Fields
-
int order
- order of finite elements
-
int dimension
- how many components
-
bool iscomplex
- complex space
-
bool eliminate_internal
- eliminate element-internal dofs ?
-
ngmg::Prolongation* prol
- prolongation operators between multigrid levels
-
ARRAY<int> definedon
- on which subdomains is the space defined ?
-
ARRAY<int> definedonbound
- on which boundaries is the space defined ?
-
BitArray dirichlet_boundaries
- prototype: what are the (homogeneous) Dirichlet boundaries ?
-
BitArray dirichlet_dofs
- dofs on Dirichlet boundary
-
FiniteElement* tet
-
FiniteElement* prism
-
FiniteElement* pyramid
-
FiniteElement* hex
-
FiniteElement* trig
-
FiniteElement* quad
-
FiniteElement* segm
-
BilinearFormIntegrator* evaluator
- Evaluator for visualization
-
BilinearFormIntegrator* boundary_evaluator
- Evaluator for visualization of boundary data
-
FESpace* low_order_space
- if non-zero, pointer to low order space
-
ARRAY<bool> directsolverclustered
- if directsolverclustered[i] is true, then the unknowns of domain i are clustered
int order
- order of finite elements
int dimension
- how many components
bool iscomplex
- complex space
bool eliminate_internal
- eliminate element-internal dofs ?
ngmg::Prolongation* prol
- prolongation operators between multigrid levels
ARRAY<int> definedon
- on which subdomains is the space defined ?
ARRAY<int> definedonbound
- on which boundaries is the space defined ?
BitArray dirichlet_boundaries
- prototype: what are the (homogeneous) Dirichlet boundaries ?
BitArray dirichlet_dofs
- dofs on Dirichlet boundary
FiniteElement* tet
FiniteElement* prism
FiniteElement* pyramid
FiniteElement* hex
FiniteElement* trig
FiniteElement* quad
FiniteElement* segm
BilinearFormIntegrator* evaluator
- Evaluator for visualization
BilinearFormIntegrator* boundary_evaluator
- Evaluator for visualization of boundary data
FESpace* low_order_space
- if non-zero, pointer to low order space
ARRAY<bool> directsolverclustered
- if directsolverclustered[i] is true, then the unknowns of domain i are clustered
FESpace(const MeshAccess & ama, int aorder, int adim, bool acomplex, bool parseflags=false)
- old constructor type
FESpace(const MeshAccess & ama, const Flags & flags, bool parseflags=false)
-
Constructor.
Used flags are:
-order=<int>: finite element order
-dim=<int>: number of components
-complex: complex space
-eliminate_internal: eliminate internal dofs
-dirichlet=<int-list>: dirichlet boundaries, 1-based
virtual ~FESpace()
virtual void Update()
- update dof-tables
virtual void Update(LocalHeap & lh)
- update dof-table, preferred form
virtual void PrintReport(ostream & ost)
- print report to stream
int GetOrder() const
- order of finite elements
int GetDimension() const
- how many components
bool IsComplex() const
- complex space ?
virtual const char* GetType()
- will be replaced by getclassname !
virtual int GetNDof() const = 0
- number of global dofs
virtual int GetNDofLevel(int level) const
- number of global dofs on the level
virtual const FiniteElement& GetFE(int elnr, LocalHeap & lh) const
- returs finite element. attention: should be thread-safe, but is not always
virtual void GetDofNrs(int elnr, ARRAY<int> & dnums) const = 0
- get dof-nrs of the element
virtual void GetExternalDofNrs(int elnr, ARRAY<int> & dnums) const
- get remaining dofs after static condensation
virtual void GetWireBasketDofNrs(int vnr, ARRAY<int> & dnums) const
- experiments with new preconditioners
virtual void GetVertexDofNrs(int vnr, ARRAY<int> & dnums) const
- experiments with new preconditioners
virtual void GetEdgeDofNrs(int ednr, ARRAY<int> & dnums) const
- experiments with new preconditioners
virtual void GetFaceDofNrs(int fanr, ARRAY<int> & dnums) const
- experiments with new preconditioners
virtual void GetInnerDofNrs(int elnr, ARRAY<int> & dnums) const
- experiments with new preconditioners
virtual const FiniteElement& GetSFE(int selnr, LocalHeap & lh) const
- returns surface element for boundary interals
virtual void GetSDofNrs(int selnr, ARRAY<int> & dnums) const = 0
- returns dofs of sourface element
bool DefinedOn(int elnr) const
- is the FESpace defined for this sub-domain nr ?
bool DefinedOnBoundary(int belnr) const
void SetDefinedOn(const BitArray & defon)
void SetDefinedOnBoundary(const BitArray & defon)
void SetDirichletBoundaries(const BitArray & dirbnds)
const FiniteElement& GetFE(ELEMENT_TYPE type) const
- Get reference element for tet, prism, trig, etc
FESpace& LowOrderFESpace()
- according low-order FESpace (if available)
const FESpace& LowOrderFESpace() const
- according low-order FESpace (if available)
virtual void LockSomeDofs(BaseMatrix & mat) const
virtual Table<int> * CreateSmoothingBlocks(int type = 0) const
virtual Table<int> * CreateSmoothingBlocks(const Flags & flags) const
virtual BitArray* CreateIntermediatePlanes(int type = 0) const
- for anisotropic plane smoothing
virtual const ngmg::Prolongation* GetProlongation() const
- Returns multigrid-prolongation
void SetProlongation(ngmg::Prolongation* aprol)
- Set multigrid prolongation
const BilinearFormIntegrator* GetEvaluator() const
- returns function-evaluator
const BilinearFormIntegrator* GetBoundaryEvaluator() const
- returns function-evaluator for boundary values
virtual MatrixGraph* GetGraph(int level, bool symmetric)
- generates matrix graph
mutable ARRAY<SpecialElement*> specialelements
- special elements for hacks (used for contact, periodic-boundary-penalty-constraints,
- Direct child classes:
- SurfaceElementFESpace
RaviartThomasFESpace
NonconformingFESpace
NonConformingFESpace
NodalFESpaceAlt
NodalFESpace
NedelecFESpace2
NedelecFESpace
HDivSymFESpace
HCurlHighOrderFESpaceAlt
HCurlHighOrderFESpace
H1HighOrderFESpace
ElementFESpace
CompoundFESpace
Alphabetic index HTML hierarchy of classes or Java
This page was generated with the help of DOC++.