GetFEM  5.5
getfem_integration.h
Go to the documentation of this file.
1 /* -*- c++ -*- (enables emacs c++ mode) */
2 /*===========================================================================
3 
4  Copyright (C) 2000-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_integration.h
32  @author Yves Renard <Yves.Renard@insa-lyon.fr>
33  @date December 17, 2000.
34  @brief Integration methods (exact and approximated) on convexes
35 
36  @section im_list List of integration methods for getfem::int_method_descriptor
37 
38  - "IM_EXACT_SIMPLEX(n)" : exact integration on simplexes.
39  - "IM_PRODUCT(IM1, IM2)" : product of two integration methods
40  - "IM_EXACT_PARALLELEPIPED(n)" : exact integration on parallelepipeds
41  - "IM_EXACT_PRISM(n)" : exact integration on prisms
42  - "IM_GAUSS1D(K)" : Gauss method on the segment,
43  order K=1, 3, 5, 7, ..., 99.
44  - "IM_GAUSSLOBATTO1D(K)" : Gauss-Lobatto method on the segment,
45  order K=1, 3, 5, 7, ..., 99.
46  - "IM_NC(N,K)" : Newton-Cotes approximative
47  integration on simplexes, order K
48  - "IM_NC_PARALLELEPIPED(N,K)" : product of Newton-Cotes integration
49  on parallelepipeds
50  - "IM_NC_PRISM(N,K)" : product of Newton-Cotes integration
51  on prisms
52  - "IM_GAUSS_PARALLELEPIPED(N,K)" : product of Gauss1D integration
53  on parallelepipeds (K=1,3,..99)
54  - "IM_TRIANGLE(K)" : Gauss methods on triangles
55  (K=1..10,13,17,19)
56  - "IM_QUAD(K)" : Gauss methods on quadrilaterons
57  (K=2, 3, 5, 7, 9, 17)
58  - "IM_HEXAHEDRON(K)" : Gauss methods on hexahedrons
59  (K=5,9,11)
60  - "IM_TETRAHEDRON(K)" : Gauss methods on tetrahedrons
61  (K=1, 2, 3, 5, 6, or 8)
62  - "IM_SIMPLEX4D(3)" : Gauss method on a 4-dimensional
63  simplex.
64  - "IM_CUBE4D(K)" : Gauss method on a 4-dimensional cube
65  (K=5,9).
66  - "IM_STRUCTURED_COMPOSITE(IM1, K)" : Composite method on a grid with
67  K divisions
68  - "IM_HCT_COMPOSITE(IM1)" : Composite integration suited to the
69  HCT composite finite element.
70  - "IM_QUAQC1_COMPOSITE(IM1)" : Composite integration suited to the
71  QUADC1 composite finite element.
72 
73  - "IM_QUASI_POLAR(IM1, IP1)" : if IM1 is an integration method on a
74  square, gives an integration method on a triangle which is
75  close to a polar integration with respect to vertex IP1.
76  if IM1 is an integration method on a tetrahedron, gives an
77  integration method on a tetrahedron which is close to a
78  cylindrical integration with respect to vertex IP1
79  (does not work very well).
80  if IM1 is an integration method on a prism. Gives an integration
81  method on a tetrahedron which is close to a
82  cylindrical integration with respect to vertex IP1.
83 
84  - "IM_QUASI_POLAR(IM1, IP1, IP2)" : IM1 should be an integration method
85  on a prism. Gives an integration method on a tetrahedron which
86  is close to a cylindrical integration with respect to IP1-IP2
87  axis.
88 
89  - "IM_PYRAMID_COMPOSITE(IM1)" : Composite integration for a pyramid
90  decomposed into two tetrahedrons
91 */
92 #ifndef GETFEM_INTEGRATION_H__
93 #define GETFEM_INTEGRATION_H__
94 
95 #include "getfem_config.h"
96 #include "bgeot_convex_ref.h"
97 #include "bgeot_geometric_trans.h"
98 #include "bgeot_node_tab.h"
99 #include "bgeot_poly_composite.h"
101 
102 
103 
104 
105 namespace getfem
106 {
107  /** Description of an exact integration of polynomials.
108  * This class is not to be manipulate by itself. Use ppoly_integration and
109  * the functions written to produce the basic descriptions.
110  */
112  protected :
113 
115 
116  mutable std::vector<long_scalar_type> int_monomials;
117  mutable std::vector< std::vector<long_scalar_type> > int_face_monomials;
118 
119  public :
120 
121  /// Dimension of convex of reference.
122  dim_type dim(void) const { return cvs->dim(); }
123  /// {Structure of convex of reference.
124  bgeot::pconvex_structure structure(void) const { return cvs; }
125 
126  virtual long_scalar_type int_monomial(const bgeot::power_index &power)
127  const = 0;
128  virtual long_scalar_type int_monomial_on_face(const bgeot::power_index
129  &power, short_type f) const = 0;
130 
131  /// Evaluate the integral of the polynomial P on the reference element.
132  long_scalar_type int_poly(const base_poly &P) const;
133  /** Evaluate the integral of the polynomial P on the face f of the
134  * reference element.
135  */
136  long_scalar_type int_poly_on_face(const base_poly &P,
137  short_type f) const;
138 
139  virtual ~poly_integration() {}
140 
141  };
142 
143  typedef std::shared_ptr<const poly_integration> ppoly_integration;
144 
145  /** Description of an approximate integration of polynomials of
146  * several variables on a reference element.
147  * This class is not to be manipulate by itself. Use papprox\_integration
148  * and the functions written to produce the basic descriptions.
149  */
150 
151  class integration_method;
152  typedef std::shared_ptr<const integration_method> pintegration_method;
153 
154 
155  class approx_integration {
156  protected :
157 
158  typedef bgeot::node_tab PT_TAB;
159  bgeot::pconvex_ref cvr;
160  bgeot::pstored_point_tab pint_points;
161  std::vector<scalar_type> int_coeffs;
162  std::vector<size_type> repartition;
163 
164  // index 0 : points for volumic integration, index > 0 : points for faces
165  std::vector<PT_TAB> pt_to_store;
166  bool valid;
167  bool built_on_the_fly;
168 
169  public :
170 
171  /// Dimension of reference convex.
172  dim_type dim(void) const { return cvr->structure()->dim(); }
173  size_type nb_points(void) const { return int_coeffs.size(); }
174  /// Number of integration nodes on the reference element.
175  size_type nb_points_on_convex(void) const { return repartition[0]; }
176  /// Number of integration nodes on the face f of the reference element.
177  size_type nb_points_on_face(short_type f) const
178  { return repartition[f+1] - repartition[f]; }
179  size_type ind_first_point_on_face(short_type f) const
180  { return repartition[f]; }
181  /// Convenience method returning the number of faces of the reference convex
182  short_type nb_convex_faces() const
183  { return cvr->structure()->nb_faces(); /* == repartition.size() - 1;*/ }
184  bool is_built_on_the_fly(void) const { return built_on_the_fly; }
185  void set_built_on_the_fly(void) { built_on_the_fly = true; }
186  /// Structure of the reference element.
187  bgeot::pconvex_structure structure(void) const
188  { return basic_structure(cvr->structure()); }
189  bgeot::pconvex_ref ref_convex(void) const { return cvr; }
190 
191  const std::vector<size_type> &repart(void) const { return repartition; }
192 
193  /// Gives an array of integration nodes.
194  // const bgeot::stored_point_tab &integration_points(void) const
195  // { return *(pint_points); }
196  bgeot::pstored_point_tab pintegration_points(void) const
197  { return pint_points; }
198  /// Gives the integration node i on the reference element.
199  const base_node &point(size_type i) const
200  { return (*pint_points)[i]; }
201  /// Gives the integration node i of the face f.
202  const base_node &
203  point_on_face(short_type f, size_type i) const
204  { return (*pint_points)[repartition[f] + i]; }
205  /// Gives an array of the integration coefficients.
206  const std::vector<scalar_type> &integration_coefficients(void) const
207  { return int_coeffs; }
208  /// Gives the integration coefficient corresponding to node i.
209  scalar_type coeff(size_type i) const { return int_coeffs[i]; }
210  /// Gives the integration coefficient corresponding to node i of face f.
211  scalar_type coeff_on_face(short_type f, size_type i) const
212  { return int_coeffs[repartition[f] + i]; }
213 
214  void add_point(const base_node &pt, scalar_type w,
215  short_type f = short_type(-1),
216  bool include_empty = false);
217  void add_point_norepeat(const base_node &pt, scalar_type w,
218  short_type f=short_type(-1));
219  void add_point_full_symmetric(base_node pt, scalar_type w);
220  void add_method_on_face(pintegration_method ppi, short_type f);
221  void valid_method(void);
222 
223  approx_integration(void) : valid(false), built_on_the_fly(false) { }
224  approx_integration(bgeot::pconvex_ref cr)
225  : cvr(cr), repartition(cr->structure()->nb_faces()+1),
226  pt_to_store(cr->structure()->nb_faces()+1), valid(false),
227  built_on_the_fly(false)
228  { std::fill(repartition.begin(), repartition.end(), 0); }
229  virtual ~approx_integration() {}
230  };
231 
232  typedef std::shared_ptr<const approx_integration> papprox_integration;
233 
234  /**
235  the list of main integration method types
236  */
237  typedef enum { IM_APPROX, IM_EXACT, IM_NONE } integration_method_type;
238 
239  /**
240  this structure is not intended to be used directly. It is built via
241  the int_method_descriptor() function
242  */
244  ppoly_integration ppi; /* for exact integrations */
245  papprox_integration pai; /* for approximate integrations (i.e. cubatures) */
246  integration_method_type im_type;
247  void remove() { pai = papprox_integration(); ppi = ppoly_integration(); }
248 
249  public:
250  integration_method_type type(void) const { return im_type; }
251  const papprox_integration &approx_method(void) const { return pai; }
252  const ppoly_integration &exact_method(void) const { return ppi; }
253 
254  void set_approx_method(papprox_integration paii)
255  { remove(); pai = paii; im_type = IM_APPROX; }
256  void set_exact_method(ppoly_integration ppii)
257  { remove(); ppi = ppii; im_type = IM_EXACT; }
258 
259  bgeot::pstored_point_tab pintegration_points(void) const {
260  if (type() == IM_EXACT) {
261  size_type n = ppi->structure()->dim();
262  std::vector<base_node> spt(1); spt[0] = base_node(n);
263  return store_point_tab(spt);
264  }
265  else if (type() == IM_APPROX)
266  return pai->pintegration_points();
267  else GMM_ASSERT1(false, "IM_NONE has no points");
268  }
269 
270  // const bgeot::stored_point_tab &integration_points(void) const
271  // { return *(pintegration_points()); }
272 
273  bgeot::pconvex_structure structure(void) const {
274  switch (type()) {
275  case IM_EXACT: return ppi->structure();
276  case IM_APPROX: return pai->structure();
277  case IM_NONE: GMM_ASSERT1(false, "IM_NONE has no structure");
278  default: GMM_ASSERT3(false, "");
279  }
280  return 0;
281  }
282 
283  integration_method(ppoly_integration p) {
284  DAL_STORED_OBJECT_DEBUG_CREATED(this, "Exact integration method");
285  ppi = p; im_type = IM_EXACT;
286  }
287 
288  integration_method(papprox_integration p) {
289  DAL_STORED_OBJECT_DEBUG_CREATED(this, "Approximate integration method");
290  pai = p; im_type = IM_APPROX;
291  }
292 
294  DAL_STORED_OBJECT_DEBUG_CREATED(this, "Integration method");
295  im_type = IM_NONE;
296  }
297 
298  virtual ~integration_method()
299  { DAL_STORED_OBJECT_DEBUG_DESTROYED(this, "Integration method"); }
300  };
301 
302 
303  /** Get an integration method from its name .
304  @see @ref im_list
305  @param name the integration method name, for example "IM_TRIANGLE(6)"
306  @param throw_if_not_found choose if an exception must be thrown
307  when the integration method does not exist (if no exception, a
308  null pointer is returned).
309  */
310  pintegration_method int_method_descriptor(std::string name,
311  bool throw_if_not_found = true);
312 
313  /** Get the string name of an integration method .
314  */
315  std::string name_of_int_method(pintegration_method p);
316 
317  /**
318  return an exact integration method for convex type handled by pgt.
319  If pgt is not linear, classical_exact_im will fail.
320  */
321  pintegration_method classical_exact_im(bgeot::pgeometric_trans pgt);
322  /**
323  try to find an approximate integration method for the geometric
324  transformation pgt which is able to integrate exactly polynomials
325  of degree <= "degree". It may return a higher order integration
326  method if no method match the exact degree.
327  */
328  pintegration_method classical_approx_im(bgeot::pgeometric_trans pgt, dim_type degree);
329 
330  /// return IM_EXACT_SIMPLEX(n)
331  pintegration_method exact_simplex_im(size_type n);
332  /// return IM_EXACT_PARALLELEPIPED(n)
333  pintegration_method exact_parallelepiped_im(size_type n);
334  /// return IM_EXACT_PRISM(n)
335  pintegration_method exact_prism_im(size_type n);
336  /// use classical_exact_im instead.
337  pintegration_method exact_classical_im(bgeot::pgeometric_trans pgt) IS_DEPRECATED;
338  /// return IM_NONE
339  pintegration_method im_none(void);
340 
341 
342  class mesh_im;
343  papprox_integration composite_approx_int_method(const bgeot::mesh_precomposite &mp,
344  const mesh_im &mf,
345  bgeot::pconvex_ref cr);
346 
347  /* try to integrate all monomials up to order 'order' and return the
348  max. error */
349  scalar_type test_integration_error(papprox_integration pim, dim_type order);
350 
351  papprox_integration get_approx_im_or_fail(pintegration_method pim);
352 
353  /* Function allowing the add of an integration method outwards
354  of getfem_integration.cc */
355 
356  typedef dal::naming_system<integration_method>::param_list im_param_list;
357 
358  void add_integration_name(std::string name,
360 
361 } /* end of namespace getfem. */
362 
363 
364 #endif /* BGEOT_INTEGRATION_H__ */
Reference convexes.
Geometric transformations on convexes.
Structure which dynamically collects points identifying points that are nearer than a certain very sm...
Handle composite polynomials.
Store a set of points, identifying points that are nearer than a certain very small distance.
Vector of integer (16 bits type) which represent the powers of a monomial.
Definition: bgeot_poly.h:63
Associate a name to a method descriptor and store method descriptors.
base class for static stored objects
this structure is not intended to be used directly.
Describe an integration method linked to a mesh.
Description of an exact integration of polynomials.
long_scalar_type int_poly(const base_poly &P) const
Evaluate the integral of the polynomial P on the reference element.
bgeot::pconvex_structure structure(void) const
{Structure of convex of reference.
long_scalar_type int_poly_on_face(const base_poly &P, short_type f) const
Evaluate the integral of the polynomial P on the face f of the reference element.
dim_type dim(void) const
Dimension of convex of reference.
Naming system.
defines and typedefs for namespace getfem
gmm::uint16_type short_type
used as the common short type integer in the library
Definition: bgeot_config.h:72
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
pconvex_structure basic_structure(pconvex_structure cv)
Original structure (if concerned)
std::shared_ptr< const bgeot::geometric_trans > pgeometric_trans
pointer type for a geometric transformation
GEneric Tool for Finite Element Methods.
pintegration_method exact_parallelepiped_im(size_type n)
return IM_EXACT_PARALLELEPIPED(n)
std::string name_of_int_method(pintegration_method p)
Get the string name of an integration method .
pintegration_method exact_simplex_im(size_type n)
return IM_EXACT_SIMPLEX(n)
pintegration_method exact_prism_im(size_type n)
return IM_EXACT_PRISM(n)
pintegration_method 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 i...
pintegration_method classical_exact_im(bgeot::pgeometric_trans pgt)
return an exact integration method for convex type handled by pgt.
pintegration_method exact_classical_im(bgeot::pgeometric_trans pgt) IS_DEPRECATED
use classical_exact_im instead.
pintegration_method int_method_descriptor(std::string name, bool throw_if_not_found=true)
Get an integration method from its name .
integration_method_type
the list of main integration method types
pintegration_method im_none(void)
return IM_NONE