Selecting finite element methods
The description of a finite element method on a whole mesh is done thanks to the structure getfem::mesh_fem, defined in the file getfem/getfem_mesh_fem.h. Basically, this structure describes the finite element method on each element of the mesh. It is possible to have an arbitrary number of finite element description for a single mesh. This is particularly necessary for mixed methods, but also to describe different data on the same mesh. One can instantiate a getfem::mesh_fem object as follows
getfem::mesh_fem mf(mymesh);
where mymesh is an already
existing mesh. The structure will be linked to this mesh and will
react when modifications will be done on it.
It is possible to specify element by element the finite element
method, so that element of mixed types can be treated, even if the
dimensions are different. For usual elements, the connection
between two elements is done when the two elements are compatibles
(same degrees of freedom on the common face). A numeration of the
degrees of freedom is automatically done with a Cuthill Mc Kee like
algorithm. You have to keep in mind that there is absolutely no
connection between the numeration of vertices of the mesh and the
numeration of the degrees of freedom. Every getfem::mesh_fem object has its own
numeration.
To select a particular finite element method on a given element,
one can use
mf.set_finite_element(i, pf);
where i is the index of the element and pf is the descriptor (of type getfem::pfem, basically a pointer to an object which inherits from getfem::virtual_fem) of the finite element method. Alternative forms of this member function are:
void mesh_fem::set_finite_element(const dal::bit_vector &cvs,
getfem::pfem pf);
void mesh_fem::set_finite_element(getfem::pfem pf);
which set the finite elements for either the convexes listed in the bit_vector cvs, or all the convexes of the mesh.
Descriptors for finite element methods and integration methods
are available thanks to the following function
getfem::pfem pf = getfem::fem_descriptor("name of method");
where "name of method" is to be
chosen among the existing methods. A name of a method can be
retrieved thanks to the following functions
std::string femname = getfem::name_of_fem(pf);
A non exhautive list (see Appendix A or getfem/getfem_fem.h for exhaustive lists) of finite element methods is given by
| "FEM_PK(n,k)" | Classical PK methods on simplexes of dimension n with degree k polynomials. |
| "FEM_QK(n,k)" | Classical QK methods on parallelepiped of dimension n. Tensorial product of degree k PK method on the segment. |
| "FEM_PK_PRISM(n,k)" | Classical methods on prism of dimension n. Tensorial product of two degree k PK method. |
| "FEM_PRODUCT(a,b)" | Tensorial product of the two polynomial finite element method a and b. |
| "FEM_PK_DISCONTINUOUS(n,k)" | discontinuous PK methods on simplexes of dimension n with degree k polynomials. |
getfem::pfem getfem::classical_fem(bgeot::pgeometric_trans pg, short_type degree); pfem getfem::classical_discontinuous_fem(bgeot::pgeometric_trans pg, short_type degree);
The mesh_fem can call directly these functions via:
void mesh_fem::set_classical_finite_element(const dal::bit_vector &cvs,
dim_type fem_degree);
void mesh_fem::set_classical_discontinuous_finite_element
(const dal::bit_vector &cvs, dim_type fem_degree);
void mesh_fem::set_classical_finite_element(dim_type fem_degree);
void mesh_fem::set_classical_discontinuous_finite_element(dim_type fem_degree);
Examples
For instance if one needs to have a description of a P1 finite element method on a triangle, the way to set it is
mf.set_finite_element(i, getfem::fem_descriptor("FEM_PK(2, 1)"));
where i is still the index of the triangle. It is also possible to select a particular method directly on a set of element, passing to mf.set_finite_element a dal::bit_vector instead of a single index. For instance
mf.set_finite_element(mymesh.convex_index(),
getfem::fem_descriptor("FEM_PK(2, 1)"));
selects the method on all the elements of the mesh.
IMPORTANT: If the finite element represents an unknown for a
vectorial problem (such as the displacement for linear elasticity),
one should use mf.set_qdim(Q) to set
the target dimension (see [] for the
definition of the target dimension Q).
If the target dimension Q is set to a value different of
1, the scalar FEMs (such as Pk fems etc.)
are automatically "vectorized", i.e. each scalar degree of freedom
is duplicated Q times, from the mesh_fem object point of vue. To sum it
up,
- if the fem of the ith element is intrinsically a vector
FEM,
mf.get_qdim() == mf.fem_of_element(i)->target_dim()
and
mf.nb_dof_of_element(i) == mf.fem_of_element(i).nb_dof(). - if the fem has a target_dim
equal to 1,
mf.nb_dof_of_element(i) == mf.get_qdim()*mf.fem_of_element(i).nb_dof().
Other methods of the mesh_fem object
Once a finite element method is defined on a mesh, it is
possible to obtain information on it with the following methods
(the list is not exhaustive).
| mf.convex_index() | Set of indexes (a dal::bit_vector) on which a finite element method is defined. |
| mf.linked_mesh() | gives a reference to the linked mesh. |
| mf.fem_of_element(i) | gives a descriptor on the finite element method defined on element of index i. |
| mf.nb_dof_of_element(i) | gives the number of degrees of freedom on the element of index i. |
| mf.ind_dof_of_element(i) | gives a container (an array) with all the global indexes of the degrees of freedom of element of index i. |
| mf.point_of_dof(i, j) | gives a bgeot::base_node which represents the point associated with the dof of local index j on element of index i. |
| mf.point_of_dof(j) | gives a bgeot::base_node which represents the point associated with the dof of global index j. |
| mf.reference_point_of_dof(i, j) | gives a bgeot::base_node which represents the point associated with the dof of local index j on element of index i in the coordinates of the reference element. |
| mf.first_convex_of_dof(j) | gives the index of the first element on which the degree of freedom of global index j is defined. |
| mf.nb_dof() | gives the total number of different degrees of freedom. |
| mf.get_qdim() | gives the target dimension Q. |
| mf.clear() | Clear the structure, no finite element method is still defined. |
| mf.dof_on_set(i) | Return a dal::bit_vector which represent the indices of dof which are in the set of convexes or the set of faces of index i (see the getfem::mesh object). |
Obtaining generic mesh_fems
It is possible to use the function
const mesh_fem &getfem::classical_mesh_fem(const getfem::mesh &mymesh, dim_type K);
to get a classical polynomial mesh_fem of order K on the given mymesh. The returned mesh_fem will be destroyed automatically when its linked_mesh is destroyed. All the mesh_fem built by this function are stored in a cache, which means that calling this function twice with the same arguments will return the same mesh_fem object. A consequence is that you should NEVER modify this mesh_fem!
