In file fespace.hpp:

class FESpace

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

[more]mutable ARRAY<SpecialElement*> specialelements
special elements for hacks (used for contact, periodic-boundary-penalty-constraints,

Public Methods

[more] FESpace(const MeshAccess & ama, int aorder, int adim, bool acomplex, bool parseflags=false)
old constructor type
[more] FESpace(const MeshAccess & ama, const Flags & flags, bool parseflags=false)
Constructor.
[more]virtual ~FESpace()
[more]virtual void Update()
update dof-tables
[more]virtual void Update(LocalHeap & lh)
update dof-table, preferred form
[more]virtual void PrintReport(ostream & ost)
print report to stream
[more]int GetOrder() const
order of finite elements
[more]int GetDimension() const
how many components
[more]bool IsComplex() const
complex space ?
[more]virtual const char* GetType()
will be replaced by getclassname !
[more]virtual int GetNDof() const = 0
number of global dofs
[more]virtual int GetNDofLevel(int level) const
number of global dofs on the level
[more]virtual const FiniteElement& GetFE(int elnr, LocalHeap & lh) const
returs finite element.
[more]virtual void GetDofNrs(int elnr, ARRAY<int> & dnums) const = 0
get dof-nrs of the element
[more]virtual void GetExternalDofNrs(int elnr, ARRAY<int> & dnums) const
get remaining dofs after static condensation
[more]virtual void GetWireBasketDofNrs(int vnr, ARRAY<int> & dnums) const
experiments with new preconditioners
[more]virtual void GetVertexDofNrs(int vnr, ARRAY<int> & dnums) const
experiments with new preconditioners
[more]virtual void GetEdgeDofNrs(int ednr, ARRAY<int> & dnums) const
experiments with new preconditioners
[more]virtual void GetFaceDofNrs(int fanr, ARRAY<int> & dnums) const
experiments with new preconditioners
[more]virtual void GetInnerDofNrs(int elnr, ARRAY<int> & dnums) const
experiments with new preconditioners
[more]virtual const FiniteElement& GetSFE(int selnr, LocalHeap & lh) const
returns surface element for boundary interals
[more]virtual void GetSDofNrs(int selnr, ARRAY<int> & dnums) const = 0
returns dofs of sourface element
[more]bool DefinedOn(int elnr) const
is the FESpace defined for this sub-domain nr ?
[more]bool DefinedOnBoundary(int belnr) const
[more]void SetDefinedOn(const BitArray & defon)
[more]void SetDefinedOnBoundary(const BitArray & defon)
[more]void SetDirichletBoundaries(const BitArray & dirbnds)
[more]const FiniteElement& GetFE(ELEMENT_TYPE type) const
Get reference element for tet, prism, trig, etc
[more]FESpace& LowOrderFESpace()
according low-order FESpace (if available)
[more]const FESpace& LowOrderFESpace() const
according low-order FESpace (if available)
[more]virtual void LockSomeDofs(BaseMatrix & mat) const
[more]virtual Table<int> * CreateSmoothingBlocks(int type = 0) const
[more]virtual Table<int> * CreateSmoothingBlocks(const Flags & flags) const
[more]virtual BitArray* CreateIntermediatePlanes(int type = 0) const
for anisotropic plane smoothing
[more]virtual const ngmg::Prolongation* GetProlongation() const
Returns multigrid-prolongation
[more]void SetProlongation(ngmg::Prolongation* aprol)
Set multigrid prolongation
[more]const BilinearFormIntegrator* GetEvaluator() const
returns function-evaluator
[more]const BilinearFormIntegrator* GetBoundaryEvaluator() const
returns function-evaluator for boundary values
[more]virtual MatrixGraph* GetGraph(int level, bool symmetric)
generates matrix graph

Protected Fields

[more]int order
order of finite elements
[more]int dimension
how many components
[more]bool iscomplex
complex space
[more]bool eliminate_internal
eliminate element-internal dofs ?
[more]ngmg::Prolongation* prol
prolongation operators between multigrid levels
[more]ARRAY<int> definedon
on which subdomains is the space defined ?
[more]ARRAY<int> definedonbound
on which boundaries is the space defined ?
[more]BitArray dirichlet_boundaries
prototype: what are the (homogeneous) Dirichlet boundaries ?
[more]BitArray dirichlet_dofs
dofs on Dirichlet boundary
[more]FiniteElement* tet
[more]FiniteElement* prism
[more]FiniteElement* pyramid
[more]FiniteElement* hex
[more]FiniteElement* trig
[more]FiniteElement* quad
[more]FiniteElement* segm
[more]BilinearFormIntegrator* evaluator
Evaluator for visualization
[more]BilinearFormIntegrator* boundary_evaluator
Evaluator for visualization of boundary data
[more]FESpace* low_order_space
if non-zero, pointer to low order space
[more]ARRAY<bool> directsolverclustered
if directsolverclustered[i] is true, then the unknowns of domain i are clustered

oint order
order of finite elements

oint dimension
how many components

obool iscomplex
complex space

obool eliminate_internal
eliminate element-internal dofs ?

ongmg::Prolongation* prol
prolongation operators between multigrid levels

oARRAY<int> definedon
on which subdomains is the space defined ?

oARRAY<int> definedonbound
on which boundaries is the space defined ?

oBitArray dirichlet_boundaries
prototype: what are the (homogeneous) Dirichlet boundaries ?

oBitArray dirichlet_dofs
dofs on Dirichlet boundary

oFiniteElement* tet

oFiniteElement* prism

oFiniteElement* pyramid

oFiniteElement* hex

oFiniteElement* trig

oFiniteElement* quad

oFiniteElement* segm

oBilinearFormIntegrator* evaluator
Evaluator for visualization

oBilinearFormIntegrator* boundary_evaluator
Evaluator for visualization of boundary data

oFESpace* low_order_space
if non-zero, pointer to low order space

oARRAY<bool> directsolverclustered
if directsolverclustered[i] is true, then the unknowns of domain i are clustered

o FESpace(const MeshAccess & ama, int aorder, int adim, bool acomplex, bool parseflags=false)
old constructor type

o 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

ovirtual ~FESpace()

ovirtual void Update()
update dof-tables

ovirtual void Update(LocalHeap & lh)
update dof-table, preferred form

ovirtual void PrintReport(ostream & ost)
print report to stream

oint GetOrder() const
order of finite elements

oint GetDimension() const
how many components

obool IsComplex() const
complex space ?

ovirtual const char* GetType()
will be replaced by getclassname !

ovirtual int GetNDof() const = 0
number of global dofs

ovirtual int GetNDofLevel(int level) const
number of global dofs on the level

ovirtual const FiniteElement& GetFE(int elnr, LocalHeap & lh) const
returs finite element. attention: should be thread-safe, but is not always

ovirtual void GetDofNrs(int elnr, ARRAY<int> & dnums) const = 0
get dof-nrs of the element

ovirtual void GetExternalDofNrs(int elnr, ARRAY<int> & dnums) const
get remaining dofs after static condensation

ovirtual void GetWireBasketDofNrs(int vnr, ARRAY<int> & dnums) const
experiments with new preconditioners

ovirtual void GetVertexDofNrs(int vnr, ARRAY<int> & dnums) const
experiments with new preconditioners

ovirtual void GetEdgeDofNrs(int ednr, ARRAY<int> & dnums) const
experiments with new preconditioners

ovirtual void GetFaceDofNrs(int fanr, ARRAY<int> & dnums) const
experiments with new preconditioners

ovirtual void GetInnerDofNrs(int elnr, ARRAY<int> & dnums) const
experiments with new preconditioners

ovirtual const FiniteElement& GetSFE(int selnr, LocalHeap & lh) const
returns surface element for boundary interals

ovirtual void GetSDofNrs(int selnr, ARRAY<int> & dnums) const = 0
returns dofs of sourface element

obool DefinedOn(int elnr) const
is the FESpace defined for this sub-domain nr ?

obool DefinedOnBoundary(int belnr) const

ovoid SetDefinedOn(const BitArray & defon)

ovoid SetDefinedOnBoundary(const BitArray & defon)

ovoid SetDirichletBoundaries(const BitArray & dirbnds)

oconst FiniteElement& GetFE(ELEMENT_TYPE type) const
Get reference element for tet, prism, trig, etc

oFESpace& LowOrderFESpace()
according low-order FESpace (if available)

oconst FESpace& LowOrderFESpace() const
according low-order FESpace (if available)

ovirtual void LockSomeDofs(BaseMatrix & mat) const

ovirtual Table<int> * CreateSmoothingBlocks(int type = 0) const

ovirtual Table<int> * CreateSmoothingBlocks(const Flags & flags) const

ovirtual BitArray* CreateIntermediatePlanes(int type = 0) const
for anisotropic plane smoothing

ovirtual const ngmg::Prolongation* GetProlongation() const
Returns multigrid-prolongation

ovoid SetProlongation(ngmg::Prolongation* aprol)
Set multigrid prolongation

oconst BilinearFormIntegrator* GetEvaluator() const
returns function-evaluator

oconst BilinearFormIntegrator* GetBoundaryEvaluator() const
returns function-evaluator for boundary values

ovirtual MatrixGraph* GetGraph(int level, bool symmetric)
generates matrix graph

omutable 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++.