GetFEM  5.5
getfem_fem.h
Go to the documentation of this file.
1 /* -*- c++ -*- (enables emacs c++ mode) */
2 /*===========================================================================
3 
4  Copyright (C) 1999-2026 Yves Renard
5 
6  This file is a part of GetFEM
7 
8  GetFEM is free software; you can redistribute it and/or modify it
9  under the terms of the GNU Lesser General Public License as published
10  by the Free Software Foundation; either version 3 of the License, or
11  (at your option) any later version along with the GCC Runtime Library
12  Exception either version 3.1 or (at your option) any later version.
13  This program is distributed in the hope that it will be useful, but
14  WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
15  or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
16  License and GCC Runtime Library Exception for more details.
17  You should have received a copy of the GNU Lesser General Public License
18  along with this program. If not, see https://www.gnu.org/licenses/.
19 
20  As a special exception, you may use this file as it is a part of a free
21  software library without restriction. Specifically, if other files
22  instantiate templates or use macros or inline functions from this file,
23  or you compile this file and link it with other files to produce an
24  executable, this file does not by itself cause the resulting executable
25  to be covered by the GNU Lesser General Public License. This exception
26  does not however invalidate any other reasons why the executable file
27  might be covered by the GNU Lesser General Public License.
28 
29 ===========================================================================*/
30 
31 /**@file getfem_fem.h
32  @author Yves Renard <Yves.Renard@insa-lyon.fr>
33  @date December 21, 1999.
34  @brief Definition of the finite element methods.
35 
36  This file defines the getfem::virtual_fem class, which is the common
37  base class of all FEM.
38 
39  @section fem_list List of FEM known by getfem::fem_descriptor :
40 
41  - "FEM_PK(N,K)" : classical Lagrange element PK on a simplex.
42 
43  - "FEM_PK_DISCONTINUOUS(N,K,alpha)" : discontinuous Lagrange
44  element PK on a simplex.
45 
46  - "FEM_QK(N,K)" : classical Lagrange element QK on a parellepiped.
47 
48  - "FEM_QK_DISCONTINUOUS(N,K,alpha)" : discontinuous Lagrange
49  element QK on a parallelepiped.
50 
51  - "FEM_Q2_INCOMPLETE(N)" : incomplete Q2 elements with 8 and 20 dof
52  (serendipity Quad 8 and Hexa 20 elements)
53 
54  - "FEM_Q2_INCOMPLETE_DISCONTINUOUS(N)" : discontinuous incomplete Q2
55  elements with 8 and 20 dof (serendipity Quad 8 and Hexa 20 elements)
56 
57  - "FEM_PRISM_PK(N,K)" : classical Lagrange element PK on a prism.
58 
59  - "FEM_PRISM_PK_DISCONTINUOUS(N,K,alpha)" : classical discontinuous
60  Lagrange element PK on a prism.
61 
62  - "FEM_PRISM_INCOMPLETE_P2" : Incomplete Lagrange element on a
63  quadratic 3D prism (serendipity, 15-node wedge element). Can be connected
64  toa standard P2 Lagrange on its triangular faces and a Q2_INCOMPLETE
65  Lagrange element on its quadrangular faces.
66 
67  - "FEM_PK_WITH_CUBIC_BUBBLE(N,K)" : classical Lagrange element PK
68  on a simplex with an additional volumic bubble function.
69 
70  - "FEM_PRODUCT(FEM1,FEM2)" : tensorial product of two polynomial
71  elements
72 
73  - "FEM_P1_NONCONFORMING" : Nonconforming P1 method on a
74  triangle.
75 
76  - "FEM_P1_BUBBLE_FACE(N)" : P1 method on a simplex with an
77  additional bubble function on face 0.
78 
79  - "FEM_P1_BUBBLE_FACE_LAG" : P1 method on a simplex with an
80  additional lagrange dof on face 0.
81 
82  - "FEM_HERMITE(N)" : Hermite element P3 on dimension N (1, 2 por 3).
83 
84  - "FEM_ARGYRIS" : Argyris element on the triangle.
85 
86  - "FEM_HCT_TRIANGLE" : Hsieh-Clough-Tocher element on the triangle
87  (composite P3 element which is C^1, 12 dof).
88 
89  - "FEM_REDUCED_HCT_TRIANGLE" : Hsieh-Clough-Tocher element on the triangle
90  (composite P3 element which is C^1, 9 dof).
91 
92  - "FEM_QUADC1_COMPOSITE" : quadrilateral element, composite P3 element
93  and C^1 (16 dof).
94 
95  - "FEM_REDUCED_QUADC1_COMPOSITE" : quadrilateral element, composite
96  P3 element and C^1 (12 dof).
97 
98  - "FEM_PK_HIERARCHICAL(N,K)" : PK element with a hierarchical basis.
99 
100  - "FEM_QK_HIERARCHICAL(N,K)" : QK element with a hierarchical basis.
101 
102  - "FEM_PRISM_PK_HIERARCHICAL(N,K)" : PK element on a prism with a
103  hierarchical basis.
104 
105  - "FEM_STRUCTURED_COMPOSITE(FEM, K)" : Composite fem on a grid with
106  K divisions.
107 
108  - "FEM_PK_HIERARCHICAL_COMPOSITE(N,K,S)" : PK composite element on
109  a grid with S subdivisions and with a hierarchical basis.
110 
111  - "FEM_PK_FULL_HIERARCHICAL_COMPOSITE(N,K,S)" : PK composite
112  element with S subdivisions and a hierarchical basis on both degree
113  and subdivision.
114 
115  - "FEM_PYRAMID_QK(K)" : Lagrange element on a 3D pyramid of degree
116  K=0, 1 or 2. Can be connected to a standard P1/P2 Lagrange element on its
117  triangular faces and a standard Q1/Q2 Lagrange element on its quadrangular
118  face.
119 
120  - "FEM_PYRAMID_QK_DISCONTINUOUS(K)" : Discontinuous Lagrange element
121  on a 3D pyramid of degree K = 0, 1 or 2.
122 
123  - "FEM_PYRAMID_Q2_INCOMPLETE" : Incomplete Lagrange element on a
124  quadratic 3D pyramid (serendipity, 13-node element). Can be connected to
125  a standard P2 Lagrange element on its triangular faces and a Q2_INCOMPLETE
126  Lagrange element on its quadrangular face.
127 
128  - "HHO(fem_interior, fem_face_1, ..., fem_face_n)" : Build a hybrid method
129  with "fem_interior" on the element itself and "fem_face_1", ...,
130  "fem_face_n" on each face. If only one method is given for the faces, it
131  is duplicated on each face.
132 
133 */
134 
135 #ifndef GETFEM_FEM_H__
136 #define GETFEM_FEM_H__
137 
139 #include "bgeot_geometric_trans.h"
140 #include "bgeot_poly_composite.h"
141 #include "getfem_integration.h"
142 #include "dal_naming_system.h"
143 #include <deque>
144 
145 namespace getfem {
146 
147  /** @defgroup dofdescr Dof description
148  * This set of functions gives a pointer to dof descriptions and
149  * tests the compatibility between two dof descriptions.
150  * A dof description describes a type of dof (Lagrange type,
151  * Hermite type ...) in order to be able to say wether two dofs are
152  * to be identified or not. The construction of dof type with
153  * the tensorial product is taken into
154  * account, making the dof_description structure a little bit complex
155  * (this structure will probably evoluate in the future).
156  * @{
157  */
158 
159  struct dof_description;
160 
161  /// Type representing a pointer on a dof_description
162  typedef dof_description *pdof_description;
163 
164  /** @brief Description of a unique dof of lagrange type (value at the node).
165  * @param d the dimension of the reference element (2 for triangles, 3 for tetrahedrons ...)
166  */
167  pdof_description lagrange_dof(dim_type d);
168 
169  /** Description of a unique dof of derivative type.
170  * (a derivative at the node).
171  * @param d the dimension of the reference element.
172  * @param r corresponds to the variable number for which the derivative is taken (0 <= r < d)
173  */
174  pdof_description derivative_dof(dim_type d, dim_type r);
175 
176  /** Description of a unique dof of second derivative type.
177  * @param d the dimension of the reference element.
178  * @param num_der1 corresponds to the variable number for which the first derivative is taken (0 <= r < d)
179  * @param num_der2 corresponds to the variable number for which the second derivative is taken (0 <= r < d)
180  */
182  dim_type num_der1,
183  dim_type num_der2);
184 
185  /** Description of a unique dof of normal derivative type
186  * (normal derivative at the node, regarding a face).
187  * @param d the dimension of the reference element.
188  */
190 
191  /** Description of a unique dof of mean value type.
192  * @param d the dimension of the reference element.
193  */
194  pdof_description mean_value_dof(dim_type d);
195 
196  pdof_description bubble1_dof(dim_type ct);
197 
198  /** Description of a global dof, i.e. a numbered dof having a global scope.
199  * @param d the dimension of the reference element.
200  */
201  pdof_description global_dof(dim_type d);
202 
203  /// Product description of the descriptions *pnd1 and *pnd2.
205 
206  pdof_description to_coord_dof(pdof_description p, dim_type ct);
207 
208  /// Description of a special dof for Xfem
210 
211  /// Returns the xfem_index of dof (0 for normal dof)
213 
214  size_type reserve_xfem_index();
215  dim_type coord_index_of_dof(pdof_description a);
216 
217  int dof_weak_compatibility(pdof_description a, pdof_description b);
218 
219  /** Gives a total order on the dof description compatible with the
220  * identification.
221  */
223 
224  /// Says if the dof is linkable.
226 
227  /// Says if the two dofs can be identified.
229 
230 
231  /* @} */
232 
233  /* ******************************************************************** */
234  /* Classes for description of a finite element. */
235  /* ******************************************************************** */
236 
237  /** @defgroup pfem Finite element description
238  *
239  * @{
240  */
241 
242 
243  class virtual_fem;
244 
245  /** type of pointer on a fem description @see getfem::virtual_fem */
246  typedef std::shared_ptr<const getfem::virtual_fem> pfem;
247 
248  class fem_precomp_;
249  typedef std::shared_ptr<const getfem::fem_precomp_> pfem_precomp;
250 
252 
253  /** @brief Base class for finite element description */
254  class virtual_fem : virtual public dal::static_stored_object,
255  public std::enable_shared_from_this<const virtual_fem> {
256  public :
257  enum vec_type { VECTORIAL_NOTRANSFORM_TYPE, VECTORIAL_PRIMAL_TYPE,
258  VECTORIAL_DUAL_TYPE };
259 
260  protected :
261 
262  mutable std::vector<pdof_description> dof_types_;
263  /* this convex structure is "owned" by the virtual_fem
264  (a deep copy is made when virtual_fems are copied)
265  But cvs_node has to be a pointer, as bgeot::convex_structure
266  inherits from dal::static_stored_object
267  */
268  std::shared_ptr<bgeot::convex_structure> cvs_node;
269  std::vector<std::vector<short_type>> face_tab; // face list for each dof
270  bgeot::convex<base_node> cv_node;
271  mutable bgeot::pstored_point_tab pspt;
272  mutable bool pspt_valid;
273  bgeot::pconvex_ref cvr; // reference element.
274  dim_type ntarget_dim; // dimension of the target space
275  mutable dim_type dim_; // dimension of the reference element
276  bool is_equiv; // true if the FEM is equivalent
277  bool is_lag; // true if the FEM is of Lagrange type
278  bool is_pol; // true if the FEM is polynomial
279  bool is_polycomp; // true if the FEM is polynomial composite
280  bool real_element_defined;
281  bool is_standard_fem; // is_equiv && !real_element_defined && scalar
282  short_type es_degree; // estimated polynomial degree of the FEM
283  short_type hier_raff; // hierarchical refinement of the FEM
284  vec_type vtype; // for vectorial elements, type of transformation
285  // from the reference element.
286  std::string debug_name_;
287 
288  public :
289  /** Number of degrees of freedom.
290 
291  @param cv the convex number for this FEM. This information is
292  rarely used, but is needed by some "special" FEMs, such as
293  getfem::interpolated_fem.
294  */
295  virtual size_type nb_dof(size_type /*cv*/) const
296  { return dof_types_.size(); }
297  /// Number of basis functions.
298  virtual size_type nb_base(size_type cv) const
299  { return nb_dof(cv); }
300  /// Number of components (nb_dof() * dimension of the target space).
302  { return nb_base(cv) * ntarget_dim; }
303  size_type nb_components(size_type cv) const
304  { return nb_dof(cv) * ntarget_dim; }
305  /// Get the array of pointer on dof description.
306  const std::vector<pdof_description> &dof_types() const
307  { return dof_types_; }
308  short_type hierarchical_raff() const { return hier_raff; }
309  /// dimension of the reference element.
310  dim_type dim() const { return dim_; }
311  dim_type &dim() { return dim_; }
312  /// dimension of the target space.
313  dim_type target_dim() const { return ntarget_dim; }
314  /// Type of vectorial element.
315  vec_type vectorial_type() const { return vtype; }
316  /// Return the convex of the reference element.
317  virtual bgeot::pconvex_ref ref_convex(size_type) const { return cvr; }
318  /// @internal
319  bgeot::pconvex_ref &mref_convex() { return cvr; }
320  /// Gives the convex of the reference element.
322  { return ref_convex(cv)->structure(); }
323  /// Gives the convex representing the nodes on the reference element.
325  { return cv_node; }
326  /// Gives the convex structure of the reference element nodes.
328  { return node_convex(cv).structure(); }
329  const std::string &debug_name() const { return debug_name_; }
330  std::string &debug_name() { return debug_name_; }
331  virtual bgeot::pstored_point_tab node_tab(size_type) const {
332  if (!pspt_valid) {
333  pspt = bgeot::store_point_tab(cv_node.points());
334  pspt_valid = true;
335  }
336  return pspt;
337  }
338  /** Gives the node corresponding to the dof i.
339  @param cv the convex number for this FEM. This information is
340  rarely used, by is needed by some "special" FEMs, such as
341  getfem::interpolated_fem.
342  @param i the local dof number (<tt>i < nb_dof(cv)</tt>)
343  */
344  const base_node &node_of_dof(size_type cv, size_type i) const
345  { return (*(node_tab(cv)))[i];}
346  virtual const std::vector<short_type> &
347  faces_of_dof(size_type /*cv*/, size_type i) const;
348  bool is_on_real_element() const { return real_element_defined; }
349  bool is_equivalent() const { return is_equiv; }
350  bool need_G() const
351  { return !(is_equivalent()) || real_element_defined; }
352  /// true if the base functions are such that @f$ \varphi_i(\textrm{node\_of\_dof(j)}) = \delta_{ij} @f$
353  bool is_lagrange() const { return is_lag; }
354  /// true if the base functions are polynomials
355  bool is_polynomial() const { return is_pol; }
356  bool is_polynomialcomp() const { return is_polycomp; }
357  bool is_standard() const { return is_standard_fem; }
358  bool &is_polynomialcomp() { return is_polycomp; }
359  bool &is_equivalent() { return is_equiv; }
360  bool &is_lagrange() { return is_lag; }
361  bool &is_polynomial() { return is_pol; }
362  bool &is_standard() { return is_standard_fem; }
363  short_type estimated_degree() const { return es_degree; }
364  short_type &estimated_degree() { return es_degree; }
365 
366  virtual void mat_trans(base_matrix &, const base_matrix &,
368  { GMM_ASSERT1(false, "This function should not be called."); }
369  /** Interpolate at an arbitrary point x given on the reference
370  element.
371 
372  @param c the fem_interpolation_context, should have been
373  suitably initialized for the point of evaluation.
374 
375  @param coeff is the vector of coefficient relatively to
376  the base functions, its length should be @c Qdim*this->nb_dof().
377 
378  @param val contains the interpolated value, on output (its
379  size should be @c Qdim*this->target_dim()).
380 
381  @param Qdim is the optional Q dimension, if the FEM is
382  considered as a "vectorized" one.
383  */
384  template<typename CVEC, typename VVEC>
385  void interpolation(const fem_interpolation_context& c,
386  const CVEC& coeff, VVEC &val, dim_type Qdim) const;
387 
388  /** Build the interpolation matrix for the interpolation at a
389  fixed point x, given on the reference element.
390 
391  The matrix @c M is filled, such that for a given @c coeff
392  vector, the interpolation is given by @c M*coeff.
393  */
394  template <typename MAT>
395  void interpolation(const fem_interpolation_context& c,
396  MAT &M, dim_type Qdim) const;
397 
398  /** Interpolation of the gradient. The output is stored in the @f$
399  Q\times N@f$ matrix @c val.
400  */
401  template<typename CVEC, typename VMAT>
402  void interpolation_grad(const fem_interpolation_context& c,
403  const CVEC& coeff, VMAT &val,
404  dim_type Qdim=1) const;
405 
406  /** Interpolation of the hessian. The output is stored in the @f$
407  Q\times (N^2)@f$ matrix @c val.
408  */
409  template<typename CVEC, typename VMAT>
410  void interpolation_hess(const fem_interpolation_context& c,
411  const CVEC& coeff, VMAT &val,
412  dim_type Qdim) const;
413 
414  /** Interpolation of the divergence. The output is stored in the
415  scalar @c val.
416  */
417  template<typename CVEC>
419  (const fem_interpolation_context& c, const CVEC& coeff,
420  typename gmm::linalg_traits<CVEC>::value_type &val) const;
421 
422  /** Give the value of all components of the base functions at the
423  * point x of the reference element. Basic function used essentially
424  * by fem_precomp.
425  */
426  virtual void base_value(const base_node &x, base_tensor &t) const = 0;
427 
428  /** Give the value of all gradients (on ref. element) of the components
429  * of the base functions at the point x of the reference element.
430  * Basic function used essentially by fem_precomp.
431  */
432  virtual void grad_base_value(const base_node &x, base_tensor &t) const = 0;
433 
434  /** Give the value of all hessians (on ref. element) of the components
435  * of the base functions at the point x of the reference element.
436  * Basic function used essentially by fem_precomp.
437  */
438  virtual void hess_base_value(const base_node &x, base_tensor &t) const = 0;
439 
440  /** Give the value of all components of the base functions at the
441  current point of the fem_interpolation_context. Used by
442  elementary computations. if withM is false the matrix M for
443  non tau-equivalent elements is not taken into account.
444  */
445  virtual void real_base_value(const fem_interpolation_context &c,
446  base_tensor &t, bool withM = true) const;
447 
448  /** Give the gradient of all components of the base functions at the
449  current point of the fem_interpolation_context. Used by
450  elementary computations. if withM is false the matrix M for
451  non tau-equivalent elements is not taken into account.
452  */
453  virtual void real_grad_base_value(const fem_interpolation_context &c,
454  base_tensor &t, bool withM = true) const;
455 
456  /** Give the hessian of all components of the base functions at the
457  current point of the fem_interpolation_context. Used by
458  elementary computations. if withM is false the matrix M for
459  non tau-equivalent elements is not taken into account.
460  */
461  virtual void real_hess_base_value(const fem_interpolation_context &c,
462  base_tensor &t, bool withM = true) const;
463 
464  virtual size_type index_of_global_dof(size_type, size_type) const
465  { GMM_ASSERT1(false, "internal error."); }
466 
467  /** internal function adding a node to an element for the creation
468  * of a finite element method. Important : the faces should be the faces
469  * on which the corresponding base function is non zero.
470  */
471  void add_node(const pdof_description &d, const base_node &pt,
472  const dal::bit_vector &faces);
473  void add_node(const pdof_description &d, const base_node &pt);
474  void init_cvs_node();
475  void unfreeze_cvs_node();
476 
477  virtual_fem &operator =(const virtual_fem &f) {
478  copy(f); return *this;
479  }
480 
481  virtual_fem() {
482  DAL_STORED_OBJECT_DEBUG_CREATED(this, "Fem");
483  ntarget_dim = 1; dim_ = 1;
484  is_equiv = is_pol = is_polycomp = is_lag = is_standard_fem = false;
485  pspt_valid = false; hier_raff = 0; real_element_defined = false;
486  es_degree = 5;
487  vtype = VECTORIAL_NOTRANSFORM_TYPE;
488  cvs_node = bgeot::new_convex_structure();
489  }
490  virtual_fem(const virtual_fem& f)
491  : dal::static_stored_object(),
492  std::enable_shared_from_this<const virtual_fem>()
493  { copy(f); DAL_STORED_OBJECT_DEBUG_CREATED(this, "Fem"); }
494  virtual ~virtual_fem() { DAL_STORED_OBJECT_DEBUG_DESTROYED(this, "Fem"); }
495  private:
496  void copy(const virtual_fem &f);
497  };
498 
499  /**
500  virtual_fem implementation as a vector of generic functions. The
501  class FUNC should provide "derivative" and "eval" member
502  functions (this is the case for bgeot::polynomial<T>).
503  */
504  template <class FUNC> class fem : public virtual_fem {
505  protected :
506  std::vector<FUNC> base_;
507  mutable std::vector<std::vector<FUNC>> grad_, hess_;
508  mutable bool grad_computed_ = false;
509  mutable bool hess_computed_ = false;
510 
511  void compute_grad_() const {
512  if (grad_computed_) return;
513  GLOBAL_OMP_GUARD
514  if (grad_computed_) return;
516  dim_type n = dim();
517  grad_.resize(R);
518  for (size_type i = 0; i < R; ++i) {
519  grad_[i].resize(n);
520  for (dim_type j = 0; j < n; ++j) {
521  grad_[i][j] = base_[i]; grad_[i][j].derivative(j);
522  }
523  }
524  grad_computed_ = true;
525  }
526 
527  void compute_hess_() const {
528  if (hess_computed_) return;
529  GLOBAL_OMP_GUARD
530  if (hess_computed_) return;
532  dim_type n = dim();
533  hess_.resize(R);
534  for (size_type i = 0; i < R; ++i) {
535  hess_[i].resize(n*n);
536  for (dim_type j = 0; j < n; ++j) {
537  for (dim_type k = 0; k < n; ++k) {
538  hess_[i][j+k*n] = base_[i];
539  hess_[i][j+k*n].derivative(j); hess_[i][j+k*n].derivative(k);
540  }
541  }
542  }
543  hess_computed_ = true;
544  }
545 
546  public :
547 
548  /// Gives the array of basic functions (components).
549  const std::vector<FUNC> &base() const { return base_; }
550  std::vector<FUNC> &base() { return base_; }
551  /** Evaluates at point x, all base functions and returns the result in
552  t(nb_base,target_dim) */
553  void base_value(const base_node &x, base_tensor &t) const {
554  bgeot::multi_index mi(2);
555  mi[1] = target_dim(); mi[0] = short_type(nb_base(0));
556  t.adjust_sizes(mi);
558  base_tensor::iterator it = t.begin();
559  for (size_type i = 0; i < R; ++i, ++it)
560  *it = bgeot::to_scalar(base_[i].eval(x.begin()));
561  }
562  /** Evaluates at point x, the gradient of all base functions w.r.t. the
563  reference element directions 0,..,dim-1 and returns the result in
564  t(nb_base,target_dim,dim) */
565  void grad_base_value(const base_node &x, base_tensor &t) const {
566  if (!grad_computed_) compute_grad_();
567  bgeot::multi_index mi(3);
568  dim_type n = dim();
569  mi[2] = n; mi[1] = target_dim(); mi[0] = short_type(nb_base(0));
570  t.adjust_sizes(mi);
572  base_tensor::iterator it = t.begin();
573  for (dim_type j = 0; j < n; ++j)
574  for (size_type i = 0; i < R; ++i, ++it)
575  *it = bgeot::to_scalar(grad_[i][j].eval(x.begin()));
576  }
577  /** Evaluates at point x, the hessian of all base functions w.r.t. the
578  reference element directions 0,..,dim-1 and returns the result in
579  t(nb_base,target_dim,dim,dim) */
580  void hess_base_value(const base_node &x, base_tensor &t) const {
581  if (!hess_computed_) compute_hess_();
582  bgeot::multi_index mi(4);
583  dim_type n = dim();
584  mi[3] = n; mi[2] = n; mi[1] = target_dim();
585  mi[0] = short_type(nb_base(0));
586  t.adjust_sizes(mi);
588  base_tensor::iterator it = t.begin();
589  for (dim_type k = 0; k < n; ++k)
590  for (dim_type j = 0; j < n; ++j)
591  for (size_type i = 0; i < R; ++i, ++it)
592  *it = bgeot::to_scalar(hess_[i][j+k*n].eval(x.begin()));
593  }
594 
595  };
596 
597  /** Classical polynomial FEM. */
599  /** Polynomial composite FEM */
601  /** Rational fration FEM */
603 
604  /** Give a pointer on the structures describing the classical
605  polynomial fem of degree k on a given convex type.
606 
607  @param pgt the geometric transformation (which defines the convex type).
608  @param k the degree of the fem.
609  @param complete a flag which requests complete Langrange polynomial
610  elements even if the provided pgt is an incomplete one (e.g. 8-node
611  quadrilateral or 20-node hexahedral).
612  @return a ppolyfem.
613  */
615  bool complete=false);
616 
617  /** Give a pointer on the structures describing the classical
618  polynomial discontinuous fem of degree k on a given convex type.
619 
620  @param pgt the geometric transformation (which defines the convex type).
621 
622  @param k the degree of the fem.
623 
624  @param alpha the "inset" factor for the dof nodes: with alpha =
625  0, the nodes are located as usual (i.e. with node on the convex border),
626  and for 0 < alpha < 1, they converge to the center of gravity of the
627  convex.
628 
629  @param complete a flag which requests complete Langrange polynomial
630  elements even if the provided pgt is an incomplete one (e.g. 8-node
631  quadrilateral or 20-node hexahedral).
632 
633  @return a ppolyfem.
634  */
636  scalar_type alpha=0, bool complete=false);
637 
638  /** get a fem descriptor from its string name. */
639  pfem fem_descriptor(const std::string &name);
640 
641  /** get the string name of a fem descriptor. */
642  std::string name_of_fem(pfem p);
643 
644  pfem PK_fem(size_type n, short_type k);
645  pfem QK_fem(size_type n, short_type k);
646  pfem PK_prism_fem(size_type n, short_type k);
647 
648  /**
649  Pre-computations on a fem (given a fixed set of points on the
650  reference convex, this object computes the value/gradient/hessian
651  of all base functions on this set of points and stores them.
652  */
653  class fem_precomp_ : virtual public dal::static_stored_object {
654  protected:
655  const pfem pf;
656  const bgeot::pstored_point_tab pspt;
657  mutable std::vector<base_tensor> c; // stored values of base functions
658  mutable std::vector<base_tensor> pc; // stored gradients of base functions
659  mutable std::vector<base_tensor> hpc; // stored hessians of base functions
660  public:
661  /// returns values of the base functions
662  inline const base_tensor &val(size_type i) const
663  { if (c.empty()) init_val(); return c[i]; }
664  /// returns gradients of the base functions
665  inline const base_tensor &grad(size_type i) const
666  { if (pc.empty()) init_grad(); return pc[i]; }
667  /// returns hessians of the base functions
668  inline const base_tensor &hess(size_type i) const
669  { if (hpc.empty()) init_hess(); return hpc[i]; }
670  inline pfem get_pfem() const { return pf; }
671  // inline const bgeot::stored_point_tab& get_point_tab() const
672  // { return *pspt; }
673  inline bgeot::pstored_point_tab get_ppoint_tab() const
674  { return pspt; }
675  fem_precomp_(const pfem, const bgeot::pstored_point_tab);
676  ~fem_precomp_() { DAL_STORED_OBJECT_DEBUG_DESTROYED(this, "Fem_precomp"); }
677  private:
678  void init_val() const;
679  void init_grad() const;
680  void init_hess() const;
681  };
682 
683 
684  /** @brief Handles precomputations for FEM. statically allocates a
685  fem-precomputation object, and returns a pointer to it. The
686  fem_precomp_ objects are "cached", i.e. they are stored in a
687  global pool and if this function is called two times with the
688  same arguments, a pointer to the same object will be returned.
689 
690  @param pf a pointer to the fem object. @param pspt a pointer to
691  a list of points in the reference convex.CAUTION: this array must
692  not be destroyed as long as the fem_precomp is used!!.
693 
694  Moreover pspt is supposed to identify uniquely the set of
695  points. This means that you should NOT alter its content at any
696  time after using this function.
697 
698  If you need a set of "temporary" getfem::fem_precomp_, create
699  them via a getfem::fem_precomp_pool structure. All memory will be
700  freed when this structure will be destroyed. */
701  pfem_precomp fem_precomp(pfem pf, bgeot::pstored_point_tab pspt,
702  dal::pstatic_stored_object dep);
703 
704  /** Request for the removal of a pfem_precomp */
705  inline void delete_fem_precomp(pfem_precomp pfp)
706  { dal::del_stored_object(pfp); }
707 
708 
709  /**
710  handle a pool (i.e. a set) of fem_precomp. The difference with
711  the global fem_precomp function is that these fem_precomp objects
712  are freed when the fem_precomp_pool is destroyed (they can eat
713  much memory). An example of use can be found in the
714  getfem::interpolation_solution functions of getfem_export.h
715  */
717  std::set<pfem_precomp> precomps;
718 
719  public :
720 
721  /** Request a pfem_precomp. If not already in the pool, the
722  pfem_precomp is computed, and added to the pool.
723 
724  @param pf a pointer to the fem object.
725  @param pspt a pointer to a list of points in the reference convex.
726 
727  CAUTION:
728  this array must not be destroyed as long as the fem_precomp is used!!
729 
730  Moreover pspt is supposed to identify uniquely the set of
731  points. This means that you should NOT alter its content until
732  the fem_precomp_pool is destroyed.
733  */
734  pfem_precomp operator()(pfem pf, bgeot::pstored_point_tab pspt) {
735  pfem_precomp p = fem_precomp(pf, pspt, 0);
736  precomps.insert(p);
737  return p;
738  }
739  void clear();
740  ~fem_precomp_pool() { clear(); }
741  };
742 
743 
744  /** structure passed as the argument of fem interpolation
745  functions. This structure can be partially filled (for example
746  the xreal will be computed if needed as long as pgp+ii is known).
747  */
750 
751  mutable base_matrix M_; // optional transformation matrix (for non tau-equivalent fems)
752  pfem pf_; // current fem
753  pfem_precomp pfp_; // optional fem_precomp_ (speed up the computations)
754  size_type convex_num_; // The element (convex) number
755  short_type face_num_; // Face number for boundary integration
756  int xfem_side_; // For the computation of a jump with fem_level_set only
757  public:
758  /// true if a fem_precomp_ has been supplied.
759  bool have_pfp() const { return pfp_ != 0; }
760  /// true if the pfem is available.
761  bool have_pf() const { return pf_ != 0; }
762  /// non tau-equivalent transformation matrix.
763  const base_matrix& M() const;
764  /** fill the tensor with the values of the base functions (taken
765  at point @c this->xref())
766  */
767  void base_value(base_tensor& t, bool withM = true) const;
768  // Optimized function for high level generic assembly
769  void pfp_base_value(base_tensor& t, const pfem_precomp &pfp__);
770  /** fill the tensor with the gradient of the base functions (taken
771  at point @c this->xref())
772  */
773  void grad_base_value(base_tensor& t, bool withM = true) const;
774  // Optimized function for high level generic assembly
775  void pfp_grad_base_value(base_tensor& t, const pfem_precomp &pfp__);
776  /** fill the tensor with the hessian of the base functions (taken
777  at point @c this->xref())
778  */
779  void hess_base_value(base_tensor& t, bool withM = true) const;
780  /** get the current FEM descriptor */
781  const pfem pf() const { return pf_; }
782  /** get the current convex number */
783  size_type convex_num() const;
784  bool is_convex_num_valid() const;
785  void invalid_convex_num() { convex_num_ = size_type(-1); }
786  /** set the current face number */
787  void set_face_num(short_type f) { face_num_ = f; }
788  /** get the current face number */
789  short_type face_num() const;
790  /** On a face ? */
791  bool is_on_face() const;
792  /** get the current fem_precomp_ */
793  pfem_precomp pfp() const { return pfp_; }
794  void set_pfp(pfem_precomp newpfp);
795  void set_pf(pfem newpf);
796  int xfem_side() const { return xfem_side_; }
797  void set_xfem_side(int side) { xfem_side_ = side; }
798  void change(bgeot::pgeotrans_precomp pgp__,
799  pfem_precomp pfp__, size_type ii__,
800  const base_matrix& G__, size_type convex_num__,
801  short_type face_num__ = short_type(-1)) {
802  bgeot::geotrans_interpolation_context::change(pgp__,ii__,G__);
803  convex_num_ = convex_num__; face_num_ = face_num__; xfem_side_ = 0;
804  set_pfp(pfp__);
805  }
806  void change(bgeot::pgeometric_trans pgt__,
807  pfem_precomp pfp__, size_type ii__,
808  const base_matrix& G__, size_type convex_num__,
809  short_type face_num__ = short_type(-1)) {
810  bgeot::geotrans_interpolation_context::change
811  (pgt__, pfp__->get_ppoint_tab(), ii__, G__);
812  convex_num_ = convex_num__; face_num_ = face_num__; xfem_side_ = 0;
813  set_pfp(pfp__);
814  }
815  void change(bgeot::pgeometric_trans pgt__,
816  pfem pf__, const base_node& xref__, const base_matrix& G__,
817  size_type convex_num__, short_type face_num__=short_type(-1)) {
818  bgeot::geotrans_interpolation_context::change(pgt__,xref__,G__);
819  pf_ = pf__; pfp_ = 0; convex_num_ = convex_num__; face_num_ = face_num__;
820  xfem_side_ = 0;
821  }
822  fem_interpolation_context()
823  : bgeot::geotrans_interpolation_context(),
824  convex_num_(size_type(-1)), face_num_(short_type(-1)), xfem_side_(0) {}
825  fem_interpolation_context(bgeot::pgeotrans_precomp pgp__,
826  pfem_precomp pfp__, size_type ii__,
827  const base_matrix& G__,
828  size_type convex_num__,
829  short_type face_num__ = short_type(-1))
830  : bgeot::geotrans_interpolation_context(pgp__,ii__,G__),
831  convex_num_(convex_num__), face_num_(face_num__), xfem_side_(0)
832  { set_pfp(pfp__); }
833  fem_interpolation_context(bgeot::pgeometric_trans pgt__,
834  pfem_precomp pfp__, size_type ii__,
835  const base_matrix& G__,
836  size_type convex_num__,
837  short_type face_num__ = short_type(-1))
838  : bgeot::geotrans_interpolation_context(pgt__,pfp__->get_ppoint_tab(),
839  ii__, G__),
840  convex_num_(convex_num__), face_num_(face_num__), xfem_side_(0)
841  { set_pfp(pfp__); }
842  fem_interpolation_context(bgeot::pgeometric_trans pgt__,
843  pfem pf__,
844  const base_node& xref__,
845  const base_matrix& G__,
846  size_type convex_num__,
847  short_type face_num__ = short_type(-1))
848  : bgeot::geotrans_interpolation_context(pgt__,xref__,G__),
849  pf_(pf__), pfp_(0), convex_num_(convex_num__), face_num_(face_num__),
850  xfem_side_(0) {}
851  };
852 
853  // IN : coeff(Qmult,nb_dof)
854  // OUT: val(Qdim), Qdim=target_dim*Qmult
855  // AUX: Z(nb_dof,target_dim)
856  template <typename CVEC, typename VVEC>
858  const CVEC& coeff, VVEC &val,
859  dim_type Qdim) const {
860  size_type Qmult = size_type(Qdim) / target_dim();
861  size_type nbdof = nb_dof(c.convex_num());
862  GMM_ASSERT1(gmm::vect_size(val) == Qdim, "dimensions mismatch");
863  GMM_ASSERT1(gmm::vect_size(coeff) == nbdof*Qmult,
864  "Wrong size for coeff vector");
865 
866  gmm::clear(val);
867  base_tensor Z; real_base_value(c, Z);
868 
869  for (size_type j = 0; j < nbdof; ++j) {
870  for (size_type q = 0; q < Qmult; ++q) {
871  typename gmm::linalg_traits<CVEC>::value_type co = coeff[j*Qmult+q];
872  for (size_type r = 0; r < target_dim(); ++r)
873  val[r + q*target_dim()] += co * Z[j + r*nbdof];
874  }
875  }
876  }
877 
878  template <typename MAT>
880  MAT &M, dim_type Qdim) const {
881  size_type Qmult = size_type(Qdim) / target_dim();
882  size_type nbdof = nb_dof(c.convex_num());
883  GMM_ASSERT1(gmm::mat_nrows(M) == Qdim && gmm::mat_ncols(M) == nbdof*Qmult,
884  "dimensions mismatch");
885 
886  gmm::clear(M);
887  base_tensor Z; real_base_value(c, Z);
888  for (size_type j = 0; j < nbdof; ++j) {
889  for (size_type q = 0; q < Qmult; ++q) {
890  for (size_type r = 0; r < target_dim(); ++r)
891  M(r+q*target_dim(), j*Qmult+q) = Z[j + r*nbdof];
892  }
893  }
894  }
895 
896 
897  // IN : coeff(Qmult,nb_dof)
898  // OUT: val(Qdim,N), Qdim=target_dim*Qmult
899  // AUX: t(nb_dof,target_dim,N)
900  template<typename CVEC, typename VMAT>
902  const CVEC& coeff, VMAT &val,
903  dim_type Qdim) const {
904  size_type N = c.N();
905  size_type nbdof = nb_dof(c.convex_num());
906  size_type Qmult = gmm::vect_size(coeff) / nbdof;
907  GMM_ASSERT1(gmm::mat_ncols(val) == N &&
908  gmm::mat_nrows(val) == target_dim()*Qmult &&
909  gmm::vect_size(coeff) == nbdof*Qmult,
910  "dimensions mismatch");
911  GMM_ASSERT1(Qdim == target_dim()*Qmult, // Qdim seems to be superfluous input, could be removed in the future
912  "dimensions mismatch");
913  base_tensor t;
914  real_grad_base_value(c, t); // t(nbdof,target_dim,N)
915 
916  gmm::clear(val);
917  for (size_type q = 0; q < Qmult; ++q) {
918  base_tensor::const_iterator it = t.begin();
919  for (size_type k = 0; k < N; ++k)
920  for (size_type r = 0; r < target_dim(); ++r)
921  for (size_type j = 0; j < nbdof; ++j, ++it)
922  val(r + q*target_dim(), k) += coeff[j*Qmult+q] * (*it);
923  }
924  }
925 
926 
927  template<typename CVEC, typename VMAT>
929  const CVEC& coeff, VMAT &val,
930  dim_type Qdim) const {
931  size_type Qmult = size_type(Qdim) / target_dim();
932  size_type N = c.N();
933  GMM_ASSERT1(gmm::mat_ncols(val) == N*N &&
934  gmm::mat_nrows(val) == Qdim, "dimensions mismatch");
935 
936  base_tensor t;
937  size_type nbdof = nb_dof(c.convex_num());
938 
939  gmm::clear(val);
940  real_hess_base_value(c, t);
941  for (size_type q = 0; q < Qmult; ++q) {
942  base_tensor::const_iterator it = t.begin();
943  for (size_type k = 0; k < N*N; ++k)
944  for (size_type r = 0; r < target_dim(); ++r)
945  for (size_type j = 0; j < nbdof; ++j, ++it)
946  val(r + q*target_dim(), k) += coeff[j*Qmult+q] * (*it);
947  }
948  }
949 
950 
951  // IN : coeff(Qmult,nb_dof)
952  // OUT: val
953  // AUX: t(nb_dof,target_dim,N), Qmult*target_dim == N
954  template<typename CVEC>
956  (const fem_interpolation_context& c, const CVEC& coeff,
957  typename gmm::linalg_traits<CVEC>::value_type &val) const {
958  size_type N = c.N();
959  size_type nbdof = nb_dof(c.convex_num());
960  size_type Qmult = gmm::vect_size(coeff) / nbdof;
961  GMM_ASSERT1(gmm::vect_size(coeff) == nbdof*Qmult , "dimensions mismatch");
962  GMM_ASSERT1(target_dim()*Qmult == N &&
963  (Qmult == 1 || target_dim() == 1),
964  "Dimensions mismatch. Divergence operator requires fem qdim equal to dim.");
965  base_tensor t;
966  real_grad_base_value(c, t); // t(nbdof,target_dim,N)
967  // for Qmult == 1 this is sub-optimal since it evaluates all (:,i,j)
968  // gradients instead of only the diagonal ones(:,i,i)
969 
970  val = scalar_type(0);
971  base_tensor::const_iterator it = t.begin();
972  if (Qmult == 1)
973  for (size_type k = 0; k < N; ++k) {
974  if (k) it += (N*nbdof + 1);
975  for (size_type j = 0; j < nbdof; ++j) {
976  if (j) ++it;
977  val += coeff[j] * (*it);
978  }
979  }
980  else // if (target_dim() == 1)
981  for (size_type k = 0; k < N; ++k) {
982  if (k) ++it;
983  for (size_type j = 0; j < nbdof; ++j) {
984  if (j) ++it;
985  val += coeff[j*N+k] * (*it);
986  }
987  }
988  }
989 
990  /**
991  Specific function for a HHO method to obtain the method in the interior.
992  If the method is not of composite type, return the argument.
993  */
995 
996 
997  /* Functions allowing the add of a finite element method outside
998  of getfem_fem.cc */
999 
1000  typedef dal::naming_system<virtual_fem>::param_list fem_param_list;
1001 
1002  void inline read_poly(bgeot::base_poly &p, int d, const char *s)
1003  { p = bgeot::read_base_poly(short_type(d), s); }
1004 
1005  void add_fem_name(std::string name,
1007 
1008 
1009  /* @} */
1010 
1011 } /* end of namespace getfem. */
1012 
1013 
1014 #endif
Geometric transformations on convexes.
Handle composite polynomials.
the geotrans_interpolation_context structure is passed as the argument of geometric transformation in...
Associate a name to a method descriptor and store method descriptors.
base class for static stored objects
structure passed as the argument of fem interpolation functions.
Definition: getfem_fem.h:749
Pre-computations on a fem (given a fixed set of points on the reference convex, this object computes ...
Definition: getfem_fem.h:653
handle a pool (i.e.
Definition: getfem_fem.h:716
virtual_fem implementation as a vector of generic functions.
Definition: getfem_fem.h:504
Base class for finite element description.
Definition: getfem_fem.h:255
Naming system.
Stores interdependent getfem objects.
Integration methods (exact and approximated) on convexes.
pdof_description derivative_dof(dim_type d, dim_type r)
Description of a unique dof of derivative type.
Definition: getfem_fem.cc:487
int dof_description_compare(pdof_description a, pdof_description b)
Gives a total order on the dof description compatible with the identification.
Definition: getfem_fem.cc:606
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.
Definition: getfem_fem.cc:496
pdof_description mean_value_dof(dim_type d)
Description of a unique dof of mean value type.
Definition: getfem_fem.cc:534
bool dof_linkable(pdof_description)
Says if the dof is linkable.
Definition: getfem_fem.cc:615
bool dof_compatibility(pdof_description, pdof_description)
Says if the two dofs can be identified.
Definition: getfem_fem.cc:618
pdof_description xfem_dof(pdof_description p, size_type ind)
Description of a special dof for Xfem.
Definition: getfem_fem.cc:433
pdof_description lagrange_dof(dim_type d)
Description of a unique dof of lagrange type (value at the node).
Definition: getfem_fem.cc:390
size_type dof_xfem_index(pdof_description)
Returns the xfem_index of dof (0 for normal dof)
Definition: getfem_fem.cc:621
pdof_description global_dof(dim_type d)
Description of a global dof, i.e.
Definition: getfem_fem.cc:542
pdof_description normal_derivative_dof(dim_type d)
Description of a unique dof of normal derivative type (normal derivative at the node,...
Definition: getfem_fem.cc:507
dof_description * pdof_description
Type representing a pointer on a dof_description.
Definition: getfem_fem.h:159
pdof_description product_dof(pdof_description pnd1, pdof_description pnd2)
Product description of the descriptions *pnd1 and *pnd2.
Definition: getfem_fem.cc:552
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.
Definition: getfem_fem.h:580
bool is_lagrange() const
true if the base functions are such that
Definition: getfem_fem.h:353
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.
Definition: getfem_fem.cc:648
virtual void hess_base_value(const base_node &x, base_tensor &t) const =0
Give the value of all hessians (on ref.
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...
Definition: getfem_fem.cc:317
dim_type target_dim() const
dimension of the target space.
Definition: getfem_fem.h:313
size_type convex_num() const
get the current convex number
Definition: getfem_fem.cc:52
std::shared_ptr< const getfem::virtual_fem > pfem
type of pointer on a fem description
Definition: getfem_fem.h:243
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...
Definition: getfem_fem.cc:313
void delete_fem_precomp(pfem_precomp pfp)
Request for the removal of a pfem_precomp.
Definition: getfem_fem.h:705
const fem< bgeot::base_rational_fraction > * prationalfracfem
Rational fration FEM.
Definition: getfem_fem.h:602
const base_node & node_of_dof(size_type cv, size_type i) const
Gives the node corresponding to the dof i.
Definition: getfem_fem.h:344
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.
Definition: getfem_fem.h:565
vec_type vectorial_type() const
Type of vectorial element.
Definition: getfem_fem.h:315
bool have_pfp() const
true if a fem_precomp_ has been supplied.
Definition: getfem_fem.h:759
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)
Definition: getfem_fem.h:553
const base_matrix & M() const
non tau-equivalent transformation matrix.
Definition: getfem_fem.cc:43
void interpolation_diverg(const fem_interpolation_context &c, const CVEC &coeff, typename gmm::linalg_traits< CVEC >::value_type &val) const
Interpolation of the divergence.
Definition: getfem_fem.h:956
void interpolation_hess(const fem_interpolation_context &c, const CVEC &coeff, VMAT &val, dim_type Qdim) const
Interpolation of the hessian.
Definition: getfem_fem.h:928
virtual size_type nb_base(size_type cv) const
Number of basis functions.
Definition: getfem_fem.h:298
const fem< bgeot::polynomial_composite > * ppolycompfem
Polynomial composite FEM.
Definition: getfem_fem.h:600
pfem fem_descriptor(const std::string &name)
get a fem descriptor from its string name.
Definition: getfem_fem.cc:4659
const base_tensor & grad(size_type i) const
returns gradients of the base functions
Definition: getfem_fem.h:665
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())
Definition: getfem_fem.cc:255
bool have_pf() const
true if the pfem is available.
Definition: getfem_fem.h:761
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...
Definition: getfem_fem.cc:4571
bgeot::pconvex_structure structure(size_type cv) const
Gives the convex structure of the reference element nodes.
Definition: getfem_fem.h:327
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())
Definition: getfem_fem.cc:202
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...
Definition: getfem_fem.cc:4566
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())
Definition: getfem_fem.cc:120
short_type face_num() const
get the current face number
Definition: getfem_fem.cc:61
bgeot::pconvex_structure basic_structure(size_type cv) const
Gives the convex of the reference element.
Definition: getfem_fem.h:321
const std::vector< pdof_description > & dof_types() const
Get the array of pointer on dof description.
Definition: getfem_fem.h:306
virtual const bgeot::convex< base_node > & node_convex(size_type) const
Gives the convex representing the nodes on the reference element.
Definition: getfem_fem.h:324
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.
Definition: getfem_fem.h:857
void interpolation_grad(const fem_interpolation_context &c, const CVEC &coeff, VMAT &val, dim_type Qdim=1) const
Interpolation of the gradient.
Definition: getfem_fem.h:901
const std::vector< FUNC > & base() const
Gives the array of basic functions (components).
Definition: getfem_fem.h:549
pfem interior_fem_of_hho_method(pfem hho_method)
Specific function for a HHO method to obtain the method in the interior.
bool is_polynomial() const
true if the base functions are polynomials
Definition: getfem_fem.h:355
const pfem pf() const
get the current FEM descriptor
Definition: getfem_fem.h:781
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.
pfem_precomp fem_precomp(pfem pf, bgeot::pstored_point_tab pspt, dal::pstatic_stored_object dep)
Handles precomputations for FEM.
Definition: getfem_fem.cc:4760
size_type nb_base_components(size_type cv) const
Number of components (nb_dof() * dimension of the target space).
Definition: getfem_fem.h:301
void set_face_num(short_type f)
set the current face number
Definition: getfem_fem.h:787
virtual void grad_base_value(const base_node &x, base_tensor &t) const =0
Give the value of all gradients (on ref.
bool is_on_face() const
On a face ?
Definition: getfem_fem.cc:67
pfem_precomp pfp() const
get the current fem_precomp_
Definition: getfem_fem.h:793
pfem_precomp operator()(pfem pf, bgeot::pstored_point_tab pspt)
Request a pfem_precomp.
Definition: getfem_fem.h:734
dim_type dim() const
dimension of the reference element.
Definition: getfem_fem.h:310
virtual bgeot::pconvex_ref ref_convex(size_type) const
Return the convex of the reference element.
Definition: getfem_fem.h:317
const base_tensor & val(size_type i) const
returns values of the base functions
Definition: getfem_fem.h:662
const fem< bgeot::base_poly > * ppolyfem
Classical polynomial FEM.
Definition: getfem_fem.h:598
virtual size_type nb_dof(size_type) const
Number of degrees of freedom.
Definition: getfem_fem.h:295
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_...
Definition: getfem_fem.cc:309
const base_tensor & hess(size_type i) const
returns hessians of the base functions
Definition: getfem_fem.h:668
std::string name_of_fem(pfem p)
get the string name of a fem descriptor.
Definition: getfem_fem.cc:4668
Basic Geometric Tools.
gmm::uint16_type short_type
used as the common short type integer in the library
Definition: bgeot_config.h:72
base_poly read_base_poly(short_type n, std::istream &f)
read a base_poly on the stream ist.
Definition: bgeot_poly.cc:210
std::shared_ptr< const convex_structure > pconvex_structure
Pointer on a convex structure description.
size_t size_type
used as the common size type in the library
Definition: bgeot_poly.h:48
std::shared_ptr< const bgeot::geometric_trans > pgeometric_trans
pointer type for a geometric transformation
Dynamic Array Library.
void del_stored_object(const pstatic_stored_object &o, bool ignore_unstored)
Delete an object and the object which depend on it.
GEneric Tool for Finite Element Methods.