| 
    GetFEM
    5.4.4
    
   | 
 
Describe a finite element method linked to a mesh. More...
#include <getfem_mesh_fem.h>
Inherits getfem::context_dependencies, and dal::static_stored_object.
Inherited by getfem::mesh_fem_global_function, getfem::mesh_fem_level_set, getfem::mesh_fem_product, getfem::mesh_fem_sum, getfem::partial_mesh_fem, and getfem::torus_mesh_fem.
Public Member Functions | |
| void | update_from_context () const | 
| this function has to be defined and should update the object when the context is modified.  | |
| const dal::bit_vector & | convex_index () const | 
| Get the set of convexes where a finite element has been assigned.  | |
| bool | is_reduced () const | 
| Return true if a reduction matrix is applied to the dofs.  | |
| const REDUCTION_MATRIX & | reduction_matrix () const | 
| Return the reduction matrix applied to the dofs.  | |
| const EXTENSION_MATRIX & | extension_matrix () const | 
| Return the extension matrix corresponding to reduction applied (RE=I).  | |
| template<typename MATR , typename MATE > | |
| void | set_reduction_matrices (const MATR &RR, const MATE &EE) | 
| Allows to set the reduction and the extension matrices.  More... | |
| void | reduce_to_basic_dof (const dal::bit_vector &kept_basic_dof) | 
| Allows to set the reduction and the extension matrices in order to keep only a certain number of dof.  | |
| void | set_reduction (bool r) | 
| Validate or invalidate the reduction (keeping the reduction matrices).  | |
| const mesh & | linked_mesh () const | 
| Return a reference to the underlying mesh.  | |
| void | set_auto_add (dim_type K, bool disc=false, scalar_type alpha=scalar_type(0), bool complete=false) | 
| Set the degree of the fem for automatic addition of element option.  More... | |
| void | set_auto_add (pfem pf) | 
| Set the fem for automatic addition of element option.  More... | |
| virtual dim_type | get_qdim () const | 
| Return the Q dimension.  More... | |
| virtual void | set_qdim (dim_type q) | 
| Change the Q dimension.  | |
| virtual void | set_qdim (dim_type M, dim_type N) | 
| Set the dimension for a matrix field.  | |
| virtual void | set_qdim (dim_type M, dim_type N, dim_type O, dim_type P) | 
| Set the dimension for a fourth order tensor field.  | |
| virtual void | set_qdim (const bgeot::multi_index &mii) | 
| Set the dimension for an arbitrary order tensor field.  | |
| dim_type | get_qdim_m () const IS_DEPRECATED | 
| for matrix fields, return the number of rows.  | |
| dim_type | get_qdim_n () const IS_DEPRECATED | 
| for matrix fields, return the number of columns.  | |
| void | set_finite_element (size_type cv, pfem pf) | 
| Set the finite element method of a convex.  More... | |
| void | set_finite_element (const dal::bit_vector &cvs, pfem pf) | 
| Set the finite element on a set of convexes.  More... | |
| void | set_finite_element (pfem pf) | 
| shortcut for set_finite_element(linked_mesh().convex_index(),pf); and set_auto_add(pf).  | |
| void | set_classical_finite_element (size_type cv, dim_type fem_degree, bool complete=false) | 
| Set a classical (i.e.  More... | |
| void | set_classical_finite_element (const dal::bit_vector &cvs, dim_type fem_degree, bool complete=false) | 
| Set a classical (i.e.  More... | |
| void | set_classical_discontinuous_finite_element (size_type cv, dim_type fem_degree, scalar_type alpha=0, bool complete=false) | 
| Similar to set_classical_finite_element, but uses discontinuous lagrange elements.  More... | |
| void | set_classical_discontinuous_finite_element (const dal::bit_vector &cvs, dim_type fem_degree, scalar_type alpha=0, bool complete=false) | 
| Similar to set_classical_finite_element, but uses discontinuous lagrange elements.  More... | |
| void | set_classical_finite_element (dim_type fem_degree, bool complete=false) | 
| Shortcut for set_classical_finite_element(linked_mesh().convex_index(),...)  | |
| void | set_classical_discontinuous_finite_element (dim_type fem_degree, scalar_type alpha=0, bool complete=false) | 
| Shortcut for set_classical_discontinuous_finite_element(linked_mesh().convex_index() ,...)  | |
| virtual pfem | fem_of_element (size_type cv) const | 
| Return the basic fem associated with an element (if no fem is associated, the function will crash! use the convex_index() of the mesh_fem to check that a fem is associated to a given convex).  More... | |
| virtual ind_dof_ct | ind_basic_dof_of_element (size_type cv) const | 
| Give an array of the dof numbers a of convex.  More... | |
| virtual ind_dof_face_ct | ind_basic_dof_of_face_of_element (size_type cv, short_type f) const | 
| Give an array of the dof numbers lying of a convex face (all degrees of freedom whose associated base function is non-zero on the convex face).  More... | |
| virtual size_type | nb_basic_dof_of_face_of_element (size_type cv, short_type f) const | 
| Return the number of dof lying on the given convex face.  More... | |
| virtual size_type | nb_basic_dof_of_element (size_type cv) const | 
| Return the number of degrees of freedom attached to a given convex.  More... | |
| virtual base_node | point_of_basic_dof (size_type cv, size_type i) const | 
| Return the geometrical location of a degree of freedom.  More... | |
| virtual base_node | point_of_basic_dof (size_type d) const | 
| Return the geometrical location of a degree of freedom.  More... | |
| virtual dim_type | basic_dof_qdim (size_type d) const | 
| Return the dof component number (0<= x <Qdim)  | |
| virtual size_type | first_convex_of_basic_dof (size_type d) const | 
| Shortcut for convex_to_dof(d)[0].  More... | |
| virtual const mesh::ind_cv_ct & | convex_to_basic_dof (size_type d) const | 
| Return the list of convexes attached to the specified dof.  More... | |
| virtual void | get_global_dof_index (std::vector< size_type > &ind) const | 
| Give an array that contains the global dof indices corresponding to the mesh_fem dofs or size_type(-1) if a dof is not global.  More... | |
| virtual void | enumerate_dof () const | 
| Renumber the degrees of freedom.  More... | |
| virtual size_type | nb_basic_dof () const | 
| Return the total number of basic degrees of freedom (before the optional reduction).  | |
| virtual size_type | nb_dof () const | 
| Return the total number of degrees of freedom.  | |
| virtual dal::bit_vector | basic_dof_on_region (const mesh_region &b) const | 
| Get a list of basic dof lying on a given mesh_region.  More... | |
| dal::bit_vector | dof_on_region (const mesh_region &b) const | 
| Get a list of dof lying on a given mesh_region.  More... | |
| mesh_fem (const mesh &me, dim_type Q=1) | |
| Build a new mesh_fem.  More... | |
| virtual void | read_from_file (std::istream &ist) | 
| Read the mesh_fem from a stream.  More... | |
| void | read_from_file (const std::string &name) | 
| Read the mesh_fem from a file.  More... | |
| virtual void | write_to_file (std::ostream &ost) const | 
| Write the mesh_fem to a stream.  | |
| void | write_to_file (const std::string &name, bool with_mesh=false) const | 
| Write the mesh_fem to a file.  More... | |
  Public Member Functions inherited from getfem::context_dependencies | |
| bool | context_check () const | 
| return true if update_from_context was called  | |
Describe a finite element method linked to a mesh.
Definition at line 152 of file getfem_mesh_fem.h.
      
  | 
  explicit | 
Build a new mesh_fem.
A mesh object must be supplied.
| me | the linked mesh. | 
| Q | the Q dimension (see mesh_fem::get_qdim). | 
Definition at line 527 of file getfem_mesh_fem.cc.
      
  | 
  inline | 
Allows to set the reduction and the extension matrices.
Should satify (RR*EE=I).
Definition at line 210 of file getfem_mesh_fem.h.
      
  | 
  inline | 
Set the degree of the fem for automatic addition of element option.
K=-1 disables the automatic addition.
Definition at line 291 of file getfem_mesh_fem.h.
      
  | 
  inline | 
Set the fem for automatic addition of element option.
pf=0 disables the automatic addition.
Definition at line 304 of file getfem_mesh_fem.h.
      
  | 
  inlinevirtual | 
Return the Q dimension.
A mesh_fem used for scalar fields has Q=1, for vector fields, Q is typically equal to linked_mesh().dim().
Reimplemented in getfem::partial_mesh_fem.
Definition at line 316 of file getfem_mesh_fem.h.
| void getfem::mesh_fem::set_finite_element | ( | size_type | cv, | 
| pfem | pf | ||
| ) | 
Set the finite element method of a convex.
| cv | the convex number. | 
| pf | the finite element. | 
Definition at line 127 of file getfem_mesh_fem.cc.
| void getfem::mesh_fem::set_finite_element | ( | const dal::bit_vector & | cvs, | 
| pfem | pf | ||
| ) | 
Set the finite element on a set of convexes.
| cvs | the set of convex indexes, as a dal::bit_vector. | 
| pf | the finite element, typically obtained with  getfem::fem_descriptor("FEM_SOMETHING(..)") 
pfem fem_descriptor(const std::string &name) get a fem descriptor from its string name. Definition: getfem_fem.cc:4660  | 
Definition at line 166 of file getfem_mesh_fem.cc.
| void getfem::mesh_fem::set_classical_finite_element | ( | size_type | cv, | 
| dim_type | fem_degree, | ||
| bool | complete = false  | 
        ||
| ) | 
Set a classical (i.e.
lagrange polynomial) finite element on a convex.
| cv | is the convex number. | 
| fem_degree | the polynomial degree of the finite element. | 
| complete | is a flag for excluding incomplete element irrespective of the element geometric transformation. | 
Definition at line 174 of file getfem_mesh_fem.cc.
| void getfem::mesh_fem::set_classical_finite_element | ( | const dal::bit_vector & | cvs, | 
| dim_type | fem_degree, | ||
| bool | complete = false  | 
        ||
| ) | 
Set a classical (i.e.
lagrange polynomial) finite element on a set of convexes.
| cvs | the set of convexes, as a dal::bit_vector. | 
| fem_degree | the polynomial degree of the finite element. | 
Definition at line 182 of file getfem_mesh_fem.cc.
| void getfem::mesh_fem::set_classical_discontinuous_finite_element | ( | size_type | cv, | 
| dim_type | fem_degree, | ||
| scalar_type | alpha = 0,  | 
        ||
| bool | complete = false  | 
        ||
| ) | 
Similar to set_classical_finite_element, but uses discontinuous lagrange elements.
| cv | the convex index. | 
| fem_degree | the polynomial degree of the finite element. | 
| alpha | the node inset, 0 <= alpha < 1, where 0 implies usual dof nodes, greater values move the nodes toward the center of gravity, and 1 means that all degrees of freedom collapse on the center of gravity. | 
Definition at line 199 of file getfem_mesh_fem.cc.
| void getfem::mesh_fem::set_classical_discontinuous_finite_element | ( | const dal::bit_vector & | cvs, | 
| dim_type | fem_degree, | ||
| scalar_type | alpha = 0,  | 
        ||
| bool | complete = false  | 
        ||
| ) | 
Similar to set_classical_finite_element, but uses discontinuous lagrange elements.
| cvs | the set of convexes, as a dal::bit_vector. | 
| fem_degree | the polynomial degree of the finite element. | 
| alpha | the node inset, 0 <= alpha < 1, where 0 implies usual dof nodes, greater values move the nodes toward the center of gravity, and 1 means that all degrees of freedom collapse on the center of gravity. | 
Definition at line 206 of file getfem_mesh_fem.cc.
      
  | 
  inlinevirtual | 
Return the basic fem associated with an element (if no fem is associated, the function will crash! use the convex_index() of the mesh_fem to check that a fem is associated to a given convex).
This fem does not take into account the optional vectorization due to qdim nor the optional reduction.
Reimplemented in getfem::partial_mesh_fem.
Definition at line 448 of file getfem_mesh_fem.h.
      
  | 
  inlinevirtual | 
Give an array of the dof numbers a of convex.
| cv | the convex number. | 
Reimplemented in getfem::partial_mesh_fem.
Definition at line 454 of file getfem_mesh_fem.h.
      
  | 
  inlinevirtual | 
Give an array of the dof numbers lying of a convex face (all degrees of freedom whose associated base function is non-zero on the convex face).
| cv | the convex number. | 
| f | the face number. | 
Reimplemented in getfem::partial_mesh_fem.
Definition at line 472 of file getfem_mesh_fem.h.
      
  | 
  inlinevirtual | 
Return the number of dof lying on the given convex face.
| cv | the convex number. | 
| f | the face number. | 
Reimplemented in getfem::partial_mesh_fem.
Definition at line 485 of file getfem_mesh_fem.h.
      
  | 
  inlinevirtual | 
Return the number of degrees of freedom attached to a given convex.
| cv | the convex number. | 
Reimplemented in getfem::partial_mesh_fem.
Definition at line 498 of file getfem_mesh_fem.h.
      
  | 
  virtual | 
Return the geometrical location of a degree of freedom.
| cv | the convex number. | 
| i | the local dof number. | 
Reimplemented in getfem::partial_mesh_fem.
Definition at line 223 of file getfem_mesh_fem.cc.
      
  | 
  virtual | 
Return the geometrical location of a degree of freedom.
| d | the global dof number. | 
Reimplemented in getfem::partial_mesh_fem.
Definition at line 231 of file getfem_mesh_fem.cc.
      
  | 
  virtual | 
Shortcut for convex_to_dof(d)[0].
| d | the global dof number. | 
Reimplemented in getfem::partial_mesh_fem.
Definition at line 258 of file getfem_mesh_fem.cc.
      
  | 
  virtual | 
Return the list of convexes attached to the specified dof.
| d | the global dof number. | 
Reimplemented in getfem::partial_mesh_fem.
Definition at line 267 of file getfem_mesh_fem.cc.
      
  | 
  virtual | 
Give an array that contains the global dof indices corresponding to the mesh_fem dofs or size_type(-1) if a dof is not global.
| ind | the returned global dof indices array. | 
Definition at line 293 of file getfem_mesh_fem.cc.
      
  | 
  virtual | 
Renumber the degrees of freedom.
Enumeration of dofs.
You should not have to call this function, as it is done automatically
Reimplemented in getfem::torus_mesh_fem.
Definition at line 321 of file getfem_mesh_fem.cc.
      
  | 
  virtual | 
Get a list of basic dof lying on a given mesh_region.
| b | the mesh_region. | 
Reimplemented in getfem::partial_mesh_fem.
Definition at line 77 of file getfem_mesh_fem.cc.
| dal::bit_vector getfem::mesh_fem::dof_on_region | ( | const mesh_region & | b | ) | const | 
Get a list of dof lying on a given mesh_region.
For a reduced mesh_fem a dof is lying on a region if its potential corresponding shape function is nonzero on this region. The extension matrix is used to make the correspondence between basic and reduced dofs.
| b | the mesh_region. | 
Definition at line 114 of file getfem_mesh_fem.cc.
      
  | 
  virtual | 
Read the mesh_fem from a stream.
| ist | the stream. | 
Reimplemented in getfem::partial_mesh_fem.
Definition at line 539 of file getfem_mesh_fem.cc.
| void getfem::mesh_fem::read_from_file | ( | const std::string & | name | ) | 
Read the mesh_fem from a file.
| name | the file name. | 
Definition at line 703 of file getfem_mesh_fem.cc.
| void getfem::mesh_fem::write_to_file | ( | const std::string & | name, | 
| bool | with_mesh = false  | 
        ||
| ) | const | 
Write the mesh_fem to a file.
| name | the file name | 
| with_mesh | if set, then the linked_mesh() will also be saved to the file. | 
Definition at line 785 of file getfem_mesh_fem.cc.