Hermes2D
2.0
|
#include <space_h1.h>
Classes | |
struct | EdgeInfo |
struct | FixedVertex |
Public Member Functions | |
H1Space (const Mesh *mesh, EssentialBCs< Scalar > *boundary_conditions, int p_init=1, Shapeset *shapeset=NULL) | |
H1Space (const Mesh *mesh, int p_init=1, Shapeset *shapeset=NULL) | |
virtual void | set_shapeset (Shapeset *shapeset) |
Sets the shapeset. | |
void | fix_vertex (int id, Scalar value=0.0) |
virtual Scalar * | get_bc_projection (SurfPos *surf_pos, int order, EssentialBoundaryCondition< Scalar > *bc) |
virtual void | copy (const Space< Scalar > *space, Mesh *new_mesh) |
Copy from Space instance 'space'. | |
![]() | |
Space (const Mesh *mesh, Shapeset *shapeset, EssentialBCs< Scalar > *essential_bcs) | |
void | init () |
Common code for constructors. | |
virtual bool | isOkay () const |
State querying helpers. | |
std::string | getClassName () const |
Get class name, for the purpose of messaging. | |
virtual | ~Space () |
Destructor. | |
virtual void | set_element_order (int id, int order) |
virtual void | set_element_orders (int *elem_orders) |
Sets polynomial order to all elements. | |
int | get_element_order (int id) const |
Returns element polynomial order. | |
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 | unrefine_all_mesh_elements (bool keep_initial_refinements=true) |
void | update_element_orders_after_refinement () |
Updates element orders when the underlying mesh has been refined. | |
int | get_num_dofs () const |
Returns the number of basis functions contained in the space. | |
Mesh * | get_mesh () const |
void | set_mesh (Mesh *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. | |
EssentialBCs< Scalar > * | get_essential_bcs () const |
Obtains an boundary conditions. | |
void | update_essential_bc_values () |
Shapeset * | get_shapeset () const |
bool | save (const char *filename) const |
Saves this space into a file. | |
virtual void | get_element_assembly_list (Element *e, AsmList< Scalar > *al, unsigned int first_dof=0) const |
Obtains an assembly list for the given element. | |
virtual void | set_element_order_internal (int id, int order) |
virtual int | assign_dofs (int first_dof=0, int stride=1) |
Builds basis functions and assigns DOF numbers to them. More... | |
void | distribute_orders (Mesh *mesh, int *parents) |
Sets polynomial orders to elements created by Mesh::regularize() using "parents". | |
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. | |
bool | is_up_to_date () const |
Returns true if the space is ready for computation, false otherwise. | |
void | get_boundary_assembly_list (Element *e, int surf_num, AsmList< Scalar > *al, unsigned int first_dof=0) const |
Obtains an edge assembly list (contains shape functions that are nonzero on the specified edge). | |
void | set_uniform_order_internal (int order, int marker) |
void | free () |
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. | |
template<> | |
void | init () |
template<> | |
void | init () |
template<> | |
void | free () |
template<> | |
void | free () |
template<> | |
Space (const Mesh *mesh, Shapeset *shapeset, EssentialBCs< double > *essential_bcs) | |
template<> | |
Space (const Mesh *mesh, Shapeset *shapeset, EssentialBCs< std::complex< double > > *essential_bcs) | |
![]() | |
void | check () const |
Method to handle the state. | |
![]() | |
XMLParsing () | |
Constructor. | |
void | set_validation (bool to_set) |
Set to validate / not to validate. | |
Protected Member Functions | |
virtual SpaceType | get_type () const |
void | init (Shapeset *shapeset, int p_init) |
Common code for the constructors. | |
virtual void | assign_vertex_dofs () |
virtual void | assign_edge_dofs () |
virtual void | assign_bubble_dofs () |
virtual void | get_vertex_assembly_list (Element *e, int iv, AsmList< Scalar > *al) const |
virtual void | get_boundary_assembly_list_internal (Element *e, int ie, AsmList< Scalar > *al) const |
void | output_component (typename Space< Scalar >::BaseComponent *¤t, typename Space< Scalar >::BaseComponent *&last, typename Space< Scalar >::BaseComponent *min, Node *&edge, typename Space< Scalar >::BaseComponent *&edge_dofs) |
Space< Scalar >::BaseComponent * | merge_baselists (typename Space< Scalar >::BaseComponent *l1, int n1, typename Space< Scalar >::BaseComponent *l2, int n2, Node *edge, typename Space< Scalar >::BaseComponent *&edge_dofs, int &ncomponents) |
void | update_constrained_nodes (Element *e, EdgeInfo *ei0, EdgeInfo *ei1, EdgeInfo *ei2, EdgeInfo *ei3) |
virtual void | update_constraints () |
bool | is_fixed_vertex (int id) const |
virtual void | post_assign () |
![]() | |
virtual int | get_edge_order_internal (Node *en) const |
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 | 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) |
void | free_bc_data () |
int | get_seq () const |
Internal. Used by DiscreteProblem to detect changes in the space. | |
Protected Attributes | |
Hermes::vector< FixedVertex > | fixed_vertices |
![]() | |
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. | |
const Mesh * | mesh |
FE mesh. | |
int | default_tri_order |
int | default_quad_order |
int | vertex_functions_count |
int | edge_functions_count |
int | bubble_functions_count |
int | first_dof |
int | next_dof |
int | stride |
int | seq |
int | mesh_seq |
int | was_assigned |
NodeData * | ndata |
node data table | |
ElementData * | edata |
element data table | |
int | nsize |
int | ndata_allocated |
number of items in ndata, allocated space | |
int | esize |
double ** | proj_mat |
double * | chol_p |
Hermes::vector< void * > | bc_data |
Used for bc projection. | |
![]() | |
bool | validate |
Internal. | |
Friends | |
class | Space< Scalar > |
Additional Inherited Members | |
![]() | |
static int | get_num_dofs (Hermes::vector< const Space< Scalar > * > spaces) |
Returns the number of basis functions contained in the spaces. | |
static int | get_num_dofs (Hermes::vector< Space< Scalar > * > spaces) |
static int | get_num_dofs (const Space< Scalar > *space) |
Returns the number of basis functions contained in the space. | |
static int | get_num_dofs (Space< Scalar > *space) |
static Space< Scalar > * | load (const char *filename, Mesh *mesh, bool validate, EssentialBCs< Scalar > *essential_bcs=NULL, Shapeset *shapeset=NULL) |
Loads a space from a file. | |
static int | assign_dofs (Hermes::vector< Space< Scalar > * > spaces) |
Assings the degrees of freedom to all Spaces in the Hermes::vector. | |
static void | update_essential_bc_values (Hermes::vector< Space< Scalar > * > spaces, double time) |
static void | update_essential_bc_values (Space< Scalar > *s, double time) |
static Node * | get_mid_edge_vertex_node (Element *e, int i, int j) |
![]() | |
static const int | H2D_UNASSIGNED_DOF = -2 |
DOF which was not assigned yet. | |
static const int | H2D_CONSTRAINED_DOF = -1 |
DOF which is constrained. | |
H1Space represents a space of continuous Scalar functions over a domain (mesh).
Typical usage:
...
Hermes::Hermes2D::EssentialBCs<double> bcs(&bc_essential1, &bc_essential2, ...);
// Initialize space.
int globalPolynomialOrder = 4;
Hermes::Hermes2D::H1Space<double> space(&mesh, &bcs, globalPolynomialOrder);
void Hermes::Hermes2D::H1Space< Scalar >::fix_vertex | ( | int | id, |
Scalar | value = 0.0 |
||
) |
Removes the degree of freedom from a vertex node with the given id (i.e., its number in the mesh file) and makes it part of the Dirichlet lift with the given value. This is a special-purpose function which normally should not be needed. It is intended for fixing the solution of a system which would otherwise be singular and for some reason a standard Dirichlet condition (with non-zero measure on the boundary) is not suitable.
Definition at line 575 of file space_h1.cpp.
|
inlineprotectedvirtual |
Internal. Return type of this space (H1 = HERMES_H1_SPACE, Hcurl = HERMES_HCURL_SPACE, Hdiv = HERMES_HDIV_SPACE, L2 = HERMES_L2_SPACE)
Implements Hermes::Hermes2D::Space< Scalar >.
Definition at line 65 of file space_h1.h.
|
protectedvirtual |
Auxiliary function the descendants may implement to perform additional tasks after the DOFs have been assigned.
Reimplemented from Hermes::Hermes2D::Space< Scalar >.
Definition at line 592 of file space_h1.cpp.
|
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 from Hermes::Hermes2D::Space< Scalar >.
Definition at line 567 of file space_h1.cpp.