Hermes2D
3.0
|
Represents a finite element space over a domain. More...
#include <space.h>
Classes | |
struct | BaseComponent |
struct | EdgeInfo |
class | ElementData |
union | NodeData |
class | ReferenceSpaceCreator |
Class for creating reference space. More... | |
Public Member Functions | |
Space (MeshSharedPtr mesh, Shapeset *shapeset, EssentialBCs< Scalar > *essential_bcs) | |
void | init () |
Common code for constructors. | |
virtual void | init (Shapeset *shapeset, int p_init, bool assign_dofs_init=true)=0 |
Common code for constructors of subclasses. | |
virtual | ~Space () |
Destructor. | |
void | free () |
Freeing the data. | |
void | free_bc_data () |
Free BC data. | |
int | get_element_order (int id) const |
Returns element polynomial order. | |
int | get_num_dofs () const |
Returns the number of basis functions contained in the space. | |
MeshSharedPtr | get_mesh () const |
virtual SpaceType | get_type () const =0 |
Internal. Return type of this space. See enum SpaceType. | |
int | get_vertex_functions_count () |
int | get_edge_functions_count () |
Returns the total (global) number of edge functions. | |
int | get_bubble_functions_count () |
Returns the total (global) number of bubble functions. | |
int | get_seq () const |
Internal. Used by DiscreteProblem to detect changes in the space. | |
EssentialBCs< Scalar > * | get_essential_bcs () const |
Obtains an boundary conditions. | |
virtual void | get_element_assembly_list (Element *e, AsmList< Scalar > *al) const |
Obtains an assembly list for the given element. | |
virtual int | get_edge_order (Element *e, int edge) const |
Internal. Obtains the order of an edge, according to the minimum rule. | |
int | get_max_dof () const |
Returns the DOF number of the last basis function. | |
void | get_boundary_assembly_list (Element *e, int surf_num, AsmList< Scalar > *al) const |
Obtains an edge assembly list (contains shape functions that are nonzero on the specified edge). | |
Shapeset * | get_shapeset () const |
virtual void | set_shapeset (Shapeset *shapeset, bool clone=false)=0 |
Sets the shapeset. | |
void | set_mesh (MeshSharedPtr mesh) |
Sets a (new) mesh and calls assign_dofs(). | |
void | set_mesh_seq (int seq) |
Sets a (new) mesh seq, and mesh_seq. | |
void | set_essential_bcs (EssentialBCs< Scalar > *essential_bcs) |
Sets the boundary condition. | |
virtual void | set_element_order (int id, int order, int order_v=-1) |
void | set_uniform_order (int order, std::string marker=HERMES_ANY) |
void | adjust_element_order (int order_change, int min_order) |
void | adjust_element_order (int horizontal_order_change, int vertical_order_change, unsigned int horizontal_min_order, unsigned int vertical_min_order) |
Version for quads. | |
void | set_uniform_order_internal (int order, int marker) |
virtual int | assign_dofs (int first_dof=0) |
Builds basis functions and assigns DOF numbers to them. More... | |
void | unrefine_all_mesh_elements (bool keep_initial_refinements=true) |
virtual Scalar * | get_bc_projection (SurfPos *surf_pos, int order, EssentialBoundaryCondition< Scalar > *bc)=0 |
void | update_essential_bc_values () |
void | save (const char *filename) const |
Saves this space into a file. | |
void | load (const char *filename) |
This method is here for rapid re-loading. | |
virtual void | copy (SpaceSharedPtr< Scalar > space, MeshSharedPtr new_mesh) |
![]() | |
XMLParsing () | |
Constructor. | |
void | set_validation (bool to_set) |
Set to validate / not to validate. | |
Static Public Member Functions | |
static int | get_num_dofs (std::vector< SpaceSharedPtr< Scalar > > spaces) |
Returns the number of basis functions contained in the spaces. | |
static int | get_num_dofs (SpaceSharedPtr< Scalar > space) |
Returns the number of basis functions contained in the space. | |
static int | assign_dofs (std::vector< SpaceSharedPtr< Scalar > > spaces) |
Assings the degrees of freedom to all Spaces in the std::vector. | |
static void | unrefine_all_mesh_elements (std::vector< SpaceSharedPtr< Scalar > > spaces, bool keep_initial_refinements=true) |
static void | update_essential_bc_values (std::vector< SpaceSharedPtr< Scalar > > spaces, double time) |
static void | update_essential_bc_values (SpaceSharedPtr< Scalar > space, double time) |
static SpaceSharedPtr< Scalar > | load (const char *filename, MeshSharedPtr mesh, bool validate=false, EssentialBCs< Scalar > *essential_bcs=nullptr, Shapeset *shapeset=nullptr) |
Loads a space from a file in XML format. | |
Protected Member Functions | |
virtual void | set_element_order_internal (int id, int order, int order_v=-1) |
void | distribute_orders (MeshSharedPtr mesh, int *parents) |
Sets polynomial orders to elements created by Mesh::regularize() using "parents". | |
bool | is_up_to_date () const |
Returns true if the space is ready for computation, false otherwise. | |
virtual bool | isOkay () const |
State querying helpers. | |
std::string | getClassName () const |
virtual int | get_edge_order_internal (Node *en) const |
Internal. | |
void | unrefine_all_mesh_elements_internal (bool keep_initial_refinements, bool only_unrefine_space_data) |
virtual void | resize_tables () |
Updates internal node and element tables. More... | |
void | update_orders_recurrent (Element *e, int order) |
virtual void | reset_dof_assignment () |
Resets assignment of DOF to an unassigned state. | |
virtual void | assign_vertex_dofs ()=0 |
virtual void | assign_edge_dofs ()=0 |
virtual void | assign_bubble_dofs ()=0 |
virtual void | get_vertex_assembly_list (Element *e, int iv, AsmList< Scalar > *al) const =0 |
virtual void | get_boundary_assembly_list_internal (Element *e, int surf_num, AsmList< Scalar > *al) const =0 |
virtual void | get_bubble_assembly_list (Element *e, AsmList< Scalar > *al) const |
void | precalculate_projection_matrix (int nv, double **&mat, double *&p) |
void | update_edge_bc (Element *e, SurfPos *surf_pos) |
virtual void | update_constraints () |
virtual void | post_assign () |
Static Protected Member Functions | |
static Node * | get_mid_edge_vertex_node (Element *e, unsigned char i, unsigned char j) |
static SpaceSharedPtr< Scalar > | init_empty_space (SpaceType spaceType, MeshSharedPtr mesh, Shapeset *shapeset) |
Protected Attributes | |
int | ndof |
Number of degrees of freedom (dimension of the space). | |
Shapeset * | shapeset |
bool | own_shapeset |
true if default shapeset is created in the constructor, false if shapeset is supplied by user. | |
EssentialBCs< Scalar > * | essential_bcs |
Boundary conditions. | |
MeshSharedPtr | mesh |
FE mesh. | |
int | vertex_functions_count |
For statistics. | |
int | edge_functions_count |
int | bubble_functions_count |
int | first_dof |
For equation systems. | |
int | next_dof |
unsigned int | seq |
Tracking changes. | |
unsigned int | seq_assigned |
Tracking changes - mark call to assign_dofs(). | |
int | mesh_seq |
Tracking changes - mesh. | |
NodeData * | ndata |
node data table | |
int | nsize |
node data table size | |
ElementData * | edata |
element data table | |
int | esize |
element data table size | |
double ** | proj_mat |
double * | chol_p |
std::vector< Scalar * > | bc_data_projections |
Used for bc projection. | |
std::vector< typename Space < Scalar >::BaseComponent * > | bc_data_base_components |
![]() | |
bool | validate |
Internal. | |
Static Protected Attributes | |
static const int | H2D_UNASSIGNED_DOF = -2 |
DOF which was not assigned yet. | |
static const int | H2D_CONSTRAINED_DOF = -1 |
DOF which is constrained. | |
Represents a finite element space over a domain.
The Space class represents a finite element space over a domain defined by 'mesh', spanned by basis functions constructed using 'shapeset'. It serves as a base class for H1Space, HcurlSpace, HdivSpace and L2Space, since most of the functionality is common for all these spaces.
There are four main functions the Space class provides:
It handles the Dirichlet (essential) boundary conditions. The user provides a pointer to an instance of the EssentialBCs class that determines which markers represent the Dirichlet part of the boundary, and which markers represent the Neumann and Newton parts (all other than those with an Essential BC specified on them). Handling of these conditions is done naturally - through weak formulation.
It stores element polynomial degrees, or 'orders'. All active elements need to have an order set for the Space to be valid. Individual orders can be set by calling set_element_order(). You can also set the same order for all elements using set_uniform_order(). Quadrilateral elements can have different orders in the vertical and horizontal directions. It is therefore necessary to form the order using the macro H2D_MAKE_QUAD_ORDER() when calling the aforementioned functions.
It builds and enumerates the basis functions. After all element orders have been set, you must call the function assign_dofs(). This function assigns the DOF (degree-of- freedom) numbers to basis functions, starting with 'first_dof' (optional parameter). It also determines constraining relationships in the mesh due to hanging nodes and builds constrained basis functions. The total number of basis functions can then be obtained by calling (). Standard basis functions are assigned positive numbers from 'first_dof' to ('first_dof' + (get_num_dofs() - 1) * 'stride'). All shape functions belonging to the Dirichlet lift are assigned DOF number of -1. This way the Dirichlet lift becomes a (virtual) basis function. This simplifies assembling.
Finally, and most importantly, the Space is able to deliver a list of shape functions existing on each element. Such a list is called an "assembly list" and is represented by the class AsmList<Scalar>. The assembly list contains the triplets (idx, dof, coef). 'idx' is the shape function index, as understood by the Shapeset. 'dof' is the number of the basis function, whose part the shape function forms. 'coef' is a Real constant, with which the shape function must be multiplied in order to fit into the basis function. This is typically 1, but can be less than that for constrained functions. Constrained vertex functions can belong to more than one basis functions. This results in more triplets with the same 'idx'. However, the assembling procedure or the Solution class do not have to worry about that. For instance, the Solution class simply obtains the values of all shape functions contained in the list, multiplies them by 'coef', and forms a linear combination of them by taking the solution vector values at the 'dof' positions. This way the solution to the PDE is obtained.
Space is an abstract class and cannot be instatiated. Use one of the specializations H1Space, HcurlSpace or L2Space instead.
The handling of irregular meshes is desribed in H1Space and HcurlSpace.
void Hermes::Hermes2D::Space< Scalar >::adjust_element_order | ( | int | order_change, |
int | min_order | ||
) |
|
virtual |
Builds basis functions and assigns DOF numbers to them.
This functions must be called after assigning element orders, and before using the space in a computation, otherwise an error will occur.
first_dof[in] | The DOF number of the first basis function. |
stride[in] | The difference between the DOF numbers of successive basis functions. |
Definition at line 864 of file space.cpp.
Referenced by Hermes::Hermes2D::H1Space< Scalar >::init(), Hermes::Hermes2D::PicardSolver< Scalar >::init_solving(), Hermes::Hermes2D::NewtonSolver< Scalar >::init_solving(), and Hermes::Hermes2D::Space< Scalar >::load().
|
virtual |
Copy from Space instance 'space'
[in] | new_mesh | Mesh where data will be copied to. |
Reimplemented in Hermes::Hermes2D::H1Space< Scalar >, Hermes::Hermes2D::HdivSpace< Scalar >, Hermes::Hermes2D::HcurlSpace< Scalar >, and Hermes::Hermes2D::L2Space< Scalar >.
Definition at line 243 of file space.cpp.
Referenced by Hermes::Hermes2D::L2Space< Scalar >::copy(), Hermes::Hermes2D::H1Space< Scalar >::copy(), and Hermes::Hermes2D::HdivSpace< Scalar >::copy().
int Hermes::Hermes2D::Space< Scalar >::get_vertex_functions_count | ( | ) |
|
staticprotected |
|
protectedvirtual |
Auxiliary function the descendants may implement to perform additional tasks after the DOFs have been assigned.
Definition at line 835 of file space.cpp.
Referenced by Hermes::Hermes2D::Space< Scalar >::assign_dofs().
|
protectedvirtual |
Updates internal node and element tables.
Since meshes only contain geometric information, the Space class keeps two tables with FEM-related information. The first one, 'ndata', contains DOF numbers and other things for each node. The second table, 'edata', holds element orders and bubble DOF numbers. Both tables are directly indexed by the node and element IDs. The function resize_tables() is called to check whether the tables are large enough to contain all node and element id's, and to reallocate them if not.
Definition at line 206 of file space.cpp.
Referenced by Hermes::Hermes2D::Space< Scalar >::assign_dofs(), and Hermes::Hermes2D::Space< Scalar >::load().
|
virtual |
Sets element polynomial order. Can be called by the user. Should not be called for many elements at once, since assign_dofs() is called at the end of this function.
|
protectedvirtual |
Sets element polynomial order. This version does not call assign_dofs() and is intended primarily for internal use.
Definition at line 325 of file space.cpp.
Referenced by Hermes::Hermes2D::Space< Scalar >::distribute_orders().
void Hermes::Hermes2D::Space< Scalar >::set_uniform_order | ( | int | order, |
std::string | marker = HERMES_ANY |
||
) |
Sets the same polynomial order for all elements in the mesh. Intended for the user and thus assign_dofs() is called at the end of this function.
void Hermes::Hermes2D::Space< Scalar >::set_uniform_order_internal | ( | int | order, |
int | marker | ||
) |
Sets the same polynomial order for all elements in the mesh. Does not call assign_dofs(). For internal use.
Definition at line 401 of file space.cpp.
Referenced by Hermes::Hermes2D::H1Space< Scalar >::init().
void Hermes::Hermes2D::Space< Scalar >::unrefine_all_mesh_elements | ( | bool | keep_initial_refinements = true | ) |
Recursively removes all son elements of the given element and makes it active. Also handles element orders.
[in] | keep_initial_refinements | Refinements in Mesh can be marked as initial (to prevent taking them back), this parameter serves to prevent taking them back with this method. |
|
static |
|
protected |
Recursively removes all son elements of the given element and makes it active. Also handles element orders.
[in] | keep_initial_refinements | Refinements in Mesh can be marked as initial (to prevent taking them back), this parameter serves to prevent taking them back with this method. |
[in] | only_unrefine_space_data | Useful when more spaces share the mesh if one wants to unrefine the underlying Mesh only once, but wants other spaces know about the change. |
|
protectedvirtual |
Called by Space to update constraining relationships between shape functions due to hanging nodes in the mesh. As this is space-specific, this function is reimplemented in H1Space and HcurlSpace.
Reimplemented in Hermes::Hermes2D::H1Space< Scalar >, Hermes::Hermes2D::HdivSpace< Scalar >, and Hermes::Hermes2D::HcurlSpace< Scalar >.
Definition at line 830 of file space.cpp.
Referenced by Hermes::Hermes2D::Space< Scalar >::assign_dofs().
void Hermes::Hermes2D::Space< Scalar >::update_essential_bc_values | ( | ) |
Updates essential BC values. Typically used for time-dependent essential boundary conditions.
Definition at line 1068 of file space.cpp.
Referenced by Hermes::Hermes2D::Space< Scalar >::assign_dofs().