Provided by: librheolef-dev_7.2-2_amd64 bug

NAME

       field - piecewise polynomial finite element function (rheolef-7.2)

DESCRIPTION

       This class represents a piecewise polynomial finite element function. Since this function
       spans onto a finite dimensional basis, it simply stores its coefficients on this basis:
       these coefficients are called the degrees of freedom (see space(2)).

INTERPOLATION

       For any function u, its interpolation on the finite element space Xh as a field uh in Xh
       expresses simply via the interpolate(3) function:

           Float u (const point& x) { return exp(x[0]*x[1]); }
           ...
           field uh = interpolate (Xh, u);

LINEAR ALGEBRA

       Linear algebra, such as uh+vh, uh-vh and lambda*uh + mu*vh, where lambda and mu are of
       type Float, are supported. The duality product between two fields lh and uh writes simply
       dual(lh,uh). As we consider finite dimensional spaces, this duality product coincides with
       the usual Euclidean dot product in IR^dim(Xh). The application of the a bilinear form(2)
       writes a(uh,vh) and is equivalent to dual(m*uh,vh).

COMMON ACCESSORS

       For convenience, uh.max(), uh.min() and uh.max_abs() returns respectively the maximum,
       minimum and maximum of the absolute value of the degrees of freedom.

NON-LINEAR ALGEBRA

       Non-linear operations, such as sqrt(uh) or 1/uh are also available. Note that non-linear
       operations do not returns in general piecewise polynomials: the value returned by sqrt(uh)
       may be filtered by interpolate(3) function:

           field vh = interpolate (Xh, sqrt(uh));

        Also, the multiplication uh*vh and the division uh/vh returns a result that is not in the
       same discrete finite element space: its result also may be filtered by interpolate(3)`:

           field wh = interpolate(Xh, uh*vh);

        All standard unary and binary math functions abs, cos, sin... are extended to scalar
       fields. Also sqr(uh), the square of a field, and min(uh,vh), max(uh,vh) are provided.
       Binary functions can be used also with a scalar, as in

           field vh = interpolate (Xh, max (abs(uh), 0));
           field wh = interpolate (Xh, pow (abs(uh), 1./3));

        For applying a user-provided function to a field, please see the compose(3) function.

ACCESS BY DOMAIN

       The restriction of a field to a geometric domain, says 'boundary' writes uh['boundary']:
       it represents the trace of the field on the boundary:

           space Xh (omega, "P1");
           uh["boundary"] = 0;

        A possible alternative uses a geo(2) domain as index:

           geo boundary = omega["boundary"];
           uh[boundary] = 0;

MULTI-VALUED FIELDS

       A vector-valued field contains several components, as:

           space Xh (omega, "P2", "vector");
           field uh (Xh);

        Conversely, for a tensor-valued field:

           space Th (omega, "P1d", "tensor");
           field sigma_h (Xh);

GENERAL SPACE PRODUCT INTERFACE

       A general multi-component field writes:

           space Th (omega, "P1d", "tensor");
           space Vh (omega, "P2", "vector");
           space Qh (omega, "P1");
           space Xh = Th*Vh*Qh;
           field xh (Xh);
           field tau_h = xh[0]; // tensor-valued
           field uh    = xh[1]; // vector-valued
           field qh    = xh[2]; // scalar

        Remark the hierarchical multi-component field structure: the first-component is tensor-
       valued and the second-one is vector-valued. There is no limitation upon the hierarchical
       number of levels in use: the hierarchy is not flattened.

       The xh.size() returns the number of field components. When the field is scalar, it returns
       zero by convention, and xh[0] is undefined.

DIRECT ACCESS TO DEGREES-OF-FREEDOM

       The field class provides a STL-like container interface for accessing the degrees-of-
       freedom (dofs) of a finite element field uh. The number of dofs is uh.ndof() and any dof
       can be accessed via uh.dof(idof). In a distributed memory environment, a non-local dof at
       the partition interface can be obtain via uh.dis_dof(dis_idof) where dis_idof is the
       (global) distributed index assoiated to the distribution uh.ownership(). See
       distributor(4).

       For better performances, a STL-like iterator interface is available, with uh.begin_dof()
       and uh.end_dof() returns iterators to the dofs on the current processor.

REPRESENTATION

       The degrees of freedom (see space(2)) are splited between unknowns and blocked, i.e.
       uh=[uh.u,uh.b] for any field uh in Xh. Access to these vectors is allowed via some
       accessors: a read-only one, as uh.u() and uh.b(), and a read-and-write one, as uh.set_u()
       and uh.set_b(), see vec(4).

       Note that blocked and unknown degrees of freedom could also be elegantly set by using a
       domain name indexation (see geo(2)):

           geo omega ("circle");
           space Xh (omega, "P1");
           Xh.block ("boundary");
           field uh (Xh);
           uh ["boundary"] = 0;

IMPLEMENTATION

       This documentation has been generated from file main/lib/field.h

       The field class is simply an alias to the field_basic class

       typedef field_basic<Float> field;

       The field_basic class provides an interface to a vector data container:

       template <class T, class M = rheo_default_memory_model>
       class field_basic: public details::field_wdof_base<field_basic<T,M>> {
       public :
       // typedefs:

         using wdof_base     = details::field_wdof_base<field_basic<T,M>>;
         using rdof_base     = details::field_rdof_base<field_basic<T,M>>;
         using size_type     = std::size_t;
         using scalar_type   = T;
         using memory_type   = M;
         using float_type    = typename float_traits<T>::type;
         using geo_type      = geo_basic  <float_type,memory_type>;
         using space_type    = space_basic<float_type,memory_type>;
         using dis_reference = typename vec<scalar_type,memory_type>::dis_reference;
         using valued_type   = space_constant::valued_type;
         using value_type    = T; // TODO: run-time dependent, set it to undeterminated<T>
         using result_type   = T; // TODO: run-time dependent, set it to undeterminated<T>
         class iterator;
         class const_iterator;

       // allocator/deallocator:

         field_basic();

         explicit field_basic (
           const space_type& V,
           const T& init_value = std::numeric_limits<T>::max());

         void resize (
           const space_type& V,
           const T& init_value = std::numeric_limits<T>::max());

         // expressions: field_expr or field_lazy are accepted here
         // TODO: merge with FieldLazy piecewise_poly
         template <class Expr, class Sfinae =
           typename std::enable_if<
             (      details::is_field_expr_affine_homogeneous<Expr>::value
               && ! details::is_field<Expr>::value
             )
             ||   details::is_field_lazy<Expr>::value
           >::type>
         field_basic (const Expr& expr);

         field_basic (const std::initializer_list<details::field_concat_value<T,M> >& init_list);

       // assignments:

         field_basic<T,M>& operator= (const field_basic<T,M>&);
         field_basic<T,M>& operator= (const std::initializer_list<details::field_concat_value<T,M> >& init_list);

       // accessors:

         const space_type&  get_space()  const { return _V; }
         const geo_type&    get_geo()    const { return _V.get_geo(); }
         std::string        get_approx() const { return _V.get_approx(); }
         valued_type        valued_tag() const { return _V.valued_tag(); }
         const std::string& valued()     const { return _V.valued(); }
       #ifdef TO_CLEAN
         std::string        name()      const { return _V.name(); }
         bool have_homogeneous_space (space_basic<T,M>& Xh) const { Xh = get_space(); return true; }
       #endif // TO_CLEAN

       // accessors & modifiers to unknown & blocked parts:

         const vec<T,M>&     u() const { return _u; }
         const vec<T,M>&     b() const { return _b; }
               vec<T,M>& set_u()       { dis_dof_indexes_requires_update(); return _u; }
               vec<T,M>& set_b()       { dis_dof_indexes_requires_update(); return _b; }

       // accessors to extremas:

         T min() const;
         T max() const;
         T max_abs() const;
         T min_abs() const;

       // accessors by degrees-of-freedom (dof):

         const distributor& ownership() const { return get_space().ownership(); }
         const communicator& comm() const { return ownership().comm(); }
         size_type     ndof() const { return ownership().size(); }
         size_type dis_ndof() const { return ownership().dis_size(); }
               T& dof (size_type idof);
         const T& dof (size_type idof) const;
         const T& dis_dof (size_type dis_idof) const;
         // write access to non-local dis_idof changes to others procs
         dis_reference dis_dof_entry (size_type dis_idof);

         iterator begin_dof();
         iterator end_dof();
         const_iterator begin_dof() const;
         const_iterator end_dof() const;

       // input/output:

         idiststream& get (idiststream& ips);
         odiststream& put (odiststream& ops) const;
         odiststream& put_field (odiststream& ops) const;

       // evaluate uh(x) where x is given locally as hat_x in K:

         T dis_evaluate (const point_basic<T>& x, size_type i_comp = 0) const;
         T operator()   (const point_basic<T>& x) const { return dis_evaluate (x,0); }
         point_basic<T> dis_vector_evaluate (const point_basic<T>& x) const;

       };

       namespace details {

       // concepts:
       template<class T, class M>
       struct is_field<field_basic<T,M>> : std::true_type {};

       // field class is not an ordinary function/functor : for compose(field,args...) filtering
       template <typename T, typename M>  struct function_traits <field_basic<T,M> > {};

       template<class T, class M>
       struct field_traits<field_basic<T,M>> {
         using size_type     = std::size_t;
         using scalar_type   = T;
         using memory_type   = M;
       };

       } // namespace details

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.