GetFEM  5.4.2
getfem_fem.h File Reference

Definition of the finite element methods. More...

#include "dal_static_stored_objects.h"
#include "bgeot_geometric_trans.h"
#include "bgeot_poly_composite.h"
#include "getfem_integration.h"
#include "dal_naming_system.h"
#include <deque>

Go to the source code of this file.

Classes

class  getfem::virtual_fem
 Base class for finite element description. More...
 
class  getfem::fem< FUNC >
 virtual_fem implementation as a vector of generic functions. More...
 
class  getfem::fem_precomp_
 Pre-computations on a fem (given a fixed set of points on the reference convex, this object computes the value/gradient/hessian of all base functions on this set of points and stores them. More...
 
class  getfem::fem_precomp_pool
 handle a pool (i.e. More...
 
class  getfem::fem_interpolation_context
 structure passed as the argument of fem interpolation functions. More...
 

Namespaces

 getfem
 GEneric Tool for Finite Element Methods.
 

Typedefs

typedef dof_description * getfem::pdof_description
 Type representing a pointer on a dof_description.
 
typedef std::shared_ptr< const getfem::virtual_femgetfem::pfem
 type of pointer on a fem description More...
 

Functions

pdof_description getfem::lagrange_dof (dim_type d)
 Description of a unique dof of lagrange type (value at the node). More...
 
pdof_description getfem::derivative_dof (dim_type d, dim_type r)
 Description of a unique dof of derivative type. More...
 
pdof_description getfem::second_derivative_dof (dim_type d, dim_type num_der1, dim_type num_der2)
 Description of a unique dof of second derivative type. More...
 
pdof_description getfem::normal_derivative_dof (dim_type d)
 Description of a unique dof of normal derivative type (normal derivative at the node, regarding a face). More...
 
pdof_description getfem::mean_value_dof (dim_type d)
 Description of a unique dof of mean value type. More...
 
pdof_description getfem::global_dof (dim_type d)
 Description of a global dof, i.e. More...
 
pdof_description getfem::product_dof (pdof_description pnd1, pdof_description pnd2)
 Product description of the descriptions *pnd1 and *pnd2.
 
pdof_description getfem::xfem_dof (pdof_description p, size_type ind)
 Description of a special dof for Xfem.
 
size_type getfem::dof_xfem_index (pdof_description)
 Returns the xfem_index of dof (0 for normal dof)
 
int getfem::dof_description_compare (pdof_description a, pdof_description b)
 Gives a total order on the dof description compatible with the identification.
 
bool getfem::dof_linkable (pdof_description)
 Says if the dof is linkable.
 
bool getfem::dof_compatibility (pdof_description, pdof_description)
 Says if the two dofs can be identified.
 
pfem getfem::classical_fem (bgeot::pgeometric_trans pgt, short_type k, bool complete=false)
 Give a pointer on the structures describing the classical polynomial fem of degree k on a given convex type. More...
 
pfem getfem::classical_discontinuous_fem (bgeot::pgeometric_trans pg, short_type k, scalar_type alpha=0, bool complete=false)
 Give a pointer on the structures describing the classical polynomial discontinuous fem of degree k on a given convex type. More...
 
pfem getfem::fem_descriptor (const std::string &name)
 get a fem descriptor from its string name.
 
std::string getfem::name_of_fem (pfem p)
 get the string name of a fem descriptor.
 
pfem_precomp getfem::fem_precomp (pfem pf, bgeot::pstored_point_tab pspt, dal::pstatic_stored_object dep)
 Handles precomputations for FEM. More...
 
void getfem::delete_fem_precomp (pfem_precomp pfp)
 Request for the removal of a pfem_precomp.
 
pfem getfem::interior_fem_of_hho_method (pfem hho_method)
 Specific function for a HHO method to obtain the method in the interior. More...
 

Variables

const typedef fem< bgeot::base_poly > * getfem::ppolyfem
 Classical polynomial FEM.
 
const typedef fem< bgeot::polynomial_composite > * getfem::ppolycompfem
 Polynomial composite FEM.
 
const typedef fem< bgeot::base_rational_fraction > * getfem::prationalfracfem
 Rational fration FEM.
 

Detailed Description

Definition of the finite element methods.

Author
Yves Renard Yves..nosp@m.Rena.nosp@m.rd@in.nosp@m.sa-l.nosp@m.yon.f.nosp@m.r
Date
December 21, 1999. This file defines the getfem::virtual_fem class, which is the common base class of all FEM.

List of FEM known by getfem::fem_descriptor :

  • "FEM_PK(N,K)" : classical Lagrange element PK on a simplex.
  • "FEM_PK_DISCONTINUOUS(N,K,alpha)" : discontinuous Lagrange element PK on a simplex.
  • "FEM_QK(N,K)" : classical Lagrange element QK on a parellepiped.
  • "FEM_QK_DISCONTINUOUS(N,K,alpha)" : discontinuous Lagrange element QK on a parallelepiped.
  • "FEM_Q2_INCOMPLETE(N)" : incomplete Q2 elements with 8 and 20 dof (serendipity Quad 8 and Hexa 20 elements)
  • "FEM_Q2_INCOMPLETE_DISCONTINUOUS(N)" : discontinuous incomplete Q2 elements with 8 and 20 dof (serendipity Quad 8 and Hexa 20 elements)
  • "FEM_PRISM_PK(N,K)" : classical Lagrange element PK on a prism.
  • "FEM_PRISM_PK_DISCONTINUOUS(N,K,alpha)" : classical discontinuous Lagrange element PK on a prism.
  • "FEM_PRISM_INCOMPLETE_P2" : Incomplete Lagrange element on a quadratic 3D prism (serendipity, 15-node wedge element). Can be connected toa standard P2 Lagrange on its triangular faces and a Q2_INCOMPLETE Lagrange element on its quadrangular faces.
  • "FEM_PK_WITH_CUBIC_BUBBLE(N,K)" : classical Lagrange element PK on a simplex with an additional volumic bubble function.
  • "FEM_PRODUCT(FEM1,FEM2)" : tensorial product of two polynomial elements
  • "FEM_P1_NONCONFORMING" : Nonconforming P1 method on a triangle.
  • "FEM_P1_BUBBLE_FACE(N)" : P1 method on a simplex with an additional bubble function on face 0.
  • "FEM_P1_BUBBLE_FACE_LAG" : P1 method on a simplex with an additional lagrange dof on face 0.
  • "FEM_HERMITE(N)" : Hermite element P3 on dimension N (1, 2 por 3).
  • "FEM_ARGYRIS" : Argyris element on the triangle.
  • "FEM_HCT_TRIANGLE" : Hsieh-Clough-Tocher element on the triangle (composite P3 element which is C^1, 12 dof).
  • "FEM_REDUCED_HCT_TRIANGLE" : Hsieh-Clough-Tocher element on the triangle (composite P3 element which is C^1, 9 dof).
  • "FEM_QUADC1_COMPOSITE" : quadrilateral element, composite P3 element and C^1 (16 dof).
  • "FEM_REDUCED_QUADC1_COMPOSITE" : quadrilateral element, composite P3 element and C^1 (12 dof).
  • "FEM_PK_HIERARCHICAL(N,K)" : PK element with a hierarchical basis.
  • "FEM_QK_HIERARCHICAL(N,K)" : QK element with a hierarchical basis.
  • "FEM_PRISM_PK_HIERARCHICAL(N,K)" : PK element on a prism with a hierarchical basis.
  • "FEM_STRUCTURED_COMPOSITE(FEM, K)" : Composite fem on a grid with K divisions.
  • "FEM_PK_HIERARCHICAL_COMPOSITE(N,K,S)" : PK composite element on a grid with S subdivisions and with a hierarchical basis.
  • "FEM_PK_FULL_HIERARCHICAL_COMPOSITE(N,K,S)" : PK composite element with S subdivisions and a hierarchical basis on both degree and subdivision.
  • "FEM_PYRAMID_QK(K)" : Lagrange element on a 3D pyramid of degree K=0, 1 or 2. Can be connected to a standard P1/P2 Lagrange element on its triangular faces and a standard Q1/Q2 Lagrange element on its quadrangular face.
  • "FEM_PYRAMID_QK_DISCONTINUOUS(K)" : Discontinuous Lagrange element on a 3D pyramid of degree K = 0, 1 or 2.
  • "FEM_PYRAMID_Q2_INCOMPLETE" : Incomplete Lagrange element on a quadratic 3D pyramid (serendipity, 13-node element). Can be connected to a standard P2 Lagrange element on its triangular faces and a Q2_INCOMPLETE Lagrange element on its quadrangular face.
  • "HHO(fem_interior, fem_face_1, ..., fem_face_n)" : Build a hybrid method with "fem_interior" on the element itself and "fem_face_1", ..., "fem_face_n" on each face. If only one method is given for the faces, it is duplicated on each face.

Definition in file getfem_fem.h.