Provided by: librheolef-dev_7.2-3build5_amd64 bug

NAME

       basis - finite element method (rheolef-7.2)

DESCRIPTION

       The class constructor takes a string that characterizes the finite element method, e.g.
       'P0', 'P1' or 'P2'. The letters characterize the finite element family and the number is
       the polynomial degree in this family. The finite element also depends upon the reference
       element, e.g. triangle, square, tetrahedron, etc. See reference_element(6) for more
       details. For instance, on the square reference element, the 'P1' string designates the
       common Q1 four-nodes finite element.

AVAILABLE FINITE ELEMENTS

       Pk

            The classical Lagrange finite element family, where k is any polynomial degree
           greater or equal to 1.

       Pkd

            The discontinuous Galerkin finite element family, where k is any polynomial degree
           greater or equal to 0. For convenience, P0d could also be written as P0.

       bubble

            The product of baycentric coordinates. Only simplicial elements are supported:
           edge(6), triangle(6) and tetrahedron(6).

       RTk
       RTkd

            The vector-valued Raviart-Thomas family, where k is any polynomial degree greater or
           equal to 0. The RTkd piecewise discontinuous version is fully implemented for the
           triangle(6). Its RTk continuous counterpart is still under development.

       Bk

            The Bernstein finite element family, where k is any polynomial degree greater or
           equal to 1. This basis was initially introduced by Bernstein (Comm. Soc. Math.
           Kharkov, 2th series, 1912) and more recently used in the context of finite elements.
           This feature is still under development.

       Sk

            The spectral finite element family, as proposed by Sherwin and Karniadakis (Oxford
           University Press, 2nd ed., 2005). Here, k is any polynomial degree greater or equal to
           1. This feature is still under development.

CONTINUOUS FEATURE

       Some finite element families could support either continuous or discontinuous junction at
       inter-element boundaries. For instance, the Pk Lagrange finite element family, with k >=
       1, has a continuous interelements junction: it defines a piecewise polynomial and globally
       continuous finite element functional space(2). Conversely, the Pkd Lagrange finite element
       family, with k >= 0, has a discontinuous interelements junction: it defines a piecewise
       polynomial and globally discontinuous finite element functional space(2).

OPTIONS

       The basis class supports some options associated to each finite element method. These
       options are appended to the string bewteen square braces, and are separated by a coma,
       e.g. 'P5[feteke,bernstein]'. See basis(1) for some examples.

LAGRANGE NODES OPTION

       equispaced

            Nodes are equispaced. This is the default.

       fekete

            Node are non-equispaced: for high-order polynomial degree, e.g. greater or equal to
           3, their placement is fully optimized, as proposed by Taylor, Wingate and Vincent,
           2000, SIAM J. Numer. Anal. With this choice, the interpolation error is dramatically
           decreased for high order polynomials. This feature is still restricted to the triangle
           reference element and to polynomial degree lower or equal to 19. Otherwise, it fall
           back to the following warburton node set.

       warburton

            Node are non-equispaced: for high-order polynomial degree, e.g. greater or equal to
           3, their placement is optimized thought an heuristic solution, as proposed by
           Warburton, 2006, J. Eng. Math. With this choice, the interpolation error is
           dramatically decreased for high order polynomials. This feature applies to any
           reference element and polynomial degree.

RAW POLYNOMIAL BASIS

       The raw (or initial) polynomial basis is used for building the dual basis, associated to
       degrees of freedom, via the Vandermonde matrix and its inverse. Changing the raw basis do
       not change the definition of the FEM basis, but only the way it is constructed. There are
       three possible raw basis:

       monomial

            The monomial basis is the simplest choice. While it is suitable for low polynomial
           degree (less than five), for higher polynomial degree, the Vandermonde matrix becomes
           ill-conditioned and its inverse leads to errors in double precision.

       dubiner

            The Dubiner basis (see Dubiner, 1991 J. Sci. Comput.) leads to much better condition
           number. This is the default.

       bernstein

            The Bernstein basis could also be an alternative raw basis.

INTERNALS: EVALUATE AND INTERPOLATE

       The basis class defines member functions that evaluates all the polynomial basis functions
       of a finite element and their derivatives at a point or a set of point.

       The basis polynomial functions are evaluated by the evaluate member function. This member
       function is called during the assembly process of the integrate(3) function.

       The interpolation nodes used by the interpolate(3) function are available by the member
       function hat_node. Conversely, the member function compute_dofs compute the degrees of
       freedom.

IMPLEMENTATION

       This documentation has been generated from file fem/lib/basis.h

       The basis class is simply an alias to the basis_basic class

       typedef basis_basic<Float> basis;

       The basis_basic class provides an interface to a data container:

       template<class T>
       class basis_basic : public smart_pointer_nocopy<basis_rep<T> >,
                           public persistent_table<basis_basic<T> > {
       public:

       // typedefs:

         typedef basis_rep<T>              rep;
         typedef smart_pointer_nocopy<rep> base;
         typedef typename rep::size_type   size_type;
         typedef typename rep::value_type  value_type;
         typedef typename rep::valued_type valued_type;

       // allocators:

         basis_basic (std::string name = "");
         void reset (std::string& name);
         void reset_family_index (size_type k);

       // accessors:

         bool is_initialized() const { return base::operator->() != 0; }
         size_type   degree() const;
         size_type   family_index() const;
         std::string family_name() const;
         std::string name() const;
         size_type   ndof (reference_element hat_K) const;
         size_type   nnod (reference_element hat_K) const;
         bool is_continuous() const;
         bool is_discontinuous() const;
         bool is_nodal() const;
         bool have_continuous_feature() const;
         bool have_compact_support_inside_element() const;
         bool is_hierarchical() const;
         size_type size() const;
         const basis_basic<T>& operator[] (size_type i_comp) const;
         bool have_index_parameter() const;
         const basis_option& option() const;
         valued_type        valued_tag() const;
         const std::string& valued()     const;
         const piola_fem<T>& get_piola_fem() const;

       };

AUTHOR

       Pierre  Saramito  <Pierre.Saramito@imag.fr>

COPYRIGHT

       Copyright   (C)  2000-2018  Pierre  Saramito  <Pierre.Saramito@imag.fr> GPLv3+: GNU GPL
       version 3 or later  <http://gnu.org/licenses/gpl.html>.  This  is  free  software:  you
       are free to change and redistribute it.  There is NO WARRANTY, to the extent permitted by
       law.