Welcome to the NetCologne GmbH open source mirroring service!

This machine mirrors various open-source projects. 20 Gbit/s uplink.

If there are any issues or you want another project mirrored, please contact mirror-service -=AT=- netcologne DOT de !

GetFEM: src/getfem/getfem_integration.h File Reference
GetFEM  5.4.2
getfem_integration.h File Reference

Integration methods (exact and approximated) on convexes. More...

Go to the source code of this file.

Classes

class  getfem::poly_integration
 Description of an exact integration of polynomials. More...
 
class  getfem::integration_method
 this structure is not intended to be used directly. More...
 

Namespaces

 getfem
 GEneric Tool for Finite Element Methods.
 

Enumerations

enum  getfem::integration_method_type
 the list of main integration method types
 

Functions

pintegration_method getfem::int_method_descriptor (std::string name, bool throw_if_not_found=true)
 Get an integration method from its name . More...
 
std::string getfem::name_of_int_method (pintegration_method p)
 Get the string name of an integration method .
 
pintegration_method getfem::classical_exact_im (bgeot::pgeometric_trans pgt)
 return an exact integration method for convex type handled by pgt. More...
 
pintegration_method getfem::classical_approx_im (bgeot::pgeometric_trans pgt, dim_type degree)
 try to find an approximate integration method for the geometric transformation pgt which is able to integrate exactly polynomials of degree <= "degree". More...
 
pintegration_method getfem::exact_simplex_im (size_type n)
 return IM_EXACT_SIMPLEX(n)
 
pintegration_method getfem::exact_parallelepiped_im (size_type n)
 return IM_EXACT_PARALLELEPIPED(n)
 
pintegration_method getfem::exact_prism_im (size_type n)
 return IM_EXACT_PRISM(n)
 
pintegration_method getfem::exact_classical_im (bgeot::pgeometric_trans pgt) IS_DEPRECATED
 use classical_exact_im instead.
 
pintegration_method getfem::im_none (void)
 return IM_NONE
 

Detailed Description

Integration methods (exact and approximated) on convexes.

Author
Yves Renard Yves..nosp@m.Rena.nosp@m.rd@in.nosp@m.sa-l.nosp@m.yon.f.nosp@m.r
Date
December 17, 2000.

List of integration methods for getfem::int_method_descriptor

  • "IM_EXACT_SIMPLEX(n)" : exact integration on simplexes.
  • "IM_PRODUCT(IM1, IM2)" : product of two integration methods
  • "IM_EXACT_PARALLELEPIPED(n)" : exact integration on parallelepipeds
  • "IM_EXACT_PRISM(n)" : exact integration on prisms
  • "IM_GAUSS1D(K)" : Gauss method on the segment, order K=1, 3, 5, 7, ..., 99.
  • "IM_GAUSSLOBATTO1D(K)" : Gauss-Lobatto method on the segment, order K=1, 3, 5, 7, ..., 99.
  • "IM_NC(N,K)" : Newton-Cotes approximative integration on simplexes, order K
  • "IM_NC_PARALLELEPIPED(N,K)" : product of Newton-Cotes integration on parallelepipeds
  • "IM_NC_PRISM(N,K)" : product of Newton-Cotes integration on prisms
  • "IM_GAUSS_PARALLELEPIPED(N,K)" : product of Gauss1D integration on parallelepipeds (K=1,3,..99)
  • "IM_TRIANGLE(K)" : Gauss methods on triangles (K=1..10,13,17,19)
  • "IM_QUAD(K)" : Gauss methods on quadrilaterons (K=2, 3, 5, 7, 9, 17)
  • "IM_HEXAHEDRON(K)" : Gauss methods on hexahedrons (K=5,9,11)
  • "IM_TETRAHEDRON(K)" : Gauss methods on tetrahedrons (K=1, 2, 3, 5, 6, or 8)
  • "IM_SIMPLEX4D(3)" : Gauss method on a 4-dimensional simplex.
  • "IM_CUBE4D(K)" : Gauss method on a 4-dimensional cube (K=5,9).
  • "IM_STRUCTURED_COMPOSITE(IM1, K)" : Composite method on a grid with K divisions
  • "IM_HCT_COMPOSITE(IM1)" : Composite integration suited to the HCT composite finite element.
  • "IM_QUAQC1_COMPOSITE(IM1)" : Composite integration suited to the QUADC1 composite finite element.
  • "IM_QUASI_POLAR(IM1, IP1)" : if IM1 is an integration method on a square, gives an integration method on a triangle which is close to a polar integration with respect to vertex IP1. if IM1 is an integration method on a tetrahedron, gives an integration method on a tetrahedron which is close to a cylindrical integration with respect to vertex IP1 (does not work very well). if IM1 is an integration method on a prism. Gives an integration method on a tetrahedron which is close to a cylindrical integration with respect to vertex IP1.
  • "IM_QUASI_POLAR(IM1, IP1, IP2)" : IM1 should be an integration method on a prism. Gives an integration method on a tetrahedron which is close to a cylindrical integration with respect to IP1-IP2 axis.
  • "IM_PYRAMID_COMPOSITE(IM1)" : Composite integration for a pyramid decomposed into two tetrahedrons

Definition in file getfem_integration.h.