Go to the documentation of this file.
136 #ifndef GETFEM_FEM_H__
137 #define GETFEM_FEM_H__
160 struct dof_description;
247 typedef std::shared_ptr<const getfem::virtual_fem>
pfem;
250 typedef std::shared_ptr<const getfem::fem_precomp_> pfem_precomp;
256 public std::enable_shared_from_this<const virtual_fem> {
258 enum vec_type { VECTORIAL_NOTRANSFORM_TYPE, VECTORIAL_PRIMAL_TYPE,
259 VECTORIAL_DUAL_TYPE };
263 mutable std::vector<pdof_description> dof_types_;
269 std::shared_ptr<bgeot::convex_structure> cvs_node;
270 std::vector<std::vector<short_type>> face_tab;
272 mutable bgeot::pstored_point_tab pspt;
273 mutable bool pspt_valid;
274 bgeot::pconvex_ref cvr;
275 dim_type ntarget_dim;
276 mutable dim_type dim_;
281 bool real_element_defined;
282 bool is_standard_fem;
287 std::string debug_name_;
297 {
return dof_types_.size(); }
303 {
return nb_base(cv) * ntarget_dim; }
305 {
return nb_dof(cv) * ntarget_dim; }
308 {
return dof_types_; }
309 short_type hierarchical_raff()
const {
return hier_raff; }
311 dim_type
dim()
const {
return dim_; }
312 dim_type &
dim() {
return dim_; }
320 bgeot::pconvex_ref &mref_convex() {
return cvr; }
330 const std::string &debug_name()
const {
return debug_name_; }
331 std::string &debug_name() {
return debug_name_; }
332 virtual bgeot::pstored_point_tab node_tab(
size_type)
const {
334 pspt = bgeot::store_point_tab(cv_node.points());
346 {
return (*(node_tab(cv)))[i];}
347 virtual const std::vector<short_type> &
349 bool is_on_real_element()
const {
return real_element_defined; }
350 bool is_equivalent()
const {
return is_equiv; }
352 {
return !(is_equivalent()) || real_element_defined; }
357 bool is_polynomialcomp()
const {
return is_polycomp; }
358 bool is_standard()
const {
return is_standard_fem; }
359 bool &is_polynomialcomp() {
return is_polycomp; }
360 bool &is_equivalent() {
return is_equiv; }
363 bool &is_standard() {
return is_standard_fem; }
364 short_type estimated_degree()
const {
return es_degree; }
365 short_type &estimated_degree() {
return es_degree; }
367 virtual void mat_trans(base_matrix &,
const base_matrix &,
369 { GMM_ASSERT1(
false,
"This function should not be called."); }
385 template<
typename CVEC,
typename VVEC>
387 const CVEC& coeff, VVEC &val, dim_type Qdim)
const;
395 template <
typename MAT>
397 MAT &M, dim_type Qdim)
const;
402 template<
typename CVEC,
typename VMAT>
404 const CVEC& coeff, VMAT &val,
405 dim_type Qdim=1)
const;
410 template<
typename CVEC,
typename VMAT>
412 const CVEC& coeff, VMAT &val,
413 dim_type Qdim)
const;
418 template<
typename CVEC>
420 (
const fem_interpolation_context& c,
const CVEC& coeff,
421 typename gmm::linalg_traits<CVEC>::value_type &val)
const;
427 virtual void base_value(
const base_node &x, base_tensor &t)
const = 0;
433 virtual void grad_base_value(
const base_node &x, base_tensor &t)
const = 0;
439 virtual void hess_base_value(
const base_node &x, base_tensor &t)
const = 0;
447 base_tensor &t,
bool withM =
true)
const;
455 base_tensor &t,
bool withM =
true)
const;
463 base_tensor &t,
bool withM =
true)
const;
466 { GMM_ASSERT1(
false,
"internal error."); }
473 const dal::bit_vector &faces);
475 void init_cvs_node();
476 void unfreeze_cvs_node();
478 virtual_fem &operator =(
const virtual_fem &f) {
479 copy(f);
return *
this;
483 DAL_STORED_OBJECT_DEBUG_CREATED(
this,
"Fem");
484 ntarget_dim = 1; dim_ = 1;
485 is_equiv = is_pol = is_polycomp = is_lag = is_standard_fem =
false;
486 pspt_valid =
false; hier_raff = 0; real_element_defined =
false;
488 vtype = VECTORIAL_NOTRANSFORM_TYPE;
489 cvs_node = bgeot::new_convex_structure();
491 virtual_fem(
const virtual_fem& f)
492 :
dal::static_stored_object(),
493 std::enable_shared_from_this<const virtual_fem>()
494 { copy(f); DAL_STORED_OBJECT_DEBUG_CREATED(
this,
"Fem"); }
495 virtual ~virtual_fem() { DAL_STORED_OBJECT_DEBUG_DESTROYED(
this,
"Fem"); }
497 void copy(
const virtual_fem &f);
507 std::vector<FUNC> base_;
508 mutable std::vector<std::vector<FUNC>> grad_, hess_;
509 mutable bool grad_computed_ =
false;
510 mutable bool hess_computed_ =
false;
512 void compute_grad_()
const {
513 if (grad_computed_)
return;
515 if (grad_computed_)
return;
521 for (dim_type j = 0; j < n; ++j) {
522 grad_[i][j] = base_[i]; grad_[i][j].derivative(j);
525 grad_computed_ =
true;
528 void compute_hess_()
const {
529 if (hess_computed_)
return;
531 if (hess_computed_)
return;
536 hess_[i].resize(n*n);
537 for (dim_type j = 0; j < n; ++j) {
538 for (dim_type k = 0; k < n; ++k) {
539 hess_[i][j+k*n] = base_[i];
540 hess_[i][j+k*n].derivative(j); hess_[i][j+k*n].derivative(k);
544 hess_computed_ =
true;
550 const std::vector<FUNC> &
base()
const {
return base_; }
551 std::vector<FUNC> &
base() {
return base_; }
555 bgeot::multi_index mi(2);
559 base_tensor::iterator it = t.begin();
561 *it = bgeot::to_scalar(base_[i].eval(x.begin()));
567 if (!grad_computed_) compute_grad_();
568 bgeot::multi_index mi(3);
573 base_tensor::iterator it = t.begin();
574 for (dim_type j = 0; j < n; ++j)
576 *it = bgeot::to_scalar(grad_[i][j].eval(x.begin()));
582 if (!hess_computed_) compute_hess_();
583 bgeot::multi_index mi(4);
589 base_tensor::iterator it = t.begin();
590 for (dim_type k = 0; k < n; ++k)
591 for (dim_type j = 0; j < n; ++j)
593 *it = bgeot::to_scalar(hess_[i][j+k*n].eval(x.begin()));
616 bool complete=
false);
637 scalar_type alpha=0,
bool complete=
false);
657 const bgeot::pstored_point_tab pspt;
658 mutable std::vector<base_tensor> c;
659 mutable std::vector<base_tensor> pc;
660 mutable std::vector<base_tensor> hpc;
664 {
if (c.empty()) init_val();
return c[i]; }
667 {
if (pc.empty()) init_grad();
return pc[i]; }
670 {
if (hpc.empty()) init_hess();
return hpc[i]; }
671 inline pfem get_pfem()
const {
return pf; }
674 inline bgeot::pstored_point_tab get_ppoint_tab()
const
676 fem_precomp_(
const pfem,
const bgeot::pstored_point_tab);
677 ~fem_precomp_() { DAL_STORED_OBJECT_DEBUG_DESTROYED(
this,
"Fem_precomp"); }
679 void init_val()
const;
680 void init_grad()
const;
681 void init_hess()
const;
703 dal::pstatic_stored_object dep);
718 std::set<pfem_precomp> precomps;
752 mutable base_matrix M_;
764 const base_matrix&
M()
const;
768 void base_value(base_tensor& t,
bool withM =
true)
const;
770 void pfp_base_value(base_tensor& t,
const pfem_precomp &pfp__);
776 void pfp_grad_base_value(base_tensor& t,
const pfem_precomp &pfp__);
785 bool is_convex_num_valid()
const;
786 void invalid_convex_num() { convex_num_ =
size_type(-1); }
794 pfem_precomp
pfp()
const {
return pfp_; }
795 void set_pfp(pfem_precomp newpfp);
796 void set_pf(
pfem newpf);
797 int xfem_side()
const {
return xfem_side_; }
798 void set_xfem_side(
int side) { xfem_side_ = side; }
799 void change(bgeot::pgeotrans_precomp pgp__,
801 const base_matrix& G__,
size_type convex_num__,
803 bgeot::geotrans_interpolation_context::change(pgp__,ii__,G__);
804 convex_num_ = convex_num__; face_num_ = face_num__; xfem_side_ = 0;
809 const base_matrix& G__,
size_type convex_num__,
811 bgeot::geotrans_interpolation_context::change
812 (pgt__, pfp__->get_ppoint_tab(), ii__, G__);
813 convex_num_ = convex_num__; face_num_ = face_num__; xfem_side_ = 0;
817 pfem pf__,
const base_node& xref__,
const base_matrix& G__,
819 bgeot::geotrans_interpolation_context::change(pgt__,xref__,G__);
820 pf_ = pf__; pfp_ = 0; convex_num_ = convex_num__; face_num_ = face_num__;
823 fem_interpolation_context()
824 :
bgeot::geotrans_interpolation_context(),
826 fem_interpolation_context(bgeot::pgeotrans_precomp pgp__,
828 const base_matrix& G__,
831 :
bgeot::geotrans_interpolation_context(pgp__,ii__,G__),
832 convex_num_(convex_num__), face_num_(face_num__), xfem_side_(0)
836 const base_matrix& G__,
839 :
bgeot::geotrans_interpolation_context(pgt__,pfp__->get_ppoint_tab(),
841 convex_num_(convex_num__), face_num_(face_num__), xfem_side_(0)
845 const base_node& xref__,
846 const base_matrix& G__,
849 :
bgeot::geotrans_interpolation_context(pgt__,xref__,G__),
850 pf_(pf__), pfp_(0), convex_num_(convex_num__), face_num_(face_num__),
857 template <
typename CVEC,
typename VVEC>
859 const CVEC& coeff, VVEC &val,
860 dim_type Qdim)
const {
863 GMM_ASSERT1(gmm::vect_size(val) == Qdim,
"dimensions mismatch");
864 GMM_ASSERT1(gmm::vect_size(coeff) == nbdof*Qmult,
865 "Wrong size for coeff vector");
872 typename gmm::linalg_traits<CVEC>::value_type co = coeff[j*Qmult+q];
874 val[r + q*
target_dim()] += co * Z[j + r*nbdof];
879 template <
typename MAT>
881 MAT &M, dim_type Qdim)
const {
884 GMM_ASSERT1(gmm::mat_nrows(M) == Qdim && gmm::mat_ncols(M) == nbdof*Qmult,
885 "dimensions mismatch");
892 M(r+q*
target_dim(), j*Qmult+q) = Z[j + r*nbdof];
901 template<
typename CVEC,
typename VMAT>
903 const CVEC& coeff, VMAT &val,
904 dim_type Qdim)
const {
907 size_type Qmult = gmm::vect_size(coeff) / nbdof;
908 GMM_ASSERT1(gmm::mat_ncols(val) == N &&
910 gmm::vect_size(coeff) == nbdof*Qmult,
911 "dimensions mismatch");
913 "dimensions mismatch");
919 base_tensor::const_iterator it = t.begin();
922 for (
size_type j = 0; j < nbdof; ++j, ++it)
923 val(r + q*
target_dim(), k) += coeff[j*Qmult+q] * (*it);
928 template<
typename CVEC,
typename VMAT>
930 const CVEC& coeff, VMAT &val,
931 dim_type Qdim)
const {
934 GMM_ASSERT1(gmm::mat_ncols(val) == N*N &&
935 gmm::mat_nrows(val) == Qdim,
"dimensions mismatch");
943 base_tensor::const_iterator it = t.begin();
946 for (
size_type j = 0; j < nbdof; ++j, ++it)
947 val(r + q*
target_dim(), k) += coeff[j*Qmult+q] * (*it);
955 template<
typename CVEC>
958 typename gmm::linalg_traits<CVEC>::value_type &val)
const {
961 size_type Qmult = gmm::vect_size(coeff) / nbdof;
962 GMM_ASSERT1(gmm::vect_size(coeff) == nbdof*Qmult ,
"dimensions mismatch");
963 GMM_ASSERT1(target_dim()*Qmult == N &&
964 (Qmult == 1 || target_dim() == 1),
965 "Dimensions mismatch. Divergence operator requires fem qdim equal to dim.");
967 real_grad_base_value(c, t);
971 val = scalar_type(0);
972 base_tensor::const_iterator it = t.begin();
975 if (k) it += (N*nbdof + 1);
978 val += coeff[j] * (*it);
986 val += coeff[j*N+k] * (*it);
1001 typedef dal::naming_system<virtual_fem>::param_list fem_param_list;
1006 void add_fem_name(std::string name,
virtual void base_value(const base_node &x, base_tensor &t) const =0
Give the value of all components of the base functions at the point x of the reference element.
void interpolation_grad(const fem_interpolation_context &c, const CVEC &coeff, VMAT &val, dim_type Qdim=1) const
Interpolation of the gradient.
const base_matrix & M() const
non tau-equivalent transformation matrix.
dim_type dim() const
dimension of the reference element.
pfem interior_fem_of_hho_method(pfem hho_method)
Specific function for a HHO method to obtain the method in the interior.
const typedef fem< bgeot::polynomial_composite > * ppolycompfem
Polynomial composite FEM.
const base_tensor & val(size_type i) const
returns values of the base functions
virtual size_type nb_dof(size_type) const
Number of degrees of freedom.
dof_description * pdof_description
Type representing a pointer on a dof_description.
const std::vector< pdof_description > & dof_types() const
Get the array of pointer on dof description.
size_t size_type
used as the common size type in the library
bool is_lagrange() const
true if the base functions are such that
void hess_base_value(const base_node &x, base_tensor &t) const
Evaluates at point x, the hessian of all base functions w.r.t.
bool is_on_face() const
On a face ?
virtual bgeot::pconvex_ref ref_convex(size_type) const
Return the convex of the reference element.
virtual const bgeot::convex< base_node > & node_convex(size_type) const
Gives the convex representing the nodes on the reference element.
virtual void real_hess_base_value(const fem_interpolation_context &c, base_tensor &t, bool withM=true) const
Give the hessian of all components of the base functions at the current point of the fem_interpolatio...
pfem_precomp fem_precomp(pfem pf, bgeot::pstored_point_tab pspt, dal::pstatic_stored_object dep)
Handles precomputations for FEM.
void hess_base_value(base_tensor &t, bool withM=true) const
fill the tensor with the hessian of the base functions (taken at point this->xref())
Integration methods (exact and approximated) on convexes.
Associate a name to a method descriptor and store method descriptors.
void interpolation(const fem_interpolation_context &c, const CVEC &coeff, VVEC &val, dim_type Qdim) const
Interpolate at an arbitrary point x given on the reference element.
void interpolation_hess(const fem_interpolation_context &c, const CVEC &coeff, VMAT &val, dim_type Qdim) const
Interpolation of the hessian.
the geotrans_interpolation_context structure is passed as the argument of geometric transformation in...
void del_stored_object(const pstatic_stored_object &o, bool ignore_unstored)
Delete an object and the object which depend on it.
pfem 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 conve...
const base_tensor & grad(size_type i) const
returns gradients of the base functions
Base class for finite element description.
size_type convex_num() const
get the current convex number
base_poly read_base_poly(short_type n, std::istream &f)
read a base_poly on the stream ist.
const typedef fem< bgeot::base_poly > * ppolyfem
Classical polynomial FEM.
Pre-computations on a fem (given a fixed set of points on the reference convex, this object computes ...
gmm::uint16_type short_type
used as the common short type integer in the library
bgeot::pconvex_structure structure(size_type cv) const
Gives the convex structure of the reference element nodes.
bool dof_compatibility(pdof_description, pdof_description)
Says if the two dofs can be identified.
GEneric Tool for Finite Element Methods.
pdof_description product_dof(pdof_description pnd1, pdof_description pnd2)
Product description of the descriptions *pnd1 and *pnd2.
bgeot::pconvex_structure basic_structure(size_type cv) const
Gives the convex of the reference element.
structure passed as the argument of fem interpolation functions.
void add_node(const pdof_description &d, const base_node &pt, const dal::bit_vector &faces)
internal function adding a node to an element for the creation of a finite element method.
short_type face_num() const
get the current face number
virtual_fem implementation as a vector of generic functions.
pdof_description second_derivative_dof(dim_type d, dim_type num_der1, dim_type num_der2)
Description of a unique dof of second derivative type.
dim_type target_dim() const
dimension of the target space.
This class deals with plain polynomials with several variables.
vec_type vectorial_type() const
Type of vectorial element.
pdof_description derivative_dof(dim_type d, dim_type r)
Description of a unique dof of derivative type.
virtual void hess_base_value(const base_node &x, base_tensor &t) const =0
Give the value of all hessians (on ref.
void delete_fem_precomp(pfem_precomp pfp)
Request for the removal of a pfem_precomp.
virtual size_type nb_base(size_type cv) const
Number of basis functions.
Stores interdependent getfem objects.
bool dof_linkable(pdof_description)
Says if the dof is linkable.
std::shared_ptr< const getfem::virtual_fem > pfem
type of pointer on a fem description
bool have_pf() const
true if the pfem is available.
bool is_polynomial() const
true if the base functions are polynomials
pdof_description mean_value_dof(dim_type d)
Description of a unique dof of mean value type.
Geometric transformations on convexes.
const base_node & node_of_dof(size_type cv, size_type i) const
Gives the node corresponding to the dof i.
pdof_description normal_derivative_dof(dim_type d)
Description of a unique dof of normal derivative type (normal derivative at the node,...
pfem_precomp operator()(pfem pf, bgeot::pstored_point_tab pspt)
Request a pfem_precomp.
size_type nb_base_components(size_type cv) const
Number of components (nb_dof() * dimension of the target space).
bool have_pfp() const
true if a fem_precomp_ has been supplied.
pdof_description lagrange_dof(dim_type d)
Description of a unique dof of lagrange type (value at the node).
const typedef fem< bgeot::base_rational_fraction > * prationalfracfem
Rational fration FEM.
virtual void grad_base_value(const base_node &x, base_tensor &t) const =0
Give the value of all gradients (on ref.
std::shared_ptr< const convex_structure > pconvex_structure
Pointer on a convex structure description.
pfem_precomp pfp() const
get the current fem_precomp_
base class for static stored objects
void grad_base_value(const base_node &x, base_tensor &t) const
Evaluates at point x, the gradient of all base functions w.r.t.
size_type dof_xfem_index(pdof_description)
Returns the xfem_index of dof (0 for normal dof)
std::shared_ptr< const bgeot::geometric_trans > pgeometric_trans
pointer type for a geometric transformation
std::string name_of_fem(pfem p)
get the string name of a fem descriptor.
void set_face_num(short_type f)
set the current face number
const pfem pf() const
get the current FEM descriptor
void base_value(const base_node &x, base_tensor &t) const
Evaluates at point x, all base functions and returns the result in t(nb_base,target_dim)
void grad_base_value(base_tensor &t, bool withM=true) const
fill the tensor with the gradient of the base functions (taken at point this->xref())
pdof_description xfem_dof(pdof_description p, size_type ind)
Description of a special dof for Xfem.
void interpolation_diverg(const fem_interpolation_context &c, const CVEC &coeff, typename gmm::linalg_traits< CVEC >::value_type &val) const
Interpolation of the divergence.
const std::vector< FUNC > & base() const
Gives the array of basic functions (components).
void base_value(base_tensor &t, bool withM=true) const
fill the tensor with the values of the base functions (taken at point this->xref())
Handle composite polynomials.
pfem 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...
virtual void real_grad_base_value(const fem_interpolation_context &c, base_tensor &t, bool withM=true) const
Give the gradient of all components of the base functions at the current point of the fem_interpolati...
const base_tensor & hess(size_type i) const
returns hessians of the base functions
pfem fem_descriptor(const std::string &name)
get a fem descriptor from its string name.
int dof_description_compare(pdof_description a, pdof_description b)
Gives a total order on the dof description compatible with the identification.
pdof_description global_dof(dim_type d)
Description of a global dof, i.e.
virtual void real_base_value(const fem_interpolation_context &c, base_tensor &t, bool withM=true) const
Give the value of all components of the base functions at the current point of the fem_interpolation_...